MsSpec-DFM/New_libraries/DFM_library/VARIOUS_FUNCTIONS_LIBRARY/2F1_real.f90

375 lines
19 KiB
Fortran
Raw Normal View History

2022-02-02 16:19:10 +01:00
!
!=======================================================================
!
MODULE CONFLUENT_HYPGEOM_REAL
!
! This module provides several subroutines to compute
! the Gauss hypergeometric function
!
! 2F1(a,b,c;x)
!
!
! --> <--
! --> x real <--
! --> <--
!
!
!
! 1) SUBROUTINE HYGFX(A,B,C,X,HF)
!
! 2) SUBROUTINE HYP(Z,A,B,C,RE,IM)
!
!
!
USE ACCURACY_REAL
!
CONTAINS
!
!=======================================================================
!
SUBROUTINE HYGFX(A,B,C,X,HF)
!
! HYGFX computes the hypergeometric function F(a,b,c,x).
!
! Licensing:
!
! This routine is copyrighted by Shanjie Zhang and Jianming Jin. However,
! they give permission to incorporate this routine into a user program
! provided that the copyright is acknowledged.
!
! Modified:
!
! 04 April 2012
!
! Author:
!
! Shanjie Zhang, Jianming Jin
!
! Reference:
!
! Shanjie Zhang, Jianming Jin,
! Computation of Special Functions,
! Wiley, 1996,
! ISBN: 0-471-11963-6,
! LC: QA351.C45.
!
! ==========================================================
!
! Purpose: Compute hypergeometric function F(a,b,c,x)
!
! Input : A --- Parameter
! B --- Parameter
! C --- Parameter, C <> 0,-1,-2,...
! X --- Argument (X < 1)
!
! Output: HF --- F(A,B,C,X)
!
! Routines called:
! (1) GAMMA for computing gamma function
! (2) PSI for computing psi function
!
! ==========================================================
!
!
! Last modified (DS) : 1 Sep 2020
!
!
USE REAL_NUMBERS, ONLY : ZERO,ONE,TWO,HALF
USE PI_ETC, ONLY : PI
USE GAMMA_FUNCTION, ONLY : GAMMA
USE DIGAMMA_FUNCTION, ONLY : PSI
!
IMPLICIT NONE
!
INTEGER :: J,K,M,NM
!
REAL (WP), INTENT(IN) :: C
REAL (WP) :: A,B,X
REAL (WP), INTENT(OUT) :: HF
!
REAL (WP) :: C0
REAL (WP) :: EPS
REAL (WP) :: G0,G1,G2,G3
REAL (WP) :: GC,R
REAL (WP) :: A0,AA,BB,C1
REAL (WP) :: F0,F1
REAL (WP) :: GA,GABC,GAM,GB,GBM,GCA,GCAB,GCB,GM
REAL (WP) :: HW
REAL (WP) :: PA,PB,R1,RM,RP
REAL (WP) :: SM,SP,SP0
REAL (WP) :: R0,X1
!
REAL, PARAMETER :: EL = 0.5772156649015329E+00_WP
!
LOGICAL :: L0,L1,L2,L3,L4,L5
!
L0 = C == INT(C) .AND. C < ZERO !
L1 = ONE - X < 1.0E-15_WP .AND. C - A - B <= ZERO !
L2 = A == INT(A) .AND. A < ZERO !
L3 = B == INT(B) .AND. B < ZERO !
L4 = C - A == INT(C - A) .AND. C - A <= ZERO !
L5 = C - B == INT(C - B) .AND. C - B <= ZERO !
!
IF(L0 .OR. L1) THEN !
WRITE(6,110) !
WRITE(6,111) !
WRITE(6,112) !
STOP !
END IF !
!
IF(0.95E+00_WP < X) THEN !
EPS = 1.0E-08_WP !
ELSE !
EPS = 1.0E-15_WP !
END IF !
!
IF(X == ZERO .OR. A == ZERO .OR. B == ZERO) THEN !
HF = ONE !
RETURN !
ELSE IF(ONE - X == EPS .AND. ZERO < C - A - B) THEN !
CALL GAMMA(C, GC) !
CALL GAMMA(C - A - B, GCAB) !
CALL GAMMA(C - A, GCA) !
CALL GAMMA(C - B, GCB) !
HF = GC * GCAB / (GCA * GCB) !
RETURN !
ELSE IF(ONE + X <= EPS .AND. & !
DABS(C - A + B -ONE) <= EPS) THEN !
G0 = DSQRT(PI) * TWO ** (- A) !
CALL GAMMA(C, G1) !
CALL GAMMA(ONE + A / TWO - B, G2) !
CALL GAMMA(HALF + HALF * A, G3) !
HF = G0 * G1 / (G2 * G3) !
RETURN !
ELSE IF(L2 .OR. L3) THEN !
IF(L2) THEN !
NM = INT(DABS(A)) !
END IF !
IF(L3) THEN !
NM = INT(DABS(B)) !
END IF !
HF = ONE !
R = ONE !
DO K = 1, NM !
R = R * (A + K - ONE) * (B + K - ONE) / & !
(K * (C + K - ONE)) * X !
HF = HF + R !
END DO !
RETURN !
ELSE IF(L4 .OR. L5) THEN !
IF(L4) THEN !
NM = INT(DABS(C - A)) !
END IF !
IF(L5) THEN !
NM = INT(DABS(C - B)) !
END IF !
HF = ONE !
R = ONE !
DO K = 1, NM !
R = R * (C - A + K - ONE) * (C - B + K - ONE) / & !
(K * (C + K - ONE)) * X !
HF = HF + R !
END DO !
HF = (ONE - X) ** (C - A - B) * HF !
RETURN !
END IF !
!
AA = A !
BB = B !
X1 = X !
!
IF(X < ZERO) THEN !
X = X / (X - ONE) !
IF(A < C .AND. B < A .AND. ZERO < B) THEN !
A = BB !
B = AA !
END IF !
B = C - B !
END IF !
!
IF(0.75E+00_WP <= X) THEN !
GM = ZERO !
IF(DABS(C - A - B - INT(C - A - B)) < 1.0E-15_WP) THEN !
M = INT(C - A - B) !
CALL GAMMA(A, GA) !
CALL GAMMA(B, GB) !
CALL GAMMA(C, GC) !
CALL GAMMA(A + M, GAM) !
CALL GAMMA(B + M, GBM) !
CALL PSI(A, PA) !
CALL PSI(B, PB) !
IF(M .NE. 0) THEN !
GM = ONE !
END IF !
DO J = 1, ABS (M) - 1 !
GM = GM * J !
END DO !
RM = ONE !
DO J = 1, ABS (M) !
RM = RM * J !
END DO !
F0 = ONE !
R0 = ONE !
R1 = ONE !
SP0 = ZERO !
SP = ZERO !
IF(0 <= M) THEN !
C0 = GM * GC / (GAM * GBM) !
C1 = - GC * (X - ONE) ** M / (GA * GB * RM) !
DO K = 1, M - 1 !
R0 = R0 * (A + K - ONE) * (B + K - ONE) / & !
(K * (K - M)) * (ONE - X) !
F0 = F0 + R0 !
END DO !
!
DO K = 1, M !
SP0 = SP0 + ONE / (A + K - ONE) + & !
ONE / (B + K - ONE) - ONE / K !
END DO !
!
F1 = PA + PB + SP0 + TWO * EL + DLOG(ONE - X) !
!
DO K = 1, 250 !
SP = SP + (ONE - A) / (K * (A + K - ONE)) + & !
(ONE - B) / (K * (B + K - ONE)) !
SM = ZERO !
DO J = 1,M !
SM = SM + (ONE - A) / ((J + K) * & !
(A + J + K - ONE)) + ONE / & !
(B + J + K - ONE) !
END DO !
RP = PA + PB + TWO * EL + SP + SM + DLOG(ONE - X) !
R1 = R1 * (A + M + K - ONE) * (B + M + K - ONE) / & !
(K * (M + K)) * (ONE - X) !
F1 = F1 + R1 * RP !
IF(DABS(F1 - HW) < DABS(F1) * EPS) THEN !
GO TO 10 !
END IF !
HW = F1 !
END DO !
!
10 CONTINUE
!
HF = F0 *C0 + F1 * C1 !
!
ELSE IF(M < 0) THEN !
!
M = - M !
C0 = GM * GC / (GA * GB * (ONE - X)**M) !
C1 = - (-1)**M * GC / (GAM * GBM * RM) !
!
DO K = 1,M-1 !
R0 = R0 * (A - M + K - ONE) * & !
(B - M + K - ONE) / & !
(K * (K - M)) * (ONE - X) !
F0 = F0 + R0 !
END DO !
!
DO K = 1,M !
SP0 = SP0 + ONE / K !
END DO !
F1 = PA + PB - SP0 + TWO * EL + DLOG(ONE - X) !
DO K = 1,250 !
SP = SP + (ONE - A) / (K * (A + K - ONE)) + & !
(ONE - B) / (K * (B + K - ONE)) !
SM = ZERO !
DO J = 1,M !
SM = SM+ ONE /(J + K) !
END DO !
RP = PA + PB + TWO * EL + SP - SM + & !
DLOG(ONE - X) !
R1 = R1 * (A + K - ONE) * (B + K - ONE) / & !
(K * (M + K)) * (ONE - X) !
F1 = F1 + R1 * RP !
IF(DABS(F1 - HW) < DABS(F1) * EPS) THEN !
GO TO 20 !
END IF !
HW = F1 !
END DO !
!
20 CONTINUE !
!
HF = F0 * C0 + F1 * C1 !
END IF !
ELSE !
CALL GAMMA(A,GA) !
CALL GAMMA(B,GB) !
CALL GAMMA(C,GC) !
CALL GAMMA(C-A,GCA) !
CALL GAMMA(C-B,GCB) !
CALL GAMMA(C-A-B,GCAB) !
CALL GAMMA(A+B-C,GABC) !
C0 = GC * GCAB / (GCA * GCB) !
C1 = GC * GABC / (GA * GB) * (ONE - X)**(C - A - B) !
HF = ZERO !
R0 = C0 !
R1 = C1 !
DO K = 1,250 !
R0 = R0 * (A + K - ONE) * (B + K - ONE) / & !
(K * (A + B - C + K)) * (ONE - X) !
R1 = R1 * (C - A + K - ONE) * (C - B + K - ONE) /& !
(K * (C - A - B + K)) * (ONE - X) !
HF = HF + R0 + R1 !
IF(DABS(HF - HW) < DABS(HF) * EPS) THEN !
GO TO 30 !
END IF !
HW = HF !
END DO !
!
30 CONTINUE !
!
HF = HF + C0 + C1 !
!
END IF !
ELSE !
A0 = ONE !
IF(C > A .AND. C < TWO * A .AND. & !
C > B .AND. C < TWO*B) THEN !
A0 = (ONE - X) ** (C - A - B) !
A = C - A !
B = C - B !
END IF !
HF = ONE !
R = ONE !
DO K = 1,250 !
R = R *(A + K - ONE) * (B + K - ONE) / & !
(K * (C + K - ONE)) * X !
HF = HF + R !
IF(DABS(HF - HW) <= DABS(HF) * EPS) THEN !
GO TO 40 !
END IF !
HW = HF !
END DO !
!
40 CONTINUE !
!
HF = A0 * HF !
!
END IF !
!
IF(X1 < ZERO) THEN !
X = X1 !
C0 = ONE / (ONE - X) ** AA !
HF = C0 * HF !
END IF
!
A = AA !
B = BB !
!
IF(120 < K) THEN !
WRITE(*,115) !
END IF !
!
! Formats:
!
110 FORMAT(' ')
111 FORMAT('HYGFX - Fatal error!')
112 FORMAT('The hypergeometric series is divergent.')
115 FORMAT('Warning! You should check the accuracy')
!
RETURN !
!
END SUBROUTINE HYGFX
!
END MODULE CONFLUENT_HYPGEOM_REAL