single integration class: for 1d integration of f(x)
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
procedure(func_1d), | public, | pointer | :: | fun | => | null() | function |
procedure(gauss_func), | public, | pointer | :: | g | => | null() | the guass quadrature formula to use |
real(kind=wp), | public | :: | a | = | zero | lower limit of integration |
|
real(kind=wp), | public | :: | b | = | zero | upper limit of integration (may be less than a) |
|
real(kind=wp), | public | :: | tol | = | zero | the requested relative error tolerance. |
|
real(kind=wp), | public | :: | val | = | zero | the value of |
core integration routine. refactored from SLATEC with selectable quadrature method
Integrate a real function of one variable over a finite interval using the specified adaptive algorithm. Intended primarily for high accuracy integration or integration of smooth functions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integration_class_1d), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | lb | lower bound of the integration |
||
real(kind=wp), | intent(in) | :: | ub | upper bound of the integration |
||
real(kind=wp), | intent(in) | :: | error_tol | is a requested pseudorelative error tolerance. normally pick a value of abs(error_tol) so that dtol < abs(error_tol) <= 1.0e-3 where dtol is the larger of 1.0e-18 and the real(wp) unit roundoff d1mach(4). ans will normally have no more error than abs(error_tol) times the integral of the absolute value of fun(x). usually, smaller values of error_tol yield more accuracy and require more function evaluations. |
||
real(kind=wp), | intent(out) | :: | ans | computed value of integral |
||
integer, | intent(out) | :: | ierr | status code:
|
||
real(kind=wp), | intent(out) | :: | err | an estimate of the absolute error in |
to set up the class
Initialize the 1D integration class. Must be called before integration is performed.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integration_class_1d), | intent(inout) | :: | me | |||
procedure(func_1d) | :: | fx | 1d function: f(x) |
|||
real(kind=wp), | intent(in) | :: | xl | x integration lower bound |
||
real(kind=wp), | intent(in) | :: | xu | x integration upper bound |
||
real(kind=wp), | intent(in) | :: | tolx | error tolerance for dx integration |
||
integer, | intent(in) | :: | methodx | quadrature method to use for x |
to integrate the function fun
Perform the 1D integration.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integration_class_1d), | intent(inout) | :: | me | |||
real(kind=wp), | intent(out) | :: | ans | |||
integer, | intent(out) | :: | ierr | |||
real(kind=wp), | intent(out) | :: | err |
type,extends(integration_class),public :: integration_class_1d
!! single integration class: for 1d integration of `f(x)`
private
procedure(func_1d),pointer :: fun => null() !! function `f(x)` to be integrated
procedure(gauss_func),pointer :: g => null() !! the guass quadrature formula to use
real(wp) :: a = zero !! lower limit of integration
real(wp) :: b = zero !! upper limit of integration (may be less than a)
real(wp) :: tol = zero !! the requested relative error tolerance.
real(wp) :: val = zero !! the value of `x`. Only used for multiple
!! integration to pass to the inner integrals
contains
private
procedure :: dgauss_generic !! core integration routine. refactored from
!! SLATEC with selectable quadrature method
procedure,public :: initialize => initialize_integration_class !! to set up the class
procedure,public :: integrate => integrate_1d !! to integrate the function `fun`
end type integration_class_1d