dqags Subroutine

public subroutine dqags(f, a, b, Epsabs, Epsrel, Result, Abserr, Neval, Ier, Limit, Lenw, Last, Iwork, Work)

1D globally adaptive integrator using interval subdivision and extrapolation

the routine calculates an approximation result to a given definite integral i = integral of f over (a,b), hopefully satisfying following claim for accuracy abs(i-result)<=max(epsabs,epsrel*abs(i)).

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

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 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 (epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28) or limit<1 or lenw<limit*4. result, abserr, neval, last are set to zero. except when limit or lenw is invalid, iwork(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.
integer, intent(in) :: Limit

dimensioning parameter for iwork. limit determines the maximum number of subintervals in the partition of the given integration interval (a,b), limit>=1. if limit<1, the routine will end with ier = 6.

integer, intent(in) :: Lenw

dimensioning parameter for work. lenw must be at least limit*4. if lenw<limit*4, the routine will end with ier = 6.

integer, intent(out) :: Last

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

integer :: Iwork(Limit)

vector of dimension at least limit, 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

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 error estimates.

Calls

proc~~dqags~~CallsGraph proc~dqags quadpack_generic::dqags proc~dqagse quadpack_generic::dqagse proc~dqags->proc~dqagse proc~dqk21 quadpack_generic::dqk21 proc~dqagse->proc~dqk21

Source Code

    subroutine dqags(f, a, b, Epsabs, Epsrel, Result, Abserr, Neval, Ier, &
                     Limit, 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
        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 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
                                    !!   `(epsabs<=0` and
                                    !!   `epsrel<max(50*rel.mach.acc.,0.5e-28)`
                                    !!   or `limit<1` or `lenw<limit*4`.
                                    !!   `result`, `abserr`, `neval`, `last` are set to
                                    !!   zero. except when limit or lenw is invalid,
                                    !!   `iwork(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`.
        integer, intent(in) :: Limit !! dimensioning parameter for `iwork`.
                                     !! `limit` determines the maximum number of subintervals
                                     !! in the partition of the given integration interval
                                     !! `(a,b)`, `limit>=1`.
                                     !! if `limit<1`, the routine will end with ier = 6.
        integer, intent(in) :: Lenw !! dimensioning parameter for `work`.
                                    !! `lenw` must be at least `limit*4`.
                                    !! if `lenw<limit*4`, the routine will end
                                    !! with ier = 6.
        integer, intent(out) :: Last !! on return, `last` equals the number of subintervals
                                     !! produced in the subdivision process, determines the
                                     !! number of significant elements actually in the `work`
                                     !! arrays.
        integer :: Iwork(Limit) !! vector of dimension at least `limit`, 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
        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 error estimates.

        integer :: lvl, l1, l2, l3

        ! check validity of limit and lenw.

        Ier = 6
        Neval = 0
        Last = 0
        Result = 0.0_wp
        Abserr = 0.0_wp
        if (Limit >= 1 .and. Lenw >= Limit*4) then

            ! prepare call for dqagse.
            l1 = Limit + 1
            l2 = Limit + l1
            l3 = Limit + l2

            call dqagse(f, a, b, Epsabs, Epsrel, Limit, Result, Abserr, Neval, Ier, &
                        Work(1), Work(l1), Work(l2), Work(l3), Iwork, Last)

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

    end subroutine dqags