10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-09-27 03:51:01 +02:00

Fixed fci

This commit is contained in:
Anthony Scemama 2019-01-08 19:02:04 +01:00
parent bca983cd84
commit a82cbdb22e
3 changed files with 16 additions and 19 deletions

View File

@ -257,8 +257,6 @@ subroutine ZMQ_pt2(E, pt2,relative_error, error, variance, norm, N_in)
call make_selection_buffer_s2(b) call make_selection_buffer_s2(b)
endif endif
call fill_H_apply_buffer_no_selection(b%cur,b%det,N_int,0) call fill_H_apply_buffer_no_selection(b%cur,b%det,N_int,0)
call copy_H_apply_buffer_to_wf()
call save_wavefunction
endif endif
call delete_selection_buffer(b) call delete_selection_buffer(b)

View File

@ -5,7 +5,7 @@ subroutine run_stochastic_cipsi
END_DOC END_DOC
integer :: i,j,k integer :: i,j,k
double precision, allocatable :: pt2(:), variance(:), norm(:), rpt2(:) double precision, allocatable :: pt2(:), variance(:), norm(:), rpt2(:)
integer :: N_det_before, N_occ_pattern_before, to_select integer :: to_select
double precision :: rss double precision :: rss
double precision, external :: memory_of_double double precision, external :: memory_of_double
@ -52,9 +52,6 @@ subroutine run_stochastic_cipsi
call save_wavefunction call save_wavefunction
endif endif
N_det_before = 0
N_occ_pattern_before = 0
double precision :: correlation_energy_ratio double precision :: correlation_energy_ratio
double precision :: error(N_states) double precision :: error(N_states)
@ -68,8 +65,6 @@ subroutine run_stochastic_cipsi
write(*,'(A)') '--------------------------------------------------------------------------------' write(*,'(A)') '--------------------------------------------------------------------------------'
N_det_before = N_det
N_occ_pattern_before = N_occ_pattern
to_select = N_det to_select = N_det
to_select = max(N_states_diag, to_select) to_select = max(N_states_diag, to_select)
@ -77,7 +72,7 @@ subroutine run_stochastic_cipsi
variance = 0.d0 variance = 0.d0
norm = 0.d0 norm = 0.d0
call ZMQ_pt2(psi_energy_with_nucl_rep,pt2,relative_error,error, variance, & call ZMQ_pt2(psi_energy_with_nucl_rep,pt2,relative_error,error, variance, &
norm, to_select) ! Stochastic PT2 norm, to_select) ! Stochastic PT2 and selection
correlation_energy_ratio = (psi_energy_with_nucl_rep(1) - hf_energy_ref) / & correlation_energy_ratio = (psi_energy_with_nucl_rep(1) - hf_energy_ref) / &
(psi_energy_with_nucl_rep(1) + pt2(1) - hf_energy_ref) (psi_energy_with_nucl_rep(1) + pt2(1) - hf_energy_ref)
@ -85,16 +80,20 @@ subroutine run_stochastic_cipsi
call save_energy(psi_energy_with_nucl_rep, pt2) call save_energy(psi_energy_with_nucl_rep, pt2)
call write_double(6,correlation_energy_ratio, 'Correlation ratio') call write_double(6,correlation_energy_ratio, 'Correlation ratio')
call print_summary(psi_energy_with_nucl_rep(1:N_states),pt2,error,variance,norm,N_det_before,N_occ_pattern_before) call print_summary(psi_energy_with_nucl_rep,pt2,error,variance,norm,N_det,N_occ_pattern,N_states)
do k=1,N_states do k=1,N_states
rpt2(:) = pt2(:)/(1.d0 + norm(k)) rpt2(:) = pt2(:)/(1.d0 + norm(k))
enddo enddo
call save_iterations(psi_energy_with_nucl_rep(1:N_states),rpt2,N_det_before) call save_iterations(psi_energy_with_nucl_rep(1:N_states),rpt2,N_det)
call print_extrapolated_energy() call print_extrapolated_energy()
N_iter += 1 N_iter += 1
! Add selected determinants
call copy_H_apply_buffer_to_wf()
call save_wavefunction
PROVIDE psi_coef PROVIDE psi_coef
PROVIDE psi_det PROVIDE psi_det
PROVIDE psi_det_sorted PROVIDE psi_det_sorted
@ -123,7 +122,7 @@ subroutine run_stochastic_cipsi
rpt2(:) = pt2(:)/(1.d0 + norm(k)) rpt2(:) = pt2(:)/(1.d0 + norm(k))
enddo enddo
call print_summary(psi_energy_with_nucl_rep(1:N_states),pt2,error,variance,norm,N_det,N_occ_pattern) call print_summary(psi_energy_with_nucl_rep(1:N_states),pt2,error,variance,norm,N_det,N_occ_pattern,N_states)
call save_iterations(psi_energy_with_nucl_rep(1:N_states),rpt2,N_det) call save_iterations(psi_energy_with_nucl_rep(1:N_states),rpt2,N_det)
call print_extrapolated_energy() call print_extrapolated_energy()

View File

@ -1,16 +1,16 @@
subroutine print_summary(e_,pt2_,error_,variance_,norm_,n_det_,n_occ_pattern_) subroutine print_summary(e_,pt2_,error_,variance_,norm_,n_det_,n_occ_pattern_,n_st)
implicit none implicit none
BEGIN_DOC BEGIN_DOC
! Print the extrapolated energy in the output ! Print the extrapolated energy in the output
END_DOC END_DOC
double precision, intent(in) :: e_(N_states), pt2_(N_states), variance_(N_states), norm_(N_states), error_(N_states) integer, intent(in) :: n_det_, n_occ_pattern_, n_st
integer, intent(in) :: n_det_, n_occ_pattern_ double precision, intent(in) :: e_(n_st), pt2_(n_st), variance_(n_st), norm_(n_st), error_(n_st)
integer :: i, k integer :: i, k
integer :: N_states_p integer :: N_states_p
character*(9) :: pt2_string character*(9) :: pt2_string
character*(512) :: fmt character*(512) :: fmt
double precision :: f(N_states) double precision :: f(n_st)
if (do_pt2) then if (do_pt2) then
pt2_string = ' ' pt2_string = ' '
@ -18,7 +18,7 @@ subroutine print_summary(e_,pt2_,error_,variance_,norm_,n_det_,n_occ_pattern_)
pt2_string = '(approx)' pt2_string = '(approx)'
endif endif
N_states_p = min(N_det_,N_states) N_states_p = min(N_det_,n_st)
do i=1,N_states_p do i=1,N_states_p
f(i) = 1.d0/(1.d0+norm_(i)) f(i) = 1.d0/(1.d0+norm_(i))
@ -57,7 +57,7 @@ subroutine print_summary(e_,pt2_,error_,variance_,norm_,n_det_,n_occ_pattern_)
print *, '' print *, ''
print *, 'N_det = ', N_det_ print *, 'N_det = ', N_det_
print *, 'N_states = ', N_states print *, 'N_states = ', n_st
if (s2_eig) then if (s2_eig) then
print *, 'N_sop = ', N_occ_pattern_ print *, 'N_sop = ', N_occ_pattern_
endif endif
@ -76,7 +76,7 @@ subroutine print_summary(e_,pt2_,error_,variance_,norm_,n_det_,n_occ_pattern_)
enddo enddo
print *, '-----' print *, '-----'
if(N_states.gt.1)then if(n_st.gt.1)then
print *, 'Variational Energy difference (au | eV)' print *, 'Variational Energy difference (au | eV)'
do i=2, N_states_p do i=2, N_states_p
print*,'Delta E = ', (e_(i) - e_(1)), & print*,'Delta E = ', (e_(i) - e_(1)), &