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

Generate and return the natural logarithm of the Beta Function \(\mathrm{B}(\alpha, \beta)\) as defined in the details section of pm_mathBeta. More...

Detailed Description

Generate and return the natural logarithm of the Beta Function \(\mathrm{B}(\alpha, \beta)\) as defined in the details section of pm_mathBeta.

Parameters
[in]alpha: 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]beta: The input scalar (or array of the same shape as other array-like arguments) of the same type and kind as alpha.
Returns
logFuncBeta : The output object of the same type, kind, and rank as highest-rank input argument containing the natural logarithm of the Beta Function.


Possible calling interfaces

logFuncBeta = getLogBeta(alpha, beta)
Generate and return the natural logarithm of the Beta Function as defined in the details section of ...
This module contains classes and procedures for computing the mathematical Beta Function and its inve...
Definition: pm_mathBeta.F90:84
Warning
The conditions 0 < alpha .and. 0 < beta must hold for the corresponding input arguments.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
Remarks
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: getLogBeta
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("getLogBeta(alpha = 2._RK, beta = 3._RK) ! -2.484906649788000")
21 call disp%show( getLogBeta(alpha = 2._RK, beta = 3._RK) )
22 call disp%skip()
23
24 call disp%skip()
25 call disp%show("getLogBeta(alpha = 5._RK, beta = 3._RK) ! -4.653960350157524")
26 call disp%show( getLogBeta(alpha = 5._RK, beta = 3._RK) )
27 call disp%skip()
28
29 call disp%skip()
30 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%show("! Compute the regularized Incomplete Beta Function for a vector of input parameters.")
32 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 call disp%skip()
34
35 call disp%skip()
36 call disp%show("getLogBeta(alpha = [0.1_RK, 1._RK, 10._RK], beta = 3._RK) ! 2.158484749020289 -1.098612288668110 -6.492239835020470")
37 call disp%show( getLogBeta(alpha = [0.1_RK, 1._RK, 10._RK], beta = 3._RK) )
38 call disp%skip()
39
40 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
41 ! Output an example array of the regularized Incomplete Beta function for visualization.
42 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
43
44 block
45
46 use pm_arraySpace, only: getLinSpace
47 integer(IK) , parameter :: NP = 1000_IK
48 integer :: fileUnit, i
49 real(RK) :: A_RK(NP)
50
51 A_RK = exp(getLinSpace(-2._RK, 2._RK, NP))
52 open(newunit = fileUnit, file = "getLogBeta.RK.txt")
53 do i = 1, NP
54 write(fileUnit, "(*(g0,:,' '))") A_RK(i), exp(getLogBeta(A_RK(i), [0.1_RK, 1.0_RK, 10._RK]))
55 end do
56 close(fileUnit)
57
58 end block
59
60end program example
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
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
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
7getLogBeta(alpha = 2._RK, beta = 3._RK) ! -2.484906649788000
8-2.4849066497880004
9
10
11getLogBeta(alpha = 5._RK, beta = 3._RK) ! -4.653960350157524
12-4.6539603501575222
13
14
15!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
16! Compute the regularized Incomplete Beta Function for a vector of input parameters.
17!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
18
19
20getLogBeta(alpha = [0.1_RK, 1._RK, 10._RK], beta = 3._RK) ! 2.158484749020289 -1.098612288668110 -6.492239835020470
21+2.1584847490202885, -1.0986122886681096, -6.4922398350204720
22
23

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"
12label = [ r"$\beta = 0.1$"
13 , r"$\beta = 1.0$"
14 , r"$\beta = 10.$"
15 ]
16
17pattern = "*." + kind + ".txt"
18fileList = glob.glob(pattern)
19if len(fileList) == 1:
20
21 df = pd.read_csv(fileList[0], delimiter = " ")
22
23 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
24 ax = plt.subplot()
25
26 for i in range(1,len(df.values[0,:]+1)):
27
28 plt.plot( df.values[:, 0]
29 , df.values[:,i]
30 , linewidth = 2
31 )
32
33 plt.xticks(fontsize = fontsize - 2)
34 plt.yticks(fontsize = fontsize - 2)
35 ax.set_xlabel("a", fontsize = fontsize)
36 ax.set_ylabel("Beta Function", fontsize = fontsize)
37 #ax.set_xscale("log")
38 #ax.set_yscale("log")
39
40 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
41 ax.tick_params(axis = "y", which = "minor")
42 ax.tick_params(axis = "x", which = "minor")
43
44 ax.legend ( label
45 , fontsize = fontsize
46 #, loc = "center left"
47 #, bbox_to_anchor = (1, 0.5)
48 )
49
50 plt.savefig(fileList[0].replace(".txt",".png"))
51
52else:
53
54 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 144 of file pm_mathBeta.F90.


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