10
1
mirror of https://github.com/pfloos/quack synced 2025-05-07 07:35:02 +02:00

start complex matrix diag

This commit is contained in:
Loris Burth 2025-03-04 14:04:03 +01:00
parent 685b0c5824
commit 6069b68bdf
6 changed files with 53 additions and 10 deletions

View File

@ -3,6 +3,6 @@
cp ./methods.test ../input/methods
cp ./options.test ../input/options
cd ..
python3 PyDuck.py -x N2 -b aug-cc-pVTZ -c -1 -m 2
python3 PyDuck.py -x N2 -b sto-3g -c -1 -m 2
cp input/methods.default input/methods
cp input/options.default input/options

View File

@ -14,7 +14,7 @@ subroutine complex_Hartree_matrix_AO_basis(nBas,P,ERI,H)
! Local variables
integer :: mu,nu,la,si
! Output variables
complex*16,intent(out) :: H(nBas,nBas)

View File

@ -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),transpose(c(:,1:nO)))
P(:,:) = 2d0*matmul(c(:,1:nO),conjg(transpose(c(:,1:nO))))
! Initialization
n_diis = 0
@ -135,7 +135,6 @@ subroutine cRHF(dotest,maxSCF,thresh,max_diis,guess_type,level_shift,nNuc,ZNuc,r
! Kinetic energy
ET = trace_matrix(nBas,matmul(P,T))
! Potential energy
EV = trace_matrix(nBas,matmul(P,V))
@ -167,13 +166,13 @@ subroutine cRHF(dotest,maxSCF,thresh,max_diis,guess_type,level_shift,nNuc,ZNuc,r
! Diagonalize Fock matrix
Fp = matmul(transpose(X),matmul(F,X))
Fp = matmul(conjg(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),transpose(c(:,1:nO)))
P(:,:) = 2d0*matmul(c(:,1:nO),conjg(transpose(c(:,1:nO))))
! Dump results

View File

@ -7,7 +7,7 @@ subroutine complex_core_guess(nBas, nOrb, Hc, X, c)
! Input variables
integer,intent(in) :: nBas, nOrb
double precision,intent(in) :: Hc(nBas,nBas)
complex*16,intent(in) :: Hc(nBas,nBas)
double precision,intent(in) :: X(nBas,nOrb)
! Local variables
@ -28,9 +28,9 @@ subroutine complex_core_guess(nBas, nOrb, Hc, X, c)
cp(:,:) = matmul(transpose(X(:,:)), matmul(Hc(:,:), X(:,:)))
call diagonalize_matrix(nOrb, cp, e)
call diagonalize_complex_matrix(nOrb, cp, e)
c(:,:) = matmul(X(:,:), cp(:,:))
write(*,*) c
deallocate(cp, e)
end subroutine

View File

@ -28,7 +28,7 @@ subroutine complex_mo_guess(nBas, nOrb, guess_type, S, Hc, X, c)
write(*,*) 'Core guess...'
call complex_core_guess(nBas, nOrb, Hc, X, c)
else
print*,'Wrong guess option'

View File

@ -128,6 +128,50 @@ subroutine diagonalize_matrix(N,A,e)
end subroutine
subroutine diagonalize_complex_matrix(N,A,e)
! Diagonalize a general complex matrix
implicit none
! Input variables
integer,intent(in) :: N
complex*16,intent(inout) :: A(N,N)
complex*16,intent(out) :: e(N)
! Local variables
integer :: lwork,info
integer :: sdim
double precision :: rwork(N)
logical :: bwork(N)
double precision,allocatable :: work(:)
complex*16 :: VS(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)
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)
deallocate(work)
A = VS
if(info /= 0) then
print*,'Problem in diagonalize_matrix (zgees)!!'
end if
end subroutine
subroutine svd(N,A,U,D,Vt)
! Compute A = U.D.Vt