ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_distBinom.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
123
124!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
125
127
128 use pm_kind, only: SK, IK, LK, RKB
130
131 implicit none
132
133 character(*, SK), parameter :: MODULE_NAME = "@pm_distBinom"
134
135!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
136
171 end type
172
173!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
174
238
239 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
240
241#if RK5_ENABLED
242 PURE elemental module function getBinomLogPMF_RK5(nsuc, ntry, psuc) result(logPMF)
243#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
244 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomLogPMF_RK5
245#endif
246 use pm_kind, only: RKG => RK5
247 integer(IK) , intent(in) :: nsuc, ntry
248 real(RKG) , intent(in) :: psuc
249 real(RKG) :: logPMF
250 end function
251#endif
252
253#if RK4_ENABLED
254 PURE elemental module function getBinomLogPMF_RK4(nsuc, ntry, psuc) result(logPMF)
255#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
256 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomLogPMF_RK4
257#endif
258 use pm_kind, only: RKG => RK4
259 integer(IK) , intent(in) :: nsuc, ntry
260 real(RKG) , intent(in) :: psuc
261 real(RKG) :: logPMF
262 end function
263#endif
264
265#if RK3_ENABLED
266 PURE elemental module function getBinomLogPMF_RK3(nsuc, ntry, psuc) result(logPMF)
267#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
268 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomLogPMF_RK3
269#endif
270 use pm_kind, only: RKG => RK3
271 integer(IK) , intent(in) :: nsuc, ntry
272 real(RKG) , intent(in) :: psuc
273 real(RKG) :: logPMF
274 end function
275#endif
276
277#if RK2_ENABLED
278 PURE elemental module function getBinomLogPMF_RK2(nsuc, ntry, psuc) result(logPMF)
279#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
280 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomLogPMF_RK2
281#endif
282 use pm_kind, only: RKG => RK2
283 integer(IK) , intent(in) :: nsuc, ntry
284 real(RKG) , intent(in) :: psuc
285 real(RKG) :: logPMF
286 end function
287#endif
288
289#if RK1_ENABLED
290 PURE elemental module function getBinomLogPMF_RK1(nsuc, ntry, psuc) result(logPMF)
291#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
292 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomLogPMF_RK1
293#endif
294 use pm_kind, only: RKG => RK1
295 integer(IK) , intent(in) :: nsuc, ntry
296 real(RKG) , intent(in) :: psuc
297 real(RKG) :: logPMF
298 end function
299#endif
300
301 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
302
303 end interface
304
305!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
306
376
377 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
378
379#if RK5_ENABLED
380 PURE elemental module subroutine setBinomLogPMF_RK5(logPMF, nsuc, ntry, logp, logq)
381#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
382 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomLogPMF_RK5
383#endif
384 use pm_kind, only: RKG => RK5
385 real(RKG) , intent(out) :: logPMF
386 real(RKG) , intent(in) :: logp, logq
387 integer(IK) , intent(in) :: nsuc, ntry
388 end subroutine
389#endif
390
391#if RK4_ENABLED
392 PURE elemental module subroutine setBinomLogPMF_RK4(logPMF, nsuc, ntry, logp, logq)
393#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
394 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomLogPMF_RK4
395#endif
396 use pm_kind, only: RKG => RK4
397 real(RKG) , intent(out) :: logPMF
398 real(RKG) , intent(in) :: logp, logq
399 integer(IK) , intent(in) :: nsuc, ntry
400 end subroutine
401#endif
402
403#if RK3_ENABLED
404 PURE elemental module subroutine setBinomLogPMF_RK3(logPMF, nsuc, ntry, logp, logq)
405#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
406 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomLogPMF_RK3
407#endif
408 use pm_kind, only: RKG => RK3
409 real(RKG) , intent(out) :: logPMF
410 real(RKG) , intent(in) :: logp, logq
411 integer(IK) , intent(in) :: nsuc, ntry
412 end subroutine
413#endif
414
415#if RK2_ENABLED
416 PURE elemental module subroutine setBinomLogPMF_RK2(logPMF, nsuc, ntry, logp, logq)
417#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
418 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomLogPMF_RK2
419#endif
420 use pm_kind, only: RKG => RK2
421 real(RKG) , intent(out) :: logPMF
422 real(RKG) , intent(in) :: logp, logq
423 integer(IK) , intent(in) :: nsuc, ntry
424 end subroutine
425#endif
426
427#if RK1_ENABLED
428 PURE elemental module subroutine setBinomLogPMF_RK1(logPMF, nsuc, ntry, logp, logq)
429#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
430 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomLogPMF_RK1
431#endif
432 use pm_kind, only: RKG => RK1
433 real(RKG) , intent(out) :: logPMF
434 real(RKG) , intent(in) :: logp, logq
435 integer(IK) , intent(in) :: nsuc, ntry
436 end subroutine
437#endif
438
439 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
440
441 end interface
442
443!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
444
508 interface getBinomCDF
509
510 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
511
512#if RK5_ENABLED
513 PURE elemental module function getBinomCDF_RK5(nsuc, ntry, psuc) result(cdf)
514#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
515 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomCDF_RK5
516#endif
517 use pm_kind, only: RKG => RK5
518 integer(IK) , intent(in) :: nsuc, ntry
519 real(RKG) , intent(in) :: psuc
520 real(RKG) :: cdf
521 end function
522#endif
523
524#if RK4_ENABLED
525 PURE elemental module function getBinomCDF_RK4(nsuc, ntry, psuc) result(cdf)
526#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
527 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomCDF_RK4
528#endif
529 use pm_kind, only: RKG => RK4
530 integer(IK) , intent(in) :: nsuc, ntry
531 real(RKG) , intent(in) :: psuc
532 real(RKG) :: cdf
533 end function
534#endif
535
536#if RK3_ENABLED
537 PURE elemental module function getBinomCDF_RK3(nsuc, ntry, psuc) result(cdf)
538#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
539 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomCDF_RK3
540#endif
541 use pm_kind, only: RKG => RK3
542 integer(IK) , intent(in) :: nsuc, ntry
543 real(RKG) , intent(in) :: psuc
544 real(RKG) :: cdf
545 end function
546#endif
547
548#if RK2_ENABLED
549 PURE elemental module function getBinomCDF_RK2(nsuc, ntry, psuc) result(cdf)
550#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
551 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomCDF_RK2
552#endif
553 use pm_kind, only: RKG => RK2
554 integer(IK) , intent(in) :: nsuc, ntry
555 real(RKG) , intent(in) :: psuc
556 real(RKG) :: cdf
557 end function
558#endif
559
560#if RK1_ENABLED
561 PURE elemental module function getBinomCDF_RK1(nsuc, ntry, psuc) result(cdf)
562#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
563 !DEC$ ATTRIBUTES DLLEXPORT :: getBinomCDF_RK1
564#endif
565 use pm_kind, only: RKG => RK1
566 integer(IK) , intent(in) :: nsuc, ntry
567 real(RKG) , intent(in) :: psuc
568 real(RKG) :: cdf
569 end function
570#endif
571
572 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
573
574 end interface
575
576!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
577
643 interface setBinomCDF
644
645 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
646
647#if RK5_ENABLED
648 PURE elemental module subroutine setBinomCDF_RK5(cdf, nsuc, ntry, psuc, info)
649#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
650 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomCDF_RK5
651#endif
652 use pm_kind, only: RKG => RK5
653 integer(IK) , intent(in) :: nsuc, ntry
654 real(RKG) , intent(in) :: psuc
655 real(RKG) , intent(out) :: cdf
656 integer(IK) , intent(out) :: info
657 end subroutine
658#endif
659
660#if RK4_ENABLED
661 PURE elemental module subroutine setBinomCDF_RK4(cdf, nsuc, ntry, psuc, info)
662#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
663 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomCDF_RK4
664#endif
665 use pm_kind, only: RKG => RK4
666 integer(IK) , intent(in) :: nsuc, ntry
667 real(RKG) , intent(in) :: psuc
668 real(RKG) , intent(out) :: cdf
669 integer(IK) , intent(out) :: info
670 end subroutine
671#endif
672
673#if RK3_ENABLED
674 PURE elemental module subroutine setBinomCDF_RK3(cdf, nsuc, ntry, psuc, info)
675#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
676 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomCDF_RK3
677#endif
678 use pm_kind, only: RKG => RK3
679 integer(IK) , intent(in) :: nsuc, ntry
680 real(RKG) , intent(in) :: psuc
681 real(RKG) , intent(out) :: cdf
682 integer(IK) , intent(out) :: info
683 end subroutine
684#endif
685
686#if RK2_ENABLED
687 PURE elemental module subroutine setBinomCDF_RK2(cdf, nsuc, ntry, psuc, info)
688#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
689 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomCDF_RK2
690#endif
691 use pm_kind, only: RKG => RK2
692 integer(IK) , intent(in) :: nsuc, ntry
693 real(RKG) , intent(in) :: psuc
694 real(RKG) , intent(out) :: cdf
695 integer(IK) , intent(out) :: info
696 end subroutine
697#endif
698
699#if RK1_ENABLED
700 PURE elemental module subroutine setBinomCDF_RK1(cdf, nsuc, ntry, psuc, info)
701#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
702 !DEC$ ATTRIBUTES DLLEXPORT :: setBinomCDF_RK1
703#endif
704 use pm_kind, only: RKG => RK1
705 integer(IK) , intent(in) :: nsuc, ntry
706 real(RKG) , intent(in) :: psuc
707 real(RKG) , intent(out) :: cdf
708 integer(IK) , intent(out) :: info
709 end subroutine
710#endif
711
712 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
713
714 end interface
715
716!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
717
718end module pm_distBinom
Generate and return the Cumulative Distribution Function (CDF) of the Binomial distribution for an in...
Generate and return the natural logarithm of the Probability Mass Function (PMF) of the Binomial dist...
Return the Cumulative Distribution Function (CDF) of the Binomial distribution.
Return the natural logarithm of the Probability Mass Function (PMF) of the Binomial distribution for ...
This module contains classes and procedures for computing various statistical quantities related to t...
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 RKB
The scalar integer constant of intrinsic default kind, representing the Best-precision real kind supp...
Definition: pm_kind.F90:1371
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 Binomial as defined in the des...
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...