10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-11-08 15:13:52 +01:00
quantum_package/src/Determinants/diagonalize_CI.irp.f

119 lines
3.4 KiB
Fortran
Raw Normal View History

2015-04-20 16:45:06 +02:00
BEGIN_PROVIDER [ character*(64), diag_algorithm ]
implicit none
BEGIN_DOC
! Diagonalization algorithm (Davidson or Lapack)
END_DOC
if (N_det > N_det_max_jacobi) then
diag_algorithm = "Davidson"
else
diag_algorithm = "Lapack"
endif
if (N_det < N_states_diag) then
diag_algorithm = "Lapack"
endif
END_PROVIDER
BEGIN_PROVIDER [ double precision, CI_energy, (N_states_diag) ]
implicit none
BEGIN_DOC
! N_states lowest eigenvalues of the CI matrix
END_DOC
integer :: j
character*(8) :: st
call write_time(output_determinants)
do j=1,N_states_diag
CI_energy(j) = CI_electronic_energy(j) + nuclear_repulsion
write(st,'(I4)') j
call write_double(output_determinants,CI_energy(j),'Energy of state '//trim(st))
call write_double(output_determinants,CI_eigenvectors_s2(j),'S^2 of state '//trim(st))
enddo
END_PROVIDER
BEGIN_PROVIDER [ double precision, CI_electronic_energy, (N_states_diag) ]
&BEGIN_PROVIDER [ double precision, CI_eigenvectors, (N_det,N_states_diag) ]
&BEGIN_PROVIDER [ double precision, CI_eigenvectors_s2, (N_states_diag) ]
implicit none
BEGIN_DOC
! Eigenvectors/values of the CI matrix
END_DOC
integer :: i,j
do j=1,N_states_diag
do i=1,N_det
CI_eigenvectors(i,j) = psi_coef(i,j)
enddo
enddo
if (diag_algorithm == "Davidson") then
call davidson_diag(psi_det,CI_eigenvectors,CI_electronic_energy, &
size(CI_eigenvectors,1),N_det,N_states_diag,N_int,output_determinants)
2015-06-29 10:35:29 +02:00
do j=1,N_states_diag
call get_s2_u0(psi_det,CI_eigenvectors(1,j),N_det,size(CI_eigenvectors,1),CI_eigenvectors_s2(j))
enddo
2015-04-20 16:45:06 +02:00
else if (diag_algorithm == "Lapack") then
double precision, allocatable :: eigenvectors(:,:), eigenvalues(:)
allocate (eigenvectors(size(H_matrix_all_dets,1),N_det))
allocate (eigenvalues(N_det))
call lapack_diag(eigenvalues,eigenvectors, &
H_matrix_all_dets,size(H_matrix_all_dets,1),N_det)
CI_electronic_energy(:) = 0.d0
do i=1,N_det
CI_eigenvectors(i,1) = eigenvectors(i,1)
enddo
integer :: i_state
double precision :: s2
if (s2_eig) then
i_state = 0
do j=1,N_det
2015-06-29 10:35:29 +02:00
call get_s2_u0(psi_det,eigenvectors(1,j),N_det,size(eigenvectors,1),s2)
print*,'s2 = ',s2
if(dabs(s2-expected_s2).le.0.3d0)then
i_state += 1
do i=1,N_det
CI_eigenvectors(i,i_state) = eigenvectors(i,j)
enddo
CI_electronic_energy(i_state) = eigenvalues(j)
CI_eigenvectors_s2(i_state) = s2
endif
if (i_state.ge.N_states_diag) then
exit
endif
enddo
else
2015-05-06 18:26:36 +02:00
do j=1,N_states_diag
call get_s2_u0(psi_det,eigenvectors(1,j),N_det,N_det,s2)
do i=1,N_det
CI_eigenvectors(i,j) = eigenvectors(i,j)
enddo
CI_electronic_energy(j) = eigenvalues(j)
CI_eigenvectors_s2(j) = s2
enddo
endif
2015-04-20 16:45:06 +02:00
deallocate(eigenvectors,eigenvalues)
endif
END_PROVIDER
subroutine diagonalize_CI
implicit none
BEGIN_DOC
! Replace the coefficients of the CI states by the coefficients of the
! eigenstates of the CI matrix
END_DOC
integer :: i,j
do j=1,N_states_diag
do i=1,N_det
psi_coef(i,j) = CI_eigenvectors(i,j)
enddo
enddo
SOFT_TOUCH psi_coef CI_electronic_energy CI_energy CI_eigenvectors CI_eigenvectors_s2
end