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

57 lines
1.4 KiB
Fortran
Raw Normal View History

2020-01-17 17:35:40 +01:00
subroutine mo_guess(nBas,nO,guess_type,S,Hc,ERI,J,K,X,cp,F,Fp,e,c,P)
2019-10-05 22:06:25 +02:00
! Guess of the molecular orbitals for HF calculation
implicit none
! Input variables
integer,intent(in) :: nBas
integer,intent(in) :: nO
integer,intent(in) :: guess_type
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)
double precision,intent(inout):: P(nBas,nBas)
! Local variables
integer :: nSCF
! Output variables
double precision,intent(out) :: c(nBas,nBas)
if(guess_type == 1) then
Fp = matmul(transpose(X),matmul(Hc,X))
cp(:,:) = Fp(:,:)
call diagonalize_matrix(nBas,cp,e)
c = matmul(X,cp)
elseif(guess_type == 2) then
2020-01-17 17:35:40 +01:00
call huckel_guess(nBas,nO,S,Hc,ERI,J,K,X,cp,F,Fp,e,c,P)
2019-10-05 22:06:25 +02:00
elseif(guess_type == 3) then
call random_number(c)
2020-01-08 10:17:19 +01:00
else
print*,'Wrong guess option'
stop
2019-10-05 22:06:25 +02:00
endif
P(:,:) = 2d0*matmul(c(:,1:nO),transpose(c(:,1:nO)))
end subroutine