10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-12-23 12:56:14 +01:00
quantum_package/plugins/Full_CI_ZMQ/fci_zmq.irp.f

128 lines
3.8 KiB
Fortran
Raw Normal View History

program fci_zmq
implicit none
2016-09-08 10:12:28 +02:00
integer :: i,j,k
2016-09-25 22:14:17 +02:00
logical, external :: detEq
2016-09-23 16:16:48 +02:00
double precision, allocatable :: pt2(:)
integer :: degree
integer :: n_det_before, to_select
2016-11-22 13:00:02 +01:00
double precision :: threshold_davidson_in
2016-09-25 22:14:17 +02:00
2016-09-23 16:16:48 +02:00
allocate (pt2(N_states))
pt2 = 1.d0
2016-11-22 13:00:02 +01:00
threshold_davidson_in = threshold_davidson
2016-12-05 09:28:04 +01:00
threshold_davidson = threshold_davidson_in * 100.d0
2016-11-22 13:00:02 +01:00
SOFT_TOUCH threshold_davidson
if (N_det > N_det_max) then
call diagonalize_CI
call save_wavefunction
psi_det = psi_det_sorted
psi_coef = psi_coef_sorted
N_det = N_det_max
soft_touch N_det psi_det psi_coef
call diagonalize_CI
call save_wavefunction
print *, 'N_det = ', N_det
print *, 'N_states = ', N_states
2016-09-23 16:16:48 +02:00
do k=1,N_states
print*,'State ',k
print *, 'PT2 = ', pt2(k)
print *, 'E = ', CI_energy(k)
print *, 'E+PT2 = ', CI_energy(k) + pt2(k)
print *, '-----'
enddo
endif
2016-09-25 22:14:17 +02:00
double precision :: E_CI_before(N_states)
2016-07-21 13:24:25 +02:00
2016-09-25 22:14:17 +02:00
print*,'Beginning the selection ...'
2016-09-23 16:16:48 +02:00
E_CI_before(1:N_states) = CI_energy(1:N_states)
n_det_before = 0
2016-08-01 16:05:40 +02:00
2016-09-23 16:16:48 +02:00
do while ( (N_det < N_det_max) .and. (maxval(abs(pt2(1:N_states))) > pt2_max) )
2016-09-25 22:14:17 +02:00
print *, 'N_det = ', N_det
print *, 'N_states = ', N_states
2016-09-23 16:16:48 +02:00
do k=1, N_states
print*,'State ',k
print *, 'PT2 = ', pt2(k)
print *, 'E = ', CI_energy(k)
print *, 'E(before)+PT2 = ', E_CI_before(k)+pt2(k)
enddo
print *, '-----'
if(N_states.gt.1)then
2016-09-25 22:14:17 +02:00
print*,'Variational Energy difference'
do i = 2, N_states
print*,'Delta E = ',CI_energy(i) - CI_energy(1)
enddo
endif
if(N_states.gt.1)then
2016-09-25 22:14:17 +02:00
print*,'Variational + perturbative Energy difference'
do i = 2, N_states
print*,'Delta E = ',E_CI_before(i)+ pt2(i) - (E_CI_before(1) + pt2(1))
enddo
endif
2016-09-23 16:16:48 +02:00
E_CI_before(1:N_states) = CI_energy(1:N_states)
2016-11-30 17:33:34 +01:00
call ezfio_set_full_ci_zmq_energy(CI_energy(1))
n_det_before = N_det
2017-01-16 16:31:49 +01:00
to_select = N_det
2017-03-03 12:02:21 +01:00
to_select = max(N_det, to_select)
2016-11-30 17:33:34 +01:00
to_select = min(to_select, N_det_max-n_det_before)
call ZMQ_selection(to_select, pt2)
PROVIDE psi_coef
PROVIDE psi_det
PROVIDE psi_det_sorted
if (N_det == N_det_max) then
threshold_davidson = threshold_davidson_in
2016-12-01 16:28:56 +01:00
SOFT_TOUCH threshold_davidson
2016-11-30 17:33:34 +01:00
endif
call diagonalize_CI
call save_wavefunction
call ezfio_set_full_ci_zmq_energy(CI_energy(1))
enddo
2016-10-10 16:03:56 +02:00
2016-11-30 17:33:34 +01:00
if (N_det < N_det_max) then
threshold_davidson = threshold_davidson_in
SOFT_TOUCH threshold_davidson
call diagonalize_CI
call save_wavefunction
call ezfio_set_full_ci_zmq_energy(CI_energy(1))
endif
2016-09-25 22:14:17 +02:00
if(do_pt2_end)then
print*,'Last iteration only to compute the PT2'
2017-01-12 16:48:34 +01:00
!threshold_selectors = max(threshold_selectors,threshold_selectors_pt2)
!threshold_generators = max(threshold_generators,threshold_generators_pt2)
!TOUCH threshold_selectors threshold_generators
2016-09-25 22:14:17 +02:00
threshold_selectors = 1.d0
2017-01-12 16:48:34 +01:00
threshold_generators = 1d0
2016-09-25 22:14:17 +02:00
E_CI_before(1:N_states) = CI_energy(1:N_states)
2017-01-16 18:19:00 +01:00
double precision :: relative_error
2017-01-16 18:51:01 +01:00
relative_error=1.d-3
2017-03-13 00:26:21 +01:00
pt2 = 0.d0
call ZMQ_pt2(pt2,relative_error)
!call ZMQ_selection(0, pt2)! pour non-stochastic
2016-09-25 22:14:17 +02:00
print *, 'Final step'
print *, 'N_det = ', N_det
print *, 'N_states = ', N_states
do k=1,N_states
print *, 'State', k
print *, 'PT2 = ', pt2
print *, 'E = ', E_CI_before
print *, 'E+PT2 = ', E_CI_before+pt2
print *, '-----'
enddo
2016-11-30 17:33:34 +01:00
call ezfio_set_full_ci_zmq_energy_pt2(E_CI_before(1)+pt2(1))
2016-09-25 22:14:17 +02:00
endif
call save_wavefunction
2016-11-30 17:33:34 +01:00
call ezfio_set_full_ci_zmq_energy(CI_energy(1))
call ezfio_set_full_ci_zmq_energy_pt2(E_CI_before(1)+pt2(1))
end