ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_mathBeta.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
80
81!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
82
84
85 use pm_kind, only: SK, IK, LK
86
87 implicit none
88
89 character(*, SK), parameter :: MODULE_NAME = "@pm_mathBeta"
90
91!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
92
143 interface getLogBeta
144
145 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
146
147#if RK5_ENABLED
148 PURE elemental module function getLogBeta_RK5(alpha, beta) result(logFuncBeta)
149#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
150 !DEC$ ATTRIBUTES DLLEXPORT :: getLogBeta_RK5
151#endif
152 use pm_kind, only: RKG => RK5
153 real(RKG) , intent(in) :: alpha, beta
154 real(RKG) :: logFuncBeta
155 end function
156#endif
157
158#if RK4_ENABLED
159 PURE elemental module function getLogBeta_RK4(alpha, beta) result(logFuncBeta)
160#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
161 !DEC$ ATTRIBUTES DLLEXPORT :: getLogBeta_RK4
162#endif
163 use pm_kind, only: RKG => RK4
164 real(RKG) , intent(in) :: alpha, beta
165 real(RKG) :: logFuncBeta
166 end function
167#endif
168
169#if RK3_ENABLED
170 PURE elemental module function getLogBeta_RK3(alpha, beta) result(logFuncBeta)
171#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
172 !DEC$ ATTRIBUTES DLLEXPORT :: getLogBeta_RK3
173#endif
174 use pm_kind, only: RKG => RK3
175 real(RKG) , intent(in) :: alpha, beta
176 real(RKG) :: logFuncBeta
177 end function
178#endif
179
180#if RK2_ENABLED
181 PURE elemental module function getLogBeta_RK2(alpha, beta) result(logFuncBeta)
182#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
183 !DEC$ ATTRIBUTES DLLEXPORT :: getLogBeta_RK2
184#endif
185 use pm_kind, only: RKG => RK2
186 real(RKG) , intent(in) :: alpha, beta
187 real(RKG) :: logFuncBeta
188 end function
189#endif
190
191#if RK1_ENABLED
192 PURE elemental module function getLogBeta_RK1(alpha, beta) result(logFuncBeta)
193#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
194 !DEC$ ATTRIBUTES DLLEXPORT :: getLogBeta_RK1
195#endif
196 use pm_kind, only: RKG => RK1
197 real(RKG) , intent(in) :: alpha, beta
198 real(RKG) :: logFuncBeta
199 end function
200#endif
201
202 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
203
204 end interface
205
206!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
207
279 interface getBetaInc
280
281 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
282
283#if RK5_ENABLED
284 impure elemental module function getBetaInc_RK5(x, alpha, beta, signed) result(betaInc)
285#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
286 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInc_RK5
287#endif
288 use pm_kind, only: RKG => RK5
289 real(RKG) :: betaInc
290 real(RKG) , intent(in) :: x, alpha, beta
291 logical(LK) , intent(in) , optional :: signed
292 end function
293#endif
294
295#if RK4_ENABLED
296 impure elemental module function getBetaInc_RK4(x, alpha, beta, signed) result(betaInc)
297#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
298 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInc_RK4
299#endif
300 use pm_kind, only: RKG => RK4
301 real(RKG) :: betaInc
302 real(RKG) , intent(in) :: x, alpha, beta
303 logical(LK) , intent(in) , optional :: signed
304 end function
305#endif
306
307#if RK3_ENABLED
308 impure elemental module function getBetaInc_RK3(x, alpha, beta, signed) result(betaInc)
309#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
310 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInc_RK3
311#endif
312 use pm_kind, only: RKG => RK3
313 real(RKG) :: betaInc
314 real(RKG) , intent(in) :: x, alpha, beta
315 logical(LK) , intent(in) , optional :: signed
316 end function
317#endif
318
319#if RK2_ENABLED
320 impure elemental module function getBetaInc_RK2(x, alpha, beta, signed) result(betaInc)
321#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
322 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInc_RK2
323#endif
324 use pm_kind, only: RKG => RK2
325 real(RKG) :: betaInc
326 real(RKG) , intent(in) :: x, alpha, beta
327 logical(LK) , intent(in) , optional :: signed
328 end function
329#endif
330
331#if RK1_ENABLED
332 impure elemental module function getBetaInc_RK1(x, alpha, beta, signed) result(betaInc)
333#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
334 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInc_RK1
335#endif
336 use pm_kind, only: RKG => RK1
337 real(RKG) :: betaInc
338 real(RKG) , intent(in) :: x, alpha, beta
339 logical(LK) , intent(in) , optional :: signed
340 end function
341#endif
342
343 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
344
345 end interface
346
347!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
348
442 interface setBetaInc
443
444 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
445 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
446 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
447
448 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
449
450#if RK5_ENABLED
451 PURE elemental module subroutine setBetaIncDef_RK5(betaInc, x, alpha, beta, logFuncBeta, signed, info)
452#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
453 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncDef_RK5
454#endif
455 use pm_kind, only: RKG => RK5
456 real(RKG) , intent(out) :: betaInc
457 real(RKG) , intent(in) :: logFuncBeta
458 real(RKG) , value :: x, alpha, beta
459 logical(LK) , intent(in) :: signed
460 integer(IK) , intent(out) :: info
461 end subroutine
462#endif
463
464#if RK4_ENABLED
465 PURE elemental module subroutine setBetaIncDef_RK4(betaInc, x, alpha, beta, logFuncBeta, signed, info)
466#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
467 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncDef_RK4
468#endif
469 use pm_kind, only: RKG => RK4
470 real(RKG) , intent(out) :: betaInc
471 real(RKG) , intent(in) :: logFuncBeta
472 real(RKG) , value :: x, alpha, beta
473 logical(LK) , intent(in) :: signed
474 integer(IK) , intent(out) :: info
475 end subroutine
476#endif
477
478#if RK3_ENABLED
479 PURE elemental module subroutine setBetaIncDef_RK3(betaInc, x, alpha, beta, logFuncBeta, signed, info)
480#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
481 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncDef_RK3
482#endif
483 use pm_kind, only: RKG => RK3
484 real(RKG) , intent(out) :: betaInc
485 real(RKG) , intent(in) :: logFuncBeta
486 real(RKG) , value :: x, alpha, beta
487 logical(LK) , intent(in) :: signed
488 integer(IK) , intent(out) :: info
489 end subroutine
490#endif
491
492#if RK2_ENABLED
493 PURE elemental module subroutine setBetaIncDef_RK2(betaInc, x, alpha, beta, logFuncBeta, signed, info)
494#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
495 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncDef_RK2
496#endif
497 use pm_kind, only: RKG => RK2
498 real(RKG) , intent(out) :: betaInc
499 real(RKG) , intent(in) :: logFuncBeta
500 real(RKG) , value :: x, alpha, beta
501 logical(LK) , intent(in) :: signed
502 integer(IK) , intent(out) :: info
503 end subroutine
504#endif
505
506#if RK1_ENABLED
507 PURE elemental module subroutine setBetaIncDef_RK1(betaInc, x, alpha, beta, logFuncBeta, signed, info)
508#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
509 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncDef_RK1
510#endif
511 use pm_kind, only: RKG => RK1
512 real(RKG) , intent(out) :: betaInc
513 real(RKG) , intent(in) :: logFuncBeta
514 real(RKG) , value :: x, alpha, beta
515 logical(LK) , intent(in) :: signed
516 integer(IK) , intent(out) :: info
517 end subroutine
518#endif
519
520 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
521
522#if RK5_ENABLED
523 impure elemental module subroutine setBetaIncGK21_RK5(betaInc, x, alpha, beta, logFuncBeta, reltol, signed, info)
524#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
525 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncGK21_RK5
526#endif
527 use pm_kind, only: RKG => RK5
528 real(RKG) , intent(out) :: betaInc
529 real(RKG) , intent(in) :: logFuncBeta
530 real(RKG) , value :: x, alpha, beta
531 real(RKG) , intent(in) :: reltol
532 logical(LK) , intent(in) :: signed
533 integer(IK) , intent(out) :: info
534 end subroutine
535#endif
536
537#if RK4_ENABLED
538 impure elemental module subroutine setBetaIncGK21_RK4(betaInc, x, alpha, beta, logFuncBeta, reltol, signed, info)
539#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
540 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncGK21_RK4
541#endif
542 use pm_kind, only: RKG => RK4
543 real(RKG) , intent(out) :: betaInc
544 real(RKG) , intent(in) :: logFuncBeta
545 real(RKG) , value :: x, alpha, beta
546 real(RKG) , intent(in) :: reltol
547 logical(LK) , intent(in) :: signed
548 integer(IK) , intent(out) :: info
549 end subroutine
550#endif
551
552#if RK3_ENABLED
553 impure elemental module subroutine setBetaIncGK21_RK3(betaInc, x, alpha, beta, logFuncBeta, reltol, signed, info)
554#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
555 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncGK21_RK3
556#endif
557 use pm_kind, only: RKG => RK3
558 real(RKG) , intent(out) :: betaInc
559 real(RKG) , intent(in) :: logFuncBeta
560 real(RKG) , value :: x, alpha, beta
561 real(RKG) , intent(in) :: reltol
562 logical(LK) , intent(in) :: signed
563 integer(IK) , intent(out) :: info
564 end subroutine
565#endif
566
567#if RK2_ENABLED
568 impure elemental module subroutine setBetaIncGK21_RK2(betaInc, x, alpha, beta, logFuncBeta, reltol, signed, info)
569#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
570 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncGK21_RK2
571#endif
572 use pm_kind, only: RKG => RK2
573 real(RKG) , intent(out) :: betaInc
574 real(RKG) , intent(in) :: logFuncBeta
575 real(RKG) , value :: x, alpha, beta
576 real(RKG) , intent(in) :: reltol
577 logical(LK) , intent(in) :: signed
578 integer(IK) , intent(out) :: info
579 end subroutine
580#endif
581
582#if RK1_ENABLED
583 impure elemental module subroutine setBetaIncGK21_RK1(betaInc, x, alpha, beta, logFuncBeta, reltol, signed, info)
584#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
585 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaIncGK21_RK1
586#endif
587 use pm_kind, only: RKG => RK1
588 real(RKG) , intent(out) :: betaInc
589 real(RKG) , intent(in) :: logFuncBeta
590 real(RKG) , value :: x, alpha, beta
591 real(RKG) , intent(in) :: reltol
592 logical(LK) , intent(in) :: signed
593 integer(IK) , intent(out) :: info
594 end subroutine
595#endif
596
597 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
598
599 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
600 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
601 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
602
603 end interface
604
605!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
606
677 interface getBetaInv
678
679 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
680
681#if RK5_ENABLED
682 PURE elemental module function getBetaInv_RK5(betaInc, alpha, beta, signed) result(betaInv)
683#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
684 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInv_RK5
685#endif
686 use pm_kind, only: RKG => RK5
687 real(RKG) :: betaInv
688 real(RKG) , intent(in) :: betaInc, alpha, beta
689 logical(LK) , intent(in) , optional :: signed
690 end function
691#endif
692
693#if RK4_ENABLED
694 PURE elemental module function getBetaInv_RK4(betaInc, alpha, beta, signed) result(betaInv)
695#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
696 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInv_RK4
697#endif
698 use pm_kind, only: RKG => RK4
699 real(RKG) :: betaInv
700 real(RKG) , intent(in) :: betaInc, alpha, beta
701 logical(LK) , intent(in) , optional :: signed
702 end function
703#endif
704
705#if RK3_ENABLED
706 PURE elemental module function getBetaInv_RK3(betaInc, alpha, beta, signed) result(betaInv)
707#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
708 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInv_RK3
709#endif
710 use pm_kind, only: RKG => RK3
711 real(RKG) :: betaInv
712 real(RKG) , intent(in) :: betaInc, alpha, beta
713 logical(LK) , intent(in) , optional :: signed
714 end function
715#endif
716
717#if RK2_ENABLED
718 PURE elemental module function getBetaInv_RK2(betaInc, alpha, beta, signed) result(betaInv)
719#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
720 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInv_RK2
721#endif
722 use pm_kind, only: RKG => RK2
723 real(RKG) :: betaInv
724 real(RKG) , intent(in) :: betaInc, alpha, beta
725 logical(LK) , intent(in) , optional :: signed
726 end function
727#endif
728
729#if RK1_ENABLED
730 PURE elemental module function getBetaInv_RK1(betaInc, alpha, beta, signed) result(betaInv)
731#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
732 !DEC$ ATTRIBUTES DLLEXPORT :: getBetaInv_RK1
733#endif
734 use pm_kind, only: RKG => RK1
735 real(RKG) :: betaInv
736 real(RKG) , intent(in) :: betaInc, alpha, beta
737 logical(LK) , intent(in) , optional :: signed
738 end function
739#endif
740
741 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
742
743 end interface
744
745!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
746
826 interface setBetaInv
827
828 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
829
830#if RK5_ENABLED
831 PURE elemental module subroutine setBetaInv_RK5(betaInv, betaInc, alpha, beta, logFuncBeta, signed, info)
832#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
833 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaInv_RK5
834#endif
835 use pm_kind, only: RKG => RK5
836 real(RKG) , intent(out) :: betaInv
837 real(RKG) , value :: betaInc, alpha, beta, logFuncBeta
838 logical(LK) , intent(in) :: signed
839 integer(IK) , intent(out) :: info
840 end subroutine
841#endif
842
843#if RK4_ENABLED
844 PURE elemental module subroutine setBetaInv_RK4(betaInv, betaInc, alpha, beta, logFuncBeta, signed, info)
845#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
846 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaInv_RK4
847#endif
848 use pm_kind, only: RKG => RK4
849 real(RKG) , intent(out) :: betaInv
850 real(RKG) , value :: betaInc, alpha, beta, logFuncBeta
851 logical(LK) , intent(in) :: signed
852 integer(IK) , intent(out) :: info
853 end subroutine
854#endif
855
856#if RK3_ENABLED
857 PURE elemental module subroutine setBetaInv_RK3(betaInv, betaInc, alpha, beta, logFuncBeta, signed, info)
858#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
859 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaInv_RK3
860#endif
861 use pm_kind, only: RKG => RK3
862 real(RKG) , intent(out) :: betaInv
863 real(RKG) , value :: betaInc, alpha, beta, logFuncBeta
864 logical(LK) , intent(in) :: signed
865 integer(IK) , intent(out) :: info
866 end subroutine
867#endif
868
869#if RK2_ENABLED
870 PURE elemental module subroutine setBetaInv_RK2(betaInv, betaInc, alpha, beta, logFuncBeta, signed, info)
871#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
872 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaInv_RK2
873#endif
874 use pm_kind, only: RKG => RK2
875 real(RKG) , intent(out) :: betaInv
876 real(RKG) , value :: betaInc, alpha, beta, logFuncBeta
877 logical(LK) , intent(in) :: signed
878 integer(IK) , intent(out) :: info
879 end subroutine
880#endif
881
882#if RK1_ENABLED
883 PURE elemental module subroutine setBetaInv_RK1(betaInv, betaInc, alpha, beta, logFuncBeta, signed, info)
884#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
885 !DEC$ ATTRIBUTES DLLEXPORT :: setBetaInv_RK1
886#endif
887 use pm_kind, only: RKG => RK1
888 real(RKG) , intent(out) :: betaInv
889 real(RKG) , value :: betaInc, alpha, beta, logFuncBeta
890 logical(LK) , intent(in) :: signed
891 integer(IK) , intent(out) :: info
892 end subroutine
893#endif
894
895 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
896
897 end interface
898
899!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
900
901end module pm_mathBeta
Generate and return the regularized Incomplete Beta Function as defined in the details section of pm...
Generate and return the regularized Inverse Incomplete Beta Function as defined in the details secti...
Generate and return the natural logarithm of the Beta Function as defined in the details section of ...
Return the regularized Incomplete Beta Function as defined in the details section of pm_mathBeta.
Return the regularized Inverse Incomplete Beta Function as defined in the details section of pm_math...
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 classes and procedures for computing the mathematical Beta Function and its inve...
Definition: pm_mathBeta.F90:83
character(*, SK), parameter MODULE_NAME
Definition: pm_mathBeta.F90:89