4
1
mirror of https://github.com/pfloos/quack synced 2024-11-09 07:33:55 +01:00
quack/src/HF/UHF.f90

260 lines
7.2 KiB
Fortran
Raw Normal View History

2021-02-15 17:27:06 +01:00
subroutine UHF(maxSCF,thresh,max_diis,guess_type,mix,nNuc,ZNuc,rNuc,ENuc,nBas,nO,S,T,V,Hc,ERI,dipole_int,X,EUHF,e,c,P,Vx)
2019-03-19 10:13:33 +01:00
! Perform unrestricted Hartree-Fock calculation
implicit none
include 'parameters.h'
! Input variables
integer,intent(in) :: maxSCF
integer,intent(in) :: max_diis
integer,intent(in) :: guess_type
2020-10-27 17:23:48 +01:00
logical,intent(in) :: mix
2019-03-19 10:13:33 +01:00
double precision,intent(in) :: thresh
integer,intent(in) :: nBas
2020-10-04 12:40:44 +02:00
integer,intent(in) :: nNuc
double precision,intent(in) :: ZNuc(nNuc)
double precision,intent(in) :: rNuc(nNuc,ncart)
double precision,intent(in) :: ENuc
2019-03-19 10:13:33 +01:00
integer,intent(in) :: nO(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(nBas,nBas,nBas,nBas)
2020-10-04 12:40:44 +02:00
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
2019-03-19 10:13:33 +01:00
! Local variables
integer :: nSCF
integer :: nBasSq
integer :: n_diis
double precision :: conv
double precision :: rcond(nspin)
double precision :: ET(nspin)
double precision :: EV(nspin)
double precision :: EJ(nsp)
double precision :: Ex(nspin)
2020-10-04 12:40:44 +02:00
double precision :: dipole(ncart)
2019-03-19 10:13:33 +01:00
double precision,allocatable :: cp(:,:,:)
double precision,allocatable :: J(:,:,:)
double precision,allocatable :: F(:,:,:)
double precision,allocatable :: Fp(:,:,:)
2019-03-19 11:21:34 +01:00
double precision,allocatable :: K(:,:,:)
2019-03-19 10:13:33 +01:00
double precision,allocatable :: err(:,:,:)
double precision,allocatable :: err_diis(:,:,:)
double precision,allocatable :: F_diis(:,:,:)
double precision,external :: trace_matrix
integer :: ispin
2019-03-19 11:21:34 +01:00
! 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)
2021-02-15 17:27:06 +01:00
double precision,intent(out) :: Vx(nBas,nspin)
2019-03-19 11:21:34 +01:00
2019-03-19 10:13:33 +01:00
! Hello world
write(*,*)
write(*,*)'************************************************'
write(*,*)'* Unrestricted Hartree-Fock calculation *'
write(*,*)'************************************************'
write(*,*)
! Useful stuff
nBasSq = nBas*nBas
! Memory allocation
2019-03-19 11:21:34 +01:00
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), &
2019-03-19 10:13:33 +01:00
err_diis(nBasSq,max_diis,nspin),F_diis(nBasSq,max_diis,nspin))
! Guess coefficients and eigenvalues
if(guess_type == 1) then
2020-10-06 15:24:24 +02:00
F(:,:,:) = 0d0
2019-03-19 10:13:33 +01:00
do ispin=1,nspin
F(:,:,ispin) = Hc(:,:)
end do
else if(guess_type == 2) then
do ispin=1,nspin
call random_number(F(:,:,ispin))
end do
end if
! Initialization
nSCF = 0
conv = 1d0
n_diis = 0
F_diis(:,:,:) = 0d0
err_diis(:,:,:) = 0d0
!------------------------------------------------------------------------
! Main SCF loop
!------------------------------------------------------------------------
write(*,*)
2019-03-19 11:21:34 +01:00
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(*,*)'----------------------------------------------------------'
2019-03-19 10:13:33 +01:00
do while(conv > thresh .and. nSCF < maxSCF)
! Increment
nSCF = nSCF + 1
! 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
2019-03-19 11:21:34 +01:00
call diagonalize_matrix(nBas,cp(:,:,ispin),e(:,ispin))
2019-03-19 10:13:33 +01:00
end do
! Back-transform eigenvectors in non-orthogonal basis
do ispin=1,nspin
c(:,:,ispin) = matmul(X(:,:),cp(:,:,ispin))
end do
2020-10-27 17:23:48 +01:00
! Mix guess for UHF solution in singlet states
if(mix .and. nSCF == 1) call mix_guess(nBas,nO,c)
2019-03-19 10:13:33 +01:00
! Compute density matrix
do ispin=1,nspin
P(:,:,ispin) = matmul(c(:,1:nO(ispin),ispin),transpose(c(:,1:nO(ispin),ispin)))
end do
! Build Coulomb repulsion
do ispin=1,nspin
call Coulomb_matrix_AO_basis(nBas,P(:,:,ispin),ERI(:,:,:,:),J(:,:,ispin))
end do
! Compute exchange potential
do ispin=1,nspin
2019-03-19 11:21:34 +01:00
call exchange_matrix_AO_basis(nBas,P(:,:,ispin),ERI(:,:,:,:),K(:,:,ispin))
2019-03-19 10:13:33 +01:00
end do
2020-10-06 15:24:24 +02:00
2019-03-19 10:13:33 +01:00
! Build Fock operator
2020-10-06 15:24:24 +02:00
2019-03-19 10:13:33 +01:00
do ispin=1,nspin
2019-03-19 11:21:34 +01:00
F(:,:,ispin) = Hc(:,:) + J(:,:,ispin) + J(:,:,mod(ispin,2)+1) + K(:,:,ispin)
2019-03-19 10:13:33 +01:00
end do
! Check convergence
do ispin=1,nspin
err(:,:,ispin) = matmul(F(:,:,ispin),matmul(P(:,:,ispin),S(:,:))) - matmul(matmul(S(:,:),P(:,:,ispin)),F(:,:,ispin))
end do
2020-10-06 15:24:24 +02:00
if(nSCF > 1) conv = maxval(abs(err(:,:,:)))
2019-03-19 10:13:33 +01:00
! DIIS extrapolation
n_diis = min(n_diis+1,max_diis)
2021-03-27 15:03:54 +01:00
if(minval(rcond(:)) > 1d-7) then
do ispin=1,nspin
if(nO(ispin) > 1) call DIIS_extrapolation(rcond(ispin),nBasSq,nBasSq,n_diis,err_diis(:,1:n_diis,ispin), &
2020-10-21 12:09:18 +02:00
F_diis(:,1:n_diis,ispin),err(:,:,ispin),F(:,:,ispin))
2021-03-27 15:03:54 +01:00
end do
else
n_diis = 0
end if
2019-03-19 10:13:33 +01:00
!------------------------------------------------------------------------
! Compute UHF 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) = trace_matrix(nBas,matmul(P(:,:,1),J(:,:,2)))
EJ(3) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,2),J(:,:,2)))
! Exchange energy
do ispin=1,nspin
2019-03-19 11:21:34 +01:00
Ex(ispin) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,ispin),K(:,:,ispin)))
2019-03-19 10:13:33 +01:00
end do
! Total energy
2019-03-19 11:21:34 +01:00
EUHF = sum(ET(:)) + sum(EV(:)) + sum(EJ(:)) + sum(Ex(:))
2019-03-19 10:13:33 +01:00
! Dump results
2019-03-19 11:21:34 +01:00
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,'|'
2019-03-19 10:13:33 +01:00
end do
2019-03-19 11:21:34 +01:00
write(*,*)'----------------------------------------------------------'
2019-03-19 10:13:33 +01:00
!------------------------------------------------------------------------
! End of SCF loop
!------------------------------------------------------------------------
! Did it actually converge?
if(nSCF == maxSCF) then
write(*,*)
write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!'
write(*,*)' Convergence failed '
write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!'
write(*,*)
stop
end if
! Compute final UHF energy
2020-10-04 12:40:44 +02:00
call dipole_moment(nBas,P(:,:,1)+P(:,:,2),nNuc,ZNuc,rNuc,dipole_int,dipole)
call print_UHF(nBas,nO,S,e,c,ENuc,ET,EV,EJ,Ex,EUHF,dipole)
2019-03-19 10:13:33 +01:00
2021-02-15 17:27:06 +01:00
! Compute Vx for post-HF calculations
do ispin=1,nspin
call exchange_potential(nBas,c(:,:,ispin),K(:,:,ispin),Vx(:,ispin))
end do
2019-03-19 10:13:33 +01:00
end subroutine UHF