123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353 |
- *DECK QAGE
- SUBROUTINE QAGE (F, A, B, EPSABS, EPSREL, KEY, LIMIT, RESULT,
- + ABSERR, NEVAL, IER, ALIST, BLIST, RLIST, ELIST, IORD, LAST)
- C***BEGIN PROLOGUE QAGE
- C***PURPOSE The routine calculates an approximation result to a given
- C definite integral I = Integral of F over (A,B),
- C hopefully satisfying following claim for accuracy
- C ABS(I-RESLT).LE.MAX(EPSABS,EPSREL*ABS(I)).
- C***LIBRARY SLATEC (QUADPACK)
- C***CATEGORY H2A1A1
- C***TYPE SINGLE PRECISION (QAGE-S, DQAGE-D)
- C***KEYWORDS AUTOMATIC INTEGRATOR, GAUSS-KRONROD RULES,
- C GENERAL-PURPOSE, GLOBALLY ADAPTIVE, INTEGRAND EXAMINATOR,
- C QUADPACK, QUADRATURE
- C***AUTHOR Piessens, Robert
- C Applied Mathematics and Programming Division
- C K. U. Leuven
- C de Doncker, Elise
- C Applied Mathematics and Programming Division
- C K. U. Leuven
- C***DESCRIPTION
- C
- C Computation of a definite integral
- C Standard fortran subroutine
- C Real version
- C
- C PARAMETERS
- C ON ENTRY
- C F - Real
- C Function subprogram defining the integrand
- C function F(X). The actual name for F needs to be
- C declared E X T E R N A L in the driver program.
- C
- C A - Real
- C Lower limit of integration
- C
- C B - Real
- C Upper limit of integration
- C
- C EPSABS - Real
- C Absolute accuracy requested
- C EPSREL - Real
- C Relative accuracy requested
- C If EPSABS.LE.0
- C and EPSREL.LT.MAX(50*REL.MACH.ACC.,0.5D-28),
- C the routine will end with IER = 6.
- C
- C KEY - Integer
- C Key for choice of local integration rule
- C A Gauss-Kronrod pair is used with
- C 7 - 15 points if KEY.LT.2,
- C 10 - 21 points if KEY = 2,
- C 15 - 31 points if KEY = 3,
- C 20 - 41 points if KEY = 4,
- C 25 - 51 points if KEY = 5,
- C 30 - 61 points if KEY.GT.5.
- C
- C LIMIT - Integer
- C Gives an upper bound on the number of subintervals
- C in the partition of (A,B), LIMIT.GE.1.
- C
- C ON RETURN
- C RESULT - Real
- C Approximation to the integral
- C
- C ABSERR - Real
- C Estimate of the modulus of the absolute error,
- C which should equal or exceed ABS(I-RESULT)
- C
- C NEVAL - Integer
- C Number of integrand evaluations
- C
- C IER - Integer
- C IER = 0 Normal and reliable termination of the
- C routine. It is assumed that the requested
- C accuracy has been achieved.
- C IER.GT.0 Abnormal termination of the routine
- C The estimates for result and error are
- C less reliable. It is assumed that the
- C requested accuracy has not been achieved.
- C ERROR MESSAGES
- C IER = 1 Maximum number of subdivisions allowed
- C has been achieved. One can allow more
- C subdivisions by increasing the value
- C of LIMIT.
- C However, if this yields no improvement it
- C is rather advised to analyze the integrand
- C in order to determine the integration
- C difficulties. If the position of a local
- C difficulty can be determined(e.g.
- C SINGULARITY, DISCONTINUITY within the
- C interval) one will probably gain from
- C splitting up the interval at this point
- C and calling the integrator on the
- C subranges. If possible, an appropriate
- C special-purpose integrator should be used
- C which is designed for handling the type of
- C difficulty involved.
- C = 2 The occurrence of roundoff error is
- C detected, which prevents the requested
- C tolerance from being achieved.
- C = 3 Extremely bad integrand behaviour occurs
- C at some points of the integration
- C interval.
- C = 6 The input is invalid, because
- C (EPSABS.LE.0 and
- C EPSREL.LT.MAX(50*REL.MACH.ACC.,0.5D-28),
- C RESULT, ABSERR, NEVAL, LAST, RLIST(1) ,
- C ELIST(1) and IORD(1) are set to zero.
- C ALIST(1) and BLIST(1) are set to A and B
- C respectively.
- C
- C ALIST - Real
- C Vector of dimension at least LIMIT, the first
- C LAST elements of which are the left
- C end points of the subintervals in the partition
- C of the given integration range (A,B)
- C
- C BLIST - Real
- C Vector of dimension at least LIMIT, the first
- C LAST elements of which are the right
- C end points of the subintervals in the partition
- C of the given integration range (A,B)
- C
- C RLIST - Real
- C Vector of dimension at least LIMIT, the first
- C LAST elements of which are the
- C integral approximations on the subintervals
- C
- C ELIST - Real
- C Vector of dimension at least LIMIT, the first
- C LAST elements of which are the moduli of the
- C absolute error estimates on the subintervals
- C
- C IORD - Integer
- C Vector of dimension at least LIMIT, the first K
- C elements of which are pointers to the
- C error estimates over the subintervals,
- C such that ELIST(IORD(1)), ...,
- C ELIST(IORD(K)) form a decreasing sequence,
- C with K = LAST if LAST.LE.(LIMIT/2+2), and
- C K = LIMIT+1-LAST otherwise
- C
- C LAST - Integer
- C Number of subintervals actually produced in the
- C subdivision process
- C
- C***REFERENCES (NONE)
- C***ROUTINES CALLED QK15, QK21, QK31, QK41, QK51, QK61, QPSRT, R1MACH
- C***REVISION HISTORY (YYMMDD)
- C 800101 DATE WRITTEN
- C 890531 Changed all specific intrinsics to generic. (WRB)
- C 890831 Modified array declarations. (WRB)
- C 890831 REVISION DATE from Version 3.2
- C 891214 Prologue converted to Version 4.0 format. (BAB)
- C***END PROLOGUE QAGE
- C
- REAL A,ABSERR,ALIST,AREA,AREA1,AREA12,AREA2,A1,A2,B,BLIST,
- 1 B1,B2,DEFABS,DEFAB1,DEFAB2,R1MACH,ELIST,EPMACH,
- 2 EPSABS,EPSREL,ERRBND,ERRMAX,ERROR1,ERROR2,ERRO12,ERRSUM,F,
- 3 RESABS,RESULT,RLIST,UFLOW
- INTEGER IER,IORD,IROFF1,IROFF2,K,KEY,KEYF,LAST,
- 1 LIMIT,MAXERR,NEVAL,NRMAX
- C
- DIMENSION ALIST(*),BLIST(*),ELIST(*),IORD(*),
- 1 RLIST(*)
- C
- EXTERNAL F
- C
- C LIST OF MAJOR VARIABLES
- C -----------------------
- C
- C ALIST - LIST OF LEFT END POINTS OF ALL SUBINTERVALS
- C CONSIDERED UP TO NOW
- C BLIST - LIST OF RIGHT END POINTS OF ALL SUBINTERVALS
- C CONSIDERED UP TO NOW
- C RLIST(I) - APPROXIMATION TO THE INTEGRAL OVER
- C (ALIST(I),BLIST(I))
- C ELIST(I) - ERROR ESTIMATE APPLYING TO RLIST(I)
- C MAXERR - POINTER TO THE INTERVAL WITH LARGEST
- C ERROR ESTIMATE
- C ERRMAX - ELIST(MAXERR)
- C AREA - SUM OF THE INTEGRALS OVER THE SUBINTERVALS
- C ERRSUM - SUM OF THE ERRORS OVER THE SUBINTERVALS
- C ERRBND - REQUESTED ACCURACY MAX(EPSABS,EPSREL*
- C ABS(RESULT))
- C *****1 - VARIABLE FOR THE LEFT SUBINTERVAL
- C *****2 - VARIABLE FOR THE RIGHT SUBINTERVAL
- C LAST - INDEX FOR SUBDIVISION
- C
- C
- C MACHINE DEPENDENT CONSTANTS
- C ---------------------------
- C
- C EPMACH IS THE LARGEST RELATIVE SPACING.
- C UFLOW IS THE SMALLEST POSITIVE MAGNITUDE.
- C
- C***FIRST EXECUTABLE STATEMENT QAGE
- EPMACH = R1MACH(4)
- UFLOW = R1MACH(1)
- C
- C TEST ON VALIDITY OF PARAMETERS
- C ------------------------------
- C
- 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.LE.0.0E+00.AND.
- 1 EPSREL.LT.MAX(0.5E+02*EPMACH,0.5E-14)) IER = 6
- IF(IER.EQ.6) GO TO 999
- C
- C FIRST APPROXIMATION TO THE INTEGRAL
- C -----------------------------------
- C
- KEYF = KEY
- IF(KEY.LE.0) KEYF = 1
- IF(KEY.GE.7) KEYF = 6
- NEVAL = 0
- IF(KEYF.EQ.1) CALL QK15(F,A,B,RESULT,ABSERR,DEFABS,RESABS)
- IF(KEYF.EQ.2) CALL QK21(F,A,B,RESULT,ABSERR,DEFABS,RESABS)
- IF(KEYF.EQ.3) CALL QK31(F,A,B,RESULT,ABSERR,DEFABS,RESABS)
- IF(KEYF.EQ.4) CALL QK41(F,A,B,RESULT,ABSERR,DEFABS,RESABS)
- IF(KEYF.EQ.5) CALL QK51(F,A,B,RESULT,ABSERR,DEFABS,RESABS)
- IF(KEYF.EQ.6) CALL QK61(F,A,B,RESULT,ABSERR,DEFABS,RESABS)
- LAST = 1
- RLIST(1) = RESULT
- ELIST(1) = ABSERR
- IORD(1) = 1
- C
- C TEST ON ACCURACY.
- C
- ERRBND = MAX(EPSABS,EPSREL*ABS(RESULT))
- IF(ABSERR.LE.0.5E+02*EPMACH*DEFABS.AND.ABSERR.GT.
- 1 ERRBND) IER = 2
- IF(LIMIT.EQ.1) IER = 1
- IF(IER.NE.0.OR.(ABSERR.LE.ERRBND.AND.ABSERR.NE.RESABS)
- 1 .OR.ABSERR.EQ.0.0E+00) GO TO 60
- C
- C INITIALIZATION
- C --------------
- C
- C
- ERRMAX = ABSERR
- MAXERR = 1
- AREA = RESULT
- ERRSUM = ABSERR
- NRMAX = 1
- IROFF1 = 0
- IROFF2 = 0
- C
- C MAIN DO-LOOP
- C ------------
- C
- DO 30 LAST = 2,LIMIT
- C
- C BISECT THE SUBINTERVAL WITH THE LARGEST ERROR ESTIMATE.
- C
- A1 = ALIST(MAXERR)
- B1 = 0.5E+00*(ALIST(MAXERR)+BLIST(MAXERR))
- A2 = B1
- B2 = BLIST(MAXERR)
- IF(KEYF.EQ.1) CALL QK15(F,A1,B1,AREA1,ERROR1,RESABS,DEFAB1)
- IF(KEYF.EQ.2) CALL QK21(F,A1,B1,AREA1,ERROR1,RESABS,DEFAB1)
- IF(KEYF.EQ.3) CALL QK31(F,A1,B1,AREA1,ERROR1,RESABS,DEFAB1)
- IF(KEYF.EQ.4) CALL QK41(F,A1,B1,AREA1,ERROR1,RESABS,DEFAB1)
- IF(KEYF.EQ.5) CALL QK51(F,A1,B1,AREA1,ERROR1,RESABS,DEFAB1)
- IF(KEYF.EQ.6) CALL QK61(F,A1,B1,AREA1,ERROR1,RESABS,DEFAB1)
- IF(KEYF.EQ.1) CALL QK15(F,A2,B2,AREA2,ERROR2,RESABS,DEFAB2)
- IF(KEYF.EQ.2) CALL QK21(F,A2,B2,AREA2,ERROR2,RESABS,DEFAB2)
- IF(KEYF.EQ.3) CALL QK31(F,A2,B2,AREA2,ERROR2,RESABS,DEFAB2)
- IF(KEYF.EQ.4) CALL QK41(F,A2,B2,AREA2,ERROR2,RESABS,DEFAB2)
- IF(KEYF.EQ.5) CALL QK51(F,A2,B2,AREA2,ERROR2,RESABS,DEFAB2)
- IF(KEYF.EQ.6) CALL QK61(F,A2,B2,AREA2,ERROR2,RESABS,DEFAB2)
- C
- C IMPROVE PREVIOUS APPROXIMATIONS TO INTEGRAL
- C AND ERROR AND TEST FOR ACCURACY.
- C
- NEVAL = NEVAL+1
- AREA12 = AREA1+AREA2
- ERRO12 = ERROR1+ERROR2
- ERRSUM = ERRSUM+ERRO12-ERRMAX
- AREA = AREA+AREA12-RLIST(MAXERR)
- IF(DEFAB1.EQ.ERROR1.OR.DEFAB2.EQ.ERROR2) GO TO 5
- IF(ABS(RLIST(MAXERR)-AREA12).LE.0.1E-04*ABS(AREA12)
- 1 .AND.ERRO12.GE.0.99E+00*ERRMAX) IROFF1 = IROFF1+1
- IF(LAST.GT.10.AND.ERRO12.GT.ERRMAX) IROFF2 = IROFF2+1
- 5 RLIST(MAXERR) = AREA1
- RLIST(LAST) = AREA2
- ERRBND = MAX(EPSABS,EPSREL*ABS(AREA))
- IF(ERRSUM.LE.ERRBND) GO TO 8
- C
- C TEST FOR ROUNDOFF ERROR AND EVENTUALLY
- C SET ERROR FLAG.
- C
- IF(IROFF1.GE.6.OR.IROFF2.GE.20) IER = 2
- C
- C SET ERROR FLAG IN THE CASE THAT THE NUMBER OF
- C SUBINTERVALS EQUALS LIMIT.
- C
- IF(LAST.EQ.LIMIT) IER = 1
- C
- C SET ERROR FLAG IN THE CASE OF BAD INTEGRAND BEHAVIOUR
- C AT A POINT OF THE INTEGRATION RANGE.
- C
- IF(MAX(ABS(A1),ABS(B2)).LE.(0.1E+01+0.1E+03*
- 1 EPMACH)*(ABS(A2)+0.1E+04*UFLOW)) IER = 3
- C
- C APPEND THE NEWLY-CREATED INTERVALS TO THE LIST.
- C
- 8 IF(ERROR2.GT.ERROR1) GO TO 10
- ALIST(LAST) = A2
- BLIST(MAXERR) = B1
- BLIST(LAST) = B2
- ELIST(MAXERR) = ERROR1
- ELIST(LAST) = ERROR2
- GO TO 20
- 10 ALIST(MAXERR) = A2
- ALIST(LAST) = A1
- BLIST(LAST) = B1
- RLIST(MAXERR) = AREA2
- RLIST(LAST) = AREA1
- ELIST(MAXERR) = ERROR2
- ELIST(LAST) = ERROR1
- C
- C CALL SUBROUTINE QPSRT TO MAINTAIN THE DESCENDING ORDERING
- C IN THE LIST OF ERROR ESTIMATES AND SELECT THE
- C SUBINTERVAL WITH THE LARGEST ERROR ESTIMATE (TO BE
- C BISECTED NEXT).
- C
- 20 CALL QPSRT(LIMIT,LAST,MAXERR,ERRMAX,ELIST,IORD,NRMAX)
- C ***JUMP OUT OF DO-LOOP
- IF(IER.NE.0.OR.ERRSUM.LE.ERRBND) GO TO 40
- 30 CONTINUE
- C
- C COMPUTE FINAL RESULT.
- C ---------------------
- C
- 40 RESULT = 0.0E+00
- DO 50 K=1,LAST
- RESULT = RESULT+RLIST(K)
- 50 CONTINUE
- ABSERR = ERRSUM
- 60 IF(KEYF.NE.1) NEVAL = (10*KEYF+1)*(2*NEVAL+1)
- IF(KEYF.EQ.1) NEVAL = 30*NEVAL+15
- 999 RETURN
- END
|