fci_tc_bi_ortho works for multi state ninja
continuous-integration/drone/push Build is failing Details

This commit is contained in:
eginer 2024-03-12 18:23:09 +01:00
parent fdc418d72a
commit a56488e3a8
5 changed files with 20 additions and 91 deletions

View File

@ -11,15 +11,13 @@ subroutine run_stochastic_cipsi
implicit none
integer :: i, j, k, ndet
integer :: to_select
logical :: print_pt2
logical :: has
type(pt2_type) :: pt2_data, pt2_data_err
double precision :: rss
double precision :: correlation_energy_ratio, E_denom, E_tc, norm
double precision :: correlation_energy_ratio
double precision :: hf_energy_ref
double precision :: relative_error
double precision, allocatable :: ept2(:), pt1(:), extrap_energy(:)
double precision, allocatable :: zeros(:)
double precision, allocatable :: zeros(:),E_tc(:), norm(:)
logical, external :: qp_stop
double precision, external :: memory_of_double
@ -32,14 +30,13 @@ subroutine run_stochastic_cipsi
write(*,*) i, Fock_matrix_tc_mo_tot(i,i)
enddo
N_iter = 1
threshold_generators = 1.d0
SOFT_TOUCH threshold_generators
rss = memory_of_double(N_states)*4.d0
call check_mem(rss, irp_here)
allocate(zeros(N_states))
allocate(zeros(N_states),E_tc(N_states), norm(N_states))
call pt2_alloc(pt2_data, N_states)
call pt2_alloc(pt2_data_err, N_states)
@ -55,8 +52,7 @@ subroutine run_stochastic_cipsi
! if (s2_eig) then
! call make_s2_eigenfunction
! endif
print_pt2 = .False.
call diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm, pt2_data, print_pt2)
call diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm)
! if (N_det > N_det_max) then
@ -67,19 +63,16 @@ subroutine run_stochastic_cipsi
! if (s2_eig) then
! call make_s2_eigenfunction
! endif
! print_pt2 = .False.
! call diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm,pt2_data,print_pt2)
! call diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm)
! call routine_save_right
! endif
allocate(ept2(1000),pt1(1000),extrap_energy(100))
correlation_energy_ratio = 0.d0
! thresh_it_dav = 5.d-5
! soft_touch thresh_it_dav
print_pt2 = .True.
do while( (N_det < N_det_max) .and. &
(maxval(abs(pt2_data % pt2(1:N_states))) > pt2_max))
@ -90,13 +83,12 @@ subroutine run_stochastic_cipsi
to_select = int(sqrt(dble(N_states))*dble(N_det)*selection_factor)
to_select = max(N_states_diag, to_select)
E_denom = E_tc ! TC Energy of the current wave function
print*,'E_tc = ',E_tc
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_denom, pt2_data, pt2_data_err, relative_error,to_select) ! Stochastic PT2 and selection
call ZMQ_pt2(E_tc, pt2_data, pt2_data_err, relative_error,to_select) ! Stochastic PT2 and selection
! stop
call print_summary_tc(psi_energy_with_nucl_rep, pt2_data, pt2_data_err, N_det, N_configuration, N_states, psi_s2)
@ -117,48 +109,19 @@ subroutine run_stochastic_cipsi
PROVIDE psi_det
PROVIDE psi_det_sorted_tc
ept2(N_iter-1) = E_tc + nuclear_repulsion + (pt2_data % pt2(1))/norm
pt1(N_iter-1) = dsqrt(pt2_data % overlap(1,1))
call diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm, pt2_data, print_pt2)
call diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm)
! stop
if (qp_stop()) exit
enddo
! print*,'data to extrapolate '
! do i = 2, N_iter
! print*,'iteration ',i
! print*,'pt1,Ept2',pt1(i),ept2(i)
! call get_extrapolated_energy(i-1,ept2(i),pt1(i),extrap_energy(i))
! do j = 2, i
! print*,'j,e,energy',j,extrap_energy(j)
! enddo
! enddo
! thresh_it_dav = 5.d-6
! soft_touch thresh_it_dav
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,pt2_data,print_pt2)
! if (.not.qp_stop()) then
! if (N_det < N_det_max) then
! thresh_it_dav = 5.d-7
! soft_touch thresh_it_dav
! call diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm,pt2_data,print_pt2)
! endif
!
! 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_denom, pt2_data, pt2_data_err, relative_error, 0) ! Stochastic PT2
! call diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm,pt2_data,print_pt2)
! endif
! call pt2_dealloc(pt2_data)
! call pt2_dealloc(pt2_data_err)
! call routine_save_right
call diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm)
call pt2_dealloc(pt2_data)
call pt2_dealloc(pt2_data_err)
end

View File

@ -1,7 +1,7 @@
! ---
subroutine diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm, pt2_data, print_pt2)
subroutine diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm )
BEGIN_DOC
! Replace the coefficients of the CI states by the coefficients of the
@ -12,50 +12,10 @@ subroutine diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm, pt2_data, print_pt2)
implicit none
integer, intent(inout) :: ndet ! number of determinants from before
double precision, intent(inout) :: E_tc(N_states), norm(N_states) ! E and norm from previous wave function
type(pt2_type) , intent(in) :: pt2_data ! PT2 from previous wave function
logical, intent(in) :: print_pt2
integer :: i, j,k
double precision:: pt2_minus,pt2_plus,pt2_tot, pt2_abs,pt1_norm,rpt2_tot
double precision :: error_pt2_minus, error_pt2_plus, error_pt2_tot, error_pt2_abs
PROVIDE mo_l_coef mo_r_coef
! print*,'*****'
! print*,'New wave function information'
! print*,'N_det tc = ',N_det
! do k = 1, N_states
! print*,'************'
! print*,'State ',k
! pt2_plus = pt2_data % variance(k)
! pt2_minus = pt2_data % pt2(k)
! pt2_abs = pt2_plus - pt2_minus
! pt2_tot = pt2_plus + pt2_minus
!
! pt1_norm = pt2_data % overlap(k,k)
! rpt2_tot = pt2_tot / (1.d0 + pt1_norm)
!
!
! print*,'norm_ground_left_right_bi_orth = ',norm_ground_left_right_bi_orth(k)
! print*,'eigval_right_tc = ',eigval_right_tc_bi_orth(k)
! print*,'*****'
!
! if(print_pt2) then
! print*,'*****'
! print*,'previous wave function info'
! print*,'norm(before) = ',norm
! print*,'E(before) = ',E_tc
! print*,'PT1 norm = ',dsqrt(pt1_norm)
! print*,'PT2 = ',pt2_tot
! print*,'rPT2 = ',rpt2_tot
! print*,'|PT2| = ',pt2_abs
! print*,'Positive PT2 = ',pt2_plus
! print*,'Negative PT2 = ',pt2_minus
! print*,'E(before) + PT2 = ',E_tc + pt2_tot/norm
! print*,'E(before) +rPT2 = ',E_tc + rpt2_tot/norm
! write(*,'(A28,X,I10,X,100(F16.8,X))')'Ndet,E,E+PT2,E+RPT2,|PT2|=',ndet,E_tc ,E_tc + pt2_tot/norm,E_tc + rpt2_tot/norm,pt2_minus, pt2_plus
! print*,'*****'
! endif
! enddo
do k = 1, N_states
E_tc(k) = eigval_right_tc_bi_orth(k)
norm(k) = norm_ground_left_right_bi_orth(k)

View File

@ -11,10 +11,16 @@ BEGIN_PROVIDER [ double precision, psi_average_norm_contrib_tc, (psi_det_size) ]
psi_average_norm_contrib_tc(:) = 0.d0
do k=1,N_states
do i=1,N_det
psi_average_norm_contrib_tc(i) = &
! print*,dabs(psi_l_coef_bi_ortho(i,k)*psi_r_coef_bi_ortho(i,k)),psi_l_coef_bi_ortho(i,k),psi_r_coef_bi_ortho(i,k)
psi_average_norm_contrib_tc(i) += &
dabs(psi_l_coef_bi_ortho(i,k)*psi_r_coef_bi_ortho(i,k))*state_average_weight(k)
enddo
enddo
! print*,'***'
! do i = 1, N_det
! print*,psi_average_norm_contrib_tc(i)
! enddo
print*,'sum(psi_average_norm_contrib_tc(1:N_det))',sum(psi_average_norm_contrib_tc(1:N_det))
f = 1.d0/sum(psi_average_norm_contrib_tc(1:N_det))
do i=1,N_det
psi_average_norm_contrib_tc(i) = psi_average_norm_contrib_tc(i)*f

View File

@ -1,3 +1,3 @@
subroutine provide_for_zmq_pt2
PROVIDE psi_selectors_coef_transp psi_det_sorted psi_det_sorted_order
PROVIDE psi_selectors_coef_transp psi_det_sorted psi_det_sorted_order psi_det_hii
end

View File

@ -306,7 +306,7 @@ subroutine run_slave_main
PROVIDE psi_bilinear_matrix_rows psi_bilinear_matrix_order
PROVIDE psi_bilinear_matrix_transp_rows_loc psi_bilinear_matrix_transp_columns
PROVIDE psi_bilinear_matrix_transp_order psi_selectors_coef_transp
PROVIDE psi_det_hii selection_weight pseudo_sym pt2_min_parallel_tasks
PROVIDE selection_weight pseudo_sym pt2_min_parallel_tasks
call provide_for_zmq_pt2
if (mpi_master) then