10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-08-16 01:28:39 +02:00

bias at long runtime

This commit is contained in:
Yann Garniron 2017-01-04 13:41:37 +01:00
parent e1d014aa41
commit d1e52144d3

View File

@ -204,34 +204,46 @@ subroutine ZMQ_pt2(pt2)
integer :: tbc(0:N_det_generators) integer :: tbc(0:N_det_generators)
integer :: i, Ncomb, generator_per_task, i_generator_end integer :: i, Ncomb, generator_per_task, i_generator_end
integer, external :: pt2_find integer, external :: pt2_find
double precision :: sumabove(comb_teeth), sum2above(comb_teeth), Nabove(comb_teeth)
double precision, external :: omp_get_wtime
double precision :: time0, time
provide nproc provide nproc
call new_parallel_job(zmq_to_qp_run_socket,"pt2")
call zmq_put_psi(zmq_to_qp_run_socket,1,ci_electronic_energy,size(ci_electronic_energy))
call zmq_set_running(zmq_to_qp_run_socket)
call create_selection_buffer(1, 1*2, b)
call random_seed() call random_seed()
computed = .false. computed = .false.
tbc(0) = first_det_of_comb - 1 tbc(0) = first_det_of_comb - 1
do i=1, tbc(0) do i=1, tbc(0)
tbc(i) = i tbc(i) = i
computed(i) = .true. computed(i) = .true.
end do end do
print *, "detererminist initial ", tbc(0)+1 pt2_detail = 0d0
!LOOP?
time0 = omp_get_wtime()
print *, "grep - time - avg - err - n_combs"
do while(.true.)
call new_parallel_job(zmq_to_qp_run_socket,"pt2")
call zmq_put_psi(zmq_to_qp_run_socket,1,ci_electronic_energy,size(ci_electronic_energy))
call zmq_set_running(zmq_to_qp_run_socket)
call create_selection_buffer(1, 1*2, b)
call get_carlo_workbatch(1d-3, computed, comb, Ncomb, tbc) call get_carlo_workbatch(1d-3, computed, comb, Ncomb, tbc)
generator_per_task = tbc(0)/1000 + 1 generator_per_task = tbc(0)/1000 + 1
print *, "TASK", tbc(0), tbc(1:tbc(0))
do i=1,tbc(0),generator_per_task do i=1,tbc(0),generator_per_task
i_generator_end = min(i+generator_per_task-1, tbc(0)) i_generator_end = min(i+generator_per_task-1, tbc(0))
!print *, "TASK", (i_generator_end-i+1), tbc(i:i_generator_end) !print *, "TASK", (i_generator_end-i+1), tbc(i:i_generator_end)
write(task,*) (i_generator_end-i+1), tbc(i:i_generator_end) write(task,*) (i_generator_end-i+1), tbc(i:i_generator_end)
call add_task_to_taskserver(zmq_to_qp_run_socket,task) call add_task_to_taskserver(zmq_to_qp_run_socket,task)
end do end do
pt2_detail = 0d0
!$OMP PARALLEL DEFAULT(shared) SHARED(b, pt2) PRIVATE(i) NUM_THREADS(nproc+1) !$OMP PARALLEL DEFAULT(shared) SHARED(b, pt2) PRIVATE(i) NUM_THREADS(nproc+1)
i = omp_get_thread_num() i = omp_get_thread_num()
@ -243,17 +255,37 @@ subroutine ZMQ_pt2(pt2)
!$OMP END PARALLEL !$OMP END PARALLEL
call end_parallel_job(zmq_to_qp_run_socket, 'pt2') call end_parallel_job(zmq_to_qp_run_socket, 'pt2')
double precision :: E0, avg, eqt double precision :: E0, avg, eqt
double precision :: sumabove(comb_teeth), sum2above(comb_teeth), Nabove(comb_teeth)
call do_carlo(tbc, Ncomb, comb, pt2_detail, sumabove, sum2above, Nabove) call do_carlo(tbc, Ncomb, comb, pt2_detail, sumabove, sum2above, Nabove)
tbc(0) = 0
!END LOOP? !END LOOP?
integer :: tooth integer :: tooth
!-8.091550677158776E-003 !-8.091550677158776E-003
call get_first_tooth(computed, tooth) call get_first_tooth(computed, tooth)
E0 = sum(pt2_detail(1,:first_det_of_teeth(tooth)-1)) print *, "TOOTH ", tooth
avg = E0 + (sumabove(tooth) / Nabove(tooth))
eqt = sqrt(1d0 / (Nabove(tooth)-1) * abs(sum2above(tooth) / Nabove(tooth) - (sumabove(tooth)/Nabove(tooth))**2)) !!! ASSERT
print *, "PT2 ", avg, "+/-", eqt do i=1,first_det_of_teeth(tooth)-1
if(not(computed(i))) stop "deter non calc"
end do
logical :: ok
ok = .false.
do i=first_det_of_teeth(tooth), first_det_of_teeth(tooth+1)-1
if(not(computed(i))) ok = .true.
end do
if(not(ok)) stop "not OK..."
!!!!!
if(Nabove(tooth) >= 30) then
E0 = sum(pt2_detail(1,:first_det_of_teeth(tooth)-1))
avg = E0 + (sumabove(tooth) / Nabove(tooth))
eqt = sqrt(1d0 / (Nabove(tooth)-1) * abs(sum2above(tooth) / Nabove(tooth) - (sumabove(tooth)/Nabove(tooth))**2))
time = omp_get_wtime()
print "(A, 5(E15.7))", "PT2stoch ", time - time0, avg, eqt, Nabove(tooth)
else
print *, Nabove(tooth), "< 30 combs"
end if
tbc(0) = 0
end do
pt2 = 0d0 pt2 = 0d0
end subroutine end subroutine
@ -271,7 +303,7 @@ subroutine do_carlo(tbc, Ncomb, comb, pt2_detail, sumabove, sum2above, Nabove)
myVal = 0d0 myVal = 0d0
myVal2 = 0d0 myVal2 = 0d0
do j=comb_teeth,1,-1 do j=comb_teeth,1,-1
if(pt2_detail(1, dets(j)) == -1d0) print *, "uncalculatedidified", dets(j), pt2_detail(1, dets(j)-1:dets(j)+1) !if(pt2_detail(1, dets(j)) == -1d0) print *, "uncalculatedidified", dets(j), pt2_detail(1, dets(j)-1:dets(j)+1)
myVal += pt2_detail(1, dets(j)) / weight(dets(j)) * comb_step myVal += pt2_detail(1, dets(j)) / weight(dets(j)) * comb_step
sumabove(j) += myVal sumabove(j) += myVal
sum2above(j) += myVal**2 sum2above(j) += myVal**2
@ -354,7 +386,7 @@ subroutine pt2_collector(b, pt2_detail)
type(selection_buffer), intent(inout) :: b type(selection_buffer), intent(inout) :: b
double precision, intent(out) :: pt2_detail(N_states, N_det) double precision, intent(inout) :: pt2_detail(N_states, N_det)
double precision :: pt2_mwen(N_states, N_det) double precision :: pt2_mwen(N_states, N_det)
integer(ZMQ_PTR),external :: new_zmq_to_qp_run_socket integer(ZMQ_PTR),external :: new_zmq_to_qp_run_socket
integer(ZMQ_PTR) :: zmq_to_qp_run_socket integer(ZMQ_PTR) :: zmq_to_qp_run_socket
@ -375,7 +407,7 @@ subroutine pt2_collector(b, pt2_detail)
allocate(val(b%N), det(N_int, 2, b%N), task_id(N_det), index(N_det)) allocate(val(b%N), det(N_int, 2, b%N), task_id(N_det), index(N_det))
done = 0 done = 0
more = 1 more = 1
pt2_detail = -1d0 !pt2_detail = 0d0
call CPU_TIME(time0) call CPU_TIME(time0)
do while (more == 1) do while (more == 1)
call pull_pt2_results(zmq_socket_pull, Nindex, index, pt2_mwen, task_id, ntask) call pull_pt2_results(zmq_socket_pull, Nindex, index, pt2_mwen, task_id, ntask)
@ -483,7 +515,7 @@ end function
BEGIN_PROVIDER [ integer, comb_teeth ] BEGIN_PROVIDER [ integer, comb_teeth ]
implicit none implicit none
comb_teeth = 10 comb_teeth = 20
END_PROVIDER END_PROVIDER
@ -503,13 +535,14 @@ subroutine get_first_tooth(computed, first_teeth)
end do end do
do i=comb_teeth, 1, -1 do i=comb_teeth, 1, -1
if(first_det_of_teeth(i) < first_teeth) then if(first_det_of_teeth(i) <= first_teeth) then
first_teeth = i first_teeth = i
exit exit
end if end if
end do end do
end subroutine end subroutine
subroutine get_carlo_workbatch(maxWorkload, computed, comb, Ncomb, tbc) subroutine get_carlo_workbatch(maxWorkload, computed, comb, Ncomb, tbc)
implicit none implicit none
double precision, intent(in) :: maxWorkload double precision, intent(in) :: maxWorkload
@ -574,7 +607,7 @@ end subroutine
&BEGIN_PROVIDER [ double precision, cweight, (N_det_generators) ] &BEGIN_PROVIDER [ double precision, cweight, (N_det_generators) ]
&BEGIN_PROVIDER [ double precision, comb_workload, (N_det_generators) ] &BEGIN_PROVIDER [ double precision, comb_workload, (N_det_generators) ]
&BEGIN_PROVIDER [ double precision, comb_step ] &BEGIN_PROVIDER [ double precision, comb_step ]
&BEGIN_PROVIDER [ integer, first_det_of_teeth, (comb_teeth) ] &BEGIN_PROVIDER [ integer, first_det_of_teeth, (comb_teeth+1) ]
&BEGIN_PROVIDER [ integer, first_det_of_comb ] &BEGIN_PROVIDER [ integer, first_det_of_comb ]
implicit none implicit none
integer :: i integer :: i
@ -597,7 +630,7 @@ end subroutine
comb_step = 1d0/dfloat(comb_teeth) comb_step = 1d0/dfloat(comb_teeth)
do i=1,N_det_generators do i=1,N_det_generators
if(weight(i)/norm_left < comb_step/2d0) then if(weight(i)/norm_left < comb_step/1d1) then
first_det_of_comb = i first_det_of_comb = i
exit exit
end if end if
@ -611,8 +644,9 @@ end subroutine
first_det_of_teeth(i) = pt2_find(stato, cweight) first_det_of_teeth(i) = pt2_find(stato, cweight)
stato -= comb_step stato -= comb_step
end do end do
print *, first_det_of_teeth(1), first_det_of_comb first_det_of_teeth(comb_teeth+1) = N_det_generators + 1
if(first_det_of_teeth(1) /= first_det_of_comb) stop "comb provider" if(first_det_of_teeth(1) /= first_det_of_comb) stop "comb provider"
END_PROVIDER END_PROVIDER