1133 lines
		
	
	
		
			48 KiB
		
	
	
	
		
			Fortran
		
	
	
	
			
		
		
	
	
			1133 lines
		
	
	
		
			48 KiB
		
	
	
	
		
			Fortran
		
	
	
	
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
MODULE INTEGRATION
 | 
						||
!
 | 
						||
!  This module contains integration routines in order to integrate
 | 
						||
!    a function F over the interval [A,B].
 | 
						||
!
 | 
						||
!  These routines are:
 | 
						||
!
 | 
						||
!
 | 
						||
!   * Lagrange                              :  INTEGR_L(F,DR,NSIZE,NMAX,A,ID)
 | 
						||
!
 | 
						||
!   * N-point Gauss-Legendre                :  GAUSS_LEG(FCT,A,B,NGL,RES)
 | 
						||
!
 | 
						||
!   * CERNLIB adaptive Gauss quadrature     :  DGAUSS1(OM,KK,A,B,EPS)
 | 
						||
!
 | 
						||
!   * double exponential transformation     :  INTDE(F,A,B,EPS,I,ERR)
 | 
						||
!
 | 
						||
!   * fast double exponential transformation:  INTDE_F(F,A,B,AW,I,ERR)
 | 
						||
!
 | 
						||
!   * Romberg                               :  RBI1(FCT,A,B,PREC,OBTPREC,NITER,ITERMIN,ITERMAX)
 | 
						||
!
 | 
						||
!   * Simpson                               :  SIMPSON(FCT,A,B,N,RES)
 | 
						||
!
 | 
						||
!   *                                       :  QANC8(FCT,A,B,AERR,RERR,RES,ERR,NBF,FLG)
 | 
						||
!
 | 
						||
!
 | 
						||
      USE ACCURACY_REAL
 | 
						||
!
 | 
						||
CONTAINS
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      SUBROUTINE INTEGR_L(F,DR,NSIZE,NMAX,A,ID)
 | 
						||
!
 | 
						||
!.....Based on Lagrange integration formula 25.4.12
 | 
						||
!
 | 
						||
!     (See Table 25.3 for numerical coefficients) - Chapter 25 of
 | 
						||
!     Abramowitz & Stegun, "Handbook of mathematical functions",
 | 
						||
!     page 886 (Dover)
 | 
						||
!
 | 
						||
!
 | 
						||
!  Input parameters:
 | 
						||
!
 | 
						||
!       * F        : function to be integrated
 | 
						||
!       * DR       : constant grid step
 | 
						||
!       * NSIZE    : dimensioning of the arrays
 | 
						||
!       * NMAX     : index of upper limit of integration on the r mesh
 | 
						||
!       * ID       : integer parameter
 | 
						||
!                       ID  = 1  --> F0 = 0 at the origin
 | 
						||
!                       ID  > 1  --> F0 not 0 at the origin
 | 
						||
!
 | 
						||
!
 | 
						||
!  Output parameters:
 | 
						||
!
 | 
						||
!       * A        : integral result
 | 
						||
!
 | 
						||
!
 | 
						||
!                      --> Real function F case <--
 | 
						||
!
 | 
						||
!   Author :  C. R. Natoli
 | 
						||
!
 | 
						||
!                                    Last modified (DS) :  2 Nov 2020
 | 
						||
!
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : ZERO,FIVE,NINE
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP), INTENT(IN)  ::  F(NSIZE),DR
 | 
						||
      REAL (WP), INTENT(OUT) ::  A
 | 
						||
!
 | 
						||
      REAL (WP)              ::  H,A0,F0
 | 
						||
      REAL (WP)              ::  S720,S251,S646,S264
 | 
						||
      REAL (WP)              ::  S106,S19,S346,S456,S74,S11
 | 
						||
!
 | 
						||
      INTEGER, INTENT(IN)    ::  NSIZE,NMAX,ID
 | 
						||
!
 | 
						||
      INTEGER                ::  K0,KX,K
 | 
						||
!
 | 
						||
!  Coefficients given by table 25.3 p. 915:
 | 
						||
!
 | 
						||
      S720 = 720.0E0_WP                                             !
 | 
						||
      S251 = 251.0E0_WP                                             !
 | 
						||
      S646 = 646.0E0_WP                                             !
 | 
						||
      S264 = 264.0E0_WP                                             !
 | 
						||
      S106 = 106.0E0_WP                                             !
 | 
						||
      S19  =  19.0E0_WP                                             !
 | 
						||
      S346 = 346.0E0_WP                                             !
 | 
						||
      S456 = 456.0E0_WP                                             !
 | 
						||
      S74  =  74.0E0_WP                                             !
 | 
						||
      S11  =  11.0E0_WP                                             !
 | 
						||
!
 | 
						||
      H  = DR                                                       !
 | 
						||
      A0 = ZERO                                                     !
 | 
						||
!
 | 
						||
      IF(ID == 1) THEN                                              !
 | 
						||
        F0 = ZERO                                                   !
 | 
						||
        K0 = 0                                                      !
 | 
						||
      ELSE                                                          !
 | 
						||
        F0 = F(1)                                                   !
 | 
						||
        K0 = 1                                                      !
 | 
						||
      END IF                                                        !
 | 
						||
!
 | 
						||
      KX = NMAX                                                     !
 | 
						||
!
 | 
						||
      A = A0 + H * ( S251 * F0      + S646 * F(K0+1) -            & !
 | 
						||
                     S264 * F(K0+2) + S106 * F(K0+3) -            & !
 | 
						||
                     S19  * F(K0+4)                               & !
 | 
						||
                   ) / S720                                         !
 | 
						||
      A = A  + H * ( -S19 * F0      + S346 * F(K0+1) +            & !
 | 
						||
                     S456 * F(K0+2) - S74  * F(K0+3) +            & !
 | 
						||
                     S11  * F(K0+4)                               & !
 | 
						||
                    ) / S720                                        !
 | 
						||
      A = A  + H * ( S11  * F0      - S74  * F(K0+1) +            & !
 | 
						||
                     S456 * F(K0+2) + S346 * F(K0+3) -            & !
 | 
						||
                     S19  * F(K0+4)                               & !
 | 
						||
                   ) / S720                                         !
 | 
						||
!
 | 
						||
      K0 = K0 + 4                                                   !
 | 
						||
!
 | 
						||
      DO K = K0, KX                                                 !
 | 
						||
        A = A + H * ( NINE * F(K)   + 19.0E0_WP * F(K-1) -        & !
 | 
						||
                      FIVE * F(K-2) +  F(K-3)                     & !
 | 
						||
                    ) / 24.0E0_WP                                   !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      END SUBROUTINE INTEGR_L
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      SUBROUTINE GAUSS_LEG(FCT,A,B,NGL,RES)
 | 
						||
!
 | 
						||
!  This subroutine performs a Gauss-Legendre integration of
 | 
						||
!    the external function FCT
 | 
						||
!
 | 
						||
!
 | 
						||
!
 | 
						||
!   Author :  D. Sébilleau
 | 
						||
!
 | 
						||
!                                           Last modified :  5 Jun 2020
 | 
						||
!
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : ZERO
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      INTEGER, INTENT(IN)    :: NGL
 | 
						||
!
 | 
						||
      INTEGER                :: J
 | 
						||
!
 | 
						||
      REAL (WP), INTENT(IN)  ::  A,B
 | 
						||
      REAL (WP), INTENT(OUT) ::  RES
 | 
						||
!
 | 
						||
      REAL (WP)              ::  XGL(NGL),WGT(NGL)
 | 
						||
!
 | 
						||
      REAL (WP), EXTERNAL   ::  FCT
 | 
						||
!
 | 
						||
!  Construct Gauss-Legendre points from Numerical Recipes subroutine
 | 
						||
!
 | 
						||
      CALL GAULEG(A,B,XGL,WGT,NGL)                                  !
 | 
						||
!
 | 
						||
!  Performing the integral
 | 
						||
!
 | 
						||
      RES = ZERO                                                    !
 | 
						||
      DO J = 1, NGL                                                 !
 | 
						||
        RES = RES + WGT(J) * FCT(XGL(J))                            !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      END SUBROUTINE GAUSS_LEG
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      SUBROUTINE GAULEG(X1,X2,X,W,N)
 | 
						||
!
 | 
						||
!  Given the lower and upper limits of integration X1 and X2,
 | 
						||
!    and given N, this routine returns arrays X[1..N] and W[1..N]
 | 
						||
!    of length N, containing the abscissas and weights
 | 
						||
!    of the Gauss-Legendre N-point quadrature formula
 | 
						||
!
 | 
						||
!            This subroutine is taken from the book :
 | 
						||
!
 | 
						||
!          "Numerical Recipes : The Art of Scientific
 | 
						||
!           Computing" par W.H. Press, B.P. Flannery,
 | 
						||
!               S.A. Teukolsky et W.T. Vetterling
 | 
						||
!               (Cambridge University Press 1992)
 | 
						||
!
 | 
						||
!                       p. 145
 | 
						||
!
 | 
						||
!  Input parameters:
 | 
						||
!
 | 
						||
!       X1       :  lower limit of integration
 | 
						||
!       X2       :  upper limit of integration
 | 
						||
!       N        :  order of the Gauss-Legendre quadrature formula
 | 
						||
!
 | 
						||
!
 | 
						||
!  Output parameters:
 | 
						||
!
 | 
						||
!       X        :  abscissas for Gauss-Legendre N-point quadrature formula
 | 
						||
!       W        :  weights for Gauss-Legendre N-point quadrature formula
 | 
						||
!
 | 
						||
!
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : ZERO,ONE,TWO,HALF,FOURTH
 | 
						||
      USE PI_ETC,           ONLY : PI
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP)             ::  X1,X2,X(N),W(N)
 | 
						||
      REAL (WP)             ::  EPS
 | 
						||
      REAL (WP)             ::  XM,XL,Z,Z1,P1,P2,P3,PP
 | 
						||
!
 | 
						||
      REAL (WP)             ::  DCOS,DFLOAT,DABS
 | 
						||
!
 | 
						||
      INTEGER N,M,I,J
 | 
						||
!
 | 
						||
      EPS=3.0E-14_WP                                                !
 | 
						||
!
 | 
						||
      M=(N+1)/2                                                     ! the roots are symmetric
 | 
						||
      XM=HALF*(X2+X1)                                               ! in the interval, so we only
 | 
						||
      XL=HALF*(X2-X1)                                               ! have to find half of them
 | 
						||
!
 | 
						||
!  Loop over the desired roots
 | 
						||
!
 | 
						||
      DO I=1,M                                                      !
 | 
						||
!
 | 
						||
!  Starting with the approximation to the ith root,
 | 
						||
!    we enter the main loop of refinement by Newton’s method
 | 
						||
!
 | 
						||
        Z=DCOS(PI*(DFLOAT(I)-FOURTH)/(DFLOAT(N)+HALF))              ! approx for ith root
 | 
						||
!
 | 
						||
   1    CONTINUE                                                    !
 | 
						||
!
 | 
						||
        P1=ONE                                                      !
 | 
						||
        P2=ZERO                                                     !
 | 
						||
!
 | 
						||
!  Loop up the recurrence relation to get the
 | 
						||
!    Legendre polynomial evaluated at Z
 | 
						||
!
 | 
						||
        DO J=1,N                                                    !
 | 
						||
          P3=P2                                                     !
 | 
						||
          P2=P1                                                     !
 | 
						||
          P1=((TWO*DFLOAT(J)-ONE)*Z*P2-(DFLOAT(J)-ONE)*P3)    &     !
 | 
						||
               /DFLOAT(J)                                           !
 | 
						||
        END DO                                                      !
 | 
						||
!
 | 
						||
!  P1 is now the desired Legendre polynomial. We next compute PP,
 | 
						||
!    its derivative,by a standard relation involving also P2,
 | 
						||
!    the polynomial of one lower order
 | 
						||
!
 | 
						||
        PP=DFLOAT(N)*(Z*P1-P2)/(Z*Z-ONE)                            !
 | 
						||
        Z1=Z                                                        !
 | 
						||
        Z=Z1-P1/PP                                                  ! Newton’s method
 | 
						||
!
 | 
						||
        IF(DABS(Z-Z1) > EPS) GO TO 1                                !
 | 
						||
!
 | 
						||
!  Scale the root to the desired interval and put in
 | 
						||
!    its symmetric counterpart
 | 
						||
!
 | 
						||
        X(I)=XM-XL*Z                                                !
 | 
						||
        X(N+1-I)=XM+XL*Z                                            !
 | 
						||
!
 | 
						||
!  Compute the weight and its symmetric counterpart
 | 
						||
!
 | 
						||
        W(I)=TWO*XL/((ONE-Z*Z)*PP*PP)                               !
 | 
						||
        W(N+1-I)=W(I)                                               !
 | 
						||
!
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      END SUBROUTINE GAULEG
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      FUNCTION DGAUSS1(FCT,OM,KK,A,B,EPS)
 | 
						||
!
 | 
						||
!     ******************************************************************
 | 
						||
!
 | 
						||
!     ADAPTIVE GAUSSIAN QUADRATURE.
 | 
						||
!
 | 
						||
!     GAUSS IS SET EQUAL TO THE APPROXIMATE VALUE OF THE INTEGRAL OF
 | 
						||
!     THE FUNCTION F OVER THE INTERVAL (A,B), WITH ACCURACY PARAMETER
 | 
						||
!     EPS.
 | 
						||
!
 | 
						||
!     ******************************************************************
 | 
						||
!
 | 
						||
!  Originally written by K.S. Kölbig for CERNLIB
 | 
						||
!
 | 
						||
!                                          First version: 12 May 1966
 | 
						||
!                                          Revised      : 15 Mar 1993
 | 
						||
!
 | 
						||
! $Id: imp64.inc,v 1.1.1.1 1996/04/01 15:02:59 mclareni Exp $
 | 
						||
!
 | 
						||
! $Log: imp64.inc,v $
 | 
						||
! Revision 1.1.1.1  1996/04/01 15:02:59  mclareni
 | 
						||
! Mathlib gen
 | 
						||
!
 | 
						||
!
 | 
						||
! imp64.inc
 | 
						||
!
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : ZERO,ONE,TWO,FIVE
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP)             ::  W(12),X(12)
 | 
						||
      REAL (WP)             ::  CONST,AA,BB,U,S8,S16,C1,C2,H
 | 
						||
      REAL (WP)             ::  Z1,HF,CST
 | 
						||
      REAL (WP)             ::  A,B,EPS
 | 
						||
      REAL (WP)             ::  DGAUSS1
 | 
						||
      REAL (WP)             ::  OM
 | 
						||
!
 | 
						||
      REAL (WP)             ::  FCT
 | 
						||
!
 | 
						||
      REAL (WP)             ::  DABS
 | 
						||
!
 | 
						||
      INTEGER               ::  I,KK
 | 
						||
      INTEGER               ::  LOGF
 | 
						||
!
 | 
						||
      PARAMETER (Z1 = ONE, HF = Z1/TWO, CST = FIVE*Z1/1000.0E0_WP)  !
 | 
						||
!
 | 
						||
      DATA X( 1) /9.6028985649753623E-01_WP/                        !
 | 
						||
      DATA X( 2) /7.9666647741362674E-01_WP/                        !
 | 
						||
      DATA X( 3) /5.2553240991632899E-01_WP/                        !
 | 
						||
      DATA X( 4) /1.8343464249564980E-01_WP/                        !
 | 
						||
      DATA X( 5) /9.8940093499164993E-01_WP/                        !
 | 
						||
      DATA X( 6) /9.4457502307323258E-01_WP/                        !
 | 
						||
      DATA X( 7) /8.6563120238783174E-01_WP/                        !
 | 
						||
      DATA X( 8) /7.5540440835500303E-01_WP/                        !
 | 
						||
      DATA X( 9) /6.1787624440264375E-01_WP/                        !
 | 
						||
      DATA X(10) /4.5801677765722739E-01_WP/                        !
 | 
						||
      DATA X(11) /2.8160355077925891E-01_WP/                        !
 | 
						||
      DATA X(12) /9.5012509837637440E-02_WP/                        !
 | 
						||
!
 | 
						||
      DATA W( 1) /1.0122853629037626E-01_WP/                        !
 | 
						||
      DATA W( 2) /2.2238103445337447E-01_WP/                        !
 | 
						||
      DATA W( 3) /3.1370664587788729E-01_WP/                        !
 | 
						||
      DATA W( 4) /3.6268378337836198E-01_WP/                        !
 | 
						||
      DATA W( 5) /2.7152459411754095E-02_WP/                        !
 | 
						||
      DATA W( 6) /6.2253523938647893E-02_WP/                        !
 | 
						||
      DATA W( 7) /9.5158511682492785E-02_WP/                        !
 | 
						||
      DATA W( 8) /1.2462897125553387E-01_WP/                        !
 | 
						||
      DATA W( 9) /1.4959598881657673E-01_WP/                        !
 | 
						||
      DATA W(10) /1.6915651939500254E-01_WP/                        !
 | 
						||
      DATA W(11) /1.8260341504492359E-01_WP/                        !
 | 
						||
      DATA W(12) /1.8945061045506850E-01_WP/                        !
 | 
						||
!
 | 
						||
      H=ZERO                                                        !
 | 
						||
!
 | 
						||
      LOGF=6                                                        !
 | 
						||
!
 | 
						||
      IF(B == A) GO TO 99                                           !
 | 
						||
!
 | 
						||
      CONST=CST/DABS(B-A)                                           !
 | 
						||
      BB=A                                                          !
 | 
						||
!
 | 
						||
    1 AA=BB                                                         !
 | 
						||
      BB=B                                                          !
 | 
						||
!
 | 
						||
    2 C1=HF*(BB+AA)                                                 !
 | 
						||
      C2=HF*(BB-AA)                                                 !
 | 
						||
!
 | 
						||
      S8=ZERO                                                       !
 | 
						||
      DO I = 1,4                                                    !
 | 
						||
        U=C2*X(I)                                                   !
 | 
						||
        S8=S8+W(I)*(FCT(C1+U,OM,KK)+FCT(C1-U,OM,KK))                !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      S16=ZERO                                                      !
 | 
						||
      DO I = 5,12                                                   !
 | 
						||
        U=C2*X(I)                                                   !
 | 
						||
        S16=S16+W(I)*(FCT(C1+U,OM,KK)+FCT(C1-U,OM,KK))              !
 | 
						||
      END DO                                                        !
 | 
						||
      S16=C2*S16                                                    !
 | 
						||
!
 | 
						||
      IF(DABS(S16-C2*S8) <= EPS*(ONE+DABS(S16))) THEN               !
 | 
						||
        H=H+S16                                                     !
 | 
						||
        IF(BB /= B) GO TO 1                                         !
 | 
						||
      ELSE                                                          !
 | 
						||
        BB=C1                                                       !
 | 
						||
        IF(ONE+CONST*DABS(C2) /= ONE) GO TO 2                       !
 | 
						||
        H=ZERO                                                      !
 | 
						||
!
 | 
						||
        WRITE(LOGF,*)' DGAUSS: D103.1, too high accuracy required'  !
 | 
						||
        STOP                                                        !
 | 
						||
!
 | 
						||
      END IF                                                        !
 | 
						||
!
 | 
						||
   99 DGAUSS1=H                                                     !
 | 
						||
!
 | 
						||
      END FUNCTION DGAUSS1
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      SUBROUTINE INTDE(F,A,B,EPS,I,ERR)
 | 
						||
!
 | 
						||
!  This subroutine is the integrator of f(x) over (a,b)
 | 
						||
!
 | 
						||
!
 | 
						||
!  Input parameters:
 | 
						||
!
 | 
						||
!       * F        : integrand f(x)
 | 
						||
!       * A        : lower limit of integration
 | 
						||
!       * B        : upper limit of integration
 | 
						||
!       * EPS      : relative error requested
 | 
						||
!
 | 
						||
!
 | 
						||
!  Output variables:
 | 
						||
!
 | 
						||
!       * I        : approximation to the integral
 | 
						||
!       * ERR      : estimate of the absolute error
 | 
						||
!
 | 
						||
!
 | 
						||
!  Remarks:
 | 
						||
!
 | 
						||
!         function
 | 
						||
!             f(x) needs to be analytic over (a,b).
 | 
						||
!         relative error
 | 
						||
!             EPS is relative error requested excluding
 | 
						||
!             cancellation of significant digits.
 | 
						||
!             i.e. EPS means : (absolute error) /
 | 
						||
!                              (integral_a^b |f(x)| dx).
 | 
						||
!             EPS does not mean : (absolute error) / I.
 | 
						||
!         error message
 | 
						||
!             ERR >= 0 : normal termination.
 | 
						||
!             ERR < 0  : abnormal termination (M >= MMAX).
 | 
						||
!                        i.e. convergent error is detected :
 | 
						||
!                            1. f(x) or (d/dx)^n f(x) has
 | 
						||
!                               discontinuous points or sharp
 | 
						||
!                               peaks over (a,b).
 | 
						||
!                               you must divide the interval
 | 
						||
!                               (a,b) at this points.
 | 
						||
!                            2. relative error of f(x) is
 | 
						||
!                               greater than eps.
 | 
						||
!                            3. f(x) has oscillatory factor
 | 
						||
!                               and frequency of the oscillation
 | 
						||
!                               is very high.
 | 
						||
!
 | 
						||
!    Copyright(C) 1996 Takuya OOURA (email: ooura@mmm.t.u-tokyo.ac.jp).
 | 
						||
!    You may use, copy, modify this code for any purpose and
 | 
						||
!    without fee. You may distribute this ORIGINAL package.
 | 
						||
!
 | 
						||
!
 | 
						||
!   Modified: D. Sébilleau   4 Jun 2020
 | 
						||
!
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : ONE,TWO,FOUR,HALF,FOURTH
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP)             ::  A,B,EPS,I,ERR
 | 
						||
      REAL (WP)             ::  EFS,HOFF
 | 
						||
      REAL (WP)             ::  PI2,EPSLN,EPSH,H0,EHP,EHM,EPST,BA,IR,H
 | 
						||
      REAL (WP)             ::  IBACK,IRBACK,T,EP,EM,XW,XA,WG,FA,FB,ERRT
 | 
						||
      REAL (WP)             ::  ERRH,ERRD
 | 
						||
!
 | 
						||
      REAL (WP), EXTERNAL   ::  F
 | 
						||
!
 | 
						||
      REAL (WP)             ::  DATAN,DLOG,DSQRT,DEXP,DABS
 | 
						||
!
 | 
						||
      INTEGER               ::  MMAX
 | 
						||
      INTEGER               ::  M
 | 
						||
!
 | 
						||
! ---- adjustable parameter ----
 | 
						||
!
 | 
						||
      MMAX = 256                                                    !
 | 
						||
      EFS  = 0.1E0_WP                                               !
 | 
						||
      HOFF = 8.5E0_WP                                               !
 | 
						||
!
 | 
						||
! ------------------------------
 | 
						||
!
 | 
						||
      PI2   = TWO*DATAN(ONE)                                        !
 | 
						||
      EPSLN = ONE-DLOG(EFS*EPS)                                     !
 | 
						||
      EPSH  = DSQRT(EFS*EPS)                                        !
 | 
						||
      H0    = HOFF/EPSLN                                            !
 | 
						||
      EHP   = DEXP(H0)                                              !
 | 
						||
      EHM   = ONE/EHP                                               !
 | 
						||
      EPST  = DEXP(-EHM*EPSLN)                                      !
 | 
						||
      BA    = B-A                                                   !
 | 
						||
      IR    = F((A+B)*HALF)*(BA*FOURTH)                             !
 | 
						||
      I     = IR*(TWO*PI2)                                          !
 | 
						||
      ERR   = DABS(I)*EPST                                          !
 | 
						||
      H     = TWO  *H0                                              !
 | 
						||
      M     = 1                                                     !
 | 
						||
!
 | 
						||
   10 CONTINUE                                                      !
 | 
						||
!
 | 
						||
      IBACK  = I                                                    !
 | 
						||
      IRBACK = IR                                                   !
 | 
						||
      T      = H*HALF                                               !
 | 
						||
!
 | 
						||
   20 CONTINUE                                                      !
 | 
						||
!
 | 
						||
      EM = DEXP(T)                                                  !
 | 
						||
      EP = PI2*EM                                                   !
 | 
						||
      EM = PI2/EM                                                   !
 | 
						||
!
 | 
						||
   30 CONTINUE                                                      !
 | 
						||
!
 | 
						||
      XW   = ONE/(ONE+DEXP(EP-EM))                                  !
 | 
						||
      XA   = BA*XW                                                  !
 | 
						||
      WG   = XA*(ONE-XW)                                            !
 | 
						||
      FA   = F(A+XA)*WG                                             !
 | 
						||
      FB   = F(B-XA)*WG                                             !
 | 
						||
      IR   = IR+(FA+FB)                                             !
 | 
						||
      I    = I+(FA+FB)*(EP+EM)                                      !
 | 
						||
      ERRT = (DABS(FA)+DABS(FB))*(EP+EM)                            !
 | 
						||
!
 | 
						||
      IF(M == 1) ERR = ERR+ERRT*EPST                                !
 | 
						||
!
 | 
						||
      EP = EP*EHP                                                   !
 | 
						||
      EM = EM*EHM                                                   !
 | 
						||
!
 | 
						||
      IF(ERRT > ERR .OR. XW > EPSH) GO TO 30                        !
 | 
						||
!
 | 
						||
      T = T+H                                                       !
 | 
						||
!
 | 
						||
      IF(T < H0) GO TO 20                                           !
 | 
						||
!
 | 
						||
      IF(M == 1) THEN                                               !
 | 
						||
        ERRH = (ERR/EPST)*EPSH*H0                                   !
 | 
						||
        ERRD = ONE + TWO*ERRH                                       !
 | 
						||
      ELSE                                                          !
 | 
						||
        ERRD = H*(DABS(I-TWO*IBACK) + FOUR*DABS(IR-TWO*IRBACK))     !
 | 
						||
      END IF                                                        !
 | 
						||
!
 | 
						||
      H = H*HALF                                                    !
 | 
						||
      M = M*2                                                       !
 | 
						||
!
 | 
						||
      IF(ERRD > ERRH .AND. M < MMAX) GO TO 10                       !
 | 
						||
!
 | 
						||
      I = I*H                                                       !
 | 
						||
!
 | 
						||
      IF(ERRD > ERRH) THEN                                          !
 | 
						||
        ERR = -ERRD*M                                               !
 | 
						||
      ELSE                                                          !
 | 
						||
        ERR = ERRH*EPSH*M / (TWO*EFS)                               !
 | 
						||
      END IF                                                        !
 | 
						||
!
 | 
						||
      END SUBROUTINE INTDE
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      SUBROUTINE INTDE_F(F,A,B,AW,I,ERR)
 | 
						||
!
 | 
						||
!  This subroutine is the integrator of f(x) over (a,b)
 | 
						||
!
 | 
						||
!
 | 
						||
!              -->                           <--
 | 
						||
!              --> This is the fast version  <--
 | 
						||
!              -->                           <--
 | 
						||
!
 | 
						||
!
 | 
						||
!  Usage:
 | 
						||
!
 | 
						||
!         CALL INTDEINI(LENAW,TINY,EPS,AW)        ! initialization of AW
 | 
						||
!         ...
 | 
						||
!         CALL INTDE_F(F,A,B,AW,I,ERR)
 | 
						||
!
 | 
						||
!
 | 
						||
!  Input parameters:
 | 
						||
!
 | 
						||
!       * F        : integrand f(x)
 | 
						||
!       * A        : lower limit of integration
 | 
						||
!       * B        : upper limit of integration
 | 
						||
!       * AW       : points and weights of the quadrature
 | 
						||
!                      formula, AW(0...LENAW-1)
 | 
						||
!
 | 
						||
!
 | 
						||
!  Output variables:
 | 
						||
!
 | 
						||
!       * I        : approximation to the integral
 | 
						||
!       * ERR      : estimate of the absolute error
 | 
						||
!
 | 
						||
!
 | 
						||
!  Remarks:
 | 
						||
!
 | 
						||
!         initial parameters
 | 
						||
!             LENAW > 1000,
 | 
						||
!             IEEE double :
 | 
						||
!                 LENAW = 8000
 | 
						||
!                 TINY  = 1.0D-307
 | 
						||
!         function
 | 
						||
!             f(x) needs to be analytic over (a,b).
 | 
						||
!         relative error
 | 
						||
!             EPS is relative error requested excluding
 | 
						||
!             cancellation of significant digits.
 | 
						||
!             i.e. EPS means : (absolute error) /
 | 
						||
!                              (integral_a^b |f(x)| dx).
 | 
						||
!             EPS does not mean : (absolute error) / I.
 | 
						||
!         error message
 | 
						||
!             ERR >= 0 : normal termination.
 | 
						||
!             ERR < 0  : abnormal termination (M >= MMAX).
 | 
						||
!                        i.e. convergent error is detected :
 | 
						||
!                            1. f(x) or (d/dx)^n f(x) has
 | 
						||
!                               discontinuous points or sharp
 | 
						||
!                               peaks over (a,b).
 | 
						||
!                               you must divide the interval
 | 
						||
!                               (a,b) at this points.
 | 
						||
!                            2. relative error of f(x) is
 | 
						||
!                               greater than eps.
 | 
						||
!                            3. f(x) has oscillatory factor
 | 
						||
!                               and frequency of the oscillation
 | 
						||
!                               is very high.
 | 
						||
!
 | 
						||
!    Copyright(C) 1996 Takuya OOURA (email: ooura@mmm.t.u-tokyo.ac.jp).
 | 
						||
!    You may use, copy, modify this code for any purpose and
 | 
						||
!    without fee. You may distribute this ORIGINAL package.
 | 
						||
!
 | 
						||
!
 | 
						||
!   Modified: D. Sébilleau   4 Jun 2020
 | 
						||
!
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : ONE,TWO,HALF
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP)             ::  A,B,AW(0 : *),I,ERR
 | 
						||
      REAL (WP)             ::  EPSH,BA,IR,XA,FA,FB,ERRT,ERRH,ERRD,H,IBACK,IRBACK
 | 
						||
!
 | 
						||
      REAL (WP), EXTERNAL   ::  F
 | 
						||
!
 | 
						||
      REAL (WP)             ::  DABS
 | 
						||
!
 | 
						||
      INTEGER               ::  NOFF,LENAWM,NK,K,J,JTMP,JM,M,KLIM
 | 
						||
!
 | 
						||
      INTEGER               ::  INT
 | 
						||
!
 | 
						||
      NOFF   = 5                                                    !
 | 
						||
      LENAWM = INT(AW(0)+HALF)                                      !
 | 
						||
      NK     = INT(AW(1)+HALF)                                      !
 | 
						||
      EPSH   = AW(4)                                                !
 | 
						||
      BA     = B - A                                                !
 | 
						||
      I      = F((A+B) * AW(NOFF))                                  !
 | 
						||
      IR     = I * AW(NOFF+1)                                       !
 | 
						||
      I      = I * AW(NOFF+2)                                       !
 | 
						||
      ERR    = DABS(I)                                              !
 | 
						||
      K      = NK + NOFF                                            !
 | 
						||
      J      = NOFF                                                 !
 | 
						||
!
 | 
						||
   10 CONTINUE                                                      !
 | 
						||
!
 | 
						||
      J   = J + 3                                                   !
 | 
						||
      XA  = BA * AW(J)                                              !
 | 
						||
      FA  = F(A+XA)                                                 !
 | 
						||
      FB  = F(B-XA)                                                 !
 | 
						||
      IR  = IR + (FA+FB) * AW(J+1)                                  !
 | 
						||
      FA  = FA * AW(J+2)                                            !
 | 
						||
      FB  = FB * AW(J+2)                                            !
 | 
						||
      I   = I + (FA+FB)                                             !
 | 
						||
      ERR = ERR + (DABS(FA)+DABS(FB))                               !
 | 
						||
!
 | 
						||
      IF (AW(J) > EPSH .AND. J < K) GO TO 10                        !
 | 
						||
!
 | 
						||
      ERRT = ERR * AW(3)                                            !
 | 
						||
      ERRH = ERR * EPSH                                             !
 | 
						||
      ERRD = ONE + TWO*ERRH                                         !
 | 
						||
      JTMP = J                                                      !
 | 
						||
!
 | 
						||
      DO WHILE (DABS(FA) > ERRT .AND. J < K)                        !
 | 
						||
        J  = J + 3                                                  !
 | 
						||
        FA = F(A + BA*AW(J))                                        !
 | 
						||
        IR = IR + FA*AW(J+1)                                        !
 | 
						||
        FA = FA * AW(J+2)                                           !
 | 
						||
        I  = I + FA                                                 !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      JM = J                                                        !
 | 
						||
      J  = JTMP                                                     !
 | 
						||
!
 | 
						||
      DO WHILE (DABS(FB) > ERRT .AND. J < K)                        !
 | 
						||
        J  = J + 3                                                  !
 | 
						||
        FB = F(B - BA*AW(J))                                        !
 | 
						||
        IR = IR + FB*AW(J+1)                                        !
 | 
						||
        FB = FB * AW(J+2)                                           !
 | 
						||
        I  = I + FB                                                 !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      IF(J < JM) JM = J                                             !
 | 
						||
!
 | 
						||
      JM   = JM - (NOFF+3)                                          !
 | 
						||
      H    = ONE                                                    !
 | 
						||
      M    = 1                                                      !
 | 
						||
      KLIM = K + NK                                                 !
 | 
						||
!
 | 
						||
      DO WHILE (ERRD > ERRH .AND. KLIM <= LENAWM)                   !
 | 
						||
        IBACK  = I                                                  !
 | 
						||
        IRBACK = IR                                                 !
 | 
						||
!
 | 
						||
   20   CONTINUE                                                    !
 | 
						||
!
 | 
						||
        JTMP = K + JM                                               !
 | 
						||
        DO J = K + 3, JTMP, 3                                       !
 | 
						||
          XA = BA*AW(J)                                             !
 | 
						||
          FA = F(A + XA)                                            !
 | 
						||
          FB = F(B - XA)                                            !
 | 
						||
          IR = IR + (FA + FB)*AW(J+1)                               !
 | 
						||
          I  = I + (FA + FB)*AW(J+2)                                !
 | 
						||
        END DO                                                      !
 | 
						||
!
 | 
						||
        K = K + NK                                                  !
 | 
						||
        J = JTMP                                                    !
 | 
						||
!
 | 
						||
   30   CONTINUE                                                    !
 | 
						||
!
 | 
						||
        J  = J + 3                                                  !
 | 
						||
        FA = F(A + BA*AW(J))                                        !
 | 
						||
        IR = IR + FA*AW(J+1)                                        !
 | 
						||
        FA = FA * AW(J+2)                                           !
 | 
						||
        I  = I + FA                                                 !
 | 
						||
!
 | 
						||
        IF(DABS(FA) > ERRT .AND. J < K) GO TO 30                    !
 | 
						||
!
 | 
						||
        J  = JTMP                                                   !
 | 
						||
!
 | 
						||
   40   CONTINUE                                                    !
 | 
						||
!
 | 
						||
        J  = J + 3                                                  !
 | 
						||
        FB = F(B - BA*AW(J))                                        !
 | 
						||
        IR = IR + FB*AW(J+1)                                        !
 | 
						||
        FB = FB * AW(J+2)                                           !
 | 
						||
        I  = I + FB                                                 !
 | 
						||
!
 | 
						||
        IF(DABS(FB) > ERRT .AND. J < K) GO TO 40                    !
 | 
						||
!
 | 
						||
        IF(K < KLIM) GO TO 20                                       !
 | 
						||
!
 | 
						||
        ERRD = H * (DABS(I-2*IBACK) + DABS(IR-2*IRBACK))            !
 | 
						||
        H    = H * HALF                                             !
 | 
						||
        M    = M * 2                                                !
 | 
						||
        KLIM = 2*KLIM - NOFF                                        !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      I = I * (H*BA)                                                !
 | 
						||
      IF(ERRD > ERRH) THEN                                          !
 | 
						||
          ERR = -ERRD * (M * DABS(BA))                              !
 | 
						||
      ELSE                                                          !
 | 
						||
          ERR = ERR * AW(2)*(M * DABS(BA))                          !
 | 
						||
      END IF                                                        !
 | 
						||
!
 | 
						||
      END SUBROUTINE INTDE_F
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      SUBROUTINE INTDEINI_F(LENAW,TINY,EPS,AW)
 | 
						||
!
 | 
						||
!  This subroutine calculates the points and weights of the quadrature
 | 
						||
!                      formula
 | 
						||
!
 | 
						||
!    Copyright(C) 1996 Takuya OOURA (email: ooura@mmm.t.u-tokyo.ac.jp).
 | 
						||
!    You may use, copy, modify this code for any purpose and
 | 
						||
!    without fee. You may distribute this ORIGINAL package.
 | 
						||
!
 | 
						||
!
 | 
						||
!   Modified: D. Sébilleau   4 Jun 2020
 | 
						||
!
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : ZERO,ONE,TWO,FOUR,HALF
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP)             ::  TINY,EPS,AW(0 : LENAW - 1)
 | 
						||
      REAL (WP)             ::  EFS,HOFF
 | 
						||
      REAL (WP)             ::  PI2,TINYLN,EPSLN,H0,EHP,EHM,H,T,EP,EM,XW,WG
 | 
						||
!
 | 
						||
      REAL (WP)             ::  DATAN,DLOG,DEXP,DSQRT
 | 
						||
!
 | 
						||
      INTEGER               ::  LENAW
 | 
						||
      INTEGER               ::  NOFF,NK,K,J
 | 
						||
!
 | 
						||
! ---- adjustable parameter ----
 | 
						||
!
 | 
						||
      EFS  = 0.1E0_WP                                               !
 | 
						||
      HOFF = 8.5E0_WP                                               !
 | 
						||
!
 | 
						||
! ------------------------------
 | 
						||
!
 | 
						||
      PI2    = TWO * DATAN(ONE)                                     !
 | 
						||
      TINYLN = -DLOG(TINY)                                          !
 | 
						||
      EPSLN  = ZERO - DLOG(EFS*EPS)                                 !
 | 
						||
      H0     = HOFF / EPSLN                                         !
 | 
						||
      EHP    = DEXP(H0)                                             !
 | 
						||
      EHM    = ONE / EHP                                            !
 | 
						||
      AW(2)  = EPS                                                  !
 | 
						||
      AW(3)  = DEXP(-EHM*EPSLN)                                     !
 | 
						||
      AW(4)  = DSQRT(EFS*EPS)                                       !
 | 
						||
      NOFF   = 5                                                    !
 | 
						||
      AW(NOFF)   = HALF                                             !
 | 
						||
      AW(NOFF+1) = H0                                               !
 | 
						||
      AW(NOFF+2) = PI2 * H0 * HALF                                  !
 | 
						||
      H  = TWO                                                      !
 | 
						||
      NK = 0                                                        !
 | 
						||
      K  = NOFF + 3                                                 !
 | 
						||
!
 | 
						||
   10 CONTINUE                                                      !
 | 
						||
!
 | 
						||
      T = H * HALF                                                  !
 | 
						||
!
 | 
						||
   20 CONTINUE                                                      !
 | 
						||
!
 | 
						||
      EM = DEXP(H0*T)                                               !
 | 
						||
      EP = PI2 * EM                                                 !
 | 
						||
      EM = PI2 / EM                                                 !
 | 
						||
      J  = K                                                        !
 | 
						||
!
 | 
						||
   30 CONTINUE                                                      !
 | 
						||
!
 | 
						||
      XW      = ONE / (ONE + DEXP(EP-EM))                           !
 | 
						||
      WG      = XW * (ONE-XW) * H0                                  !
 | 
						||
      AW(J)   = XW                                                  !
 | 
						||
      AW(J+1) = WG * FOUR                                           !
 | 
						||
      AW(J+2) = WG * (EP+EM)                                        !
 | 
						||
      EP      = EP * EHP                                            !
 | 
						||
      EM      = EM * EHM                                            !
 | 
						||
      J       = J + 3                                               !
 | 
						||
!
 | 
						||
      IF (EP < TINYLN .AND. J <= (LENAW-3)) GO TO 30                !
 | 
						||
!
 | 
						||
      T = T + H                                                     !
 | 
						||
      K = K + NK                                                    !
 | 
						||
!
 | 
						||
      IF(T < ONE) GO TO 20                                          !
 | 
						||
!
 | 
						||
      H = H * HALF                                                  !
 | 
						||
!
 | 
						||
      IF(NK == 0) THEN                                              !
 | 
						||
        IF(J > (LENAW-6)) J = J - 3                                 !
 | 
						||
        NK    = J - NOFF                                            !
 | 
						||
        K     = K + NK                                              !
 | 
						||
        AW(1) = NK                                                  !
 | 
						||
      END IF                                                        !
 | 
						||
!
 | 
						||
      IF((2*K - NOFF - 3) <= LENAW) GO TO 10                        !
 | 
						||
!
 | 
						||
      AW(0) = DFLOAT(K-3)                                           !
 | 
						||
!
 | 
						||
      END SUBROUTINE INTDEINI_F
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      FUNCTION RBI1(FCT,A,B,PREC,OBTPREC,NITER,ITERMIN,ITERMAX)
 | 
						||
!
 | 
						||
!*******************************************************
 | 
						||
!* Integral of a function FCT(X) by Romberg's method   *
 | 
						||
!* --------------------------------------------------- *
 | 
						||
!* INPUTS:                                             *
 | 
						||
!*          A      begin value of x variable           *
 | 
						||
!*          B      end value of x variable             *
 | 
						||
!*       PREC      desired precision                   *
 | 
						||
!*    ITERMIN      minimum number of iterations        *
 | 
						||
!*    ITERMAX      maximum number of iterations        *
 | 
						||
!*                                                     *
 | 
						||
!* OUTPUTS:                                            *
 | 
						||
!*    OBTPREC      obtained precision for integral     *
 | 
						||
!*      NITER      number of iterations done           *
 | 
						||
!*   INTEGRAL      the integral of FCT(X) from a to b  *
 | 
						||
!*                                                     *
 | 
						||
!*******************************************************
 | 
						||
!
 | 
						||
!   Last modified: D. Sébilleau   5 June 2020
 | 
						||
!
 | 
						||
!
 | 
						||
      USE DIMENSION_CODE,   ONLY : MAXITER
 | 
						||
      USE REAL_NUMBERS,     ONLY : ONE,TWO,FOUR
 | 
						||
!
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP)             ::  RBI1
 | 
						||
      REAL (WP)             ::  A,B,PREC,OBTPREC
 | 
						||
      REAL (WP)             ::  T(0:MAXITER,0:MAXITER)
 | 
						||
      REAL (WP)             ::  PAS,R,S,TA
 | 
						||
!
 | 
						||
      REAL (WP), EXTERNAL   ::  FCT
 | 
						||
!
 | 
						||
      REAL (WP)             ::  DABS
 | 
						||
!
 | 
						||
      INTEGER               ::  NITER,ITERMIN,ITERMAX,I,J
 | 
						||
!
 | 
						||
      IF (ITERMAX > MAXITER)  ITERMAX=MAXITER                       !
 | 
						||
!
 | 
						||
      R = FCT(A)                                                    !
 | 
						||
      TA = (R + FCT(B) ) / TWO                                      !
 | 
						||
      NITER=0                                                       !
 | 
						||
      PAS=B-A                                                       !
 | 
						||
      T(0,0)=TA*PAS                                                 !
 | 
						||
 100  NITER=NITER+1                                                 !
 | 
						||
      PAS=PAS/TWO                                                   !
 | 
						||
      S=TA                                                          !
 | 
						||
!
 | 
						||
      DO I=1, 2**NITER-1                                            !
 | 
						||
        S = S + FCT(A+PAS*I)                                        !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      T(0,NITER)=S*PAS                                              !
 | 
						||
      R=ONE                                                         !
 | 
						||
      DO I=1, NITER                                                 !
 | 
						||
        R=R*FOUR                                                    !
 | 
						||
        J=NITER-I                                                   !
 | 
						||
        T(I,J)=(R*T(I-1,J+1) - T(I-1,J))/(R-ONE)                    !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      OBTPREC = DABS(T(NITER,0) - T(NITER-1,0))                     !
 | 
						||
!
 | 
						||
      IF (NITER   > ITERMAX)  GO TO 200                             !
 | 
						||
      IF (NITER   < ITERMIN)  GO TO 100                             !
 | 
						||
      IF (OBTPREC > PREC)     GO TO 100                             !
 | 
						||
!
 | 
						||
 200  RBI1 = T(NITER,0)                                             !
 | 
						||
!
 | 
						||
      END FUNCTION RBI1
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      SUBROUTINE SIMPSON(FCT,A,B,N,RES)
 | 
						||
!
 | 
						||
!*******************************************************
 | 
						||
!* Integral of a function FCT(X) by Simpson's method   *
 | 
						||
!* --------------------------------------------------- *
 | 
						||
!* INPUTS:                                             *
 | 
						||
!*          A      begin value of x variable           *
 | 
						||
!*          B      end value of x variable             *
 | 
						||
!*          N      number of integration steps         *
 | 
						||
!*                                                     *
 | 
						||
!* OUTPUT:                                             *
 | 
						||
!*          RES    the integral of FCT(X) from a to b  *
 | 
						||
!*                                                     *
 | 
						||
!*******************************************************
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : TWO,THREE
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP)             ::  A,B,RES
 | 
						||
      REAL (WP)             ::  STEP,R
 | 
						||
!
 | 
						||
      REAL (WP), EXTERNAL   ::  FCT
 | 
						||
!
 | 
						||
      INTEGER               ::   N,I
 | 
						||
!
 | 
						||
      STEP = (B-A)/TWO/N                                            !
 | 
						||
      R = FCT(A)                                                    !
 | 
						||
      RES = (R+FCT(B))/TWO                                          !
 | 
						||
!
 | 
						||
      DO I=1, 2*N-1                                                 !
 | 
						||
        R = FCT(A+I*STEP)                                           !
 | 
						||
        IF(MOD(I,2) /= 0) THEN                                      !
 | 
						||
          RES = RES + R + R                                         !
 | 
						||
        ELSE                                                        !
 | 
						||
          RES = RES + R                                             !
 | 
						||
        END IF                                                      !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      RES = RES * STEP*TWO/THREE                                    !
 | 
						||
!
 | 
						||
      END SUBROUTINE SIMPSON
 | 
						||
!
 | 
						||
!=======================================================================
 | 
						||
!
 | 
						||
      SUBROUTINE QANC8 (FCT,A,B,AERR,RERR,RES,ERR,NBF,FLG)
 | 
						||
!
 | 
						||
!     Integrate a real function FCT(X) from X = A to X = B,
 | 
						||
!     with given absolute and relative precisions, AERR, RERR.
 | 
						||
!
 | 
						||
!     Inputs:
 | 
						||
!
 | 
						||
!     FCT       :  external user-defined function for any X value
 | 
						||
!                  in interval [A,B]
 | 
						||
!     A,B       :  limits of interval
 | 
						||
!     AERR,RERR :  respectively absolute error and relative error
 | 
						||
!                  required by user
 | 
						||
!
 | 
						||
!     Outputs:
 | 
						||
!
 | 
						||
!     RES       :  value of integral
 | 
						||
!     ERR       :  estimated error
 | 
						||
!     NBF       :  number of necessary FCT(X) evaluations
 | 
						||
!     FLG       :  indicator
 | 
						||
!                    = 0.0       correct result
 | 
						||
!                    = NNN.RRR   no convergence du to a singularity
 | 
						||
!                      the singular point abcissa is given by formula:
 | 
						||
!                                 XS = B-.RRR*(B-A)
 | 
						||
 | 
						||
!     Reference :
 | 
						||
!
 | 
						||
!     G.E. Forsythe, Computer Methods for Mathematical
 | 
						||
!                    Computations, Prentice-Hall, Inc. (1977)
 | 
						||
!
 | 
						||
! -----------------------------------------------------------------------
 | 
						||
!
 | 
						||
      USE REAL_NUMBERS,     ONLY : ZERO,ONE,TWO,HALF
 | 
						||
!
 | 
						||
      IMPLICIT NONE
 | 
						||
!
 | 
						||
      REAL (WP), INTENT(IN)   ::    A,B,AERR,RERR
 | 
						||
      REAL (WP), INTENT(OUT)  ::    RES,FLG
 | 
						||
      REAL (WP)               ::    ERR
 | 
						||
      REAL (WP)               ::    QR(31),F(16),X(16),FS(8,30),XS(8,30)
 | 
						||
      REAL (WP)               ::    W0,W1,W2,W3,W4
 | 
						||
      REAL (WP)               ::    COR,SUM
 | 
						||
      REAL (WP)               ::    X0,QP,PAS1,PAS,QL,QN,QD,ERR1,TOL1
 | 
						||
      REAL (WP)               ::    F0,TEMP
 | 
						||
      REAL (WP)               ::    DABS,MAX
 | 
						||
!
 | 
						||
      REAL (WP), EXTERNAL     ::    FCT
 | 
						||
!
 | 
						||
      INTEGER, INTENT(OUT)    ::    NBF
 | 
						||
      INTEGER                 ::    LMIN,LMAX,LOUT,NMAX,NFIN
 | 
						||
      INTEGER                 ::    L,NIM,J,I
 | 
						||
!
 | 
						||
      LMIN = 1                                                      !
 | 
						||
      LMAX = 30                                                     !
 | 
						||
      LOUT = 6                                                      !
 | 
						||
      NMAX = 5000                                                   !
 | 
						||
      NFIN = NMAX-8*(LMAX-LOUT+2**(LOUT+1))                         !
 | 
						||
      W0   =   3956.E0_WP/14175.E0_WP                               !
 | 
						||
      W1   =  23552.E0_WP/14175.E0_WP                               !
 | 
						||
      W2   =  -3712.E0_WP/14175.E0_WP                               !
 | 
						||
      W3   =  41984.E0_WP/14175.E0_WP                               !
 | 
						||
      W4   = -18160.E0_WP/14175.E0_WP                               !
 | 
						||
      FLG  = ZERO                                                   !
 | 
						||
      RES  = ZERO                                                   !
 | 
						||
      COR  = ZERO                                                   !
 | 
						||
      ERR  = ZERO                                                   !
 | 
						||
      SUM  = ZERO                                                   !
 | 
						||
      NBF  = 0                                                      !
 | 
						||
!
 | 
						||
      IF (A == B) RETURN                                            !
 | 
						||
!
 | 
						||
      L     = 0                                                     !
 | 
						||
      NIM   = 1                                                     !
 | 
						||
      X0    = A                                                     !
 | 
						||
      X(16) = B                                                     !
 | 
						||
      QP    = ZERO                                                  !
 | 
						||
      F0    = FCT(X0)                                               !
 | 
						||
      PAS1  = (B-A)/16.E0_WP                                        !
 | 
						||
      X(8)  = (X0+X(16))*HALF                                       !
 | 
						||
      X(4)  = (X0+X(8))*HALF                                        !
 | 
						||
      X(12) = (X(8)+X(16))*HALF                                     !
 | 
						||
      X(2)  = (X0+X(4))*HALF                                        !
 | 
						||
      X(6)  = (X(4)+X(8))*HALF                                      !
 | 
						||
      X(10) = (X(8)+X(12))*HALF                                     !
 | 
						||
      X(14) = (X(12)+X(16))*HALF                                    !
 | 
						||
!
 | 
						||
      DO J = 2,16,2                                                 !
 | 
						||
        F(J) = FCT(X(J))                                            !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      NBF  = 9                                                      !
 | 
						||
   30 X(1) = (X0+X(2))*HALF                                         !
 | 
						||
      F(1) = FCT(X(1))                                              !
 | 
						||
!
 | 
						||
      DO J = 3,15,2                                                 !
 | 
						||
        X(J) = (X(J-1)+X(J+1))*HALF                                 !
 | 
						||
        F(J) = FCT(X(J))                                            !
 | 
						||
      END DO
 | 
						||
!
 | 
						||
      NBF  = NBF+8                                                  !
 | 
						||
      PAS  = (X(16)-X0)/16.E0_WP                                    !
 | 
						||
      QL   = (W0*(F0+F(8))+W1*(F(1)+F(7))+W2*(F(2)+F(6))        &   !
 | 
						||
                                 +W3*(F(3)+F(5))+W4*F(4))*PAS       !
 | 
						||
      QR(L+1) = (W0*(F(8)+F(16))+W1*(F(9)+F(15))                &   !
 | 
						||
             +W2*(F(10)+F(14))+W3*(F(11)+F(13))+W4*F(12))*PAS       !
 | 
						||
      QN   = QL + QR(L+1)                                           !
 | 
						||
      QD   = QN - QP                                                !
 | 
						||
      SUM  = SUM + QD                                               !
 | 
						||
      ERR1 = DABS(QD)/1023.E0_WP                                    !
 | 
						||
      TOL1 = MAX(AERR,RERR*DABS(SUM))*(PAS/PAS1)                    !
 | 
						||
!
 | 
						||
      IF (L    <   LMIN)    GO TO 50                                !
 | 
						||
      IF (L    >=  LMAX)    GO TO 62                                !
 | 
						||
      IF (NBF  >   NFIN)    GO TO 60                                !
 | 
						||
      IF (ERR1 <=  TOL1)    GO TO 70                                !
 | 
						||
!
 | 
						||
   50 NIM = 2*NIM                                                   !
 | 
						||
      L   = L+1                                                     !
 | 
						||
!
 | 
						||
      DO I = 1,8                                                    !
 | 
						||
        FS(I,L) = F(I+8)                                            !
 | 
						||
        XS(I,L) = X(I+8)                                            !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      QP = QL                                                       !
 | 
						||
!
 | 
						||
      DO I = 1,8                                                    !
 | 
						||
        F(18-2*I) = F(9-I)                                          !
 | 
						||
        X(18-2*I) = X(9-I)                                          !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      GO TO 30                                                      !
 | 
						||
!
 | 
						||
   60 NFIN = 2*NFIN                                                 !
 | 
						||
      LMAX = LOUT                                                   !
 | 
						||
      FLG  = FLG + (B-X0)/(B-A)                                     !
 | 
						||
!
 | 
						||
      GO TO 70                                                      !
 | 
						||
!
 | 
						||
   62 FLG = FLG + ONE                                               !
 | 
						||
   70 RES = RES + QN                                                !
 | 
						||
      ERR = ERR + ERR1                                              !
 | 
						||
      COR = COR + QD/1023.E0_WP                                     !
 | 
						||
!
 | 
						||
   72 IF (NIM == 2*(NIM/2)) GO TO 75                                !
 | 
						||
      NIM = NIM/2                                                   !
 | 
						||
      L   = L-1                                                     !
 | 
						||
!
 | 
						||
      GO TO 72                                                      !
 | 
						||
!
 | 
						||
   75 NIM = NIM+1                                                   !
 | 
						||
      IF (L <= 0) GO TO 80                                          !
 | 
						||
      QP  = QR(L)                                                   !
 | 
						||
      X0  = X(16)                                                   !
 | 
						||
      F0  = F(16)                                                   !
 | 
						||
!
 | 
						||
      DO I = 1,8                                                    !
 | 
						||
        F(2*I) = FS(I,L)                                            !
 | 
						||
        X(2*I) = XS(I,L)                                            !
 | 
						||
      END DO                                                        !
 | 
						||
!
 | 
						||
      GO TO 30                                                      !
 | 
						||
!
 | 
						||
   80 RES = RES + COR                                               !
 | 
						||
      IF (ERR == ZERO)   RETURN                                     !
 | 
						||
   82 TEMP = DABS(RES) + ERR                                        !
 | 
						||
      IF (TEMP /= DABS(RES)) RETURN                                 !
 | 
						||
      ERR = TWO*ERR                                                 !
 | 
						||
!
 | 
						||
      GO TO 82                                                      !
 | 
						||
!
 | 
						||
      END SUBROUTINE QANC8
 | 
						||
!
 | 
						||
END MODULE INTEGRATION
 |