Modernized QUADPACK: a Fortran subroutine package for the numerical computation of definite one-dimensional integrals
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public, | parameter | :: | quadpack_RK | = | wp |
the real kind used in this module |
1D globally adaptive integrator using Gauss-Kronrod quadrature, oscillating integrand
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested if epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28), the routine will end with ier = 6. |
||
integer, | intent(in) | :: | Key |
key for choice of local integration rule. a gauss-kronrod pair is used with: |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
integer, | intent(in) | :: | Limit |
dimensioning parameter for |
||
integer, | intent(in) | :: | Lenw |
dimensioning parameter for |
||
integer, | intent(out) | :: | Last |
on return, |
||
integer | :: | Iwork(Limit) |
vector of dimension at least |
|||
real(kind=wp) | :: | Work(Lenw) |
vector of dimension at least |
same as dqag but provides more information and control
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
integer, | intent(in) | :: | Key |
key for choice of local integration rule a gauss-kronrod pair is used with |
||
integer, | intent(in) | :: | Limit |
gives an upperbound on the number of subintervals
in the partition of |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Last |
number of subintervals actually produced in the subdivision process |
1D globally adaptive integrator, infinite intervals
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | Bound |
finite bound of integration range (has no meaning if interval is doubly-infinite) |
||
integer, | intent(in) | :: | Inf |
indicating the kind of integration range involved: |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
integer, | intent(in) | :: | Limit |
dimensioning parameter for |
||
integer, | intent(in) | :: | Lenw |
dimensioning parameter for |
||
integer, | intent(out) | :: | Last |
on return, |
||
integer | :: | Iwork(Limit) |
vector of dimension at least |
|||
real(kind=wp) | :: | Work(Lenw) |
vector of dimension at least |
same as dqagi but provides more information and control
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | Bound |
finite bound of integration range (has no meaning if interval is doubly-infinite) |
||
integer, | intent(in) | :: | Inf |
indicating the kind of integration range involved
* inf = 1 corresponds to |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
integer, | intent(in) | :: | Limit |
gives an upper bound on the number of subintervals
in the partition of |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension |
||
integer, | intent(out) | :: | Last |
number of subintervals actually produced in the subdivision process |
1D globally adaptive integrator, singularities or discontinuities
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
integer, | intent(in) | :: | Npts2 |
number equal to two more than the number of
user-supplied break points within the integration
range, |
||
real(kind=wp), | intent(in) | :: | Points(Npts2) |
vector of dimension npts2, the first |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
integer, | intent(in) | :: | Leniw |
dimensioning parameter for |
||
integer, | intent(in) | :: | Lenw |
dimensioning parameter for |
||
integer, | intent(out) | :: | Last |
on return, |
||
integer | :: | Iwork(Leniw) |
vector of dimension at least |
|||
real(kind=wp) | :: | Work(Lenw) |
vector of dimension at least |
same as dqagp but provides more information and control
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f | ||||
real(kind=wp) | :: | a | ||||
real(kind=wp) | :: | b | ||||
integer, | intent(in) | :: | Npts2 |
number equal to two more than the number of
user-supplied break points within the integration
range, |
||
real(kind=wp), | intent(in) | :: | Points(Npts2) |
vector of dimension npts2, the first (npts2-2) elements of which are the user provided break points. if these points do not constitute an ascending sequence there will be an automatic sorting. |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
integer, | intent(in) | :: | Limit |
gives an upper bound on the number of subintervals
in the partition of |
||
real(kind=wp) | :: | Result | ||||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Pts(Npts2) |
vector of dimension at least npts2, containing the integration limits and the break points of the interval in ascending sequence. |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Level(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Ndin(Npts2) |
vector of dimension at least npts2, after first
integration over the intervals |
||
integer, | intent(out) | :: | Last |
number of subintervals actually produced in the subdivisions process |
1D globally adaptive integrator using interval subdivision and extrapolation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand
function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
integer, | intent(in) | :: | Limit |
dimensioning parameter for |
||
integer, | intent(in) | :: | Lenw |
dimensioning parameter for |
||
integer, | intent(out) | :: | Last |
on return, |
||
integer | :: | Iwork(Limit) |
vector of dimension at least |
|||
real(kind=wp) | :: | Work(Lenw) |
vector of dimension at least |
same as dqags but provides more information and control
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand
function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
integer, | intent(in) | :: | Limit |
gives an upperbound on the number of subintervals
in the partition of |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages:
* ier = 1 maximum number of subdivisions allowed
has been achieved. one can allow more sub-
divisions by increasing the value of limit
(and taking the according dimension
adjustments into account). however, if
this yields no improvement it is advised
to analyze the integrand in order to
determine the integration difficulties. if
the position of a local difficulty can be
determined (e.g. singularity,
discontinuity within the interval) one
will probably gain from splitting up the
interval at this point and calling the
integrator on the subranges. if possible,
an appropriate special-purpose integrator
should be used, which is designed for
handling the type of difficulty involved.
* ier = 2 the occurrence of roundoff error is
detected, which prevents the requested
tolerance from being achieved.
the error may be under-estimated.
* ier = 3 extremely bad integrand behaviour
occurs at some points of the integration
interval.
* ier = 4 the algorithm does not converge.
roundoff error is detected in the
extrapolation table.
it is presumed that the requested
tolerance cannot be achieved, and that the
returned result is the best which can be
obtained.
* ier = 5 the integral is probably divergent, or
slowly convergent. it must be noted that
divergence can occur with any other value
of ier.
* ier = 6 the input is invalid, because
|
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Last |
number of subintervals actually produced in the subdivision process |
compute Cauchy principal value of f(x)/(x-c)
over a finite interval
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
under limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | c |
parameter in the weight function, |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate or the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages:
* ier = 1 maximum number of subdivisions allowed
has been achieved. one can allow more sub-
divisions by increasing the value of limit
(and taking the according dimension
adjustments into account). however, if
this yields no improvement it is advised
to analyze the integrand in order to
determine the integration difficulties.
if the position of a local difficulty
can be determined (e.g. singularity,
discontinuity within the interval) one
will probably gain from splitting up the
interval at this point and calling
appropriate integrators on the subranges.
* ier = 2 the occurrence of roundoff error is
detected, which prevents the requested
tolerance from being achieved.
* ier = 3 extremely bad integrand behaviour occurs
at some points of the integration
interval.
* ier = 6 the input is invalid, because
|
||
integer, | intent(in) | :: | Limit |
dimensioning parameter for |
||
integer, | intent(in) | :: | Lenw |
dimensioning parameter for |
||
integer, | intent(out) | :: | Last |
on return, |
||
integer | :: | Iwork(Limit) |
vector of dimension at least |
|||
real(kind=wp) | :: | Work(Lenw) |
vector of dimension at least |
same as dqawc but provides more information and control
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp) | :: | c | ||||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
integer, | intent(in) | :: | Limit |
gives an upper bound on the number of subintervals
in the partition of |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Last |
number of subintervals actually produced in the subdivision process |
Fourier sine/cosine transform for user supplied interval a
to infinity
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | Omega |
parameter in the integrand weight function |
||
integer, | intent(in) | :: | Integr |
indicates which of the weight functions is used: |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested, |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
integer, | intent(in) | :: | Limlst |
limlst gives an upper bound on the number of
cycles, |
||
integer, | intent(out) | :: | Lst |
on return, lst indicates the number of cycles
actually needed for the integration.
if |
||
integer, | intent(in) | :: | Leniw |
dimensioning parameter for |
||
integer, | intent(in) | :: | Maxp1 |
maxp1 gives an upper bound on the number of
chebyshev moments which can be stored, i.e. for
the intervals of lengths |
||
integer, | intent(in) | :: | Lenw |
dimensioning parameter for |
||
integer | :: | Iwork(Leniw) |
vector of dimension at least |
|||
real(kind=wp) | :: | Work(Lenw) |
vector of dimension at least |
same as dqawf but provides more information and control
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | Omega |
parameter in the weight function |
||
integer, | intent(in) | :: | Integr |
indicates which weight function is used: |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested, |
||
integer, | intent(in) | :: | Limlst |
limlst gives an upper bound on the number of
cycles, |
||
integer, | intent(in) | :: | Limit |
gives an upper bound on the number of subintervals
allowed in the partition of each cycle, |
||
integer, | intent(in) | :: | Maxp1 |
gives an upper bound on the number of
chebyshev moments which can be stored, i.e.
for the intervals of lengths
|
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
real(kind=wp), | intent(out) | :: | Rslst(Limlst) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Erlst(Limlst) |
vector of dimension at least |
||
integer, | intent(out) | :: | Ierlst(Limlst) |
vector of dimension at least |
||
integer, | intent(out) | :: | Lst |
number of subintervals needed for the integration
if |
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Nnlog(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Chebmo(Maxp1,25) |
array of dimension at least |
1D integration of cos(omega*x)*f(x)
or sin(omega*x)*f(x)
over a finite interval, adaptive subdivision with extrapolation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | Omega |
parameter in the integrand weight function |
||
integer, | intent(in) | :: | Integr |
indicates which of the weight functions is used |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
integer, | intent(in) | :: | Leniw |
dimensioning parameter for |
||
integer, | intent(in) | :: | Maxp1 |
gives an upper bound on the number of chebyshev
moments which can be stored, i.e. for the
intervals of lengths |
||
integer, | intent(in) | :: | Lenw |
dimensioning parameter for |
||
integer, | intent(out) | :: | Last |
on return, |
||
integer | :: | Iwork(Leniw) |
vector of dimension at least leniw
on return, the first |
|||
real(kind=wp) | :: | Work(Lenw) |
vector of dimension at least |
same as dqawo but provides more information and control
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | Omega |
parameter in the integrand weight function |
||
integer, | intent(in) | :: | Integr |
indicates which of the weight functions is to be used: |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested.
if |
||
integer, | intent(in) | :: | Limit |
gives an upper bound on the number of subdivisions
in the partition of |
||
integer, | intent(in) | :: | Icall |
if dqawoe is to be used only once, icall must
be set to 1. assume that during this call, the
chebyshev moments (for clenshaw-curtis integration
of degree 24) have been computed for intervals of
lengths |
||
integer, | intent(in) | :: | Maxp1 |
gives an upper bound on the number of chebyshev
moments which can be stored, i.e. for the
intervals of lengths |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
||
integer, | intent(out) | :: | Last |
on return, |
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Nnlog(Limit) |
vector of dimension at least |
||
integer, | intent(inout) | :: | Momcom |
indicating that the chebyshev moments
have been computed for intervals of lengths
|
||
real(kind=wp), | intent(inout) | :: | Chebmo(Maxp1,25) |
array of dimension |
1D integration of functions with powers and or logs over a finite interval
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration, b>a if b<=a, the routine will end with ier = 6. |
||
real(kind=wp), | intent(in) | :: | alfa |
parameter in the integrand function, |
||
real(kind=wp), | intent(in) | :: | beta |
parameter in the integrand function, |
||
integer, | intent(in) | :: | integr |
indicates which weight function is to be used: |
||
real(kind=wp), | intent(in) | :: | epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | epsrel |
relative accuracy requested.
if |
||
real(kind=wp), | intent(out) | :: | result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | ier |
error messages: |
||
integer, | intent(in) | :: | limit |
dimensioning parameter for |
||
integer, | intent(in) | :: | lenw |
dimensioning parameter for |
||
integer, | intent(out) | :: | last |
on return, |
||
integer | :: | iwork(limit) |
vector of dimension limit, the first |
|||
real(kind=wp) | :: | work(lenw) |
on return: |
same as dqaws but provides more information and control
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration, |
||
real(kind=wp), | intent(in) | :: | Alfa |
parameter in the weight function, |
||
real(kind=wp), | intent(in) | :: | Beta |
parameter in the weight function, |
||
integer, | intent(in) | :: | Integr |
indicates which weight function is to be used: |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested.
if |
||
integer, | intent(in) | :: | Limit |
gives an upper bound on the number of subintervals
in the partition of |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
|
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Last |
number of subintervals actually produced in the subdivision process |
1D integral for Cauchy principal values using a 25 point quadrature rule
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
left end point of the integration interval |
||
real(kind=wp), | intent(in) | :: | b |
right end point of the integration interval, |
||
real(kind=wp), | intent(in) | :: | c |
parameter in the weight function |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral.
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(inout) | :: | Krul |
key which is decreased by 1 if the 15-point gauss-kronrod scheme has been used |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
1D integral for sin/cos integrand using a 25 point quadrature rule
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand
function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | Omega |
parameter in the weight function |
||
integer, | intent(in) | :: | Integr |
indicates which weight function is to be used |
||
integer, | intent(in) | :: | Nrmom |
the length of interval |
||
integer, | intent(in) | :: | Maxp1 |
gives an upper bound on the number of chebyshev
moments which can be stored, i.e. for the
intervals of lengths |
||
integer, | intent(in) | :: | Ksave |
key which is one when the moments for the current interval have been computed |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute
error, which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral j |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of |
||
integer, | intent(inout) | :: | Momcom |
for each interval length we need to compute the
chebyshev moments. momcom counts the number of
intervals for which these moments have already been
computed. if |
||
real(kind=wp), | intent(inout) | :: | Chebmo(Maxp1,25) |
array of dimension at least |
25-point clenshaw-curtis integration
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand f(x). |
|||
real(kind=wp), | intent(in) | :: | a |
left end point of the original interval |
||
real(kind=wp), | intent(in) | :: | b |
right end point of the original interval, |
||
real(kind=wp), | intent(in) | :: | Bl |
lower limit of integration, |
||
real(kind=wp), | intent(in) | :: | Br |
upper limit of integration, |
||
real(kind=wp), | intent(in) | :: | Alfa |
parameter in the weight function |
||
real(kind=wp), | intent(in) | :: | Beta |
parameter in the weight function |
||
real(kind=wp), | intent(in) | :: | Ri(25) |
modified chebyshev moments for the application of the generalized clenshaw-curtis method (computed in subroutine dqmomo) |
||
real(kind=wp), | intent(in) | :: | Rj(25) |
modified chebyshev moments for the application of the generalized clenshaw-curtis method (computed in subroutine dqmomo) |
||
real(kind=wp), | intent(in) | :: | Rg(25) |
modified chebyshev moments for the application of the generalized clenshaw-curtis method (computed in subroutine dqmomo) |
||
real(kind=wp), | intent(in) | :: | Rh(25) |
modified chebyshev moments for the application of the generalized clenshaw-curtis method (computed in subroutine dqmomo) |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of abs(f*w-i/(b-a)) |
||
integer, | intent(in) | :: | Integr |
which determines the weight function
* = 1 |
||
integer, | intent(out) | :: | Nev |
number of integrand evaluations |
chebyshev series expansion
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in) | :: | x(11) |
vector of dimension 11 containing the
values |
||
real(kind=wp), | intent(inout) | :: | Fval(25) |
vector of dimension 25 containing the
function values at the points
|
||
real(kind=wp), | intent(out) | :: | Cheb12(13) |
vector of dimension 13 containing the chebyshev coefficients for degree 12 |
||
real(kind=wp), | intent(out) | :: | Cheb24(25) |
vector of dimension 25 containing the chebyshev coefficients for degree 24 |
estimate 1D integral on finite interval using a 15 point gauss-kronrod rule and give error estimate, non-automatic
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should not exceed |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral j |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of |
estimate 1D integral on (semi)infinite interval using a 15 point gauss-kronrod quadrature rule, non-automatic
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | Boun |
finite bound of original integration range (set to zero if inf = +2) |
||
integer, | intent(in) | :: | Inf |
the integral is computed as the sum of two
integrals, one over |
||
real(kind=wp), | intent(in) | :: | a |
lower limit for integration over subrange of |
||
real(kind=wp), | intent(in) | :: | b |
upper limit for integration over subrange of |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i.
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral j |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of
|
estimate 1D integral with special singular weight functions using a 15 point gauss-kronrod quadrature rule
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
procedure(weight_func) | :: | w |
function subprogram defining the integrand weight function |
|||
real(kind=wp), | intent(in) | :: | p1 |
parameter in the weight function |
||
real(kind=wp), | intent(in) | :: | p2 |
parameter in the weight function |
||
real(kind=wp), | intent(in) | :: | p3 |
parameter in the weight function |
||
real(kind=wp), | intent(in) | :: | p4 |
parameter in the weight function |
||
integer, | intent(in) | :: | Kp |
key for indicating the type of weight function |
||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral of |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of |
estimate 1D integral on finite interval using a 21 point gauss-kronrod rule and give error estimate, non-automatic
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should not exceed |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral j |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of |
estimate 1D integral on finite interval using a 31 point gauss-kronrod rule and give error estimate, non-automatic
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the modulus,
which should not exceed |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral j |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of |
estimate 1D integral on finite interval using a 41 point gauss-kronrod rule and give error estimate, non-automatic
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should not exceed |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral j |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of abs(f-i/(b-a))
over |
estimate 1D integral on finite interval using a 51 point gauss-kronrod rule and give error estimate, non-automatic
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subroutine defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i.
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should not exceed |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral j |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of |
estimate 1D integral on finite interval using a 61 point gauss-kronrod rule and give error estimate, non-automatic
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i
|
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
real(kind=wp), | intent(out) | :: | Resabs |
approximation to the integral j |
||
real(kind=wp), | intent(out) | :: | Resasc |
approximation to the integral of |
1D integration of k
-th degree Chebyshev polynomial times a function with singularities
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in) | :: | Alfa |
parameter in the weight function |
||
real(kind=wp), | intent(in) | :: | Beta |
parameter in the weight function |
||
real(kind=wp), | intent(out) | :: | Ri(25) |
|
||
real(kind=wp), | intent(out) | :: | Rj(25) |
|
||
real(kind=wp), | intent(out) | :: | Rg(25) |
|
||
real(kind=wp), | intent(out) | :: | Rh(25) |
|
||
integer, | intent(in) | :: | Integr |
input parameter indicating the modified moments to be computed: |
1D non-adaptive automatic integrator
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
real(kind=wp), | intent(out) | :: | Result |
approximation to the integral i result is obtained by applying the 21-point gauss-kronrod rule (res21) obtained by optimal addition of abscissae to the 10-point gauss rule (res10), or by applying the 43-point rule (res43) obtained by optimal addition of abscissae to the 21-point gauss-kronrod rule, or by applying the 87-point rule (res87) obtained by optimal addition of abscissae to the 43-point rule. |
||
real(kind=wp), | intent(out) | :: | Abserr |
estimate of the modulus of the absolute error,
which should equal or exceed |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages: |
This subroutine attempts to calculate the integral of f(x)
over the interval a
to b
with relative error not
exceeding epsil
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration. |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration. |
||
real(kind=wp), | intent(out) | :: | result |
the value of the integral to the specified relative accuracy. |
||
real(kind=wp), | intent(in) | :: | epsil |
relative accuracy required. when the relative
difference of two successive formulae does not
exceed |
||
integer, | intent(out) | :: | npts |
number integrand evaluations. |
||
integer, | intent(out) | :: | icheck |
on exit normally |
Integrate a function tabulated at arbitrarily spaced abscissas using overlapping parabolas.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in), | dimension(:) | :: | x |
array of abscissas, which must be in increasing order. |
|
real(kind=wp), | intent(in), | dimension(:) | :: | y |
array of function values. i.e., |
|
integer, | intent(in) | :: | n |
The integer number of function values supplied.
|
||
real(kind=wp), | intent(in) | :: | xlo |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | xup |
upper limit of integration. Must have |
||
real(kind=wp), | intent(out) | :: | ans |
computed approximate value of integral |
||
integer, | intent(out) | :: | ierr |
A status code: |
Integrate a function using a 7-point adaptive Newton-Cotes quadrature rule.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | fun |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration (may be less than |
||
real(kind=wp), | intent(in) | :: | err |
a requested error tolerance. Normally, pick a value
|
||
real(kind=wp), | intent(out) | :: | ans |
computed value of the integral. Hopefully, |
||
integer, | intent(out) | :: | ierr |
a status code: |
||
integer, | intent(out) | :: | k |
the number of function evaluations actually used to do
the integration. A value of |
Integrate a real function of one variable over a finite interval using an adaptive 8-point Legendre-Gauss algorithm.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower bound of the integration |
||
real(kind=wp), | intent(in) | :: | b |
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
|
||
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 |
Numerically evaluate integral using adaptive Simpson rule.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower bound of the integration |
||
real(kind=wp), | intent(in) | :: | b |
upper bound of the integration |
||
real(kind=wp), | intent(in) | :: | error_tol |
relative error tolerance |
||
real(kind=wp), | intent(out) | :: | ans |
computed value of integral |
||
integer, | intent(out) | :: | ierr |
status code: |
Numerically evaluate integral using adaptive Lobatto rule
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower bound of the integration |
||
real(kind=wp), | intent(in) | :: | b |
upper bound of the integration |
||
real(kind=wp), | intent(in) | :: | error_tol |
relative error tolerance |
||
real(kind=wp), | intent(out) | :: | ans |
computed value of integral |
||
integer, | intent(out) | :: | ierr |
status code: |