10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-07-03 18:05:59 +02:00

debugging

This commit is contained in:
Anthony Scemama 2018-05-22 19:11:10 +02:00
parent 72b0c84cb6
commit 9245c090e2
10 changed files with 46 additions and 34 deletions

View File

@ -63,11 +63,11 @@ subroutine run_dressing(N_st,energy)
enddo
print *, 'Variational energy <Psi|H|Psi>'
do i=1,N_st
print *, i, psi_energy(i)
print *, i, psi_energy(i)+nuclear_repulsion
enddo
print *, 'Dressed energy <Psi|H+Delta|Psi>'
do i=1,N_st
print *, i, ci_energy_dressed(i)
print *, i, ci_energy_dressed(i)+nuclear_repulsion
enddo
endif

View File

@ -59,10 +59,8 @@ integer, external :: zmq_get_dvector, zmq_get_N_det_generators
if (zmq_get_dvector(zmq_to_qp_run_socket,1,'energy',energy,N_states) == -1) cycle
if (zmq_get_dvector(zmq_to_qp_run_socket,1,'dress_stoch_istate',tmp,1) == -1) cycle
dress_stoch_istate = int(tmp)
TOUCH dress_stoch_istate
TOUCH state_average_weight
psi_energy(1:N_states) = energy(1:N_states)
TOUCH psi_energy dress_stoch_istate state_average_weight
PROVIDE psi_bilinear_matrix_columns_loc psi_det_alpha_unique psi_det_beta_unique
PROVIDE psi_bilinear_matrix_rows psi_det_sorted_order psi_bilinear_matrix_order

View File

@ -133,16 +133,18 @@ subroutine ZMQ_dress(E, dress, delta_out, delta_s2_out, relative_error, lndet)
print *, irp_here, ': Failed in zmq_set_running'
endif
!!$OMP PARALLEL DEFAULT(shared) NUM_THREADS(nproc) &
! !$OMP PRIVATE(i)
!i = omp_get_thread_num()
!if (i==0) then
call omp_set_nested(.true.)
!$OMP PARALLEL DEFAULT(shared) NUM_THREADS(2) &
!$OMP PRIVATE(i)
i = omp_get_thread_num()
if (i==0) then
call dress_collector(zmq_socket_pull,E, relative_error, delta, delta_s2, dress,&
dress_stoch_istate)
!else
! call dress_slave_inproc(i)
!endif
!!$OMP END PARALLEL
else
call dress_slave_inproc(i)
endif
!$OMP END PARALLEL
call omp_set_nested(.false.)
delta_out(dress_stoch_istate,1:N_det) = delta(dress_stoch_istate,1:N_det)
delta_s2_out(dress_stoch_istate,1:N_det) = delta_s2(dress_stoch_istate,1:N_det)
call end_parallel_job(zmq_to_qp_run_socket, zmq_socket_pull, 'dress')
@ -237,6 +239,7 @@ subroutine dress_collector(zmq_socket_pull, E, relative_error, delta, delta_s2,
call wall_time(time)
end if
print *, cur_cp, ind
if(cur_cp == -1) then
call dress_pulled(ind, int_buf, double_buf, det_buf, N_buf)
@ -348,7 +351,7 @@ end function
!
! gen_per_cp : number of generators per checkpoint
END_DOC
comb_teeth = 64
comb_teeth = min(1+N_det/10,64)
N_cps_max = 16
gen_per_cp = (N_det_generators / N_cps_max) + 1
END_PROVIDER
@ -505,6 +508,12 @@ END_PROVIDER
done_cp_at_det(dress_jobs(i)) = cur_cp
needed_by_cp(cur_cp) += 1
end do
print *, 'needed_by_cp'
do i=1,cur_cp
print *, i, needed_by_cp(i)
enddo
under_det = 0

View File

@ -65,12 +65,14 @@ END_PROVIDER
BEGIN_PROVIDER [ integer , N_det_delta_ij ]
implicit none
!N_det_delta_ij = 0!N_det
N_det_delta_ij = 1
END_PROVIDER
BEGIN_PROVIDER [ double precision, delta_ij, (N_states, N_det, 2) ]
implicit none
if(.true.) delta_ij(:,:N_det_delta_ij, :) = delta_ij_tmp(:,:,:)
if(.true.) then
delta_ij(:,:N_det_delta_ij, :) = delta_ij_tmp(:,:,:)
endif
delta_ij(:,N_det_delta_ij+1:,:) = 0d0
END_PROVIDER

View File

@ -48,7 +48,7 @@ subroutine run_dress_slave(thread,iproce,energy)
double precision :: fac
if(iproce /= 0) stop "RUN DRESS SLAVE is OMP"
! if(iproce /= 0) stop "RUN DRESS SLAVE is OMP"
allocate(delta_det(N_states, N_det, 0:comb_teeth+1, 2))
allocate(cp(N_states, N_det, N_cp, 2))

View File

@ -581,7 +581,7 @@ END_PROVIDER
double precision, allocatable :: mrcc(:)
double precision :: E_CI_before!, relative_error
double precision, save :: target_error = 0d0
double precision, save :: target_error = 2d-2
allocate(mrcc(N_states))
@ -594,11 +594,10 @@ END_PROVIDER
threshold_selectors = 1.d0
threshold_generators = 1d0
if(target_error /= 0d0) then
target_error = target_error / 2d0 ! (-mrcc_E0_denominator(1) + mrcc_previous_E(1)) / 1d1
target_error = target_error * 0.5d0 ! (-mrcc_E0_denominator(1) + mrcc_previous_E(1)) / 1d1
else
target_error = 1d-4
end if
target_error = 0d0
call ZMQ_mrcc(E_CI_before, mrcc, delta_ij_mrcc_zmq, delta_ij_s2_mrcc_zmq, abs(target_error))
mrcc_previous_E(:) = mrcc_E0_denominator(:)

View File

@ -310,7 +310,6 @@ subroutine mrcc_collector(zmq_socket_pull, E, relative_error, delta, delta_s2, m
end if
end do
if(cur_cp == 0) then
print *, "no checkpoint reached so far..."
cycle pullLoop
end if
!!!!!!!!!!!!

View File

@ -51,6 +51,7 @@ subroutine generator_done(i_gen, int_buf, double_buf, det_buf, N_buf, iproc)
double precision, intent(out) :: double_buf(N_dress_double_buffer)
integer(bit_kind), intent(out) :: det_buf(N_int, 2, N_dress_det_buffer)
integer :: i
int_buf(:) = 0
call sort_selection_buffer(sb(iproc))
det_buf(:,:,:sb(iproc)%cur) = sb(iproc)%det(:,:,:sb(iproc)%cur)
@ -115,15 +116,18 @@ subroutine delta_ij_done()
old_det_gen = N_det_generators
call sort_selection_buffer(global_sb)
call fill_H_apply_buffer_no_selection(global_sb%cur,global_sb%det,N_int,0)
call copy_H_apply_buffer_to_wf()
if (s2_eig.or.(N_states > 1) ) then
call make_s2_eigenfunction
if (dress_stoch_istate == N_states) then
! Add buffer only when the last state is computed
call sort_selection_buffer(global_sb)
call fill_H_apply_buffer_no_selection(global_sb%cur,global_sb%det,N_int,0)
call copy_H_apply_buffer_to_wf()
if (s2_eig.or.(N_states > 1) ) then
call make_s2_eigenfunction
endif
call undress_with_alpha(old_generators, old_det_gen, psi_det(1,1,N_det_delta_ij+1), N_det-N_det_delta_ij)
call save_wavefunction
endif
call undress_with_alpha(old_generators, old_det_gen, psi_det(1,1,N_det_delta_ij+1), N_det-N_det_delta_ij)
call save_wavefunction
end subroutine

View File

@ -192,8 +192,9 @@ subroutine copy_H_apply_buffer_to_wf
call normalize(psi_coef,N_det)
SOFT_TOUCH N_det psi_det psi_coef
! logical :: found_duplicates
! call remove_duplicates_in_psi_det(found_duplicates)
logical :: found_duplicates
call remove_duplicates_in_psi_det(found_duplicates)
end
subroutine remove_duplicates_in_psi_det(found_duplicates)

View File

@ -101,8 +101,8 @@ subroutine map_load_from_disk(filename,map)
k = map % consolidated_idx (i+2_8)
l = map % consolidated_idx (i+1_8)
n_elements = int(k - l, 4)
key_p => map % consolidated_key (l:l+n_elements)
value_p => map % consolidated_value ( l:l+n_elements )
key_p => map % consolidated_key (l:l+n_elements-1)
value_p => map % consolidated_value ( l:l+n_elements-1 )
map % map(i) % key => key_p
map % map(i) % value => value_p
map % map(i) % sorted = .True.