9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-12-22 03:23:29 +01:00
qp2/plugins/local/cipsi_tc_bi_ortho/stochastic_cipsi.irp.f

128 lines
3.8 KiB
Fortran
Raw Normal View History

2023-07-02 21:49:25 +02:00
! ---
2023-02-07 17:28:11 +01:00
subroutine run_stochastic_cipsi
2023-07-02 21:49:25 +02:00
2023-02-07 17:28:11 +01:00
BEGIN_DOC
2023-07-02 21:49:25 +02:00
! Selected Full Configuration Interaction with Stochastic selection and PT2.
2023-02-07 17:28:11 +01:00
END_DOC
2023-07-02 21:49:25 +02:00
use selection_types
implicit none
integer :: i, j, k, ndet
integer :: to_select
logical :: has
type(pt2_type) :: pt2_data, pt2_data_err
double precision :: rss
double precision :: correlation_energy_ratio
2023-07-02 21:49:25 +02:00
double precision :: hf_energy_ref
double precision :: relative_error
double precision, allocatable :: zeros(:),E_tc(:), norm(:)
2023-07-02 21:49:25 +02:00
logical, external :: qp_stop
double precision, external :: memory_of_double
PROVIDE mo_l_coef mo_r_coef
2023-02-07 17:28:11 +01:00
PROVIDE H_apply_buffer_allocated distributed_davidson
2023-07-02 21:49:25 +02:00
print*, ' Diagonal elements of the Fock matrix '
2023-02-07 17:28:11 +01:00
do i = 1, mo_num
2023-07-02 21:49:25 +02:00
write(*,*) i, Fock_matrix_tc_mo_tot(i,i)
2023-02-07 17:28:11 +01:00
enddo
2023-07-02 21:49:25 +02:00
2023-02-07 17:28:11 +01:00
threshold_generators = 1.d0
SOFT_TOUCH threshold_generators
rss = memory_of_double(N_states)*4.d0
2023-07-02 21:49:25 +02:00
call check_mem(rss, irp_here)
2023-02-07 17:28:11 +01:00
allocate(zeros(N_states),E_tc(N_states), norm(N_states))
2023-02-07 17:28:11 +01:00
call pt2_alloc(pt2_data, N_states)
call pt2_alloc(pt2_data_err, N_states)
2023-07-02 21:49:25 +02:00
relative_error = PT2_relative_error
2023-02-07 17:28:11 +01:00
2023-07-02 21:49:25 +02:00
zeros = 0.d0
pt2_data % pt2 = -huge(1.e0)
pt2_data % rpt2 = -huge(1.e0)
pt2_data % overlap = 0.d0
2023-02-07 17:28:11 +01:00
pt2_data % variance = huge(1.e0)
!!!! WARNING !!!! SEEMS TO BE PROBLEM WTH make_s2_eigenfunction !!!! THE DETERMINANTS CAN APPEAR TWICE IN THE WFT DURING SELECTION
! if (s2_eig) then
! call make_s2_eigenfunction
! endif
call diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm)
2023-02-07 17:28:11 +01:00
! if (N_det > N_det_max) then
2023-03-15 11:55:03 +01:00
! psi_det(1:N_int,1:2,1:N_det) = psi_det_generators(1:N_int,1:2,1:N_det)
! psi_coef(1:N_det,1:N_states) = psi_coef_sorted_gen(1:N_det,1:N_states)
! N_det = N_det_max
! soft_touch N_det psi_det psi_coef
! if (s2_eig) then
! call make_s2_eigenfunction
! endif
! call diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm)
2023-02-07 17:28:11 +01:00
! call routine_save_right
! endif
2023-02-07 17:28:11 +01:00
correlation_energy_ratio = 0.d0
! thresh_it_dav = 5.d-5
! soft_touch thresh_it_dav
2023-07-02 21:49:25 +02:00
do while( (N_det < N_det_max) .and. &
(maxval(abs(pt2_data % pt2(1:N_states))) > pt2_max))
2023-02-07 17:28:11 +01:00
2023-07-02 21:49:25 +02:00
print*,'maxval(abs(pt2_data % pt2(1:N_states)))',maxval(abs(pt2_data % pt2(1:N_states)))
print*,pt2_max
write(*,'(A)') '--------------------------------------------------------------------------------'
2023-02-07 17:28:11 +01:00
to_select = int(sqrt(dble(N_states))*dble(N_det)*selection_factor)
to_select = max(N_states_diag, to_select)
2024-03-12 17:45:50 +01:00
print*,'E_tc = ',E_tc
2023-02-07 17:28:11 +01:00
call pt2_dealloc(pt2_data)
call pt2_dealloc(pt2_data_err)
call pt2_alloc(pt2_data, N_states)
call pt2_alloc(pt2_data_err, N_states)
call ZMQ_pt2(E_tc, pt2_data, pt2_data_err, relative_error,to_select) ! Stochastic PT2 and selection
2023-02-07 17:28:11 +01:00
! stop
2024-03-12 15:30:52 +01:00
call print_summary_tc(psi_energy_with_nucl_rep, pt2_data, pt2_data_err, N_det, N_configuration, N_states, psi_s2)
2023-04-24 01:22:24 +02:00
call save_energy(psi_energy_with_nucl_rep, pt2_data % pt2)
call increment_n_iter(psi_energy_with_nucl_rep, pt2_data)
call print_extrapolated_energy()
2023-04-24 01:32:05 +02:00
! call print_mol_properties()
2023-04-24 01:22:24 +02:00
call write_cipsi_json(pt2_data,pt2_data_err)
2023-02-07 17:28:11 +01:00
if (qp_stop()) exit
! Add selected determinants
call copy_H_apply_buffer_to_wf_tc()
2023-07-02 21:49:25 +02:00
PROVIDE psi_l_coef_bi_ortho psi_r_coef_bi_ortho
PROVIDE psi_det
PROVIDE psi_det_sorted_tc
2023-02-07 17:28:11 +01:00
call diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm)
! stop
2023-02-07 17:28:11 +01:00
if (qp_stop()) exit
enddo
call pt2_dealloc(pt2_data)
call pt2_dealloc(pt2_data_err)
call pt2_alloc(pt2_data, N_states)
call pt2_alloc(pt2_data_err, N_states)
call ZMQ_pt2(E_tc, pt2_data, pt2_data_err, relative_error,0) ! Stochastic PT2 and selection
call diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm)
call pt2_dealloc(pt2_data)
call pt2_dealloc(pt2_data_err)
2023-02-07 17:28:11 +01:00
end