mirror of
https://github.com/pfloos/quack
synced 2024-11-04 13:13:51 +01:00
UGW under progress
This commit is contained in:
parent
11949b84bc
commit
1d3e156838
@ -1,4 +1,4 @@
|
||||
# nAt nEla nElb nCore nRyd
|
||||
1 1 1 0 0
|
||||
1 2 1 0 0
|
||||
# Znuc x y z
|
||||
Li 0.0 0.0 0.0
|
||||
|
@ -112,7 +112,7 @@ subroutine G0W0(doACFDT,exchange_kernel,doXBS,COHSEX,SOSEX,BSE,TDA_W,TDA,
|
||||
|
||||
if(SOSEX) call excitation_density_SOSEX(nBas,nC,nO,nR,nS,ERI,XpY(:,:,ispin),rhox(:,:,:,ispin))
|
||||
|
||||
call self_energy_correlation_diag(.false.,COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eHF, &
|
||||
call self_energy_correlation_diag(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eHF, &
|
||||
Omega(:,ispin),rho(:,:,:,ispin),rhox(:,:,:,ispin),EcGM,SigC)
|
||||
|
||||
! Compute renormalization factor
|
||||
|
256
src/QuAcK/UG0W0.f90
Normal file
256
src/QuAcK/UG0W0.f90
Normal file
@ -0,0 +1,256 @@
|
||||
subroutine UG0W0(doACFDT,exchange_kernel,doXBS,COHSEX,BSE,TDA_W,TDA,dBSE,dTDA,evDyn, &
|
||||
singlet_manifold,triplet_manifold,linearize,eta,nBas,nC,nO,nV,nR,nS, &
|
||||
ENuc,EUHF,Hc,ERI_aa,ERI_ab,ERI_bb,PHF,cHF,eHF,eGW)
|
||||
|
||||
! Perform unrestricted G0W0 calculation
|
||||
|
||||
implicit none
|
||||
include 'parameters.h'
|
||||
include 'quadrature.h'
|
||||
|
||||
! Input variables
|
||||
|
||||
logical,intent(in) :: doACFDT
|
||||
logical,intent(in) :: exchange_kernel
|
||||
logical,intent(in) :: doXBS
|
||||
logical,intent(in) :: COHSEX
|
||||
logical,intent(in) :: BSE
|
||||
logical,intent(in) :: TDA_W
|
||||
logical,intent(in) :: TDA
|
||||
logical,intent(in) :: dBSE
|
||||
logical,intent(in) :: dTDA
|
||||
logical,intent(in) :: evDyn
|
||||
logical,intent(in) :: singlet_manifold
|
||||
logical,intent(in) :: triplet_manifold
|
||||
logical,intent(in) :: linearize
|
||||
double precision,intent(in) :: eta
|
||||
|
||||
integer,intent(in) :: nBas
|
||||
integer,intent(in) :: nC(nspin)
|
||||
integer,intent(in) :: nO(nspin)
|
||||
integer,intent(in) :: nV(nspin)
|
||||
integer,intent(in) :: nR(nspin)
|
||||
integer,intent(in) :: nS(nspin)
|
||||
double precision,intent(in) :: ENuc
|
||||
double precision,intent(in) :: EUHF
|
||||
double precision,intent(in) :: eHF(nBas,nspin)
|
||||
double precision,intent(in) :: cHF(nBas,nBas,nspin)
|
||||
double precision,intent(in) :: PHF(nBas,nBas,nspin)
|
||||
double precision,intent(in) :: Hc(nBas,nBas,nspin)
|
||||
double precision,intent(in) :: ERI_aa(nBas,nBas,nBas,nBas)
|
||||
double precision,intent(in) :: ERI_ab(nBas,nBas,nBas,nBas)
|
||||
double precision,intent(in) :: ERI_bb(nBas,nBas,nBas,nBas)
|
||||
|
||||
! Local variables
|
||||
|
||||
logical :: print_W = .true.
|
||||
integer :: ispin
|
||||
integer :: iblock
|
||||
integer :: bra
|
||||
integer :: ket
|
||||
integer :: nSa
|
||||
integer :: nSb
|
||||
integer :: nSt
|
||||
double precision :: EcRPA(nspin)
|
||||
double precision :: EcBSE(nspin)
|
||||
double precision :: EcAC(nspin)
|
||||
double precision,allocatable :: SigC(:,:)
|
||||
double precision,allocatable :: Z(:,:)
|
||||
double precision,allocatable :: Omega(:)
|
||||
double precision,allocatable :: XpY_a(:,:)
|
||||
double precision,allocatable :: XpY_b(:,:)
|
||||
double precision,allocatable :: XmY_a(:,:)
|
||||
double precision,allocatable :: XmY_b(:,:)
|
||||
double precision,allocatable :: rho(:,:,:,:)
|
||||
|
||||
double precision,allocatable :: eGWlin(:,:)
|
||||
|
||||
! Output variables
|
||||
|
||||
double precision :: eGW(nBas,nspin)
|
||||
|
||||
! Hello world
|
||||
|
||||
write(*,*)
|
||||
write(*,*)'************************************************'
|
||||
write(*,*)'| One-shot G0W0 calculation |'
|
||||
write(*,*)'| *** Unrestricted version *** |'
|
||||
write(*,*)'************************************************'
|
||||
write(*,*)
|
||||
|
||||
! Initialization
|
||||
|
||||
EcRPA(:) = 0d0
|
||||
|
||||
! COHSEX approximation
|
||||
|
||||
if(COHSEX) write(*,*) 'COHSEX approximation activated!'
|
||||
write(*,*)
|
||||
|
||||
! TDA for W
|
||||
|
||||
if(TDA_W) write(*,*) 'Tamm-Dancoff approximation for dynamic screening!'
|
||||
write(*,*)
|
||||
|
||||
! TDA
|
||||
|
||||
if(TDA) write(*,*) 'Tamm-Dancoff approximation activated!'
|
||||
write(*,*)
|
||||
|
||||
! Memory allocation
|
||||
|
||||
nSa = nS(1)
|
||||
nSb = nS(2)
|
||||
nSt = nSa + nSb
|
||||
|
||||
allocate(SigC(nBas,nspin),Z(nBas,nspin),Omega(nSt),XpY_a(nSa,nSa),XpY_b(nSb,nSb),XmY_a(nSa,nSa),XmY_b(nSb,nSb), &
|
||||
rho(nBas,nBas,nSt,nspin),eGWlin(nBas,nspin))
|
||||
|
||||
! Compute linear response
|
||||
|
||||
!----------------------------------------------
|
||||
! alpha-alpha block
|
||||
!----------------------------------------------
|
||||
|
||||
ispin = 1
|
||||
iblock = 3
|
||||
|
||||
call linear_response(iblock,.true.,TDA_W,.false.,eta,nBas,nC(ispin),nO(ispin),nV(ispin),nR(ispin),nSa,1d0, &
|
||||
eHF(:,ispin),ERI_aa,rho(:,:,1:nSa,ispin),EcRPA(ispin),Omega(1:nSa),XpY_a,XmY_a)
|
||||
|
||||
if(print_W) call print_excitation('RPA@HF alpha',iblock,nSa,Omega(1:nSa))
|
||||
|
||||
!----------------------------------------------
|
||||
! alpha-beta block
|
||||
!----------------------------------------------
|
||||
|
||||
ispin = 2
|
||||
iblock = 3
|
||||
|
||||
call linear_response(iblock,.true.,TDA_W,.false.,eta,nBas,nC(ispin),nO(ispin),nV(ispin),nR(ispin),nSb,1d0, &
|
||||
eHF(:,ispin),ERI_bb,rho(:,:,nSa+1:nSt,ispin),EcRPA(ispin),Omega(nSa+1:nSt),XpY_b,XmY_b)
|
||||
|
||||
if(print_W) call print_excitation('RPA@HF beta ',iblock,nSb,Omega(nSa+1:nSt))
|
||||
|
||||
!----------------------------------------------
|
||||
! Excitation densities for alpha self-energy
|
||||
!----------------------------------------------
|
||||
|
||||
call unrestricted_excitation_density(nBas,nC,nO,nR,nSa,nSb,nSt,ERI_aa,ERI_ab,ERI_bb,XpY_a,XpY_b,rho)
|
||||
|
||||
!----------------------
|
||||
! Compute self-energy
|
||||
!----------------------
|
||||
|
||||
call unrestricted_self_energy_correlation_diag(eta,nBas,nC,nO,nV,nR,nSa,nSb,nSt,eHF,Omega,rho,SigC)
|
||||
|
||||
! Compute renormalization factor
|
||||
|
||||
! call renormalization_factor(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eHF, &
|
||||
! Omega(:,ispin),rho(:,:,:,ispin),rhox(:,:,:,ispin),Z(:))
|
||||
|
||||
! Solve the quasi-particle equation
|
||||
Z(:,:) = 1d0
|
||||
eGWlin(:,:) = eHF(:,:) + Z(:,:)*SigC(:,:)
|
||||
|
||||
if(linearize) then
|
||||
|
||||
write(*,*) ' *** Quasiparticle energies obtained by linearization *** '
|
||||
write(*,*)
|
||||
|
||||
eGW(:,:) = eGWlin(:,:)
|
||||
|
||||
else
|
||||
|
||||
! Find graphical solution of the QP equation
|
||||
|
||||
! do is=1,nspin
|
||||
! call QP_graph(nBas,nC(:,is),nO(:,is),nV(:,is),nR(:,is),nS(:,is),eta,eHF(:,is),Omega(:,is), &
|
||||
! rho(:,:,:,ispin),eGWlin(:,is),eGW(:,is))
|
||||
! end do
|
||||
|
||||
end if
|
||||
|
||||
! Dump results
|
||||
|
||||
do ispin=1,nspin
|
||||
call print_G0W0(nBas,nO(ispin),eHF(:,ispin),ENuc,EUHF,SigC(:,ispin),Z(:,ispin),eGW(:,ispin),EcRPA(ispin))
|
||||
end do
|
||||
|
||||
! Compute the RPA correlation energy
|
||||
|
||||
! call linear_response(ispin,.true.,TDA_W,.false.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGW,ERI, &
|
||||
! rho(:,:,:,ispin),EcRPA(ispin),Omega(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
|
||||
|
||||
write(*,*)
|
||||
write(*,*)'-------------------------------------------------------------------------------'
|
||||
write(*,'(2X,A50,F20.10)') 'Tr@RPA@G0W0 correlation energy (singlet) =',EcRPA(1)
|
||||
write(*,'(2X,A50,F20.10)') 'Tr@RPA@G0W0 correlation energy (triplet) =',EcRPA(2)
|
||||
write(*,'(2X,A50,F20.10)') 'Tr@RPA@G0W0 correlation energy =',EcRPA(1) + EcRPA(2)
|
||||
write(*,'(2X,A50,F20.10)') 'Tr@RPA@G0W0 total energy =',ENuc + EUHF + EcRPA(1) + EcRPA(2)
|
||||
write(*,*)'-------------------------------------------------------------------------------'
|
||||
write(*,*)
|
||||
|
||||
! Perform BSE calculation
|
||||
|
||||
! if(BSE) then
|
||||
|
||||
! call Bethe_Salpeter(TDA_W,TDA,dBSE,dTDA,evDyn,singlet_manifold,triplet_manifold,eta, &
|
||||
! nBas,nC,nO,nV,nR,nS,ERI,eHF,eGW,Omega,XpY,XmY,rho,EcRPA,EcBSE)
|
||||
|
||||
! if(exchange_kernel) then
|
||||
!
|
||||
! EcRPA(1) = 0.5d0*EcRPA(1)
|
||||
! EcRPA(2) = 1.5d0*EcRPA(1)
|
||||
!
|
||||
! end if
|
||||
|
||||
! write(*,*)
|
||||
! write(*,*)'-------------------------------------------------------------------------------'
|
||||
! write(*,'(2X,A50,F20.10)') 'Tr@BSE@G0W0 correlation energy (singlet) =',EcBSE(1)
|
||||
! write(*,'(2X,A50,F20.10)') 'Tr@BSE@G0W0 correlation energy (triplet) =',EcBSE(2)
|
||||
! write(*,'(2X,A50,F20.10)') 'Tr@BSE@G0W0 correlation energy =',EcBSE(1) + EcBSE(2)
|
||||
! write(*,'(2X,A50,F20.10)') 'Tr@BSE@G0W0 total energy =',ENuc + EUHF + EcBSE(1) + EcBSE(2)
|
||||
! write(*,*)'-------------------------------------------------------------------------------'
|
||||
! write(*,*)
|
||||
|
||||
! Compute the BSE correlation energy via the adiabatic connection
|
||||
|
||||
! if(doACFDT) then
|
||||
|
||||
! write(*,*) '------------------------------------------------------'
|
||||
! write(*,*) 'Adiabatic connection version of BSE correlation energy'
|
||||
! write(*,*) '------------------------------------------------------'
|
||||
! write(*,*)
|
||||
|
||||
! if(doXBS) then
|
||||
|
||||
! write(*,*) '*** scaled screening version (XBS) ***'
|
||||
! write(*,*)
|
||||
|
||||
! end if
|
||||
|
||||
! call ACFDT(exchange_kernel,doXBS,.true.,TDA_W,TDA,BSE,singlet_manifold,triplet_manifold,eta, &
|
||||
! nBas,nC,nO,nV,nR,nS,ERI,eHF,eGW,Omega,XpY,XmY,rho,EcAC)
|
||||
|
||||
! if(exchange_kernel) then
|
||||
!
|
||||
! EcAC(1) = 0.5d0*EcAC(1)
|
||||
! EcAC(2) = 1.5d0*EcAC(1)
|
||||
!
|
||||
! end if
|
||||
|
||||
! write(*,*)
|
||||
! write(*,*)'-------------------------------------------------------------------------------'
|
||||
! write(*,'(2X,A50,F20.10)') 'AC@BSE@G0W0 correlation energy (singlet) =',EcAC(1)
|
||||
! write(*,'(2X,A50,F20.10)') 'AC@BSE@G0W0 correlation energy (triplet) =',EcAC(2)
|
||||
! write(*,'(2X,A50,F20.10)') 'AC@BSE@G0W0 correlation energy =',EcAC(1) + EcAC(2)
|
||||
! write(*,'(2X,A50,F20.10)') 'AC@BSE@G0W0 total energy =',ENuc + EUHF + EcAC(1) + EcAC(2)
|
||||
! write(*,*)'-------------------------------------------------------------------------------'
|
||||
! write(*,*)
|
||||
|
||||
! end if
|
||||
|
||||
! end if
|
||||
|
||||
end subroutine UG0W0
|
@ -140,14 +140,14 @@ subroutine evGW(maxSCF,thresh,max_diis,doACFDT,exchange_kernel,doXBS,COHSEX,SOSE
|
||||
|
||||
if(G0W) then
|
||||
|
||||
call self_energy_correlation_diag(.false.,COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eHF, &
|
||||
call self_energy_correlation_diag(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eHF, &
|
||||
Omega(:,ispin),rho(:,:,:,ispin),rhox(:,:,:,ispin),EcGM,SigC)
|
||||
call renormalization_factor(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eHF, &
|
||||
Omega(:,ispin),rho(:,:,:,ispin),rhox(:,:,:,ispin),Z(:))
|
||||
|
||||
else
|
||||
|
||||
call self_energy_correlation_diag(.false.,COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eGW, &
|
||||
call self_energy_correlation_diag(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eGW, &
|
||||
Omega(:,ispin),rho(:,:,:,ispin),rhox(:,:,:,ispin),EcGM,SigC)
|
||||
call renormalization_factor(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,eGW, &
|
||||
Omega(:,ispin),rho(:,:,:,ispin),rhox(:,:,:,ispin),Z(:))
|
||||
|
46
src/QuAcK/print_UG0W0.f90
Normal file
46
src/QuAcK/print_UG0W0.f90
Normal file
@ -0,0 +1,46 @@
|
||||
subroutine print_UG0W0(nBas,nO,e,ENuc,EHF,SigmaC,Z,eGW,EcRPA)
|
||||
|
||||
! Print one-electron energies and other stuff for G0W0
|
||||
|
||||
implicit none
|
||||
include 'parameters.h'
|
||||
|
||||
integer,intent(in) :: nBas,nO
|
||||
double precision,intent(in) :: ENuc
|
||||
double precision,intent(in) :: EHF
|
||||
double precision,intent(in) :: EcRPA
|
||||
double precision,intent(in) :: e(nBas),SigmaC(nBas),Z(nBas),eGW(nBas)
|
||||
|
||||
integer :: x,HOMO,LUMO
|
||||
double precision :: Gap
|
||||
|
||||
! HOMO and LUMO
|
||||
|
||||
HOMO = nO
|
||||
LUMO = HOMO + 1
|
||||
Gap = eGW(LUMO)-eGW(HOMO)
|
||||
|
||||
! Dump results
|
||||
|
||||
write(*,*)'-------------------------------------------------------------------------------'
|
||||
write(*,*)' One-shot G0W0 calculation'
|
||||
write(*,*)'-------------------------------------------------------------------------------'
|
||||
write(*,'(1X,A1,1X,A3,1X,A1,1X,A15,1X,A1,1X,A15,1X,A1,1X,A15,1X,A1,1X,A15,1X,A1,1X)') &
|
||||
'|','#','|','e_HF (eV)','|','Sigma_c (eV)','|','Z','|','e_QP (eV)','|'
|
||||
write(*,*)'-------------------------------------------------------------------------------'
|
||||
|
||||
do x=1,nBas
|
||||
write(*,'(1X,A1,1X,I3,1X,A1,1X,F15.6,1X,A1,1X,F15.6,1X,A1,1X,F15.6,1X,A1,1X,F15.6,1X,A1,1X)') &
|
||||
'|',x,'|',e(x)*HaToeV,'|',SigmaC(x)*HaToeV,'|',Z(x),'|',eGW(x)*HaToeV,'|'
|
||||
enddo
|
||||
|
||||
write(*,*)'-------------------------------------------------------------------------------'
|
||||
write(*,'(2X,A30,F15.6)') 'G0W0 HOMO energy (eV):',eGW(HOMO)*HaToeV
|
||||
write(*,'(2X,A30,F15.6)') 'G0W0 LUMO energy (eV):',eGW(LUMO)*HaToeV
|
||||
write(*,'(2X,A30,F15.6)') 'G0W0 HOMO-LUMO gap (eV):',Gap*HaToeV
|
||||
write(*,*)'-------------------------------------------------------------------------------'
|
||||
write(*,*)
|
||||
|
||||
end subroutine print_UG0W0
|
||||
|
||||
|
@ -39,12 +39,8 @@ subroutine self_energy_correlation(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,e,Omega,
|
||||
do x=nC+1,nBas-nR
|
||||
do y=nC+1,nBas-nR
|
||||
do i=nC+1,nO
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
SigC(x,y) = SigC(x,y) + 4d0*rho(x,i,jb)*rho(y,i,jb)/Omega(jb)
|
||||
enddo
|
||||
do jb=1,nS
|
||||
SigC(x,y) = SigC(x,y) + 4d0*rho(x,i,jb)*rho(y,i,jb)/Omega(jb)
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
@ -55,12 +51,8 @@ subroutine self_energy_correlation(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,e,Omega,
|
||||
do x=nC+1,nBas-nR
|
||||
do y=nC+1,nBas-nR
|
||||
do p=nC+1,nBas-nR
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
SigC(x,y) = SigC(x,y) - 2d0*rho(x,p,jb)*rho(y,p,jb)/Omega(jb)
|
||||
enddo
|
||||
do jb=1,nS
|
||||
SigC(x,y) = SigC(x,y) - 2d0*rho(x,p,jb)*rho(y,p,jb)/Omega(jb)
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
@ -78,13 +70,9 @@ subroutine self_energy_correlation(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,e,Omega,
|
||||
do x=nC+1,nBas-nR
|
||||
do y=nC+1,nBas-nR
|
||||
do i=nC+1,nO
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(x) - e(i) + Omega(jb)
|
||||
SigC(x,y) = SigC(x,y) + 2d0*rho(x,i,jb)*rho(y,i,jb)*eps/(eps**2 + eta**2)
|
||||
enddo
|
||||
do jb=1,nS
|
||||
eps = e(x) - e(i) + Omega(jb)
|
||||
SigC(x,y) = SigC(x,y) + 2d0*rho(x,i,jb)*rho(y,i,jb)*eps/(eps**2 + eta**2)
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
@ -95,13 +83,9 @@ subroutine self_energy_correlation(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,e,Omega,
|
||||
do x=nC+1,nBas-nR
|
||||
do y=nC+1,nBas-nR
|
||||
do a=nO+1,nBas-nR
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(x) - e(a) - Omega(jb)
|
||||
SigC(x,y) = SigC(x,y) + 2d0*rho(x,a,jb)*rho(y,a,jb)*eps/(eps**2 + eta**2)
|
||||
enddo
|
||||
do jb=1,nS
|
||||
eps = e(x) - e(a) - Omega(jb)
|
||||
SigC(x,y) = SigC(x,y) + 2d0*rho(x,a,jb)*rho(y,a,jb)*eps/(eps**2 + eta**2)
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
@ -114,13 +98,9 @@ subroutine self_energy_correlation(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,e,Omega,
|
||||
do x=nC+1,nBas-nR
|
||||
do y=nC+1,nBas-nR
|
||||
do i=nC+1,nO
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(x) - e(i) + Omega(jb)
|
||||
SigC(x,y) = SigC(x,y) - rho(x,i,jb)*rhox(y,i,jb)*eps/(eps**2 + eta**2)
|
||||
enddo
|
||||
do jb=1,nS
|
||||
eps = e(x) - e(i) + Omega(jb)
|
||||
SigC(x,y) = SigC(x,y) - rho(x,i,jb)*rhox(y,i,jb)*eps/(eps**2 + eta**2)
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
@ -131,13 +111,9 @@ subroutine self_energy_correlation(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,e,Omega,
|
||||
do x=nC+1,nBas-nR
|
||||
do y=nC+1,nBas-nR
|
||||
do a=nO+1,nBas-nR
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(x) - e(a) - Omega(jb)
|
||||
SigC(x,y) = SigC(x,y) - rho(x,a,jb)*rhox(y,a,jb)*eps/(eps**2 + eta**2)
|
||||
enddo
|
||||
do jb=1,nS
|
||||
eps = e(x) - e(a) - Omega(jb)
|
||||
SigC(x,y) = SigC(x,y) - rho(x,a,jb)*rhox(y,a,jb)*eps/(eps**2 + eta**2)
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
|
@ -1,4 +1,4 @@
|
||||
subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,e,Omega,rho,rhox,EcGM,SigC)
|
||||
subroutine self_energy_correlation_diag(COHSEX,SOSEX,eta,nBas,nC,nO,nV,nR,nS,e,Omega,rho,rhox,EcGM,SigC)
|
||||
|
||||
! Compute diagonal of the correlation part of the self-energy
|
||||
|
||||
@ -7,7 +7,6 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
! Input variables
|
||||
|
||||
logical,intent(in) :: unrestricted
|
||||
logical,intent(in) :: COHSEX
|
||||
logical,intent(in) :: SOSEX
|
||||
double precision,intent(in) :: eta
|
||||
@ -24,7 +23,7 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
! Local variables
|
||||
|
||||
integer :: i,j,a,b,p,q,jb
|
||||
integer :: i,a,p,q,jb
|
||||
double precision :: eps
|
||||
double precision,external :: SigC_dcgw
|
||||
|
||||
@ -47,12 +46,8 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
do p=nC+1,nBas-nR
|
||||
do i=nC+1,nO
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
SigC(p) = SigC(p) + 4d0*rho(p,i,jb)**2/Omega(jb)
|
||||
end do
|
||||
do jb=1,nS
|
||||
SigC(p) = SigC(p) + 4d0*rho(p,i,jb)**2/Omega(jb)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
@ -61,12 +56,8 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
do p=nC+1,nBas-nR
|
||||
do q=nC+1,nBas-nR
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
SigC(p) = SigC(p) - 2d0*rho(p,q,jb)**2/Omega(jb)
|
||||
end do
|
||||
do jb=1,nS
|
||||
SigC(p) = SigC(p) - 2d0*rho(p,q,jb)**2/Omega(jb)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
@ -88,13 +79,9 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
do p=nC+1,nBas-nR
|
||||
do i=nC+1,nO
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(p) - e(i) + Omega(jb)
|
||||
SigC(p) = SigC(p) - rho(p,i,jb)*rhox(p,i,jb)*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
do jb=1,nS
|
||||
eps = e(p) - e(i) + Omega(jb)
|
||||
SigC(p) = SigC(p) - rho(p,i,jb)*rhox(p,i,jb)*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
@ -103,13 +90,9 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
do p=nC+1,nBas-nR
|
||||
do a=nO+1,nBas-nR
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(p) - e(a) - Omega(jb)
|
||||
SigC(p) = SigC(p) - rho(p,a,jb)*rhox(p,a,jb)*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
do jb=1,nS
|
||||
eps = e(p) - e(a) - Omega(jb)
|
||||
SigC(p) = SigC(p) - rho(p,a,jb)*rhox(p,a,jb)*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
@ -118,13 +101,9 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
do i=nC+1,nO
|
||||
do a=nO+1,nBas-nR
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(a) - e(i) + Omega(jb)
|
||||
EcGM = EcGM + rho(a,i,jb)*rhox(a,i,jb)*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
do jb=1,nS
|
||||
eps = e(a) - e(i) + Omega(jb)
|
||||
EcGM = EcGM + rho(a,i,jb)*rhox(a,i,jb)*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
@ -139,13 +118,9 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
do p=nC+1,nBas-nR
|
||||
do i=nC+1,nO
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(p) - e(i) + Omega(jb)
|
||||
SigC(p) = SigC(p) + 2d0*rho(p,i,jb)**2*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
do jb=1,nS
|
||||
eps = e(p) - e(i) + Omega(jb)
|
||||
SigC(p) = SigC(p) + 2d0*rho(p,i,jb)**2*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
@ -154,13 +129,9 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
|
||||
do p=nC+1,nBas-nR
|
||||
do a=nO+1,nBas-nR
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(p) - e(a) - Omega(jb)
|
||||
SigC(p) = SigC(p) + 2d0*rho(p,a,jb)**2*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
do jb=1,nS
|
||||
eps = e(p) - e(a) - Omega(jb)
|
||||
SigC(p) = SigC(p) + 2d0*rho(p,a,jb)**2*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
@ -170,26 +141,13 @@ subroutine self_energy_correlation_diag(unrestricted,COHSEX,SOSEX,eta,nBas,nC,nO
|
||||
EcGM = 0d0
|
||||
do i=nC+1,nO
|
||||
do a=nO+1,nBas-nR
|
||||
jb = 0
|
||||
do j=nC+1,nO
|
||||
do b=nO+1,nBas-nR
|
||||
jb = jb + 1
|
||||
eps = e(a) - e(i) + Omega(jb)
|
||||
EcGM = EcGM - 2d0*rho(a,i,jb)*rho(a,i,jb)*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
do jb=1,nS
|
||||
eps = e(a) - e(i) + Omega(jb)
|
||||
EcGM = EcGM - 2d0*rho(a,i,jb)*rho(a,i,jb)*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
|
||||
end if
|
||||
|
||||
! Unrestricted reference
|
||||
|
||||
if(unrestricted) then
|
||||
|
||||
SigC(:) = 0.5d0*SigC(:)
|
||||
EcGM = 0.5d0*EcGM
|
||||
|
||||
end if
|
||||
|
||||
end subroutine self_energy_correlation_diag
|
||||
|
107
src/QuAcK/unrestricted_excitation_density.f90
Normal file
107
src/QuAcK/unrestricted_excitation_density.f90
Normal file
@ -0,0 +1,107 @@
|
||||
subroutine unrestricted_excitation_density(nBas,nC,nO,nR,nSa,nSb,nSt,ERI_aa,ERI_ab,ERI_bb,XpY_a,XpY_b,rho)
|
||||
|
||||
! Compute excitation densities for unrestricted reference
|
||||
|
||||
implicit none
|
||||
include 'parameters.h'
|
||||
|
||||
! Input variables
|
||||
|
||||
integer,intent(in) :: nBas
|
||||
integer,intent(in) :: nC(nspin)
|
||||
integer,intent(in) :: nO(nspin)
|
||||
integer,intent(in) :: nR(nspin)
|
||||
integer,intent(in) :: nSa
|
||||
integer,intent(in) :: nSb
|
||||
integer,intent(in) :: nSt
|
||||
double precision,intent(in) :: ERI_aa(nBas,nBas,nBas,nBas)
|
||||
double precision,intent(in) :: ERI_ab(nBas,nBas,nBas,nBas)
|
||||
double precision,intent(in) :: ERI_bb(nBas,nBas,nBas,nBas)
|
||||
double precision,intent(in) :: XpY_a(nSa,nSa)
|
||||
double precision,intent(in) :: XpY_b(nSb,nSb)
|
||||
|
||||
! Local variables
|
||||
|
||||
integer :: ia,jb,p,q,j,b
|
||||
|
||||
! Output variables
|
||||
|
||||
double precision,intent(out) :: rho(nBas,nBas,nSt,nspin)
|
||||
|
||||
! Initialization
|
||||
|
||||
rho(:,:,:,:) = 0d0
|
||||
|
||||
!-------------
|
||||
! alpha block
|
||||
!-------------
|
||||
|
||||
do p=nC(1)+1,nBas-nR(1)
|
||||
do q=nC(1)+1,nBas-nR(1)
|
||||
|
||||
! Same-spin contribution
|
||||
do ia=1,nSa
|
||||
jb = 0
|
||||
do j=nC(1)+1,nO(1)
|
||||
do b=nO(1)+1,nBas-nR(1)
|
||||
jb = jb + 1
|
||||
|
||||
rho(p,q,ia,1) = rho(p,q,ia,1) + ERI_aa(p,j,q,b)*XpY_a(ia,jb)
|
||||
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
|
||||
! Opposite-spin contribution
|
||||
do ia=1,nSb
|
||||
jb = 0
|
||||
do j=nC(2)+1,nO(2)
|
||||
do b=nO(2)+1,nBas-nR(2)
|
||||
jb = jb + 1
|
||||
|
||||
rho(p,q,nSa+ia,1) = rho(p,q,nSa+ia,1) + ERI_ab(p,j,q,b)*XpY_b(ia,jb)
|
||||
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
|
||||
enddo
|
||||
enddo
|
||||
|
||||
!------------
|
||||
! Beta block
|
||||
!------------
|
||||
|
||||
do p=nC(2)+1,nBas-nR(2)
|
||||
do q=nC(2)+1,nBas-nR(2)
|
||||
|
||||
! Same-spin contribution
|
||||
do ia=1,nSb
|
||||
jb = 0
|
||||
do j=nC(2)+1,nO(2)
|
||||
do b=nO(2)+1,nBas-nR(2)
|
||||
jb = jb + 1
|
||||
|
||||
rho(p,q,ia,2) = rho(p,q,ia,2) + ERI_bb(p,j,q,b)*XpY_b(ia,jb)
|
||||
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
|
||||
! Opposite-spin contribution
|
||||
do ia=1,nSa
|
||||
jb = 0
|
||||
do j=nC(1)+1,nO(1)
|
||||
do b=nO(1)+1,nBas-nR(1)
|
||||
jb = jb + 1
|
||||
|
||||
rho(p,q,nSb+ia,2) = rho(p,q,nSb+ia,2) + ERI_ab(j,p,b,q)*XpY_a(ia,jb)
|
||||
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
|
||||
enddo
|
||||
enddo
|
||||
|
||||
end subroutine unrestricted_excitation_density
|
88
src/QuAcK/unrestricted_self_energy_correlation_diag.f90
Normal file
88
src/QuAcK/unrestricted_self_energy_correlation_diag.f90
Normal file
@ -0,0 +1,88 @@
|
||||
subroutine unrestricted_self_energy_correlation_diag(eta,nBas,nC,nO,nV,nR,nSa,nSb,nSt,e,Omega,rho,SigC)
|
||||
|
||||
! Compute diagonal of the correlation part of the self-energy
|
||||
|
||||
implicit none
|
||||
include 'parameters.h'
|
||||
|
||||
! Input variables
|
||||
|
||||
double precision,intent(in) :: eta
|
||||
integer,intent(in) :: nBas
|
||||
integer,intent(in) :: nC(nspin)
|
||||
integer,intent(in) :: nO(nspin)
|
||||
integer,intent(in) :: nV(nspin)
|
||||
integer,intent(in) :: nR(nspin)
|
||||
integer,intent(in) :: nSa
|
||||
integer,intent(in) :: nSb
|
||||
integer,intent(in) :: nSt
|
||||
double precision,intent(in) :: e(nBas,nspin)
|
||||
double precision,intent(in) :: Omega(nSt)
|
||||
double precision,intent(in) :: rho(nBas,nBas,nSt,nspin)
|
||||
|
||||
! Local variables
|
||||
|
||||
integer :: i,j,a,b,p,q,jb
|
||||
double precision :: eps
|
||||
|
||||
! Output variables
|
||||
|
||||
double precision,intent(out) :: SigC(nBas,nspin)
|
||||
|
||||
! Initialize
|
||||
|
||||
SigC(:,:) = 0d0
|
||||
|
||||
!--------------
|
||||
! Spin-up part
|
||||
!--------------
|
||||
|
||||
! Occupied part of the correlation self-energy
|
||||
|
||||
do p=nC(1)+1,nBas-nR(1)
|
||||
do i=nC(1)+1,nO(1)
|
||||
do jb=1,nSt
|
||||
eps = e(p,1) - e(i,1) + Omega(jb)
|
||||
SigC(p,1) = SigC(p,1) + rho(p,i,jb,1)**2*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
|
||||
! Virtual part of the correlation self-energy
|
||||
|
||||
do p=nC(1)+1,nBas-nR(1)
|
||||
do a=nO(1)+1,nBas-nR(1)
|
||||
do jb=1,nSt
|
||||
eps = e(p,1) - e(a,1) - Omega(jb)
|
||||
SigC(p,1) = SigC(p,1) + rho(p,a,jb,1)**2*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
|
||||
!----------------
|
||||
! Spin-down part
|
||||
!----------------
|
||||
|
||||
! Occupied part of the correlation self-energy
|
||||
|
||||
do p=nC(2)+1,nBas-nR(2)
|
||||
do i=nC(2)+1,nO(2)
|
||||
do jb=1,nSt
|
||||
eps = e(p,2) - e(i,2) + Omega(jb)
|
||||
SigC(p,2) = SigC(p,2) + rho(p,i,jb,2)**2*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
|
||||
! Virtual part of the correlation self-energy
|
||||
|
||||
do p=nC(2)+1,nBas-nR(2)
|
||||
do a=nO(2)+1,nBas-nR(2)
|
||||
do jb=1,nSt
|
||||
eps = e(p,2) - e(a,2) - Omega(jb)
|
||||
SigC(p,2) = SigC(p,2) + rho(p,a,jb,2)**2*eps/(eps**2 + eta**2)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
|
||||
end subroutine unrestricted_self_energy_correlation_diag
|
Loading…
Reference in New Issue
Block a user