dqawse Subroutine

public subroutine dqawse(f, a, b, Alfa, Beta, Integr, Epsabs, Epsrel, Limit, Result, Abserr, Neval, Ier, Alist, Blist, Rlist, Elist, Iord, Last)

same as dqaws but provides more information and control

the routine calculates an approximation result to a given definite integral i = integral of f*w over (a,b), (where w shows a singular behaviour at the end points, see parameter integr). 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, b>a. if b<=a, the routine will end with ier = 6.

real(kind=wp), intent(in) :: Alfa

parameter in the weight function, alfa>(-1) if alfa<=(-1), the routine will end with ier = 6.

real(kind=wp), intent(in) :: Beta

parameter in the weight function, beta>(-1) if beta<=(-1), the routine will end with ier = 6.

integer, intent(in) :: Integr

indicates which weight function is to be used:

  • = 1 (x-a)**alfa*(b-x)**beta
  • = 2 (x-a)**alfa*(b-x)**beta*log(x-a)
  • = 3 (x-a)**alfa*(b-x)**beta*log(b-x)
  • = 4 (x-a)**alfa*(b-x)**beta*log(x-a)*log(b-x)

if integr<1 or integr>4, the routine will end with ier = 6.

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) :: Limit

gives an upper bound on the number of subintervals in the partition of (a,b), limit>=2 if limit<2, 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 the 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. however, if this yields no improvement, it is advised to analyze the integrand in order to determine the integration difficulties which prevent the requested tolerance from being achieved. in case of a jump discontinuity or a local singularity of algebraico-logarithmic type at one or more interior points of the integration range, one should proceed by splitting up the interval at these points and calling the integrator 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 b<=a or alfa<=(-1) or beta<=(-1), or integr<1 or integr>4, or epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28), or limit<2. result, abserr, neval, rlist(1), elist(1), iord(1) and last 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 of which are pointers to the error estimates over the subintervals, so that elist(iord(1)), ..., elist(iord(k)) with k = last if last<=(limit/2+2), and k = limit+1-last otherwise form a decreasing sequence

integer, intent(out) :: Last

number of subintervals actually produced in the subdivision process


Calls

proc~~dqawse~~CallsGraph proc~dqawse quadpack_generic::dqawse proc~dqc25s quadpack_generic::dqc25s proc~dqawse->proc~dqc25s proc~dqmomo quadpack_generic::dqmomo proc~dqawse->proc~dqmomo proc~dqcheb quadpack_generic::dqcheb proc~dqc25s->proc~dqcheb proc~dqk15w quadpack_generic::dqk15w proc~dqc25s->proc~dqk15w

Called by

proc~~dqawse~~CalledByGraph proc~dqawse quadpack_generic::dqawse proc~dqaws quadpack_generic::dqaws proc~dqaws->proc~dqawse

Source Code

    subroutine dqawse(f, a, b, Alfa, Beta, Integr, Epsabs, Epsrel, 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, `b>a`.
                                  !! if `b<=a`, the routine will end with ier = 6.
        real(wp), intent(in) :: Alfa !! parameter in the weight function, `alfa>(-1)`
                                     !! if `alfa<=(-1)`, the routine will end with
                                     !! ier = 6.
        real(wp), intent(in) :: Beta !! parameter in the weight function, `beta>(-1)`
                                     !! if `beta<=(-1)`, the routine will end with
                                     !! ier = 6.
        integer, intent(in) :: Integr !! indicates which weight function is to be used:
                                      !!
                                      !! * = 1  `(x-a)**alfa*(b-x)**beta`
                                      !! * = 2  `(x-a)**alfa*(b-x)**beta*log(x-a)`
                                      !! * = 3  `(x-a)**alfa*(b-x)**beta*log(b-x)`
                                      !! * = 4  `(x-a)**alfa*(b-x)**beta*log(x-a)*log(b-x)`
                                      !!
                                      !! if `integr<1` or `integr>4`, the routine
                                      !! will end with ier = 6.
        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) :: Limit !! gives an upper bound on the number of subintervals
                                     !! in the partition of `(a,b)`, `limit>=2`
                                     !! if `limit<2`, 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 the 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. however, if this yields no
                                    !!   improvement, it is advised to analyze the
                                    !!   integrand in order to determine the
                                    !!   integration difficulties which prevent the
                                    !!   requested tolerance from being achieved.
                                    !!   in case of a jump discontinuity or a local
                                    !!   singularity of algebraico-logarithmic type
                                    !!   at one or more interior points of the
                                    !!   integration range, one should proceed by
                                    !!   splitting up the interval at these
                                    !!   points and calling the integrator 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
                                    !!   `b<=a` or `alfa<=(-1)` or `beta<=(-1)`, or
                                    !!   `integr<1` or `integr>4`, or
                                    !!   `epsabs<=0` and
                                    !!   `epsrel<max(50*rel.mach.acc.,0.5e-28)`,
                                    !!   or `limit<2`.
                                    !!   `result`, `abserr`, `neval`, `rlist(1)`, `elist(1)`,
                                    !!   `iord(1)` and `last` 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) :: Rlist(Limit) !! vector of dimension at least `limit`,the first
                                              !! `last` elements of which are the integral
                                              !! approximations on the subintervals
        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
        integer, intent(out) :: Iord(Limit) !! vector of dimension at least `limit`, the first `k`
                                            !! of which are pointers to the error
                                            !! estimates over the subintervals, so that
                                            !! `elist(iord(1)), ..., elist(iord(k))` with `k = last`
                                            !! if `last<=(limit/2+2)`, and `k = limit+1-last`
                                            !! otherwise form a decreasing sequence
        integer, intent(out) :: Last !! number of subintervals actually produced in
                                     !! the subdivision process

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

        ! test on validity of parameters

        Ier = 6
        Neval = 0
        Last = 0
        Rlist(1) = 0.0_wp
        Elist(1) = 0.0_wp
        Iord(1) = 0
        Result = 0.0_wp
        Abserr = 0.0_wp
        if (.not. (b <= a .or. (Epsabs == 0.0_wp .and. Epsrel < max(50.0_wp* &
            epmach, 0.5e-28_wp)) .or. Alfa <= (-1.0_wp) .or. Beta <= (-1.0_wp) &
            .or. Integr < 1 .or. Integr > 4 .or. Limit < 2)) then
            Ier = 0

            ! compute the modified chebyshev moments.

            call dqmomo(Alfa, Beta, ri, rj, rg, rh, Integr)

            ! integrate over the intervals (a,(a+b)/2) and ((a+b)/2,b).

            centre = 0.5_wp*(b + a)
            call dqc25s(f, a, b, a, centre, Alfa, Beta, ri, rj, rg, rh, area1, error1, &
                        resas1, Integr, nev)
            Neval = nev
            call dqc25s(f, a, b, centre, b, Alfa, Beta, ri, rj, rg, rh, area2, error2, &
                        resas2, Integr, nev)
            Last = 2
            Neval = Neval + nev
            Result = area1 + area2
            Abserr = error1 + error2

            ! test on accuracy.

            errbnd = max(Epsabs, Epsrel*abs(Result))

            ! initialization

            if (error2 > error1) then
                Alist(1) = centre
                Alist(2) = a
                Blist(1) = b
                Blist(2) = centre
                Rlist(1) = area2
                Rlist(2) = area1
                Elist(1) = error2
                Elist(2) = error1
            else
                Alist(1) = a
                Alist(2) = centre
                Blist(1) = centre
                Blist(2) = b
                Rlist(1) = area1
                Rlist(2) = area2
                Elist(1) = error1
                Elist(2) = error2
            end if
            Iord(1) = 1
            Iord(2) = 2
            if (Limit == 2) Ier = 1
            if (Abserr > errbnd .and. Ier /= 1) then
                errmax = Elist(1)
                maxerr = 1
                nrmax = 1
                area = Result
                errsum = Abserr
                iroff1 = 0
                iroff2 = 0

                ! main do-loop

                do Last = 3, Limit

                    ! bisect the subinterval with largest error estimate.

                    a1 = Alist(maxerr)
                    b1 = 0.5_wp*(Alist(maxerr) + Blist(maxerr))
                    a2 = b1
                    b2 = Blist(maxerr)

                    call dqc25s(f, a, b, a1, b1, Alfa, Beta, ri, rj, rg, rh, area1, &
                                error1, resas1, Integr, nev)
                    Neval = Neval + nev
                    call dqc25s(f, a, b, a2, b2, Alfa, Beta, ri, rj, rg, rh, area2, &
                                error2, resas2, Integr, nev)
                    Neval = Neval + nev

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

                    area12 = area1 + area2
                    erro12 = error1 + error2
                    errsum = errsum + erro12 - errmax
                    area = area + area12 - Rlist(maxerr)
                    if (a /= a1 .and. b /= b2) then
                        if (resas1 /= error1 .and. resas2 /= error2) then
                            ! test for roundoff error.
                            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
                    end if
                    Rlist(maxerr) = area1
                    Rlist(Last) = area2

                    ! test on accuracy.

                    errbnd = max(Epsabs, Epsrel*abs(area))
                    if (errsum > errbnd) then

                        ! set error flag in the case that the number of interval
                        ! bisections exceeds limit.

                        if (Last == Limit) Ier = 1

                        ! set error flag in the case of roundoff error.

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

                        ! set error flag in the case of bad integrand behaviour
                        ! at interior points of 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 largest error estimate (to be bisected next).

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

                ! compute final result.
                Result = 0.0_wp
                do k = 1, Last
                    Result = Result + Rlist(k)
                end do
                Abserr = errsum
            end if
        end if

    end subroutine dqawse