same as dqawc but provides more information and control
the routine calculates an approximation result to a
cauchy principal value i = integral of f*w
over (a,b)
(w(x) = 1/(x-c), (c/=a, c/=b)
, hopefully satisfying
following claim for accuracy
abs(i-result)<=max(epsabs,epsrel*abs(i))
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(func) | :: | f |
function subprogram defining the integrand function |
|||
real(kind=wp), | intent(in) | :: | a |
lower limit of integration |
||
real(kind=wp), | intent(in) | :: | b |
upper limit of integration |
||
real(kind=wp) | :: | c | ||||
real(kind=wp), | intent(in) | :: | Epsabs |
absolute accuracy requested |
||
real(kind=wp), | intent(in) | :: | Epsrel |
relative accuracy requested
if |
||
integer, | intent(in) | :: | Limit |
gives an upper bound on the number of subintervals
in the partition of |
||
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 |
||
integer, | intent(out) | :: | Neval |
number of integrand evaluations |
||
integer, | intent(out) | :: | Ier |
error messages:
|
||
real(kind=wp), | intent(out) | :: | Alist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Blist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Rlist(Limit) |
vector of dimension at least |
||
real(kind=wp), | intent(out) | :: | Elist(Limit) |
vector of dimension |
||
integer, | intent(out) | :: | Iord(Limit) |
vector of dimension at least |
||
integer, | intent(out) | :: | Last |
number of subintervals actually produced in the subdivision process |
subroutine dqawce(f, a, b, c, 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 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>=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 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. however, if this yields no !! improvement it is advised to analyze the !! the integrand, in order to determine the !! 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 !! appropriate integrators 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 interior points of !! the integration interval. !! * ier = 6 the input is invalid, because !! `c = a` or `c = b` or !! `(epsabs<=0 and epsrel<max(50*rel.mach.acc.,0.5e-28))` !! or `limit<1`. !! `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 `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, 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) :: aa, bb, c integer :: iroff1, iroff2, k, krule, nev, nrmax real(wp) :: area1, a1, b1, error1 !! variable for the left subinterval real(wp) :: area2, a2, b2, error2 !! variable for the right subinterval real(wp) :: area12 !! `area1 + area2` real(wp) :: erro12 !! `error1 + error2` real(wp) :: errmax !! elist(maxerr) real(wp) :: area !! sum of the integrals over the subintervals real(wp) :: errsum !! sum of the errors over the subintervals real(wp) :: errbnd !! requested accuracy `max(epsabs,epsrel*abs(result))` integer :: maxerr !! pointer to the interval with largest error estimate ! test on validity of parameters Ier = 6 Neval = 0 Last = 0 Alist(1) = a Blist(1) = b Rlist(1) = 0.0_wp Elist(1) = 0.0_wp Iord(1) = 0 Result = 0.0_wp Abserr = 0.0_wp if (.not. (c == a .or. c == b .or. (Epsabs <= 0.0_wp .and. Epsrel < max & (50.0_wp*epmach, 0.5e-28_wp)))) then ! first approximation to the integral aa = a bb = b if (a > b) then aa = b bb = a end if Ier = 0 krule = 1 call dqc25c(f, aa, bb, c, Result, Abserr, krule, Neval) Last = 1 Rlist(1) = Result Elist(1) = Abserr Iord(1) = 1 Alist(1) = a Blist(1) = b ! test on accuracy errbnd = max(Epsabs, Epsrel*abs(Result)) if (Limit == 1) Ier = 1 if (Abserr >= min(0.01_wp*abs(Result), errbnd) .and. Ier /= 1) then ! initialization Alist(1) = aa Blist(1) = bb Rlist(1) = Result 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 nrmax-th largest ! error estimate. a1 = Alist(maxerr) b1 = 0.5_wp*(Alist(maxerr) + Blist(maxerr)) b2 = Blist(maxerr) if (c <= b1 .and. c > a1) b1 = 0.5_wp*(c + b2) if (c > b1 .and. c < b2) b1 = 0.5_wp*(a1 + c) a2 = b1 krule = 2 call dqc25c(f, a1, b1, c, area1, error1, krule, nev) Neval = Neval + nev call dqc25c(f, a2, b2, c, area2, error2, krule, nev) Neval = Neval + nev ! improve previous approximations to integral ! and error and test for accuracy. area12 = area1 + area2 erro12 = error1 + error2 errsum = errsum + erro12 - errmax area = area + area12 - Rlist(maxerr) if (abs(Rlist(maxerr) - area12) < 0.1e-4_wp*abs(area12) & .and. erro12 >= 0.99_wp*errmax .and. krule == 0) & iroff1 = iroff1 + 1 if (Last > 10 .and. erro12 > errmax .and. krule == 0) & iroff2 = iroff2 + 1 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 .and. iroff2 > 20) Ier = 2 ! set error flag in the case that number of interval ! bisections exceeds 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 nrmax-th 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 if (aa == b) Result = -Result end if end subroutine dqawce