! !======================================================================= ! MODULE CALCULATORS_3 ! USE ACCURACY_REAL ! ! This module contains the subroutines allowing to compute ! various properties of the electron/plasma liquids: ! ! * local field corrections : CALC_LFC ! * I(q) function : CALC_IQF ! * structure factor : CALC_SFC ! * pair correlation function : CALC_PCF ! ! * pair distribution function : CALC_PDF ! * vertex function : CALC_VTX ! * plasmon damping coefficient: CALC_DMP ! ! * q bounds : CALC_QBD ! * relaxation time : CALC_RLX ! * screening wave number : CALC_SCR ! * omega = q * v_F : CALC_QVF ! ! * moments of Im[epsilon] : CALC_MEP ! * moments of S(q,omega) : CALC_MSF ! * moments of loss function : CALC_MLO ! ! * zeros of Re [epsilon] : CALC_RE0 ! ! * inelastic mean free path : CALC_MFP ! ! * Fourier-space Nevalinna : CALC_NEV ! * time-space memory function : CALC_MEM ! CONTAINS ! !======================================================================= ! SUBROUTINE CALC_LFC(X) ! ! This subroutine computes the local field correction G(q, omega) ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Intermediate parameters: ! ! * RS : Wigner-Seitz radius (in units of a_0) ! * T : temperature in SI ! ! ! Output parameters: ! ! * E : energy array ! * LFCR : real part of local field correction ! * LFCI : imaginary part of local field correction ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 3 Dec 2020 ! USE DIMENSION_CODE, ONLY : NSIZE USE REAL_NUMBERS, ONLY : ZERO USE MATERIAL_PROP, ONLY : RS USE EXT_FIELDS, ONLY : T USE LF_VALUES, ONLY : GSTDY,GQ_TYPE,GQO_TYPE USE SF_VALUES, ONLY : SQ_TYPE ! USE E_GRID ! USE LOCAL_FIELD_STATIC USE LOCAL_FIELD_STATIC_2 USE LOCAL_FIELD_DYNAMIC ! USE OUT_VALUES_3, ONLY : I_LF USE PRINT_FILES, ONLY : IO_LF ! IMPLICIT NONE ! REAL (WP) :: X REAL (WP) :: LFCR(NSIZE),LFCI(NSIZE),E(NSIZE) REAL (WP) :: RLFC,ILFC REAL (WP) :: EN,ETA REAL (WP) :: Y ! INTEGER :: IE,LOGF ! LOGF = 6 ! ! Y = X + X ! q/k_F ! IF(GSTDY == ' STATIC') THEN ! ! ! Static local field correction: ! ! 1) G(q) not based on S(q) ! IF(GQ_TYPE /= 'HNCA' .AND. GQ_TYPE /= 'IKPA') THEN ! CALL LFIELD_STATIC(X,RS,T,GQ_TYPE,RLFC) ! ELSE ! ! 1) G(q) based on S(q) ! IF(SQ_TYPE /= 'GEA' .AND. SQ_TYPE /= 'ICH' .AND. & ! SQ_TYPE /= 'PKA' .AND. SQ_TYPE /= 'SIN' .AND. & ! SQ_TYPE /= 'SPA') THEN ! CALL LFIELD_STATIC_2(X,RS,T,GQ_TYPE,RLFC) ! ELSE ! WRITE(LOGF,10) ! STOP ! END IF ! ! END IF ! ! LFCR(1) = RLFC ! LFCI(1) = ZERO ! ! IF(I_LF == 1) THEN ! writing to WRITE(IO_LF,*) Y,ZERO,RLFC,ZERO ! file END IF ! ! ELSE ! ! ! Dynamic local field correction ! CALL LFIELD_DYNAMIC(X,RS,E_MIN,E_MAX,N_E,T,ETA, & ! GQO_TYPE,E,LFCR,LFCI) ! ! IF(I_LF == 1) THEN ! DO IE = 1, N_E ! writing to WRITE(IO_LF,*) Y,E(IE),LFCR(IE),LFCI(IE) ! file END DO ! END IF ! ! END IF ! ! ! Format: ! 10 FORMAT(//,10X,'<<<<< ERROR IN CALCULATOR_3 :: CALC_LFC >>>>>', & /,10X,'<<<<< WRONG CHOICE OF GQ_TYPE :: SQ_TYPE >>>>>', & //,10X,' --> G(q) choice uses S(q) which uses G(q) <--',//) ! END SUBROUTINE CALC_LFC ! !======================================================================= ! SUBROUTINE CALC_IQF(X) ! ! This subroutine computes the local field correction I(q) = G(q,inf) ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Intermediate parameters: ! ! * RS : Wigner-Seitz radius (in units of a_0) ! * T : temperature in SI ! ! ! Output parameters: ! ! * IQR : real part of I(q) ! * IQI : imaginary of I(q) ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 15 Sep 2020 ! USE MATERIAL_PROP, ONLY : RS,DMN USE EXT_FIELDS, ONLY : T ! USE REAL_NUMBERS, ONLY : ZERO USE LF_VALUES, ONLY : GQ_TYPE,IQ_TYPE USE SF_VALUES, ONLY : SQ_TYPE USE ENERGIES, ONLY : EC_TYPE USE IQ_FUNCTIONS_1 ! USE OUT_VALUES_3, ONLY : I_IQ USE PRINT_FILES, ONLY : IO_IQ ! IMPLICIT NONE ! REAL (WP), INTENT(IN) :: X REAL (WP) :: IQ REAL (WP) :: Y ! Y = X + X ! q/k_F ! IF(DMN == '3D') THEN ! CALL IQ_3D(X,RS,IQ_TYPE,IQ) ! ELSE IF(DMN == '2D') THEN ! CONTINUE ! not implemented yet ELSE IF(DMN == '1D') THEN ! CONTINUE ! not implemented yet END IF ! ! IF(I_IQ == 1) THEN ! writing to WRITE(IO_IQ,*) Y,IQ ! file END IF ! ! END SUBROUTINE CALC_IQF ! !======================================================================= ! SUBROUTINE CALC_SFC(X) ! ! This subroutine computes the structure factor S(q, omega) ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Intermediate parameters: ! ! * RS : Wigner-Seitz radius (in units of a_0) ! * T : temperature in SI ! ! ! Output parameters: ! ! * E : energy array ! * SFCR : real part of the structure factor ! * SFCI : imaginary part of structure factor ! ! ! Note: as S(q,omega) is proportional to the inverse of a ! frequency, we renormalize it by multiplicating it ! by E_F / h_bar to obtain a dimensionless quantity ! ! ! ! Author : D. Sébilleau ! ! Last modified : 21 Dec 2020 ! USE DIMENSION_CODE, ONLY : NSIZE ! USE MATERIAL_PROP, ONLY : DMN,RS USE EXT_FIELDS, ONLY : T ! USE REAL_NUMBERS, ONLY : ZERO,ONE,FOURTH USE CONSTANTS_P1, ONLY : H_BAR USE FERMI_SI, ONLY : EF_SI ! USE E_GRID ! USE LF_VALUES, ONLY : GQ_TYPE,IQ_TYPE USE SF_VALUES USE STRUCTURE_FACTOR_STATIC USE STRUCTURE_FACTOR_STATIC_2 USE STRUCTURE_FACTOR_DYNAMIC USE STRUCTURE_FACTOR_DYNAMIC_2 ! USE UTIC_VALUES USE RELAXATION_TIME_STATIC USE PLASMON_ENE_SI USE PLASMON_DISP_REAL USE DECAY_RATE USE UTIC_PARAMETERS ! USE OUT_VALUES_3, ONLY : I_SF USE PRINT_FILES, ONLY : IO_SF ! IMPLICIT NONE ! REAL (WP) :: X REAL (WP) :: RSFC,ISFC REAL (WP) :: Z,EN,SQ,SQ_RN REAL (WP) :: Y REAL (WP) :: HOM_Q REAL (WP) :: E1,E2,E3,E4,E5,E6 ! REAL (WP) :: FLOAT ! INTEGER :: IE,LOGF ! LOGF = 6 ! ! Y = X + X ! q/k_F ! ! Storing the omega-independent UTIC parameters whenever necessary ! IF(SQO_TYPE == 'UTI') THEN ! IF(DMN == '3D') THEN ! TAU_Q = UTIC_RT_3D(X,RS,T,SQ_TYPE,GQ_TYPE) ! tau(q) CALL PLASMON_DISP_3D_2(X,RS,T,'UTI_MOD',HOM_Q) ! OM_Q = HOM_Q / H_BAR ! omega(q) GAM_Q = UTIC_DR_3D(X,RS,T,SQ_TYPE,GQ_TYPE,IQ_TYPE) ! gamma(q) CALL UTIC_PARAM(X,RS,T,MO_Q,MO_0) ! Omega(q), Omega(0) ! E1 = H_BAR / (TAU_Q * EF_SI) ! E2 = ENE_P_SI / EF_SI ! E3 = H_BAR * OM_Q / EF_SI ! E4 = ABS(H_BAR * GAM_Q / EF_SI) ! E5 = H_BAR * MO_Q / EF_SI ! E6 = H_BAR * MO_0 / EF_SI ! ! WRITE(1,*) Y,E1,E2,E3,E4,E5,E6 ! END IF END IF ! ! IF(SSTDY == ' STATIC') THEN ! ! ! Static structure factor ! ! 1) S(q) not based on G(q) ! IF(SQ_TYPE /= 'GEA' .AND. SQ_TYPE /= 'ICH' .AND. & ! SQ_TYPE /= 'PKA' .AND. SQ_TYPE /= 'SIN' .AND. & ! SQ_TYPE /= 'SPA') THEN ! CALL STFACT_STATIC(X,RS,T,SQ_TYPE,RSFC) ! ! ! 2) S(q) based on G(q) ! ELSE ! ! IF(GQ_TYPE /= 'HNCA' .AND. GQ_TYPE /= 'IKPA') THEN ! CALL STFACT_STATIC_2(X,RS,T,SQ_TYPE,RSFC) ! ELSE WRITE(LOGF,10) ! STOP ! END IF ! ! END IF ! IF(I_SF == 1) THEN ! writing to WRITE(IO_SF,*) Y,ZERO,RSFC,ZERO ! file END IF ! ! ELSE ! ! ! Dynamic structure factor ! DO IE = 1, N_E ! energy loop ! EN = E_MIN + FLOAT(IE - 1) * E_STEP ! E = hbar omega / E_F ! Z = FOURTH * EN / (X * X) ! Z = omega / omega_q ! IF(SQO_TYPE /= 'EPS') THEN ! CALL STFACT_DYNAMIC(X,Z,RS,T,SQO_TYPE,SQ_TYPE,SQ) ! ELSE CALL STFACT_DYNAMIC_FROM_EPS(X,Z,RS,T,SQ) ! END IF ! ! SQ_RN = SQ * EF_SI / H_BAR ! renormalized S(q,omega) ! IF(I_SF == 1) THEN ! writing to WRITE(IO_SF,*) Y,EN,SQ_RN,ZERO ! file END IF ! ! END DO ! ! END IF ! ! ! Format: ! 10 FORMAT(//,10X,'<<<<< ERROR IN CALCULATOR_3 :: CALC_SFC >>>>>', & /,10X,'<<<<< WRONG CHOICE OF SQ_TYPE :: GQ_TYPE >>>>>', & //,10X,' --> S(q) choice uses G(q) which uses S(q) <--',//) ! END SUBROUTINE CALC_SFC ! !======================================================================= ! SUBROUTINE CALC_PCF ! ! This subroutine computes the pair correlation function g(r) ! ! ! ! Intermediate parameters: ! ! * RS : Wigner-Seitz radius (in units of a_0) ! * T : temperature in SI ! * DMN : problem dimension ! ! ! Output parameters: ! ! * GR : g(r) ar r ! ! ! ! Author : D. Sébilleau ! ! Last modified : 29 Jul 2020 ! ! USE MATERIAL_PROP, ONLY : RS,DMN USE EXT_FIELDS, ONLY : T USE PC_VALUES, ONLY : GR_TYPE USE PD_VALUES, ONLY : RH_TYPE USE PAIR_CORRELATION ! USE R_GRID ! USE OUT_VALUES_3, ONLY : I_PC USE PRINT_FILES, ONLY : IO_PC ! IMPLICIT NONE ! REAL (WP) :: R,GR ! REAL (WP) :: FLOAT ! INTEGER :: IR ! DO IR = 1, N_R ! r loop ! R = R_MIN + FLOAT(IR - 1) * R_STEP ! r/a0 point ! IF(DMN == '3D') THEN ! CALL PAIR_CORRELATION_3D(R,RS,T,GR_TYPE,RH_TYPE,GR) ! ELSE IF(DMN == '2D') THEN ! CONTINUE ! ELSE IF(DMN == '1D') THEN ! CONTINUE ! END IF ! ! IF(I_PC == 1) THEN ! WRITE(IO_PC,*) R,GR ! END IF ! ! ENDDO ! ! END SUBROUTINE CALC_PCF ! !======================================================================= ! SUBROUTINE CALC_PDF ! ! This subroutine computes the pair distribution function rho2(r) ! ! ! ! Intermediate parameters: ! ! * RS : Wigner-Seitz radius (in units of a_0) ! * T : temperature in SI ! * DMN : problem dimension ! ! ! Output parameters: ! ! * R2 : rho2(r) ar r ! ! ! ! Author : D. Sébilleau ! ! Last modified : 29 Jul 2020 ! ! USE MATERIAL_PROP, ONLY : RS,DMN USE EXT_FIELDS, ONLY : T USE PD_VALUES, ONLY : RH_TYPE USE PAIR_DISTRIBUTION ! USE R_GRID ! USE OUT_VALUES_3, ONLY : I_P2 USE PRINT_FILES, ONLY : IO_P2 ! IMPLICIT NONE ! REAL (WP) :: R,R2 ! REAL (WP) :: FLOAT ! INTEGER :: IR ! DO IR = 1, N_R ! r loop ! R = R_MIN + FLOAT(IR - 1) * R_STEP ! r/a0 point ! IF(DMN == '3D') THEN ! CALL PAIR_DISTRIBUTION_3D(R,RS,T,RH_TYPE,R2) ! ELSE IF(DMN == '2D') THEN ! CONTINUE ! ELSE IF(DMN == '1D') THEN ! CONTINUE ! END IF ! ! IF(I_P2 == 1) THEN ! WRITE(IO_P2,*) R,R2 ! END IF ! ! END DO ! ! END SUBROUTINE CALC_PDF ! !======================================================================= ! SUBROUTINE CALC_VTX(X) ! ! This subroutine computes the vertex function Gamma(q, omega) ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Intermediate parameters: ! ! * RS : Wigner-Seitz radius (in units of a_0) ! * T : temperature in SI ! ! ! Output parameters: ! ! * E : energy array ! * VTXR : real part of vertex function ! * VTXI : imaginary part vertex function ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 30 Apr 2020 ! USE DIMENSION_CODE, ONLY : NSIZE USE REAL_NUMBERS, ONLY : ZERO,ONE,FOURTH,SMALL,INF USE COMPLEX_NUMBERS, ONLY : IC USE MATERIAL_PROP, ONLY : DMN,RS USE EXT_FIELDS, ONLY : T,H USE LF_VALUES, ONLY : GSTDY,GQ_TYPE,GQO_TYPE USE DF_VALUES, ONLY : ESTDY,EPS_T,D_FUNC ! USE E_GRID ! USE SCREENING_TYPE USE SCREENING_VEC ! USE LOCAL_FIELD_STATIC USE LOCAL_FIELD_DYNAMIC USE DFUNC_STATIC USE DFUNCT_STAN_DYNAMIC USE DFUNCL_STAN_DYNAMIC USE DFUNCL_MAGN_DYNAMIC ! USE OUT_VALUES_3, ONLY : I_VX USE PRINT_FILES, ONLY : IO_VX ! IMPLICIT NONE ! CHARACTER (LEN = 4) :: D_FUNCL,D_FUNCT ! REAL (WP), INTENT(IN) :: X REAL (WP) :: E(NSIZE) REAL (WP) :: LFCR(NSIZE),LFCI(NSIZE) REAL (WP) :: EPSR(NSIZE),EPSI(NSIZE) REAL (WP) :: VTXR,VTXI REAL (WP) :: RLFC,ILFC REAL (WP) :: REPS,IEPS REAL (WP) :: ETA REAL (WP) :: Z,EN REAL (WP) :: Y ! REAL (WP) :: A,NU,KS_SI ! REAL (WP) :: FLOAT ! COMPLEX (WP) :: EPS,LFC,VTX ! INTEGER :: IE ! A = ZERO ! temporary NU = ZERO ! ! Y = X + X ! q / k_F ! ! Computing the screening wave vector ! CALL SCREENING_VECTOR(SC_TYPE,DMN,X,RS,T,KS_SI) ! ! 1) Computing the local field correction ! IF(GSTDY == ' STATIC') THEN ! ! ! Static local field correction ! CALL LFIELD_STATIC(X,RS,T,GQ_TYPE,RLFC) ! ! LFCR(1) = RLFC ! LFCI(1) = ZERO ! ! ELSE ! ! ! Dynamic local field correction ! CALL LFIELD_DYNAMIC(X,RS,E_MIN,E_MAX,N_E,T,ETA, & ! GQO_TYPE,E,LFCR,LFCI) ! ! END IF ! ! ! 2) Computing the dielectric function ! IF(ESTDY == ' STATIC') THEN ! ! ! Static dielectric function ! IF(EPS_T == 'LONG') THEN ! ! D_FUNCL = 'RPA1' ! CALL DFUNCL_STATIC(X,D_FUNCL,REPS,IEPS) ! longitudinal eps EPSR(1) = REPS ! EPSI(1) = IEPS ! ! ELSE ! D_FUNCT = D_FUNC ! CONTINUE ! transverse eps END IF ! ! ELSE ! ! ! Dynamic dielectric function ! DO IE = 1, N_E ! energy loop ! EN = E_MIN + FLOAT(IE - 1) * E_STEP ! EN = hbar omega / E_F ! Z = FOURTH * EN / (X * X) ! Z = omega / omega_q ! IF(EPS_T == 'LONG') THEN ! longitudinal eps ! D_FUNCL = 'RPA1' ! IF(H < SMALL) THEN ! CALL DFUNCL_DYNAMIC(X,Z,RS,T,D_FUNCL,IE,REPS,IEPS) ! no magnetic field ELSE ! CALL DFUNCL_DYNAMIC_M(X,Z,KS_SI,A,NU,D_FUNCL,REPS,IEPS) ! magnetic field END IF ! ELSE ! transverse eps D_FUNCT = D_FUNC ! IF(H < SMALL) THEN ! CALL DFUNCT_DYNAMIC(X,Z,D_FUNCT,REPS,IEPS) ! no magnetic field ELSE ! CONTINUE ! magnetic field END IF ! END IF ! ! EPSR(IE) = REPS ! EPSI(IE) = IEPS ! E(IE) = EN ! ! END DO ! end of energy loop ! END IF ! ! ! 3) Computing the vertex function and writing to file ! ! As we have chosen D_FUNCL='RPA1', PI = PI0 and the ! vertex function is: ! ! 1 ! --------------------- ! 1 + LFC * (EPS - 1) ! ! IF(I_VX == 1) THEN ! DO IE = 1, N_E ! ! LFC = LFCR(IE) + IC * LFCI(IE) ! EPS = EPSR(IE) + IC * EPSI(IE) ! VTX = ONE / (ONE + LFC * (EPS - ONE)) ! VTXR = REAL(VTX,KIND=WP) ! VTXI = AIMAG(VTX) ! WRITE(IO_VX,*) Y,E(IE),VTXR,VTXI ! ! END DO ! END IF ! ! END SUBROUTINE CALC_VTX ! !======================================================================= ! SUBROUTINE CALC_DMP(IX,X) ! ! This subroutine computes the plasmon damping coefficient gamma_q ! ! ! Im[ epsilon ] | ! gamma_q = - _______________ | ! | ! d Re[ epsilon ]/d omega | omega=Omega_q ! ! where epsilon is the dielectric function. ! ! ! Input parameters: ! ! * IX : index of X-point ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Intermediate parameters: ! ! * RS : Wigner-Seitz radius (in units of a_0) ! * T : temperature in SI ! ! ! Output parameters: ! ! * E : energy array ! * VTXR : real part of vertex function ! * VTXI : imaginary part vertex function ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 30 Apr 2021 ! ! USE DIMENSION_CODE, ONLY : NSIZE USE MATERIAL_PROP, ONLY : RS USE EXT_FIELDS, ONLY : T ! USE REAL_NUMBERS, ONLY : ZERO,TWO,FOURTH,TTINY,INF USE FERMI_SI, ONLY : KF_SI ! USE E_GRID ! USE DF_VALUES, ONLY : D_FUNC ! USE INTEGRATION, ONLY : INTEGR_L USE DFUNCL_STAN_DYNAMIC ! USE FIND_ZERO, ONLY : FIND_ZERO_FUNC USE PLASMON_DAMPING, ONLY : EXACT_DAMPING ! USE OUT_VALUES_3, ONLY : I_DC USE PRINT_FILES, ONLY : IO_DC ! IMPLICIT NONE ! INTEGER,INTENT(IN) :: IX INTEGER :: IE INTEGER :: IDERIV ! REAL (WP), INTENT(IN) :: X REAL (WP) :: D,TAU REAL (WP) :: Y,E,V,Z REAL (WP) :: REPS,IEPS REAL (WP) :: EN(NSIZE) REAL (WP) :: EPSR(NSIZE),EPSI(NSIZE) REAL (WP) :: ZEROF,GAMMA_Q ! REAL (WP) :: FLOAT ! IDERIV = 5 ! ! Y = X + X ! q/k_F ! ! Constructing the e-grid ! DO IE = 1, N_E ! E_F ! E = E_MIN + FLOAT(IE - 1) * E_STEP ! in units of V = E ! hbar * omega / E_F Z = FOURTH * V / (X * X) ! omega / omega_q ! ! Computing the dielectric function epsilon(q,E) ! CALL DFUNCL_DYNAMIC(X,Z,RS,T,D_FUNC,IE,REPS,IEPS) ! ! EN(IE) = E ! EPSR(IE) = REPS ! EPSI(IE) = IEPS ! ! END DO ! ! ! Computing the plasmon energy at q: ! Only the highest energy zero is kept ! CALL FIND_ZERO_FUNC(EN,EPSR,N_E,ZEROF) ! ! ! Computing the damping coefficient ! CALL EXACT_DAMPING(IX,IDERIV,N_E,EN,EPSR,EPSI,ZEROF,GAMMA_Q) ! ! IF(I_DC == 1) THEN ! WRITE(IO_DC,*) Y,GAMMA_Q ! END IF ! ! END SUBROUTINE CALC_DMP ! !======================================================================= ! SUBROUTINE CALC_QBD ! ! This subroutine computes the plasmon q-bounds ! ! ! Output parameters: ! ! * Q_MIN : plasmon lower q-bound ! * Q_MAX : plasmon upper q-bound ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 6 Oct 2020 ! ! USE Q_BOUNDS ! USE OUT_VALUES_3, ONLY : I_QC USE PRINT_FILES, ONLY : IO_QC ! IMPLICIT NONE ! REAL (WP) :: Q_MIN,Q_MAX ! CALL QBOUNDS(Q_MIN,Q_MAX) ! IF(I_QC == 1) THEN ! WRITE(IO_QC,*) Q_MIN,Q_MAX ! END IF ! ! END SUBROUTINE CALC_QBD ! !======================================================================= ! SUBROUTINE CALC_RLX(X) ! ! This subroutine computes the relaxation time ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Output parameters: ! ! * TAU : relaxation time ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 7 Oct 2020 ! ! USE MATERIAL_PROP, ONLY : DMN,RS USE EXT_FIELDS, ONLY : T USE PLASMA, ONLY : ZION ! USE SCREENING_VEC, ONLY : DEBYE_VECTOR USE DAMPING_VALUES, ONLY : VI_TYPE USE RELAXATION_TIME_STATIC USE VISCOSITY ! USE EL_PHO_INTER ! USE REAL_NUMBERS, ONLY : ZERO ! USE OUT_VALUES_3, ONLY : I_RL USE PRINT_FILES, ONLY : IO_RL ! IMPLICIT NONE ! REAL (WP), INTENT(IN) :: X REAL (WP) :: Y,ETA,KD_SI,TAU REAL (WP) :: LR,S_L ! Y = X + X ! q/k_F ! LR = ZERO ! residual mfp (temporary) S_L = ZERO ! scattering length (temporary) ! ! Computing the Debye momentum ! CALL DEBYE_VECTOR('3D',T,RS,KD_SI) ! ! ! Computation of the viscosity ! IF(DMN == '3D') THEN ! CALL VISCOSITY_3D(RS,T,ZION,KD_SI,X,ZERO,NA,MA,RA, & ! DEBYE_T,EP_C,LR,VI_TYPE,ETA) ! ELSE IF(DMN == '2D') THEN ! CALL VISCOSITY_2D(T,S_L,VI_TYPE,ETA) ! ELSE IF(DMN == '1D') THEN ! ETA = ZERO ! not yet implemented END IF ! ! ! Computing the relaxation time ! CALL RELAXATION_TIME(X,TAU) ! ! IF(I_RL == 1) THEN ! WRITE(IO_RL,*) Y,TAU ! END IF ! ! END SUBROUTINE CALC_RLX ! !======================================================================= ! SUBROUTINE CALC_SCR(X) ! ! This subroutine computes the screening wave number ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Output parameters: ! ! * KS : screening wave number ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 13 Oct 2020 ! ! USE MATERIAL_PROP, ONLY : DMN,RS USE EXT_FIELDS, ONLY : T ! USE CONSTANTS_P1, ONLY : BOHR USE SCREENING_TYPE USE SCREENING_VEC USE SCREENING_VEC2 ! USE OUT_VALUES_3, ONLY : I_KS USE PRINT_FILES, ONLY : IO_KS ! IMPLICIT NONE ! REAL (WP), INTENT(IN) :: X REAL (WP) :: Y REAL (WP) :: KS_SI,KS ! Y = X + X ! q / k_F ! IF(SC_TYPE == 'DH' .OR. SC_TYPE == 'TF') THEN ! CALL SCREENING_VECTOR(SC_TYPE,DMN,X,RS,T,KS_SI) ! ELSE ! CALL SCREENING_VECTOR2(SC_TYPE,DMN,X,RS,T,KS_SI) ! END IF ! ! KS = KS_SI * BOHR ! KS in units of a_0 ! IF(I_KS == 1) THEN ! WRITE(IO_KS,*) Y,KS ! END IF ! ! END SUBROUTINE CALC_SCR ! !======================================================================= ! SUBROUTINE CALC_QVF(X) ! ! This subroutine computes the omega = q * v_F (U = 1) equation, as well ! as U + X = 1 and U - X = 1 ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Output parameters: ! ! * OM : omega ! ! ! hbar omega ! Note: We use here the fact that V = ------------ = 4 * U * X ! E_F ! ! ! Therefore U = 1 <==> V = 4 * X ! U + X = 1 <==> V = 4 * (1 - X) * X ! U - X = 1 <==> V = 4 * (1 + X) * X ! ! ! Author : D. Sébilleau ! ! Last modified : 25 Nox 2020 ! ! USE REAL_NUMBERS, ONLY : ONE,FOUR ! USE OUT_VALUES_3, ONLY : I_OQ USE PRINT_FILES, ONLY : IO_OQ ! IMPLICIT NONE ! REAL (WP), INTENT(IN) :: X REAL (WP) :: Y,OM0,OM1,OM2 ! Y = X + X ! q / k_F ! OM0 = FOUR * X ! OM1 = FOUR * (ONE - X) * X ! OM2 = FOUR * (ONE + X) * X ! ! IF(I_OQ == 1) THEN ! WRITE(IO_OQ,*) Y,OM0,OM1,OM2 ! END IF ! ! END SUBROUTINE CALC_QVF ! !======================================================================= ! SUBROUTINE CALC_MEP(X) ! ! This subroutine computes the moments of the ! imaginary part of the dielectric function ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Output parameters: ! ! * MEP : moment of S(q,omega) function in reduced units ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 21 Oct 2020 ! ! USE MOMENTS USE MOMENTS_CALC ! USE OUT_VALUES_3, ONLY : I_ME USE PRINT_FILES, ONLY : IO_ME ! IMPLICIT NONE ! REAL (WP), INTENT(IN) :: X REAL (WP) :: Y,MEP ! Y = X + X ! q/k_F ! CALL MOMENTS_EPSILON(X,N_M,MEP) ! ! IF(I_ME == 1) THEN ! writing to WRITE(IO_ME,*) Y,MEP ! file END IF ! ! END SUBROUTINE CALC_MEP ! !======================================================================= ! SUBROUTINE CALC_MSF(X) ! ! This subroutine computes the moments of the ! dynamical structure factor ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Output parameters: ! ! * MSF : moment of S(q,omega) function in reduced units ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 20 Oct 2020 ! ! USE MOMENTS USE MOMENTS_CALC ! USE OUT_VALUES_3, ONLY : I_MS USE PRINT_FILES, ONLY : IO_MS ! IMPLICIT NONE ! REAL (WP), INTENT(IN) :: X REAL (WP) :: Y,MSF ! Y = X + X ! q/k_F ! CALL MOMENTS_STRUCT_FACTOR(X,N_M,MSF) ! ! IF(I_MS == 1) THEN ! writing to WRITE(IO_MS,*) Y,MSF ! file END IF ! ! END SUBROUTINE CALC_MSF ! !======================================================================= ! SUBROUTINE CALC_MLO(X) ! ! This subroutine computes the moments of the loss function ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Output parameters: ! ! * MLO : moment of loss function in reduced units ! ! ! ! ! Author : D. Sébilleau ! ! Last modified : 20 Oct 2020 ! ! USE MOMENTS USE MOMENTS_CALC ! USE OUT_VALUES_3, ONLY : I_ML USE PRINT_FILES, ONLY : IO_ML ! IMPLICIT NONE ! REAL (WP), INTENT(IN) :: X REAL (WP) :: Y,MLO ! Y = X + X ! q/k_F ! CALL MOMENTS_LOSS_FUNCTION(X,N_M,MLO) ! ! IF(I_ML == 1) THEN ! writing to WRITE(IO_ML,*) Y,MLO ! file END IF ! ! END SUBROUTINE CALC_MLO ! !======================================================================= ! SUBROUTINE CALC_RE0(X) ! ! This subroutine computes the zeros of the real part of ! the dielectric function: Re[ epsilon(q,omega) ] ! ! ! Input parameters: ! ! * X : dimensionless factor --> X = q / (2 * k_F) ! ! ! Intermediate parameters: ! ! * RS : Wigner-Seitz radius (in units of a_0) ! * T : temperature in SI ! ! ! Output parameters: ! ! * Y : q / k_F ! * RE0 : zeros of real part dielectric function ! ! ! Note: By setting a non-zero value to variable SHF, the subroutine ! will solve Re[ epsilon(q,omega) ] = - SHF ! ! ! Author : D. Sébilleau ! ! Last modified : 30 Apr 2021 ! ! USE MATERIAL_PROP, ONLY : DMN,RS USE EXT_FIELDS, ONLY : T,H USE DF_VALUES, ONLY : ESTDY,EPS_T,D_FUNC ! USE DIMENSION_CODE, ONLY : NSIZE USE REAL_NUMBERS, ONLY : ZERO,FOURTH,SMALL,INF ! USE E_GRID ! USE SCREENING_TYPE USE SCREENING_VEC ! USE DFUNC_STATIC USE DFUNCT_STAN_DYNAMIC USE DFUNCL_STAN_DYNAMIC USE DFUNCL_MAGN_DYNAMIC ! USE FIND_ZERO ! USE OUT_VALUES_3, ONLY : I_ZE USE PRINT_FILES, ONLY : IO_ZE ! IMPLICIT NONE ! CHARACTER (LEN = 4) :: D_FUNCL,D_FUNCT ! REAL (WP), INTENT(IN) :: X REAL (WP) :: EN REAL (WP) :: Y,Z REAL (WP) :: REPS,IEPS REAL (WP) :: E(NSIZE) REAL (WP) :: EPSR(NSIZE) REAL (WP) :: ZEROF ! REAL (WP) :: A,NU,KS_SI REAL (WP) :: LFT,TAU,DR,D,ETA ! REAL (WP), PARAMETER :: SHF = 0.0E0_WP ! REAL (WP) :: FLOAT ! INTEGER :: IE ! A = ZERO ! temporary NU = ZERO ! ! Y = X + X ! q / k_F ! ! Computing the screening wave vector ! CALL SCREENING_VECTOR(SC_TYPE,DMN,X,RS,T,KS_SI) ! ! Computing the dynamic dielectric function ! DO IE = 1, N_E ! energy loop ! EN = E_MIN + FLOAT(IE - 1) * E_STEP ! EN = hbar omega / E_F ! Z = FOURTH * EN / (X * X) ! Z = omega / omega_q ! IF(EPS_T == 'LONG') THEN ! longitudinal eps ! D_FUNCL = D_FUNC ! IF(H < SMALL) THEN ! CALL DFUNCL_DYNAMIC(X,Z,RS,T,D_FUNCL,IE,REPS,IEPS) ! no magnetic field ELSE ! CALL DFUNCL_DYNAMIC_M(X,Z,KS_SI,A,NU,D_FUNCL,REPS,IEPS) ! magnetic field END IF ! ELSE ! transverse eps D_FUNCT = D_FUNC ! IF(H < SMALL) THEN ! CALL DFUNCT_DYNAMIC(X,Z,D_FUNCT,REPS,IEPS) ! no magnetic field ELSE ! CONTINUE ! magnetic field END IF ! END IF ! ! EPSR(IE) = REPS + SHF ! E(IE) = EN ! ! END DO ! end of energy loop ! ! Computing the zeros of EPSR ! IF(I_ZE == 1) THEN ! CALL PRINT_ZERO_FUNC(Y,E,EPSR,N_E) ! END IF ! ! END SUBROUTINE CALC_RE0 ! !======================================================================= ! SUBROUTINE CALC_MFP ! ! This subroutine computes the inelastic mean free path ! ! ! ! Author : D. Sébilleau ! ! Last modified : 19 Oct 2020 ! ! USE REAL_NUMBERS, ONLY : SMALL USE FERMI_SI, ONLY : EF_SI USE ENE_CHANGE, ONLY : EV,ANG ! USE ELECTRON_MEAN_FREE_PATH USE IMFP ! USE OUT_VALUES_3, ONLY : I_MF USE PRINT_FILES, ONLY : IO_MF ! IMPLICIT NONE ! INTEGER :: IE INTEGER, PARAMETER :: NE_MAX = 1480 ! max. number of energy points ! REAL (WP) :: LAMBDA REAL (WP) :: E_STEP REAL (WP) :: E,EK REAL (WP) :: E_MIN,E_MAX ! REAL (WP) :: FLOAT ! E_MIN = EK_INI ! lower value in eV E_MAX = EK_FIN ! upper value in eV ! E_STEP = (E_MAX - E_MIN) / FLOAT(NE_MAX - 1) ! e-step in eV ! DO IE = 1, NE_MAX ! ! E = E_MIN + FLOAT(IE - 1) * E_STEP ! E in eV EK = E * EV ! E in SI ! CALL MEAN_FREE_PATH(EK,LAMBDA) ! IMFP in SI LAMBDA = LAMBDA / ANG ! IMFP in Angström ! IF(I_MF == 1) THEN ! WRITE(IO_MF,*) E,LAMBDA ! END IF ! ! END DO ! end of energy loop ! END SUBROUTINE CALC_MFP ! !======================================================================= ! SUBROUTINE CALC_NEV(X) ! ! This subroutine computes the Fourier space Nevalinna/memory function ! ! ! Author : D. Sébilleau ! ! Last modified : 28 Jan 2021 ! ! USE MATERIAL_PROP, ONLY : RS USE EXT_FIELDS, ONLY : T ! USE E_GRID ! USE REAL_NUMBERS, ONLY : FOURTH ! USE RELAXATION_TIME_STATIC USE MEMORY_FUNCTIONS_F USE NEVALINNA_FUNCTIONS ! USE DF_VALUES, ONLY : D_FUNC,NEV_TYPE,MEM_TYPE,ALPHA,BETA USE DAMPING_SI USE DAMPING_VALUES, ONLY : PCT ! USE OUT_VALUES_3, ONLY : I_NV USE PRINT_FILES, ONLY : IO_NV ! IMPLICIT NONE ! INTEGER :: IE,I_F ! REAL (WP), INTENT(IN) :: X ! REAL (WP) :: Y,Z REAL (WP) :: E REAL (WP) :: NEVR,NEVI ! REAL (WP) :: FLOAT,REAL,AIMAG ! COMPLEX (WP) :: FUNC ! Y = X + X ! q / k_F ! ! Check for Nevanlinna or memory function --> switch I_F ! IF(D_FUNC(1:3) == 'NEV') THEN ! I_F = 1 ! ELSE IF(D_FUNC(1:3) == 'MEM') THEN ! I_F = 2 ! ELSE ! I_F = 0 ! END IF ! ! DO IE = 1, N_E ! energy loop ! E = E_MIN + FLOAT(IE - 1) * E_STEP ! E = hbar omega / E_F ! Z = FOURTH * E / (X * X) ! Z = omega / omega_q ! IF(I_F == 1) THEN ! FUNC = NEVAN2(X,Z,RS,T,TAU,NEV_TYPE) ! ELSE IF(I_F == 2) THEN ! FUNC = MEMORY_F(E,TAU,TAU2,PCT,ALPHA,BETA,MEM_TYPE) ! END IF ! ! NEVR = REAL(FUNC,KIND=WP) ! NEVI = AIMAG(FUNC) ! ! IF(I_NV == 1) THEN ! WRITE(IO_NV,*) Y,E,NEVR,NEVI ! END IF ! ! END DO ! end of energy loop ! END SUBROUTINE CALC_NEV ! !======================================================================= ! SUBROUTINE CALC_MEM ! ! This subroutine computes the time-domain memory function ! as a function of t / tau ! ! When tau is q-dependent, the only vamue considered is the first one ! ! Author : D. Sébilleau ! ! Last modified : 28 Jan 2021 ! ! USE REAL_NUMBERS, ONLY : ZERO,TEN,HALF ! USE Q_GRID, ONLY : Q_MIN ! USE DF_VALUES, ONLY : ALPHA,BETA,MEM_TYPE USE DAMPING_SI USE DAMPING_VALUES, ONLY : PCT ! USE RELAXATION_TIME_STATIC USE MEMORY_FUNCTIONS_T ! USE OUT_VALUES_3, ONLY : I_MT USE PRINT_FILES, ONLY : IO_MT ! IMPLICIT NONE ! INTEGER :: IT INTEGER, PARAMETER :: NT_MAX = 200 ! max. number of time points ! REAL (WP) :: MEMR REAL (WP) :: T ! t / tau REAL (WP) :: T1 ! t REAL (WP) :: T_MIN,T_MAX,T_STEP REAL (WP) :: X ! REAL (WP) :: FLOAT ! T_MIN = ZERO ! lower time value T_MAX = TEN ! upper time value X = Q_MIN * HALF ! initial value of q / 2 k_F ! T_STEP = (T_MAX - T_MIN) / FLOAT(NT_MAX - 1) ! t-step in units of tau ! DO IT = 1, NT_MAX ! ! T = T_MIN + FLOAT(IT - 1) * T_STEP ! t in units of tau ! T1 = T * TAU ! time is SI ! MEMR = MEMORY_T(T1,TAU,TAU2,PCT,ALPHA,BETA,MEM_TYPE) * & ! TAU * TAU ! in units of 1 / tau^2 ! IF(I_MT == 1) THEN ! WRITE(IO_MT,*) T,MEMR ! END IF ! ! END DO ! end of time loop ! END SUBROUTINE CALC_MEM ! END MODULE CALCULATORS_3