ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_distPoweto.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
64
65!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
66
68
69 use pm_kind, only: SK, IK
70
71 implicit none
72
73 character(*, SK), parameter :: MODULE_NAME = "@pm_distPoweto"
74
75!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
76
111 end type
112
113!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
114
187
188 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
189
190#if RK5_ENABLED
191 PURE elemental module function getPowetoLogPDFNF_D0_RK5(alpha, logMinX, logMaxX) result(logPDFNF)
192#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
193 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDFNF_D0_RK5
194#endif
195 use pm_kind, only: RKG => RK5
196 real(RKG) , intent(in) :: alpha
197 real(RKG) , intent(in) , optional :: logMinX, logMaxX
198 real(RKG) :: logPDFNF
199 end function
200#endif
201
202#if RK4_ENABLED
203 PURE elemental module function getPowetoLogPDFNF_D0_RK4(alpha, logMinX, logMaxX) result(logPDFNF)
204#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
205 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDFNF_D0_RK4
206#endif
207 use pm_kind, only: RKG => RK4
208 real(RKG) , intent(in) :: alpha
209 real(RKG) , intent(in) , optional :: logMinX, logMaxX
210 real(RKG) :: logPDFNF
211 end function
212#endif
213
214#if RK3_ENABLED
215 PURE elemental module function getPowetoLogPDFNF_D0_RK3(alpha, logMinX, logMaxX) result(logPDFNF)
216#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
217 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDFNF_D0_RK3
218#endif
219 use pm_kind, only: RKG => RK3
220 real(RKG) , intent(in) :: alpha
221 real(RKG) , intent(in) , optional :: logMinX, logMaxX
222 real(RKG) :: logPDFNF
223 end function
224#endif
225
226#if RK2_ENABLED
227 PURE elemental module function getPowetoLogPDFNF_D0_RK2(alpha, logMinX, logMaxX) result(logPDFNF)
228#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
229 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDFNF_D0_RK2
230#endif
231 use pm_kind, only: RKG => RK2
232 real(RKG) , intent(in) :: alpha
233 real(RKG) , intent(in) , optional :: logMinX, logMaxX
234 real(RKG) :: logPDFNF
235 end function
236#endif
237
238#if RK1_ENABLED
239 PURE elemental module function getPowetoLogPDFNF_D0_RK1(alpha, logMinX, logMaxX) result(logPDFNF)
240#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
241 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDFNF_D0_RK1
242#endif
243 use pm_kind, only: RKG => RK1
244 real(RKG) , intent(in) :: alpha
245 real(RKG) , intent(in) , optional :: logMinX, logMaxX
246 real(RKG) :: logPDFNF
247 end function
248#endif
249
250 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
251
252 end interface
253
254!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
255
323
324 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
325
326#if RK5_ENABLED
327 PURE elemental module function getPowetoLogPDF_D0_RK5(logx, alpha, logMinX, logMaxX) result(logPDF)
328#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
329 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDF_D0_RK5
330#endif
331 use pm_kind, only: RKG => RK5
332 real(RKG) , intent(in) :: logx
333 real(RKG) , intent(in) :: alpha
334 real(RKG) , intent(in) , optional :: logMinX, logMaxX
335 real(RKG) :: logPDF
336 end function
337#endif
338
339#if RK4_ENABLED
340 PURE elemental module function getPowetoLogPDF_D0_RK4(logx, alpha, logMinX, logMaxX) result(logPDF)
341#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
342 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDF_D0_RK4
343#endif
344 use pm_kind, only: RKG => RK4
345 real(RKG) , intent(in) :: logx
346 real(RKG) , intent(in) :: alpha
347 real(RKG) , intent(in) , optional :: logMinX, logMaxX
348 real(RKG) :: logPDF
349 end function
350#endif
351
352#if RK3_ENABLED
353 PURE elemental module function getPowetoLogPDF_D0_RK3(logx, alpha, logMinX, logMaxX) result(logPDF)
354#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
355 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDF_D0_RK3
356#endif
357 use pm_kind, only: RKG => RK3
358 real(RKG) , intent(in) :: logx
359 real(RKG) , intent(in) :: alpha
360 real(RKG) , intent(in) , optional :: logMinX, logMaxX
361 real(RKG) :: logPDF
362 end function
363#endif
364
365#if RK2_ENABLED
366 PURE elemental module function getPowetoLogPDF_D0_RK2(logx, alpha, logMinX, logMaxX) result(logPDF)
367#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
368 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDF_D0_RK2
369#endif
370 use pm_kind, only: RKG => RK2
371 real(RKG) , intent(in) :: logx
372 real(RKG) , intent(in) :: alpha
373 real(RKG) , intent(in) , optional :: logMinX, logMaxX
374 real(RKG) :: logPDF
375 end function
376#endif
377
378#if RK1_ENABLED
379 PURE elemental module function getPowetoLogPDF_D0_RK1(logx, alpha, logMinX, logMaxX) result(logPDF)
380#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
381 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogPDF_D0_RK1
382#endif
383 use pm_kind, only: RKG => RK1
384 real(RKG) , intent(in) :: logx
385 real(RKG) , intent(in) :: alpha
386 real(RKG) , intent(in) , optional :: logMinX, logMaxX
387 real(RKG) :: logPDF
388 end function
389#endif
390
391 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
392
393 end interface
394
395!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
396
454
455 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
456
457#if RK5_ENABLED
458 PURE elemental module subroutine setPowetoLogPDF_D0_RK5(logPDF, logx, alpha, logPDFNF)
459#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
460 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogPDF_D0_RK5
461#endif
462 use pm_kind, only: RKG => RK5
463 real(RKG) , intent(out) :: logPDF
464 real(RKG) , intent(in) :: logx
465 real(RKG) , intent(in) :: alpha
466 real(RKG) , intent(in) :: logPDFNF
467 end subroutine
468#endif
469
470#if RK4_ENABLED
471 PURE elemental module subroutine setPowetoLogPDF_D0_RK4(logPDF, logx, alpha, logPDFNF)
472#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
473 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogPDF_D0_RK4
474#endif
475 use pm_kind, only: RKG => RK4
476 real(RKG) , intent(out) :: logPDF
477 real(RKG) , intent(in) :: logx
478 real(RKG) , intent(in) :: alpha
479 real(RKG) , intent(in) :: logPDFNF
480 end subroutine
481#endif
482
483#if RK3_ENABLED
484 PURE elemental module subroutine setPowetoLogPDF_D0_RK3(logPDF, logx, alpha, logPDFNF)
485#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
486 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogPDF_D0_RK3
487#endif
488 use pm_kind, only: RKG => RK3
489 real(RKG) , intent(out) :: logPDF
490 real(RKG) , intent(in) :: logx
491 real(RKG) , intent(in) :: alpha
492 real(RKG) , intent(in) :: logPDFNF
493 end subroutine
494#endif
495
496#if RK2_ENABLED
497 PURE elemental module subroutine setPowetoLogPDF_D0_RK2(logPDF, logx, alpha, logPDFNF)
498#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
499 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogPDF_D0_RK2
500#endif
501 use pm_kind, only: RKG => RK2
502 real(RKG) , intent(out) :: logPDF
503 real(RKG) , intent(in) :: logx
504 real(RKG) , intent(in) :: alpha
505 real(RKG) , intent(in) :: logPDFNF
506 end subroutine
507#endif
508
509#if RK1_ENABLED
510 PURE elemental module subroutine setPowetoLogPDF_D0_RK1(logPDF, logx, alpha, logPDFNF)
511#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
512 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogPDF_D0_RK1
513#endif
514 use pm_kind, only: RKG => RK1
515 real(RKG) , intent(out) :: logPDF
516 real(RKG) , intent(in) :: logx
517 real(RKG) , intent(in) :: alpha
518 real(RKG) , intent(in) :: logPDFNF
519 end subroutine
520#endif
521
522 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
523
524 end interface
525
526!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
527
601
602 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
603
604#if RK5_ENABLED
605 PURE elemental module function getPowetoLogCDFNF_D0_RK5(alpha, logMinX, logMaxX) result(logCDFNF)
606#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
607 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDFNF_D0_RK5
608#endif
609 use pm_kind, only: RKG => RK5
610 real(RKG) , intent(in) :: alpha
611 real(RKG) , intent(in) , optional :: logMinX, logMaxX
612 real(RKG) :: logCDFNF
613 end function
614#endif
615
616#if RK4_ENABLED
617 PURE elemental module function getPowetoLogCDFNF_D0_RK4(alpha, logMinX, logMaxX) result(logCDFNF)
618#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
619 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDFNF_D0_RK4
620#endif
621 use pm_kind, only: RKG => RK4
622 real(RKG) , intent(in) :: alpha
623 real(RKG) , intent(in) , optional :: logMinX, logMaxX
624 real(RKG) :: logCDFNF
625 end function
626#endif
627
628#if RK3_ENABLED
629 PURE elemental module function getPowetoLogCDFNF_D0_RK3(alpha, logMinX, logMaxX) result(logCDFNF)
630#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
631 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDFNF_D0_RK3
632#endif
633 use pm_kind, only: RKG => RK3
634 real(RKG) , intent(in) :: alpha
635 real(RKG) , intent(in) , optional :: logMinX, logMaxX
636 real(RKG) :: logCDFNF
637 end function
638#endif
639
640#if RK2_ENABLED
641 PURE elemental module function getPowetoLogCDFNF_D0_RK2(alpha, logMinX, logMaxX) result(logCDFNF)
642#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
643 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDFNF_D0_RK2
644#endif
645 use pm_kind, only: RKG => RK2
646 real(RKG) , intent(in) :: alpha
647 real(RKG) , intent(in) , optional :: logMinX, logMaxX
648 real(RKG) :: logCDFNF
649 end function
650#endif
651
652#if RK1_ENABLED
653 PURE elemental module function getPowetoLogCDFNF_D0_RK1(alpha, logMinX, logMaxX) result(logCDFNF)
654#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
655 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDFNF_D0_RK1
656#endif
657 use pm_kind, only: RKG => RK1
658 real(RKG) , intent(in) :: alpha
659 real(RKG) , intent(in) , optional :: logMinX, logMaxX
660 real(RKG) :: logCDFNF
661 end function
662#endif
663
664 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
665
666 end interface
667
668!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
669
737
738 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
739
740#if RK5_ENABLED
741 PURE elemental module function getPowetoLogCDF_D0_RK5(logx, alpha, logMinX, logMaxX) result(logCDF)
742#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
743 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDF_D0_RK5
744#endif
745 use pm_kind, only: RKG => RK5
746 real(RKG) , intent(in) :: logx
747 real(RKG) , intent(in) :: alpha
748 real(RKG) , intent(in) , optional :: logMinX, logMaxX
749 real(RKG) :: logCDF
750 end function
751#endif
752
753#if RK4_ENABLED
754 PURE elemental module function getPowetoLogCDF_D0_RK4(logx, alpha, logMinX, logMaxX) result(logCDF)
755#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
756 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDF_D0_RK4
757#endif
758 use pm_kind, only: RKG => RK4
759 real(RKG) , intent(in) :: logx
760 real(RKG) , intent(in) :: alpha
761 real(RKG) , intent(in) , optional :: logMinX, logMaxX
762 real(RKG) :: logCDF
763 end function
764#endif
765
766#if RK3_ENABLED
767 PURE elemental module function getPowetoLogCDF_D0_RK3(logx, alpha, logMinX, logMaxX) result(logCDF)
768#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
769 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDF_D0_RK3
770#endif
771 use pm_kind, only: RKG => RK3
772 real(RKG) , intent(in) :: logx
773 real(RKG) , intent(in) :: alpha
774 real(RKG) , intent(in) , optional :: logMinX, logMaxX
775 real(RKG) :: logCDF
776 end function
777#endif
778
779#if RK2_ENABLED
780 PURE elemental module function getPowetoLogCDF_D0_RK2(logx, alpha, logMinX, logMaxX) result(logCDF)
781#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
782 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDF_D0_RK2
783#endif
784 use pm_kind, only: RKG => RK2
785 real(RKG) , intent(in) :: logx
786 real(RKG) , intent(in) :: alpha
787 real(RKG) , intent(in) , optional :: logMinX, logMaxX
788 real(RKG) :: logCDF
789 end function
790#endif
791
792#if RK1_ENABLED
793 PURE elemental module function getPowetoLogCDF_D0_RK1(logx, alpha, logMinX, logMaxX) result(logCDF)
794#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
795 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogCDF_D0_RK1
796#endif
797 use pm_kind, only: RKG => RK1
798 real(RKG) , intent(in) :: logx
799 real(RKG) , intent(in) :: alpha
800 real(RKG) , intent(in) , optional :: logMinX, logMaxX
801 real(RKG) :: logCDF
802 end function
803#endif
804
805 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
806
807 end interface
808
809!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
810
878
879 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
880
881#if RK5_ENABLED
882 PURE elemental module subroutine setPowetoLogCDF_D0_RK5(logCDF, logx, alpha, logMinX, logCDFNF)
883#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
884 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogCDF_D0_RK5
885#endif
886 use pm_kind, only: RKG => RK5
887 real(RKG) , intent(out) :: logCDF
888 real(RKG) , intent(in) :: logx
889 real(RKG) , intent(in) :: alpha
890 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
891 end subroutine
892#endif
893
894#if RK4_ENABLED
895 PURE elemental module subroutine setPowetoLogCDF_D0_RK4(logCDF, logx, alpha, logMinX, logCDFNF)
896#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
897 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogCDF_D0_RK4
898#endif
899 use pm_kind, only: RKG => RK4
900 real(RKG) , intent(out) :: logCDF
901 real(RKG) , intent(in) :: logx
902 real(RKG) , intent(in) :: alpha
903 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
904 end subroutine
905#endif
906
907#if RK3_ENABLED
908 PURE elemental module subroutine setPowetoLogCDF_D0_RK3(logCDF, logx, alpha, logMinX, logCDFNF)
909#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
910 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogCDF_D0_RK3
911#endif
912 use pm_kind, only: RKG => RK3
913 real(RKG) , intent(out) :: logCDF
914 real(RKG) , intent(in) :: logx
915 real(RKG) , intent(in) :: alpha
916 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
917 end subroutine
918#endif
919
920#if RK2_ENABLED
921 PURE elemental module subroutine setPowetoLogCDF_D0_RK2(logCDF, logx, alpha, logMinX, logCDFNF)
922#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
923 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogCDF_D0_RK2
924#endif
925 use pm_kind, only: RKG => RK2
926 real(RKG) , intent(out) :: logCDF
927 real(RKG) , intent(in) :: logx
928 real(RKG) , intent(in) :: alpha
929 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
930 end subroutine
931#endif
932
933#if RK1_ENABLED
934 PURE elemental module subroutine setPowetoLogCDF_D0_RK1(logCDF, logx, alpha, logMinX, logCDFNF)
935#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
936 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogCDF_D0_RK1
937#endif
938 use pm_kind, only: RKG => RK1
939 real(RKG) , intent(out) :: logCDF
940 real(RKG) , intent(in) :: logx
941 real(RKG) , intent(in) :: alpha
942 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
943 end subroutine
944#endif
945
946 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
947
948 end interface
949
950!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
951
1019
1020 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1021
1022#if RK5_ENABLED
1023 PURE elemental module function getPowetoLogQuan_D0_RK5(logCDF, alpha, logMinX, logMaxX) result(logx)
1024#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1025 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogQuan_D0_RK5
1026#endif
1027 use pm_kind, only: RKG => RK5
1028 real(RKG) , intent(in) :: logCDF
1029 real(RKG) , intent(in) :: alpha
1030 real(RKG) , intent(in) , optional :: logMinX, logMaxX
1031 real(RKG) :: logx
1032 end function
1033#endif
1034
1035#if RK4_ENABLED
1036 PURE elemental module function getPowetoLogQuan_D0_RK4(logCDF, alpha, logMinX, logMaxX) result(logx)
1037#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1038 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogQuan_D0_RK4
1039#endif
1040 use pm_kind, only: RKG => RK4
1041 real(RKG) , intent(in) :: logCDF
1042 real(RKG) , intent(in) :: alpha
1043 real(RKG) , intent(in) , optional :: logMinX, logMaxX
1044 real(RKG) :: logx
1045 end function
1046#endif
1047
1048#if RK3_ENABLED
1049 PURE elemental module function getPowetoLogQuan_D0_RK3(logCDF, alpha, logMinX, logMaxX) result(logx)
1050#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1051 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogQuan_D0_RK3
1052#endif
1053 use pm_kind, only: RKG => RK3
1054 real(RKG) , intent(in) :: logCDF
1055 real(RKG) , intent(in) :: alpha
1056 real(RKG) , intent(in) , optional :: logMinX, logMaxX
1057 real(RKG) :: logx
1058 end function
1059#endif
1060
1061#if RK2_ENABLED
1062 PURE elemental module function getPowetoLogQuan_D0_RK2(logCDF, alpha, logMinX, logMaxX) result(logx)
1063#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1064 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogQuan_D0_RK2
1065#endif
1066 use pm_kind, only: RKG => RK2
1067 real(RKG) , intent(in) :: logCDF
1068 real(RKG) , intent(in) :: alpha
1069 real(RKG) , intent(in) , optional :: logMinX, logMaxX
1070 real(RKG) :: logx
1071 end function
1072#endif
1073
1074#if RK1_ENABLED
1075 PURE elemental module function getPowetoLogQuan_D0_RK1(logCDF, alpha, logMinX, logMaxX) result(logx)
1076#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1077 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogQuan_D0_RK1
1078#endif
1079 use pm_kind, only: RKG => RK1
1080 real(RKG) , intent(in) :: logCDF
1081 real(RKG) , intent(in) :: alpha
1082 real(RKG) , intent(in) , optional :: logMinX, logMaxX
1083 real(RKG) :: logx
1084 end function
1085#endif
1086
1087 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1088
1089 end interface
1090
1091!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1092
1162
1163 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1164
1165#if RK5_ENABLED
1166 PURE elemental module subroutine setPowetoLogQuan_D0_RK5(logx, logCDF, alpha, logMinX, logCDFNF)
1167#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1168 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogQuan_D0_RK5
1169#endif
1170 use pm_kind, only: RKG => RK5
1171 real(RKG) , intent(out) :: logx
1172 real(RKG) , intent(in) :: logCDF
1173 real(RKG) , intent(in) :: alpha
1174 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1175 end subroutine
1176#endif
1177
1178#if RK4_ENABLED
1179 PURE elemental module subroutine setPowetoLogQuan_D0_RK4(logx, logCDF, alpha, logMinX, logCDFNF)
1180#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1181 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogQuan_D0_RK4
1182#endif
1183 use pm_kind, only: RKG => RK4
1184 real(RKG) , intent(out) :: logx
1185 real(RKG) , intent(in) :: logCDF
1186 real(RKG) , intent(in) :: alpha
1187 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1188 end subroutine
1189#endif
1190
1191#if RK3_ENABLED
1192 PURE elemental module subroutine setPowetoLogQuan_D0_RK3(logx, logCDF, alpha, logMinX, logCDFNF)
1193#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1194 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogQuan_D0_RK3
1195#endif
1196 use pm_kind, only: RKG => RK3
1197 real(RKG) , intent(out) :: logx
1198 real(RKG) , intent(in) :: logCDF
1199 real(RKG) , intent(in) :: alpha
1200 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1201 end subroutine
1202#endif
1203
1204#if RK2_ENABLED
1205 PURE elemental module subroutine setPowetoLogQuan_D0_RK2(logx, logCDF, alpha, logMinX, logCDFNF)
1206#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1207 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogQuan_D0_RK2
1208#endif
1209 use pm_kind, only: RKG => RK2
1210 real(RKG) , intent(out) :: logx
1211 real(RKG) , intent(in) :: logCDF
1212 real(RKG) , intent(in) :: alpha
1213 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1214 end subroutine
1215#endif
1216
1217#if RK1_ENABLED
1218 PURE elemental module subroutine setPowetoLogQuan_D0_RK1(logx, logCDF, alpha, logMinX, logCDFNF)
1219#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1220 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogQuan_D0_RK1
1221#endif
1222 use pm_kind, only: RKG => RK1
1223 real(RKG) , intent(out) :: logx
1224 real(RKG) , intent(in) :: logCDF
1225 real(RKG) , intent(in) :: alpha
1226 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1227 end subroutine
1228#endif
1229
1230 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1231
1232 end interface
1233
1234!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1235
1304
1305 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1306
1307#if RK5_ENABLED
1308 impure elemental module function getPowetoLogRand_D0_RK5(alpha, logMinX, logMaxX, logCDFNF) result(logRand)
1309#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1310 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogRand_D0_RK5
1311#endif
1312 use pm_kind, only: RKG => RK5
1313 real(RKG) , intent(in) :: alpha
1314 real(RKG) , intent(in) , optional :: logMinX, logMaxX, logCDFNF
1315 real(RKG) :: logRand
1316 end function
1317#endif
1318
1319#if RK4_ENABLED
1320 impure elemental module function getPowetoLogRand_D0_RK4(alpha, logMinX, logMaxX, logCDFNF) result(logRand)
1321#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1322 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogRand_D0_RK4
1323#endif
1324 use pm_kind, only: RKG => RK4
1325 real(RKG) , intent(in) :: alpha
1326 real(RKG) , intent(in) , optional :: logMinX, logMaxX, logCDFNF
1327 real(RKG) :: logRand
1328 end function
1329#endif
1330
1331#if RK3_ENABLED
1332 impure elemental module function getPowetoLogRand_D0_RK3(alpha, logMinX, logMaxX, logCDFNF) result(logRand)
1333#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1334 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogRand_D0_RK3
1335#endif
1336 use pm_kind, only: RKG => RK3
1337 real(RKG) , intent(in) :: alpha
1338 real(RKG) , intent(in) , optional :: logMinX, logMaxX, logCDFNF
1339 real(RKG) :: logRand
1340 end function
1341#endif
1342
1343#if RK2_ENABLED
1344 impure elemental module function getPowetoLogRand_D0_RK2(alpha, logMinX, logMaxX, logCDFNF) result(logRand)
1345#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1346 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogRand_D0_RK2
1347#endif
1348 use pm_kind, only: RKG => RK2
1349 real(RKG) , intent(in) :: alpha
1350 real(RKG) , intent(in) , optional :: logMinX, logMaxX, logCDFNF
1351 real(RKG) :: logRand
1352 end function
1353#endif
1354
1355#if RK1_ENABLED
1356 impure elemental module function getPowetoLogRand_D0_RK1(alpha, logMinX, logMaxX, logCDFNF) result(logRand)
1357#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1358 !DEC$ ATTRIBUTES DLLEXPORT :: getPowetoLogRand_D0_RK1
1359#endif
1360 use pm_kind, only: RKG => RK1
1361 real(RKG) , intent(in) :: alpha
1362 real(RKG) , intent(in) , optional :: logMinX, logMaxX, logCDFNF
1363 real(RKG) :: logRand
1364 end function
1365#endif
1366
1367 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1368
1369 end interface
1370
1371!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1372
1446
1447 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1448
1449#if RK5_ENABLED
1450 PURE elemental module subroutine setPowetoLogRand_D0_RK5(logRand, negExpRand, alpha, logMinX, logCDFNF)
1451#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1452 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogRand_D0_RK5
1453#endif
1454 use pm_kind, only: RKG => RK5
1455 real(RKG) , intent(out) :: logRand
1456 real(RKG) , intent(in) :: negExpRand, alpha
1457 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1458 end subroutine
1459#endif
1460
1461#if RK4_ENABLED
1462 PURE elemental module subroutine setPowetoLogRand_D0_RK4(logRand, negExpRand, alpha, logMinX, logCDFNF)
1463#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1464 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogRand_D0_RK4
1465#endif
1466 use pm_kind, only: RKG => RK4
1467 real(RKG) , intent(out) :: logRand
1468 real(RKG) , intent(in) :: negExpRand, alpha
1469 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1470 end subroutine
1471#endif
1472
1473#if RK3_ENABLED
1474 PURE elemental module subroutine setPowetoLogRand_D0_RK3(logRand, negExpRand, alpha, logMinX, logCDFNF)
1475#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1476 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogRand_D0_RK3
1477#endif
1478 use pm_kind, only: RKG => RK3
1479 real(RKG) , intent(out) :: logRand
1480 real(RKG) , intent(in) :: negExpRand, alpha
1481 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1482 end subroutine
1483#endif
1484
1485#if RK2_ENABLED
1486 PURE elemental module subroutine setPowetoLogRand_D0_RK2(logRand, negExpRand, alpha, logMinX, logCDFNF)
1487#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1488 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogRand_D0_RK2
1489#endif
1490 use pm_kind, only: RKG => RK2
1491 real(RKG) , intent(out) :: logRand
1492 real(RKG) , intent(in) :: negExpRand, alpha
1493 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1494 end subroutine
1495#endif
1496
1497#if RK1_ENABLED
1498 PURE elemental module subroutine setPowetoLogRand_D0_RK1(logRand, negExpRand, alpha, logMinX, logCDFNF)
1499#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1500 !DEC$ ATTRIBUTES DLLEXPORT :: setPowetoLogRand_D0_RK1
1501#endif
1502 use pm_kind, only: RKG => RK1
1503 real(RKG) , intent(out) :: logRand
1504 real(RKG) , intent(in) :: negExpRand, alpha
1505 real(RKG) , intent(in) , optional :: logMinX, logCDFNF
1506 end subroutine
1507#endif
1508
1509 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1510
1511 end interface
1512
1513!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1514
1515end module pm_distPoweto
Generate and return the natural logarithm of the normalization factor of the Cumulative Distribution ...
Generate and return the natural logarithm of the Cumulative Distribution Function (CDF) of the (Trunc...
Generate and return the natural logarithm of the normalization factor of the Probability Density Func...
Generate and return the natural logarithm of the Probability Density Function (PDF) of the (Truncated...
Generate and return the natural logarithm of the Inverse Cumulative Distribution (Quantile) Function ...
Generate and return a scalar (or array of arbitrary rank) of the natural logarithm(s) of random value...
Return the natural logarithm of the Cumulative Distribution Function (CDF) of the (Truncated) Poweto ...
Return the natural logarithm of the Probability Density Function (PDF) of the (Truncated) Poweto dist...
Return the natural logarithm of the Inverse Cumulative Distribution (Quantile) Function of the (Trunc...
Return a scalar (or array of arbitrary rank) of the natural logarithm(s) of random value(s) from the ...
This module contains classes and procedures for computing various statistical quantities related to t...
character(*, SK), parameter MODULE_NAME
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 Poweto as defined in the descr...