MsSpec-DFM/New_libraries/DFM_library/CONFINEMENT_LIBRARY/confinement_ff.f90

649 lines
23 KiB
Fortran
Raw Blame History

This file contains ambiguous Unicode characters

This file contains Unicode characters that might be confused with other characters. If you think that this is intentional, you can safely ignore this warning. Use the Escape button to reveal them.

!
!=======================================================================
!
MODULE CONFINEMENT_FF
!
USE ACCURACY_REAL
!
CONTAINS
!
!=======================================================================
!
FUNCTION CONFIN_FF(X)
!
! This function computes the confinement form factor
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
!
!
! Intermediate parameters:
!
! * CONFIN : type of confinement
! CONFIN = 'NO-CONF'
! CONFIN = 'INF_QWW'
! CONFIN = 'CC-1111' cylindrical within subband 1
! CONFIN = 'CC-1122' cylindrical between subbands 1 and 2
! CONFIN = 'CC-1221' cylindrical between subbands 1 and 2
! CONFIN = 'CC-2222' cylindrical within subband 2
! CONFIN = 'HC-1111' harmonic within subband 1
! CONFIN = 'HC-1122' harmonic between subbands 1 and 2
! CONFIN = 'HC-1221' harmonic between subbands 1 and 2
! CONFIN = 'HC-2222' harmonic within subband 2
! CONFIN = 'INVLAYE'
! CONFIN = 'IQWE_LB'
! CONFIN = 'PC1_QWI'
! CONFIN = 'PC2_QWI'
! CONFIN = 'SWC_QWI'
! * R0 : radius of the quantum wire in SI --> m
! * L : length of the quantum well (SI)
! * DL : distance between the stacked layers (SI)
!
!
! Output parameters:
!
! * CONFIN_FF : form factor
!
!
!
! Author : D. Sébilleau
!
! Last modified : 3 Aug 2020
!
!
USE REAL_NUMBERS, ONLY : ONE
USE CONFIN_VAL
USE MULTILAYER
!
IMPLICIT NONE
!
REAL (WP) :: CONFIN_FF,X
!
IF(CONFIN == 'NO-CONF') THEN !
CONFIN_FF=ONE !
ELSE IF(CONFIN == 'DSEPLAY') THEN !
CONFIN_FF=DSEPLAY_FF(X,DL) !
ELSE IF(CONFIN == 'CC-1111') THEN !
CONFIN_FF=INF_QWW_FF(X,R0,CONFIN) !
ELSE IF(CONFIN == 'CC-1122') THEN !
CONFIN_FF=INF_QWW_FF(X,R0,CONFIN) !
ELSE IF(CONFIN == 'CC-1221') THEN !
CONFIN_FF=INF_QWW_FF(X,R0,CONFIN) !
ELSE IF(CONFIN == 'CC-2222') THEN !
CONFIN_FF=INF_QWW_FF(X,R0,CONFIN) !
ELSE IF(CONFIN == 'HC-1111') THEN !
CONFIN_FF=HCM_QWW_FF(X,R0,CONFIN) !
ELSE IF(CONFIN == 'HC-1122') THEN !
CONFIN_FF=HCM_QWW_FF(X,R0,CONFIN) !
ELSE IF(CONFIN == 'HC-1221') THEN !
CONFIN_FF=HCM_QWW_FF(X,R0,CONFIN) !
ELSE IF(CONFIN == 'HC-2222') THEN !
CONFIN_FF=HCM_QWW_FF(X,R0,CONFIN) !
ELSE IF(CONFIN == 'INVLAYE') THEN !
CONFIN_FF=INVLAYE_FF(X,DL,N_DEP,N_INV,EPS_2,EPS_1) !
ELSE IF(CONFIN == 'IQWE_LB') THEN !
CONFIN_FF=IQWE_LB_FF(X,L) !
ELSE IF(CONFIN == 'PC1_QWI') THEN !
CONFIN_FF=PC1_QWI_FF(X,OM0) !
ELSE IF(CONFIN == 'PC2_QWI') THEN !
CONFIN_FF=PC2_QWI_FF(X,OM0) !
ELSE IF(CONFIN == 'SOF_COR') THEN !
CONFIN_FF=SOF_COR_FF(X,R0) !
ELSE IF(CONFIN == 'SWC_QWI') THEN !
CONFIN_FF=SWC_QWI_FF(X,L) !
END IF !
!
END FUNCTION CONFIN_FF
!
!=======================================================================
!
FUNCTION DSEPLAY_FF(X,D)
!
! This function computes the form factor of a 2D layer within a
! stacking of layers separated by a distance D along the z axis
!
! Reference: (1) A. C. Sharma, Solid State Comm. 70, 1171-1174 (1989)
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * D : distance between the stacked layers (SI)
!
! Output parameters:
!
! * DSEPLAY_FF : form factor
!
!
! Method: from eq. (7a) and (7b) adapted for SI, we have
!
! e^2
! EPS(q,omega,qz) = 1 - --------- S(q,qz) * PI(q,omega)
! q EPS_0
!
! Applying eq. (10), we find that the Coulomb potential
! for intralayer interaction must be
!
! / + pi/d
! e^2 d |
! VC(q) = --------- * ------ | S(q,qz) dqz
! q EPS_0 2 pi |
! / - pi/d
!
! Author : D. Sébilleau
!
! Last modified : 10 Jun 2020
!
!
USE REAL_NUMBERS, ONLY : TWO
USE FERMI_SI, ONLY : KF_SI
USE SPECIFIC_INT_1, ONLY : SQQZ_INT
!
IMPLICIT NONE
!
REAL (WP) :: X,D
REAL (WP) :: DSEPLAY_FF
REAL (WP) :: Q
!
Q=TWO*X*KF_SI !
!
DSEPLAY_FF=SQQZ_INT(Q,D) !
!
END FUNCTION DSEPLAY_FF
!
!=======================================================================
!
FUNCTION INF_QWW_FF(X,R0,CONFIN)
!
! This function computes the form factor of the a quantum-well wire
! with an infinite barrier when only the two lower subbands are filled
!
! Reference: (1) A. Gold and A. Ghazali, Phys. Rev. B 41, 7626 (1990)
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * R0 : radius of the quantum wire in SI --> m
! * CONFIN : type of confinement
! CONFIN = 'CC-1111' cylindrical within subband 1
! CONFIN = 'CC-1122' cylindrical between subbands 1 and 2
! CONFIN = 'CC-1221' cylindrical between subbands 1 and 2
! CONFIN = 'CC-2222' cylindrical within subband 2
!
! Output parameters:
!
! * INF_QWW_FF : form factor
!
! Author : D. Sébilleau
!
! Last modified : 10 Jun 2020
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,SIX,EIGHT,THIRD
USE FERMI_SI, ONLY : KF_SI
USE BESSEL
!
IMPLICIT NONE
!
CHARACTER*7 CONFIN
!
REAL (WP) :: X,R0
REAL (WP) :: INF_QWW_FF
REAL (WP) :: Q,QR0,QR02,QR04,QR06
REAL (WP) :: I3,K3,I4,K4
!
Q=TWO*X*KF_SI !
QR0=Q*R0 !
QR02=QR0*QR0 !
QR04=QR02*QR02 !
QR06=QR04*QR02 !
!
! Computing the Bessel functions
!
I3=BESSI(3,QR0) !
K3=BESSK(3,QR0) !
I4=BESSI(4,QR0) !
K4=BESSK(4,QR0) !
!
IF(CONFIN == 'CC-1111') THEN !
INF_QWW_FF=72.0E0_WP*( 0.10E0_WP - TWO*THIRD/QR02 + & !
32.0E0_WP*THIRD/QR04 - 64.0E0_WP*I3*K3/QR04 & ! ref. 1 eq. (11a)
) / QR02 !
ELSE IF(CONFIN == 'CC-1122') THEN !
INF_QWW_FF=288.0E0_WP*( 0.050E0_WP*THIRD - & !
0.20E0_WP*THIRD/QR02 + & !
EIGHT*THIRD/QR04 - & ! ref. 1 eq. (11b)
64.0E0_WP*K3/QR04 * & !
(I3 - SIX*I4/QR0) & !
) / QR02 !
ELSE IF(CONFIN == 'CC-1221') THEN !
INF_QWW_FF=576.0E0_WP*( 0.050E0_WP*THIRD - & !
0.80E0_WP*THIRD/QR02 + & !
EIGHT/QR04 - 64.0E0_WP*I4*K4/QR04 & ! ref. 1 eq. (11c)
) / QR02 !
ELSE IF(CONFIN == 'CC-2222') THEN !
INF_QWW_FF=1152.0E0_WP*( ONE/210.0E0_WP - & !
0.20E0_WP*THIRD/QR02 + & !
12.80E0_WP*THIRD/QR04 + & !
96.0E0_WP/QR06 - & !
64.0E0_WP/QR04 * & !
(I3 - SIX*I4/QR0) * & ! ref. 1 eq. (11d)
(K3 + SIX*K4/QR0) & !
) / QR02 !
END IF !
!
END FUNCTION INF_QWW_FF
!
!=======================================================================
!
FUNCTION HCM_QWW_FF(X,B,CONFIN)
!
! This function computes the form factor of the a quantum-well wire
! with a harmonic confinement when only the two lower subbands are filled
!
! Reference: (1) G. Y. Hu and R. F. O'Connell,
! J. Phys.: Condens. Matter 2 9381 (1990)
!
! Note: Here, the two sub-bands are labelled 1 and 2 instead of 0 and 1
! in ref. (1). This is for consistency with function INF_QWW_FF
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * B : confinement parameter (b = sqrt(h_bar / m* omega_0)
! * CONFIN : type of confinement
! CONFIN = 'HC-1111' harmonic within subband 1
! CONFIN = 'HC-1122' harmonic between subbands 1 and 2
! CONFIN = 'HC-1221' harmonic between subbands 1 and 2
! CONFIN = 'HC-2222' harmonic within subband 2
!
! Output parameters:
!
! * HCM_QWW_FF : form factor
!
! Author : D. Sébilleau
!
! Last modified : 3 Aug 2020
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,HALF,FOURTH,EIGHTH
USE FERMI_SI, ONLY : KF_SI
USE BESSEL, ONLY : BESSK0,BESSK1
!
IMPLICIT NONE
!
CHARACTER*7 CONFIN
!
REAL (WP) :: X,B
REAL (WP) :: HCM_QWW_FF
REAL (WP) :: Q,Q2B2,COEF
REAL (WP) :: K0,K1
!
REAL (WP) :: DEXP
!
Q=TWO*X*KF_SI !
Q2B2=(Q*B)*Q*B !
!
COEF=DEXP(FOURTH*Q2B2) !
!
! Computing the Bessel functions
!
K0=BESSK0(FOURTH*Q2B2) !
K1=BESSK1(FOURTH*Q2B2) !
!
IF(CONFIN == 'HC-1111') THEN !
HCM_QWW_FF=COEF*K0 ! ref. (1) eq. (2.5a)
ELSE IF(CONFIN == 'HC-1122') THEN !
HCM_QWW_FF=COEF*( K0 + FOURTH*Q2B2*(K0-K1) ) ! ref. (1) eq. (2.5c)
ELSE IF(CONFIN == 'HC-1221') THEN !
HCM_QWW_FF=COEF*FOURTH*Q2B2*(K1-K0) ! ref. (1) eq. (2.5d)
ELSE IF(CONFIN == 'HC-2222') THEN !
HCM_QWW_FF=COEF*( (ONE+HALF*Q2B2+EIGHTH*Q2B2*Q2B2)*K0 - & !
FOURTH*Q2B2*(ONE+HALF*Q2B2)*K1 & ! ref. (1) eq. (2.5b)
) !
END IF !
!
END FUNCTION HCM_QWW_FF
!
!=======================================================================
!
FUNCTION INVLAYE_FF(X,D,N_DEP,N_INV,EPSO,EPSS)
!
! This function computes the form factor of the surface inversion layer
! of a semiconductor
!
! Reference: (1) T. K. Lee, C. S. Tin and J. J. Quinn,
! Solid. State Comm. 16, 1309-1312 (1975)
! (2) M. Jonson, J. Phys. C: Solid State Phys. 9, 3055-3071 (1976)
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * D : width of the insulation layer (SI)
! * N_DEP : electron concentration in depletion layer (SI)
! * N_INV : electron concentration in inversion layer (SI)
! * EPSO : dielectric constant of the insulating layer (oxide)
! * EPSS : dielectric constant of the semiconductor
!
! Output parameters:
!
! * INVLAYE_FF : form factor
!
! Author : D. Sébilleau
!
! Last modified : 10 Jun 2020
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,THREE,THIRD
USE CONSTANTS_P1, ONLY : H_BAR,E,M_E
USE FERMI_SI, ONLY : KF_SI
USE PI_ETC, ONLY : PI
!
IMPLICIT NONE
!
REAL (WP) :: X,D,N_DEP,N_INV,EPSO,EPSS
REAL (WP) :: INVLAYE_FF
REAL (WP) :: B
REAL (WP) :: Q,NU1,NU2,NUM,DEN
REAL (WP) :: Y,Y2,Y3,Y4
!
REAL (WP) :: DTANH
!
! Computation of parameter B
!
NUM=48.0E0_WP*PI*E*E*M_E !
DEN=EPSS*H_BAR*H_BAR !
!
B=( NUM*(N_DEP+11.0E0_WP*N_INV/32.0E0_WP)/DEN )**THIRD ! ref. 2 eq. (21)
!
Q=TWO*X*KF_SI !
Y=Q/B !
Y2=Y*Y !
Y3=Y2*Y !
Y4=Y3*Y !
!
NU1=0.125E0_WP*Y*( 33.0E0_WP + 54.0E0_WP*Y + 44.0E0_WP*Y2 + & !
18.0E0_WP*Y3+THREE*Y4 & !
) !
NU2=TWO*EPSS/( EPSS+EPSO/DTANH(Q*D) ) !
DEN=(ONE+Y)**6 !
!
INVLAYE_FF=(NU1+NU2)/DEN ! ref. 1 eq. (5)
!
END FUNCTION INVLAYE_FF
!
!=======================================================================
!
FUNCTION IQWE_LB_FF(X,L)
!
! This function computes the form factor of the a quantum well
! with an infinite barrier when only the lower subband is filled
!
! Reference: (1) T. Vazifehshenas and T. Salavati-fard,
! Physica E 41, 12971300 (2009)
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * L : length of the quantum well (SI)
!
! Output parameters:
!
! * IQWE_LB_FF : form factor
!
! Author : D. Sébilleau
!
! Last modified : 10 Jun 2020
!
!
USE REAL_NUMBERS, ONLY : ONE,TWO,THREE,FOUR,EIGHT
USE FERMI_SI, ONLY : KF_SI
USE PI_ETC, ONLY : PI2
!
IMPLICIT NONE
!
REAL (WP) :: X,L
REAL (WP) :: IQWE_LB_FF
REAL (WP) :: Q,Y,Y2
REAL (WP) :: NU1,NU2,DE1,DE2
!
REAL (WP) :: DEXP
!
Q=TWO*X*KF_SI !
Y=Q*L !
Y2=Y*Y !
!
NU1=THREE*Y + EIGHT*PI2/Y !
NU2=32.0E0_WP*PI2*PI2*(ONE-DEXP(-Y)) !
DE1=Y2 + FOUR*PI2 !
DE2=Y2 * DE1*DE1 !
!
IQWE_LB_FF=NU1/DE1 - NU2/DE2 !
!
END FUNCTION IQWE_LB_FF
!
!=======================================================================
!
FUNCTION PC1_QWI_FF(X,OM0)
!
! This function computes the form factor of a quantum wire under
! an harmonic confinement potential of the form 1/2 m omega_0^2 y^2
! in the y direction
!
! Reference: (1) M. Tas, PhD thesis, Middle East Technical University (2004)
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * OM0 : frequency of the confinement potential (SI)
!
! Output parameters:
!
! * PC1_QWI_FF : form factor
!
! Author : D. Sébilleau
!
! Last modified : 10 Jun 2020
!
!
USE REAL_NUMBERS, ONLY : TWO,FOURTH
USE CONSTANTS_P1, ONLY : H_BAR,M_E
USE FERMI_SI, ONLY : KF_SI
USE BESSEL
!
IMPLICIT NONE
!
REAL (WP) :: X,OM0
REAL (WP) :: PC1_QWI_FF
REAL (WP) :: Q,B,ZZ,K0
!
REAL (WP) :: DSQRT,DEXP
!
Q=TWO*X*KF_SI !
!
! Characteristic length of the harmonic potential
! (serves as effective diameter of quantum wire)
!
B=DSQRT(H_BAR/(M_E*OM0)) !
!
ZZ=FOURTH*Q*Q*B*B !
!
! Computing the Bessel function
!
K0=BESSK(0,ZZ) !
!
PC1_QWI_FF=DEXP(ZZ)*K0 ! ref. 1 eq. (3.70)
!
END FUNCTION PC1_QWI_FF
!
!=======================================================================
!
FUNCTION PC2_QWI_FF(X,OM0)
!
! This function computes the form factor of a quantum wire under
! an harmonic confinement potential of the form 1/8 m omega_0^2 (x^2 + y^2)
! in the y direction
!
! Reference: (1) M. Tas, PhD thesis, Middle East Technical University (2004)
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * OM0 : frequency of the confinement potential (SI)
!
! Output parameters:
!
! * PC2_QWI_FF : form factor
!
! Author : D. Sébilleau
!
! Last modified : 10 Jun 2020
!
!
USE REAL_NUMBERS, ONLY : TWO
USE CONSTANTS_P1, ONLY : H_BAR,M_E
USE FERMI_SI, ONLY : KF_SI
USE BESSEL, ONLY : EXPINT
!
IMPLICIT NONE
!
REAL (WP) :: X,OM0
REAL (WP) :: PC2_QWI_FF
REAL (WP) :: Q,B,ZZ
!
REAL (WP) :: DSQRT,DEXP
!
Q=TWO*X*KF_SI !
!
! Characteristic length of the harmonic potential
! (serves as effective diameter of quantum wire)
!
B=DSQRT(H_BAR/(M_E*OM0)) !
!
ZZ=Q*Q*B*B !
!
PC2_QWI_FF=DEXP(ZZ) * EXPINT(1,ZZ) ! ref. 1 eq. (3.75)
!
END FUNCTION PC2_QWI_FF
!
!=======================================================================
!
FUNCTION SOF_COR_FF(X,R0)
!
! This function computes the form factor of a quantum wire under
! the soft core potential
!
! Reference: (1) N. Nessi and A. Iucci, Phys. Rev. B 87, 085137 (2013)
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * R0 : quantum wire radius (SI)
!
! Output parameters:
!
! * SOF_COR_FF : form factor
!
! Author : D. Sébilleau
!
! Last modified : 10 Jun 2020
!
!
USE REAL_NUMBERS, ONLY : TWO
USE FERMI_SI, ONLY : KF_SI
USE BESSEL
!
IMPLICIT NONE
!
REAL (WP) :: X,R0
REAL (WP) :: Q,SOF_COR_FF
!
Q=TWO*X*KF_SI !
!
SOF_COR_FF=TWO*BESSK(0,Q*R0) !
!
END FUNCTION SOF_COR_FF
!
!=======================================================================
!
FUNCTION SWC_QWI_FF(X,A)
!
! This function computes the form factor of the a quantum-well wire
! modelled as a square well with an infinite barrier when only
! the lower subband is filled
!
! The barrier is from -a/2 to a/2
!
!
! Reference: (1) M. Tas, PhD thesis, Middle East Technical University (2004)
!
!
! Input parameters:
!
! * X : dimensionless factor --> X = q / (2 * k_F)
! * A : length of the quantum well (SI)
!
! Output parameters:
!
! * SWC_QWI_FF : form factor
!
! Author : D. Sébilleau
!
! Last modified : 10 Jun 2020
!
!
USE DIMENSION_CODE, ONLY : NZ_MAX
USE REAL_NUMBERS, ONLY : ONE,TWO
USE FERMI_SI, ONLY : KF_SI
USE PI_ETC, ONLY : PI,PI_INV
USE INTEGRATION, ONLY : INTEGR_L
USE BESSEL
!
IMPLICIT NONE
! ! points
REAL (WP) :: X,A
REAL (WP) :: SWC_QWI_FF
REAL (WP) :: Q,XX,F(NZ_MAX),H,AA
!
REAL (WP) :: DFLOAT,DCOS,DSIN
!
INTEGER :: K,ID
!
ID=1 !
!
Q=TWO*X*KF_SI !
!
! Constructing the integrand function
!
DO K=1,NZ_MAX !
!
XX=DFLOAT(K-1)/DFLOAT(NZ_MAX-1) !
F(K)=BESSK(0,Q*A*XX) * ( & !
TWO-(ONE-XX)*DCOS(TWO*PI*XX) + & !
1.5E0_WP*PI_INV*DSIN(TWO*PI*XX) & !
) !
!
ENDDO !
!
H=ONE/DFLOAT(NZ_MAX-1) ! step
!
! Performing the integration
!
CALL INTEGR_L(F,H,NZ_MAX,NZ_MAX,AA,ID) !
SWC_QWI_FF=TWO*AA ! ref. 1 eq. (3.87)
!
END FUNCTION SWC_QWI_FF
!
END MODULE CONFINEMENT_FF