4
1
mirror of https://github.com/pfloos/quack synced 2024-06-23 13:42:19 +02:00
quack/src/HF/huckel_guess.f90

59 lines
1.5 KiB
Fortran
Raw Normal View History

2020-01-17 17:35:40 +01:00
subroutine huckel_guess(nBas,nO,S,Hc,ERI,J,K,X,cp,F,Fp,e,c,P)
2019-10-05 22:06:25 +02:00
! Hickel guess of the molecular orbitals for HF calculation
implicit none
! Input variables
integer,intent(in) :: nBas
integer,intent(in) :: nO
double precision,intent(in) :: S(nBas,nBas)
double precision,intent(in) :: Hc(nBas,nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
2020-01-17 17:35:40 +01:00
double precision,intent(inout):: J(nBas,nBas)
double precision,intent(inout):: K(nBas,nBas)
2019-10-05 22:06:25 +02:00
double precision,intent(in) :: X(nBas,nBas)
double precision,intent(inout):: cp(nBas,nBas)
2020-01-15 22:29:43 +01:00
double precision,intent(inout):: F(nBas,nBas)
2019-10-05 22:06:25 +02:00
double precision,intent(inout):: Fp(nBas,nBas)
double precision,intent(inout):: e(nBas)
2020-01-17 17:35:40 +01:00
double precision,intent(inout):: P(nBas,nBas)
2019-10-05 22:06:25 +02:00
! Local variables
integer :: mu,nu
double precision :: a
! Output variables
double precision,intent(out) :: c(nBas,nBas)
a = 1.75d0
2020-01-17 17:35:40 +01:00
Fp(:,:) = matmul(transpose(X(:,:)),matmul(Hc(:,:),X(:,:)))
cp(:,:) = Fp(:,:)
call diagonalize_matrix(nBas,cp,e)
c(:,:) = matmul(X(:,:),cp(:,:))
call Coulomb_matrix_AO_basis(nBas,P,ERI,J)
call exchange_matrix_AO_basis(nBas,P,ERI,K)
F(:,:) = Hc(:,:) + J(:,:) + 0.5d0*K(:,:)
2020-01-15 09:22:54 +01:00
do mu=1,nBas
do nu=mu+1,nBas
2020-01-15 22:29:43 +01:00
F(mu,nu) = 0.5d0*a*S(mu,nu)*(Hc(mu,mu) + Hc(nu,nu))
F(nu,mu) = F(mu,nu)
2020-01-15 09:22:54 +01:00
enddo
enddo
2020-01-15 22:29:43 +01:00
Fp(:,:) = matmul(transpose(X(:,:)),matmul(F(:,:),X(:,:)))
2019-10-05 22:06:25 +02:00
cp(:,:) = Fp(:,:)
call diagonalize_matrix(nBas,cp,e)
2020-01-15 09:22:54 +01:00
c(:,:) = matmul(X(:,:),cp(:,:))
2019-10-05 22:06:25 +02:00
end subroutine