From 9f2bb3ba7dacbdbb2c432fa767b2cab7aaff4d91 Mon Sep 17 00:00:00 2001 From: EnzoMonino Date: Fri, 25 Feb 2022 09:40:39 +0100 Subject: [PATCH] qsUGT --- input/methods | 6 +- input/options | 6 +- mol/h2.xyz | 2 +- src/GT/print_qsUGT.f90 | 82 ++++ src/GT/qsUGT.f90 | 443 ++++++++++++++++++++ src/GT/unrestricted_self_energy_Tmatrix.f90 | 207 +++++++++ src/QuAcK/QuAcK.f90 | 7 +- 7 files changed, 744 insertions(+), 9 deletions(-) create mode 100644 src/GT/print_qsUGT.f90 create mode 100644 src/GT/qsUGT.f90 create mode 100644 src/GT/unrestricted_self_energy_Tmatrix.f90 diff --git a/input/methods b/input/methods index a13b74a..b4f1b5b 100644 --- a/input/methods +++ b/input/methods @@ -1,5 +1,5 @@ # RHF UHF KS MOM - T F F F + F T F F # MP2* MP3 MP2-F12 F F F # CCD pCCD DCD CCSD CCSD(T) @@ -13,9 +13,9 @@ # G0F2* evGF2* qsGF2* G0F3 evGF3 F F F F F # G0W0* evGW* qsGW* ufG0W0 ufGW - F F F T F + F F F F F # G0T0 evGT qsGT - F F F + F F T # MCMP2 F # * unrestricted version available diff --git a/input/options b/input/options index e9b8028..5ef344f 100644 --- a/input/options +++ b/input/options @@ -1,15 +1,15 @@ # HF: maxSCF thresh DIIS n_diis guess_type ortho_type mix_guess level_shift stability - 256 0.0000001 T 5 2 1 F 0.0 F + 256 0.00001 T 5 2 1 T 0.0 F # MP: # CC: maxSCF thresh DIIS n_diis 64 0.00001 T 5 # spin: TDA singlet triplet spin_conserved spin_flip - F T T T T + T T T T T # GF: maxSCF thresh DIIS n_diis lin eta renorm reg 256 0.00001 T 5 T 0.0 3 F # GW: maxSCF thresh DIIS n_diis lin eta COHSEX SOSEX TDA_W G0W GW0 reg - 256 0.00001 T 5 T 0.0 F F F F F F + 256 0.00001 T 5 T 0.0 F F T F F F # GT: maxSCF thresh DIIS n_diis lin eta TDA_T reg 256 0.00001 T 5 T 0.0 F F # ACFDT: AC Kx XBS diff --git a/mol/h2.xyz b/mol/h2.xyz index 5bca03f..fc52308 100644 --- a/mol/h2.xyz +++ b/mol/h2.xyz @@ -1,4 +1,4 @@ 2 H 0. 0. 0. -H 0. 0. 0.30 +H 0. 0. 1.0 diff --git a/src/GT/print_qsUGT.f90 b/src/GT/print_qsUGT.f90 new file mode 100644 index 0000000..5170806 --- /dev/null +++ b/src/GT/print_qsUGT.f90 @@ -0,0 +1,82 @@ +subroutine print_qsUGT(nBas,nO,nSCF,Conv,eHF,ENuc,EUHF,SigT,Z,eGT,ET,EV,EJ,Ex,EcGM,EcRPA,EqsGT) + +! Print one-electron energies and other stuff for UG0T0 + + implicit none + include 'parameters.h' + + integer,intent(in) :: nBas + integer,intent(in) :: nO(nspin) + integer,intent(in) :: nSCF + double precision,intent(in) :: Conv + double precision,intent(in) :: ENuc + double precision,intent(in) :: EUHF + 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) :: EcGM(nspin) + double precision,intent(in) :: EcRPA(nspin) + double precision,intent(in) :: EqsGT + double precision,intent(in) :: eHF(nBas,nspin) + double precision,intent(in) :: SigT(nBas,nBas,nspin) + double precision,intent(in) :: Z(nBas,nspin) + double precision,intent(in) :: eGT(nBas,nspin) + + integer :: p + integer :: ispin + double precision :: HOMO(nspin) + double precision :: LUMO(nspin) + double precision :: Gap(nspin) + +! HOMO and LUMO + do ispin=1,nspin + if(nO(ispin) > 0) then + HOMO(ispin) = eGT(nO(ispin),ispin) + LUMO(ispin) = eGT(nO(ispin)+1,ispin) + Gap(ispin) = LUMO(ispin) - HOMO(ispin) + else + HOMO(ispin) = 0d0 + LUMO(ispin) = eGT(1,ispin) + Gap(ispin) = 0d0 + end if + end do + +! Dump results + + write(*,*)'-------------------------------------------------------------------------------' + if(nSCF < 10) then + write(*,'(1X,A21,I1,A1,I1,A12)')' Self-consistent qsG',nSCF,'T',nSCF,' calculation' + else + write(*,'(1X,A21,I2,A1,I2,A12)')' Self-consistent qsG',nSCF,'T',nSCF,' calculation' + endif + 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_T (eV)','|','Z','|','e_QP (eV)','|' + write(*,*)'-------------------------------------------------------------------------------' + + do p=1,nBas + write(*,'(A1,I3,A1,2F15.6,A1,2F15.6,A1,2F15.6,A1,2F15.6,A1)') & + '|',p,'|',eHF(p,1)*HaToeV,eHF(p,2)*HaToeV,'|',SigT(p,p,1)*HaToeV,SigT(p,p,2)*HaToeV,'|', & + Z(p,1),Z(p,2),'|',eGT(p,1)*HaToeV,eGT(p,2)*HaToeV,'|' + enddo + + write(*,*)'-------------------------------------------------------------------------------' + write(*,'(2X,A10,I3)') 'Iteration ',nSCF + write(*,'(2X,A14,F15.5)')'Convergence = ',Conv + write(*,*)'-------------------------------------------------------------------------------' + write(*,*)'-------------------------------------------------------------------------------' + write(*,'(2X,A50,F15.6,A3)') 'qsUGT HOMO energy (eV) =',maxval(HOMO(:))*HaToeV,' eV' + write(*,'(2X,A50,F15.6,A3)') 'qsUGT LUMO energy (eV) =',minval(LUMO(:))*HaToeV,' eV' + write(*,'(2X,A50,F15.6,A3)') 'qsUGT HOMO-LUMO gap (eV) =',(minval(LUMO(:))-maxval(HOMO(:)))*HaToeV,' eV' + write(*,*)'-------------------------------------------------------------------------------' +write(*,'(2X,A50,F15.6,A3)') ' qsGT total energy:',ENuc + EqsGT,' au' +write(*,'(2X,A50,F15.6,A3)') ' qsGT exchange energy:',sum(Ex(:)),' au' +write(*,'(2X,A50,F15.6,A3)') ' GM@qsGT correlation energy:',sum(EcGM(:)),' au' +write(*,'(2X,A50,F15.6,A3)') 'ppRPA@qsGT correlation energy:',sum(EcRPA(:)),' au' +write(*,*)'-------------------------------------------' +write(*,*) + +end subroutine print_qsUGT + + diff --git a/src/GT/qsUGT.f90 b/src/GT/qsUGT.f90 new file mode 100644 index 0000000..1b57c7e --- /dev/null +++ b/src/GT/qsUGT.f90 @@ -0,0 +1,443 @@ +subroutine qsUGT(maxSCF,thresh,max_diis,doACFDT,exchange_kernel,doXBS,BSE, & + TDA_T,TDA,dBSE,dTDA,evDyn,spin_conserved,spin_flip,& + eta,regularize,nBas,nC,nO,nV,nR,nS,nNuc,ZNuc,rNuc,ENuc,EUHF,S,X,T,V,Hc,ERI_AO,ERI_aaaa,& + ERI_aabb,ERI_bbbb,dipole_int_AO,dipole_int_aa,dipole_int_bb,PHF,cHF,eHF) + +! Perform a quasiparticle self-consistent GT calculation + + implicit none + include 'parameters.h' + +! Input variables + integer,intent(in) :: maxSCF + integer,intent(in) :: max_diis + double precision,intent(in) :: thresh + logical,intent(in) :: doACFDT + logical,intent(in) :: exchange_kernel + logical,intent(in) :: doXBS + logical,intent(in) :: BSE + logical,intent(in) :: TDA_T + logical,intent(in) :: TDA + logical,intent(in) :: dBSE + logical,intent(in) :: dTDA + logical,intent(in) :: evDyn + logical,intent(in) :: spin_conserved + logical,intent(in) :: spin_flip + double precision,intent(in) :: eta + logical,intent(in) :: regularize + + integer,intent(in) :: nNuc + double precision,intent(in) :: ZNuc(nNuc) + double precision,intent(in) :: rNuc(nNuc,ncart) + + 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) :: S(nBas,nBas) + double precision,intent(in) :: T(nBas,nBas) + double precision,intent(in) :: V(nBas,nBas) + double precision,intent(in) :: Hc(nBas,nBas) + double precision,intent(in) :: X(nBas,nBas) + double precision,intent(in) :: ERI_AO(nBas,nBas,nBas,nBas) + double precision,intent(in) :: ERI_aaaa(nBas,nBas,nBas,nBas) + double precision,intent(in) :: ERI_aabb(nBas,nBas,nBas,nBas) + double precision,intent(in) :: ERI_bbbb(nBas,nBas,nBas,nBas) + double precision,intent(in) :: dipole_int_AO(nBas,nBas,ncart) + double precision,intent(in) :: dipole_int_aa(nBas,nBas,ncart) + double precision,intent(in) :: dipole_int_bb(nBas,nBas,ncart) + +! Local variables + integer :: nSCF + integer :: nBasSq + double precision :: dipole(ncart) + integer :: n_diis + double precision :: rcond(nspin) + double precision,external :: trace_matrix + double precision :: Conv + double precision :: ET(nspin) + double precision :: EV(nspin) + double precision :: EJ(nsp) + double precision :: Ex(nspin) + double precision :: EqsGT + integer :: ispin,is + integer :: iblock + integer :: nH_sc,nH_sf,nHaa,nHab,nHbb + integer :: nP_sc,nP_sf,nPaa,nPab,nPbb + double precision :: EcRPA(nspin),Ecaa,Ecbb + double precision :: EcBSE(nspin) + double precision :: EcAC(nspin) + double precision :: EcGM(nspin) + double precision,allocatable :: Omega1ab(:),Omega1aa(:),Omega1bb(:) + double precision,allocatable :: X1ab(:,:),X1aa(:,:),X1bb(:,:) + double precision,allocatable :: Y1ab(:,:),Y1aa(:,:),Y1bb(:,:) + double precision,allocatable :: rho1ab(:,:,:),rho1aa(:,:,:),rho1bb(:,:,:) + double precision,allocatable :: Omega2ab(:),Omega2aa(:),Omega2bb(:) + double precision,allocatable :: X2ab(:,:),X2aa(:,:),X2bb(:,:) + double precision,allocatable :: Y2ab(:,:),Y2aa(:,:),Y2bb(:,:) + double precision,allocatable :: rho2ab(:,:,:),rho2aa(:,:,:),rho2bb(:,:,:) + double precision,allocatable :: c(:,:,:) + double precision,allocatable :: cp(:,:,:) + double precision,allocatable :: P(:,:,:) + double precision,allocatable :: F(:,:,:) + double precision,allocatable :: Fp(:,:,:) + double precision,allocatable :: J(:,:,:) + double precision,allocatable :: K(:,:,:) + double precision,allocatable :: SigT(:,:,:) + double precision,allocatable :: SigTp(:,:,:) + double precision,allocatable :: SigTm(:,:,:) + double precision,allocatable :: Z(:,:) + double precision,allocatable :: eGT(:,:) + double precision,allocatable :: eOld(:,:) + double precision,allocatable :: error_diis(:,:,:) + double precision,allocatable :: e_diis(:,:,:) + double precision,allocatable :: F_diis(:,:,:) + double precision,allocatable :: error(:,:,:) + + +! Hello world + + write(*,*) + write(*,*)'************************************************' + write(*,*)'| Self-consistent qsUGT calculation |' + write(*,*)'************************************************' + write(*,*) + +! Dimensions of the pp-URPA linear reponse matrices + + nPaa = nV(1)*(nV(1)-1)/2 + nPbb = nV(2)*(nV(2)-1)/2 + + nHaa = nO(1)*(nO(1)-1)/2; + nHbb = nO(2)*(nO(2)-1)/2; + + nPab = nV(1)*nV(2) + nHab = nO(1)*nO(2) + + nP_sc = nPab + nH_sc = nHab + + nP_sf = nPaa + nPbb + nH_sf = nHaa + nHbb + + nBasSq = nBas*nBas + + +! Memory allocation + + allocate(SigT(nBas,nbas,nspin),SigTp(nBas,nbas,nspin),SigTm(nBas,nbas,nspin), & + Z(nBas,nspin),eGT(nBas,nspin),eOld(nBas,nspin), & + error_diis(nBas,max_diis,nspin),e_diis(nBasSq,max_diis,nspin), & + F_diis(nBasSq,max_diis,nspin),error(nBas,nBas,nspin),& + c(nBas,nBas,nspin),cp(nBas,nBas,nspin),P(nBas,nBas,nspin),F(nBas,nBas,nspin), & + Fp(nBas,nBas,nspin),J(nBas,nBas,nspin),K(nBas,nBas,nspin)) + + allocate(Omega1ab(nPab),X1ab(nPab,nPab),Y1ab(nHab,nPab), & + Omega2ab(nHab),X2ab(nPab,nHab),Y2ab(nHab,nHab), & + rho1ab(nBas,nBas,nPab),rho2ab(nBas,nBas,nHab), & + Omega1aa(nPaa),X1aa(nPaa,nPaa),Y1aa(nHaa,nPaa), & + Omega2aa(nHaa),X2aa(nPaa,nHaa),Y2aa(nHaa,nHaa), & + rho1aa(nBas,nBas,nPaa),rho2aa(nBas,nBas,nHaa), & + Omega1bb(nPbb),X1bb(nPbb,nPbb),Y1bb(nHbb,nPbb), & + Omega2bb(nPbb),X2bb(nPbb,nPbb),Y2bb(nHbb,nPbb), & + rho1bb(nBas,nBas,nPbb),rho2bb(nBas,nBas,nHbb)) + +!Initialization + + nSCF = -1 + n_diis = 0 + Conv = 1d0 + P(:,:,:) = PHF(:,:,:) + e_diis(:,:,:) = 0d0 + error_diis(:,:,:) = 0d0 + eGT(:,:) = eHF(:,:) + eOld(:,:) = eHF(:,:) + c(:,:,:) = cHF(:,:,:) + Z(:,:) = 1d0 + rcond(:) = 0d0 + +!------------------------------------------------------------------------ +! Main loop +!------------------------------------------------------------------------ + + do while(Conv > thresh .and. nSCF <= maxSCF) + +! Increment + + nSCF = nSCF + 1 + +! Buid Coulomb matrix + do ispin=1,nspin + call Coulomb_matrix_AO_basis(nBas,P(:,:,ispin),ERI_AO(:,:,:,:), & + J(:,:,ispin)) + end do + +! Compute exchange part of the self-energy + do ispin=1,nspin + call exchange_matrix_AO_basis(nBas,P(:,:,ispin),ERI_AO(:,:,:,:), & + K(:,:,ispin)) + end do + +! AO to MO transformation of two-electron integrals + + ! 4-index transform for (aa|aa) block + + call AOtoMO_integral_transform(1,1,1,1,nBas,c,ERI_AO,ERI_aaaa) + + ! 4-index transform for (aa|bb) block + + call AOtoMO_integral_transform(1,1,2,2,nBas,c,ERI_AO,ERI_aabb) + + ! 4-index transform for (bb|bb) block + + call AOtoMO_integral_transform(2,2,2,2,nBas,c,ERI_AO,ERI_bbbb) + +!---------------------------------------------- +! alpha-beta block +!---------------------------------------------- + + ispin = 1 + iblock = 3 +! iblock = 1 + +! Compute linear response + + call unrestricted_linear_response_pp(iblock,TDA,nBas,nC,nO,nV,nR,nPaa,nPab,nPbb, & + nPab,nHaa,nHab,nHbb,nHab,1d0,eGT,ERI_aaaa, & + ERI_aabb,ERI_bbbb,Omega1ab,X1ab,Y1ab, & + Omega2ab,X2ab,Y2ab,EcRPA(ispin)) + +! EcRPA(ispin) = 1d0*EcRPA(ispin) + +!---------------------------------------------- +! alpha-alpha block +!---------------------------------------------- + + ispin = 2 + iblock = 4 + +! Compute linear response + + call unrestricted_linear_response_pp(iblock,TDA,nBas,nC,nO,nV,nR,nPaa,nPab,nPbb, & + nPaa,nHaa,nHab,nHbb,nHaa,1d0,eGT,ERI_aaaa, & + ERI_aabb,ERI_bbbb,Omega1aa,X1aa,Y1aa, & + Omega2aa,X2aa,Y2aa,EcRPA(ispin)) + +! EcRPA(ispin) = 2d0*EcRPA(ispin) +! EcRPA(ispin) = 3d0*EcRPA(ispin) + Ecaa = EcRPA(2) + +!---------------------------------------------- +! beta-beta block +!---------------------------------------------- + + ispin = 2 + iblock = 7 + +! Compute linear response + + call unrestricted_linear_response_pp(iblock,TDA,nBas,nC,nO,nV,nR,nPaa,nPab,nPbb, & + nPbb,nHaa,nHab,nHbb,nHbb,1d0,eGT,ERI_aaaa, & + ERI_aabb,ERI_bbbb,Omega1bb,X1bb,Y1bb, & + Omega2bb,X2bb,Y2bb,EcRPA(ispin)) + +! EcRPA(ispin) = 2d0*EcRPA(ispin) +! EcRPA(ispin) = 3d0*EcRPA(ispin) + Ecbb = EcRPA(2) + EcRPA(2) = Ecaa + Ecbb + EcRPA(1) = EcRPA(1) - EcRPA(2) + EcRPA(2) = 3d0*EcRPA(2) + +!---------------------------------------------- +! Compute T-matrix version of the self-energy +!---------------------------------------------- + + EcGM = 0d0 + SigT(:,:,:) = 0d0 + Z(:,:) = 0d0 + +!alpha-beta block + + iblock = 3 + + call unrestricted_excitation_density_Tmatrix(iblock,nBas,nC,nO,nV,nR,nHab,nPab, & + ERI_aaaa,ERI_aabb,ERI_bbbb,X1ab,Y1ab, & + rho1ab,X2ab,Y2ab,rho2ab) +!alpha-alpha block + + iblock = 4 + + call unrestricted_excitation_density_Tmatrix(iblock,nBas,nC,nO,nV,nR,nHaa,nPaa, & + ERI_aaaa,ERI_aabb,ERI_bbbb,X1aa,Y1aa, & + rho1aa,X2aa,Y2aa,rho2aa) + +!beta-beta block + + iblock = 7 + + call unrestricted_excitation_density_Tmatrix(iblock,nBas,nC,nO,nV,nR,nHbb,nPbb, & + ERI_aaaa,ERI_aabb,ERI_bbbb,X1bb,Y1bb, & + rho1bb,X2bb,Y2bb,rho2bb) + + call unrestricted_self_energy_Tmatrix(eta,nBas,nC,nO,nV,nR,nHaa,nHab,nHbb,nPaa,& + nPab,nPbb,eGT,Omega1aa,Omega1ab,Omega1bb,& + rho1aa,rho1ab,rho1bb,Omega2aa,Omega2ab,& + Omega2bb,rho2aa,rho2ab,rho2bb,EcGM,SigT) + + call unrestricted_renormalization_factor_Tmatrix(eta,nBas,nC,nO,nV,nR,nHaa,nHab,nHbb,& + nPaa,nPab,nPbb,eGT,Omega1aa,Omega1ab,& + Omega1bb,rho1aa,rho1ab,rho1bb, & + Omega2aa,Omega2ab,Omega2bb,rho2aa, & + rho2ab,rho2bb,Z) + + + Z(:,:) = 1d0/(1d0 - Z(:,:)) + +! Make correlation self-energy Hermitian and transform it back to AO basis + + do ispin=1,nspin + SigTp(:,:,ispin) = 0.5d0*(SigT(:,:,ispin) + transpose(SigT(:,:,ispin))) + SigTm(:,:,ispin) = 0.5d0*(SigT(:,:,ispin) - transpose(SigT(:,:,ispin))) + end do + + do ispin=1,nspin + call MOtoAO_transform(nBas,S,c(:,:,ispin),SigTp(:,:,ispin)) + end do + +! Solve the quasi-particle equation + + do ispin=1,nspin + F(:,:,ispin) = Hc(:,:) + J(:,:,ispin) + J(:,:,mod(ispin,2)+1) + K(:,:,ispin) & + + SigT(:,:,ispin) + end do + +! Compute commutator and convergence criteria + + do ispin=1,nspin + error_diis(:,:,ispin) = matmul(F(:,:,ispin),matmul(P(:,:,ispin),S(:,:))) & + - matmul(matmul(S(:,:),P(:,:,ispin)),F(:,:,ispin)) + end do + +! DIIS extrapolation + + n_diis = min(n_diis+1,max_diis) + if(minval(rcond(:)) > 1d-7) then + do ispin=1,nspin + if(nO(ispin) > 1) call DIIS_extrapolation(rcond(ispin),nBasSq,nBasSq,n_diis, & + error_diis(:,1:n_diis,ispin), & + F_diis(:,1:n_diis,ispin),& + error_diis(:,:,ispin),F(:,:,ispin)) + end do + else + n_diis = 0 + end if + +! Transform Fock matrix in orthogonal basis + + do ispin=1,nspin + Fp(:,:,ispin) = matmul(transpose(X(:,:)),matmul(F(:,:,ispin),X(:,:))) + end do + +! Diagonalize Fock matrix to get eigenvectors and eigenvalues + + cp(:,:,:) = Fp(:,:,:) + do ispin=1,nspin + call diagonalize_matrix(nBas,cp(:,:,ispin),eGT(:,ispin)) + end do + +! Back-transform eigenvectors in non-orthogonal basis + + do ispin=1,nspin + c(:,:,ispin) = matmul(X(:,:),cp(:,:,ispin)) + end do + +! Back-transform self-energy + + do ispin=1,nspin + SigTp(:,:,ispin) = matmul(transpose(c(:,:,ispin)),matmul(SigTp(:,:,ispin),c(:,:,ispin))) + end do + +! Compute density matrix + + do ispin=1,nspin + P(:,:,ispin) = matmul(c(:,1:nO(ispin),ispin),transpose(c(:,1:nO(ispin),ispin))) + end do + +! Save quasiparticles energy for next cycle + + Conv = maxval(abs(eGT(:,:) - eOld(:,:))) + eOld(:,:) = eGT(:,:) + +!------------------------------------------------------------------------ +! Compute total energy +!------------------------------------------------------------------------ + + ! Kinetic energy + + do ispin=1,nspin + ET(ispin) = trace_matrix(nBas,matmul(P(:,:,ispin),T(:,:))) + end do + + ! Potential energy + + do ispin=1,nspin + EV(ispin) = trace_matrix(nBas,matmul(P(:,:,ispin),V(:,:))) + end do + +! Coulomb energy + + EJ(1) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,1),J(:,:,1))) + EJ(2) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,1),J(:,:,2))) & + + 0.5d0*trace_matrix(nBas,matmul(P(:,:,2),J(:,:,1))) + EJ(3) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,2),J(:,:,2))) + + ! Exchange energy + + do ispin=1,nspin + Ex(ispin) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,ispin),K(:,:,ispin))) + end do + + ! Total energy + + EqsGT = sum(ET(:)) + sum(EV(:)) + sum(EJ(:)) + sum(Ex(:)) + +! Print results + + call dipole_moment(nBas,P,nNuc,ZNuc,rNuc,dipole_int_AO,dipole) + call print_qsUGT(nBas,nO,nSCF,Conv,thresh,eHF,eGT,c,SigTp,Z,ENuc,ET,EV,EJ,Ex,EcGM,EcRPA,EqsGT,dipole) + + enddo +!------------------------------------------------------------------------ +! End main loop +!------------------------------------------------------------------------ + +! Did it actually converge? + + if(nSCF == maxSCF+1) then + + write(*,*) + write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!' + write(*,*)' Convergence failed ' + write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!' + write(*,*) + + stop + + endif + +! Free memory + + deallocate(Omega1ab,X1ab,Y1ab,Omega2ab,X2ab,Y2ab,rho1ab,rho2ab, & + Omega1aa,X1aa,Y1aa,Omega2aa,X2aa,Y2aa,rho1aa,rho2aa, & + Omega1bb,X1bb,Y1bb,Omega2bb,X2bb,Y2bb,rho1bb,rho2bb) + + deallocate(c,cp,P,F,Fp,J,K,SigT,SigTp,SigTm,Z,error,error_diis,F_diis) + +end subroutine qsUGT diff --git a/src/GT/unrestricted_self_energy_Tmatrix.f90 b/src/GT/unrestricted_self_energy_Tmatrix.f90 new file mode 100644 index 0000000..3f386dd --- /dev/null +++ b/src/GT/unrestricted_self_energy_Tmatrix.f90 @@ -0,0 +1,207 @@ +subroutine unrestricted_self_energy_Tmatrix(eta,nBas,nC,nO,nV,nR,nHaa,nHab,nHbb,nPaa,& + nPab,nPbb,e,Omega1aa,Omega1ab,Omega1bb,& + rho1aa,rho1ab,rho1bb,Omega2aa,Omega2ab,& + Omega2bb,rho2aa,rho2ab,rho2bb,EcGM,SigT) + +! Compute the correlation part of the T-matrix 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) :: nHaa,nHab,nHbb + integer,intent(in) :: nPaa,nPab,nPbb + double precision,intent(in) :: e(nBas,nspin) + double precision,intent(in) :: Omega1aa(nPaa),Omega1ab(nPab),Omega1bb(nPbb) + double precision,intent(in) :: rho1aa(nBas,nBas,nPaa),rho1ab(nBas,nBas,nPab) + double precision,intent(in) :: rho1bb(nBas,nBas,nPbb) + double precision,intent(in) :: Omega2aa(nHaa),Omega2ab(nHab),Omega2bb(nHbb) + double precision,intent(in) :: rho2aa(nBas,nBas,nHaa),rho2ab(nBas,nBas,nHab) + double precision,intent(in) :: rho2bb(nBas,nBas,nHbb) + +! Local variables + + integer :: i,j,a,b,p,q,cd,kl + double precision :: eps + +! Output variables + + double precision,intent(inout) :: EcGM(nspin) + double precision,intent(inout) :: SigT(nBas,nBas,nspin) + +!---------------------------------------------- +! Occupied part of the T-matrix self-energy +!---------------------------------------------- + +!spin up part + + do p=nC(1)+1,nBas-nR(1) + do q=nC(1)+1,nBas-nR(1) + do i=nC(1)+1,nO(1) + do cd=1,nPaa + eps = e(p,1) + e(i,1) - Omega1aa(cd) + SigT(p,q,1) = SigT(p,q,1) + rho1aa(p,i,cd)*rho1aa(q,i,cd)*eps/(eps**2 + eta**2) + enddo + enddo + + do i=nC(2)+1,nO(2) + do cd=1,nPab + eps = e(p,1) + e(i,1) - Omega1ab(cd) + SigT(p,q,1) = SigT(p,q,1) + rho1ab(p,i,cd)*rho1ab(q,i,cd)*eps/(eps**2 + eta**2) + end do + end do + enddo + enddo + +!spin down part + + do p=nC(2)+1,nBas-nR(2) + do q=nC(2)+1,nBas-nR(2) + do i=nC(2)+1,nO(2) + do cd=1,nPbb + eps = e(p,2) + e(i,2) - Omega1bb(cd) + SigT(p,q,2) = SigT(p,q,2) + rho1bb(p,i,cd)*rho1bb(q,i,cd)*eps/(eps**2 + eta**2) + enddo + enddo + + do i=nC(2)+1,nO(2) + do cd=1,nPab + eps = e(p,2) + e(i,2) - Omega1ab(cd) + SigT(p,q,2) = SigT(p,q,2) + rho1ab(p,i,cd)*rho1ab(q,i,cd)*eps/(eps**2 + eta**2) + end do + end do + enddo + enddo + +!---------------------------------------------- +! Virtual part of the T-matrix self-energy +!---------------------------------------------- + +! spin up part + + do p=nC(1)+1,nBas-nR(1) + do q=nC(1)+1,nBas-nR(1) + do a=nO(1)+1,nBas-nR(1) + do kl=1,nHaa + eps = e(p,1) + e(a,1) - Omega2aa(kl) + SigT(p,q,1) = SigT(p,q,1) + rho2aa(p,a,kl)*rho2aa(q,a,kl)*eps/(eps**2 + eta**2) + enddo + end do + + do a=nO(1)+1,nBas-nR(1) + do kl=1,nHab + eps = e(p,1) + e(a,1) - Omega2ab(kl) + SigT(p,q,1) = SigT(p,q,1) + rho2ab(p,a,kl)*rho2ab(q,a,kl)*eps/(eps**2 + eta**2) + end do + end do + enddo + enddo + +!spin down part + + do p=nC(2)+1,nBas-nR(2) + do q=nC(2)+1,nBas-nR(2) + do a=nO(2)+1,nBas-nR(2) + do kl=1,nHbb + eps = e(p,2) + e(a,2) - Omega2bb(kl) + SigT(p,q,2) = SigT(p,q,2) + rho2bb(p,a,kl)*rho2bb(q,a,kl)*eps/(eps**2 + eta**2) + enddo + enddo + + do a=nO(2)+1,nBas-nR(2) + do kl=1,nHab + eps = e(p,2) + e(a,2) - Omega2ab(kl) + SigT(p,q,2) = SigT(p,q,2) + rho2ab(p,a,kl)*rho2ab(q,a,kl)*eps/(eps**2 + eta**2) + end do + end do + enddo + enddo + +!---------------------------------------------- +! Galitskii-Migdal correlation energy +!---------------------------------------------- + +!spin up part + + do i=nC(1)+1,nO(1) + do j=nC(1)+1,nO(1) + do cd=1,nPaa + eps = e(i,1) + e(j,1) - Omega1aa(cd) + EcGM(1) = EcGM(1) + rho1aa(i,j,cd)*rho1aa(i,j,cd)*eps/(eps**2 + eta**2) + enddo + enddo + enddo + + do i=nC(1)+1,nO(1) + do j=nC(2)+1,nO(2) + do cd=1,nPab + eps = e(i,1) + e(j,1) - Omega1ab(cd) + EcGM(1) = EcGM(1) + rho1ab(i,j,cd)*rho1ab(i,j,cd)*eps/(eps**2 + eta**2) + end do + end do + end do + + do a=nO(1)+1,nBas-nR(1) + do b=nO(1)+1,nBas-nR(1) + do kl=1,nHaa + eps = e(a,1) + e(b,1) - Omega2aa(kl) + EcGM(1) = EcGM(1) - rho2aa(a,b,kl)*rho2aa(a,b,kl)*eps/(eps**2 + eta**2) + enddo + enddo + enddo + + do a=nO(1)+1,nBas-nR(1) + do b=nO(1)+1,nBas-nR(1) + do kl=1,nHab + eps = e(a,1) + e(b,1) - Omega2ab(kl) + EcGM(1) = EcGM(1) - rho2ab(a,b,kl)*rho2ab(a,b,kl)*eps/(eps**2 + eta**2) + enddo + enddo + enddo + +! spin down part + + do i=nC(2)+1,nO(2) + do j=nC(2)+1,nO(2) + do cd=1,nPbb + eps = e(i,2) + e(j,2) - Omega1bb(cd) + EcGM(2) = EcGM(2) + rho1bb(i,j,cd)*rho1bb(i,j,cd)*eps/(eps**2 + eta**2) + enddo + enddo + enddo + + do i=nC(1)+1,nO(1) + do j=nC(2)+1,nO(2) + do cd=1,nPab + eps = e(i,2) + e(j,2) - Omega1ab(cd) + EcGM(2) = EcGM(2) + rho1ab(i,j,cd)*rho1ab(i,j,cd)*eps/(eps**2 + eta**2) + end do + end do + end do + + do a=nO(1)+1,nBas-nR(1) + do b=nO(2)+1,nBas-nR(2) + do kl=1,nHab + eps = e(a,2) + e(b,2) - Omega2ab(kl) + EcGM(2) = EcGM(2) - rho2ab(a,b,kl)*rho2ab(a,b,kl)*eps/(eps**2 + eta**2) + enddo + enddo + enddo + + do a=nO(2)+1,nBas-nR(2) + do b=nO(2)+1,nBas-nR(2) + do kl=1,nHbb + eps = e(a,2) + e(b,2) - Omega2bb(kl) + EcGM(2) = EcGM(2) - rho2bb(a,b,kl)*rho2bb(a,b,kl)*eps/(eps**2 + eta**2) + enddo + enddo + enddo + +end subroutine unrestricted_self_energy_Tmatrix diff --git a/src/QuAcK/QuAcK.f90 b/src/QuAcK/QuAcK.f90 index ff07b53..dd72207 100644 --- a/src/QuAcK/QuAcK.f90 +++ b/src/QuAcK/QuAcK.f90 @@ -1222,8 +1222,11 @@ program QuAcK if(unrestricted) then - print*,'!!! qsGT NYI at the unrestricted level !!!' - + !print*,'!!! qsGT NYI at the unrestricted level !!!' + call qsUGT(maxSCF_GT,thresh_GT,n_diis_GT,doACFDT,exchange_kernel,doXBS,BSE,TDA_T, & + TDA,dBSE,dTDA,evDyn,spin_conserved,spin_flip,eta_GT,regGT,nBas,nC,nO,nV,& + nR,nS,nNuc,ZNuc,rNuc,ENuc,EUHF,S,X,T,V,Hc,ERI_AO,ERI_MO_aaaa,ERI_MO_aabb,& + ERI_MO_bbbb,dipole_int_AO,dipole_int_aa,dipole_int_bb,PHF,cHF,eHF) else call qsGT(maxSCF_GT,thresh_GT,n_diis_GT,doACFDT,exchange_kernel,doXBS, &