4
1
mirror of https://github.com/pfloos/quack synced 2024-06-02 11:25:28 +02:00
quack/src/GW/print_UG0W0.f90

69 lines
3.3 KiB
Fortran
Raw Normal View History

2020-11-04 14:47:45 +01:00
subroutine print_UG0W0(nBas,nO,eHF,ENuc,EUHF,SigC,Z,eGW,EcRPA,EcGM)
2020-09-21 23:04:26 +02:00
! Print one-electron energies and other stuff for G0W0
implicit none
include 'parameters.h'
2020-09-22 22:13:51 +02:00
integer,intent(in) :: nBas
integer,intent(in) :: nO(nspin)
2020-09-21 23:04:26 +02:00
double precision,intent(in) :: ENuc
2020-11-03 22:07:16 +01:00
double precision,intent(in) :: EUHF
2020-09-21 23:04:26 +02:00
double precision,intent(in) :: EcRPA
2020-11-04 14:47:45 +01:00
double precision,intent(in) :: EcGM(nspin)
2020-11-03 22:07:16 +01:00
double precision,intent(in) :: eHF(nBas,nspin)
2020-09-22 22:13:51 +02:00
double precision,intent(in) :: SigC(nBas,nspin)
double precision,intent(in) :: Z(nBas,nspin)
double precision,intent(in) :: eGW(nBas,nspin)
integer :: p
integer :: ispin
2023-12-26 12:37:56 +01:00
double precision :: eHOMO(nspin)
double precision :: eLUMO(nspin)
double precision :: Gap
2020-09-21 23:04:26 +02:00
! HOMO and LUMO
do ispin=1,nspin
2023-12-26 12:37:56 +01:00
eHOMO(ispin) = maxval(eGW(1:nO(ispin),ispin))
eLUMO(ispin) = minval(eGW(nO(ispin)+1:nBas,ispin))
end do
2023-12-26 12:37:56 +01:00
Gap = minval(eLUMO) -maxval(eHOMO)
2020-09-21 23:04:26 +02:00
! Dump results
2023-11-15 00:52:12 +01:00
write(*,*)'----------------------------------------------------------------'// &
'----------------------------------------------------------------'
2023-11-27 23:25:10 +01:00
write(*,*)' G0W0@UHF calculation '
2023-11-15 00:52:12 +01:00
write(*,*)'----------------------------------------------------------------'// &
'----------------------------------------------------------------'
2020-09-24 14:39:37 +02:00
write(*,'(A1,A3,A1,A30,A1,A30,A1,A30,A1,A30,A1)') &
2023-11-27 23:25:10 +01:00
'|',' ','|','e_HF (eV) ','|','Sig_GW (eV) ','|','Z ','|','e_GW (eV) ','|'
2020-09-22 22:13:51 +02:00
write(*,'(A1,A3,A1,2A15,A1,2A15,A1,2A15,A1,2A15,A1)') &
2020-09-24 14:39:37 +02:00
'|','#','|','up ','dw ','|','up ','dw ','|','up ','dw ','|','up ','dw ','|'
2023-11-15 00:52:12 +01:00
write(*,*)'----------------------------------------------------------------'// &
'----------------------------------------------------------------'
2020-09-22 22:13:51 +02:00
do p=1,nBas
write(*,'(A1,I3,A1,2F15.6,A1,2F15.6,A1,2F15.6,A1,2F15.6,A1)') &
2020-11-03 22:07:16 +01:00
'|',p,'|',eHF(p,1)*HaToeV,eHF(p,2)*HaToeV,'|',SigC(p,1)*HaToeV,SigC(p,2)*HaToeV,'|', &
2020-09-22 22:13:51 +02:00
Z(p,1),Z(p,2),'|',eGW(p,1)*HaToeV,eGW(p,2)*HaToeV,'|'
2023-12-03 18:47:30 +01:00
end do
2020-09-21 23:04:26 +02:00
2023-11-15 00:52:12 +01:00
write(*,*)'----------------------------------------------------------------'// &
'----------------------------------------------------------------'
2023-12-26 12:37:56 +01:00
write(*,'(2X,A60,F15.6,A3)') 'G0W0@UHF HOMO energy = ',maxval(eHOMO)*HaToeV,' eV'
write(*,'(2X,A60,F15.6,A3)') 'G0W0@UHF LUMO energy = ',minval(eLUMO)*HaToeV,' eV'
write(*,'(2X,A60,F15.6,A3)') 'G0W0@UHF HOMO-LUMO gap = ',Gap*HaToeV,' eV'
2023-11-15 00:52:12 +01:00
write(*,*)'----------------------------------------------------------------'// &
'----------------------------------------------------------------'
2023-11-28 14:44:26 +01:00
write(*,'(2X,A60,F15.6,A3)') 'phRPA@G0W0@UHF total energy = ',ENuc + EUHF + EcRPA,' au'
2023-11-27 23:25:10 +01:00
write(*,'(2X,A60,F15.6,A3)') 'phRPA@G0W0@UHF correlation energy = ',EcRPA,' au'
2023-11-28 14:44:26 +01:00
write(*,'(2X,A60,F15.6,A3)') ' GM@G0W0@UHF total energy = ',ENuc + EUHF + sum(EcGM),' au'
2023-11-27 23:25:10 +01:00
write(*,'(2X,A60,F15.6,A3)') ' GM@G0W0@UHF correlation energy = ',sum(EcGM),' au'
2023-11-15 00:52:12 +01:00
write(*,*)'----------------------------------------------------------------'// &
'----------------------------------------------------------------'
2020-09-21 23:04:26 +02:00
write(*,*)
2023-07-28 15:00:17 +02:00
end subroutine