10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-12-24 13:23:41 +01:00

Fixed sBK

This commit is contained in:
Anthony Scemama 2018-09-26 15:53:54 +02:00
parent 7fc89d857e
commit 502fd9d1fd
3 changed files with 27 additions and 20 deletions

View File

@ -15,8 +15,7 @@ END_PROVIDER
integer :: i integer :: i
integer :: e integer :: e
e = elec_num - n_core_orb * 2 e = elec_num - n_core_orb * 2
pt2_n_tasks_max = 1+min((e*(e-1))/2, int(dsqrt(dble(N_det_generators)))/10) pt2_n_tasks_max = 1+min((e*(e-1))/2, int(dsqrt(dble(N_det_selectors)))/10)
pt2_n_tasks_max = 1
do i=1,N_det_generators do i=1,N_det_generators
if (maxval(dabs(psi_coef_sorted_gen(i,1:N_states))) > 0.001d0) then if (maxval(dabs(psi_coef_sorted_gen(i,1:N_states))) > 0.001d0) then
pt2_F(i) = pt2_n_tasks_max pt2_F(i) = pt2_n_tasks_max

View File

@ -30,7 +30,7 @@ END_PROVIDER
pt2_N_teeth = 1 pt2_N_teeth = 1
else else
pt2_minDetInFirstTeeth = min(5, N_det_generators) pt2_minDetInFirstTeeth = min(5, N_det_generators)
do pt2_N_teeth=20,2,-1 do pt2_N_teeth=100,2,-1
if(testTeethBuilding(pt2_minDetInFirstTeeth, pt2_N_teeth)) exit if(testTeethBuilding(pt2_minDetInFirstTeeth, pt2_N_teeth)) exit
end do end do
end if end if
@ -89,7 +89,7 @@ logical function testTeethBuilding(minF, N)
end function end function
BEGIN_PROVIDER[ integer, dress_N_cp_max ] BEGIN_PROVIDER[ integer, dress_N_cp_max ]
dress_N_cp_max = 64 dress_N_cp_max = 28
END_PROVIDER END_PROVIDER
BEGIN_PROVIDER[integer, pt2_J, (N_det_generators)] BEGIN_PROVIDER[integer, pt2_J, (N_det_generators)]
@ -102,6 +102,7 @@ END_PROVIDER
pt2_J = pt2_J_ pt2_J = pt2_J_
dress_R1 = dress_R1_ dress_R1 = dress_R1_
!return
do m=1,dress_N_cp do m=1,dress_N_cp
nmov = 0 nmov = 0
@ -209,6 +210,11 @@ END_PROVIDER
enddo enddo
dress_N_cp = m-1 dress_N_cp = m-1
if (dress_N_cp == 0) then
print *, irp_here, 'dress_N_cp = 0'
stop -1
endif
dress_R1_(dress_N_cp) = N_j dress_R1_(dress_N_cp) = N_j
dress_M_m(dress_N_cp) = N_c dress_M_m(dress_N_cp) = N_c
!!!!!!!!!!!!!! !!!!!!!!!!!!!!
@ -510,6 +516,7 @@ subroutine dress_collector(zmq_socket_pull, E, relative_error, delta, delta_s2,
more = 1 more = 1
do while (.not. found) do while (.not. found)
print *, 'm, dotfm', m, dot_f(m)
if(dot_f(m) == 0) then if(dot_f(m) == 0) then
E0 = 0 E0 = 0
do i=dress_dot_n_0(m),1,-1 do i=dress_dot_n_0(m),1,-1
@ -527,16 +534,17 @@ subroutine dress_collector(zmq_socket_pull, E, relative_error, delta, delta_s2,
end do end do
end do end do
t = dress_dot_t(m) t = dress_dot_t(m)
!print *, 'm dressncp', m, dress_N_cp
avg = E0 + S(t) / dble(c) avg = E0 + S(t) / dble(c)
if (c > 2) then if (c > 2) then
eqt = dabs((S2(t) / c) - (S(t)/c)**2) eqt = dabs((S2(t) / c) - (S(t)/c)**2)
eqt = sqrt(eqt / (dble(c)-1.5d0)) error = sqrt(eqt / (dble(c)-1.5d0))
error = eqt
time = omp_get_wtime() time = omp_get_wtime()
print '(G10.3, 2X, F16.10, 2X, G16.3, 2X, F16.4, A20)', c, avg+E(istate), eqt, time-time0, '' print '(G10.3, 2X, F16.10, 2X, G16.3, 2X, F16.4, A20)', c, avg+E(istate), error, time-time0, ''
else if ( m==dress_N_cp ) then
error = 0.d0
else else
eqt = 1.d0 error =1.d0
error = eqt
endif endif
m += 1 m += 1
if(dabs(error / avg) <= relative_error) then if(dabs(error / avg) <= relative_error) then

View File

@ -36,7 +36,7 @@ subroutine run_dress_slave(thread,iproce,energy)
! integer(kind=OMP_LOCK_KIND) :: lck_det(0:pt2_N_teeth+1) ! integer(kind=OMP_LOCK_KIND) :: lck_det(0:pt2_N_teeth+1)
! integer(kind=OMP_LOCK_KIND) :: lck_sto(dress_N_cp) ! integer(kind=OMP_LOCK_KIND) :: lck_sto(dress_N_cp)
double precision :: fac double precision :: fac
integer :: ending integer :: ending, ending_tmp
integer, external :: zmq_get_dvector, zmq_get_int integer, external :: zmq_get_dvector, zmq_get_int
! double precision, external :: omp_get_wtime ! double precision, external :: omp_get_wtime
double precision :: time, time0 double precision :: time, time0
@ -96,7 +96,8 @@ subroutine run_dress_slave(thread,iproce,energy)
call push_dress_results(zmq_socket_push, 0, 0, edI_task, edI_index, breve_delta_m, task_buf, ntask_buf) call push_dress_results(zmq_socket_push, 0, 0, edI_task, edI_index, breve_delta_m, task_buf, ntask_buf)
end if end if
!$OMP FLUSH m=0
!$OMP BARRIER
do while( (cp_done > cp_sent) .or. (m /= dress_N_cp+1) ) do while( (cp_done > cp_sent) .or. (m /= dress_N_cp+1) )
!$OMP CRITICAL (send) !$OMP CRITICAL (send)
if(ntask_tbd == 0) then if(ntask_tbd == 0) then
@ -116,13 +117,13 @@ subroutine run_dress_slave(thread,iproce,energy)
ntask_tbd -= 1 ntask_tbd -= 1
else else
m = dress_N_cp + 1 m = dress_N_cp + 1
i= zmq_get_int(zmq_to_qp_run_socket, worker_id, "ending", ending) if (zmq_get_int(zmq_to_qp_run_socket, worker_id, "ending", ending_tmp) /= -1) then
ending = ending_tmp
endif
end if end if
will_send = 0 will_send = 0
cp_max(iproc) = m cp_max(iproc) = m
! print *, cp_max(:)
! print *, ''
cp_done = minval(cp_max)-1 cp_done = minval(cp_max)-1
if(cp_done > cp_sent) then if(cp_done > cp_sent) then
will_send = cp_sent + 1 will_send = cp_sent + 1
@ -147,6 +148,7 @@ subroutine run_dress_slave(thread,iproce,energy)
edI_index(n_tasks) = i edI_index(n_tasks) = i
end if end if
end do end do
write(0,*) 'will send', will_send, n_tasks
call push_dress_results(zmq_socket_push, will_send, sum_f, edI_task, edI_index, & call push_dress_results(zmq_socket_push, will_send, sum_f, edI_task, edI_index, &
breve_delta_m, task_buf, n_tasks) breve_delta_m, task_buf, n_tasks)
end if end if
@ -160,7 +162,6 @@ subroutine run_dress_slave(thread,iproce,energy)
time0 = omp_get_wtime() time0 = omp_get_wtime()
call alpha_callback(breve_delta_m, i_generator, subset, pt2_F(i_generator), iproc) call alpha_callback(breve_delta_m, i_generator, subset, pt2_F(i_generator), iproc)
time = omp_get_wtime() time = omp_get_wtime()
!print '(I0.11, I4, A12, F12.3)', i_generator, subset, "GREPMETIME", time-time0
t = dress_T(i_generator) t = dress_T(i_generator)
!$OMP CRITICAL(t_crit) !$OMP CRITICAL(t_crit)
@ -200,7 +201,6 @@ subroutine run_dress_slave(thread,iproce,energy)
ntask_buf = 0 ntask_buf = 0
end if end if
end if end if
!$OMP FLUSH
end do end do
!$OMP BARRIER !$OMP BARRIER
@ -208,11 +208,11 @@ subroutine run_dress_slave(thread,iproce,energy)
call push_dress_results(zmq_socket_push, 0, 0, edI_task, edI_index, breve_delta_m, task_buf, ntask_buf) call push_dress_results(zmq_socket_push, 0, 0, edI_task, edI_index, breve_delta_m, task_buf, ntask_buf)
ntask_buf = 0 ntask_buf = 0
end if end if
!$OMP SINGLE
if(purge_task_id /= 0) then !$OMP SINGLE
do while(ending == dress_N_cp+1) if(purge_task_id /= 0) then
do while (zmq_get_int(zmq_to_qp_run_socket, worker_id, "ending", ending) == -1)
call sleep(1) call sleep(1)
i= zmq_get_int(zmq_to_qp_run_socket, worker_id, "ending", ending)
end do end do
will_send = ending will_send = ending