2020-07-03 11:00:49 +02:00
|
|
|
subroutine LIM_RKS(x_rung,x_DFA,c_rung,c_DFA,LDA_centered,nEns,aCC_w1,aCC_w2,nGrid,weight, &
|
2020-08-02 13:09:30 +02:00
|
|
|
maxSCF,thresh,max_diis,guess_type,nBas,AO,dAO,nO,nV,S,T,V,Hc,ERI,X,ENuc,&
|
|
|
|
c,occnum,Cx_choice)
|
2020-03-18 11:10:29 +01:00
|
|
|
|
|
|
|
! Perform restricted Kohn-Sham calculation for ensembles
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
include 'parameters.h'
|
|
|
|
|
|
|
|
! Input variables
|
|
|
|
|
|
|
|
integer,intent(in) :: x_rung,c_rung
|
|
|
|
character(len=12),intent(in) :: x_DFA,c_DFA
|
2020-04-01 22:59:52 +02:00
|
|
|
logical,intent(in) :: LDA_centered
|
2020-03-18 11:10:29 +01:00
|
|
|
integer,intent(in) :: nEns
|
2020-07-03 11:00:49 +02:00
|
|
|
double precision,intent(in) :: aCC_w1(3)
|
|
|
|
double precision,intent(in) :: aCC_w2(3)
|
2020-03-18 11:10:29 +01:00
|
|
|
integer,intent(in) :: nGrid
|
|
|
|
double precision,intent(in) :: weight(nGrid)
|
|
|
|
integer,intent(in) :: maxSCF,max_diis,guess_type
|
|
|
|
double precision,intent(in) :: thresh
|
|
|
|
integer,intent(in) :: nBas
|
|
|
|
double precision,intent(in) :: AO(nBas,nGrid)
|
|
|
|
double precision,intent(in) :: dAO(ncart,nBas,nGrid)
|
|
|
|
|
|
|
|
integer,intent(in) :: nO,nV
|
|
|
|
double precision,intent(in) :: S(nBas,nBas)
|
|
|
|
double precision,intent(in) :: T(nBas,nBas)
|
|
|
|
double precision,intent(in) :: V(nBas,nBas)
|
|
|
|
double precision,intent(in) :: Hc(nBas,nBas)
|
|
|
|
double precision,intent(in) :: X(nBas,nBas)
|
|
|
|
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
|
|
|
|
double precision,intent(in) :: ENuc
|
2020-08-02 13:09:30 +02:00
|
|
|
double precision,intent(in) :: occnum(2,2,3)
|
|
|
|
integer,intent(in) :: Cx_choice
|
2020-08-01 11:45:17 +02:00
|
|
|
|
2020-03-18 11:10:29 +01:00
|
|
|
|
2020-03-31 23:33:48 +02:00
|
|
|
double precision,intent(out) :: c(nBas,nBas)
|
2020-03-18 16:06:29 +01:00
|
|
|
|
2020-03-18 11:10:29 +01:00
|
|
|
! Local variables
|
|
|
|
|
|
|
|
integer :: iEns
|
2020-04-26 23:17:27 +02:00
|
|
|
double precision :: Ew(nEnS)
|
2020-03-18 11:10:29 +01:00
|
|
|
double precision :: wLIM(nEns)
|
2020-04-09 21:58:10 +02:00
|
|
|
double precision :: Om(nEns)
|
2020-03-18 11:10:29 +01:00
|
|
|
|
|
|
|
! Hello world
|
|
|
|
|
|
|
|
write(*,*)
|
|
|
|
write(*,*)'************************************************'
|
|
|
|
write(*,*)'* Linear-interpolation method *'
|
|
|
|
write(*,*)'* for excitation energies *'
|
|
|
|
write(*,*)'************************************************'
|
|
|
|
write(*,*)
|
|
|
|
|
2020-04-26 23:17:27 +02:00
|
|
|
! Initializatio
|
|
|
|
|
|
|
|
Ew(:) = 0d0
|
|
|
|
Om(:) = 0d0
|
|
|
|
|
2020-03-18 11:10:29 +01:00
|
|
|
!------------------------------------------------------------------------
|
|
|
|
! Zero-weight calculation
|
|
|
|
!------------------------------------------------------------------------
|
|
|
|
|
|
|
|
write(*,'(A40)') '*************************************************'
|
|
|
|
write(*,'(A40)') ' ZERO-WEIGHT CALCULATION '
|
|
|
|
write(*,'(A40)') '*************************************************'
|
|
|
|
|
2020-04-26 23:17:27 +02:00
|
|
|
wLIM(1) = 1d0
|
|
|
|
wLIM(2) = 0d0
|
|
|
|
wLIM(3) = 0d0
|
2020-03-18 11:10:29 +01:00
|
|
|
|
|
|
|
do iEns=1,nEns
|
|
|
|
write(*,'(A20,I2,A2,F16.10)') ' Weight of state ',iEns,': ',wLIM(iEns)
|
|
|
|
end do
|
|
|
|
write(*,'(A40)') '*************************************************'
|
|
|
|
write(*,*)
|
|
|
|
|
2020-07-03 11:00:49 +02:00
|
|
|
call GOK_RKS(.false.,x_rung,x_DFA,c_rung,c_DFA,LDA_centered,nEns,wLIM,aCC_w1,aCC_w2,nGrid,weight,maxSCF,thresh, &
|
2020-08-02 13:09:30 +02:00
|
|
|
max_diis,guess_type,nBas,AO,dAO,nO,nV,S,T,V,Hc,ERI,X,ENuc,Ew(1),c,occnum,Cx_choice)
|
2020-03-18 11:10:29 +01:00
|
|
|
|
|
|
|
!------------------------------------------------------------------------
|
|
|
|
! Equiensemble calculation
|
|
|
|
!------------------------------------------------------------------------
|
|
|
|
|
|
|
|
write(*,'(A40)') '*************************************************'
|
2020-04-26 23:17:27 +02:00
|
|
|
write(*,'(A40)') ' TWO-STATE EQUI-WEIGHT CALCULATION '
|
2020-03-18 11:10:29 +01:00
|
|
|
write(*,'(A40)') '*************************************************'
|
|
|
|
|
2020-04-26 23:17:27 +02:00
|
|
|
wLIM(1) = 0.5d0
|
2020-05-15 13:08:07 +02:00
|
|
|
wLIM(2) = 0.0d0
|
|
|
|
wLIM(3) = 0.5d0
|
2020-04-26 23:17:27 +02:00
|
|
|
|
2020-03-18 11:10:29 +01:00
|
|
|
do iEns=1,nEns
|
2020-04-26 23:17:27 +02:00
|
|
|
write(*,'(A20,I2,A2,F16.10)') ' Weight of state ',iEns,': ',wLIM(iEns)
|
2020-03-18 11:10:29 +01:00
|
|
|
end do
|
|
|
|
write(*,'(A40)') '*************************************************'
|
|
|
|
write(*,*)
|
|
|
|
|
2020-07-03 11:00:49 +02:00
|
|
|
call GOK_RKS(.true.,x_rung,x_DFA,c_rung,c_DFA,LDA_centered,nEns,wLIM,aCC_w1,aCC_w2,nGrid,weight,maxSCF,thresh, &
|
2020-08-02 13:09:30 +02:00
|
|
|
max_diis,guess_type,nBas,AO,dAO,nO,nV,S,T,V,Hc,ERI,X,ENuc,Ew(2),c,occnum,Cx_choice)
|
2020-04-26 23:17:27 +02:00
|
|
|
|
|
|
|
!------------------------------------------------------------------------
|
|
|
|
! Equiensemble calculation
|
|
|
|
!------------------------------------------------------------------------
|
|
|
|
|
|
|
|
write(*,'(A40)') '*************************************************'
|
|
|
|
write(*,'(A40)') ' THREE-STATE EQUI-WEIGHT CALCULATION '
|
|
|
|
write(*,'(A40)') '*************************************************'
|
|
|
|
|
|
|
|
wLIM(1) = 1d0/3d0
|
|
|
|
wLIM(2) = 1d0/3d0
|
|
|
|
wLIM(3) = 1d0/3d0
|
|
|
|
|
|
|
|
do iEns=1,nEns
|
|
|
|
write(*,'(A20,I2,A2,F16.10)') ' Weight of state ',iEns,': ',wLIM(iEns)
|
|
|
|
end do
|
|
|
|
write(*,'(A40)') '*************************************************'
|
|
|
|
write(*,*)
|
|
|
|
|
2020-07-03 11:00:49 +02:00
|
|
|
! call GOK_RKS(.true.,x_rung,x_DFA,c_rung,c_DFA,LDA_centered,nEns,wLIM,aCC_w1,aCC_w2,nGrid,weight,maxSCF,thresh, &
|
2020-08-02 13:09:30 +02:00
|
|
|
! max_diis,guess_type,nBas,AO,dAO,nO,nV,S,T,V,Hc,ERI,X,ENuc,Ew(3),c,occnum,Cx_choice)
|
2020-04-26 23:17:27 +02:00
|
|
|
|
2020-03-18 11:10:29 +01:00
|
|
|
|
|
|
|
!------------------------------------------------------------------------
|
|
|
|
! LIM excitation energies
|
|
|
|
!------------------------------------------------------------------------
|
|
|
|
|
2020-04-27 23:04:00 +02:00
|
|
|
Om(2) = 2d0*(Ew(2) - Ew(1))
|
2020-05-15 13:08:07 +02:00
|
|
|
! Om(3) = 3d0*(Ew(3) - Ew(2)) + 0.5d0*Om(2)
|
2020-03-18 11:10:29 +01:00
|
|
|
|
|
|
|
write(*,'(A60)') '-------------------------------------------------'
|
|
|
|
write(*,'(A60)') ' LINEAR INTERPOLATION METHOD EXCITATION ENERGIES '
|
|
|
|
write(*,'(A60)') '-------------------------------------------------'
|
|
|
|
|
2020-04-26 23:17:27 +02:00
|
|
|
write(*,'(A44,F16.10,A3)') ' Ensemble energy #1 ',Ew(1),' au'
|
|
|
|
write(*,'(A44,F16.10,A3)') ' Ensemble energy #2 ',Ew(2),' au'
|
|
|
|
write(*,'(A44,F16.10,A3)') ' Ensemble energy #3 ',Ew(3),' au'
|
2020-03-18 11:10:29 +01:00
|
|
|
write(*,'(A60)') '-------------------------------------------------'
|
|
|
|
do iEns=2,nEns
|
|
|
|
write(*,'(A40,I2,A2,F16.10,A3)') ' Excitation energy 1 ->',iEns,': ',Om(iEns), ' au'
|
2020-04-26 23:17:27 +02:00
|
|
|
end do
|
|
|
|
write(*,*)
|
|
|
|
do iEns=2,nEns
|
2020-03-18 11:10:29 +01:00
|
|
|
write(*,'(A40,I2,A2,F16.10,A3)') ' Excitation energy 1 ->',iEns,': ',Om(iEns)*HaToeV,' eV'
|
|
|
|
end do
|
|
|
|
write(*,'(A60)') '-------------------------------------------------'
|
|
|
|
|
|
|
|
|
|
|
|
end subroutine LIM_RKS
|