-
Notifications
You must be signed in to change notification settings - Fork 3
/
quadpack.f90
8740 lines (7991 loc) · 263 KB
/
quadpack.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!--------------------------------------------------------------------------------
!
! Copyright (C) 2017 L. J. Allen, H. G. Brown, A. J. D’Alfonso, S.D. Findlay, B. D. Forbes
!
! This program is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with this program. If not, see <http://www.gnu.org/licenses/>.
!
!--------------------------------------------------------------------------------
!This module is adapted from the quadpack library for numerical integration,
!available at <http://people.sc.fsu.edu/~jburkardt/f_src/quadpack/quadpack.html>.
!The only modifications were to make make the real variables double precision (64-bit)
!and to place the routines within a module. The original MuSTEM codes used the
!Romberg integration routine described in the book "Numerical Recipes in FORTRAN: The
!Art of Scientific Computing" by Press et al. (1986), however the copyright conditions
!of that software package prevented them being included in this open source implementation
!of MuSTEM.
module quadpack
contains
subroutine aaaa
!*****************************************************************************80
!
!! AAAA is a dummy subroutine with QUADPACK documentation in its comments.
!
! 1. introduction
!
! quadpack is a fortran subroutine package for the numerical
! computation of definite one-dimensional integrals. it originated
! from a joint project of r. piessens and e. de doncker (appl.
! math. and progr. div.- k.u.leuven, belgium), c. ueberhuber (inst.
! fuer math.- techn.u.wien, austria), and d. kahaner (nation. bur.
! of standards- washington d.c., u.s.a.).
!
! 2. survey
!
! - qags : is an integrator based on globally adaptive interval
! subdivision in connection with extrapolation (de doncker,
! 1978) by the epsilon algorithm (wynn, 1956).
!
! - qagp : serves the same purposes as qags, but also allows
! for eventual user-supplied information, i.e. the
! abscissae of internal singularities, discontinuities
! and other difficulties of the integrand function.
! the algorithm is a modification of that in qags.
!
! - qagi : handles integration over infinite intervals. the
! infinite range is mapped onto a finite interval and
! then the same strategy as in qags is applied.
!
! - qawo : is a routine for the integration of cos(omega*x)*f(x)
! or sin(omega*x)*f(x) over a finite interval (a,b).
! omega is is specified by the user
! the rule evaluation component is based on the
! modified clenshaw-curtis technique.
! an adaptive subdivision scheme is used connected with
! an extrapolation procedure, which is a modification
! of that in qags and provides the possibility to deal
! even with singularities in f.
!
! - qawf : calculates the fourier cosine or fourier sine
! transform of f(x), for user-supplied interval (a,
! infinity), omega, and f. the procedure of qawo is
! used on successive finite intervals, and convergence
! acceleration by means of the epsilon algorithm (wynn,
! 1956) is applied to the series of the integral
! contributions.
!
! - qaws : integrates w(x)*f(x) over (a,b) with a < b finite,
! and w(x) = ((x-a)**alfa)*((b-x)**beta)*v(x)
! where v(x) = 1 or log(x-a) or log(b-x)
! or log(x-a)*log(b-x)
! and -1 < alfa, -1 < beta.
! the user specifies a, b, alfa, beta and the type of
! the function v.
! a globally adaptive subdivision strategy is applied,
! with modified clenshaw-curtis integration on the
! subintervals which contain a or b.
!
! - qawc : computes the cauchy principal value of f(x)/(x-c)
! over a finite interval (a,b) and for
! user-determined c.
! the strategy is globally adaptive, and modified
! clenshaw-curtis integration is used on the subranges
! which contain the point x = c.
!
! each of the routines above also has a "more detailed" version
! with a name ending in e, as qage. these provide more
! information and control than the easier versions.
!
!
! the preceeding routines are all automatic. that is, the user
! inputs his problem and an error tolerance. the routine
! attempts to perform the integration to within the requested
! absolute or relative error.
! there are, in addition, a number of non-automatic integrators.
! these are most useful when the problem is such that the
! user knows that a fixed rule will provide the accuracy
! required. typically they return an error estimate but make
! no attempt to satisfy any particular input error request.
!
! qk15
! qk21
! qk31
! qk41
! qk51
! qk61
! estimate the integral on [a,b] using 15, 21,..., 61
! point rule and return an error estimate.
! qk15i 15 point rule for (semi)infinite interval.
! qk15w 15 point rule for special singular weight functions.
! qc25c 25 point rule for cauchy principal values
! qc25o 25 point rule for sin/cos integrand.
! qmomo integrates k-th degree chebychev polynomial times
! function with various explicit singularities.
!
! 3. guidelines for the use of quadpack
!
! here it is not our purpose to investigate the question when
! automatic quadrature should be used. we shall rather attempt
! to help the user who already made the decision to use quadpack,
! with selecting an appropriate routine or a combination of
! several routines for handling his problem.
!
! for both quadrature over finite and over infinite intervals,
! one of the first questions to be answered by the user is
! related to the amount of computer time he wants to spend,
! versus his -own- time which would be needed, for example, for
! manual subdivision of the interval or other analytic
! manipulations.
!
! (1) the user may not care about computer time, or not be
! willing to do any analysis of the problem. especially when
! only one or a few integrals must be calculated, this attitude
! can be perfectly reasonable. in this case it is clear that
! either the most sophisticated of the routines for finite
! intervals, qags, must be used, or its analogue for infinite
! intervals, qagi. these routines are able to cope with
! rather difficult, even with improper integrals.
! this way of proceeding may be expensive. but the integrator
! is supposed to give you an answer in return, with additional
! information in the case of a failure, through its error
! estimate and flag. yet it must be stressed that the programs
! cannot be totally reliable.
!
! (2) the user may want to examine the integrand function.
! if bad local difficulties occur, such as a discontinuity, a
! singularity, derivative singularity or high peak at one or
! more points within the interval, the first advice is to
! split up the interval at these points. the integrand must
! then be examinated over each of the subintervals separately,
! so that a suitable integrator can be selected for each of
! them. if this yields problems involving relative accuracies
! to be imposed on -finite- subintervals, one can make use of
! qagp, which must be provided with the positions of the local
! difficulties. however, if strong singularities are present
! and a high accuracy is requested, application of qags on the
! subintervals may yield a better result.
!
! for quadrature over finite intervals we thus dispose of qags
! and
! - qng for well-behaved integrands,
! - qag for functions with an oscillating behavior of a non
! specific type,
! - qawo for functions, eventually singular, containing a
! factor cos(omega*x) or sin(omega*x) where omega is known,
! - qaws for integrands with algebraico-logarithmic end point
! singularities of known type,
! - qawc for cauchy principal values.
!
! remark
!
! on return, the work arrays in the argument lists of the
! adaptive integrators contain information about the interval
! subdivision process and hence about the integrand behavior:
! the end points of the subintervals, the local integral
! contributions and error estimates, and eventually other
! characteristics. for this reason, and because of its simple
! globally adaptive nature, the routine qag in particular is
! well-suited for integrand examination. difficult spots can
! be located by investigating the error estimates on the
! subintervals.
!
! for infinite intervals we provide only one general-purpose
! routine, qagi. it is based on the qags algorithm applied
! after a transformation of the original interval into (0,1).
! yet it may eventuate that another type of transformation is
! more appropriate, or one might prefer to break up the
! original interval and use qagi only on the infinite part
! and so on. these kinds of actions suggest a combined use of
! different quadpack integrators. note that, when the only
! difficulty is an integrand singularity at the finite
! integration limit, it will in general not be necessary to
! break up the interval, as qagi deals with several types of
! singularity at the boundary point of the integration range.
! it also handles slowly convergent improper integrals, on
! the condition that the integrand does not oscillate over
! the entire infinite interval. if it does we would advise
! to sum succeeding positive and negative contributions to
! the integral -e.g. integrate between the zeros- with one
! or more of the finite-range integrators, and apply
! convergence acceleration eventually by means of quadpack
! subroutine qelg which implements the epsilon algorithm.
! such quadrature problems include the fourier transform as
! a special case. yet for the latter we have an automatic
! integrator available, qawf.
!
return
end subroutine
subroutine qag ( f, a, b, epsabs, epsrel, key, result, abserr, neval, ier )
!*****************************************************************************80
!
!! QAG approximates an integral over a finite interval.
!
! Discussion:
!
! The routine calculates an approximation RESULT to a definite integral
! I = integral of F over (A,B),
! hopefully satisfying
! || I - RESULT || <= max ( EPSABS, EPSREL * ||I|| ).
!
! QAG is a simple globally adaptive integrator using the strategy of
! Aind (Piessens, 1973). It is possible to choose between 6 pairs of
! Gauss-Kronrod quadrature formulae for the rule evaluation component.
! The pairs of high degree of precision are suitable for handling
! integration difficulties due to a strongly oscillating integrand.
!
! Author:
!
! Robert Piessens, Elise de Doncker-Kapenger,
! Christian Ueberhuber, David Kahaner
!
! Reference:
!
! Robert Piessens, Elise de Doncker-Kapenger,
! Christian Ueberhuber, David Kahaner,
! QUADPACK, a Subroutine Package for Automatic Integration,
! Springer Verlag, 1983
!
! Parameters:
!
! Input, external real ( kind = 8 ) F, the name of the function routine, of the form
! function f ( x )
! real ( kind = 8 ) f
! real ( kind = 8 ) x
! which evaluates the integrand function.
!
! Input, real ( kind = 8 ) A, B, the limits of integration.
!
! Input, real ( kind = 8 ) EPSABS, EPSREL, the absolute and relative accuracy requested.
!
! Input, integer ( kind = 4 ) KEY, chooses the order of the local integration rule:
! 1, 7 Gauss points, 15 Gauss-Kronrod points,
! 2, 10 Gauss points, 21 Gauss-Kronrod points,
! 3, 15 Gauss points, 31 Gauss-Kronrod points,
! 4, 20 Gauss points, 41 Gauss-Kronrod points,
! 5, 25 Gauss points, 51 Gauss-Kronrod points,
! 6, 30 Gauss points, 61 Gauss-Kronrod points.
!
! Output, real ( kind = 8 ) RESULT, the estimated value of the integral.
!
! Output, real ( kind = 8 ) ABSERR, an estimate of || I - RESULT ||.
!
! Output, integer ( kind = 4 ) NEVAL, the number of times the integral was evaluated.
!
! Output, integer ( kind = 4 ) IER, return code.
! 0, normal and reliable termination of the routine. It is assumed that the
! requested accuracy has been achieved.
! 1, maximum number of subdivisions allowed has been achieved. One can
! allow more subdivisions by increasing the value of LIMIT in QAG.
! However, if this yields no improvement it is advised to analyze the
! integrand to determine the integration difficulties. If the position
! of a local difficulty can be determined, such as a singularity or
! 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.
! 2, the occurrence of roundoff error is detected, which prevents the
! requested tolerance from being achieved.
! 3, extremely bad integrand behavior occurs at some points of the
! integration interval.
! 6, the input is invalid, because EPSABS < 0 and EPSREL < 0.
!
! Local parameters:
!
! LIMIT is the maximum number of subintervals allowed in
! the subdivision process of QAGE.
!
implicit none
integer ( kind = 4 ), parameter :: limit = 500
real ( kind = 8 ),intent(in):: a
real ( kind = 8 ) abserr
real ( kind = 8 ) alist(limit)
real ( kind = 8 ),intent(in):: b
real ( kind = 8 ) blist(limit)
real ( kind = 8 ) elist(limit)
real ( kind = 8 ) epsabs
real ( kind = 8 ) epsrel
real ( kind = 8 ), external :: f
integer ( kind = 4 ) ier
integer ( kind = 4 ) iord(limit)
integer ( kind = 4 ) key
integer ( kind = 4 ) last
integer ( kind = 4 ) neval
real ( kind = 8 ) result
real ( kind = 8 ) rlist(limit)
call qage ( f, a, b, epsabs, epsrel, key, limit, result, abserr, neval, &
ier, alist, blist, rlist, elist, iord, last )
return
end subroutine
subroutine qage ( f, a, b, epsabs, epsrel, key, limit, result, abserr, neval, &
ier, alist, blist, rlist, elist, iord, last )
!*****************************************************************************80
!
!! QAGE estimates a definite integral.
!
! Discussion:
!
! The routine calculates an approximation RESULT to a definite integral
! I = integral of F over (A,B),
! hopefully satisfying
! || I - RESULT || <= max ( EPSABS, EPSREL * ||I|| ).
!
! Author:
!
! Robert Piessens, Elise de Doncker-Kapenger,
! Christian Ueberhuber, David Kahaner
!
! Reference:
!
! Robert Piessens, Elise de Doncker-Kapenger,
! Christian Ueberhuber, David Kahaner,
! QUADPACK, a Subroutine Package for Automatic Integration,
! Springer Verlag, 1983
!
! Parameters:
!
! Input, external real ( kind = 8 ) F, the name of the function routine, of the form
! function f ( x )
! real ( kind = 8 ) f
! real ( kind = 8 ) x
! which evaluates the integrand function.
!
! Input, real ( kind = 8 ) A, B, the limits of integration.
!
! Input, real ( kind = 8 ) EPSABS, EPSREL, the absolute and relative accuracy requested.
!
! Input, integer ( kind = 4 ) KEY, chooses the order of the local integration rule:
! 1, 7 Gauss points, 15 Gauss-Kronrod points,
! 2, 10 Gauss points, 21 Gauss-Kronrod points,
! 3, 15 Gauss points, 31 Gauss-Kronrod points,
! 4, 20 Gauss points, 41 Gauss-Kronrod points,
! 5, 25 Gauss points, 51 Gauss-Kronrod points,
! 6, 30 Gauss points, 61 Gauss-Kronrod points.
!
! Input, integer ( kind = 4 ) LIMIT, the maximum number of subintervals that
! can be used.
!
! Output, real ( kind = 8 ) RESULT, the estimated value of the integral.
!
! Output, real ( kind = 8 ) ABSERR, an estimate of || I - RESULT ||.
!
! Output, integer ( kind = 4 ) NEVAL, the number of times the integral was evaluated.
!
! Output, integer ( kind = 4 ) IER, return code.
! 0, normal and reliable termination of the routine. It is assumed that the
! requested accuracy has been achieved.
! 1, maximum number of subdivisions allowed has been achieved. One can
! allow more subdivisions by increasing the value of LIMIT in QAG.
! However, if this yields no improvement it is advised to analyze the
! integrand to determine the integration difficulties. If the position
! of a local difficulty can be determined, such as a singularity or
! 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.
! 2, the occurrence of roundoff error is detected, which prevents the
! requested tolerance from being achieved.
! 3, extremely bad integrand behavior occurs at some points of the
! integration interval.
! 6, the input is invalid, because EPSABS < 0 and EPSREL < 0.
!
! Workspace, real ( kind = 8 ) ALIST(LIMIT), BLIST(LIMIT), contains in entries 1
! through LAST the left and right ends of the partition subintervals.
!
! Workspace, real ( kind = 8 ) RLIST(LIMIT), contains in entries 1 through LAST
! the integral approximations on the subintervals.
!
! Workspace, real ( kind = 8 ) ELIST(LIMIT), contains in entries 1 through LAST
! the absolute error estimates on the subintervals.
!
! Output, integer ( kind = 4 ) IORD(LIMIT), the first K elements of which are pointers
! to the error estimates over the subintervals, such that
! elist(iord(1)), ..., elist(iord(k)) form a decreasing sequence, with
! k = last if last <= (limit/2+2), and k = limit+1-last otherwise.
!
! Output, integer ( kind = 4 ) LAST, the number of subintervals actually produced
! in the subdivision process.
!
! Local parameters:
!
! alist - list of left end points of all subintervals
! considered up to now
! blist - list of right end points of all subintervals
! considered up to now
! elist(i) - error estimate applying to rlist(i)
! maxerr - pointer to the interval with largest error estimate
! errmax - elist(maxerr)
! area - sum of the integrals over the subintervals
! errsum - sum of the errors over the subintervals
! errbnd - requested accuracy max(epsabs,epsrel*abs(result))
! *****1 - variable for the left subinterval
! *****2 - variable for the right subinterval
! last - index for subdivision
!
implicit none
integer ( kind = 4 ) limit
real ( kind = 8 ),intent(in):: a
real ( kind = 8 ) abserr
real ( kind = 8 ) alist(limit)
real ( kind = 8 ) area
real ( kind = 8 ) area1
real ( kind = 8 ) area12
real ( kind = 8 ) area2
real ( kind = 8 ) a1
real ( kind = 8 ) a2
real ( kind = 8 ),intent(in):: b
real ( kind = 8 ) blist(limit)
real ( kind = 8 ) b1
real ( kind = 8 ) b2
real ( kind = 8 ) c
real ( kind = 8 ) defabs
real ( kind = 8 ) defab1
real ( kind = 8 ) defab2
real ( kind = 8 ) elist(limit)
real ( kind = 8 ) epsabs
real ( kind = 8 ) epsrel
real ( kind = 8 ) errbnd
real ( kind = 8 ) errmax
real ( kind = 8 ) error1
real ( kind = 8 ) error2
real ( kind = 8 ) erro12
real ( kind = 8 ) errsum
real ( kind = 8 ), external :: f
integer ( kind = 4 ) ier
integer ( kind = 4 ) iord(limit)
integer ( kind = 4 ) iroff1
integer ( kind = 4 ) iroff2
integer ( kind = 4 ) key
integer ( kind = 4 ) keyf
integer ( kind = 4 ) last
integer ( kind = 4 ) maxerr
integer ( kind = 4 ) neval
integer ( kind = 4 ) nrmax
real ( kind = 8 ) resabs
real ( kind = 8 ) result
real ( kind = 8 ) rlist(limit)
!
! Test on validity of parameters.
!
ier = 0
neval = 0
last = 0
result = 0.0E+00
abserr = 0.0E+00
alist(1) = a
blist(1) = b
rlist(1) = 0.0E+00
elist(1) = 0.0E+00
iord(1) = 0
if ( epsabs < 0.0E+00 .and. epsrel < 0.0E+00 ) then
ier = 6
return
end if
!
! First approximation to the integral.
!
keyf = key
keyf = max ( keyf, 1 )
keyf = min ( keyf, 6 )
c = keyf
neval = 0
if ( keyf == 1 ) then
call qk15 ( f, a, b, result, abserr, defabs, resabs )
else if ( keyf == 2 ) then
call qk21 ( f, a, b, result, abserr, defabs, resabs )
else if ( keyf == 3 ) then
call qk31 ( f, a, b, result, abserr, defabs, resabs )
else if ( keyf == 4 ) then
call qk41 ( f, a, b, result, abserr, defabs, resabs )
else if ( keyf == 5 ) then
call qk51 ( f, a, b, result, abserr, defabs, resabs )
else if ( keyf == 6 ) then
call qk61 ( f, a, b, result, abserr, defabs, resabs )
end if
last = 1
rlist(1) = result
elist(1) = abserr
iord(1) = 1
!
! Test on accuracy.
!
errbnd = max ( epsabs, epsrel * abs ( result ) )
if ( abserr <= 5.0E+01 * epsilon ( defabs ) * defabs .and. &
errbnd < abserr ) then
ier = 2
end if
if ( limit == 1 ) then
ier = 1
end if
if ( ier /= 0 .or. &
( abserr <= errbnd .and. abserr /= resabs ) .or. &
abserr == 0.0E+00 ) then
if ( keyf /= 1 ) then
neval = (10*keyf+1) * (2*neval+1)
else
neval = 30 * neval + 15
end if
return
end if
!
! Initialization.
!
errmax = abserr
maxerr = 1
area = result
errsum = abserr
nrmax = 1
iroff1 = 0
iroff2 = 0
do last = 2, limit
!
! Bisect the subinterval with the largest error estimate.
!
a1 = alist(maxerr)
b1 = 0.5E+00 * ( alist(maxerr) + blist(maxerr) )
a2 = b1
b2 = blist(maxerr)
if ( keyf == 1 ) then
call qk15 ( f, a1, b1, area1, error1, resabs, defab1 )
else if ( keyf == 2 ) then
call qk21 ( f, a1, b1, area1, error1, resabs, defab1 )
else if ( keyf == 3 ) then
call qk31 ( f, a1, b1, area1, error1, resabs, defab1 )
else if ( keyf == 4 ) then
call qk41 ( f, a1, b1, area1, error1, resabs, defab1)
else if ( keyf == 5 ) then
call qk51 ( f, a1, b1, area1, error1, resabs, defab1 )
else if ( keyf == 6 ) then
call qk61 ( f, a1, b1, area1, error1, resabs, defab1 )
end if
if ( keyf == 1 ) then
call qk15 ( f, a2, b2, area2, error2, resabs, defab2 )
else if ( keyf == 2 ) then
call qk21 ( f, a2, b2, area2, error2, resabs, defab2 )
else if ( keyf == 3 ) then
call qk31 ( f, a2, b2, area2, error2, resabs, defab2 )
else if ( keyf == 4 ) then
call qk41 ( f, a2, b2, area2, error2, resabs, defab2 )
else if ( keyf == 5 ) then
call qk51 ( f, a2, b2, area2, error2, resabs, defab2 )
else if ( keyf == 6 ) then
call qk61 ( f, a2, b2, area2, error2, resabs, defab2 )
end if
!
! Improve previous approximations to integral and error and
! test for accuracy.
!
neval = neval + 1
area12 = area1 + area2
erro12 = error1 + error2
errsum = errsum + erro12 - errmax
area = area + area12 - rlist(maxerr)
if ( defab1 /= error1 .and. defab2 /= error2 ) then
if ( abs ( rlist(maxerr) - area12 ) <= 1.0E-05 * abs ( area12 ) &
.and. 9.9E-01 * errmax <= erro12 ) then
iroff1 = iroff1 + 1
end if
if ( 10 < last .and. errmax < erro12 ) then
iroff2 = iroff2 + 1
end if
end if
rlist(maxerr) = area1
rlist(last) = area2
errbnd = max ( epsabs, epsrel * abs ( area ) )
!
! Test for roundoff error and eventually set error flag.
!
if ( errbnd < errsum ) then
if ( 6 <= iroff1 .or. 20 <= iroff2 ) then
ier = 2
end if
!
! Set error flag in the case that the number of subintervals
! equals limit.
!
if ( last == limit ) then
ier = 1
end if
!
! Set error flag in the case of bad integrand behavior
! at a point of the integration range.
!
if ( max ( abs ( a1 ), abs ( b2 ) ) <= ( 1.0E+00 + c * 1.0E+03 * &
epsilon ( a1 ) ) * ( abs ( a2 ) + 1.0E+04 * tiny ( a2 ) ) ) then
ier = 3
end if
end if
!
! Append the newly-created intervals to the list.
!
if ( error2 <= error1 ) then
alist(last) = a2
blist(maxerr) = b1
blist(last) = b2
elist(maxerr) = error1
elist(last) = error2
else
alist(maxerr) = a2
alist(last) = a1
blist(last) = b1
rlist(maxerr) = area2
rlist(last) = area1
elist(maxerr) = error2
elist(last) = error1
end if
!
! Call QSORT to maintain the descending ordering
! in the list of error estimates and select the subinterval
! with the largest error estimate (to be bisected next).
!
call qsort ( limit, last, maxerr, errmax, elist, iord, nrmax )
if ( ier /= 0 .or. errsum <= errbnd ) then
exit
end if
end do
!
! Compute final result.
!
result = sum ( rlist(1:last) )
abserr = errsum
if ( keyf /= 1 ) then
neval = ( 10 * keyf + 1 ) * ( 2 * neval + 1 )
else
neval = 30 * neval + 15
end if
return
end subroutine
subroutine qagi ( f, bound, inf, epsabs, epsrel, result, abserr, neval, ier )
!*****************************************************************************80
!
!! QAGI estimates an integral over a semi-infinite or infinite interval.
!
! Discussion:
!
! The routine calculates an approximation RESULT to a definite integral
! I = integral of F over (A, +Infinity),
! or
! I = integral of F over (-Infinity,A)
! or
! I = integral of F over (-Infinity,+Infinity),
! hopefully satisfying
! || I - RESULT || <= max ( EPSABS, EPSREL * ||I|| ).
!
! Author:
!
! Robert Piessens, Elise de Doncker-Kapenger,
! Christian Ueberhuber, David Kahaner
!
! Reference:
!
! Robert Piessens, Elise de Doncker-Kapenger,
! Christian Ueberhuber, David Kahaner,
! QUADPACK, a Subroutine Package for Automatic Integration,
! Springer Verlag, 1983
!
! Parameters:
!
! Input, external real ( kind = 8 ) F, the name of the function routine, of the form
! function f ( x )
! real ( kind = 8 ) f
! real ( kind = 8 ) x
! which evaluates the integrand function.
!
! Input, real ( kind = 8 ) BOUND, the value of the finite endpoint of the integration
! range, if any, that is, if INF is 1 or -1.
!
! Input, integer ( kind = 4 ) INF, indicates the type of integration range.
! 1: ( BOUND, +Infinity),
! -1: ( -Infinity, BOUND),
! 2: ( -Infinity, +Infinity).
!
! Input, real ( kind = 8 ) EPSABS, EPSREL, the absolute and relative accuracy requested.
!
! Output, real ( kind = 8 ) RESULT, the estimated value of the integral.
!
! Output, real ( kind = 8 ) ABSERR, an estimate of || I - RESULT ||.
!
! Output, integer ( kind = 4 ) NEVAL, the number of times the integral was evaluated.
!
! Output, integer ( kind = 4 ) IER, error indicator.
! 0, normal and reliable termination of the routine. It is assumed that
! the requested accuracy has been achieved.
! > 0, abnormal termination of the routine. The estimates for result
! and error are less reliable. It is assumed that the requested
! accuracy has not been achieved.
! 1, maximum number of subdivisions allowed has been achieved. One can
! allow more subdivisions by increasing the data value of LIMIT in QAGI
! (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.
! 2, the occurrence of roundoff error is detected, which prevents the
! requested tolerance from being achieved. The error may be
! under-estimated.
! 3, extremely bad integrand behavior occurs at some points of the
! integration interval.
! 4, the algorithm does not converge. Roundoff error is detected in the
! extrapolation table. It is assumed that the requested tolerance
! cannot be achieved, and that the returned result is the best which
! can be obtained.
! 5, the integral is probably divergent, or slowly convergent. It must
! be noted that divergence can occur with any other value of IER.
! 6, the input is invalid, because INF /= 1 and INF /= -1 and INF /= 2, or
! epsabs < 0 and epsrel < 0. result, abserr, neval are set to zero.
!
! Local parameters:
!
! the dimension of rlist2 is determined by the value of
! limexp in QEXTR.
!
! alist - list of left end points of all subintervals
! considered up to now
! blist - list of right end points of all subintervals
! considered up to now
! rlist(i) - approximation to the integral over
! (alist(i),blist(i))
! rlist2 - array of dimension at least (limexp+2),
! containing the part of the epsilon table
! which is still needed for further computations
! elist(i) - error estimate applying to rlist(i)
! maxerr - pointer to the interval with largest error
! estimate
! errmax - elist(maxerr)
! erlast - error on the interval currently subdivided
! (before that subdivision has taken place)
! area - sum of the integrals over the subintervals
! errsum - sum of the errors over the subintervals
! errbnd - requested accuracy max(epsabs,epsrel*
! abs(result))
! *****1 - variable for the left subinterval
! *****2 - variable for the right subinterval
! last - index for subdivision
! nres - number of calls to the extrapolation routine
! numrl2 - number of elements currently in rlist2. if an
! appropriate approximation to the compounded
! integral has been obtained, it is put in
! rlist2(numrl2) after numrl2 has been increased
! by one.
! small - length of the smallest interval considered up
! to now, multiplied by 1.5
! erlarg - sum of the errors over the intervals larger
! than the smallest interval considered up to now
! extrap - logical variable denoting that the routine
! is attempting to perform extrapolation. i.e.
! before subdividing the smallest interval we
! try to decrease the value of erlarg.
! noext - logical variable denoting that extrapolation
! is no longer allowed (true-value)
!
implicit none
integer ( kind = 4 ), parameter :: limit = 500
real ( kind = 8 ) abseps
real ( kind = 8 ) abserr
real ( kind = 8 ) alist(limit)
real ( kind = 8 ) area
real ( kind = 8 ) area1
real ( kind = 8 ) area12
real ( kind = 8 ) area2
real ( kind = 8 ) a1
real ( kind = 8 ) a2
real ( kind = 8 ) blist(limit)
real ( kind = 8 ) boun
real ( kind = 8 ) bound
real ( kind = 8 ) b1
real ( kind = 8 ) b2
real ( kind = 8 ) correc
real ( kind = 8 ) defabs
real ( kind = 8 ) defab1
real ( kind = 8 ) defab2
real ( kind = 8 ) dres
real ( kind = 8 ) elist(limit)
real ( kind = 8 ) epsabs
real ( kind = 8 ) epsrel
real ( kind = 8 ) erlarg
real ( kind = 8 ) erlast
real ( kind = 8 ) errbnd
real ( kind = 8 ) errmax
real ( kind = 8 ) error1
real ( kind = 8 ) error2
real ( kind = 8 ) erro12
real ( kind = 8 ) errsum
real ( kind = 8 ) ertest
logical extrap
real ( kind = 8 ), external :: f
integer ( kind = 4 ) id
integer ( kind = 4 ) ier
integer ( kind = 4 ) ierro
integer ( kind = 4 ) inf
integer ( kind = 4 ) iord(limit)
integer ( kind = 4 ) iroff1
integer ( kind = 4 ) iroff2
integer ( kind = 4 ) iroff3
integer ( kind = 4 ) jupbnd
integer ( kind = 4 ) k
integer ( kind = 4 ) ksgn
integer ( kind = 4 ) ktmin
integer ( kind = 4 ) last
integer ( kind = 4 ) maxerr
integer ( kind = 4 ) neval
logical noext
integer ( kind = 4 ) nres
integer ( kind = 4 ) nrmax
integer ( kind = 4 ) numrl2
real ( kind = 8 ) resabs
real ( kind = 8 ) reseps
real ( kind = 8 ) result
real ( kind = 8 ) res3la(3)
real ( kind = 8 ) rlist(limit)
real ( kind = 8 ) rlist2(52)
real ( kind = 8 ) small
!
! Test on validity of parameters.
!
ier = 0
neval = 0
last = 0
result = 0.0E+00
abserr = 0.0E+00
alist(1) = 0.0E+00
blist(1) = 1.0E+00
rlist(1) = 0.0E+00
elist(1) = 0.0E+00
iord(1) = 0
if ( epsabs < 0.0E+00 .and. epsrel < 0.0E+00 ) then
ier = 6
return
end if
!
! First approximation to the integral.
!
! Determine the interval to be mapped onto (0,1).
! If INF = 2 the integral is computed as i = i1+i2, where
! i1 = integral of f over (-infinity,0),
! i2 = integral of f over (0,+infinity).
!
if ( inf == 2 ) then
boun = 0.0E+00
else
boun = bound
end if
call qk15i ( f, boun, inf, 0.0_8, 1.0_8, result, abserr, defabs, resabs )
!
! Test on accuracy.
!
last = 1
rlist(1) = result
elist(1) = abserr
iord(1) = 1
dres = abs ( result )
errbnd = max ( epsabs, epsrel * dres )
if ( abserr <= 100.0E+00 * epsilon ( defabs ) * defabs .and. &
errbnd < abserr ) then
ier = 2
end if
if ( limit == 1 ) then
ier = 1
end if
if ( ier /= 0 .or. (abserr <= errbnd .and. abserr /= resabs ) .or. &
abserr == 0.0E+00 ) go to 130
!
! Initialization.
!
rlist2(1) = result
errmax = abserr
maxerr = 1
area = result
errsum = abserr
abserr = huge ( abserr )
nrmax = 1
nres = 0
ktmin = 0
numrl2 = 2
extrap = .false.
noext = .false.
ierro = 0
iroff1 = 0
iroff2 = 0
iroff3 = 0
if ( ( 1.0E+00 - 5.0E+01 * epsilon ( defabs ) ) * defabs <= dres ) then
ksgn = 1
else
ksgn = -1
end if
do last = 2, limit
!
! Bisect the subinterval with nrmax-th largest error estimate.
!
a1 = alist(maxerr)
b1 = 5.0E-01 * ( alist(maxerr) + blist(maxerr) )
a2 = b1
b2 = blist(maxerr)
erlast = errmax
call qk15i ( f, boun, inf, a1, b1, area1, error1, resabs, defab1 )
call qk15i ( f, boun, inf, a2, b2, area2, error2, resabs, defab2 )
!
! 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 ( defab1 /= error1 .and. defab2 /= error2 ) then
if ( abs ( rlist(maxerr) - area12 ) <= 1.0E-05 * abs ( area12 ) &
.and. 9.9E-01 * errmax <= erro12 ) then
if ( extrap ) then
iroff2 = iroff2 + 1
end if
if ( .not. extrap ) then