dqagp Subroutine

public subroutine dqagp(f, a, b, Npts2, Points, Epsabs, Epsrel, Result, Abserr, Neval, Ier, Leniw, Lenw, Last, Iwork, Work)

1D globally adaptive integrator, singularities or discontinuities

the routine calculates an approximation result to a given definite integral i = integral of f over (a,b), hopefully satisfying following claim for accuracy break points of the integration interval, where local difficulties of the integrand may occur (e.g. singularities, discontinuities), are provided by the user.

History

  • QUADPACK: date written 800101, revision date 830518 (yymmdd)

Arguments

Type IntentOptional Attributes Name
procedure(func) :: f

function subprogram defining the integrand function f(x).

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, npts>=2. if npts2<2, the routine will end with ier = 6.

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 epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28), the routine will end with ier = 6.

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 abs(i-result)

integer, intent(out) :: Neval

number of integrand evaluations

integer, intent(out) :: Ier
  • ier = 0 normal and reliable termination of the routine. it is assumed that the requested accuracy has been achieved.
  • ier>0 abnormal termination of the routine. the estimates for integral and error are less reliable. it is assumed that the requested accuracy has not been achieved.

error messages:

  • ier = 1 maximum number of subdivisions allowed has been achieved. one can allow more subdivisions 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 (i.e. singularity, discontinuity within the interval), it should be supplied to the routine as an element of the vector points. if necessary an appropriate special-purpose integrator must 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>0.
  • ier = 6 the input is invalid because npts2<2 or break points are specified outside the integration range or (epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28)) result, abserr, neval, last are set to zero. except when leniw or lenw or npts2 is invalid, iwork(1), iwork(limit+1), work(limit*2+1) and work(limit*3+1) are set to zero. work(1) is set to a and work(limit+1) to b (where limit = (leniw-npts2)/2).
integer, intent(in) :: Leniw

dimensioning parameter for iwork. leniw determines limit = (leniw-npts2)/2, which is the maximum number of subintervals in the partition of the given integration interval (a,b), leniw>=(3*npts2-2). if leniw<(3*npts2-2), the routine will end with ier = 6.

integer, intent(in) :: Lenw

dimensioning parameter for work. lenw must be at least leniw*2-npts2. if lenw<leniw*2-npts2, the routine will end with ier = 6.

integer, intent(out) :: Last

on return, last equals the number of subintervals produced in the subdivision process, which determines the number of significant elements actually in the work arrays.

integer :: Iwork(Leniw)

vector of dimension at least leniw. on return, the first k elements of which contain pointers to the error estimates over the subintervals, such that work(limit*3+iwork(1)),...,work(limit*3+iwork(k)) form a decreasing sequence, with k = last if last<=(limit/2+2), and k = limit+1-last otherwise iwork(limit+1), ...,iwork(limit+last) contain the subdivision levels of the subintervals, i.e. if (aa,bb) is a subinterval of (p1,p2) where p1 as well as p2 is a user-provided break point or integration limit, then (aa,bb) has level l if abs(bb-aa) = abs(p2-p1)*2**(-l), iwork(limit*2+1), ..., iwork(limit*2+npts2) have no significance for the user, note that limit = (leniw-npts2)/2.

real(kind=wp) :: Work(Lenw)

vector of dimension at least lenw. on return:

  • work(1), ..., work(last) contain the left end points of the subintervals in the partition of (a,b),
  • work(limit+1), ..., work(limit+last) contain the right end points,
  • work(limit*2+1), ..., work(limit*2+last) contain the integral approximations over the subintervals,
  • work(limit*3+1), ..., work(limit*3+last) contain the corresponding error estimates,
  • work(limit*4+1), ..., work(limit*4+npts2) contain the integration limits and the break points sorted in an ascending sequence.

note that limit = (leniw-npts2)/2.


Calls

proc~~dqagp~~CallsGraph proc~dqagp quadpack_generic::dqagp proc~dqagpe quadpack_generic::dqagpe proc~dqagp->proc~dqagpe proc~dqk21 quadpack_generic::dqk21 proc~dqagpe->proc~dqk21

Source Code

    subroutine dqagp(f, a, b, Npts2, Points, Epsabs, Epsrel, Result, Abserr, &
                     Neval, Ier, Leniw, Lenw, Last, Iwork, Work)
        implicit none

        procedure(func) :: f !! function subprogram defining the integrand function `f(x)`.
        real(wp), intent(in) :: a !! lower limit of integration
        real(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, `npts>=2`.
                                     !! if `npts2<2`, the routine will end with ier = 6.
        real(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(wp), intent(in) :: Epsabs !! absolute accuracy requested
        real(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.
        real(wp), intent(out) :: Result !! approximation to the integral
        real(wp), intent(out) :: Abserr !! estimate of the modulus of the absolute error,
                                        !! which should equal or exceed `abs(i-result)`
        integer, intent(out) :: Neval !! number of integrand evaluations
        integer, intent(out) :: Ier !! * ier = 0 normal and reliable termination of the
                                    !!   routine. it is assumed that the requested
                                    !!   accuracy has been achieved.
                                    !! * ier>0 abnormal termination of the routine.
                                    !!   the estimates for integral and error are
                                    !!   less reliable. it is assumed that the
                                    !!   requested accuracy has not been achieved.
                                    !!
                                    !! error messages:
                                    !!
                                    !! * ier = 1 maximum number of subdivisions allowed
                                    !!   has been achieved. one can allow more
                                    !!   subdivisions 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 (i.e. singularity,
                                    !!   discontinuity within the interval), it
                                    !!   should be supplied to the routine as an
                                    !!   element of the vector points. if necessary
                                    !!   an appropriate special-purpose integrator
                                    !!   must 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>0.
                                    !! * ier = 6 the input is invalid because
                                    !!   `npts2<2` or
                                    !!   break points are specified outside
                                    !!   the integration range or
                                    !!   `(epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28))`
                                    !!   `result`, `abserr`, `neval`, `last` are set to
                                    !!   zero. except when `leniw` or `lenw` or `npts2` is
                                    !!   invalid, `iwork(1)`, `iwork(limit+1)`,
                                    !!   `work(limit*2+1) and work(limit*3+1)`
                                    !!   are set to zero.
                                    !!   `work(1)` is set to a and `work(limit+1)`
                                    !!   to `b` (where `limit = (leniw-npts2)/2`).
        integer, intent(in) :: Leniw !! dimensioning parameter for `iwork`.
                                     !! `leniw` determines `limit = (leniw-npts2)/2`,
                                     !! which is the maximum number of subintervals in the
                                     !! partition of the given integration interval `(a,b)`,
                                     !! `leniw>=(3*npts2-2)`.
                                     !! if `leniw<(3*npts2-2)`, the routine will end with
                                     !! ier = 6.
        integer, intent(in) :: Lenw !! dimensioning parameter for `work`.
                                    !! `lenw` must be at least `leniw*2-npts2`.
                                    !! if `lenw<leniw*2-npts2`, the routine will end
                                    !! with ier = 6.
        integer, intent(out) :: Last !! on return, `last` equals the number of subintervals
                                     !! produced in the subdivision process, which
                                     !! determines the number of significant elements
                                     !! actually in the work arrays.
        integer :: Iwork(Leniw) !! vector of dimension at least `leniw`. on return,
                                !! the first `k` elements of which contain
                                !! pointers to the error estimates over the
                                !! subintervals, such that
                                !! `work(limit*3+iwork(1)),...,work(limit*3+iwork(k))`
                                !! form a decreasing
                                !! sequence, with `k = last` if `last<=(limit/2+2)`, and
                                !! `k = limit+1-last` otherwise
                                !! `iwork(limit+1), ...,iwork(limit+last)` contain the
                                !! subdivision levels of the subintervals, i.e.
                                !! if `(aa,bb)` is a subinterval of `(p1,p2)`
                                !! where `p1` as well as `p2` is a user-provided
                                !! break point or integration limit, then `(aa,bb)` has
                                !! level `l` if `abs(bb-aa) = abs(p2-p1)*2**(-l)`,
                                !! `iwork(limit*2+1), ..., iwork(limit*2+npts2)` have
                                !! no significance for the user,
                                !! note that `limit = (leniw-npts2)/2`.
        real(wp) :: Work(Lenw) !! vector of dimension at least `lenw`.
                               !! on return:
                               !!
                               !! * `work(1), ..., work(last)` contain the left
                               !!   end points of the subintervals in the
                               !!   partition of `(a,b)`,
                               !! * `work(limit+1), ..., work(limit+last)` contain
                               !!   the right end points,
                               !! * `work(limit*2+1), ..., work(limit*2+last)` contain
                               !!   the integral approximations over the subintervals,
                               !! * `work(limit*3+1), ..., work(limit*3+last)`
                               !!   contain the corresponding error estimates,
                               !! * `work(limit*4+1), ..., work(limit*4+npts2)`
                               !!   contain the integration limits and the
                               !!   break points sorted in an ascending sequence.
                               !!
                               !! note that `limit = (leniw-npts2)/2`.

        integer :: limit, lvl, l1, l2, l3, l4

        ! check validity of limit and lenw.
        Ier = 6
        Neval = 0
        Last = 0
        Result = 0.0_wp
        Abserr = 0.0_wp
        if (Leniw >= (3*Npts2 - 2) .and. Lenw >= (Leniw*2 - Npts2) .and. Npts2 >= 2) then

            ! prepare call for dqagpe.
            limit = (Leniw - Npts2)/2
            l1 = limit + 1
            l2 = limit + l1
            l3 = limit + l2
            l4 = limit + l3

            call dqagpe(f, a, b, Npts2, Points, Epsabs, Epsrel, limit, Result, &
                        Abserr, Neval, Ier, Work(1), Work(l1), Work(l2), Work(l3), &
                        Work(l4), Iwork(1), Iwork(l1), Iwork(l2), Last)

            ! call error handler if necessary.
            lvl = 0
        end if
        if (Ier == 6) lvl = 1
        if (Ier /= 0) call xerror('abnormal return from dqagp', Ier, lvl)

    end subroutine dqagp