1      subroutine dqng(f,a,b,epsabs,epsrel,result,abserr,neval,ier)
2c***begin prologue  dqng
3c***date written   800101   (yymmdd)
4c***revision date  810101   (yymmdd)
5c***category no.  h2a1a1
6c***keywords  automatic integrator, smooth integrand,
7c             non-adaptive, gauss-kronrod(patterson)
8c***author  piessens,robert,appl. math. & progr. div. - k.u.leuven
9c           de doncker,elise,appl math & progr. div. - k.u.leuven
10c           kahaner,david,nbs - modified (2/82)
11c***purpose  the routine calculates an approximation result to a
12c            given definite integral i = integral of f over (a,b),
13c            hopefully satisfying following claim for accuracy
14c            abs(i-result).le.max(epsabs,epsrel*abs(i)).
15c***description
16c
17c non-adaptive integration
18c standard fortran subroutine
19c double precision version
20c
21c           f      - double precision
22c                    function subprogram defining the integrand function
23c                    f(x). the actual name for f needs to be declared
24c                    e x t e r n a l in the driver program.
25c
26c           a      - double precision
27c                    lower limit of integration
28c
29c           b      - double precision
30c                    upper limit of integration
31c
32c           epsabs - double precision
33c                    absolute accuracy requested
34c           epsrel - double precision
35c                    relative accuracy requested
36c                    if  epsabs.le.0
37c                    and epsrel.lt.max(50*rel.mach.acc.,0.5d-28),
38c                    the routine will end with ier = 6.
39c
40c         on return
41c           result - double precision
42c                    approximation to the integral i
43c                    result is obtained by applying the 21-point
44c                    gauss-kronrod rule (res21) obtained by optimal
45c                    addition of abscissae to the 10-point gauss rule
46c                    (res10), or by applying the 43-point rule (res43)
47c                    obtained by optimal addition of abscissae to the
48c                    21-point gauss-kronrod rule, or by applying the
49c                    87-point rule (res87) obtained by optimal addition
50c                    of abscissae to the 43-point rule.
51c
52c           abserr - double precision
53c                    estimate of the modulus of the absolute error,
54c                    which should equal or exceed abs(i-result)
55c
56c           neval  - integer
57c                    number of integrand evaluations
58c
59c           ier    - ier = 0 normal and reliable termination of the
60c                            routine. it is assumed that the requested
61c                            accuracy has been achieved.
62c                    ier.gt.0 abnormal termination of the routine. it is
63c                            assumed that the requested accuracy has
64c                            not been achieved.
65c           error messages
66c                    ier = 1 the maximum number of steps has been
67c                            executed. the integral is probably too
68c                            difficult to be calculated by dqng.
69c                        = 6 the input is invalid, because
70c                            epsabs.le.0 and
71c                            epsrel.lt.max(50*rel.mach.acc.,0.5d-28).
72c                            result, abserr and neval are set to zero.
73c
74c***references  (none)
75c***routines called  d1mach,xerror
76c***end prologue  dqng
77c
78      double precision a,absc,abserr,b,centr,dabs,dhlgth,dmax1,dmin1,
79     *  d1mach,epmach,epsabs,epsrel,f,fcentr,fval,fval1,fval2,fv1,fv2,
80     *  fv3,fv4,hlgth,result,res10,res21,res43,res87,resabs,resasc,
81     *  reskh,savfun,uflow,w10,w21a,w21b,w43a,w43b,w87a,w87b,x1,x2,x3,x4
82      integer ier,ipx,k,l,neval
83      external f
84c
85      dimension fv1(5),fv2(5),fv3(5),fv4(5),x1(5),x2(5),x3(11),x4(22),
86     *  w10(5),w21a(5),w21b(6),w43a(10),w43b(12),w87a(21),w87b(23),
87     *  savfun(21)
88c
89c           the following data statements contain the
90c           abscissae and weights of the integration rules used.
91c
92c           x1      abscissae common to the 10-, 21-, 43- and 87-
93c                   point rule
94c           x2      abscissae common to the 21-, 43- and 87-point rule
95c           x3      abscissae common to the 43- and 87-point rule
96c           x4      abscissae of the 87-point rule
97c           w10     weights of the 10-point formula
98c           w21a    weights of the 21-point formula for abscissae x1
99c           w21b    weights of the 21-point formula for abscissae x2
100c           w43a    weights of the 43-point formula for abscissae x1, x3
101c           w43b    weights of the 43-point formula for abscissae x3
102c           w87a    weights of the 87-point formula for abscissae x1,
103c                   x2, x3
104c           w87b    weights of the 87-point formula for abscissae x4
105c
106c
107c gauss-kronrod-patterson quadrature coefficients for use in
108c quadpack routine qng.  these coefficients were calculated with
109c 101 decimal digit arithmetic by l. w. fullerton, bell labs, nov 1981.
110c
111      data x1    (  1) / 0.9739065285 1717172007 7964012084 452 d0 /
112      data x1    (  2) / 0.8650633666 8898451073 2096688423 493 d0 /
113      data x1    (  3) / 0.6794095682 9902440623 4327365114 874 d0 /
114      data x1    (  4) / 0.4333953941 2924719079 9265943165 784 d0 /
115      data x1    (  5) / 0.1488743389 8163121088 4826001129 720 d0 /
116      data w10   (  1) / 0.0666713443 0868813759 3568809893 332 d0 /
117      data w10   (  2) / 0.1494513491 5058059314 5776339657 697 d0 /
118      data w10   (  3) / 0.2190863625 1598204399 5534934228 163 d0 /
119      data w10   (  4) / 0.2692667193 0999635509 1226921569 469 d0 /
120      data w10   (  5) / 0.2955242247 1475287017 3892994651 338 d0 /
121c
122      data x2    (  1) / 0.9956571630 2580808073 5527280689 003 d0 /
123      data x2    (  2) / 0.9301574913 5570822600 1207180059 508 d0 /
124      data x2    (  3) / 0.7808177265 8641689706 3717578345 042 d0 /
125      data x2    (  4) / 0.5627571346 6860468333 9000099272 694 d0 /
126      data x2    (  5) / 0.2943928627 0146019813 1126603103 866 d0 /
127      data w21a  (  1) / 0.0325581623 0796472747 8818972459 390 d0 /
128      data w21a  (  2) / 0.0750396748 1091995276 7043140916 190 d0 /
129      data w21a  (  3) / 0.1093871588 0229764189 9210590325 805 d0 /
130      data w21a  (  4) / 0.1347092173 1147332592 8054001771 707 d0 /
131      data w21a  (  5) / 0.1477391049 0133849137 4841515972 068 d0 /
132      data w21b  (  1) / 0.0116946388 6737187427 8064396062 192 d0 /
133      data w21b  (  2) / 0.0547558965 7435199603 1381300244 580 d0 /
134      data w21b  (  3) / 0.0931254545 8369760553 5065465083 366 d0 /
135      data w21b  (  4) / 0.1234919762 6206585107 7958109831 074 d0 /
136      data w21b  (  5) / 0.1427759385 7706008079 7094273138 717 d0 /
137      data w21b  (  6) / 0.1494455540 0291690566 4936468389 821 d0 /
138c
139      data x3    (  1) / 0.9993333609 0193208139 4099323919 911 d0 /
140      data x3    (  2) / 0.9874334029 0808886979 5961478381 209 d0 /
141      data x3    (  3) / 0.9548079348 1426629925 7919200290 473 d0 /
142      data x3    (  4) / 0.9001486957 4832829362 5099494069 092 d0 /
143      data x3    (  5) / 0.8251983149 8311415084 7066732588 520 d0 /
144      data x3    (  6) / 0.7321483889 8930498261 2354848755 461 d0 /
145      data x3    (  7) / 0.6228479705 3772523864 1159120344 323 d0 /
146      data x3    (  8) / 0.4994795740 7105649995 2214885499 755 d0 /
147      data x3    (  9) / 0.3649016613 4658076804 3989548502 644 d0 /
148      data x3    ( 10) / 0.2222549197 7660129649 8260928066 212 d0 /
149      data x3    ( 11) / 0.0746506174 6138332204 3914435796 506 d0 /
150      data w43a  (  1) / 0.0162967342 8966656492 4281974617 663 d0 /
151      data w43a  (  2) / 0.0375228761 2086950146 1613795898 115 d0 /
152      data w43a  (  3) / 0.0546949020 5825544214 7212685465 005 d0 /
153      data w43a  (  4) / 0.0673554146 0947808607 5553166302 174 d0 /
154      data w43a  (  5) / 0.0738701996 3239395343 2140695251 367 d0 /
155      data w43a  (  6) / 0.0057685560 5976979618 4184327908 655 d0 /
156      data w43a  (  7) / 0.0273718905 9324884208 1276069289 151 d0 /
157      data w43a  (  8) / 0.0465608269 1042883074 3339154433 824 d0 /
158      data w43a  (  9) / 0.0617449952 0144256449 6240336030 883 d0 /
159      data w43a  ( 10) / 0.0713872672 6869339776 8559114425 516 d0 /
160      data w43b  (  1) / 0.0018444776 4021241410 0389106552 965 d0 /
161      data w43b  (  2) / 0.0107986895 8589165174 0465406741 293 d0 /
162      data w43b  (  3) / 0.0218953638 6779542810 2523123075 149 d0 /
163      data w43b  (  4) / 0.0325974639 7534568944 3882222526 137 d0 /
164      data w43b  (  5) / 0.0421631379 3519181184 7627924327 955 d0 /
165      data w43b  (  6) / 0.0507419396 0018457778 0189020092 084 d0 /
166      data w43b  (  7) / 0.0583793955 4261924837 5475369330 206 d0 /
167      data w43b  (  8) / 0.0647464049 5144588554 4689259517 511 d0 /
168      data w43b  (  9) / 0.0695661979 1235648452 8633315038 405 d0 /
169      data w43b  ( 10) / 0.0728244414 7183320815 0939535192 842 d0 /
170      data w43b  ( 11) / 0.0745077510 1417511827 3571813842 889 d0 /
171      data w43b  ( 12) / 0.0747221475 1740300559 4425168280 423 d0 /
172c
173      data x4    (  1) / 0.9999029772 6272923449 0529830591 582 d0 /
174      data x4    (  2) / 0.9979898959 8667874542 7496322365 960 d0 /
175      data x4    (  3) / 0.9921754978 6068722280 8523352251 425 d0 /
176      data x4    (  4) / 0.9813581635 7271277357 1916941623 894 d0 /
177      data x4    (  5) / 0.9650576238 5838461912 8284110607 926 d0 /
178      data x4    (  6) / 0.9431676131 3367059681 6416634507 426 d0 /
179      data x4    (  7) / 0.9158064146 8550720959 1826430720 050 d0 /
180      data x4    (  8) / 0.8832216577 7131650137 2117548744 163 d0 /
181      data x4    (  9) / 0.8457107484 6241566660 5902011504 855 d0 /
182      data x4    ( 10) / 0.8035576580 3523098278 8739474980 964 d0 /
183      data x4    ( 11) / 0.7570057306 8549555832 8942793432 020 d0 /
184      data x4    ( 12) / 0.7062732097 8732181982 4094274740 840 d0 /
185      data x4    ( 13) / 0.6515894665 0117792253 4422205016 736 d0 /
186      data x4    ( 14) / 0.5932233740 5796108887 5273770349 144 d0 /
187      data x4    ( 15) / 0.5314936059 7083193228 5268948562 671 d0 /
188      data x4    ( 16) / 0.4667636230 4202284487 1966781659 270 d0 /
189      data x4    ( 17) / 0.3994248478 5921880473 2101665817 923 d0 /
190      data x4    ( 18) / 0.3298748771 0618828826 5053371824 597 d0 /
191      data x4    ( 19) / 0.2585035592 0216155180 2280975429 025 d0 /
192      data x4    ( 20) / 0.1856953965 6834665201 5917141167 606 d0 /
193      data x4    ( 21) / 0.1118422131 7990746817 2398359241 362 d0 /
194      data x4    ( 22) / 0.0373521233 9461987081 4998165437 704 d0 /
195      data w87a  (  1) / 0.0081483773 8414917290 0002878448 190 d0 /
196      data w87a  (  2) / 0.0187614382 0156282224 3935059003 794 d0 /
197      data w87a  (  3) / 0.0273474510 5005228616 1582829741 283 d0 /
198      data w87a  (  4) / 0.0336777073 1163793004 6581056957 588 d0 /
199      data w87a  (  5) / 0.0369350998 2042790761 4589586742 499 d0 /
200      data w87a  (  6) / 0.0028848724 3021153050 1334156248 695 d0 /
201      data w87a  (  7) / 0.0136859460 2271270188 8950035273 128 d0 /
202      data w87a  (  8) / 0.0232804135 0288831112 3409291030 404 d0 /
203      data w87a  (  9) / 0.0308724976 1171335867 5466394126 442 d0 /
204      data w87a  ( 10) / 0.0356936336 3941877071 9351355457 044 d0 /
205      data w87a  ( 11) / 0.0009152833 4520224136 0843392549 948 d0 /
206      data w87a  ( 12) / 0.0053992802 1930047136 7738743391 053 d0 /
207      data w87a  ( 13) / 0.0109476796 0111893113 4327826856 808 d0 /
208      data w87a  ( 14) / 0.0162987316 9678733526 2665703223 280 d0 /
209      data w87a  ( 15) / 0.0210815688 8920383511 2433060188 190 d0 /
210      data w87a  ( 16) / 0.0253709697 6925382724 3467999831 710 d0 /
211      data w87a  ( 17) / 0.0291896977 5647575250 1446154084 920 d0 /
212      data w87a  ( 18) / 0.0323732024 6720278968 5788194889 595 d0 /
213      data w87a  ( 19) / 0.0347830989 5036514275 0781997949 596 d0 /
214      data w87a  ( 20) / 0.0364122207 3135178756 2801163687 577 d0 /
215      data w87a  ( 21) / 0.0372538755 0304770853 9592001191 226 d0 /
216      data w87b  (  1) / 0.0002741455 6376207235 0016527092 881 d0 /
217      data w87b  (  2) / 0.0018071241 5505794294 8341311753 254 d0 /
218      data w87b  (  3) / 0.0040968692 8275916486 4458070683 480 d0 /
219      data w87b  (  4) / 0.0067582900 5184737869 9816577897 424 d0 /
220      data w87b  (  5) / 0.0095499576 7220164653 6053581325 377 d0 /
221      data w87b  (  6) / 0.0123294476 5224485369 4626639963 780 d0 /
222      data w87b  (  7) / 0.0150104473 4638895237 6697286041 943 d0 /
223      data w87b  (  8) / 0.0175489679 8624319109 9665352925 900 d0 /
224      data w87b  (  9) / 0.0199380377 8644088820 2278192730 714 d0 /
225      data w87b  ( 10) / 0.0221949359 6101228679 6332102959 499 d0 /
226      data w87b  ( 11) / 0.0243391471 2600080547 0360647041 454 d0 /
227      data w87b  ( 12) / 0.0263745054 1483920724 1503786552 615 d0 /
228      data w87b  ( 13) / 0.0282869107 8877120065 9968002987 960 d0 /
229      data w87b  ( 14) / 0.0300525811 2809269532 2521110347 341 d0 /
230      data w87b  ( 15) / 0.0316467513 7143992940 4586051078 883 d0 /
231      data w87b  ( 16) / 0.0330504134 1997850329 0785944862 689 d0 /
232      data w87b  ( 17) / 0.0342550997 0422606178 7082821046 821 d0 /
233      data w87b  ( 18) / 0.0352624126 6015668103 3782717998 428 d0 /
234      data w87b  ( 19) / 0.0360769896 2288870118 5500318003 895 d0 /
235      data w87b  ( 20) / 0.0366986044 9845609449 8018047441 094 d0 /
236      data w87b  ( 21) / 0.0371205492 6983257611 4119958413 599 d0 /
237      data w87b  ( 22) / 0.0373342287 5193504032 1235449094 698 d0 /
238      data w87b  ( 23) / 0.0373610737 6267902341 0321241766 599 d0 /
239c
240c           list of major variables
241c           -----------------------
242c
243c           centr  - mid point of the integration interval
244c           hlgth  - half-length of the integration interval
245c           fcentr - function value at mid point
246c           absc   - abscissa
247c           fval   - function value
248c           savfun - array of function values which have already been
249c                    computed
250c           res10  - 10-point gauss result
251c           res21  - 21-point kronrod result
252c           res43  - 43-point result
253c           res87  - 87-point result
254c           resabs - approximation to the integral of abs(f)
255c           resasc - approximation to the integral of abs(f-i/(b-a))
256c
257c           machine dependent constants
258c           ---------------------------
259c
260c           epmach is the largest relative spacing.
261c           uflow is the smallest positive magnitude.
262c
263c***first executable statement  dqng
264      epmach = d1mach(4)
265      uflow = d1mach(1)
266c
267c           test on validity of parameters
268c           ------------------------------
269c
270      result = 0.0d+00
271      abserr = 0.0d+00
272      neval = 0
273      ier = 6
274      if(epsabs.le.0.0d+00.and.epsrel.lt.dmax1(0.5d+02*epmach,0.5d-28))
275     *  go to 80
276      hlgth = 0.5d+00*(b-a)
277      dhlgth = dabs(hlgth)
278      centr = 0.5d+00*(b+a)
279      fcentr = f(centr)
280      neval = 21
281      ier = 1
282c
283c          compute the integral using the 10- and 21-point formula.
284c
285      do 70 l = 1,3
286      go to (5,25,45),l
287    5 res10 = 0.0d+00
288      res21 = w21b(6)*fcentr
289      resabs = w21b(6)*dabs(fcentr)
290      do 10 k=1,5
291        absc = hlgth*x1(k)
292        fval1 = f(centr+absc)
293        fval2 = f(centr-absc)
294        fval = fval1+fval2
295        res10 = res10+w10(k)*fval
296        res21 = res21+w21a(k)*fval
297        resabs = resabs+w21a(k)*(dabs(fval1)+dabs(fval2))
298        savfun(k) = fval
299        fv1(k) = fval1
300        fv2(k) = fval2
301   10 continue
302      ipx = 5
303      do 15 k=1,5
304        ipx = ipx+1
305        absc = hlgth*x2(k)
306        fval1 = f(centr+absc)
307        fval2 = f(centr-absc)
308        fval = fval1+fval2
309        res21 = res21+w21b(k)*fval
310        resabs = resabs+w21b(k)*(dabs(fval1)+dabs(fval2))
311        savfun(ipx) = fval
312        fv3(k) = fval1
313        fv4(k) = fval2
314   15 continue
315c
316c          test for convergence.
317c
318      result = res21*hlgth
319      resabs = resabs*dhlgth
320      reskh = 0.5d+00*res21
321      resasc = w21b(6)*dabs(fcentr-reskh)
322      do 20 k = 1,5
323        resasc = resasc+w21a(k)*(dabs(fv1(k)-reskh)+dabs(fv2(k)-reskh))
324     *                  +w21b(k)*(dabs(fv3(k)-reskh)+dabs(fv4(k)-reskh))
325   20 continue
326      abserr = dabs((res21-res10)*hlgth)
327      resasc = resasc*dhlgth
328      go to 65
329c
330c          compute the integral using the 43-point formula.
331c
332   25 res43 = w43b(12)*fcentr
333      neval = 43
334      do 30 k=1,10
335        res43 = res43+savfun(k)*w43a(k)
336   30 continue
337      do 40 k=1,11
338        ipx = ipx+1
339        absc = hlgth*x3(k)
340        fval = f(absc+centr)+f(centr-absc)
341        res43 = res43+fval*w43b(k)
342        savfun(ipx) = fval
343   40 continue
344c
345c          test for convergence.
346c
347      result = res43*hlgth
348      abserr = dabs((res43-res21)*hlgth)
349      go to 65
350c
351c          compute the integral using the 87-point formula.
352c
353   45 res87 = w87b(23)*fcentr
354      neval = 87
355      do 50 k=1,21
356        res87 = res87+savfun(k)*w87a(k)
357   50 continue
358      do 60 k=1,22
359        absc = hlgth*x4(k)
360        res87 = res87+w87b(k)*(f(absc+centr)+f(centr-absc))
361   60 continue
362      result = res87*hlgth
363      abserr = dabs((res87-res43)*hlgth)
364   65 if(resasc.ne.0.0d+00.and.abserr.ne.0.0d+00)
365     *  abserr = resasc*dmin1(0.1d+01,(0.2d+03*abserr/resasc)**1.5d+00)
366      if (resabs.gt.uflow/(0.5d+02*epmach)) abserr = dmax1
367     *  ((epmach*0.5d+02)*resabs,abserr)
368      if (abserr.le.dmax1(epsabs,epsrel*dabs(result))) ier = 0
369c ***jump out of do-loop
370      if (ier.eq.0) go to 999
371   70 continue
372   80 call xerror('abnormal return from dqng ',26,ier,0)
373  999 return
374      end
375