ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_distKolm.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
81
82!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
83
85
86 use pm_kind, only: SK, IK
87
88 implicit none
89
90 character(*, SK), parameter :: MODULE_NAME = "@pm_distKolm"
91
92!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
93
128 end type
129
130!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
131
189 interface getKolmPDF
190
191 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
192
193#if RK5_ENABLED
194 PURE elemental module function getKolmPDF_RK5(x) result(pdf)
195#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
196 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmPDF_RK5
197#endif
198 use pm_kind, only: RKG => RK5
199 real(RKG) , intent(in) :: x
200 real(RKG) :: pdf
201 end function
202#endif
203
204#if RK4_ENABLED
205 PURE elemental module function getKolmPDF_RK4(x) result(pdf)
206#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
207 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmPDF_RK4
208#endif
209 use pm_kind, only: RKG => RK4
210 real(RKG) , intent(in) :: x
211 real(RKG) :: pdf
212 end function
213#endif
214
215#if RK3_ENABLED
216 PURE elemental module function getKolmPDF_RK3(x) result(pdf)
217#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
218 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmPDF_RK3
219#endif
220 use pm_kind, only: RKG => RK3
221 real(RKG) , intent(in) :: x
222 real(RKG) :: pdf
223 end function
224#endif
225
226#if RK2_ENABLED
227 PURE elemental module function getKolmPDF_RK2(x) result(pdf)
228#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
229 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmPDF_RK2
230#endif
231 use pm_kind, only: RKG => RK2
232 real(RKG) , intent(in) :: x
233 real(RKG) :: pdf
234 end function
235#endif
236
237#if RK1_ENABLED
238 PURE elemental module function getKolmPDF_RK1(x) result(pdf)
239#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
240 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmPDF_RK1
241#endif
242 use pm_kind, only: RKG => RK1
243 real(RKG) , intent(in) :: x
244 real(RKG) :: pdf
245 end function
246#endif
247
248 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
249
250 end interface
251
252!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
253
310 interface setKolmPDF
311
312 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
313
314#if RK5_ENABLED
315 PURE elemental module subroutine setKolmPDF_RK5(pdf, x)
316#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
317 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmPDF_RK5
318#endif
319 use pm_kind, only: RKG => RK5
320 real(RKG) , intent(out) :: pdf
321 real(RKG) , intent(in) :: x
322 end subroutine
323#endif
324
325#if RK4_ENABLED
326 PURE elemental module subroutine setKolmPDF_RK4(pdf, x)
327#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
328 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmPDF_RK4
329#endif
330 use pm_kind, only: RKG => RK4
331 real(RKG) , intent(out) :: pdf
332 real(RKG) , intent(in) :: x
333 end subroutine
334#endif
335
336#if RK3_ENABLED
337 PURE elemental module subroutine setKolmPDF_RK3(pdf, x)
338#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
339 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmPDF_RK3
340#endif
341 use pm_kind, only: RKG => RK3
342 real(RKG) , intent(out) :: pdf
343 real(RKG) , intent(in) :: x
344 end subroutine
345#endif
346
347#if RK2_ENABLED
348 PURE elemental module subroutine setKolmPDF_RK2(pdf, x)
349#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
350 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmPDF_RK2
351#endif
352 use pm_kind, only: RKG => RK2
353 real(RKG) , intent(out) :: pdf
354 real(RKG) , intent(in) :: x
355 end subroutine
356#endif
357
358#if RK1_ENABLED
359 PURE elemental module subroutine setKolmPDF_RK1(pdf, x)
360#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
361 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmPDF_RK1
362#endif
363 use pm_kind, only: RKG => RK1
364 real(RKG) , intent(out) :: pdf
365 real(RKG) , intent(in) :: x
366 end subroutine
367#endif
368
369 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
370
371 end interface
372
373!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
374
432 interface getKolmCDF
433
434 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
435
436#if RK5_ENABLED
437 PURE elemental module function getKolmCDF_RK5(x) result(cdf)
438#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
439 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmCDF_RK5
440#endif
441 use pm_kind, only: RKG => RK5
442 real(RKG) , intent(in) :: x
443 real(RKG) :: cdf
444 end function
445#endif
446
447#if RK4_ENABLED
448 PURE elemental module function getKolmCDF_RK4(x) result(cdf)
449#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
450 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmCDF_RK4
451#endif
452 use pm_kind, only: RKG => RK4
453 real(RKG) , intent(in) :: x
454 real(RKG) :: cdf
455 end function
456#endif
457
458#if RK3_ENABLED
459 PURE elemental module function getKolmCDF_RK3(x) result(cdf)
460#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
461 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmCDF_RK3
462#endif
463 use pm_kind, only: RKG => RK3
464 real(RKG) , intent(in) :: x
465 real(RKG) :: cdf
466 end function
467#endif
468
469#if RK2_ENABLED
470 PURE elemental module function getKolmCDF_RK2(x) result(cdf)
471#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
472 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmCDF_RK2
473#endif
474 use pm_kind, only: RKG => RK2
475 real(RKG) , intent(in) :: x
476 real(RKG) :: cdf
477 end function
478#endif
479
480#if RK1_ENABLED
481 PURE elemental module function getKolmCDF_RK1(x) result(cdf)
482#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
483 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmCDF_RK1
484#endif
485 use pm_kind, only: RKG => RK1
486 real(RKG) , intent(in) :: x
487 real(RKG) :: cdf
488 end function
489#endif
490
491 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
492
493 end interface
494
495!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
496
553 interface setKolmCDF
554
555 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
556
557#if RK5_ENABLED
558 PURE elemental module subroutine setKolmCDF_RK5(cdf, x)
559#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
560 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmCDF_RK5
561#endif
562 use pm_kind, only: RKG => RK5
563 real(RKG) , intent(out) :: cdf
564 real(RKG) , intent(in) :: x
565 end subroutine
566#endif
567
568#if RK4_ENABLED
569 PURE elemental module subroutine setKolmCDF_RK4(cdf, x)
570#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
571 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmCDF_RK4
572#endif
573 use pm_kind, only: RKG => RK4
574 real(RKG) , intent(out) :: cdf
575 real(RKG) , intent(in) :: x
576 end subroutine
577#endif
578
579#if RK3_ENABLED
580 PURE elemental module subroutine setKolmCDF_RK3(cdf, x)
581#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
582 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmCDF_RK3
583#endif
584 use pm_kind, only: RKG => RK3
585 real(RKG) , intent(out) :: cdf
586 real(RKG) , intent(in) :: x
587 end subroutine
588#endif
589
590#if RK2_ENABLED
591 PURE elemental module subroutine setKolmCDF_RK2(cdf, x)
592#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
593 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmCDF_RK2
594#endif
595 use pm_kind, only: RKG => RK2
596 real(RKG) , intent(out) :: cdf
597 real(RKG) , intent(in) :: x
598 end subroutine
599#endif
600
601#if RK1_ENABLED
602 PURE elemental module subroutine setKolmCDF_RK1(cdf, x)
603#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
604 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmCDF_RK1
605#endif
606 use pm_kind, only: RKG => RK1
607 real(RKG) , intent(out) :: cdf
608 real(RKG) , intent(in) :: x
609 end subroutine
610#endif
611
612 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
613
614 end interface
615
616!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
617
675 interface getKolmQuan
676
677 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
678
679#if RK5_ENABLED
680 PURE elemental module function getKolmQuan_RK5(cdf) result(quan)
681#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
682 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmQuan_RK5
683#endif
684 use pm_kind, only: RKG => RK5
685 real(RKG) , intent(in) :: cdf
686 real(RKG) :: quan
687 end function
688#endif
689
690#if RK4_ENABLED
691 PURE elemental module function getKolmQuan_RK4(cdf) result(quan)
692#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
693 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmQuan_RK4
694#endif
695 use pm_kind, only: RKG => RK4
696 real(RKG) , intent(in) :: cdf
697 real(RKG) :: quan
698 end function
699#endif
700
701#if RK3_ENABLED
702 PURE elemental module function getKolmQuan_RK3(cdf) result(quan)
703#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
704 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmQuan_RK3
705#endif
706 use pm_kind, only: RKG => RK3
707 real(RKG) , intent(in) :: cdf
708 real(RKG) :: quan
709 end function
710#endif
711
712#if RK2_ENABLED
713 PURE elemental module function getKolmQuan_RK2(cdf) result(quan)
714#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
715 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmQuan_RK2
716#endif
717 use pm_kind, only: RKG => RK2
718 real(RKG) , intent(in) :: cdf
719 real(RKG) :: quan
720 end function
721#endif
722
723#if RK1_ENABLED
724 PURE elemental module function getKolmQuan_RK1(cdf) result(quan)
725#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
726 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmQuan_RK1
727#endif
728 use pm_kind, only: RKG => RK1
729 real(RKG) , intent(in) :: cdf
730 real(RKG) :: quan
731 end function
732#endif
733
734 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
735
736 end interface
737
738!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
739
792 interface setKolmQuan
793
794 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
795
796#if RK5_ENABLED
797 PURE elemental module subroutine setKolmQuan_RK5(quan, cdf)
798#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
799 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmQuan_RK5
800#endif
801 use pm_kind, only: RKG => RK5
802 real(RKG) , intent(out) :: quan
803 real(RKG) , intent(in) :: cdf
804 end subroutine
805#endif
806
807#if RK4_ENABLED
808 PURE elemental module subroutine setKolmQuan_RK4(quan, cdf)
809#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
810 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmQuan_RK4
811#endif
812 use pm_kind, only: RKG => RK4
813 real(RKG) , intent(out) :: quan
814 real(RKG) , intent(in) :: cdf
815 end subroutine
816#endif
817
818#if RK3_ENABLED
819 PURE elemental module subroutine setKolmQuan_RK3(quan, cdf)
820#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
821 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmQuan_RK3
822#endif
823 use pm_kind, only: RKG => RK3
824 real(RKG) , intent(out) :: quan
825 real(RKG) , intent(in) :: cdf
826 end subroutine
827#endif
828
829#if RK2_ENABLED
830 PURE elemental module subroutine setKolmQuan_RK2(quan, cdf)
831#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
832 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmQuan_RK2
833#endif
834 use pm_kind, only: RKG => RK2
835 real(RKG) , intent(out) :: quan
836 real(RKG) , intent(in) :: cdf
837 end subroutine
838#endif
839
840#if RK1_ENABLED
841 PURE elemental module subroutine setKolmQuan_RK1(quan, cdf)
842#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
843 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmQuan_RK1
844#endif
845 use pm_kind, only: RKG => RK1
846 real(RKG) , intent(out) :: quan
847 real(RKG) , intent(in) :: cdf
848 end subroutine
849#endif
850
851 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
852
853 end interface
854
855!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
856
915 interface getKolmRand
916
917 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
918
919#if RK5_ENABLED
920 impure elemental module function getKolmRand_RK5(unif) result(rand)
921#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
922 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmRand_RK5
923#endif
924 use pm_kind, only: RKG => RK5
925 real(RKG) , intent(in) :: unif
926 real(RKG) :: rand
927 end function
928#endif
929
930#if RK4_ENABLED
931 impure elemental module function getKolmRand_RK4(unif) result(rand)
932#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
933 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmRand_RK4
934#endif
935 use pm_kind, only: RKG => RK4
936 real(RKG) , intent(in) :: unif
937 real(RKG) :: rand
938 end function
939#endif
940
941#if RK3_ENABLED
942 impure elemental module function getKolmRand_RK3(unif) result(rand)
943#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
944 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmRand_RK3
945#endif
946 use pm_kind, only: RKG => RK3
947 real(RKG) , intent(in) :: unif
948 real(RKG) :: rand
949 end function
950#endif
951
952#if RK2_ENABLED
953 impure elemental module function getKolmRand_RK2(unif) result(rand)
954#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
955 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmRand_RK2
956#endif
957 use pm_kind, only: RKG => RK2
958 real(RKG) , intent(in) :: unif
959 real(RKG) :: rand
960 end function
961#endif
962
963#if RK1_ENABLED
964 impure elemental module function getKolmRand_RK1(unif) result(rand)
965#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
966 !DEC$ ATTRIBUTES DLLEXPORT :: getKolmRand_RK1
967#endif
968 use pm_kind, only: RKG => RK1
969 real(RKG) , intent(in) :: unif
970 real(RKG) :: rand
971 end function
972#endif
973
974 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
975
976 end interface
977
978!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
979
1037 interface setKolmRand
1038
1039 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1040
1041#if RK5_ENABLED
1042 PURE elemental module subroutine setKolmRand_RK5(rand, unif)
1043#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1044 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmRand_RK5
1045#endif
1046 use pm_kind, only: RKG => RK5
1047 real(RKG) , intent(out) :: rand
1048 real(RKG) , intent(in) :: unif
1049 end subroutine
1050#endif
1051
1052#if RK4_ENABLED
1053 PURE elemental module subroutine setKolmRand_RK4(rand, unif)
1054#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1055 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmRand_RK4
1056#endif
1057 use pm_kind, only: RKG => RK4
1058 real(RKG) , intent(out) :: rand
1059 real(RKG) , intent(in) :: unif
1060 end subroutine
1061#endif
1062
1063#if RK3_ENABLED
1064 PURE elemental module subroutine setKolmRand_RK3(rand, unif)
1065#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1066 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmRand_RK3
1067#endif
1068 use pm_kind, only: RKG => RK3
1069 real(RKG) , intent(out) :: rand
1070 real(RKG) , intent(in) :: unif
1071 end subroutine
1072#endif
1073
1074#if RK2_ENABLED
1075 PURE elemental module subroutine setKolmRand_RK2(rand, unif)
1076#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1077 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmRand_RK2
1078#endif
1079 use pm_kind, only: RKG => RK2
1080 real(RKG) , intent(out) :: rand
1081 real(RKG) , intent(in) :: unif
1082 end subroutine
1083#endif
1084
1085#if RK1_ENABLED
1086 PURE elemental module subroutine setKolmRand_RK1(rand, unif)
1087#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1088 !DEC$ ATTRIBUTES DLLEXPORT :: setKolmRand_RK1
1089#endif
1090 use pm_kind, only: RKG => RK1
1091 real(RKG) , intent(out) :: rand
1092 real(RKG) , intent(in) :: unif
1093 end subroutine
1094#endif
1095
1096 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1097
1098 end interface
1099
1100!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1101
1102end module pm_distKolm
Generate and return the Cumulative Distribution Function (CDF) of the Kolmogorov distribution for an ...
Generate and return the Probability Density Function (PDF) of the Kolmogorov distribution for an inpu...
Generate and return a scalar (or array of arbitrary rank) of the quantile corresponding to the specif...
Generate and return a scalar (or array of arbitrary rank) of the random value(s) from the Kolmogorov ...
Return the Cumulative Distribution Function (CDF) of the Kolmogorov distribution for an input x withi...
Return the Probability Density Function (PDF) of the Kolmogorov distribution for an input x within th...
Return a scalar (or array of arbitrary rank) of the quantile corresponding to the specified CDF of Ko...
Return a scalar (or array of arbitrary rank) of the random value(s) from the Kolmogorov distribution.
This module contains classes and procedures for computing various statistical quantities related to t...
Definition: pm_distKolm.F90:84
character(*, SK), parameter MODULE_NAME
Definition: pm_distKolm.F90:90
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 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 is the derived type for signifying distributions that are of type Kolmogorov as defined in the d...