4
1
mirror of https://github.com/pfloos/quack synced 2024-06-26 15:12:17 +02:00
quack/src/QuAcK/G0F2.f90

100 lines
2.1 KiB
Fortran
Raw Normal View History

2020-03-19 10:21:18 +01:00
subroutine G0F2(linearize,nBas,nC,nO,nV,nR,V,e0)
! Perform a one-shot second-order Green function calculation
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: linearize
integer,intent(in) :: nBas
integer,intent(in) :: nO
integer,intent(in) :: nC
integer,intent(in) :: nV
integer,intent(in) :: nR
double precision,intent(in) :: e0(nBas)
double precision,intent(in) :: V(nBas,nBas,nBas,nBas)
! Local variables
double precision :: eps
2020-03-19 10:59:20 +01:00
double precision :: VV
2020-03-19 10:21:18 +01:00
double precision,allocatable :: eGF2(:)
2020-03-21 16:31:39 +01:00
double precision,allocatable :: Sig(:)
2020-03-19 10:21:18 +01:00
double precision,allocatable :: Z(:)
integer :: i,j,a,b,p
! Hello world
write(*,*)
write(*,*)'************************************************'
write(*,*)'| One-shot second-order Green function |'
write(*,*)'************************************************'
write(*,*)
! Memory allocation
2020-03-21 16:31:39 +01:00
allocate(Sig(nBas),Z(nBas),eGF2(nBas))
2020-03-19 10:21:18 +01:00
2020-03-21 16:31:39 +01:00
if(linearize) then
write(*,*) '*** Quasiparticle equation will be linearized ***'
write(*,*)
2020-03-19 10:21:18 +01:00
2020-03-21 16:31:39 +01:00
end if
2020-03-19 10:21:18 +01:00
2020-03-21 16:31:39 +01:00
! Frequency-dependent second-order contribution
Sig(:) = 0d0
Z(:) = 0d0
2020-03-19 10:21:18 +01:00
do p=nC+1,nBas-nR
do i=nC+1,nO
do j=nC+1,nO
do a=nO+1,nBas-nR
2020-03-19 10:59:20 +01:00
eps = e0(p) + e0(a) - e0(i) - e0(j)
VV = (2d0*V(p,a,i,j) - V(p,a,j,i))*V(p,a,i,j)
2020-03-21 16:31:39 +01:00
Sig(p) = Sig(p) + VV/eps
Z(p) = Z(p) + VV/eps**2
2020-03-19 10:21:18 +01:00
end do
end do
end do
end do
do p=nC+1,nBas-nR
do i=nC+1,nO
do a=nO+1,nBas-nR
do b=nO+1,nBas-nR
2020-03-19 10:59:20 +01:00
eps = e0(p) + e0(i) - e0(a) - e0(b)
VV = (2d0*V(p,i,a,b) - V(p,i,b,a))*V(p,i,a,b)
2020-03-21 16:31:39 +01:00
Sig(p) = Sig(p) + VV/eps
Z(p) = Z(p) + VV/eps**2
2020-03-19 10:21:18 +01:00
end do
end do
end do
end do
2020-03-19 10:59:20 +01:00
Z(:) = 1d0/(1d0 + Z(:))
2020-03-19 10:21:18 +01:00
if(linearize) then
2020-03-21 16:31:39 +01:00
eGF2(:) = e0(:) + Z(:)*Sig(:)
2020-03-19 10:21:18 +01:00
else
2020-03-21 16:31:39 +01:00
eGF2(:) = e0(:) + Sig(:)
2020-03-19 10:21:18 +01:00
end if
2020-03-21 16:31:39 +01:00
2020-03-19 10:21:18 +01:00
! Print results
2020-03-21 16:31:39 +01:00
call print_G0F2(nBas,nO,e0,Sig,eGF2,Z)
2020-03-19 10:21:18 +01:00
end subroutine G0F2