mirror of
https://github.com/pfloos/quack
synced 2025-05-06 23:34:42 +02:00
complex diag implemented/wrapped
This commit is contained in:
parent
6069b68bdf
commit
587f897e81
@ -25,7 +25,7 @@ subroutine complex_Hartree_matrix_AO_basis(nBas,P,ERI,H)
|
||||
do nu=1,nBas
|
||||
do la=1,nBas
|
||||
do mu=1,nBas
|
||||
H(mu,nu) = H(mu,nu) + P(la,si)*ERI(mu,la,nu,si)
|
||||
H(mu,nu) = H(mu,nu) + P(la,si)*ERI(mu,la,nu,si)
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
|
@ -93,7 +93,7 @@ subroutine cRHF(dotest,maxSCF,thresh,max_diis,guess_type,level_shift,nNuc,ZNuc,r
|
||||
! Guess coefficients and density matrix
|
||||
|
||||
call complex_mo_guess(nBas,nBas,guess_type,S,Hc,X,c)
|
||||
P(:,:) = 2d0*matmul(c(:,1:nO),conjg(transpose(c(:,1:nO))))
|
||||
P(:,:) = 2d0*matmul(c(:,1:nO),transpose(c(:,1:nO)))
|
||||
! Initialization
|
||||
|
||||
n_diis = 0
|
||||
@ -166,13 +166,13 @@ subroutine cRHF(dotest,maxSCF,thresh,max_diis,guess_type,level_shift,nNuc,ZNuc,r
|
||||
|
||||
! Diagonalize Fock matrix
|
||||
|
||||
Fp = matmul(conjg(transpose(X)),matmul(F,X))
|
||||
Fp = matmul(transpose(X),matmul(F,X))
|
||||
cp(:,:) = Fp(:,:)
|
||||
call diagonalize_matrix(nBas,cp,eHF)
|
||||
c = matmul(X,cp)
|
||||
! Density matrix
|
||||
|
||||
P(:,:) = 2d0*matmul(c(:,1:nO),conjg(transpose(c(:,1:nO))))
|
||||
P(:,:) = 2d0*matmul(c(:,1:nO),transpose(c(:,1:nO)))
|
||||
|
||||
! Dump results
|
||||
|
||||
|
@ -28,9 +28,8 @@ subroutine complex_core_guess(nBas, nOrb, Hc, X, c)
|
||||
|
||||
cp(:,:) = matmul(transpose(X(:,:)), matmul(Hc(:,:), X(:,:)))
|
||||
|
||||
call diagonalize_complex_matrix(nOrb, cp, e)
|
||||
call complex_diagonalize_matrix(nOrb, cp, e)
|
||||
c(:,:) = matmul(X(:,:), cp(:,:))
|
||||
write(*,*) c
|
||||
deallocate(cp, e)
|
||||
|
||||
end subroutine
|
||||
|
@ -128,7 +128,7 @@ subroutine diagonalize_matrix(N,A,e)
|
||||
|
||||
end subroutine
|
||||
|
||||
subroutine diagonalize_complex_matrix(N,A,e)
|
||||
subroutine complex_diagonalize_matrix(N,A,e)
|
||||
|
||||
! Diagonalize a general complex matrix
|
||||
|
||||
@ -143,28 +143,27 @@ subroutine diagonalize_complex_matrix(N,A,e)
|
||||
! Local variables
|
||||
|
||||
integer :: lwork,info
|
||||
integer :: sdim
|
||||
double precision :: rwork(N)
|
||||
logical :: bwork(N)
|
||||
double precision,allocatable :: work(:)
|
||||
complex*16 :: VS(N,N)
|
||||
double precision :: rwork(2*N)
|
||||
complex*16,allocatable :: work(:)
|
||||
complex*16 :: VL(N,N)
|
||||
complex*16 :: VR(N,N)
|
||||
|
||||
! Memory allocation
|
||||
|
||||
allocate(work(1))
|
||||
|
||||
lwork = -1
|
||||
call zgees('V','N',N,A,N,sdim,e,VS,N,work,lwork,rwork,bwork,info)
|
||||
call zgeev('N','V',N,A,N,e,VL,N,VR,N,work,lwork,rwork,info)
|
||||
lwork = int(work(1))
|
||||
|
||||
deallocate(work)
|
||||
|
||||
allocate(work(lwork))
|
||||
|
||||
call zgees('V','N',N,A,N,sdim,e,VS,N,work,lwork,rwork,bwork,info)
|
||||
call zgeev('N','V',N,A,N,e,VL,N,VR,N,work,lwork,rwork,info)
|
||||
|
||||
deallocate(work)
|
||||
A = VS
|
||||
A = VR
|
||||
|
||||
if(info /= 0) then
|
||||
print*,'Problem in diagonalize_matrix (zgees)!!'
|
||||
|
Loading…
x
Reference in New Issue
Block a user