ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_quadPack::GK41_type Type Reference

This is the indicator type for generating instances of objects that indicate the use of 20-point Gauss-Legendre quadrature with 41-point Kronrod extension formulae as the drivers in the (Adaptive) Global Gauss-Kronrod Quadrature of getQuadGK. More...

Detailed Description

This is the indicator type for generating instances of objects that indicate the use of 20-point Gauss-Legendre quadrature with 41-point Kronrod extension formulae as the drivers in the (Adaptive) Global Gauss-Kronrod Quadrature of getQuadGK.

This is an empty derived type that exists solely for generating unique objects that are distinguishable as input arguments to procedures under the generic interface getQuadGK.


Possible calling interfaces

type(GK41_type) :: GK41
This module contains classes and procedures for non-adaptive and adaptive global numerical quadrature...
This is the indicator type for generating instances of objects that indicate the use of 20-point Gaus...
See also
GK15_type
GK21_type
GK31_type
GK41_type
GK51_type
GK61_type
getQuadGK
getQuadErr


Example usage

1program example
2
3 use pm_kind, only: SK, IK, RK => RKH ! All other real kinds are also supported.
4 use pm_quadpack, only: getQuadErr, GK41
5 use pm_quadpack, only: setNodeWeightGK
7 use pm_distNorm, only: getNormCDF
8 use pm_io, only: display_type
9 use pm_val2str, only: getStr
10
11 implicit none
12
13 integer(IK) , parameter :: NINTMAX = 1000_IK
14 real(RK) , parameter :: ABSTOL = 0._RK, RELTOL = epsilon(0._RK) * 10000
15 real(RK) :: lb, ub, truth, integral, abserr, sinfo(4, NINTMAX)
16 integer(IK) :: err, neval, nint, sindex(NINTMAX)
17
18 type(display_type) :: disp
19 disp = display_type(file = "main.out.F90")
20
21 call disp%skip()
22 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%show("! Compute the Gauss-Kronrod quadrature of the probability density function of the Standard Normal distribution.")
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%skip()
26
27 call disp%skip()
28 call disp%show("[ABSTOL, RELTOL]")
29 call disp%show( [ABSTOL, RELTOL] )
30 call disp%show("lb = -100._RK; ub = +100._RK")
31 lb = -100._RK; ub = +100._RK
32 call disp%show("err = getQuadErr(getNormPDF, lb, ub, ABSTOL, RELTOL, GK41, integral, abserr, sinfo, sindex, neval, nint)")
33 err = getQuadErr(getNormPDF, lb, ub, ABSTOL, RELTOL, GK41, integral, abserr, sinfo, sindex, neval, nint)
34 call disp%show("if (err /= 0_IK) error stop 'integration failed: err = '//getStr(err)")
35 if (err /= 0_IK) error stop 'integration failed: err = '//getStr(err)
36 call disp%show("truth = getNormCDF(ub) - getNormCDF(lb)")
37 truth = getNormCDF(ub) - getNormCDF(lb)
38 call disp%show("[truth, integral, abserr, abs(integral - truth) / truth]")
39 call disp%show( [truth, integral, abserr, abs(integral - truth) / truth] )
40 call disp%show("[neval, nint]")
41 call disp%show( [neval, nint] )
42 call disp%skip()
43
44 call disp%skip()
45 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
46 call disp%show("! Compute the Gauss-Kronrod quadrature of the probability density function of the Standard LogNormal distribution.")
47 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
48 call disp%skip()
49
50 call disp%skip()
51 call disp%show("[ABSTOL, RELTOL]")
52 call disp%show( [ABSTOL, RELTOL] )
53 call disp%show("lb = exp(-10._RK); ub = exp(+10._RK)")
54 lb = exp(-10._RK); ub = exp(+10._RK)
55 call disp%show("err = getQuadErr(getLogNormPDF, lb, ub, ABSTOL, RELTOL, GK41, integral, abserr, sinfo, sindex, neval, nint)")
56 err = getQuadErr(getLogNormPDF, lb, ub, ABSTOL, RELTOL, GK41, integral, abserr, sinfo, sindex, neval, nint)
57 call disp%show("if (err /= 0_IK) error stop 'integration failed: err = '//getStr(err)")
58 if (err /= 0_IK) error stop 'integration failed: err = '//getStr(err)
59 call disp%show("truth = getLogNormCDF(ub) - getLogNormCDF(lb)")
60 truth = getLogNormCDF(ub) - getLogNormCDF(lb)
61 call disp%show("[truth, integral, abserr, abs(integral - truth) / truth]")
62 call disp%show( [truth, integral, abserr, abs(integral - truth) / truth] )
63 call disp%show("[neval, nint]")
64 call disp%show( [neval, nint] )
65 call disp%skip()
66
67contains
68
69 function getIntSinCos(x) result(integrand)
70 real(RK) , intent(in) :: x
71 real(RK) :: integrand
72 integrand = cos(100._RK * sin(x))
73 end function
74
75 function getNormPDF(x) result(pdf)
76 use pm_distNorm, only: getNormLogPDF
77 real(RK) , intent(in) :: x
78 real(RK) :: pdf
79 pdf = exp(getNormLogPDF(x))
80 end function
81
82 function getLogNormPDF(x) result(pdf)
84 real(RK) , intent(in) :: x
85 real(RK) :: pdf
86 pdf = exp(getLogNormLogPDF(x))
87 end function
88
89end program example
Generate and return the Cumulative Distribution Function (CDF) of the univariate Lognormal distributi...
Generate the natural logarithm of probability density function (PDF) of the univariate Lognormal dist...
Generate and return the Cumulative Distribution Function (CDF) of the univariate Normal distribution.
Generate the natural logarithm of probability density function (PDF) of the univariate Normal distrib...
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
Compute the 1D integral of the input scalar (potentially singular) integrand getFunc on a finite or s...
Return the Kronrod nodes and weights of the extension to the -point Gauss-Legendre quadrature,...
Generate and return the conversion of the input value to an output Fortran string,...
Definition: pm_val2str.F90:167
This module contains classes and procedures for computing various statistical quantities related to t...
This module contains classes and procedures for computing various statistical quantities related to t...
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 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
type(GK41_type), parameter GK41
The scalar constant object of type GK41_type that indicates the use of 20-point Gauss-Legendre quadra...
This module contains the generic procedures for converting values of different types and kinds to For...
Definition: pm_val2str.F90:58
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 Gauss-Kronrod quadrature of the probability density function of the Standard Normal distribution.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7[ABSTOL, RELTOL]
8+0.00000000000000000000000000000000000, +0.192592994438723585305597794258492732E-29
9lb = -100._RK; ub = +100._RK
10err = getQuadErr(getNormPDF, lb, ub, ABSTOL, RELTOL, GK41, integral, abserr, sinfo, sindex, neval, nint)
11if (err /= 0_IK) error stop 'integration failed: err = '//getStr(err)
12truth = getNormCDF(ub) - getNormCDF(lb)
13[truth, integral, abserr, abs(integral - truth) / truth]
14+1.00000000000000000000000000000000000, +1.00000000000000000000000000000000019, +0.348826834552332163342274808104495182E-30, +0.192592994438723585305597794258492732E-33
15[neval, nint]
16+943, +12
17
18
19!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
20! Compute the Gauss-Kronrod quadrature of the probability density function of the Standard LogNormal distribution.
21!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
22
23
24[ABSTOL, RELTOL]
25+0.00000000000000000000000000000000000, +0.192592994438723585305597794258492732E-29
26lb = exp(-10._RK); ub = exp(+10._RK)
27err = getQuadErr(getLogNormPDF, lb, ub, ABSTOL, RELTOL, GK41, integral, abserr, sinfo, sindex, neval, nint)
28if (err /= 0_IK) error stop 'integration failed: err = '//getStr(err)
29truth = getLogNormCDF(ub) - getLogNormCDF(lb)
30[truth, integral, abserr, abs(integral - truth) / truth]
31+0.999999999999999999999984760293951620, +0.999999999999999999999984760293951620, +0.952513811497278733282817811682810727E-32, +0.00000000000000000000000000000000000
32[neval, nint]
33+2009, +25
34
35
Test:
test_pm_quadPack


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 at Austin

Definition at line 1264 of file pm_quadPack.F90.


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