ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_mathBeta::getBetaInc Interface Reference

Generate and return the regularized Incomplete Beta Function \(I_x(\alpha, \beta)\) as defined in the details section of pm_mathBeta. More...

Detailed Description

Generate and return the regularized Incomplete Beta Function \(I_x(\alpha, \beta)\) as defined in the details section of pm_mathBeta.

Parameters
[in]x: The input scalar (or array of the same shape as other array-like arguments) of type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128).
[in]alpha: The input scalar (or array of the same shape as other array-like arguments) of the same type and kind as x.
[in]beta: The input scalar (or array of the same shape as other array-like arguments) of the same type and kind as x.
[in]signed: The input scalar (or array of the same shape as other array-like arguments) of type logical of default kind LK.
  1. If signed = .false., the input x must be in range 0 <= x <= 1 and the output betaInc will be the expected incomplete Beta function in range 0 <= betaInc <= 1.
  2. If signed = .true., then the following rules hold:
    1. If the condition x < 0 holds, then the value x = 1 - x < 0 will be used instead of x.
    2. If the output betaInc is near 1, the output will be returned as betaInc = betaInc - 1 < 0 instead of betaInc.
      Therefore, the user is expected to be aware of the convention and apply the necessary correction (addition by 1) before using the output value.
Specifying signed = .true. can lead to considerably more accurate calculations (by orders of magnitudes) for values of x and betaInc that are near 1.
The loss of precision near 1 occurs because of inadequate resolution of real number representations in digital computers near 1 which is orders of magnitude worse than the precision near 0.
(optional, default = .false., following the principle of least surprise.)
Returns
betaInc : The output object of the same type, kind, and rank as highest-rank input argument containing the regularized Incomplete Beta Function.


Possible calling interfaces

betaInc = getBetaInc(x, alpha, beta, signed = signed)
Generate and return the regularized Incomplete Beta Function as defined in the details section of pm...
This module contains classes and procedures for computing the mathematical Beta Function and its inve...
Definition: pm_mathBeta.F90:84
Warning
The conditions 0 <= x .and. .not. signed .or. 0 <= x must hold for the corresponding input arguments.
The conditions 0 < alpha .and. 0 < beta must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The procedures under this generic interface will abort the program if the computation of the continued fraction representation of the regularized beta function fails to converge.
If error control is necessary, use the generic interface seBetaInc.
Remarks
The procedures under discussion are impure.
The procedures under discussion are elemental.
See also
getLogBeta
getBetaInc
getBetaLogPDF


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_kind, only: RK ! all processor kinds are supported.
5 use pm_io, only: display_type
6 use pm_mathBeta, only: getBetaInc
7
8 implicit none
9
10 type(display_type) :: disp
11 disp = display_type(file = "main.out.F90")
12
13 call disp%skip()
14 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
15 call disp%show("! Compute the regularized (lower) Incomplete Beta Function using its series representation.")
16 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
17 call disp%skip()
18
19 call disp%skip()
20 call disp%show("getBetaInc(x = 0._RK, alpha = 2._RK, beta = 3._RK)")
21 call disp%show( getBetaInc(x = 0._RK, alpha = 2._RK, beta = 3._RK) )
22 call disp%skip()
23
24 call disp%skip()
25 call disp%show("getBetaInc(x = .5_RK, alpha = 2._RK, beta = 3._RK) ! 0.6875")
26 call disp%show( getBetaInc(x = .5_RK, alpha = 2._RK, beta = 3._RK) )
27 call disp%skip()
28
29 call disp%skip()
30 call disp%show("getBetaInc(x = 1._RK, alpha = 2._RK, beta = 3._RK)")
31 call disp%show( getBetaInc(x = 1._RK, alpha = 2._RK, beta = 3._RK) )
32 call disp%skip()
33
34 call disp%skip()
35 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
36 call disp%show("! Compute the regularized (lower) Incomplete Beta Function for a vector of points.")
37 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
38 call disp%skip()
39
40 call disp%skip()
41 call disp%show("getBetaInc(x = [0._RK, 0.5_RK, 1._RK], alpha = 2._RK, beta = 3._RK)")
42 call disp%show( getBetaInc(x = [0._RK, 0.5_RK, 1._RK], alpha = 2._RK, beta = 3._RK) )
43 call disp%skip()
44
45 call disp%skip()
46 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
47 call disp%show("! Compute the regularized Incomplete Beta Function for a vector of input parameters.")
48 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
49 call disp%skip()
50
51 call disp%skip()
52 call disp%show("getBetaInc(x = [0._RK, 0.5_RK, 1._RK], alpha = [0.1_RK, 1._RK, 10._RK], beta = 3._RK) ! 0, 0.875, 1")
53 call disp%show( getBetaInc(x = [0._RK, 0.5_RK, 1._RK], alpha = [0.1_RK, 1._RK, 10._RK], beta = 3._RK) )
54 call disp%skip()
55
56 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
57 ! Output an example array of the regularized Incomplete Beta function for visualization.
58 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59
60 !block
61 ! use pm_arraySpace, only: setLinSpace
62 ! integer(IK) , parameter :: nx = 1000_IK
63 ! real(RK) :: x(nx)
64 ! integer :: fileUnit, i
65 ! call setLinSpace(x, 0._RK, 1._RK)
66 ! open(newunit = fileUnit, file = "getBetaInc.RK.txt")
67 ! do i = 1, nx
68 ! write(fileUnit, "(*(g0,:,' '))") x_RK(i), getBetaInc(x(i), 5._RK, [1.0_RK, 5.0_RK, 10._RK])
69 ! end do
70 ! close(fileUnit)
71 !end block
72
73 block
74 use pm_kind, only: RKG => RKH, RKH
75 use pm_arraySpace, only: setLinSpace
76 integer(IK) , parameter :: nx = 1000
77 real(RKG) , parameter :: alpha(*) = [0.1_RKG, 10._RKG, 1._RKG, 0.1_RKG, 10._RKG], beta(*) = [0.1_RKG, 0.1_RKG, 1._RKG, 10._RKG, 10._RKG]
78 real(RKG) :: betaInc(max(size(alpha), size(beta)))
79 real(RKG) :: x(nx)
80 integer :: fileUnit, i
81 call setLinSpace(x, 0._RKG, 1._RKG, fopen = .true._LK, lopen = .true._LK)
82 open(newunit = fileUnit, file = "getBetaInc.RK.txt")
83 do i = 1, nx
84 betaInc = getBetaInc(x(i), alpha, beta)
85 write(fileUnit, "(*(g0,:,','))") x(i), merge(1 + betaInc, betaInc, betaInc < 0)
86 end do
87 close(fileUnit)
88 end block
89
90end program example
Return the linSpace output argument with size(linSpace) elements of evenly-spaced values over the int...
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
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 RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3! Compute the regularized (lower) Incomplete Beta Function using its series representation.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7getBetaInc(x = 0._RK, alpha = 2._RK, beta = 3._RK)
8+0.0000000000000000
9
10
11getBetaInc(x = .5_RK, alpha = 2._RK, beta = 3._RK) ! 0.6875
12+0.68750000000000000
13
14
15getBetaInc(x = 1._RK, alpha = 2._RK, beta = 3._RK)
16+1.0000000000000000
17
18
19!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
20! Compute the regularized (lower) Incomplete Beta Function for a vector of points.
21!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
22
23
24getBetaInc(x = [0._RK, 0.5_RK, 1._RK], alpha = 2._RK, beta = 3._RK)
25+0.0000000000000000, +0.68750000000000000, +1.0000000000000000
26
27
28!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
29! Compute the regularized Incomplete Beta Function for a vector of input parameters.
30!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
31
32
33getBetaInc(x = [0._RK, 0.5_RK, 1._RK], alpha = [0.1_RK, 1._RK, 10._RK], beta = 3._RK) ! 0, 0.875, 1
34+0.0000000000000000, +0.87500000000000000, +1.0000000000000000
35
36

Postprocessing of the example output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6import glob
7import sys
8
9fontsize = 17
10
11kind = "RK"
12#label = [ r"$\alpha, \beta = 5., 1.0$"
13# , r"$\alpha, \beta = 5., 5.0$"
14# , r"$\alpha, \beta = 5., 10.$"
15# ]
16label = [ r"$\alpha, \beta = .1, .1$"
17 , r"$\alpha, \beta = 10, .1$"
18 , r"$\alpha, \beta = 1., 1.$"
19 , r"$\alpha, \beta = .1, 10$"
20 , r"$\alpha, \beta = 10, 10$"
21 ]
22
23pattern = "*." + kind + ".txt"
24fileList = glob.glob(pattern)
25if len(fileList) == 1:
26
27 df = pd.read_csv(fileList[0], delimiter = ",")
28
29 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
30 ax = plt.subplot()
31
32 for i in range(1,len(df.values[0,:]+1)):
33
34 plt.plot( df.values[:, 0]
35 , df.values[:,i]
36 , linewidth = 2
37 )
38
39 plt.xticks(fontsize = fontsize - 2)
40 plt.yticks(fontsize = fontsize - 2)
41 ax.set_xlabel("x", fontsize = fontsize)
42 ax.set_ylabel("Regularized Lower\nIncomplete Beta Function", fontsize = fontsize)
43 #ax.set_xscale("log")
44 #ax.set_yscale("log")
45
46 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
47 ax.tick_params(axis = "y", which = "minor")
48 ax.tick_params(axis = "x", which = "minor")
49
50 ax.legend ( label
51 , fontsize = fontsize
52 #, loc = "center left"
53 #, bbox_to_anchor = (1, 0.5)
54 )
55
56 plt.savefig(fileList[0].replace(".txt",".png"))
57
58else:
59
60 sys.exit("Ambiguous file list exists.")

Visualization of the example output
Test:
test_pm_mathBeta


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, Oct 16, 2009, 11:14 AM, Michigan

Definition at line 280 of file pm_mathBeta.F90.


The documentation for this interface was generated from the following file: