10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-11-14 01:53:55 +01:00
quantum_package/plugins/Hartree_Fock/diagonalize_fock.irp.f

120 lines
3.4 KiB
Fortran
Raw Normal View History

2015-06-17 18:22:08 +02:00
BEGIN_PROVIDER [ double precision, diagonal_Fock_matrix_mo, (ao_num) ]
&BEGIN_PROVIDER [ double precision, eigenvectors_Fock_matrix_mo, (ao_num_align,mo_tot_num) ]
implicit none
BEGIN_DOC
! Diagonal Fock matrix in the MO basis
END_DOC
integer :: i,j
integer :: liwork, lwork, n, info
integer, allocatable :: iwork(:)
double precision, allocatable :: work(:), F(:,:), S(:,:)
2015-11-27 19:14:36 +01:00
allocate( F(mo_tot_num_align,mo_tot_num) )
do j=1,mo_tot_num
do i=1,mo_tot_num
F(i,j) = Fock_matrix_mo(i,j)
enddo
enddo
2016-02-17 10:52:57 +01:00
if(no_oa_or_av_opt)then
integer :: iorb,jorb
do i = 1, n_act_orb
iorb = list_act(i)
do j = 1, n_inact_orb
jorb = list_inact(j)
F(iorb,jorb) = 0.d0
F(jorb,iorb) = 0.d0
enddo
do j = 1, n_virt_orb
jorb = list_virt(j)
F(iorb,jorb) = 0.d0
F(jorb,iorb) = 0.d0
enddo
do j = 1, n_core_orb
jorb = list_core(j)
F(iorb,jorb) = 0.d0
F(jorb,iorb) = 0.d0
enddo
enddo
endif
2015-11-27 19:14:36 +01:00
2016-01-02 21:45:28 +01:00
! Insert level shift here
2016-01-18 23:11:55 +01:00
do i = elec_beta_num+1, elec_alpha_num
F(i,i) += 0.5d0*level_shift
enddo
2016-01-02 21:45:28 +01:00
do i = elec_alpha_num+1, mo_tot_num
2016-01-18 23:11:55 +01:00
F(i,i) += level_shift
2016-01-02 21:45:28 +01:00
enddo
2015-11-27 19:14:36 +01:00
n = mo_tot_num
lwork = 1+6*n + 2*n*n
liwork = 3 + 5*n
allocate(work(lwork), iwork(liwork) )
lwork = -1
liwork = -1
call dsyevd( 'V', 'U', mo_tot_num, F, &
size(F,1), diagonal_Fock_matrix_mo, &
work, lwork, iwork, liwork, info)
if (info /= 0) then
print *, irp_here//' failed : ', info
stop 1
endif
lwork = int(work(1))
liwork = iwork(1)
deallocate(work,iwork)
allocate(work(lwork), iwork(liwork) )
call dsyevd( 'V', 'U', mo_tot_num, F, &
size(F,1), diagonal_Fock_matrix_mo, &
work, lwork, iwork, liwork, info)
if (info /= 0) then
print *, irp_here//' failed : ', info
stop 1
endif
2015-06-17 18:22:08 +02:00
2015-11-27 19:14:36 +01:00
call dgemm('N','N',ao_num,mo_tot_num,mo_tot_num, 1.d0, &
mo_coef, size(mo_coef,1), F, size(F,1), &
0.d0, eigenvectors_Fock_matrix_mo, size(eigenvectors_Fock_matrix_mo,1))
deallocate(work, iwork, F)
2015-06-17 18:22:08 +02:00
2016-01-02 21:45:28 +01:00
! endif
2015-06-17 18:22:08 +02:00
END_PROVIDER
BEGIN_PROVIDER [double precision, diagonal_Fock_matrix_mo_sum, (mo_tot_num)]
implicit none
BEGIN_DOC
! diagonal element of the fock matrix calculated as the sum over all the interactions
! with all the electrons in the RHF determinant
! diagonal_Fock_matrix_mo_sum(i) = sum_{j=1, N_elec} 2 J_ij -K_ij
END_DOC
integer :: i,j
double precision :: accu
do j = 1,elec_alpha_num
2015-06-17 18:22:08 +02:00
accu = 0.d0
do i = 1, elec_alpha_num
2015-06-17 18:22:08 +02:00
accu += 2.d0 * mo_bielec_integral_jj_from_ao(i,j) - mo_bielec_integral_jj_exchange_from_ao(i,j)
enddo
diagonal_Fock_matrix_mo_sum(j) = accu + mo_mono_elec_integral(j,j)
2015-06-17 18:22:08 +02:00
enddo
do j = elec_alpha_num+1,mo_tot_num
2015-06-17 18:22:08 +02:00
accu = 0.d0
do i = 1, elec_alpha_num
2015-06-17 18:22:08 +02:00
accu += 2.d0 * mo_bielec_integral_jj_from_ao(i,j) - mo_bielec_integral_jj_exchange_from_ao(i,j)
enddo
diagonal_Fock_matrix_mo_sum(j) = accu + mo_mono_elec_integral(j,j)
2015-06-17 18:22:08 +02:00
enddo
END_PROVIDER