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

Generate and return the integral of the input function getFunc() in the closed range [lb, ub] using the Adaptive Romberg extrapolation method. More...

Detailed Description

Generate and return the integral of the input function getFunc() in the closed range [lb, ub] using the Adaptive Romberg extrapolation method.

This Romberg integration method is quite powerful for sufficiently smooth (e.g., analytic) integrands getFunc(), integrated over bounded intervals which contain no singularities, and where the end points are also nonsingular.

Parameters
getFunc: The input function to be integrated (i.e., the integrand).
  1. On entry, it must take an input scalar of the same type and kind as quadRomb.
  2. On exit, it must generate an input scalar of the same type and kind as quadRomb, representing the corresponding function value.
    The following illustrates the general interface of getFunc:
    function getFunc(x) result(func)
    use pm_kind, only: RK => RKG
    real(RK) , intent(in) :: x
    real(RK) :: func
    end function
    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
    where RKG refers to any desired real kind supported by the processor.
[in]lb: The input scalar of the same type and kind as the output quadRomb, containing the lower bound of the integration.
[in]ub: The input scalar of the same type and kind as the output quadRomb, containing the upper bound of the integration.
[in]tol: The input scalar of the same type and kind as the output quadRomb, containing the relative error the integration.
The algorithm converges if \(\ms{relerr} ~(~\equiv |\Delta\ms{quadRomb}|~) ~\leq~ \ms{tol} \times | \ms{quadRomb} |\).
Note that tol > epsilon(0._RKG) must hold at all times for integration to converge. Here RKG is the desired real kind of the output.
Ideally, set tol to a value such that tol < epsilon(0._RKG) * 100 holds to ensure convergence.
[in]nref: The input scalar integer of default kind IK, representing the number of refinements to be used in the Romberg method.
Think of nref as the maximum possible degree of the polynomial extrapolation used for approximating the integral at any stage.
  • The smaller values of nref can delay an accurate estimation of the integral via the Romberg method.
  • The larger values of nref can delay the first estimation of the integral by requiring more function evaluations.
  • The computational precision of the real kind used in this procedure imposes an upper limit on the value of nref.
    The maximum value for nref is roughly equal to int(log(epsilon(1._RKG)) / log(0.25)) with RKG representing the real kind used for the integration.
  • The maximum nref for real32, real64, real128 are respectively 12, 26, 56.
  • If the specified nref is larger than the maximum possible value, the integration will fail to converge.
  • The number nref = 2 corresponds to the famous Simpson integration rule.
  • A number between 4-6 is frequently a reasonable choice.
[out]relerr: The output scalar of the same type and kind as the output quadRomb containing the final estimated relative error in the result.
By definition, this is always a positive value if the integration converges.
Specify the relerr optional output argument to monitor convergence.
If relerr < 0., then the integration has failed to converge.
(optional. If missing and the integration fails to converge, the program will halt by calling error stop.)
[out]neval: The output scalar integer of default kind IK, representing the number of function evaluations made within the integrator.
(optional. It can be present if and only if relerr argument is also present.)
Returns
quadRomb : The output scalar real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128), containing the integration result.


Possible calling interfaces

quadRomb = getQuadRomb(getFunc, lb, ub, tol, nref)
quadRomb = getQuadRomb(getFunc, lb, ub, tol, nref, relerr)
quadRomb = getQuadRomb(getFunc, lb, ub, tol, nref, relerr, neval)
quadRomb = getQuadRomb(getFunc, lb, ub, tol, nref, interval)
quadRomb = getQuadRomb(getFunc, lb, ub, tol, nref, interval, relerr)
quadRomb = getQuadRomb(getFunc, lb, ub, tol, nref, interval, relerr, neval)
Generate and return the integral of the input function getFunc() in the closed range [lb,...
This module contains classes and procedures to perform numerical integrations.
Definition: pm_quadRomb.F90:31
Warning
The procedures of this generic interface will behave differently if the integration fails to converge:
  • If the optional relerr output argument is missing, the program will halt by calling error stop upon integration convergence failure.
  • If the optional relerr output argument is present, the program will return relerr < 0. to indicate integration convergence failure.
See also
getQuadErr


Example usage

1program example
2
3 use pm_kind, only: SK, IK, QP => RKH
4 use pm_mathBeta, only: getBetaInc
5 use pm_distBeta, only: getBetaPDF
7 use pm_io, only: display_type
8
9 implicit none
10
11 integer, parameter :: SP = kind(0.e0), DP = kind(0.d0)
12
13 integer(IK) :: neval
14
15 real(SP) :: quad_sp, quadref_sp, relerr_sp, alpha_sp, beta_sp
16 real(DP) :: quad_dp, quadref_dp, relerr_dp, alpha_dp, beta_dp
17 real(QP) :: quad_qp, quadref_qp, relerr_qp, alpha_qp, beta_qp
18
19 type(display_type) :: disp
20 disp = display_type(file = "main.out.F90")
21
22 call disp%skip()
23 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%show("! Compute the Cumulative Distribution Function (CDF) over a closed interval of the Beta distribution by numerical integration.")
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%skip()
29
30 call disp%skip()
31 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
32 call disp%show("! Compute the numerical integration with single precision.")
33 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
34 call disp%skip()
35
36 call disp%skip()
37 call disp%show("alpha_sp = 2.")
38 alpha_sp = 2.
39 call disp%show("beta_sp = 5.")
40 beta_sp = 5.
41 call disp%show("quadref_sp = getBetaInc(1., alpha = alpha_sp, beta = beta_sp) - getBetaInc(0., alpha = alpha_sp, beta = beta_sp)")
42 quadref_sp = getBetaInc(1., alpha = alpha_sp, beta = beta_sp) - getBetaInc(0., alpha = alpha_sp, beta = beta_sp)
43 call disp%show("quadref_sp")
44 call disp%show( quadref_sp )
45 call disp%show("quad_sp = getQuadRomb(getFunc = getBetaPDF_SP, lb = 0., ub = 1., tol = epsilon(1.) * 100, nref = 4_IK, relerr = relerr_sp, neval = neval)")
46 quad_sp = getQuadRomb(getFunc = getBetaPDF_SP, lb = 0., ub = 1., tol = epsilon(1.) * 100, nref = 4_IK, relerr = relerr_sp, neval = neval)
47 call disp%show("if (relerr_sp < 0.) error stop 'Integration failed to converge.'")
48 if (relerr_sp < 0.) error stop 'Integration failed to converge.'
49 call disp%show("relerr_sp ! < 0. if integration fails.")
50 call disp%show( relerr_sp )
51 call disp%show("neval ! # calls to the integrand.")
52 call disp%show( neval )
53 call disp%show("quad_sp ! integral")
54 call disp%show( quad_sp )
55 call disp%skip()
56
57 call disp%skip()
58 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
59 call disp%show("! Compute the numerical integration with double precision.")
60 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
61 call disp%skip()
62
63 call disp%skip()
64 call disp%show("alpha_dp = 2.d0")
65 alpha_dp = 2.d0
66 call disp%show("beta_dp = 5.d0")
67 beta_dp = 5.d0
68 call disp%show("quadref_dp = getBetaInc(1.d0, alpha = alpha_dp, beta = beta_dp) - getBetaInc(0.d0, alpha = alpha_dp, beta = beta_dp)")
69 quadref_dp = getBetaInc(1.d0, alpha = alpha_dp, beta = beta_dp) - getBetaInc(0.d0, alpha = alpha_dp, beta = beta_dp)
70 call disp%show("quadref_dp")
71 call disp%show( quadref_dp )
72 call disp%show("quad_dp = getQuadRomb(getFunc = getBetaPDF_DP, lb = 0.d0, ub = 1.d0, tol = epsilon(1.d0) * 100, nref = 4_IK, relerr = relerr_dp, neval = neval)")
73 quad_dp = getQuadRomb(getFunc = getBetaPDF_DP, lb = 0.d0, ub = 1.d0, tol = epsilon(1.d0) * 100, nref = 4_IK, relerr = relerr_dp, neval = neval)
74 call disp%show("if (relerr_dp < 0.) error stop 'Integration failed to converge.'")
75 if (relerr_dp < 0.) error stop 'Integration failed to converge.'
76 call disp%show("relerr_dp ! < 0. if integration fails.")
77 call disp%show( relerr_dp )
78 call disp%show("neval ! # calls to the integrand.")
79 call disp%show( neval )
80 call disp%show("quad_dp ! integral")
81 call disp%show( quad_dp )
82 call disp%skip()
83
84 call disp%skip()
85 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
86 call disp%show("! Compute the numerical integration with quadro precision.")
87 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
88 call disp%skip()
89
90 call disp%skip()
91 call disp%show("alpha_qp = 2._qp")
92 alpha_qp = 2._qp
93 call disp%show("beta_qp = 5._qp")
94 beta_qp = 5._qp
95 call disp%show("quadref_qp = getBetaInc(1._QP, alpha = alpha_qp, beta = beta_qp) - getBetaInc(0._QP, alpha = alpha_qp, beta = beta_qp)")
96 quadref_qp = getBetaInc(1._QP, alpha = alpha_qp, beta = beta_qp) - getBetaInc(0._QP, alpha = alpha_qp, beta = beta_qp)
97 call disp%show("quadref_qp")
98 call disp%show( quadref_qp )
99 call disp%show("quad_qp = getQuadRomb(getFunc = getBetaPDF_QP, lb = 0._QP, ub = 1._QP, tol = epsilon(1._QP) * 100, nref = 4_IK, relerr = relerr_qp, neval = neval)")
100 quad_qp = getQuadRomb(getFunc = getBetaPDF_QP, lb = 0._QP, ub = 1._QP, tol = epsilon(1._QP) * 100, nref = 4_IK, relerr = relerr_qp, neval = neval)
101 call disp%show("if (relerr_qp < 0.) error stop 'Integration failed to converge.'")
102 if (relerr_qp < 0.) error stop 'Integration failed to converge.'
103 call disp%show("relerr_qp ! < 0. if integration fails.")
104 call disp%show( relerr_qp )
105 call disp%show("neval ! # calls to the integrand.")
106 call disp%show( neval )
107 call disp%show("quad_qp ! integral")
108 call disp%show( quad_qp )
109 call disp%skip()
110
111contains
112
113 function getBetaPDF_SP(x) result(betaPDF)
114 real , intent(in) :: x
115 real :: betaPDF
116 betaPDF = getBetaPDF(x, alpha = alpha_sp, beta = beta_sp)
117 end function
118
119 function getBetaPDF_DP(x) result(betaPDF)
120 real(DP), intent(in) :: x
121 real(DP) :: betaPDF
122 betaPDF = getBetaPDF(x, alpha = alpha_dp, beta = beta_dp)
123 end function
124
125 function getBetaPDF_QP(x) result(betaPDF)
126 real(QP), intent(in) :: x
127 real(QP) :: betaPDF
128 betaPDF = getBetaPDF(x, alpha = alpha_qp, beta = beta_qp)
129 end function
130
131end program example
Generate and return the Probability Density Function (PDF) of the Beta distribution for an input x wi...
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
Generate and return the regularized Incomplete Beta Function as defined in the details section of pm...
This module contains classes and procedures for computing various statistical quantities related to t...
Definition: pm_distBeta.F90:99
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
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
This module contains classes and procedures for computing the mathematical Beta Function and its inve...
Definition: pm_mathBeta.F90:83
Generate and return an object of type display_type.
Definition: pm_io.F90:10282
This is the indicator type for generating instances of objects that indicate the integration interval...
This is the indicator type for generating instances of objects that indicate the integration interval...

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!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
4! Compute the Cumulative Distribution Function (CDF) over a closed interval of the Beta distribution by numerical integration.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
10! Compute the numerical integration with single precision.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13
14alpha_sp = 2.
15beta_sp = 5.
16quadref_sp = getBetaInc(1., alpha = alpha_sp, beta = beta_sp) - getBetaInc(0., alpha = alpha_sp, beta = beta_sp)
17quadref_sp
18+1.00000000
19quad_sp = getQuadRomb(getFunc = getBetaPDF_SP, lb = 0., ub = 1., tol = epsilon(1.) * 100, nref = 4_IK, relerr = relerr_sp, neval = neval)
20if (relerr_sp < 0.) error stop 'Integration failed to converge.'
21relerr_sp ! < 0. if integration fails.
22+0.946105527E-9
23neval ! # calls to the integrand.
24+9
25quad_sp ! integral
26+1.00000012
27
28
29!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
30! Compute the numerical integration with double precision.
31!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
32
33
34alpha_dp = 2.d0
35beta_dp = 5.d0
36quadref_dp = getBetaInc(1.d0, alpha = alpha_dp, beta = beta_dp) - getBetaInc(0.d0, alpha = alpha_dp, beta = beta_dp)
37quadref_dp
38+1.0000000000000000
39quad_dp = getQuadRomb(getFunc = getBetaPDF_DP, lb = 0.d0, ub = 1.d0, tol = epsilon(1.d0) * 100, nref = 4_IK, relerr = relerr_dp, neval = neval)
40if (relerr_dp < 0.) error stop 'Integration failed to converge.'
41relerr_dp ! < 0. if integration fails.
42+0.11564823173178713E-17
43neval ! # calls to the integrand.
44+9
45quad_dp ! integral
46+1.0000000000000002
47
48
49!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
50! Compute the numerical integration with quadro precision.
51!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
52
53
54alpha_qp = 2._qp
55beta_qp = 5._qp
56quadref_qp = getBetaInc(1._QP, alpha = alpha_qp, beta = beta_qp) - getBetaInc(0._QP, alpha = alpha_qp, beta = beta_qp)
57quadref_qp
58+1.00000000000000000000000000000000000
59quad_qp = getQuadRomb(getFunc = getBetaPDF_QP, lb = 0._QP, ub = 1._QP, tol = epsilon(1._QP) * 100, nref = 4_IK, relerr = relerr_qp, neval = neval)
60if (relerr_qp < 0.) error stop 'Integration failed to converge.'
61relerr_qp ! < 0. if integration fails.
62+0.152851582887875861353649043062295812E-35
63neval ! # calls to the integrand.
64+9
65quad_qp ! integral
66+1.00000000000000000000000000000000019
67
68
Test:
test_pm_quadRomb


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, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin
Joshua Alexander Osborne, May 28, 2020, 8:58 PM, Arlington, TX

Definition at line 396 of file pm_quadRomb.F90.


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