MsSpec-DFM/New_libraries/DFM_library/UTILITIES_LIBRARY/integration2.f90

878 lines
37 KiB
Fortran
Raw Normal View History

2022-02-02 16:19:10 +01:00
!
!=======================================================================
!
MODULE INTEGRATION2
!
! This module contains integration routines in order to integrate
! a function F over the interval [0,+INF].
!
! These routines are:
!
!
! * Lagrange : INTEGR_L_0_INF(F,X1,NSIZE,NMAX1,ID,F_INF,HOP,A)
!
! * double exponential transformation : INTDEI(F,A,EPS,I,ERR) <-- standard version
! INTDEI_1(F,A,X,EPS,I,ERR) <-- 1-parameter version (X)
! INTDEI_2(F,A,X,Y,EPS,I,ERR) <-- 2-parameter version (X,Y)
! INTDEI_F(F,A,AW,I,ERR) <-- fast version
!
USE ACCURACY_REAL
!
CONTAINS
!
!=======================================================================
!
SUBROUTINE INTEGR_L_0_INF(F,X1,N_SIZE,NMAX1,ID,F_INF,HOP,A)
!
! This subroutine integrates from 0 to infinity a function f(x)
! defined by the arrays X and F up to NMAX and for which
! the asymptotic value at infinity F_INF is known.
!
!
! Input parameters:
!
! * F : function array
! * X1 : mesh array (constant step)
! * N_SIZE : dimensioning of the arrays
! * NMAX1 : index of upper limit of integration on the mesh
! * ID : integer parameter
! ID = 1 --> F0 = 0 at the origin
! ID > 1 --> F0 not 0 at the origin
! * F_INF : limit of f(x) for x --> infinity
! * HOP : hopping parameter defining the mesh X(NMAX) --> infinity
!
!
! Output parameters:
!
! * A : integral result
!
!
! Method: The integral is separated as
!
! / X1(NMAX) / + INF = X1(NMAX+NMAX/HOP)
! | |
! | f(x) dx + | f(x) dx
! | |
! / 0 / X1(NMAX)
!
!
! * The first integral is computed using Lagrange integration.
! * For the second integral, a new mesh containing NMAX/HOP points
! is constructed. This is done by a 5-point Lagrange interpolation
! from X1(NMAX-3*HOP), X1(NMAX-2*HOP), X1(NMAX-HOP), X1(NMAX)
! and X1(NMAX+NMAX/HOP) whose F value is taken as F_INF
! * Once the new F and X1 arrays are constructed, the second integral
! is also computed using a Lagrange integration.
!
!
! --> Real function F case <--
!
!
! Author : D. Sébilleau
!
! Last modified : 5 Aug 2020
!
!
USE DIMENSION_CODE, ONLY : NSIZE
USE REAL_NUMBERS, ONLY : ZERO,ONE,TWO,THREE
USE INTEGRATION, ONLY : INTEGR_L
USE INTERPOLATION, ONLY : LAG_5P_INTERP
!
IMPLICIT NONE
!
REAL (WP) :: X1(N_SIZE),F(N_SIZE) ! 1st integral
REAL (WP) :: X2(NSIZE),G(NSIZE) ! 2nd integral
REAL (WP) :: XX(5),AA(5) ! interpolation points
REAL (WP) :: F_INF
REAL (WP) :: H1,H2
REAL (WP) :: A1,A2,A
!
REAL (WP) :: DFLOAT
!
INTEGER :: N_SIZE,NMAX1,ID,HOP
INTEGER :: NMAX2,IP
!
H1=X1(2)-X1(1) ! step for 1st integral
!
! Computing the first integral
!
CALL INTEGR_L(F,H1,N_SIZE,NMAX1,A1,ID) !
!
NMAX2=NMAX1/HOP ! nb of points of mesh X2
H2=H1*HOP ! step for mesh X2
!
! Defining the [ X(NMAX, + INF ] mesh
!
DO IP=1,NMAX2 !
X2(IP)=X1(NMAX1)+DFLOAT(IP-1)*H2 !
END DO !
!
! Constructing the Lagrange interpolation points
!
XX(1)=X1(NMAX1-3*HOP) !
XX(2)=X1(NMAX1-2*HOP) !
XX(3)=X1(NMAX1-HOP) !
XX(4)=X1(NMAX1) !
XX(5)=X2(NMAX2) !
!
AA(1)=F(NMAX1-3*HOP) !
AA(2)=F(NMAX1-2*HOP) !
AA(3)=F(NMAX1-HOP) !
AA(4)=F(NMAX1) !
AA(5)=F_INF !
!
! Evaluating f(x) = G over the X2 mesh
!
DO IP=1,NMAX2 !
G(IP)=LAG_5P_INTERP(XX,AA,X2(IP)) !
END DO !
!
! Computing the second integral
!
CALL INTEGR_L(G,H2,NSIZE,NMAX2,A2,1) !
!
A=A1+A2 !
!
END SUBROUTINE INTEGR_L_0_INF
!
!=======================================================================
!
SUBROUTINE INTDEI(F,A,EPS,I,ERR)
!
! This subroutine is the integrator of f(x) over (a,infinity)
!
!
! Input parameters:
!
! * F : integrand f(x)
! * A : lower 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,infinity).
! relative error
! EPS is relative error requested excluding
! cancellation of significant digits.
! i.e. EPS means : (absolute error) /
! (integral_a^infinity |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,infinity).
! you must divide the interval
! (a,infinity) at this points.
! 2. relative error of f(x) is
! greater than EPS.
! 3. f(x) has oscillatory factor
! and decay of f(x) is very slow
! as x -> infinity.
! 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 Aug 2020
!
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,FOUR,HALF
!
IMPLICIT NONE
!
REAL (WP) :: A,EPS,I,ERR
REAL (WP) :: EFS,HOFF
REAL (WP) :: PI4,EPSLN,EPSH,H0,EHP,EHM,EPST,IR,H,IBACK
REAL (WP) :: IRBACK,T,EP,EM,XP,XM,FP,FM,ERRT,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 = 11.0E0_WP !
!
! ------------------------------
!
PI4 = DATAN(ONE) !
EPSLN = ONE - DLOG(EFS*EPS) !
EPSH = DSQRT(EFS*EPS) !
H0 = HOFF /EPSLN !
EHP = DEXP(H0) !
EHM = ONE/EHP !
EPST = DEXP(-EHM*EPSLN) !
IR = F(A+1) !
I = IR*(TWO*PI4) !
ERR = DABS(I)*EPST !
H = TWO*H0 !
M = 1 !
!
10 CONTINUE !
!
IBACK = I !
IRBACK = IR !
T = H * HALF !
!
20 CONTINUE !
!
EM = DEXP(T) !
EP = PI4 * EM !
EM = PI4 / EM !
!
30 CONTINUE !
!
XP = DEXP(EP-EM) !
XM = ONE/XP !
FP = F(A+XP)*XP !
FM = F(A+XM)*XM !
IR = IR+(FP+FM) !
I = I+(FP+FM)*(EP+EM) !
ERRT = (DABS(FP)+DABS(FM))*(EP+EM) !
!
IF(M == 1) ERR = ERR+ERRT*EPST !
!
EP = EP*EHP !
EM = EM*EHM !
!
IF(ERRT > ERR .OR. XM > 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 INTDEI
!
!=======================================================================
!
SUBROUTINE INTDEI_1(F,A,X,EPS,I,ERR)
!
! This subroutine is the integrator of f(x) over (a,infinity)
! with ONE extra parameter in the definition of the integrand
!
!
! Input parameters:
!
! * F : integrand f(x)
! * A : lower limit of integration
! * X : extra parameter of F
! * 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,infinity).
! relative error
! EPS is relative error requested excluding
! cancellation of significant digits.
! i.e. EPS means : (absolute error) /
! (integral_a^infinity |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,infinity).
! you must divide the interval
! (a,infinity) at this points.
! 2. relative error of f(x) is
! greater than EPS.
! 3. f(x) has oscillatory factor
! and decay of f(x) is very slow
! as x -> infinity.
! 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.
!
!
! Last modified: D. Sébilleau 4 Aug 2020
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,FOUR,HALF
!
IMPLICIT NONE
!
REAL (WP) :: A,X,EPS,I,ERR
REAL (WP) :: EFS,HOFF
REAL (WP) :: PI4,EPSLN,EPSH,H0,EHP,EHM,EPST,IR,H,IBACK
REAL (WP) :: IRBACK,T,EP,EM,XP,XM,FP,FM,ERRT,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 = 11.0E0_WP !
!
! ------------------------------
!
PI4 = DATAN(ONE) !
EPSLN = ONE - DLOG(EFS*EPS) !
EPSH = DSQRT(EFS*EPS) !
H0 = HOFF /EPSLN !
EHP = DEXP(H0) !
EHM = ONE/EHP !
EPST = DEXP(-EHM*EPSLN) !
IR = F(A+1,X) !
I = IR*(TWO*PI4) !
ERR = DABS(I)*EPST !
H = TWO*H0 !
M = 1 !
!
10 CONTINUE !
!
IBACK = I !
IRBACK = IR !
T = H * HALF !
!
20 CONTINUE !
!
EM = DEXP(T) !
EP = PI4 * EM !
EM = PI4 / EM !
!
30 CONTINUE !
!
XP = DEXP(EP-EM) !
XM = ONE/XP !
FP = F(A+XP,X)*XP !
FM = F(A+XM,X)*XM !
IR = IR+(FP+FM) !
I = I+(FP+FM)*(EP+EM) !
ERRT = (DABS(FP)+DABS(FM))*(EP+EM) !
!
IF(M == 1) ERR = ERR+ERRT*EPST !
!
EP = EP*EHP !
EM = EM*EHM !
!
IF(ERRT > ERR .OR. XM > 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 INTDEI_1
!
!=======================================================================
!
SUBROUTINE INTDEI_2(F,A,X,Y,EPS,I,ERR)
!
! This subroutine is the integrator of f(x) over (a,infinity)
! with TWO extra parameters in the definition of the integrand
!
!
! Input parameters:
!
! * F : integrand f(x)
! * A : lower limit of integration
! * X : first extra parameter of F
! * Y : second extra parameter of F
! * 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,infinity).
! relative error
! EPS is relative error requested excluding
! cancellation of significant digits.
! i.e. EPS means : (absolute error) /
! (integral_a^infinity |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,infinity).
! you must divide the interval
! (a,infinity) at this points.
! 2. relative error of f(x) is
! greater than EPS.
! 3. f(x) has oscillatory factor
! and decay of f(x) is very slow
! as x -> infinity.
! 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.
!
!
! Last modified: D. Sébilleau 4 Aug 2020
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,FOUR,HALF
!
IMPLICIT NONE
!
REAL (WP) :: A,X,Y,EPS,I,ERR
REAL (WP) :: EFS,HOFF
REAL (WP) :: PI4,EPSLN,EPSH,H0,EHP,EHM,EPST,IR,H,IBACK
REAL (WP) :: IRBACK,T,EP,EM,XP,XM,FP,FM,ERRT,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 = 11.0E0_WP !
!
! ------------------------------
!
PI4 = DATAN(ONE) !
EPSLN = ONE - DLOG(EFS*EPS) !
EPSH = DSQRT(EFS*EPS) !
H0 = HOFF /EPSLN !
EHP = DEXP(H0) !
EHM = ONE/EHP !
EPST = DEXP(-EHM*EPSLN) !
IR = F(A+1,X,Y) !
I = IR*(TWO*PI4) !
ERR = DABS(I)*EPST !
H = TWO*H0 !
M = 1 !
!
10 CONTINUE !
!
IBACK = I !
IRBACK = IR !
T = H * HALF !
!
20 CONTINUE !
!
EM = DEXP(T) !
EP = PI4 * EM !
EM = PI4 / EM !
!
30 CONTINUE !
!
XP = DEXP(EP-EM) !
XM = ONE/XP !
FP = F(A+XP,X,Y)*XP !
FM = F(A+XM,X,Y)*XM !
IR = IR+(FP+FM) !
I = I+(FP+FM)*(EP+EM) !
ERRT = (DABS(FP)+DABS(FM))*(EP+EM) !
!
IF(M == 1) ERR = ERR+ERRT*EPST !
!
EP = EP*EHP !
EM = EM*EHM !
!
IF(ERRT > ERR .OR. XM > 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 INTDEI_2
!
!=======================================================================
!
SUBROUTINE INTDEI_F(F,A,AW,I,ERR)
!
! This subroutine is the integrator of f(x) over (a,infinity)
!
!
!
! --> <--
! --> This is the fast version <--
! --> <--
!
!
! Usage:
!
! CALL INTDEIINI(LENAW,TINY,EPS,AW) ! initialization of AW
! ...
! CALL INTDEI_F(F,A,AW,I,ERR)
!
!
! Input parameters:
!
! * F : integrand f(x)
! * A : lower 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,infinity).
! relative error
! EPS is relative error requested excluding
! cancellation of significant digits.
! i.e. EPS means : (absolute error) /
! (integral_a^infinity |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,infinity).
! you must divide the interval
! (a,infinity) at this points.
! 2. relative error of f(x) is
! greater than EPS.
! 3. f(x) has oscillatory factor
! and decay of f(x) is very slow
! as x -> infinity.
! 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 Aug 2020
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,HALF
!
IMPLICIT NONE
!
REAL (WP) :: A,AW(0 : *),I,ERR
REAL (WP) :: EPSH,IR,FP,FM,ERRT,ERRH,ERRD
REAL (WP) :: H,IBACK,IRBACK
!
REAL (WP), EXTERNAL :: F
!
REAL (WP) :: DABS
!
INTEGER :: NOFF,LENAWM,NK,K,J,JTMP,JM,M,KL,KLIM
!
INTEGER :: INT
!
NOFF = 5 !
LENAWM = INT(AW(0) + HALF) !
NK = INT(AW(1) + HALF) !
EPSH = AW(4) !
I = F(A + AW(NOFF)) !
IR = I * AW(NOFF+1) !
I = I * AW(NOFF+2) !
ERR = DABS(I) !
K = NK + NOFF !
J = NOFF !
!
10 CONTINUE !
!
J = J + 6 !
FM = F(A + AW(J)) !
FP = F(A + AW(J+1)) !
IR = IR + (FM*AW(J+2) + FP*AW(J+3)) !
FM = FM * AW(J+4) !
FP = FP * AW(J+5) !
I = I + (FM+FP) !
ERR = ERR + (DABS(FM)+DABS(FP)) !
!
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(FM) > ERRT .AND. J < K) !
J = J + 6 !
FM = F(A + AW(J)) !
IR = IR + FM*AW(J+2) !
FM = FM * AW(J+4) !
I = I + FM !
END DO !
!
JM = J !
J = JTMP !
!
DO WHILE (DABS(FP) > ERRT .AND. J < K) !
J = J + 6 !
FP = F(A + AW(J+1)) !
IR = IR + FP*AW(J+3) !
FP = FP * AW(J+5) !
I = I + FP !
END DO !
!
IF(J < JM) JM = J !
!
JM = JM - (NOFF+6) !
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+6,JTMP,6 !
FM = F(A + AW(J)) !
FP = F(A + AW(J+1)) !
IR = IR + (FM*AW(J+2) + FP*AW(J+3)) !
I = I + (FM*AW( +4) + FP*AW(J+5)) !
END DO !
!
K = K + NK !
J = JTMP !
!
30 CONTINUE !
!
J = J + 6 !
FM = F(A + AW(J)) !
IR = IR + FM*AW(J+2) !
FM = FM * AW(J+4) !
I = I + FM !
!
IF(DABS(FM) > ERRT .AND. J < K) GO TO 30 !
!
J = JTMP !
!
40 CONTINUE !
!
J = J + 6 !
FP = F(A + AW(J+1)) !
IR = IR + FP*AW(J+3) !
FP = FP * AW(J+5) !
I = I + FP !
!
IF(DABS(FP) > 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 !
!
IF(ERRD > ERRH) THEN !
ERR = -ERRD * M !
ELSE !
ERR = ERR * (AW(2)*M) !
END IF !
!
END SUBROUTINE INTDEI_F
!
!=======================================================================
!
SUBROUTINE INTDEIINI_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 5 Aug 2020
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,FOUR,HALF
!
IMPLICIT NONE
!
REAL (WP) :: TINY,EPS,AW(0 : LENAW - 1)
REAL (WP) :: EFS,HOFF
REAL (WP) :: PI4,TINYLN,EPSLN,H0,EHP,EHM
REAL (WP) :: H,T,EP,EM,XP,XM,WWP,WWM
!
REAL (WP) :: DATAN,DLOG,DSQRT,DEXP,DFLOAT
!
INTEGER :: LENAW
INTEGER :: NOFF,NK,K,J
!
! ---- adjustable parameter ----
!
EFS = 0.1E0_WP !
HOFF = 11.0E0_WP !
!
! ------------------------------
!
PI4 = DATAN(ONE) !
TINYLN = -DLOG(TINY) !
EPSLN = ONE - 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) = ONE !
AW(NOFF+1) = FOUR * H0 !
AW(NOFF+2) = TWO * PI4 * H0 !
H = TWO !
NK = 0 !
K = NOFF + 6 !
!
10 CONTINUE !
!
T = H * HALF !
!
20 CONTINUE !
!
EM = DEXP(H0*T) !
EP = PI4 * EM !
EM = PI4 / EM !
J = K !
!
30 CONTINUE !
!
XP = DEXP(EP-EM) !
XM = ONE / XP !
WWP= XP * ((EP+EM)*H0) !
WWM= XM * ((EP+EM)*H0) !
AW(J) = XM !
AW(J+1) = XP !
AW(J+2) = XM * (FOUR*H0) !
AW(J+3) = XP * (FOUR*H0) !
AW(J+4) = WWM !
AW(J+5) = WWP !
EP = EP * EHP !
EM = EM * EHM !
J = J + 6 !
!
IF(EP < TINYLN .AND. J <= (LENAW-6)) 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-12)) J = J - 6 !
NK = J - NOFF !
K = K + NK !
AW(1) = NK !
END IF !
!
IF((2*K - NOFF - 6) <= LENAW) GO TO 10 !
!
AW(0) = DFLOAT(K-6) !
!
END SUBROUTINE INTDEIINI_F
!
END MODULE INTEGRATION2