ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_distBeta.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
96
97!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
98
100
101 use pm_kind, only: SK, IK, LK
102 use pm_distUnif, only: rngf_type
104
105 implicit none
106
107 character(*, SK), parameter :: MODULE_NAME = "@pm_distBeta"
108
109!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
110
145 end type
146
147!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
148
227 interface getBetaPDF
228
229 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
230
231#if RK5_ENABLED
232 PURE elemental module function getBetaPDF_RK5(x, alpha, beta) result(pdf)
233#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
234 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaPDF_RK5
235#endif
236 use pm_kind, only: RKG => RK5
237 real(RKG) , intent(in) :: x, alpha, beta
238 real(RKG) :: pdf
239 end function
240#endif
241
242#if RK4_ENABLED
243 PURE elemental module function getBetaPDF_RK4(x, alpha, beta) result(pdf)
244#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
245 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaPDF_RK4
246#endif
247 use pm_kind, only: RKG => RK4
248 real(RKG) , intent(in) :: x, alpha, beta
249 real(RKG) :: pdf
250 end function
251#endif
252
253#if RK3_ENABLED
254 PURE elemental module function getBetaPDF_RK3(x, alpha, beta) result(pdf)
255#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
256 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaPDF_RK3
257#endif
258 use pm_kind, only: RKG => RK3
259 real(RKG) , intent(in) :: x, alpha, beta
260 real(RKG) :: pdf
261 end function
262#endif
263
264#if RK2_ENABLED
265 PURE elemental module function getBetaPDF_RK2(x, alpha, beta) result(pdf)
266#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
267 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaPDF_RK2
268#endif
269 use pm_kind, only: RKG => RK2
270 real(RKG) , intent(in) :: x, alpha, beta
271 real(RKG) :: pdf
272 end function
273#endif
274
275#if RK1_ENABLED
276 PURE elemental module function getBetaPDF_RK1(x, alpha, beta) result(pdf)
277#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
278 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaPDF_RK1
279#endif
280 use pm_kind, only: RKG => RK1
281 real(RKG) , intent(in) :: x, alpha, beta
282 real(RKG) :: pdf
283 end function
284#endif
285
286 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
287
288 end interface
289
290!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
291
358
359 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
360
361#if RK5_ENABLED
362 PURE elemental module function getBetaLogPDFD_RK5(x, alpha, beta) result(logPDF)
363#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
364 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFD_RK5
365#endif
366 use pm_kind, only: RKG => RK5
367 real(RKG) , intent(in) :: x, alpha, beta
368 real(RKG) :: logPDF
369 end function
370#endif
371
372#if RK4_ENABLED
373 PURE elemental module function getBetaLogPDFD_RK4(x, alpha, beta) result(logPDF)
374#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
375 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFD_RK4
376#endif
377 use pm_kind, only: RKG => RK4
378 real(RKG) , intent(in) :: x, alpha, beta
379 real(RKG) :: logPDF
380 end function
381#endif
382
383#if RK3_ENABLED
384 PURE elemental module function getBetaLogPDFD_RK3(x, alpha, beta) result(logPDF)
385#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
386 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFD_RK3
387#endif
388 use pm_kind, only: RKG => RK3
389 real(RKG) , intent(in) :: x, alpha, beta
390 real(RKG) :: logPDF
391 end function
392#endif
393
394#if RK2_ENABLED
395 PURE elemental module function getBetaLogPDFD_RK2(x, alpha, beta) result(logPDF)
396#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
397 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFD_RK2
398#endif
399 use pm_kind, only: RKG => RK2
400 real(RKG) , intent(in) :: x, alpha, beta
401 real(RKG) :: logPDF
402 end function
403#endif
404
405#if RK1_ENABLED
406 PURE elemental module function getBetaLogPDFD_RK1(x, alpha, beta) result(logPDF)
407#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
408 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFD_RK1
409#endif
410 use pm_kind, only: RKG => RK1
411 real(RKG) , intent(in) :: x, alpha, beta
412 real(RKG) :: logPDF
413 end function
414#endif
415
416 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
417
418#if RK5_ENABLED
419 PURE elemental module function getBetaLogPDFL_RK5(x, alpha, beta, logBeta) result(logPDF)
420#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
421 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFL_RK5
422#endif
423 use pm_kind, only: RKG => RK5
424 real(RKG) , intent(in) :: x, alpha, beta, logBeta
425 real(RKG) :: logPDF
426 end function
427#endif
428
429#if RK4_ENABLED
430 PURE elemental module function getBetaLogPDFL_RK4(x, alpha, beta, logBeta) result(logPDF)
431#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
432 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFL_RK4
433#endif
434 use pm_kind, only: RKG => RK4
435 real(RKG) , intent(in) :: x, alpha, beta, logBeta
436 real(RKG) :: logPDF
437 end function
438#endif
439
440#if RK3_ENABLED
441 PURE elemental module function getBetaLogPDFL_RK3(x, alpha, beta, logBeta) result(logPDF)
442#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
443 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFL_RK3
444#endif
445 use pm_kind, only: RKG => RK3
446 real(RKG) , intent(in) :: x, alpha, beta, logBeta
447 real(RKG) :: logPDF
448 end function
449#endif
450
451#if RK2_ENABLED
452 PURE elemental module function getBetaLogPDFL_RK2(x, alpha, beta, logBeta) result(logPDF)
453#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
454 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFL_RK2
455#endif
456 use pm_kind, only: RKG => RK2
457 real(RKG) , intent(in) :: x, alpha, beta, logBeta
458 real(RKG) :: logPDF
459 end function
460#endif
461
462#if RK1_ENABLED
463 PURE elemental module function getBetaLogPDFL_RK1(x, alpha, beta, logBeta) result(logPDF)
464#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
465 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaLogPDFL_RK1
466#endif
467 use pm_kind, only: RKG => RK1
468 real(RKG) , intent(in) :: x, alpha, beta, logBeta
469 real(RKG) :: logPDF
470 end function
471#endif
472
473 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
474
475 end interface
476
477!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
478
543
544 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
545
546#if RK5_ENABLED
547 PURE elemental module subroutine setBetaLogPDFD_RK5(logPDF, x, alpha, beta)
548#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
549 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFD_RK5
550#endif
551 use pm_kind, only: RKG => RK5
552 real(RKG) , intent(out) :: logPDF
553 real(RKG) , intent(in) :: x, alpha, beta
554 end subroutine
555#endif
556
557#if RK4_ENABLED
558 PURE elemental module subroutine setBetaLogPDFD_RK4(logPDF, x, alpha, beta)
559#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
560 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFD_RK4
561#endif
562 use pm_kind, only: RKG => RK4
563 real(RKG) , intent(out) :: logPDF
564 real(RKG) , intent(in) :: x, alpha, beta
565 end subroutine
566#endif
567
568#if RK3_ENABLED
569 PURE elemental module subroutine setBetaLogPDFD_RK3(logPDF, x, alpha, beta)
570#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
571 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFD_RK3
572#endif
573 use pm_kind, only: RKG => RK3
574 real(RKG) , intent(out) :: logPDF
575 real(RKG) , intent(in) :: x, alpha, beta
576 end subroutine
577#endif
578
579#if RK2_ENABLED
580 PURE elemental module subroutine setBetaLogPDFD_RK2(logPDF, x, alpha, beta)
581#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
582 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFD_RK2
583#endif
584 use pm_kind, only: RKG => RK2
585 real(RKG) , intent(out) :: logPDF
586 real(RKG) , intent(in) :: x, alpha, beta
587 end subroutine
588#endif
589
590#if RK1_ENABLED
591 PURE elemental module subroutine setBetaLogPDFD_RK1(logPDF, x, alpha, beta)
592#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
593 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFD_RK1
594#endif
595 use pm_kind, only: RKG => RK1
596 real(RKG) , intent(out) :: logPDF
597 real(RKG) , intent(in) :: x, alpha, beta
598 end subroutine
599#endif
600
601 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
602
603#if RK5_ENABLED
604 PURE elemental module subroutine setBetaLogPDFL_RK5(logPDF, x, alpha, beta, logBeta)
605#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
606 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFL_RK5
607#endif
608 use pm_kind, only: RKG => RK5
609 real(RKG) , intent(out) :: logPDF
610 real(RKG) , intent(in) :: x, alpha, beta, logBeta
611 end subroutine
612#endif
613
614#if RK4_ENABLED
615 PURE elemental module subroutine setBetaLogPDFL_RK4(logPDF, x, alpha, beta, logBeta)
616#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
617 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFL_RK4
618#endif
619 use pm_kind, only: RKG => RK4
620 real(RKG) , intent(out) :: logPDF
621 real(RKG) , intent(in) :: x, alpha, beta, logBeta
622 end subroutine
623#endif
624
625#if RK3_ENABLED
626 PURE elemental module subroutine setBetaLogPDFL_RK3(logPDF, x, alpha, beta, logBeta)
627#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
628 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFL_RK3
629#endif
630 use pm_kind, only: RKG => RK3
631 real(RKG) , intent(out) :: logPDF
632 real(RKG) , intent(in) :: x, alpha, beta, logBeta
633 end subroutine
634#endif
635
636#if RK2_ENABLED
637 PURE elemental module subroutine setBetaLogPDFL_RK2(logPDF, x, alpha, beta, logBeta)
638#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
639 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFL_RK2
640#endif
641 use pm_kind, only: RKG => RK2
642 real(RKG) , intent(out) :: logPDF
643 real(RKG) , intent(in) :: x, alpha, beta, logBeta
644 end subroutine
645#endif
646
647#if RK1_ENABLED
648 PURE elemental module subroutine setBetaLogPDFL_RK1(logPDF, x, alpha, beta, logBeta)
649#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
650 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaLogPDFL_RK1
651#endif
652 use pm_kind, only: RKG => RK1
653 real(RKG) , intent(out) :: logPDF
654 real(RKG) , intent(in) :: x, alpha, beta, logBeta
655 end subroutine
656#endif
657
658 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
659
660 end interface
661
662!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
663
732 interface getBetaCDF
733
734 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
735
736#if RK5_ENABLED
737 impure elemental module function getBetaCDF_RK5(x, alpha, beta, signed) result(cdf)
738#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
739 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaCDF_RK5
740#endif
741 use pm_kind, only: RKG => RK5
742 logical(LK) , intent(in), optional :: signed
743 real(RKG) , intent(in) :: x, alpha, beta
744 real(RKG) :: cdf
745 end function
746#endif
747
748#if RK4_ENABLED
749 impure elemental module function getBetaCDF_RK4(x, alpha, beta, signed) result(cdf)
750#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
751 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaCDF_RK4
752#endif
753 use pm_kind, only: RKG => RK4
754 logical(LK) , intent(in), optional :: signed
755 real(RKG) , intent(in) :: x, alpha, beta
756 real(RKG) :: cdf
757 end function
758#endif
759
760#if RK3_ENABLED
761 impure elemental module function getBetaCDF_RK3(x, alpha, beta, signed) result(cdf)
762#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
763 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaCDF_RK3
764#endif
765 use pm_kind, only: RKG => RK3
766 logical(LK) , intent(in), optional :: signed
767 real(RKG) , intent(in) :: x, alpha, beta
768 real(RKG) :: cdf
769 end function
770#endif
771
772#if RK2_ENABLED
773 impure elemental module function getBetaCDF_RK2(x, alpha, beta, signed) result(cdf)
774#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
775 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaCDF_RK2
776#endif
777 use pm_kind, only: RKG => RK2
778 logical(LK) , intent(in), optional :: signed
779 real(RKG) , intent(in) :: x, alpha, beta
780 real(RKG) :: cdf
781 end function
782#endif
783
784#if RK1_ENABLED
785 impure elemental module function getBetaCDF_RK1(x, alpha, beta, signed) result(cdf)
786#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
787 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaCDF_RK1
788#endif
789 use pm_kind, only: RKG => RK1
790 logical(LK) , intent(in), optional :: signed
791 real(RKG) , intent(in) :: x, alpha, beta
792 real(RKG) :: cdf
793 end function
794#endif
795
796 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
797
798 end interface
799
800!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
801
881 interface setBetaCDF
882
883 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
884
885#if RK5_ENABLED
886 PURE elemental module subroutine setBetaCDF_RK5(cdf, x, alpha, beta, logFuncBeta, signed, info)
887#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
888 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaCDF_RK5
889#endif
890 use pm_kind, only: RKG => RK5
891 real(RKG) , intent(out) :: cdf
892 real(RKG) , intent(in) :: x, alpha, beta, logFuncBeta
893 logical(LK) , intent(in) :: signed
894 integer(IK) , intent(out) :: info
895 end subroutine
896#endif
897
898#if RK4_ENABLED
899 PURE elemental module subroutine setBetaCDF_RK4(cdf, x, alpha, beta, logFuncBeta, signed, info)
900#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
901 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaCDF_RK4
902#endif
903 use pm_kind, only: RKG => RK4
904 real(RKG) , intent(out) :: cdf
905 real(RKG) , intent(in) :: x, alpha, beta, logFuncBeta
906 logical(LK) , intent(in) :: signed
907 integer(IK) , intent(out) :: info
908 end subroutine
909#endif
910
911#if RK3_ENABLED
912 PURE elemental module subroutine setBetaCDF_RK3(cdf, x, alpha, beta, logFuncBeta, signed, info)
913#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
914 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaCDF_RK3
915#endif
916 use pm_kind, only: RKG => RK3
917 real(RKG) , intent(out) :: cdf
918 real(RKG) , intent(in) :: x, alpha, beta, logFuncBeta
919 logical(LK) , intent(in) :: signed
920 integer(IK) , intent(out) :: info
921 end subroutine
922#endif
923
924#if RK2_ENABLED
925 PURE elemental module subroutine setBetaCDF_RK2(cdf, x, alpha, beta, logFuncBeta, signed, info)
926#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
927 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaCDF_RK2
928#endif
929 use pm_kind, only: RKG => RK2
930 real(RKG) , intent(out) :: cdf
931 real(RKG) , intent(in) :: x, alpha, beta, logFuncBeta
932 logical(LK) , intent(in) :: signed
933 integer(IK) , intent(out) :: info
934 end subroutine
935#endif
936
937#if RK1_ENABLED
938 PURE elemental module subroutine setBetaCDF_RK1(cdf, x, alpha, beta, logFuncBeta, signed, info)
939#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
940 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaCDF_RK1
941#endif
942 use pm_kind, only: RKG => RK1
943 real(RKG) , intent(out) :: cdf
944 real(RKG) , intent(in) :: x, alpha, beta, logFuncBeta
945 logical(LK) , intent(in) :: signed
946 integer(IK) , intent(out) :: info
947 end subroutine
948#endif
949
950 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
951
952 end interface
953
954!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
955
1039 interface setBetaRand
1040
1041 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1042 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1043 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1044
1045 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1046
1047#if RK5_ENABLED
1048 impure elemental module subroutine setBetaRandRNGD_D0_RK5(rand, alpha, beta)
1049#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1050 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D0_RK5
1051#endif
1052 use pm_kind, only: RKG => RK5
1053 real(RKG) , intent(out) :: rand
1054 real(RKG) , intent(in) :: alpha, beta
1055 end subroutine
1056#endif
1057
1058#if RK4_ENABLED
1059 impure elemental module subroutine setBetaRandRNGD_D0_RK4(rand, alpha, beta)
1060#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1061 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D0_RK4
1062#endif
1063 use pm_kind, only: RKG => RK4
1064 real(RKG) , intent(out) :: rand
1065 real(RKG) , intent(in) :: alpha, beta
1066 end subroutine
1067#endif
1068
1069#if RK3_ENABLED
1070 impure elemental module subroutine setBetaRandRNGD_D0_RK3(rand, alpha, beta)
1071#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1072 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D0_RK3
1073#endif
1074 use pm_kind, only: RKG => RK3
1075 real(RKG) , intent(out) :: rand
1076 real(RKG) , intent(in) :: alpha, beta
1077 end subroutine
1078#endif
1079
1080#if RK2_ENABLED
1081 impure elemental module subroutine setBetaRandRNGD_D0_RK2(rand, alpha, beta)
1082#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1083 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D0_RK2
1084#endif
1085 use pm_kind, only: RKG => RK2
1086 real(RKG) , intent(out) :: rand
1087 real(RKG) , intent(in) :: alpha, beta
1088 end subroutine
1089#endif
1090
1091#if RK1_ENABLED
1092 impure elemental module subroutine setBetaRandRNGD_D0_RK1(rand, alpha, beta)
1093#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1094 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D0_RK1
1095#endif
1096 use pm_kind, only: RKG => RK1
1097 real(RKG) , intent(out) :: rand
1098 real(RKG) , intent(in) :: alpha, beta
1099 end subroutine
1100#endif
1101
1102 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1103
1104#if RK5_ENABLED
1105 impure elemental module subroutine setBetaRandRNGF_D0_RK5(rng, rand, alpha, beta)
1106#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1107 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D0_RK5
1108#endif
1109 use pm_kind, only: RKG => RK5
1110 type(rngf_type) , intent(in) :: rng
1111 real(RKG) , intent(out) :: rand
1112 real(RKG) , intent(in) :: alpha, beta
1113 end subroutine
1114#endif
1115
1116#if RK4_ENABLED
1117 impure elemental module subroutine setBetaRandRNGF_D0_RK4(rng, rand, alpha, beta)
1118#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1119 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D0_RK4
1120#endif
1121 use pm_kind, only: RKG => RK4
1122 type(rngf_type) , intent(in) :: rng
1123 real(RKG) , intent(out) :: rand
1124 real(RKG) , intent(in) :: alpha, beta
1125 end subroutine
1126#endif
1127
1128#if RK3_ENABLED
1129 impure elemental module subroutine setBetaRandRNGF_D0_RK3(rng, rand, alpha, beta)
1130#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1131 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D0_RK3
1132#endif
1133 use pm_kind, only: RKG => RK3
1134 type(rngf_type) , intent(in) :: rng
1135 real(RKG) , intent(out) :: rand
1136 real(RKG) , intent(in) :: alpha, beta
1137 end subroutine
1138#endif
1139
1140#if RK2_ENABLED
1141 impure elemental module subroutine setBetaRandRNGF_D0_RK2(rng, rand, alpha, beta)
1142#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1143 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D0_RK2
1144#endif
1145 use pm_kind, only: RKG => RK2
1146 type(rngf_type) , intent(in) :: rng
1147 real(RKG) , intent(out) :: rand
1148 real(RKG) , intent(in) :: alpha, beta
1149 end subroutine
1150#endif
1151
1152#if RK1_ENABLED
1153 impure elemental module subroutine setBetaRandRNGF_D0_RK1(rng, rand, alpha, beta)
1154#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1155 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D0_RK1
1156#endif
1157 use pm_kind, only: RKG => RK1
1158 type(rngf_type) , intent(in) :: rng
1159 real(RKG) , intent(out) :: rand
1160 real(RKG) , intent(in) :: alpha, beta
1161 end subroutine
1162#endif
1163
1164 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1165
1166#if RK5_ENABLED
1167 PURE module subroutine setBetaRandRNGX_D0_RK5(rng, rand, alpha, beta)
1168#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1169 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D0_RK5
1170#endif
1171 use pm_kind, only: RKG => RK5
1172 type(xoshiro256ssw_type), intent(inout) :: rng
1173 real(RKG) , intent(out) :: rand
1174 real(RKG) , intent(in) :: alpha, beta
1175 end subroutine
1176#endif
1177
1178#if RK4_ENABLED
1179 PURE module subroutine setBetaRandRNGX_D0_RK4(rng, rand, alpha, beta)
1180#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1181 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D0_RK4
1182#endif
1183 use pm_kind, only: RKG => RK4
1184 type(xoshiro256ssw_type), intent(inout) :: rng
1185 real(RKG) , intent(out) :: rand
1186 real(RKG) , intent(in) :: alpha, beta
1187 end subroutine
1188#endif
1189
1190#if RK3_ENABLED
1191 PURE module subroutine setBetaRandRNGX_D0_RK3(rng, rand, alpha, beta)
1192#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1193 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D0_RK3
1194#endif
1195 use pm_kind, only: RKG => RK3
1196 type(xoshiro256ssw_type), intent(inout) :: rng
1197 real(RKG) , intent(out) :: rand
1198 real(RKG) , intent(in) :: alpha, beta
1199 end subroutine
1200#endif
1201
1202#if RK2_ENABLED
1203 PURE module subroutine setBetaRandRNGX_D0_RK2(rng, rand, alpha, beta)
1204#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1205 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D0_RK2
1206#endif
1207 use pm_kind, only: RKG => RK2
1208 type(xoshiro256ssw_type), intent(inout) :: rng
1209 real(RKG) , intent(out) :: rand
1210 real(RKG) , intent(in) :: alpha, beta
1211 end subroutine
1212#endif
1213
1214#if RK1_ENABLED
1215 PURE module subroutine setBetaRandRNGX_D0_RK1(rng, rand, alpha, beta)
1216#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1217 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D0_RK1
1218#endif
1219 use pm_kind, only: RKG => RK1
1220 type(xoshiro256ssw_type), intent(inout) :: rng
1221 real(RKG) , intent(out) :: rand
1222 real(RKG) , intent(in) :: alpha, beta
1223 end subroutine
1224#endif
1225
1226 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1227
1228 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1229 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1230 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1231
1232 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1233
1234#if RK5_ENABLED
1235 impure module subroutine setBetaRandRNGD_D1_RK5(rand, alpha, beta)
1236#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1237 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D1_RK5
1238#endif
1239 use pm_kind, only: RKG => RK5
1240 real(RKG) , intent(out) :: rand(:)
1241 real(RKG) , intent(in) :: alpha, beta
1242 end subroutine
1243#endif
1244
1245#if RK4_ENABLED
1246 impure module subroutine setBetaRandRNGD_D1_RK4(rand, alpha, beta)
1247#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1248 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D1_RK4
1249#endif
1250 use pm_kind, only: RKG => RK4
1251 real(RKG) , intent(out) :: rand(:)
1252 real(RKG) , intent(in) :: alpha, beta
1253 end subroutine
1254#endif
1255
1256#if RK3_ENABLED
1257 impure module subroutine setBetaRandRNGD_D1_RK3(rand, alpha, beta)
1258#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1259 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D1_RK3
1260#endif
1261 use pm_kind, only: RKG => RK3
1262 real(RKG) , intent(out) :: rand(:)
1263 real(RKG) , intent(in) :: alpha, beta
1264 end subroutine
1265#endif
1266
1267#if RK2_ENABLED
1268 impure module subroutine setBetaRandRNGD_D1_RK2(rand, alpha, beta)
1269#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1270 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D1_RK2
1271#endif
1272 use pm_kind, only: RKG => RK2
1273 real(RKG) , intent(out) :: rand(:)
1274 real(RKG) , intent(in) :: alpha, beta
1275 end subroutine
1276#endif
1277
1278#if RK1_ENABLED
1279 impure module subroutine setBetaRandRNGD_D1_RK1(rand, alpha, beta)
1280#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1281 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGD_D1_RK1
1282#endif
1283 use pm_kind, only: RKG => RK1
1284 real(RKG) , intent(out) :: rand(:)
1285 real(RKG) , intent(in) :: alpha, beta
1286 end subroutine
1287#endif
1288
1289 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1290
1291#if RK5_ENABLED
1292 impure module subroutine setBetaRandRNGF_D1_RK5(rng, rand, alpha, beta)
1293#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1294 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D1_RK5
1295#endif
1296 use pm_kind, only: RKG => RK5
1297 type(rngf_type) , intent(in) :: rng
1298 real(RKG) , intent(out) :: rand(:)
1299 real(RKG) , intent(in) :: alpha, beta
1300 end subroutine
1301#endif
1302
1303#if RK4_ENABLED
1304 impure module subroutine setBetaRandRNGF_D1_RK4(rng, rand, alpha, beta)
1305#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1306 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D1_RK4
1307#endif
1308 use pm_kind, only: RKG => RK4
1309 type(rngf_type) , intent(in) :: rng
1310 real(RKG) , intent(out) :: rand(:)
1311 real(RKG) , intent(in) :: alpha, beta
1312 end subroutine
1313#endif
1314
1315#if RK3_ENABLED
1316 impure module subroutine setBetaRandRNGF_D1_RK3(rng, rand, alpha, beta)
1317#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1318 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D1_RK3
1319#endif
1320 use pm_kind, only: RKG => RK3
1321 type(rngf_type) , intent(in) :: rng
1322 real(RKG) , intent(out) :: rand(:)
1323 real(RKG) , intent(in) :: alpha, beta
1324 end subroutine
1325#endif
1326
1327#if RK2_ENABLED
1328 impure module subroutine setBetaRandRNGF_D1_RK2(rng, rand, alpha, beta)
1329#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1330 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D1_RK2
1331#endif
1332 use pm_kind, only: RKG => RK2
1333 type(rngf_type) , intent(in) :: rng
1334 real(RKG) , intent(out) :: rand(:)
1335 real(RKG) , intent(in) :: alpha, beta
1336 end subroutine
1337#endif
1338
1339#if RK1_ENABLED
1340 impure module subroutine setBetaRandRNGF_D1_RK1(rng, rand, alpha, beta)
1341#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1342 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGF_D1_RK1
1343#endif
1344 use pm_kind, only: RKG => RK1
1345 type(rngf_type) , intent(in) :: rng
1346 real(RKG) , intent(out) :: rand(:)
1347 real(RKG) , intent(in) :: alpha, beta
1348 end subroutine
1349#endif
1350
1351 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1352
1353#if RK5_ENABLED
1354 PURE module subroutine setBetaRandRNGX_D1_RK5(rng, rand, alpha, beta)
1355#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1356 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D1_RK5
1357#endif
1358 use pm_kind, only: RKG => RK5
1359 type(xoshiro256ssw_type), intent(inout) :: rng
1360 real(RKG) , intent(out) :: rand(:)
1361 real(RKG) , intent(in) :: alpha, beta
1362 end subroutine
1363#endif
1364
1365#if RK4_ENABLED
1366 PURE module subroutine setBetaRandRNGX_D1_RK4(rng, rand, alpha, beta)
1367#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1368 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D1_RK4
1369#endif
1370 use pm_kind, only: RKG => RK4
1371 type(xoshiro256ssw_type), intent(inout) :: rng
1372 real(RKG) , intent(out) :: rand(:)
1373 real(RKG) , intent(in) :: alpha, beta
1374 end subroutine
1375#endif
1376
1377#if RK3_ENABLED
1378 PURE module subroutine setBetaRandRNGX_D1_RK3(rng, rand, alpha, beta)
1379#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1380 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D1_RK3
1381#endif
1382 use pm_kind, only: RKG => RK3
1383 type(xoshiro256ssw_type), intent(inout) :: rng
1384 real(RKG) , intent(out) :: rand(:)
1385 real(RKG) , intent(in) :: alpha, beta
1386 end subroutine
1387#endif
1388
1389#if RK2_ENABLED
1390 PURE module subroutine setBetaRandRNGX_D1_RK2(rng, rand, alpha, beta)
1391#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1392 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D1_RK2
1393#endif
1394 use pm_kind, only: RKG => RK2
1395 type(xoshiro256ssw_type), intent(inout) :: rng
1396 real(RKG) , intent(out) :: rand(:)
1397 real(RKG) , intent(in) :: alpha, beta
1398 end subroutine
1399#endif
1400
1401#if RK1_ENABLED
1402 PURE module subroutine setBetaRandRNGX_D1_RK1(rng, rand, alpha, beta)
1403#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1404 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaRandRNGX_D1_RK1
1405#endif
1406 use pm_kind, only: RKG => RK1
1407 type(xoshiro256ssw_type), intent(inout) :: rng
1408 real(RKG) , intent(out) :: rand(:)
1409 real(RKG) , intent(in) :: alpha, beta
1410 end subroutine
1411#endif
1412
1413 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1414
1415 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1416 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1417 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1418
1419 end interface
1420
1421!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1422
1423end module pm_distBeta
Generate and return the CDF of the Beta distribution for the given parameters as defined in the deta...
Generate and return the natural logarithm of the Probability Density Function (PDF) of the Beta distr...
Generate and return the Probability Density Function (PDF) of the Beta distribution for an input x wi...
Return the CDF of the Beta distribution for the given parameters as defined in the details section o...
Return the natural logarithm of the Probability Density Function (PDF) of the Beta distribution for a...
Return a scalar or array of arbitrary rank of Beta-distributed random values in range (or ,...
This module contains classes and procedures for computing various statistical quantities related to t...
Definition: pm_distBeta.F90:99
character(*, SK), parameter MODULE_NAME
This module contains classes and procedures for computing various statistical quantities related to t...
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 is the derived type for signifying distributions that are of type Beta as defined in the descrip...
This is a concrete derived type whose instances can be used to define/request the default uniform ran...
This is the derived type for declaring and generating objects of type xoshiro256ssw_type containing a...