4
1
mirror of https://github.com/pfloos/quack synced 2024-07-25 12:17:35 +02:00

dyn corr for BSE

This commit is contained in:
Pierre-Francois Loos 2020-04-22 09:55:58 +02:00
parent a2c2eb32b5
commit e855727d8a
4 changed files with 32 additions and 27 deletions

View File

@ -57,7 +57,7 @@ subroutine Bethe_Salpeter(TDA,singlet_manifold,triplet_manifold,eta, &
! Compute dynamic correction for BSE via perturbation theory
call BSE_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA(:,ispin),OmBSE(:,ispin), &
call Bethe_Salpeter_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA(:,ispin),OmBSE(:,ispin), &
XpY(:,:,ispin),XmY(:,:,ispin),rho(:,:,:,ispin))
end if
@ -82,7 +82,7 @@ subroutine Bethe_Salpeter(TDA,singlet_manifold,triplet_manifold,eta, &
! Compute dynamic correction for BSE via perturbation theory
call BSE_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA(:,ispin),OmBSE(:,ispin), &
call Bethe_Salpeter_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA(:,ispin),OmBSE(:,ispin), &
XpY(:,:,ispin),XmY(:,:,ispin),rho(:,:,:,ispin))

View File

@ -11,7 +11,7 @@ subroutine Bethe_Salpeter_A_matrix_dynamic(eta,nBas,nC,nO,nV,nR,nS,lambda,OmRPA,
double precision,intent(in) :: eta
double precision,intent(in) :: lambda
double precision,intent(in) :: OmRPA(nS)
double precision,intent(in) :: OmBSE(nS)
double precision,intent(in) :: OmBSE
double precision,intent(in) :: rho(nBas,nBas,nS)
! Local variables
@ -48,11 +48,11 @@ subroutine Bethe_Salpeter_A_matrix_dynamic(eta,nBas,nC,nO,nV,nR,nS,lambda,OmRPA,
chi = 0d0
do kc=1,nS
eps = (OmBSE(kc) - OmRPA(kc))**2 + eta**2
chi = chi + rho(i,j,kc)*rho(a,b,kc)*(OmBSE(kc) - OmRPA(kc))/eps
eps = (OmBSE - OmRPA(kc))**2 + eta**2
chi = chi + rho(i,j,kc)*rho(a,b,kc)*(OmBSE - OmRPA(kc))/eps
eps = (OmBSE(kc) + OmRPA(kc))**2 + eta**2
chi = chi - rho(i,j,kc)*rho(a,b,kc)*(OmBSE(kc) + OmRPA(kc))/eps
eps = (OmBSE + OmRPA(kc))**2 + eta**2
chi = chi - rho(i,j,kc)*rho(a,b,kc)*(OmBSE + OmRPA(kc))/eps
enddo

View File

@ -11,7 +11,7 @@ subroutine Bethe_Salpeter_B_matrix_dynamic(eta,nBas,nC,nO,nV,nR,nS,lambda,OmRPA,
double precision,intent(in) :: eta
double precision,intent(in) :: lambda
double precision,intent(in) :: OmRPA(nS)
double precision,intent(in) :: OmBSE(nS)
double precision,intent(in) :: OmBSE
double precision,intent(in) :: rho(nBas,nBas,nS)
! Local variables
@ -48,11 +48,11 @@ subroutine Bethe_Salpeter_B_matrix_dynamic(eta,nBas,nC,nO,nV,nR,nS,lambda,OmRPA,
chi = 0d0
do kc=1,nS
eps = (OmBSE(kc) - OmRPA(kc))**2 + eta**2
chi = chi + rho(i,b,kc)*rho(a,j,kc)*(OmBSE(kc) - OmRPA(kc))/eps
eps = (OmBSE - OmRPA(kc))**2 + eta**2
chi = chi + rho(i,b,kc)*rho(a,j,kc)*(OmBSE - OmRPA(kc))/eps
eps = (OmBSE(kc) + OmRPA(kc))**2 + eta**2
chi = chi - rho(i,b,kc)*rho(a,j,kc)*(OmBSE(kc) + OmRPA(kc))/eps
eps = (OmBSE + OmRPA(kc))**2 + eta**2
chi = chi - rho(i,b,kc)*rho(a,j,kc)*(OmBSE + OmRPA(kc))/eps
enddo

View File

@ -1,4 +1,4 @@
subroutine BSE_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA,OmBSE,XpY,XmY,rho)
subroutine Bethe_Salpeter_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA,OmBSE,XpY,XmY,rho)
! Compute dynamical effects via perturbation theory for BSE
@ -35,22 +35,27 @@ subroutine BSE_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA,OmBSE,XpY,
write(*,*) '-----------------------------------------------------------------------------------------------------------'
write(*,'(2X,A5,1X,A30,1X,A30,1X,A30)') '#','Static excitation (eV)','Dynamic correction (eV)','Dynamic excitation (eV)'
write(*,*) '-----------------------------------------------------------------------------------------------------------'
do ia=1,min(nS,maxS)
do ia=1,maxS
X(:) = 0.5d0*(XpY(ia,:) + XmY(ia,:))
Y(:) = 0.5d0*(XpY(ia,:) - XmY(ia,:))
X(:) = 0.5d0*(XpY(:,ia) + XmY(:,ia))
Y(:) = 0.5d0*(XpY(:,ia) - XmY(:,ia))
call Bethe_Salpeter_A_matrix_dynamic(eta,nBas,nC,nO,nV,nR,nS,1d0,OmRPA(:),OmBSE(:),rho(:,:,:),A_dyn(:,:))
call Bethe_Salpeter_A_matrix_dynamic(eta,nBas,nC,nO,nV,nR,nS,1d0,OmRPA(:),OmBSE(ia),rho(:,:,:),A_dyn(:,:))
if(TDA) then
B_dyn(:,:) = 0d0
else
call Bethe_Salpeter_B_matrix_dynamic(eta,nBas,nC,nO,nV,nR,nS,1d0,OmRPA(:),OmBSE(:),rho(:,:,:),B_dyn(:,:))
end if
OmDyn(ia) = dot_product(X(:),matmul(A_dyn(:,:),X(:))) - dot_product(Y(:),matmul(A_dyn(:,:),Y(:))) &
+ dot_product(X(:),matmul(B_dyn(:,:),Y(:))) - dot_product(Y(:),matmul(B_dyn(:,:),X(:)))
OmDyn(ia) = dot_product(X(:),matmul(A_dyn(:,:),X(:)))
else
call Bethe_Salpeter_B_matrix_dynamic(eta,nBas,nC,nO,nV,nR,nS,1d0,OmRPA(:),OmBSE(ia),rho(:,:,:),B_dyn(:,:))
OmDyn(ia) = dot_product(X(:),matmul(A_dyn(:,:),X(:))) &
- dot_product(Y(:),matmul(A_dyn(:,:),Y(:))) &
+ dot_product(X(:),matmul(B_dyn(:,:),Y(:))) &
- dot_product(Y(:),matmul(B_dyn(:,:),X(:)))
end if
write(*,'(2X,I5,15X,F15.6,15X,F15.6,15X,F15.6)') ia,OmBSE(ia)*HaToeV,OmDyn(ia)*HaToeV,(OmBSE(ia)+OmDyn(ia))*HaToeV
@ -58,4 +63,4 @@ subroutine BSE_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA,OmBSE,XpY,
write(*,*) '-----------------------------------------------------------------------------------------------------------'
write(*,*)
end subroutine BSE_dynamic_perturbation
end subroutine Bethe_Salpeter_dynamic_perturbation