ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation. |
This is the derived type for generating test integrand objects of algebraic form as described below. More...
Public Member Functions | |
procedure | get => getIntDoncker2 |
Public Member Functions inherited from pm_quadTest::integrand_type | |
procedure(get_proc), deferred | get |
The function member returning the value of the unweighted integrand (whether Cauchy/sin/cos/algebraically types of weights) at a specified input point x . More... | |
Additional Inherited Members | |
Public Attributes inherited from pm_quadTest::integrand_type | |
real(RKH) | lb |
The scalar of type real of the highest kind supported by the library RKH, containing the lower limit of integration. More... | |
real(RKH) | ub |
The scalar of type real of the highest kind supported by the library RKH, containing the upper limit of integration. More... | |
real(RKH) | integral |
The scalar of type real of the highest kind supported by the library RKH, containing the true result of integration. More... | |
real(RKH), dimension(:), allocatable | break |
The scalar of type real of the highest kind supported by the library RKH, containing the points of difficulties of integration. More... | |
type(wcauchy_type), allocatable | wcauchy |
The scalar of type wcauchy_type, containing the Cauchy singularity of the integrand. More... | |
character(:, SK), allocatable | desc |
The scalar allocatable character of default kind SK containing a description of the integrand and integration limits and difficulties. More... | |
This is the derived type for generating test integrand objects of algebraic form as described below.
The full integrand is defined over a finite interval as,
\begin{equation} f(x) = \frac{\exp(x)}{\sqrt{-x}} ~,~ x \in (\ms{lb}, \ms{ub} < 0.) ~. \end{equation}
The integrand has the precise value,
\begin{equation} \int_{\ms{lb} = -\infty}^{\ms{ub}} f(x) dx = \sqrt{\pi} \big(\mathrm{erf}(\sqrt{-\ms{lb}} - \mathrm{erf}(\sqrt{-\ms{ub}}) \big) ~. \end{equation}
This integrand is an extension of the example discussed in Doncker et al (1976), Automatic Computation of Integrals with Singular integrand.
[in] | lb | : The input negative scalar of type real of kind RKH, containing the lower limit of integration.(optional, default = getInfNeg(0._RKG)) |
[in] | ub | : The input positive scalar of the same type and kind as lb , containing the upper limit of integration.(optional, default = 0._RK ) |
Possible calling interfaces ⛓
lb < ub < 0._RK
must for the relevant input arguments.CHECK_ENABLED=1
.
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.
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.
Definition at line 1589 of file pm_quadTest.F90.
procedure pm_quadTest::intDoncker2_type::get |
Definition at line 1591 of file pm_quadTest.F90.
References pm_kind::RKH.