10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-06-02 11:25:26 +02:00
quantum_package/plugins/Hartree_Fock/Roothaan_Hall_SCF.irp.f

286 lines
8.4 KiB
Fortran
Raw Permalink Normal View History

subroutine Roothaan_Hall_SCF
BEGIN_DOC
! Roothaan-Hall algorithm for SCF Hartree-Fock calculation
END_DOC
implicit none
2017-06-02 14:20:59 +02:00
double precision :: energy_SCF,energy_SCF_previous,Delta_energy_SCF
double precision :: max_error_DIIS,max_error_DIIS_alpha,max_error_DIIS_beta
2017-06-02 23:53:44 +02:00
double precision, allocatable :: Fock_matrix_DIIS(:,:,:),error_matrix_DIIS(:,:,:)
integer :: iteration_SCF,dim_DIIS,index_dim_DIIS
integer :: i,j
2017-06-16 15:59:23 +02:00
double precision, allocatable :: mo_coef_save(:,:)
2018-01-10 16:54:12 +01:00
PROVIDE ao_md5 mo_occ level_shift
2017-06-16 15:59:23 +02:00
allocate(mo_coef_save(ao_num,mo_tot_num), &
2017-06-02 23:53:44 +02:00
Fock_matrix_DIIS (ao_num,ao_num,max_dim_DIIS), &
error_matrix_DIIS(ao_num,ao_num,max_dim_DIIS) &
)
2018-01-05 15:12:27 +01:00
call write_time(6)
2018-01-05 15:12:27 +01:00
write(6,'(A4, 1X, A16, 1X, A16, 1X, A16)') &
'====','================','================','================'
2018-01-05 15:12:27 +01:00
write(6,'(A4, 1X, A16, 1X, A16, 1X, A16)') &
' N ', 'Energy ', 'Energy diff ', 'DIIS error '
2018-01-05 15:12:27 +01:00
write(6,'(A4, 1X, A16, 1X, A16, 1X, A16)') &
'====','================','================','================'
! Initialize energies and density matrices
energy_SCF_previous = HF_energy
2017-06-02 23:53:44 +02:00
Delta_energy_SCF = 1.d0
iteration_SCF = 0
dim_DIIS = 0
max_error_DIIS = 1.d0
!
! Start of main SCF loop
!
2017-06-03 00:09:02 +02:00
do while(( (max_error_DIIS > threshold_DIIS_nonzero).or.(dabs(Delta_energy_SCF) > thresh_SCF) ) .and. (iteration_SCF < n_it_SCF_max))
! Increment cycle number
iteration_SCF += 1
! Current size of the DIIS space
dim_DIIS = min(dim_DIIS+1,max_dim_DIIS)
2017-06-02 23:53:44 +02:00
if (scf_algorithm == 'DIIS') then
! Store Fock and error matrices at each iteration
do j=1,ao_num
do i=1,ao_num
index_dim_DIIS = mod(dim_DIIS-1,max_dim_DIIS)+1
Fock_matrix_DIIS (i,j,index_dim_DIIS) = Fock_matrix_AO(i,j)
error_matrix_DIIS(i,j,index_dim_DIIS) = FPS_SPF_matrix_AO(i,j)
enddo
enddo
2017-06-02 23:53:44 +02:00
! Compute the extrapolated Fock matrix
call extrapolate_Fock_matrix( &
error_matrix_DIIS,Fock_matrix_DIIS, &
Fock_matrix_AO,size(Fock_matrix_AO,1), &
iteration_SCF,dim_DIIS &
)
Fock_matrix_AO_alpha = Fock_matrix_AO*0.5d0
2017-06-03 00:09:02 +02:00
Fock_matrix_AO_beta = Fock_matrix_AO*0.5d0
TOUCH Fock_matrix_AO_alpha Fock_matrix_AO_beta
2017-06-02 23:53:44 +02:00
endif
2017-06-02 14:20:59 +02:00
MO_coef = eigenvectors_Fock_matrix_MO
2017-06-03 00:09:02 +02:00
TOUCH MO_coef
2017-06-02 23:53:44 +02:00
! Calculate error vectors
2017-06-02 23:53:44 +02:00
max_error_DIIS = maxval(Abs(FPS_SPF_Matrix_MO))
! SCF energy
energy_SCF = HF_energy
Delta_Energy_SCF = energy_SCF - energy_SCF_previous
2017-06-03 00:09:02 +02:00
if ( (SCF_algorithm == 'DIIS').and.(Delta_Energy_SCF > 0.d0) ) then
Fock_matrix_AO(1:ao_num,1:ao_num) = Fock_matrix_DIIS (1:ao_num,1:ao_num,index_dim_DIIS)
Fock_matrix_AO_alpha = Fock_matrix_AO*0.5d0
Fock_matrix_AO_beta = Fock_matrix_AO*0.5d0
TOUCH Fock_matrix_AO_alpha Fock_matrix_AO_beta
endif
double precision :: level_shift_save
level_shift_save = level_shift
2017-06-16 15:59:23 +02:00
mo_coef_save(1:ao_num,1:mo_tot_num) = mo_coef(1:ao_num,1:mo_tot_num)
2017-06-03 00:09:02 +02:00
do while (Delta_Energy_SCF > 0.d0)
2017-06-16 15:59:23 +02:00
mo_coef(1:ao_num,1:mo_tot_num) = mo_coef_save
TOUCH mo_coef
2017-06-07 23:05:00 +02:00
level_shift = level_shift + 0.1d0
2017-06-16 15:59:23 +02:00
mo_coef(1:ao_num,1:mo_tot_num) = eigenvectors_Fock_matrix_MO(1:ao_num,1:mo_tot_num)
TOUCH mo_coef level_shift
2017-06-03 00:09:02 +02:00
Delta_Energy_SCF = HF_energy - energy_SCF_previous
2017-06-07 23:05:00 +02:00
energy_SCF = HF_energy
2017-06-16 15:59:23 +02:00
if (level_shift-level_shift_save > 1.d0) exit
dim_DIIS=0
2017-06-03 00:09:02 +02:00
enddo
level_shift = level_shift_save
SOFT_TOUCH level_shift
energy_SCF_previous = energy_SCF
2017-06-02 23:53:44 +02:00
! Print results at the end of each iteration
2018-01-05 15:12:27 +01:00
write(6,'(I4, 1X, F16.10, 1X, F16.10, 1X, F16.10, 1X, I3)') &
2017-06-02 14:20:59 +02:00
iteration_SCF, energy_SCF, Delta_energy_SCF, max_error_DIIS, dim_DIIS
2017-06-02 23:53:44 +02:00
if (Delta_energy_SCF < 0.d0) then
call save_mos
endif
enddo
!
! End of Main SCF loop
!
2018-01-05 15:12:27 +01:00
write(6,'(A4, 1X, A16, 1X, A16, 1X, A16)') &
'====','================','================','================'
2018-01-05 15:12:27 +01:00
write(6,*)
if(.not.no_oa_or_av_opt)then
2017-09-13 09:06:32 +02:00
call mo_as_eigvectors_of_mo_matrix(Fock_matrix_mo,size(Fock_matrix_mo,1),size(Fock_matrix_mo,2),mo_label,1,.true.)
endif
2018-01-05 15:12:27 +01:00
call write_double(6, Energy_SCF, 'Hartree-Fock energy')
call ezfio_set_hartree_fock_energy(Energy_SCF)
2018-01-05 15:12:27 +01:00
call write_time(6)
end
2017-06-02 14:20:59 +02:00
subroutine extrapolate_Fock_matrix( &
error_matrix_DIIS,Fock_matrix_DIIS, &
Fock_matrix_AO_,size_Fock_matrix_AO, &
iteration_SCF,dim_DIIS &
)
BEGIN_DOC
! Compute the extrapolated Fock matrix using the DIIS procedure
END_DOC
implicit none
2017-06-02 14:20:59 +02:00
double precision,intent(in) :: Fock_matrix_DIIS(ao_num,ao_num,*),error_matrix_DIIS(ao_num,ao_num,*)
integer,intent(in) :: iteration_SCF, size_Fock_matrix_AO
double precision,intent(inout):: Fock_matrix_AO_(size_Fock_matrix_AO,ao_num)
integer,intent(inout) :: dim_DIIS
double precision,allocatable :: B_matrix_DIIS(:,:),X_vector_DIIS(:)
2017-06-02 14:20:59 +02:00
double precision,allocatable :: C_vector_DIIS(:)
double precision,allocatable :: scratch(:,:)
integer :: i,j,k,i_DIIS,j_DIIS
allocate( &
B_matrix_DIIS(dim_DIIS+1,dim_DIIS+1), &
X_vector_DIIS(dim_DIIS+1), &
2017-06-02 14:20:59 +02:00
C_vector_DIIS(dim_DIIS+1), &
scratch(ao_num,ao_num) &
)
! Compute the matrices B and X
do j=1,dim_DIIS
do i=1,dim_DIIS
j_DIIS = mod(iteration_SCF-j,max_dim_DIIS)+1
i_DIIS = mod(iteration_SCF-i,max_dim_DIIS)+1
! Compute product of two errors vectors
2017-06-02 14:20:59 +02:00
call dgemm('N','N',ao_num,ao_num,ao_num, &
1.d0, &
error_matrix_DIIS(1,1,i_DIIS),size(error_matrix_DIIS,1), &
error_matrix_DIIS(1,1,j_DIIS),size(error_matrix_DIIS,1), &
0.d0, &
scratch,size(scratch,1))
! Compute Trace
B_matrix_DIIS(i,j) = 0.d0
2017-06-02 14:20:59 +02:00
do k=1,ao_num
2017-06-02 23:53:44 +02:00
B_matrix_DIIS(i,j) = B_matrix_DIIS(i,j) + scratch(k,k)
enddo
enddo
enddo
! Pad B matrix and build the X matrix
do i=1,dim_DIIS
B_matrix_DIIS(i,dim_DIIS+1) = -1.d0
B_matrix_DIIS(dim_DIIS+1,i) = -1.d0
2017-06-02 14:20:59 +02:00
C_vector_DIIS(i) = 0.d0
enddo
B_matrix_DIIS(dim_DIIS+1,dim_DIIS+1) = 0.d0
2017-06-02 14:20:59 +02:00
C_vector_DIIS(dim_DIIS+1) = -1.d0
! Solve the linear system C = B.X
integer :: info
integer,allocatable :: ipiv(:)
allocate( &
ipiv(dim_DIIS+1) &
)
2017-06-02 14:20:59 +02:00
double precision, allocatable :: AF(:,:)
allocate (AF(dim_DIIS+1,dim_DIIS+1))
double precision :: rcond, ferr, berr
2017-06-08 22:15:42 +02:00
integer :: iwork(dim_DIIS+1), lwork
2017-06-02 14:20:59 +02:00
call dsysvx('N','U',dim_DIIS+1,1, &
B_matrix_DIIS,size(B_matrix_DIIS,1), &
2017-06-08 22:15:42 +02:00
AF, size(AF,1), &
ipiv, &
C_vector_DIIS,size(C_vector_DIIS,1), &
X_vector_DIIS,size(X_vector_DIIS,1), &
rcond, &
ferr, &
berr, &
scratch,-1, &
iwork, &
info &
)
lwork = int(scratch(1,1))
deallocate(scratch)
allocate(scratch(lwork,1))
call dsysvx('N','U',dim_DIIS+1,1, &
B_matrix_DIIS,size(B_matrix_DIIS,1), &
2017-06-02 14:20:59 +02:00
AF, size(AF,1), &
ipiv, &
2017-06-02 14:20:59 +02:00
C_vector_DIIS,size(C_vector_DIIS,1), &
X_vector_DIIS,size(X_vector_DIIS,1), &
2017-06-02 14:20:59 +02:00
rcond, &
ferr, &
berr, &
scratch,size(scratch), &
2017-06-02 14:20:59 +02:00
iwork, &
info &
)
2017-06-02 14:20:59 +02:00
if(info < 0) then
stop 'bug in DIIS'
endif
2017-06-07 23:05:00 +02:00
if (rcond > 1.d-12) then
2017-06-02 23:53:44 +02:00
! Compute extrapolated Fock matrix
2017-06-02 23:53:44 +02:00
!$OMP PARALLEL DO PRIVATE(i,j,k) DEFAULT(SHARED)
2017-06-02 14:20:59 +02:00
do j=1,ao_num
do i=1,ao_num
2017-06-02 23:53:44 +02:00
Fock_matrix_AO_(i,j) = 0.d0
enddo
do k=1,dim_DIIS
do i=1,ao_num
Fock_matrix_AO_(i,j) = Fock_matrix_AO_(i,j) + &
X_vector_DIIS(k)*Fock_matrix_DIIS(i,j,dim_DIIS-k+1)
enddo
enddo
enddo
2017-06-02 23:53:44 +02:00
!$OMP END PARALLEL DO
else
dim_DIIS = 0
endif
end