ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_mathGammaAM.F90
Go to the documentation of this file.
1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3!!!! !!!!
4!!!! ParaMonte: Parallel Monte Carlo and Machine Learning Library. !!!!
5!!!! !!!!
6!!!! Copyright (C) 2012-present, The Computational Data Science Lab !!!!
7!!!! !!!!
8!!!! This file is part of the ParaMonte library. !!!!
9!!!! !!!!
10!!!! LICENSE !!!!
11!!!! !!!!
12!!!! https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md !!!!
13!!!! !!!!
14!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16
56
58
59 use pm_kind, only: SK, IK, LK
60
61 implicit none
62
63 character(*, SK), parameter :: MODULE_NAME = "@pm_mathGammaAM"
64
65!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
66
143
144 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
145
146#if RK5_ENABLED
147 impure elemental module function getGammaIncLowAM_RK5(x, kappa) result(gammaIncLow)
148#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
149 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncLowAM_RK5
150#endif
151 use pm_kind, only: RKG => RK5
152 real(RKG) , intent(in) :: x, kappa
153 real(RKG) :: gammaIncLow
154 end function
155#endif
156
157#if RK4_ENABLED
158 impure elemental module function getGammaIncLowAM_RK4(x, kappa) result(gammaIncLow)
159#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
160 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncLowAM_RK4
161#endif
162 use pm_kind, only: RKG => RK4
163 real(RKG) , intent(in) :: x, kappa
164 real(RKG) :: gammaIncLow
165 end function
166#endif
167
168#if RK3_ENABLED
169 impure elemental module function getGammaIncLowAM_RK3(x, kappa) result(gammaIncLow)
170#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
171 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncLowAM_RK3
172#endif
173 use pm_kind, only: RKG => RK3
174 real(RKG) , intent(in) :: x, kappa
175 real(RKG) :: gammaIncLow
176 end function
177#endif
178
179#if RK2_ENABLED
180 impure elemental module function getGammaIncLowAM_RK2(x, kappa) result(gammaIncLow)
181#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
182 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncLowAM_RK2
183#endif
184 use pm_kind, only: RKG => RK2
185 real(RKG) , intent(in) :: x, kappa
186 real(RKG) :: gammaIncLow
187 end function
188#endif
189
190#if RK1_ENABLED
191 impure elemental module function getGammaIncLowAM_RK1(x, kappa) result(gammaIncLow)
192#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
193 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncLowAM_RK1
194#endif
195 use pm_kind, only: RKG => RK1
196 real(RKG) , intent(in) :: x, kappa
197 real(RKG) :: gammaIncLow
198 end function
199#endif
200
201 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
202
203 end interface
204
205!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
206
280
281 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
282
283#if RK5_ENABLED
284 impure elemental module function getGammaIncUppAM_RK5(x, kappa) result(gammaIncUpp)
285#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
286 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncUppAM_RK5
287#endif
288 use pm_kind, only: RKG => RK5
289 real(RKG) , intent(in) :: x, kappa
290 real(RKG) :: gammaIncUpp
291 end function
292#endif
293
294#if RK4_ENABLED
295 impure elemental module function getGammaIncUppAM_RK4(x, kappa) result(gammaIncUpp)
296#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
297 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncUppAM_RK4
298#endif
299 use pm_kind, only: RKG => RK4
300 real(RKG) , intent(in) :: x, kappa
301 real(RKG) :: gammaIncUpp
302 end function
303#endif
304
305#if RK3_ENABLED
306 impure elemental module function getGammaIncUppAM_RK3(x, kappa) result(gammaIncUpp)
307#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
308 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncUppAM_RK3
309#endif
310 use pm_kind, only: RKG => RK3
311 real(RKG) , intent(in) :: x, kappa
312 real(RKG) :: gammaIncUpp
313 end function
314#endif
315
316#if RK2_ENABLED
317 impure elemental module function getGammaIncUppAM_RK2(x, kappa) result(gammaIncUpp)
318#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
319 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncUppAM_RK2
320#endif
321 use pm_kind, only: RKG => RK2
322 real(RKG) , intent(in) :: x, kappa
323 real(RKG) :: gammaIncUpp
324 end function
325#endif
326
327#if RK1_ENABLED
328 impure elemental module function getGammaIncUppAM_RK1(x, kappa) result(gammaIncUpp)
329#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
330 !DEC$ ATTRIBUTES DLLEXPORT :: getGammaIncUppAM_RK1
331#endif
332 use pm_kind, only: RKG => RK1
333 real(RKG) , intent(in) :: x, kappa
334 real(RKG) :: gammaIncUpp
335 end function
336#endif
337
338 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
339
340 end interface
341
342!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
343
417
418 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
419
420#if RK5_ENABLED
421 PURE elemental module subroutine setGammaIncLowAM_RK5(gammaIncLow, x, logGammaKappa, kappa, info)
422#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
423 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncLowAM_RK5
424#endif
425 use pm_kind, only: RKG => RK5
426 real(RKG) , intent(out) :: gammaIncLow
427 real(RKG) , intent(in) :: x, logGammaKappa, kappa
428 integer(IK) , intent(out) :: info
429 end subroutine
430#endif
431
432#if RK4_ENABLED
433 PURE elemental module subroutine setGammaIncLowAM_RK4(gammaIncLow, x, logGammaKappa, kappa, info)
434#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
435 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncLowAM_RK4
436#endif
437 use pm_kind, only: RKG => RK4
438 real(RKG) , intent(out) :: gammaIncLow
439 real(RKG) , intent(in) :: x, logGammaKappa, kappa
440 integer(IK) , intent(out) :: info
441 end subroutine
442#endif
443
444#if RK3_ENABLED
445 PURE elemental module subroutine setGammaIncLowAM_RK3(gammaIncLow, x, logGammaKappa, kappa, info)
446#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
447 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncLowAM_RK3
448#endif
449 use pm_kind, only: RKG => RK3
450 real(RKG) , intent(out) :: gammaIncLow
451 real(RKG) , intent(in) :: x, logGammaKappa, kappa
452 integer(IK) , intent(out) :: info
453 end subroutine
454#endif
455
456#if RK2_ENABLED
457 PURE elemental module subroutine setGammaIncLowAM_RK2(gammaIncLow, x, logGammaKappa, kappa, info)
458#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
459 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncLowAM_RK2
460#endif
461 use pm_kind, only: RKG => RK2
462 real(RKG) , intent(out) :: gammaIncLow
463 real(RKG) , intent(in) :: x, logGammaKappa, kappa
464 integer(IK) , intent(out) :: info
465 end subroutine
466#endif
467
468#if RK1_ENABLED
469 PURE elemental module subroutine setGammaIncLowAM_RK1(gammaIncLow, x, logGammaKappa, kappa, info)
470#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
471 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncLowAM_RK1
472#endif
473 use pm_kind, only: RKG => RK1
474 real(RKG) , intent(out) :: gammaIncLow
475 real(RKG) , intent(in) :: x, logGammaKappa, kappa
476 integer(IK) , intent(out) :: info
477 end subroutine
478#endif
479
480 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
481
482 end interface
483
484!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
485
565
566 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
567
568#if RK5_ENABLED
569 PURE elemental module subroutine setGammaIncUppAM_RK5(gammaIncUpp, x, logGammaKappa, kappa, info)
570#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
571 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncUppAM_RK5
572#endif
573 use pm_kind, only: RKG => RK5
574 real(RKG) , intent(out) :: gammaIncUpp
575 real(RKG) , intent(in) :: x, logGammaKappa, kappa
576 integer(IK) , intent(out) :: info
577 end subroutine
578#endif
579
580#if RK4_ENABLED
581 PURE elemental module subroutine setGammaIncUppAM_RK4(gammaIncUpp, x, logGammaKappa, kappa, info)
582#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
583 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncUppAM_RK4
584#endif
585 use pm_kind, only: RKG => RK4
586 real(RKG) , intent(out) :: gammaIncUpp
587 real(RKG) , intent(in) :: x, logGammaKappa, kappa
588 integer(IK) , intent(out) :: info
589 end subroutine
590#endif
591
592#if RK3_ENABLED
593 PURE elemental module subroutine setGammaIncUppAM_RK3(gammaIncUpp, x, logGammaKappa, kappa, info)
594#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
595 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncUppAM_RK3
596#endif
597 use pm_kind, only: RKG => RK3
598 real(RKG) , intent(out) :: gammaIncUpp
599 real(RKG) , intent(in) :: x, logGammaKappa, kappa
600 integer(IK) , intent(out) :: info
601 end subroutine
602#endif
603
604#if RK2_ENABLED
605 PURE elemental module subroutine setGammaIncUppAM_RK2(gammaIncUpp, x, logGammaKappa, kappa, info)
606#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
607 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncUppAM_RK2
608#endif
609 use pm_kind, only: RKG => RK2
610 real(RKG) , intent(out) :: gammaIncUpp
611 real(RKG) , intent(in) :: x, logGammaKappa, kappa
612 integer(IK) , intent(out) :: info
613 end subroutine
614#endif
615
616#if RK1_ENABLED
617 PURE elemental module subroutine setGammaIncUppAM_RK1(gammaIncUpp, x, logGammaKappa, kappa, info)
618#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
619 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncUppAM_RK1
620#endif
621 use pm_kind, only: RKG => RK1
622 real(RKG) , intent(out) :: gammaIncUpp
623 real(RKG) , intent(in) :: x, logGammaKappa, kappa
624 integer(IK) , intent(out) :: info
625 end subroutine
626#endif
627
628 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
629
630 end interface
631
632!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
633
635 ! This algorithm is taken from Abergel and Moisan, 2020, Algorithm 1006: Fast and Accurate Evaluation of a Generalized Incomplete Gamma Function
636 ! It claims superb performance and accuracy. However, I could not substantiate their claims regarding accuracy, particularly,
637 ! for very large `kappa` and `x` values `(kappa, x) > 1e7`.
638 interface setGammaIncAM
639
640 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
641
642#if RK5_ENABLED
643 PURE elemental module subroutine setGammaIncAMDef_RK5(gammaInc, x, kappa, info)
644#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
645 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMDef_RK5
646#endif
647 use pm_kind, only: RKG => RK5
648 real(RKG) , intent(out) :: gammaInc
649 real(RKG) , intent(in) :: x, kappa
650 integer(IK) , intent(out) :: info
651 end subroutine
652#endif
653
654#if RK4_ENABLED
655 PURE elemental module subroutine setGammaIncAMDef_RK4(gammaInc, x, kappa, info)
656#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
657 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMDef_RK4
658#endif
659 use pm_kind, only: RKG => RK4
660 real(RKG) , intent(out) :: gammaInc
661 real(RKG) , intent(in) :: x, kappa
662 integer(IK) , intent(out) :: info
663 end subroutine
664#endif
665
666#if RK3_ENABLED
667 PURE elemental module subroutine setGammaIncAMDef_RK3(gammaInc, x, kappa, info)
668#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
669 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMDef_RK3
670#endif
671 use pm_kind, only: RKG => RK3
672 real(RKG) , intent(out) :: gammaInc
673 real(RKG) , intent(in) :: x, kappa
674 integer(IK) , intent(out) :: info
675 end subroutine
676#endif
677
678#if RK2_ENABLED
679 PURE elemental module subroutine setGammaIncAMDef_RK2(gammaInc, x, kappa, info)
680#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
681 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMDef_RK2
682#endif
683 use pm_kind, only: RKG => RK2
684 real(RKG) , intent(out) :: gammaInc
685 real(RKG) , intent(in) :: x, kappa
686 integer(IK) , intent(out) :: info
687 end subroutine
688#endif
689
690#if RK1_ENABLED
691 PURE elemental module subroutine setGammaIncAMDef_RK1(gammaInc, x, kappa, info)
692#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
693 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMDef_RK1
694#endif
695 use pm_kind, only: RKG => RK1
696 real(RKG) , intent(out) :: gammaInc
697 real(RKG) , intent(in) :: x, kappa
698 integer(IK) , intent(out) :: info
699 end subroutine
700#endif
701
702 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
703
704#if RK5_ENABLED
705 PURE elemental module subroutine setGammaIncAMNXIK_RK5(gammaInc, x, logGammaKappa, kappa, info)
706#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
707 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMNXIK_RK5
708#endif
709 use pm_kind, only: RKG => RK5
710 real(RKG) , intent(out) :: gammaInc
711 real(RKG) , intent(in) :: x, logGammaKappa
712 integer(IK) , intent(in) :: kappa
713 integer(IK) , intent(out) :: info
714 end subroutine
715#endif
716
717#if RK4_ENABLED
718 PURE elemental module subroutine setGammaIncAMNXIK_RK4(gammaInc, x, logGammaKappa, kappa, info)
719#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
720 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMNXIK_RK4
721#endif
722 use pm_kind, only: RKG => RK4
723 real(RKG) , intent(out) :: gammaInc
724 real(RKG) , intent(in) :: x, logGammaKappa
725 integer(IK) , intent(in) :: kappa
726 integer(IK) , intent(out) :: info
727 end subroutine
728#endif
729
730#if RK3_ENABLED
731 PURE elemental module subroutine setGammaIncAMNXIK_RK3(gammaInc, x, logGammaKappa, kappa, info)
732#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
733 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMNXIK_RK3
734#endif
735 use pm_kind, only: RKG => RK3
736 real(RKG) , intent(out) :: gammaInc
737 real(RKG) , intent(in) :: x, logGammaKappa
738 integer(IK) , intent(in) :: kappa
739 integer(IK) , intent(out) :: info
740 end subroutine
741#endif
742
743#if RK2_ENABLED
744 PURE elemental module subroutine setGammaIncAMNXIK_RK2(gammaInc, x, logGammaKappa, kappa, info)
745#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
746 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMNXIK_RK2
747#endif
748 use pm_kind, only: RKG => RK2
749 real(RKG) , intent(out) :: gammaInc
750 real(RKG) , intent(in) :: x, logGammaKappa
751 integer(IK) , intent(in) :: kappa
752 integer(IK) , intent(out) :: info
753 end subroutine
754#endif
755
756#if RK1_ENABLED
757 PURE elemental module subroutine setGammaIncAMNXIK_RK1(gammaInc, x, logGammaKappa, kappa, info)
758#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
759 !DEC$ ATTRIBUTES DLLEXPORT :: setGammaIncAMNXIK_RK1
760#endif
761 use pm_kind, only: RKG => RK1
762 real(RKG) , intent(out) :: gammaInc
763 real(RKG) , intent(in) :: x, logGammaKappa
764 integer(IK) , intent(in) :: kappa
765 integer(IK) , intent(out) :: info
766 end subroutine
767#endif
768
769 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
770
771 end interface
773
774!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
775
776end module pm_mathGammaAM ! LCOV_EXCL_LINE
Generate and return the regularized Lower Incomplete Gamma function for the specified shape parameter...
Generate and return the regularized Upper Incomplete Gamma function for the specified shape parameter...
Return the regularized Lower Incomplete Gamma function for the specified shape parameter ( ) and uppe...
Return the regularized Upper Incomplete Gamma function for the specified shape parameter ( ) and lowe...
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter RK5
Definition: pm_kind.F90:478
integer, parameter RK4
Definition: pm_kind.F90:489
integer, parameter RK2
Definition: pm_kind.F90:511
integer, parameter RK3
Definition: pm_kind.F90:500
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RK1
Definition: pm_kind.F90:522
This module contains procedures and generic interfaces for the Lower and Upper Incomplete Gamma funct...
character(*, SK), parameter MODULE_NAME