dqage Subroutine

public subroutine dqage(f, a, b, Epsabs, Epsrel, Key, Limit, Result, Abserr, Neval, Ier, Alist, Blist, Rlist, Elist, Iord, Last)

same as dqag but provides more information and control

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-reslt)<=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.

integer, intent(in) :: Key

key for choice of local integration rule a gauss-kronrod pair is used with

  • 7 - 15 points if key<2,
  • 10 - 21 points if key = 2,
  • 15 - 31 points if key = 3,
  • 20 - 41 points if key = 4,
  • 25 - 51 points if key = 5,
  • 30 - 61 points if key>5.
integer, intent(in) :: Limit

gives an upperbound on the number of subintervals in the partition of (a,b), limit>=1.

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 result 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. however, if this yields no improvement it is rather 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.
  • ier = 3 extremely bad integrand behaviour occurs at some points of the integration interval.
  • ier = 6 the input is invalid, because (epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28_wp), result, abserr, neval, last, rlist(1) , elist(1) and iord(1) are set to zero. alist(1) and blist(1) are set to a and b respectively.
real(kind=wp), intent(out) :: Alist(Limit)

vector of dimension at least limit, the first last elements of which are the left end points of the subintervals in the partition of the given integration range (a,b)

real(kind=wp), intent(out) :: Blist(Limit)

vector of dimension at least limit, the first last elements of which are the right end points of the subintervals in the partition of the given integration range (a,b)

real(kind=wp), intent(out) :: Rlist(Limit)

vector of dimension at least limit, the first last elements of which are the integral approximations on the subintervals

real(kind=wp), intent(out) :: Elist(Limit)

vector of dimension at least limit, the first last elements of which are the moduli of the absolute error estimates on the subintervals

integer, intent(out) :: Iord(Limit)

vector of dimension at least limit, the first k elements of which are pointers to the error estimates over the subintervals, such that elist(iord(1)), ..., elist(iord(k)) form a decreasing sequence, with k = last if last<=(limit/2+2), and k = limit+1-last otherwise

integer, intent(out) :: Last

number of subintervals actually produced in the subdivision process


Calls

proc~~dqage~~CallsGraph proc~dqage quadpack_generic::dqage proc~dqk15 quadpack_generic::dqk15 proc~dqage->proc~dqk15 proc~dqk21 quadpack_generic::dqk21 proc~dqage->proc~dqk21 proc~dqk31 quadpack_generic::dqk31 proc~dqage->proc~dqk31 proc~dqk41 quadpack_generic::dqk41 proc~dqage->proc~dqk41 proc~dqk51 quadpack_generic::dqk51 proc~dqage->proc~dqk51 proc~dqk61 quadpack_generic::dqk61 proc~dqage->proc~dqk61

Called by

proc~~dqage~~CalledByGraph proc~dqage quadpack_generic::dqage proc~dqag quadpack_generic::dqag proc~dqag->proc~dqage

Source Code

    subroutine dqage(f, a, b, Epsabs, Epsrel, Key, Limit, Result, Abserr, &
                     Neval, Ier, Alist, Blist, Rlist, Elist, Iord, Last)
        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.
        integer, intent(in) :: Key !! key for choice of local integration rule
                                   !!  a gauss-kronrod pair is used with
                                   !!
                                   !!  * 7 - 15 points if key<2,
                                   !!  * 10 - 21 points if key = 2,
                                   !!  * 15 - 31 points if key = 3,
                                   !!  * 20 - 41 points if key = 4,
                                   !!  * 25 - 51 points if key = 5,
                                   !!  * 30 - 61 points if key>5.
        integer, intent(in) :: Limit !! gives an upperbound on the number of subintervals
                                     !! in the partition of `(a,b)`, `limit>=1`.
        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 result 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.
                                    !!          however, if this yields no improvement it
                                    !!          is rather 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.
                                    !!  * ier = 3 extremely bad integrand behaviour occurs
                                    !!          at some points of the integration
                                    !!          interval.
                                    !!  * ier = 6 the input is invalid, because
                                    !!          (epsabs<=0 and
                                    !!           epsrel<max(50*rel.mach.acc.,0.5e-28_wp),
                                    !!          result, abserr, neval, last, rlist(1) ,
                                    !!          `elist(1)` and `iord(1)` are set to zero.
                                    !!          alist(1) and blist(1) are set to a and b
                                    !!          respectively.
        real(wp), intent(out) :: Alist(Limit) !! vector of dimension at least `limit`, the first
                                              !! `last` elements of which are the left
                                              !! end points of the subintervals in the partition
                                              !! of the given integration range `(a,b)`
        real(wp), intent(out) :: Blist(Limit) !! vector of dimension at least `limit`, the first
                                              !! `last` elements of which are the right
                                              !! end points of the subintervals in the partition
                                              !! of the given integration range `(a,b)`
        real(wp), intent(out) :: Elist(Limit) !! vector of dimension at least `limit`, the first
                                              !! `last` elements of which are the moduli of the
                                              !! absolute error estimates on the subintervals
        real(wp), intent(out) :: Rlist(Limit) !! vector of dimension at least `limit`, the first
                                              !! `last` elements of which are the
                                              !! integral approximations on the subintervals
        integer, intent(out) :: Iord(Limit) !! vector of dimension at least `limit`, the first `k`
                                            !! elements of which are pointers to the
                                            !! error estimates over the subintervals,
                                            !! such that `elist(iord(1))`, ...,
                                            !! `elist(iord(k))` form a decreasing sequence,
                                            !! with `k = last` if `last<=(limit/2+2)`, and
                                            !! `k = limit+1-last` otherwise
        integer, intent(out) :: Last !! number of subintervals actually produced in the
                                     !! subdivision process

        real(wp) :: area1, a1, b1, defab1, error1 !! variable for the left subinterval
        real(wp) :: area2, a2, b2, defab2, error2 !! variable for the right subinterval
        real(wp) :: area !! sum of the integrals over the subintervals
        real(wp) :: area12 !! `area1 + area2`
        real(wp) :: erro12 !! `error1 + error2`
        real(wp) :: errsum !! sum of the errors over the subintervals
        real(wp) :: errmax !! `elist(maxerr)`
        real(wp) :: errbnd !! requested accuracy `max(epsabs,epsrel*abs(result))`
        integer :: maxerr !! pointer to the interval with largest error estimate
        real(wp) :: resabs, defabs
        integer :: iroff1, iroff2, k, keyf, nrmax

        ! test on validity of parameters

        Ier = 0
        Neval = 0
        Last = 0
        Result = 0.0_wp
        Abserr = 0.0_wp
        Alist(1) = a
        Blist(1) = b
        Rlist(1) = 0.0_wp
        Elist(1) = 0.0_wp
        Iord(1) = 0
        if (Epsabs <= 0.0_wp .and. Epsrel < max(50.0_wp*epmach, 0.5e-28_wp)) Ier = 6
        if (Ier /= 6) then

            ! first approximation to the integral

            keyf = Key
            if (Key <= 0) keyf = 1
            if (Key >= 7) keyf = 6
            Neval = 0
            select case (keyf)
            case (1); call dqk15(f, a, b, Result, Abserr, defabs, resabs)
            case (2); call dqk21(f, a, b, Result, Abserr, defabs, resabs)
            case (3); call dqk31(f, a, b, Result, Abserr, defabs, resabs)
            case (4); call dqk41(f, a, b, Result, Abserr, defabs, resabs)
            case (5); call dqk51(f, a, b, Result, Abserr, defabs, resabs)
            case (6); call dqk61(f, a, b, Result, Abserr, defabs, resabs)
            end select
            Last = 1
            Rlist(1) = Result
            Elist(1) = Abserr
            Iord(1) = 1

            ! test on accuracy.

            errbnd = max(Epsabs, Epsrel*abs(Result))
            if (Abserr <= 50.0_wp*epmach*defabs .and. Abserr > errbnd) Ier = 2
            if (Limit == 1) Ier = 1

            if (.not. (Ier /= 0 .or. (Abserr <= errbnd .and. Abserr /= resabs) &
                       .or. Abserr == 0.0_wp)) then

                ! initialization
                errmax = Abserr
                maxerr = 1
                area = Result
                errsum = Abserr
                nrmax = 1
                iroff1 = 0
                iroff2 = 0

                ! main do-loop

                do Last = 2, Limit

                    ! bisect the subinterval with the largest error estimate.

                    a1 = Alist(maxerr)
                    b1 = 0.5_wp*(Alist(maxerr) + Blist(maxerr))
                    a2 = b1
                    b2 = Blist(maxerr)
                    select case (keyf)
                    case (1)
                        call dqk15(f, a1, b1, area1, error1, resabs, defab1)
                        call dqk15(f, a2, b2, area2, error2, resabs, defab2)
                    case (2)
                        call dqk21(f, a1, b1, area1, error1, resabs, defab1)
                        call dqk21(f, a2, b2, area2, error2, resabs, defab2)
                    case (3)
                        call dqk31(f, a1, b1, area1, error1, resabs, defab1)
                        call dqk31(f, a2, b2, area2, error2, resabs, defab2)
                    case (4)
                        call dqk41(f, a1, b1, area1, error1, resabs, defab1)
                        call dqk41(f, a2, b2, area2, error2, resabs, defab2)
                    case (5)
                        call dqk51(f, a1, b1, area1, error1, resabs, defab1)
                        call dqk51(f, a2, b2, area2, error2, resabs, defab2)
                    case (6)
                        call dqk61(f, a1, b1, area1, error1, resabs, defab1)
                        call dqk61(f, a2, b2, area2, error2, resabs, defab2)
                    end select

                    ! improve previous approximations to integral
                    ! and error and test for accuracy.

                    Neval = Neval + 1
                    area12 = area1 + area2
                    erro12 = error1 + error2
                    errsum = errsum + erro12 - errmax
                    area = area + area12 - Rlist(maxerr)
                    if (defab1 /= error1 .and. defab2 /= error2) then
                        if (abs(Rlist(maxerr) - area12) <= 0.1e-4_wp*abs(area12) &
                            .and. erro12 >= 0.99_wp*errmax) iroff1 = iroff1 + 1
                        if (Last > 10 .and. erro12 > errmax) iroff2 = iroff2 + 1
                    end if
                    Rlist(maxerr) = area1
                    Rlist(Last) = area2
                    errbnd = max(Epsabs, Epsrel*abs(area))
                    if (errsum > errbnd) then

                        ! test for roundoff error and eventually set error flag.

                        if (iroff1 >= 6 .or. iroff2 >= 20) Ier = 2

                        ! set error flag in the case that the number of subintervals
                        ! equals limit.

                        if (Last == Limit) Ier = 1

                        ! set error flag in the case of bad integrand behaviour
                        ! at a point of the integration range.

                        if (max(abs(a1), abs(b2)) &
                            <= (1.0_wp + 100.0_wp*epmach) &
                            *(abs(a2) + 1000.0_wp*uflow)) Ier = 3
                    end if

                    ! append the newly-created intervals to the list.

                    if (error2 > error1) then
                        Alist(maxerr) = a2
                        Alist(Last) = a1
                        Blist(Last) = b1
                        Rlist(maxerr) = area2
                        Rlist(Last) = area1
                        Elist(maxerr) = error2
                        Elist(Last) = error1
                    else
                        Alist(Last) = a2
                        Blist(maxerr) = b1
                        Blist(Last) = b2
                        Elist(maxerr) = error1
                        Elist(Last) = error2
                    end if

                    ! call subroutine dqpsrt to maintain the descending ordering
                    ! in the list of error estimates and select the subinterval
                    ! with the largest error estimate (to be bisected next).

                    call dqpsrt(Limit, Last, maxerr, errmax, Elist, Iord, nrmax)
                    if (Ier /= 0 .or. errsum <= errbnd) exit  ! jump out of do-loop
                end do

                ! compute final result.

                Result = 0.0_wp
                do k = 1, Last
                    Result = Result + Rlist(k)
                end do
                Abserr = errsum
            end if
            if (keyf /= 1) Neval = (10*keyf + 1)*(2*Neval + 1)
            if (keyf == 1) Neval = 30*Neval + 15
        end if

    end subroutine dqage