mirror of
https://github.com/pfloos/quack
synced 2024-11-03 20:53:53 +01:00
OK with ROHF
This commit is contained in:
parent
4e057b71f5
commit
66913f9a37
@ -1,6 +1,16 @@
|
||||
subroutine ROHF_fock_matrix(nBas,nOa,nOb,S,c,Fa,Fb,F)
|
||||
|
||||
! Construct the ROHF Fock matrix
|
||||
! Construct the ROHF Fock matrix in the AO basis
|
||||
! For open shells, the ROHF Fock matrix in the MO basis reads
|
||||
!
|
||||
! | F-K | F + K/2 | F |
|
||||
! |---------------------------------|
|
||||
! | F + K/2 | F | F - K/2 |
|
||||
! |---------------------------------|
|
||||
! | F | F - K/2 | F + K |
|
||||
!
|
||||
! with F = 1/2 (Fa + Fb) and K = Fb - Fa
|
||||
!
|
||||
|
||||
implicit none
|
||||
include 'parameters.h'
|
||||
@ -32,14 +42,14 @@ subroutine ROHF_fock_matrix(nBas,nOa,nOb,S,c,Fa,Fb,F)
|
||||
|
||||
! Roothan canonicalization parameters
|
||||
|
||||
aC = -0.5d0
|
||||
bC = +1.5d0
|
||||
aC = +1.5d0
|
||||
bC = -0.5d0
|
||||
|
||||
aO = +0.5d0
|
||||
bO = +0.5d0
|
||||
|
||||
aV = +1.5d0
|
||||
bV = -0.5d0
|
||||
aV = -0.5d0
|
||||
bV = +1.5d0
|
||||
|
||||
! Number of closed, open, and virtual orbitals
|
||||
|
||||
@ -49,8 +59,8 @@ subroutine ROHF_fock_matrix(nBas,nOa,nOb,S,c,Fa,Fb,F)
|
||||
|
||||
! Block-by-block Fock matrix
|
||||
|
||||
! call AOtoMO_transform(nBas,c,Fa)
|
||||
! call AOtoMO_transform(nBas,c,Fb)
|
||||
call AOtoMO_transform(nBas,c,Fa)
|
||||
call AOtoMO_transform(nBas,c,Fb)
|
||||
|
||||
F(1:nC, 1:nC ) = aC*Fa(1:nC, 1:nC ) + bC*Fb(1:nC, 1:nC )
|
||||
F(1:nC, nC+1:nC+nO ) = Fb(1:nC, nC+1:nC+nO )
|
||||
@ -64,8 +74,8 @@ subroutine ROHF_fock_matrix(nBas,nOa,nOb,S,c,Fa,Fb,F)
|
||||
F(nO+nC+1:nC+nO+nV, nC+1:nC+nO ) = Fa(nO+nC+1:nC+nO+nV, nC+1:nC+nO )
|
||||
F(nO+nC+1:nC+nO+nV,nO+nC+1:nC+nO+nV) = aV*Fa(nO+nC+1:nC+nO+nV,nO+nC+1:nC+nO+nV) + bV*Fb(nO+nC+1:nC+nO+nV,nO+nC+1:nC+nO+nV)
|
||||
|
||||
! call MOtoAO_transform(nBas,S,c,F)
|
||||
! call MOtoAO_transform(nBas,S,c,Fa)
|
||||
! call MOtoAO_transform(nBas,S,c,Fb)
|
||||
call MOtoAO_transform(nBas,S,c,F)
|
||||
call MOtoAO_transform(nBas,S,c,Fa)
|
||||
call MOtoAO_transform(nBas,S,c,Fb)
|
||||
|
||||
end subroutine
|
||||
|
Loading…
Reference in New Issue
Block a user