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