From 599182f5e5088e4650d2982e0c49c52144fc0c40 Mon Sep 17 00:00:00 2001 From: Pierre-Francois Loos Date: Tue, 19 Mar 2019 11:21:34 +0100 Subject: [PATCH] UHF is working --- examples/molecule.Be | 4 +-- examples/molecule.He | 4 +-- examples/molecule.Ne | 4 +-- input/basis | 34 ++++++++++++++---- input/methods | 4 +-- input/molecule | 6 ++-- input/weight | 34 ++++++++++++++---- src/QuAcK/DIIS_extrapolation.f90 | 11 ++---- src/QuAcK/QuAcK.f90 | 14 ++------ src/QuAcK/RHF.f90 | 20 +++++++---- src/QuAcK/UHF.f90 | 48 +++++++++++++------------- src/QuAcK/exchange_matrix_AO_basis.f90 | 4 +-- src/QuAcK/print_UHF.f90 | 39 +++++++++------------ src/QuAcK/print_qsGW.f90 | 2 +- src/QuAcK/qsGW.f90 | 2 +- src/QuAcK/read_molecule.f90 | 33 +++++++++++++++--- 16 files changed, 157 insertions(+), 106 deletions(-) diff --git a/examples/molecule.Be b/examples/molecule.Be index 3473f29..7c4725d 100644 --- a/examples/molecule.Be +++ b/examples/molecule.Be @@ -1,4 +1,4 @@ -# nAt nEl nEla nElb nCore nRyd - 1 4 2 2 0 0 +# nAt nEla nElb nCore nRyd + 1 2 2 0 0 # Znuc x y z Be 0.0 0.0 0.0 diff --git a/examples/molecule.He b/examples/molecule.He index e84372f..c78e87e 100644 --- a/examples/molecule.He +++ b/examples/molecule.He @@ -1,4 +1,4 @@ -# nAt nEl nEla nElb nCore nRyd - 1 2 1 1 0 0 +# nAt nEla nElb nCore nRyd + 1 1 1 0 0 # Znuc x y z He 0.0 0.0 0.0 diff --git a/examples/molecule.Ne b/examples/molecule.Ne index 29ad6e0..edeba31 100644 --- a/examples/molecule.Ne +++ b/examples/molecule.Ne @@ -1,4 +1,4 @@ -# nAt nEl nEla nElb nCore nRyd - 1 10 5 5 0 0 +# nAt nEla nElb nCore nRyd + 1 5 5 0 0 # Znuc x y z Ne 0.0 0.0 0.0 diff --git a/input/basis b/input/basis index b9ca7b5..17fdf8c 100644 --- a/input/basis +++ b/input/basis @@ -1,9 +1,29 @@ -1 3 -S 3 1.00 - 38.3600000 0.0238090 - 5.7700000 0.1548910 - 1.2400000 0.4699870 +1 6 +S 8 1.00 + 1469.0000000 0.0007660 + 220.5000000 0.0058920 + 50.2600000 0.0296710 + 14.2400000 0.1091800 + 4.5810000 0.2827890 + 1.5800000 0.4531230 + 0.5640000 0.2747740 + 0.0734500 0.0097510 +S 8 1.00 + 1469.0000000 -0.0001200 + 220.5000000 -0.0009230 + 50.2600000 -0.0046890 + 14.2400000 -0.0176820 + 4.5810000 -0.0489020 + 1.5800000 -0.0960090 + 0.5640000 -0.1363800 + 0.0734500 0.5751020 S 1 1.00 - 0.2976000 1.0000000 + 0.0280500 1.0000000 +P 3 1.00 + 1.5340000 0.0227840 + 0.2749000 0.1391070 + 0.0736200 0.5003750 P 1 1.00 - 1.2750000 1.0000000 + 0.0240300 1.0000000 +D 1 1.00 + 0.1239000 1.0000000 diff --git a/input/methods b/input/methods index 3713171..31b24aa 100644 --- a/input/methods +++ b/input/methods @@ -1,9 +1,9 @@ # RHF UHF MOM F T F # MP2 MP3 MP2-F12 - T T F + F F F # CCD CCSD CCSD(T) - T F F + F F F # CIS TDHF ADC F F F # GF2 GF3 diff --git a/input/molecule b/input/molecule index 0a9db5b..058d6dd 100644 --- a/input/molecule +++ b/input/molecule @@ -1,4 +1,4 @@ -# nAt nEl nEla nElb nCore nRyd - 1 2 0 0 0 0 +# nAt nEla nElb nCore nRyd + 1 2 1 0 0 # Znuc x y z - He 0.0 0.0 0.0 + Li 0.0 0.0 0.0 diff --git a/input/weight b/input/weight index b9ca7b5..17fdf8c 100644 --- a/input/weight +++ b/input/weight @@ -1,9 +1,29 @@ -1 3 -S 3 1.00 - 38.3600000 0.0238090 - 5.7700000 0.1548910 - 1.2400000 0.4699870 +1 6 +S 8 1.00 + 1469.0000000 0.0007660 + 220.5000000 0.0058920 + 50.2600000 0.0296710 + 14.2400000 0.1091800 + 4.5810000 0.2827890 + 1.5800000 0.4531230 + 0.5640000 0.2747740 + 0.0734500 0.0097510 +S 8 1.00 + 1469.0000000 -0.0001200 + 220.5000000 -0.0009230 + 50.2600000 -0.0046890 + 14.2400000 -0.0176820 + 4.5810000 -0.0489020 + 1.5800000 -0.0960090 + 0.5640000 -0.1363800 + 0.0734500 0.5751020 S 1 1.00 - 0.2976000 1.0000000 + 0.0280500 1.0000000 +P 3 1.00 + 1.5340000 0.0227840 + 0.2749000 0.1391070 + 0.0736200 0.5003750 P 1 1.00 - 1.2750000 1.0000000 + 0.0240300 1.0000000 +D 1 1.00 + 0.1239000 1.0000000 diff --git a/src/QuAcK/DIIS_extrapolation.f90 b/src/QuAcK/DIIS_extrapolation.f90 index 4fb89dc..e53db06 100644 --- a/src/QuAcK/DIIS_extrapolation.f90 +++ b/src/QuAcK/DIIS_extrapolation.f90 @@ -1,4 +1,4 @@ -subroutine DIIS_extrapolation(n_err,n_e,n_diis,error,e,error_in,e_inout) +subroutine DIIS_extrapolation(rcond,n_err,n_e,n_diis,error,e,error_in,e_inout) ! Perform DIIS extrapolation @@ -13,11 +13,11 @@ subroutine DIIS_extrapolation(n_err,n_e,n_diis,error,e,error_in,e_inout) ! Local variables - double precision :: rcond double precision,allocatable :: A(:,:),b(:),w(:) ! Output variables + double precision,intent(out) :: rcond integer,intent(inout) :: n_diis double precision,intent(inout):: e_inout(n_e) @@ -49,13 +49,6 @@ subroutine DIIS_extrapolation(n_err,n_e,n_diis,error,e,error_in,e_inout) ! Extrapolate - if(rcond > 1d-14) then - e_inout(:) = matmul(w(1:n_diis),transpose(e(:,1:n_diis))) - else - - n_diis = 0 - - endif end subroutine DIIS_extrapolation diff --git a/src/QuAcK/QuAcK.f90 b/src/QuAcK/QuAcK.f90 index 55c8f82..e17b337 100644 --- a/src/QuAcK/QuAcK.f90 +++ b/src/QuAcK/QuAcK.f90 @@ -122,7 +122,7 @@ program QuAcK ! nS = number of single excitation ! = nO*nV - call read_molecule(nNuc,nEl,nO,nC,nR) + call read_molecule(nNuc,nEl(:),nO(:),nC(:),nR(:)) allocate(ZNuc(nNuc),rNuc(nNuc,3)) ! Read geometry @@ -178,10 +178,7 @@ program QuAcK if(doRHF) then call cpu_time(start_HF) -! call SPHF(maxSCF_HF,thresh_HF,n_diis_HF,guess_type, & -! nBas,nEl,S,T,V,Hc,ERI_AO_basis,X,ENuc,ERHF,cHF,eHF,PHF) - call RHF(maxSCF_HF,thresh_HF,n_diis_HF,guess_type, & - nBas,nO,S,T,V,Hc,ERI_AO_basis,X,ENuc,ERHF,cHF,eHF,PHF) + call RHF(maxSCF_HF,thresh_HF,n_diis_HF,guess_type,nBas,nO,S,T,V,Hc,ERI_AO_basis,X,ENuc,ERHF,eHF,cHF,PHF) call cpu_time(end_HF) t_HF = end_HF - start_HF @@ -196,11 +193,8 @@ program QuAcK if(doUHF) then - nO(2) = nO(1) - call cpu_time(start_HF) - call UHF(maxSCF_HF,thresh_HF,n_diis_HF,guess_type, & - nBas,nO,S,T,V,Hc,ERI_AO_basis,X,ENuc,EUHF) + call UHF(maxSCF_HF,thresh_HF,n_diis_HF,guess_type,nBas,nO,S,T,V,Hc,ERI_AO_basis,X,ENuc,EUHF,eHF,cHF,PHF) call cpu_time(end_HF) t_HF = end_HF - start_HF @@ -239,7 +233,6 @@ program QuAcK if(doMP2) then call cpu_time(start_MP2) -! call SPMP2(nBas,nC,nEl,nBas-nEl,nR,ERI_MO_basis,ENuc,ERHF,eHF,EcMP2) call MP2(nBas,nC,nO,nV,nR,ERI_MO_basis,ENuc,ERHF,eHF,EcMP2) call cpu_time(end_MP2) @@ -342,7 +335,6 @@ program QuAcK call cpu_time(start_TDHF) call TDHF(singlet_manifold,triplet_manifold,nBas,nC,nO,nV,nR,nS,ERI_MO_basis,eHF) -! call SPTDHF(singlet_manifold,triplet_manifold,nBas,nC,nEl,nBas-nEl,nR,nEl*(nBas-nEl),ERI_MO_basis,eHF) call cpu_time(end_TDHF) t_TDHF = end_TDHF - start_TDHF diff --git a/src/QuAcK/RHF.f90 b/src/QuAcK/RHF.f90 index 25621d0..387de91 100644 --- a/src/QuAcK/RHF.f90 +++ b/src/QuAcK/RHF.f90 @@ -1,4 +1,4 @@ -subroutine RHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,ERHF,c,e,P) +subroutine RHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,ERHF,e,c,P) ! Perform restricted Hartree-Fock calculation @@ -18,13 +18,17 @@ subroutine RHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,ERH integer :: nSCF,nBasSq,n_diis double precision :: ET,EV,EJ,EK,Conv,Gap + double precision :: rcond double precision,external :: trace_matrix double precision,allocatable :: error(:,:),error_diis(:,:),F_diis(:,:) double precision,allocatable :: J(:,:),K(:,:),cp(:,:),F(:,:),Fp(:,:) ! Output variables - double precision,intent(out) :: ERHF,c(nBas,nBas),e(nBas),P(nBas,nBas) + double precision,intent(out) :: ERHF + double precision,intent(out) :: e(nBas) + double precision,intent(out) :: c(nBas,nBas) + double precision,intent(out) :: P(nBas,nBas) ! Hello world @@ -92,7 +96,7 @@ subroutine RHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,ERH call Coulomb_matrix_AO_basis(nBas,P,ERI,J) call exchange_matrix_AO_basis(nBas,P,ERI,K) - F(:,:) = Hc(:,:) + J(:,:) + K(:,:) + F(:,:) = Hc(:,:) + J(:,:) + 0.5d0*K(:,:) ! Check convergence @@ -102,7 +106,11 @@ subroutine RHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,ERH ! DIIS extrapolation n_diis = min(n_diis+1,max_diis) - call DIIS_extrapolation(nBasSq,nBasSq,n_diis,error_diis,F_diis,error,F) + call DIIS_extrapolation(rcond,nBasSq,nBasSq,n_diis,error_diis,F_diis,error,F) + +! Reset DIIS if required + + if(abs(rcond) < 1d-15) n_diis = 0 ! Diagonalize Fock matrix @@ -119,7 +127,7 @@ subroutine RHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,ERH ERHF = trace_matrix(nBas,matmul(P,Hc)) & + 0.5d0*trace_matrix(nBas,matmul(P,J)) & - + 0.5d0*trace_matrix(nBas,matmul(P,K)) + + 0.25d0*trace_matrix(nBas,matmul(P,K)) ! Compute HOMO-LUMO gap @@ -163,7 +171,7 @@ subroutine RHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,ERH ET = trace_matrix(nBas,matmul(P,T)) EV = trace_matrix(nBas,matmul(P,V)) EJ = 0.5d0*trace_matrix(nBas,matmul(P,J)) - EK = 0.5d0*trace_matrix(nBas,matmul(P,K)) + EK = 0.25d0*trace_matrix(nBas,matmul(P,K)) ERHF = ET + EV + EJ + EK call print_RHF(nBas,nO,e,C,ENuc,ET,EV,EJ,EK,ERHF) diff --git a/src/QuAcK/UHF.f90 b/src/QuAcK/UHF.f90 index dc230eb..28fc3c8 100644 --- a/src/QuAcK/UHF.f90 +++ b/src/QuAcK/UHF.f90 @@ -1,4 +1,4 @@ -subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUHF) +subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUHF,e,c,P) ! Perform unrestricted Hartree-Fock calculation @@ -33,25 +33,26 @@ subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUH double precision :: EV(nspin) double precision :: EJ(nsp) double precision :: Ex(nspin) - double precision :: Ec(nsp) - double precision :: EUHF - double precision,allocatable :: eps(:,:) - double precision,allocatable :: c(:,:,:) double precision,allocatable :: cp(:,:,:) double precision,allocatable :: J(:,:,:) double precision,allocatable :: F(:,:,:) double precision,allocatable :: Fp(:,:,:) - double precision,allocatable :: Fx(:,:,:) + double precision,allocatable :: K(:,:,:) double precision,allocatable :: err(:,:,:) double precision,allocatable :: err_diis(:,:,:) double precision,allocatable :: F_diis(:,:,:) double precision,external :: trace_matrix - double precision,allocatable :: P(:,:,:) - integer :: ispin +! Output variables + + double precision,intent(out) :: EUHF + double precision,intent(out) :: e(nBas,nspin) + double precision,intent(out) :: c(nBas,nBas,nspin) + double precision,intent(out) :: P(nBas,nBas,nspin) + ! Hello world write(*,*) @@ -66,9 +67,8 @@ subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUH ! Memory allocation - allocate(eps(nBas,nspin),c(nBas,nBas,nspin),cp(nBas,nBas,nspin), & - J(nBas,nBas,nspin),F(nBas,nBas,nspin),Fp(nBas,nBas,nspin), & - Fx(nBas,nBas,nspin),err(nBas,nBas,nspin),P(nBas,nBas,nspin), & + allocate(J(nBas,nBas,nspin),F(nBas,nBas,nspin),Fp(nBas,nBas,nspin), & + K(nBas,nBas,nspin),err(nBas,nBas,nspin),cp(nBas,nBas,nspin), & err_diis(nBasSq,max_diis,nspin),F_diis(nBasSq,max_diis,nspin)) ! Guess coefficients and eigenvalues @@ -101,10 +101,10 @@ subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUH !------------------------------------------------------------------------ write(*,*) - write(*,*)'------------------------------------------------------------------------------------------' - write(*,'(1X,A1,1X,A3,1X,A1,1X,A16,1X,A1,1X,A16,1X,A1,1X,A16,1X,A1,1X,A10,1X,A1,1X)') & - '|','#','|','E(KS)','|','Ex(KS)','|','Ec(KS)','|','Conv','|' - write(*,*)'------------------------------------------------------------------------------------------' + write(*,*)'----------------------------------------------------------' + write(*,'(1X,A1,1X,A3,1X,A1,1X,A16,1X,A1,1X,A16,1X,A1,1X,A10,1X,A1,1X)') & + '|','#','|','E(UHF)','|','Ex(UHF)','|','Conv','|' + write(*,*)'----------------------------------------------------------' do while(conv > thresh .and. nSCF < maxSCF) @@ -122,7 +122,7 @@ subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUH cp(:,:,:) = Fp(:,:,:) do ispin=1,nspin - call diagonalize_matrix(nBas,cp(:,:,ispin),eps(:,ispin)) + call diagonalize_matrix(nBas,cp(:,:,ispin),e(:,ispin)) end do ! Back-transform eigenvectors in non-orthogonal basis @@ -146,12 +146,12 @@ subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUH ! Compute exchange potential do ispin=1,nspin - call exchange_matrix_AO_basis(nBas,P(:,:,ispin),ERI(:,:,:,:),Fx(:,:,ispin)) + call exchange_matrix_AO_basis(nBas,P(:,:,ispin),ERI(:,:,:,:),K(:,:,ispin)) end do ! Build Fock operator do ispin=1,nspin - F(:,:,ispin) = Hc(:,:) + J(:,:,ispin) + J(:,:,mod(ispin,2)+1) + Fx(:,:,ispin) + F(:,:,ispin) = Hc(:,:) + J(:,:,ispin) + J(:,:,mod(ispin,2)+1) + K(:,:,ispin) end do ! Check convergence @@ -198,20 +198,20 @@ subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUH ! Exchange energy do ispin=1,nspin - Ex(ispin) = trace_matrix(nBas,matmul(P(:,:,ispin),Fx(:,:,ispin))) + Ex(ispin) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,ispin),K(:,:,ispin))) end do ! Total energy - EUHF = sum(ET(:)) + sum(EV(:)) + sum(EJ(:)) + sum(Ex(:)) + sum(Ec(:)) + EUHF = sum(ET(:)) + sum(EV(:)) + sum(EJ(:)) + sum(Ex(:)) ! Dump results - write(*,'(1X,A1,1X,I3,1X,A1,1X,F16.10,1X,A1,1X,F16.10,1X,A1,1X,F16.10,1X,A1,1X,F10.6,1X,A1,1X)') & - '|',nSCF,'|',EUHF + ENuc,'|',sum(Ex(:)),'|',sum(Ec(:)),'|',conv,'|' + write(*,'(1X,A1,1X,I3,1X,A1,1X,F16.10,1X,A1,1X,F16.10,1X,A1,1X,F10.6,1X,A1,1X)') & + '|',nSCF,'|',EUHF + ENuc,'|',sum(Ex(:)),'|',conv,'|' end do - write(*,*)'------------------------------------------------------------------------------------------' + write(*,*)'----------------------------------------------------------' !------------------------------------------------------------------------ ! End of SCF loop !------------------------------------------------------------------------ @@ -232,6 +232,6 @@ subroutine UHF(maxSCF,thresh,max_diis,guess_type,nBas,nO,S,T,V,Hc,ERI,X,ENuc,EUH ! Compute final UHF energy - call print_UHF(nBas,nO(:),eps(:,:),c(:,:,:),ENuc,ET(:),EV(:),EJ(:),Ex(:),Ec(:),EUHF) + call print_UHF(nBas,nO(:),e(:,:),c(:,:,:),ENuc,ET(:),EV(:),EJ(:),Ex(:),EUHF) end subroutine UHF diff --git a/src/QuAcK/exchange_matrix_AO_basis.f90 b/src/QuAcK/exchange_matrix_AO_basis.f90 index b60a59a..2f0f1d2 100644 --- a/src/QuAcK/exchange_matrix_AO_basis.f90 +++ b/src/QuAcK/exchange_matrix_AO_basis.f90 @@ -24,12 +24,10 @@ subroutine exchange_matrix_AO_basis(nBas,P,G,K) do si=1,nBas do la=1,nBas do mu=1,nBas - K(mu,nu) = K(mu,nu) + P(la,si)*G(mu,la,si,nu) + K(mu,nu) = K(mu,nu) - P(la,si)*G(mu,la,si,nu) enddo enddo enddo enddo - K = -0.5d0*K - end subroutine exchange_matrix_AO_basis diff --git a/src/QuAcK/print_UHF.f90 b/src/QuAcK/print_UHF.f90 index d0a16c1..6c52e4d 100644 --- a/src/QuAcK/print_UHF.f90 +++ b/src/QuAcK/print_UHF.f90 @@ -1,4 +1,4 @@ -subroutine print_UHF(nBas,nO,eps,c,ENuc,ET,EV,EJ,Ex,Ec,Ew) +subroutine print_UHF(nBas,nO,e,c,ENuc,ET,EV,EJ,Ex,EUHF) ! Print one- and two-electron energies and other stuff for UHF calculation @@ -7,15 +7,14 @@ subroutine print_UHF(nBas,nO,eps,c,ENuc,ET,EV,EJ,Ex,Ec,Ew) integer,intent(in) :: nBas integer,intent(in) :: nO(nspin) - double precision,intent(in) :: eps(nBas,nspin) + double precision,intent(in) :: e(nBas,nspin) double precision,intent(in) :: c(nBas,nBas,nspin) double precision,intent(in) :: ENuc double precision,intent(in) :: ET(nspin) double precision,intent(in) :: EV(nspin) double precision,intent(in) :: EJ(nsp) double precision,intent(in) :: Ex(nspin) - double precision,intent(in) :: Ec(nsp) - double precision,intent(in) :: Ew + double precision,intent(in) :: EUHF integer :: HOMO(nspin) integer :: LUMO(nspin) @@ -27,8 +26,8 @@ subroutine print_UHF(nBas,nO,eps,c,ENuc,ET,EV,EJ,Ex,Ec,Ew) LUMO(:) = HOMO(:) + 1 - Gap(1) = eps(LUMO(1),1) - eps(HOMO(1),1) - Gap(2) = eps(LUMO(2),2) - eps(HOMO(2),2) + Gap(1) = e(LUMO(1),1) - e(HOMO(1),1) + Gap(2) = e(LUMO(2),2) - e(HOMO(2),2) ! Dump results @@ -47,10 +46,10 @@ subroutine print_UHF(nBas,nO,eps,c,ENuc,ET,EV,EJ,Ex,Ec,Ew) write(*,'(A40,1X,F16.10,A3)') ' Potential a energy: ',EV(1),' au' write(*,'(A40,1X,F16.10,A3)') ' Potential b energy: ',EV(2),' au' write(*,'(A60)') '-------------------------------------------------' - write(*,'(A40,1X,F16.10,A3)') ' Two-electron a energy: ',sum(EJ(:)) + sum(Ex(:)) + sum(Ec(:)),' au' - write(*,'(A40,1X,F16.10,A3)') ' Two-electron aa energy: ',EJ(1) + Ex(1) + Ec(1),' au' - write(*,'(A40,1X,F16.10,A3)') ' Two-electron ab energy: ',EJ(2) + Ec(2),' au' - write(*,'(A40,1X,F16.10,A3)') ' Two-electron bb energy: ',EJ(3) + Ex(2) + Ec(3),' au' + write(*,'(A40,1X,F16.10,A3)') ' Two-electron a energy: ',sum(EJ(:)) + sum(Ex(:)),' au' + write(*,'(A40,1X,F16.10,A3)') ' Two-electron aa energy: ',EJ(1) + Ex(1),' au' + write(*,'(A40,1X,F16.10,A3)') ' Two-electron ab energy: ',EJ(2),' au' + write(*,'(A40,1X,F16.10,A3)') ' Two-electron bb energy: ',EJ(3) + Ex(2),' au' write(*,'(A40,1X,F16.10,A3)') ' Coulomb energy: ',sum(EJ(:)),' au' write(*,'(A40,1X,F16.10,A3)') ' Coulomb aa energy: ',EJ(1),' au' write(*,'(A40,1X,F16.10,A3)') ' Coulomb ab energy: ',EJ(2),' au' @@ -58,21 +57,17 @@ subroutine print_UHF(nBas,nO,eps,c,ENuc,ET,EV,EJ,Ex,Ec,Ew) write(*,'(A40,1X,F16.10,A3)') ' Exchange energy: ',sum(Ex(:)),' au' write(*,'(A40,1X,F16.10,A3)') ' Exchange a energy: ',Ex(1),' au' write(*,'(A40,1X,F16.10,A3)') ' Exchange b energy: ',Ex(2),' au' - write(*,'(A40,1X,F16.10,A3)') ' Correlation energy: ',sum(Ec(:)),' au' - write(*,'(A40,1X,F16.10,A3)') ' Correlation aa energy: ',Ec(1),' au' - write(*,'(A40,1X,F16.10,A3)') ' Correlation ab energy: ',Ec(2),' au' - write(*,'(A40,1X,F16.10,A3)') ' Correlation bb energy: ',Ec(3),' au' write(*,'(A60)') '-------------------------------------------------' - write(*,'(A40,1X,F16.10,A3)') ' Electronic energy: ',Ew,' au' + write(*,'(A40,1X,F16.10,A3)') ' Electronic energy: ',EUHF,' au' write(*,'(A40,1X,F16.10,A3)') ' Nuclear repulsion: ',ENuc,' au' - write(*,'(A40,1X,F16.10,A3)') ' UHF energy: ',Ew + ENuc,' au' + write(*,'(A40,1X,F16.10,A3)') ' UHF energy: ',EUHF + ENuc,' au' write(*,'(A60)') '-------------------------------------------------' - write(*,'(A40,F13.6,A3)') ' UHF HOMO a energy:',eps(HOMO(1),1)*HatoeV,' eV' - write(*,'(A40,F13.6,A3)') ' UHF LUMO a energy:',eps(LUMO(1),1)*HatoeV,' eV' + write(*,'(A40,F13.6,A3)') ' UHF HOMO a energy:',e(HOMO(1),1)*HatoeV,' eV' + write(*,'(A40,F13.6,A3)') ' UHF LUMO a energy:',e(LUMO(1),1)*HatoeV,' eV' write(*,'(A40,F13.6,A3)') ' UHF HOMOa-LUMOa gap:',Gap(1)*HatoeV,' eV' write(*,'(A60)') '-------------------------------------------------' - write(*,'(A40,F13.6,A3)') ' UHF HOMO b energy:',eps(HOMO(2),2)*HatoeV,' eV' - write(*,'(A40,F13.6,A3)') ' UHF LUMO b energy:',eps(LUMO(2),2)*HatoeV,' eV' + write(*,'(A40,F13.6,A3)') ' UHF HOMO b energy:',e(HOMO(2),2)*HatoeV,' eV' + write(*,'(A40,F13.6,A3)') ' UHF LUMO b energy:',e(LUMO(2),2)*HatoeV,' eV' write(*,'(A40,F13.6,A3)') ' UHF HOMOb-LUMOb gap :',Gap(2)*HatoeV,' eV' write(*,'(A60)') '-------------------------------------------------' write(*,*) @@ -91,12 +86,12 @@ subroutine print_UHF(nBas,nO,eps,c,ENuc,ET,EV,EJ,Ex,Ec,Ew) write(*,'(A50)') '---------------------------------------' write(*,'(A50)') ' UHF spin-up orbital energies ' write(*,'(A50)') '---------------------------------------' - call matout(nBas,1,eps(:,1)) + call matout(nBas,1,e(:,1)) write(*,*) write(*,'(A50)') '---------------------------------------' write(*,'(A50)') ' UHF spin-down orbital energies ' write(*,'(A50)') '---------------------------------------' - call matout(nBas,1,eps(:,2)) + call matout(nBas,1,e(:,2)) write(*,*) end subroutine print_UHF diff --git a/src/QuAcK/print_qsGW.f90 b/src/QuAcK/print_qsGW.f90 index f92be3d..f040258 100644 --- a/src/QuAcK/print_qsGW.f90 +++ b/src/QuAcK/print_qsGW.f90 @@ -31,7 +31,7 @@ subroutine print_qsGW(nBas,nO,nSCF,Conv,thresh,eHF,eGW,c,ENuc,P,T,V,Hc,J,K,F,Sig ET = trace_matrix(nBas,matmul(P,T)) EV = trace_matrix(nBas,matmul(P,V)) EJ = 0.5d0*trace_matrix(nBas,matmul(P,J)) - Ex = 0.5d0*trace_matrix(nBas,matmul(P,K)) + Ex = 0.25d0*trace_matrix(nBas,matmul(P,K)) EqsGW = ET + EV + EJ + Ex Ec = 0d0 diff --git a/src/QuAcK/qsGW.f90 b/src/QuAcK/qsGW.f90 index 1b34714..acb1544 100644 --- a/src/QuAcK/qsGW.f90 +++ b/src/QuAcK/qsGW.f90 @@ -126,7 +126,7 @@ subroutine qsGW(maxSCF,thresh,max_diis,COHSEX,SOSEX,BSE,TDA,G0W,GW0,singlet_mani ! Solve the quasi-particle equation - F(:,:) = Hc(:,:) + J(:,:) + K(:,:) + SigCp(:,:) + F(:,:) = Hc(:,:) + J(:,:) + 0.5d0*K(:,:) + SigCp(:,:) ! Compute commutator and convergence criteria diff --git a/src/QuAcK/read_molecule.f90 b/src/QuAcK/read_molecule.f90 index 4625d03..5d8886e 100644 --- a/src/QuAcK/read_molecule.f90 +++ b/src/QuAcK/read_molecule.f90 @@ -8,7 +8,16 @@ subroutine read_molecule(nNuc,nEl,nO,nC,nR) ! Input variables - integer,intent(out) :: nNuc,nEl,nO,nC,nR + integer,intent(out) :: nNuc + integer,intent(out) :: nEl(nspin) + integer,intent(out) :: nO(nspin) + integer,intent(out) :: nC(nspin) + integer,intent(out) :: nR(nspin) + +! Local variables + + integer :: nCore + integer :: nRyd ! Open file with geometry specification @@ -17,9 +26,18 @@ subroutine read_molecule(nNuc,nEl,nO,nC,nR) ! Read number of atoms and number of electrons read(1,*) - read(1,*) nNuc,nEl,nC,nR + read(1,*) nNuc,nEl(1),nEl(2),nCore,nRyd - nO = nEl/2 + if(mod(nCore,2) /= 0 .or. mod(nRyd,2) /= 0) then + + print*, 'The number of core and Rydberg electrons must be even!' + stop + + end if + + nO(:) = nEl(:) + nC(:) = nCore/2 + nR(:) = nRyd/2 ! Print results @@ -28,7 +46,14 @@ subroutine read_molecule(nNuc,nEl,nO,nC,nR) write(*,'(A28)') '----------------------' write(*,*) write(*,'(A28)') '----------------------' - write(*,'(A28,1X,I16)') 'Number of electrons',nEl + write(*,'(A28,1X,I16)') 'Number of spin-up electrons',nEl(1) + write(*,'(A28,1X,I16)') 'Number of spin-down electrons',nEl(2) + write(*,'(A28,1X,I16)') ' Total number of electrons',sum(nEl(:)) + write(*,'(A28)') '----------------------' + write(*,*) + write(*,'(A28)') '----------------------' + write(*,'(A28,1X,I16)') 'Number of core electrons',sum(nC(:)) + write(*,'(A28,1X,I16)') 'Number of Rydberg electrons',sum(nR(:)) write(*,'(A28)') '----------------------' write(*,*)