dqaws Subroutine

public subroutine dqaws(f, a, b, alfa, beta, integr, epsabs, epsrel, result, abserr, neval, ier, limit, lenw, last, iwork, work)

1D integration of functions with powers and or logs over a finite interval

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 integrand function, alfa>(-1) if alfa<=(-1), the routine will end with ier = 6.

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

parameter in the integrand 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.

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 (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 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 or integr<1 or integr>4 or (epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28)) or limit<2 or lenw<limit*4. result, abserr, neval, last are set to zero. except when lenw or limit 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>=2. if limit<2, 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, which determines the significant number of elements actually in the work arrays.

integer :: iwork(limit)

vector of dimension 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)

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~~dqaws~~CallsGraph proc~dqaws quadpack_generic::dqaws proc~dqawse quadpack_generic::dqawse proc~dqaws->proc~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

Source Code

    subroutine dqaws(f, a, b, alfa, beta, integr, 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, b>a
                                  !! if b<=a, the routine will end with ier = 6.
        real(wp), intent(in) :: alfa !! parameter in the integrand function, `alfa>(-1)`
                                     !! if `alfa<=(-1)`, the routine will end with
                                     !! ier = 6.
        real(wp), intent(in) :: beta !! parameter in the integrand 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.
        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 (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
                                    !!   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
                                    !!   or `integr<1` or `integr>4` or
                                    !!   `(epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28))`
                                    !!   or `limit<2` or `lenw<limit*4`.
                                    !!   `result`, `abserr`, `neval`, `last` are set to
                                    !!   zero. except when `lenw` or `limit` 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>=2`.
                                    !! if `limit<2`, 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,
                                     !! which determines the significant number of
                                     !! elements actually in the work arrays.
        integer :: iwork(limit) !! vector of dimension 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) !! 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 >= 2 .and. lenw >= limit*4) then

            ! prepare call for dqawse.

            l1 = limit + 1
            l2 = limit + l1
            l3 = limit + l2

            call dqawse(f, a, b, alfa, beta, integr, 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 dqaws', ier, lvl)

    end subroutine dqaws