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

1045 lines
49 KiB
Fortran

!
!=======================================================================
!
MODULE INTEGRATION3
!
! This module contains integration routines in order to integrate
! a function F over the interval [0,+INF] when f(x) has
! an oscillatory behaviour.
!
! These routines are:
!
! * double exponential transformation : INTDEO(F,A,OMEGA,EPS,I,ERR) <-- standard version
! INTDEO_1(F,AA,A,OMEGA,EPS,I,ERR) <-- 1-parameter version (AA)
! INTDEO_2(F,AA,LL,A,OMEGA,EPS,I,ERR) <-- 2-parameter version (AA,LL)
! INTDEO_F(F,A,OMEGA,AW,I,ERR) <-- fast version ^
! |
USE ACCURACY_REAL !
! integer
CONTAINS
!
!=======================================================================
!
SUBROUTINE INTDEO(F,A,OMEGA,EPS,I,ERR)
!
! This subroutine is the integrator of f(x) over (a,infinity) when
! f(x) has an oscillatory behaviour:
!
! f(x) = g(x) * sin(omega * x + theta) as x -> infinity
!
!
! Input parameters:
!
! * F : integrand f(x)
! * A : lower limit of integration
! * OMEGA : frequency of oscillation
! * 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^R |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.
!
! 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) :: A,OMEGA,EPS,I,ERR
REAL (WP) :: EFS,ENOFF,PQOFF,PPOFF
REAL (WP) :: PI4,EPSLN,EPSH,FRQ4,PER2,PP,PQ,EHP,EHM,IR,H
REAL (WP) :: IBACK,IRBACK,T,EP,EM,TK,XW,WG,XA,FP,FM,ERRH
REAL (WP) :: TN,ERRD
!
REAL (WP), EXTERNAL :: F
!
REAL (WP) :: DATAN,DLOG,DSQRT,DEXP,DABS
!
INTEGER :: MMAX,LMAX
INTEGER :: N,M,L,K
!
INTEGER :: INT
!
! ---- adjustable parameter ----
!
MMAX = 256 !
LMAX = 5 !
EFS = 0.1E0_WP !
ENOFF = 0.40E0_WP !
PQOFF = 2.9E0_WP !
PPOFF = -0.72E0_WP !
!
! ------------------------------
!
PI4 = DATAN(ONE) !
EPSLN = ONE - DLOG(EFS*EPS) !
EPSH = DSQRT(EFS*EPS) !
N = INT(ENOFF*EPSLN) !
FRQ4 = DABS(OMEGA) / (TWO*PI4) !
PER2 = FOUR * PI4 / DABS(OMEGA) !
PQ = PQOFF / EPSLN !
PP = PPOFF - DLOG(PQ*PQ*FRQ4) !
EHP = DEXP(TWO*PQ) !
EHM = ONE / EHP !
XW = DEXP(PP-TWO*PI4) !
I = F(A + DSQRT(XW * (PER2*HALF))) !
IR = I*XW !
I = I*(PER2*HALF) !
ERR = DABS(I) !
H = 2 !
M = 1 !
!
10 CONTINUE !
!
IBACK = I !
IRBACK = IR !
T = H * HALF !
!
20 CONTINUE !
!
EM = DEXP(TWO*PQ*T) !
EP = PI4 * EM !
EM = PI4 / EM !
TK = T !
!
30 CONTINUE !
!
XW = DEXP(PP-EP-EM) !
WG = DSQRT(FRQ4*XW + TK*TK) !
XA = XW / (TK+WG) !
WG = (PQ*XW*(EP-EM)+XA) / WG !
FM = F(A + XA) !
FP = F(A + XA + PER2*TK) !
IR = IR + (FP + FM)*XW !
FM = FM * WG !
FP = FP * (PER2-WG) !
I = I + (FP+FM) !
!
IF (M == 1) ERR = ERR + (DABS(FP) + DABS(FM)) !
EP = EP * EHP !
EM = EM * EHM !
TK = TK + ONE !
!
IF (EP < EPSLN) GO TO 30 !
!
IF (M == 1) THEN !
ERRH = ERR * EPSH !
ERR = ERR * EPS !
END IF !
!
TN = TK !
!
DO WHILE(DABS(FM) > ERR) !
XW = DEXP(PP-EP-EM) !
XA = XW / TK * HALF !
WG = XA * (ONE / TK + TWO*PQ*(EP - EM)) !
FM = F(A + XA) !
IR = IR + FM*XW !
FM = FM * WG !
I = I + FM !
EP = EP * EHP !
EM = EM * EHM !
TK = TK + ONE !
END DO !
!
FM = F(A + PER2*TN) !
EM = PER2 * FM !
I = I + EM !
!
IF(DABS(FP) > ERR .OR. DABS(EM) > ERR) THEN !
L = 0 !
40 CONTINUE !
L = L + 1 !
TN = TN + N !
EM = FM !
FM = F(A + PER2*TN) !
XA = FM !
EP = FM !
EM = EM + FM !
XW = ONE !
WG = ONE !
DO K = 1, N-1 !
XW = XW * (N+1-K) / K !
WG = WG + XW !
FP = F(A + PER2*(TN-K)) !
XA = XA + FP !
EP = EP + FP*WG !
EM = EM + FP*XW !
END DO !
WG = PER2 * N / (WG*N + XW) !
EM = WG * DABS(EM) !
IF(EM <= ERR .OR. L >= LMAX) GO TO 50 !
I = I + PER2*XA !
GO TO 40 !
50 CONTINUE !
I = I + WG*EP !
IF(EM > ERR) ERR = EM !
END IF !
!
T = T + H !
!
IF(T < ONE) GO TO 20 !
!
IF(M == 1) THEN !
ERRD = ONE + TWO*ERRH !
ELSE !
ERRD = H*(DABS(I-TWO*IBACK) + PQ*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 !
ELSE !
ERR = ERR * (M*HALF) !
END IF !
!
END SUBROUTINE INTDEO
!
!=======================================================================
!
SUBROUTINE INTDEO_1(F,G,AA,A,OMEGA,EPS,I,ERR)
!
! This subroutine is the integrator of g(x) over (a,infinity) when
! g(x) has an oscillatory behaviour:
!
! g(x) = f(x) * sin(omega * x + theta) * x, as x -> infinity
!
!
! This version: f(x) has ONE extra parameters: AA (real)
!
! so that f(x) = g(a*x)
!
!
! Input parameters:
!
! * F : integrand f(x)
! * AA : parameter used in function f
! * A : lower limit of integration
! * OMEGA : frequency of oscillation
! * EPS : relative error requested
!
!
! Intermediate parameters:
!
! * G : function g(a,x) to be integrated
!
!
! 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^R |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.
!
! 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
!
!
! --> f(x) changed into g(f,a,omega,x)
! where G(f,a,omega,x) = f(a,x) * sin(omega*x) * x
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,FOUR,HALF
!
IMPLICIT NONE
!
REAL (WP) :: A,AA,OMEGA,EPS,I,ERR
REAL (WP) :: EFS,ENOFF,PQOFF,PPOFF
REAL (WP) :: PI4,EPSLN,EPSH,FRQ4,PER2,PP,PQ,EHP,EHM,IR,H
REAL (WP) :: IBACK,IRBACK,T,EP,EM,TK,XW,WG,XA,FP,FM,ERRH
REAL (WP) :: TN,ERRD
!
REAL (WP), EXTERNAL :: F,G
!
REAL (WP) :: DATAN,DLOG,DSQRT,DEXP,DABS
!
INTEGER :: MMAX,LMAX
INTEGER :: N,M,L,K
!
INTEGER :: INT
!
! ---- adjustable parameter ----
!
MMAX = 256 !
LMAX = 5 !
EFS = 0.1E0_WP !
ENOFF = 0.40E0_WP !
PQOFF = 2.9E0_WP !
PPOFF = -0.72E0_WP !
!
! ------------------------------
!
PI4 = DATAN(ONE) !
EPSLN = ONE - DLOG(EFS*EPS) !
EPSH = DSQRT(EFS*EPS) !
N = INT(ENOFF*EPSLN) !
FRQ4 = DABS(OMEGA) / (TWO*PI4) !
PER2 = FOUR * PI4 / DABS(OMEGA) !
PQ = PQOFF / EPSLN !
PP = PPOFF - DLOG(PQ*PQ*FRQ4) !
EHP = DEXP(TWO*PQ) !
EHM = ONE / EHP !
XW = DEXP(PP-TWO*PI4) !
I = G(F,AA,OMEGA,A + DSQRT(XW * (PER2*HALF))) !
IR = I*XW !
I = I*(PER2*HALF) !
ERR = DABS(I) !
H = 2 !
M = 1 !
!
10 CONTINUE !
!
IBACK = I !
IRBACK = IR !
T = H * HALF !
!
20 CONTINUE !
!
EM = DEXP(TWO*PQ*T) !
EP = PI4 * EM !
EM = PI4 / EM !
TK = T !
!
30 CONTINUE !
!
XW = DEXP(PP-EP-EM) !
WG = DSQRT(FRQ4*XW + TK*TK) !
XA = XW / (TK+WG) !
WG = (PQ*XW*(EP-EM)+XA) / WG !
FM = G(F,AA,OMEGA,A + XA) !
FP = G(F,AA,OMEGA,A + XA + PER2*TK) !
IR = IR + (FP + FM)*XW !
FM = FM * WG !
FP = FP * (PER2-WG) !
I = I + (FP+FM) !
!
IF (M == 1) ERR = ERR + (DABS(FP) + DABS(FM)) !
EP = EP * EHP !
EM = EM * EHM !
TK = TK + ONE !
!
IF (EP < EPSLN) GO TO 30 !
!
IF (M == 1) THEN !
ERRH = ERR * EPSH !
ERR = ERR * EPS !
END IF !
!
TN = TK !
!
DO WHILE(DABS(FM) > ERR) !
XW = DEXP(PP-EP-EM) !
XA = XW / TK * HALF !
WG = XA * (ONE / TK + TWO*PQ*(EP - EM)) !
FM = G(F,AA,OMEGA,A + XA) !
IR = IR + FM*XW !
FM = FM * WG !
I = I + FM !
EP = EP * EHP !
EM = EM * EHM !
TK = TK + ONE !
END DO !
!
FM = G(F,AA,OMEGA,A + PER2*TN) !
EM = PER2 * FM !
I = I + EM !
!
IF(DABS(FP) > ERR .OR. DABS(EM) > ERR) THEN !
L = 0 !
40 CONTINUE !
L = L + 1 !
TN = TN + N !
EM = FM !
FM = G(F,AA,OMEGA,A + PER2*TN) !
XA = FM !
EP = FM !
EM = EM + FM !
XW = ONE !
WG = ONE !
DO K = 1, N-1 !
XW = XW * (N+1-K) / K !
WG = WG + XW !
FP = F(A + PER2*(TN-K)) !
XA = XA + FP !
EP = EP + FP*WG !
EM = EM + FP*XW !
END DO !
WG = PER2 * N / (WG*N + XW) !
EM = WG * DABS(EM) !
IF(EM <= ERR .OR. L >= LMAX) GO TO 50 !
I = I + PER2*XA !
GO TO 40 !
50 CONTINUE !
I = I + WG*EP !
IF(EM > ERR) ERR = EM !
END IF !
!
T = T + H !
!
IF(T < ONE) GO TO 20 !
!
IF(M == 1) THEN !
ERRD = ONE + TWO*ERRH !
ELSE !
ERRD = H*(DABS(I-TWO*IBACK) + PQ*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 !
ELSE !
ERR = ERR * (M*HALF) !
END IF !
!
END SUBROUTINE INTDEO_1
!
!=======================================================================
!
SUBROUTINE INTDEO_2(F,G,AA,LL,A,OMEGA,EPS,I,ERR)
!
! This subroutine is the integrator of g(x) over (a,infinity) when
! g(x) has an oscillatory behaviour:
!
! g(x) = f(x) * sin(omega * x + theta), as x -> infinity
!
!
! This version: f(x) has TWO extra parameters: AA (real)
! LL (integer)
! so that f(x) = g_l(a*x)
!
!
! Input parameters:
!
! * F : function f(x) to be Fourier transformed
! * AA : parameter used in function F
! * LL : order of the oscillatory function F
! * A : lower limit of integration
! * OMEGA : frequency of oscillation
! * EPS : relative error requested
!
!
! Intermediate parameters:
!
! * G : function g(x) to be integrated
!
!
! 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^R |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.
!
! 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
!
!
! --> f(x) changed into g(f,a,l,omega,x)
! where G(f,a,l,omega,x) = f(a,x) * j_l(omega*x) * x^2
! --> a replaced by AA, as A is used for the lower integration bound
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,FOUR,HALF
!
IMPLICIT NONE
!
REAL (WP) :: A,AA,OMEGA,EPS,I,ERR
REAL (WP) :: EFS,ENOFF,PQOFF,PPOFF
REAL (WP) :: PI4,EPSLN,EPSH,FRQ4,PER2,PP,PQ,EHP,EHM,IR,H
REAL (WP) :: IBACK,IRBACK,T,EP,EM,TK,XW,WG,XA,FP,FM,ERRH
REAL (WP) :: TN,ERRD
!
REAL (WP), EXTERNAL :: F,G
!
REAL (WP) :: DATAN,DLOG,DSQRT,DEXP,DABS
!
INTEGER :: MMAX,LMAX
INTEGER :: N,M,L,K
INTEGER :: LL
!
INTEGER :: INT
!
! ---- adjustable parameter ----
!
MMAX = 256 !
LMAX = 5 !
EFS = 0.1E0_WP !
ENOFF = 0.40E0_WP !
PQOFF = 2.9E0_WP !
PPOFF = -0.72E0_WP !
!
! ------------------------------
!
PI4 = DATAN(ONE) !
EPSLN = ONE - DLOG(EFS*EPS) !
EPSH = DSQRT(EFS*EPS) !
N = INT(ENOFF*EPSLN) !
FRQ4 = DABS(OMEGA) / (TWO*PI4) !
PER2 = FOUR * PI4 / DABS(OMEGA) !
PQ = PQOFF / EPSLN !
PP = PPOFF - DLOG(PQ*PQ*FRQ4) !
EHP = DEXP(TWO*PQ) !
EHM = ONE / EHP !
XW = DEXP(PP-TWO*PI4) !
I = G(F,AA,OMEGA,A + DSQRT(XW * (PER2*HALF))) !
IR = I*XW !
I = I*(PER2*HALF) !
ERR = DABS(I) !
H = 2 !
M = 1 !
!
10 CONTINUE !
!
IBACK = I !
IRBACK = IR !
T = H * HALF !
!
20 CONTINUE !
!
EM = DEXP(TWO*PQ*T) !
EP = PI4 * EM !
EM = PI4 / EM !
TK = T !
!
30 CONTINUE !
!
XW = DEXP(PP-EP-EM) !
WG = DSQRT(FRQ4*XW + TK*TK) !
XA = XW / (TK+WG) !
WG = (PQ*XW*(EP-EM)+XA) / WG !
FM = G(F,AA,LL,OMEGA,A + XA) !
FP = G(F,AA,LL,OMEGA,A + XA + PER2*TK) !
IR = IR + (FP + FM)*XW !
FM = FM * WG !
FP = FP * (PER2-WG) !
I = I + (FP+FM) !
!
IF (M == 1) ERR = ERR + (DABS(FP) + DABS(FM)) !
EP = EP * EHP !
EM = EM * EHM !
TK = TK + ONE !
!
IF (EP < EPSLN) GO TO 30 !
!
IF (M == 1) THEN !
ERRH = ERR * EPSH !
ERR = ERR * EPS !
END IF !
!
TN = TK !
!
DO WHILE(DABS(FM) > ERR) !
XW = DEXP(PP-EP-EM) !
XA = XW / TK * HALF !
WG = XA * (ONE / TK + TWO*PQ*(EP - EM)) !
FM = G(F,AA,LL,OMEGA,A + XA) !
IR = IR + FM*XW !
FM = FM * WG !
I = I + FM !
EP = EP * EHP !
EM = EM * EHM !
TK = TK + ONE !
END DO !
!
FM = G(F,AA,LL,OMEGA,A + PER2*TN) !
EM = PER2 * FM !
I = I + EM !
!
IF(DABS(FP) > ERR .OR. DABS(EM) > ERR) THEN !
L = 0 !
40 CONTINUE !
L = L + 1 !
TN = TN + N !
EM = FM !
FM = G(F,AA,LL,OMEGA,A + PER2*TN) !
XA = FM !
EP = FM !
EM = EM + FM !
XW = ONE !
WG = ONE !
DO K = 1, N-1 !
XW = XW * (N+1-K) / K !
WG = WG + XW !
FP = F(A + PER2*(TN-K)) !
XA = XA + FP !
EP = EP + FP*WG !
EM = EM + FP*XW !
END DO !
WG = PER2 * N / (WG*N + XW) !
EM = WG * DABS(EM) !
IF(EM <= ERR .OR. L >= LMAX) GO TO 50 !
I = I + PER2*XA !
GO TO 40 !
50 CONTINUE !
I = I + WG*EP !
IF(EM > ERR) ERR = EM !
END IF !
!
T = T + H !
!
IF(T < ONE) GO TO 20 !
!
IF(M == 1) THEN !
ERRD = ONE + TWO*ERRH !
ELSE !
ERRD = H*(DABS(I-TWO*IBACK) + PQ*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 !
ELSE !
ERR = ERR * (M*HALF) !
END IF !
!
END SUBROUTINE INTDEO_2
!
!=======================================================================
!
SUBROUTINE INTDEO_F(F,A,OMEGA,AW,I,ERR)
!
! This subroutine is the integrator of f(x) over (a,infinity) when
! f(x) has an oscillatory behaviour:
!
! f(x) = g(x) * sin(omega * x + theta) as x -> infinity
!
!
!
! --> <--
! --> This is the fast version <--
! --> <--
!
!
! Usage:
!
! CALL INTDEOINI(LENAW,TINY,EPS,AW) ! initialization of AW
! ...
! CALL INTDEO_F(F,A,AW,I,ERR)
!
!
! Input parameters:
!
! * F : integrand f(x)
! * A : lower limit of integration
! * OMEGA : frequency of oscillation
! * EPS : relative error requested
!
!
! 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^R |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.
!
! 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 : ZERO,ONE,TWO,HALF
!
IMPLICIT NONE
!
REAL (WP) :: A,OMEGA,AW(0 : *),I,ERR
REAL (WP) :: EPS,PER,PERW,W02,IR,H,IBACK,IRBACK,T,TK
REAL (WP) :: XA,FM,FP,ERRH,S0,S1,S2,ERRD
!
REAL (WP), EXTERNAL :: F
!
REAL (WP) :: DABS
!
INTEGER :: LENAWM,NK0,NOFF0,NK,NOFF,LMAX,M,K,J,JM,L
!
INTEGER :: INT
!
LENAWM = INT(AW(0) + HALF) !
NK0 = INT(AW(1) + HALF) !
NOFF0 = 6 !
NK = INT(AW(2) + HALF) !
NOFF = 2*NK0 + NOFF0 !
LMAX = INT(AW(3) + HALF) !
EPS = AW(4) !
PER = ONE / DABS(OMEGA) !
W02 = TWO * AW(NOFF+2) !
PERW = PER * W02 !
I = F(A + AW(NOFF)*PER) !
IR = I * AW(NOFF+1) !
I = I * AW(NOFF+2) !
ERR = DABS(I) !
H = TWO !
M = 1 !
K = NOFF !
!
10 CONTINUE !
!
IBACK = I !
IRBACK = IR !
T = H * HALF !
!
20 CONTINUE !
!
IF(K == NOFF) THEN !
TK = ONE !
K = K + NK !
J = NOFF !
!
30 CONTINUE !
!
J = J + 3 !
XA = PER * AW(J) !
FM = F(A + XA) !
FP = F(A + XA + PERW*TK) !
IR = IR + (FM+FP) * AW(J+1) !
FM = FM * AW(J+2) !
FP = FP * (W02-AW(J+2)) !
I = I + (FM+FP) !
ERR = ERR + (DABS(FM)+DABS(FP)) !
TK = TK + ONE !
!
IF(AW(J) > EPS .AND. J < K) GO TO 30 !
!
ERRH = ERR * AW(5) !
ERR = ERR * EPS !
JM = J - NOFF !
ELSE !
TK = T !
DO J = K+3, K+JM, 3 !
XA = PER * AW(J) !
FM = F(A + XA) !
FP = F(A + XA + PERW*TK) !
IR = IR + (FM+FP) * AW(J+1) !
FM = FM * AW(J+2) !
FP = FP * (W02-AW(J+2)) !
I = I + (FM+FP) !
TK = TK + ONE !
END DO !
J = K + JM !
K = K + NK !
END IF !
!
DO WHILE (DABS(FM) > ERR .AND. J < K) !
J = J + 3 !
FM = F(A + PER*AW(J)) !
IR = IR + FM*AW(J+1) !
FM = FM * AW(J+2) !
I = I + FM !
END DO !
!
FM = F(A + PERW*TK) !
S2 = W02 * FM !
I = I + S2 !
!
IF(DABS(FP) > ERR .OR. DABS(S2) > ERR) THEN !
L = 0 !
!
40 CONTINUE !
!
L = L + 1 !
S0 = ZERO !
S1 = ZERO !
S2 = FM * AW(NOFF0+1) !
!
DO J = NOFF0+2, NOFF-2, 2 !
TK = TK + ONE !
FM = F(A + PERW*TK) !
S0 = S0 + FM !
S1 = S1 + FM*AW(J) !
S2 = S2 + FM*AW(J+1) !
END DO
!
IF(S2 <= ERR .OR. L >= LMAX) GO TO 50 !
!
I = I + W02*S0 !
GO TO 40 !
!
50 CONTINUE !
!
I = I + S1 !
!
IF(S2 > ERR) ERR = S2 !
!
END IF !
!
T = T + H !
!
IF(T < ONE) GO TO 20 !
!
IF(M == 1) THEN !
ERRD = ONE + TWO*ERRH !
ELSE !
ERRD = H * (DABS(I-2*IBACK) + DABS(IR- 2*IRBACK)) !
END IF !
!
H = H * HALF !
M = M * 2 !
!
IF(ERRD > ERRH .AND. (2*K - NOFF) <= LENAWM) GO TO 10 !
!
I = I * (H*PER) !
!
IF(ERRD > ERRH) THEN !
ERR = -ERRD*PER !
ELSE !
ERR = ERR * (PER*M*HALF) !
END IF !
!
END SUBROUTINE INTDEO_F
!
!=======================================================================
!
SUBROUTINE INTDEOINI_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 : ZERO,ONE,TWO,FOUR,HALF
!
IMPLICIT NONE
!
REAL (WP) :: TINY,EPS,AW(0 : LENAW - 1)
REAL (WP) :: EFS,ENOFF,PQOFF,PPOFF
REAL (WP) :: PI4,TINYLN,EPSLN,FRQ4,PER2,PP,PQ,EHP,EHM,H
REAL (WP) :: T,EP,EM,TK,XW,WG,XA
!
REAL (WP) :: DATAN,DLOG,DSQRT,DEXP,DFLOAT
!
INTEGER :: NOFF0,NK0,NOFF,K,NK,J
INTEGER :: LENAW
INTEGER :: LMAX
!
INTEGER :: INT
!
! ---- adjustable parameter ----
!
LMAX = 5 !
EFS = 0.1E0_WP !
ENOFF = 0.40E0_WP !
PQOFF = 2.9E0_WP !
PPOFF = -0.72E0_WP !
!
! ------------------------------
!
PI4 = DATAN(ONE) !
TINYLN = -DLOG(TINY) !
EPSLN = ONE - DLOG(EFS*EPS) !
FRQ4 = ONE / (TWO*PI4) !
PER2 = FOUR * PI4 !
PQ = PQOFF / EPSLN !
PP = PPOFF - DLOG(PQ*PQ*FRQ4) !
EHP = DEXP(TWO*PQ) !
EHM = ONE / EHP !
AW(3) = DFLOAT(LMAX) !
AW(4) = EPS !
AW(5) = DSQRT(EFS*EPS) !
NOFF0 = 6 !
NK0 = 1 + INT(ENOFF*EPSLN) !
AW(1) = NK0 !
NOFF = 2*NK0 + NOFF0 !
WG = ZERO !
XW = ONE !
!
DO K = 1, NK0 !
WG = WG + XW !
AW(NOFF - 2*K) = WG !
AW(NOFF - 2*K + 1) = XW !
XW = XW * (NK0-K) / K !
END DO !
!
WG = PER2 / WG !
!
DO K = NOFF0, NOFF-2, 2 !
AW(K) = AW(K)*WG !
AW(K+1) = AW(K+1)*WG !
END DO !
!
XW = DEXP(PP - TWO*PI4) !
AW(NOFF) = DSQRT(XW * (PER2*HALF)) !
AW(NOFF+1) = XW * PQ !
AW(NOFF+2) = PER2 * HALF !
H = TWO !
NK = 0 !
K = NOFF + 3 !
!
10 CONTINUE !
!
T = H * HALF !
!
20 CONTINUE !
!
EM = DEXP(2*PQ*T) !
EP = PI4 * EM !
EM = PI4 / EM !
TK = T !
J = K !
!
30 CONTINUE !
!
XW = DEXP(PP - EP - EM) !
WG = DSQRT(FRQ4*XW + TK*TK) !
XA = XW / (TK + WG) !
WG = (PQ*XW*(EP-EM) + XA) / WG !
AW(J) = XA !
AW(J+1) = XW * PQ !
AW(J+2) = WG !
EP = EP * EHP !
EM = EM * EHM !
TK = TK + ONE !
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(2) = NK !
END IF !
!
IF ((2*K - NOFF - 3) <= LENAW) GO TO 10 !
!
AW(0) = DFLOAT(K-3) !
!
END SUBROUTINE INTDEOINI_F
!
END MODULE INTEGRATION3