dqawo Subroutine

public subroutine dqawo(f, a, b, Omega, Integr, Epsabs, Epsrel, Result, Abserr, Neval, Ier, Leniw, Maxp1, Lenw, Last, Iwork, Work)

1D integration of cos(omega*x)*f(x) or sin(omega*x)*f(x) over a finite interval, adaptive subdivision with extrapolation

the routine calculates an approximation result to a given definite integral i=integral of f(x)*w(x) over (a,b) where w(x) = cos(omega*x) or w(x) = sin(omega*x), 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 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) :: Omega

parameter in the integrand weight function

integer, intent(in) :: Integr

indicates which of the weight functions is used

  • integr = 1 w(x) = cos(omega*x)
  • integr = 2 w(x) = sin(omega*x)

if integr/=1.and.integr/=2, 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 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 (= leniw/2) has been achieved. one can allow more subdivisions by increasing the value of leniw (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 interior 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 due to roundoff in the extrapolation table, 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 (integr/=1 and integr/=2), or leniw<2 or maxp1<1 or lenw<leniw*2+maxp1*25. result, abserr, neval, last are set to zero. except when leniw, maxp1 or lenw are invalid, work(limit*2+1), work(limit*3+1), iwork(1), iwork(limit+1) are set to zero, work(1) is set to a and work(limit+1) to b.
integer, intent(in) :: Leniw

dimensioning parameter for iwork. leniw/2 equals the maximum number of subintervals allowed in the partition of the given integration interval (a,b), leniw>=2. if leniw<2, the routine will end with ier = 6.

integer, intent(in) :: Maxp1

gives an upper bound on the number of chebyshev moments which can be stored, i.e. for the intervals of lengths abs(b-a)*2**(-l), l=0,1, ..., maxp1-2, maxp1>=1 if maxp1<1, the routine will end with ier = 6.

integer, intent(in) :: Lenw

dimensioning parameter for work lenw must be at least leniw*2+maxp1*25. if lenw<(leniw*2+maxp1*25), 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 limit = lenw/2 , and k = last if last<=(limit/2+2), and k = limit+1-last otherwise. furthermore, iwork(limit+1), ..., iwork(limit+last) indicate the subdivision levels of the subintervals, such that iwork(limit+i) = l means that the subinterval numbered i is of length abs(b-a)*2**(1-l).

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.
  • work(limit*4+1), ..., work(limit*4+maxp1*25) provide space for storing the chebyshev moments.

note that limit = lenw/2.


Calls

proc~~dqawo~~CallsGraph proc~dqawo quadpack_generic::dqawo proc~dqawoe quadpack_generic::dqawoe proc~dqawo->proc~dqawoe proc~dqc25f quadpack_generic::dqc25f proc~dqawoe->proc~dqc25f proc~dqcheb quadpack_generic::dqcheb proc~dqc25f->proc~dqcheb proc~dqk15w quadpack_generic::dqk15w proc~dqc25f->proc~dqk15w

Source Code

    subroutine dqawo(f, a, b, Omega, Integr, Epsabs, Epsrel, Result, Abserr, &
                     Neval, Ier, Leniw, Maxp1, Lenw, Last, Iwork, Work)
        implicit none

        procedure(func) :: f !! function subprogram defining the 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) :: Omega !! parameter in the integrand weight function
        integer, intent(in) :: Integr !! indicates which of the weight functions is used
                                      !!
                                      !! * integr = 1  `w(x) = cos(omega*x)`
                                      !! * integr = 2  `w(x) = sin(omega*x)`
                                      !!
                                      !! if `integr/=1.and.integr/=2`, 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 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
                                    !!   `(= leniw/2)` has been achieved. one can
                                    !!   allow more subdivisions by increasing the
                                    !!   value of leniw (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 interior 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
                                    !!   due to roundoff in the extrapolation
                                    !!   table, 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 `(integr/=1 and integr/=2)`,
                                    !!   or `leniw<2` or `maxp1<1` or
                                    !!   `lenw<leniw*2+maxp1*25`.
                                    !!   `result`, `abserr`, `neval`, `last` are set to
                                    !!   zero. except when `leniw`, `maxp1` or `lenw` are
                                    !!   invalid, `work(limit*2+1)`, `work(limit*3+1)`,
                                    !!   `iwork(1)`, `iwork(limit+1)` are set to zero,
                                    !!   `work(1)` is set to `a` and `work(limit+1)` to
                                    !!   `b`.
        integer, intent(in) :: Leniw !! dimensioning parameter for `iwork`.
                                     !! `leniw/2` equals the maximum number of subintervals
                                     !! allowed in the partition of the given integration
                                     !! interval `(a,b)`, `leniw>=2`.
                                     !! if `leniw<2`, the routine will end with ier = 6.
        integer, intent(in) :: Maxp1 !! gives an upper bound on the number of chebyshev
                                     !! moments which can be stored, i.e. for the
                                     !! intervals of lengths `abs(b-a)*2**(-l)`,
                                     !! `l=0,1, ..., maxp1-2, maxp1>=1`
                                     !! if `maxp1<1`, the routine will end with ier = 6.
        integer, intent(in) :: Lenw !! dimensioning parameter for `work`
                                    !! `lenw` must be at least `leniw*2+maxp1*25`.
                                    !! if `lenw<(leniw*2+maxp1*25)`, 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 `limit = lenw/2` , and `k = last`
                                !! if `last<=(limit/2+2)`, and `k = limit+1-last`
                                !! otherwise.
                                !! furthermore, `iwork(limit+1), ..., iwork(limit+last)`
                                !! indicate the subdivision levels of the
                                !! subintervals, such that `iwork(limit+i) = l` means
                                !! that the subinterval numbered `i` is of length
                                !! `abs(b-a)*2**(1-l)`.
        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.
                               !! * `work(limit*4+1), ..., work(limit*4+maxp1*25)`
                               !!   provide space for storing the chebyshev moments.
                               !!
                               !! note that `limit = lenw/2`.

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

        ! check validity of leniw, maxp1 and lenw.
        Ier = 6
        Neval = 0
        Last = 0
        Result = 0.0_wp
        Abserr = 0.0_wp
        if (Leniw >= 2 .and. Maxp1 >= 1 .and. Lenw >= (Leniw*2 + Maxp1*25)) then
            ! prepare call for dqawoe
            limit = Leniw/2
            l1 = limit + 1
            l2 = limit + l1
            l3 = limit + l2
            l4 = limit + l3
            call dqawoe(f, a, b, Omega, Integr, Epsabs, Epsrel, limit, 1, Maxp1, &
                        Result, Abserr, Neval, Ier, Last, Work(1), Work(l1), &
                        Work(l2), Work(l3), Iwork(1), Iwork(l1), momcom, &
                        Work(l4))
            ! call error handler if necessary
            lvl = 0
        end if
        if (Ier == 6) lvl = 1
        if (Ier /= 0) call xerror('abnormal return from dqawo', Ier, lvl)

    end subroutine dqawo