10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-12-23 04:43:50 +01:00
quantum_package/plugins/Hartree_Fock/DIIS.irp.f

140 lines
4.9 KiB
Fortran
Raw Normal View History

2017-06-02 23:53:44 +02:00
BEGIN_PROVIDER [ double precision, threshold_DIIS_nonzero ]
implicit none
BEGIN_DOC
! If threshold_DIIS is zero, choose sqrt(thresh_scf)
END_DOC
if (threshold_DIIS == 0.d0) then
threshold_DIIS_nonzero = dsqrt(thresh_scf)
else
threshold_DIIS_nonzero = threshold_DIIS
endif
ASSERT (threshold_DIIS_nonzero >= 0.d0)
2017-06-02 14:20:59 +02:00
2017-06-02 23:53:44 +02:00
END_PROVIDER
2017-06-02 23:53:44 +02:00
BEGIN_PROVIDER [double precision, FPS_SPF_Matrix_AO, (AO_num, AO_num)]
implicit none
BEGIN_DOC
! Commutator FPS - SPF
END_DOC
double precision, allocatable :: scratch(:,:)
allocate( &
2017-06-27 20:41:54 +02:00
scratch(AO_num, AO_num) &
2017-06-02 23:53:44 +02:00
)
! Compute FP
call dgemm('N','N',AO_num,AO_num,AO_num, &
1.d0, &
Fock_Matrix_AO,Size(Fock_Matrix_AO,1), &
HF_Density_Matrix_AO,Size(HF_Density_Matrix_AO,1), &
0.d0, &
scratch,Size(scratch,1))
! Compute FPS
call dgemm('N','N',AO_num,AO_num,AO_num, &
1.d0, &
scratch,Size(scratch,1), &
AO_Overlap,Size(AO_Overlap,1), &
0.d0, &
FPS_SPF_Matrix_AO,Size(FPS_SPF_Matrix_AO,1))
! Compute SP
call dgemm('N','N',AO_num,AO_num,AO_num, &
1.d0, &
AO_Overlap,Size(AO_Overlap,1), &
HF_Density_Matrix_AO,Size(HF_Density_Matrix_AO,1), &
0.d0, &
scratch,Size(scratch,1))
! Compute FPS - SPF
call dgemm('N','N',AO_num,AO_num,AO_num, &
-1.d0, &
scratch,Size(scratch,1), &
Fock_Matrix_AO,Size(Fock_Matrix_AO,1), &
1.d0, &
FPS_SPF_Matrix_AO,Size(FPS_SPF_Matrix_AO,1))
2017-06-02 14:20:59 +02:00
2017-06-02 23:53:44 +02:00
END_PROVIDER
2017-06-02 23:53:44 +02:00
bEGIN_PROVIDER [double precision, FPS_SPF_Matrix_MO, (mo_tot_num, mo_tot_num)]
2017-06-02 14:20:59 +02:00
implicit none
begin_doc
! Commutator FPS - SPF in MO basis
end_doc
2017-06-02 23:53:44 +02:00
call ao_to_mo(FPS_SPF_Matrix_AO, size(FPS_SPF_Matrix_AO,1), &
FPS_SPF_Matrix_MO, size(FPS_SPF_Matrix_MO,1))
END_PROVIDER
2017-06-02 14:20:59 +02:00
BEGIN_PROVIDER [ double precision, eigenvalues_Fock_matrix_AO, (AO_num) ]
2017-06-27 20:41:54 +02:00
&BEGIN_PROVIDER [ double precision, eigenvectors_Fock_matrix_AO, (AO_num,AO_num) ]
BEGIN_DOC
! Eigenvalues and eigenvectors of the Fock matrix over the AO basis
END_DOC
implicit none
double precision, allocatable :: scratch(:,:),work(:),Xt(:,:)
integer :: lwork,info
integer :: i,j
lwork = 3*AO_num - 1
2017-06-02 23:53:44 +02:00
allocate( &
2017-06-27 20:41:54 +02:00
scratch(AO_num,AO_num), &
2017-06-02 23:53:44 +02:00
work(lwork), &
Xt(AO_num,AO_num) &
)
! Calculate Xt
do i=1,AO_num
do j=1,AO_num
2017-09-13 09:06:32 +02:00
Xt(i,j) = S_half_inv(j,i)
enddo
enddo
! Calculate Fock matrix in orthogonal basis: F' = Xt.F.X
call dgemm('N','N',AO_num,AO_num,AO_num, &
1.d0, &
Fock_matrix_AO,size(Fock_matrix_AO,1), &
2017-09-13 09:06:32 +02:00
S_half_inv,size(S_half_inv,1), &
0.d0, &
eigenvectors_Fock_matrix_AO,size(eigenvectors_Fock_matrix_AO,1))
call dgemm('N','N',AO_num,AO_num,AO_num, &
1.d0, &
Xt,size(Xt,1), &
eigenvectors_Fock_matrix_AO,size(eigenvectors_Fock_matrix_AO,1), &
0.d0, &
scratch,size(scratch,1))
! Diagonalize F' to obtain eigenvectors in orthogonal basis C' and eigenvalues
call dsyev('V','U',AO_num, &
scratch,size(scratch,1), &
eigenvalues_Fock_matrix_AO, &
work,lwork,info)
if(info /= 0) then
print *, irp_here//' failed : ', info
stop 1
endif
! Back-transform eigenvectors: C =X.C'
call dgemm('N','N',AO_num,AO_num,AO_num, &
1.d0, &
2017-09-13 09:06:32 +02:00
S_half_inv,size(S_half_inv,1), &
scratch,size(scratch,1), &
0.d0, &
eigenvectors_Fock_matrix_AO,size(eigenvectors_Fock_matrix_AO,1))
END_PROVIDER