10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-07-03 18:05:59 +02:00
quantum_package/plugins/dress_zmq/run_dress_slave.irp.f
Yann Garniron 82d37fa2b2 Merge branch 'master' of https://github.com/scemama/quantum_package into stashed
Conflicts:
	plugins/mrcepa0/mrcc_stoch_routines.irp.f
2018-01-09 12:20:19 +01:00

182 lines
5.7 KiB
Fortran

BEGIN_PROVIDER [ integer, fragment_count ]
implicit none
BEGIN_DOC
! Number of fragments for the deterministic part
END_DOC
fragment_count = 1
END_PROVIDER
subroutine run_dress_slave(thread,iproc,energy)
use f77_zmq
implicit none
double precision, intent(in) :: energy(N_states_diag)
integer, intent(in) :: thread, iproc
integer :: rc, i, subset, i_generator
integer :: worker_id, task_id, ctask, ltask
character*(512) :: task
integer(ZMQ_PTR),external :: new_zmq_to_qp_run_socket
integer(ZMQ_PTR) :: zmq_to_qp_run_socket
integer(ZMQ_PTR), external :: new_zmq_push_socket
integer(ZMQ_PTR) :: zmq_socket_push
logical :: done
double precision,allocatable :: dress_detail(:)
integer :: ind
integer(bit_kind),allocatable :: abuf(:,:,:)
integer(bit_kind) :: mask(N_int,2), omask(N_int,2)
double precision,allocatable :: delta_ij_loc(:,:,:)
double precision,allocatable :: delta_ii_loc(:,:)
integer :: h,p,n
logical :: ok
double precision :: contrib(N_states)
allocate(delta_ij_loc(N_states,N_det_non_ref,2) &
,delta_ii_loc(N_states,2))
allocate(abuf(N_int, 2, N_det_non_ref))
allocate(dress_detail(N_states))
zmq_to_qp_run_socket = new_zmq_to_qp_run_socket()
zmq_socket_push = new_zmq_push_socket(thread)
call connect_to_taskserver(zmq_to_qp_run_socket,worker_id,thread)
if(worker_id == -1) then
print *, "WORKER -1"
call end_zmq_to_qp_run_socket(zmq_to_qp_run_socket)
call end_zmq_push_socket(zmq_socket_push,thread)
return
end if
dress_detail = 0d0
do
call get_task_from_taskserver(zmq_to_qp_run_socket,worker_id, task_id, task)
if(task_id /= 0) then
read (task,*) subset, i_generator
contrib = 0d0
delta_ij_loc = 0d0
delta_ii_loc = 0d0
if(do_dress_with_alpha_buffer .or. do_dress_with_alpha) then
do h=1, hh_shortcut(0)
call apply_hole_local(psi_det_generators(1,1,i_generator), hh_exists(1, h), mask, ok, N_int)
if(.not. ok) cycle
omask = 0_bit_kind
if(hh_exists(1, h) /= 0) omask = mask
n = 1
do p=hh_shortcut(h), hh_shortcut(h+1)-1
call apply_particle_local(mask, pp_exists(1, p), abuf(1,1,n), ok, N_int)
if(ok) n = n + 1
if(n > N_det_non_ref) stop "Buffer too small in dress..."
end do
n = n - 1
if(n /= 0) then
if(do_dress_with_alpha_buffer) then
call dress_with_alpha_buffer(delta_ij_loc(1,1,1), delta_ii_loc(1,1), delta_ij_loc(1,1,2), delta_ii_loc(1,2), &
i_generator,n,abuf,N_int,omask,contrib)
else
stop 'dress_with_alpha not implemented yet'
end if
endif
end do
else if(do_dress_with_generator) then
stop 'dress_with_generator not implemented yet'
else
stop 'no dressing level defined'
end if
dress_detail(:) = contrib
call task_done_to_taskserver(zmq_to_qp_run_socket,worker_id,task_id)
call push_dress_results(zmq_socket_push, i_generator, dress_detail, delta_ij_loc(1,1,1), task_id)
dress_detail(:) = 0d0
else
exit
end if
end do
call disconnect_from_taskserver(zmq_to_qp_run_socket,zmq_socket_push,worker_id)
call end_zmq_to_qp_run_socket(zmq_to_qp_run_socket)
call end_zmq_push_socket(zmq_socket_push,thread)
end subroutine
subroutine push_dress_results(zmq_socket_push, ind, dress_detail, delta_loc, task_id)
use f77_zmq
implicit none
integer(ZMQ_PTR), intent(in) :: zmq_socket_push
double precision, intent(in) :: dress_detail(N_states, N_det_generators)
double precision, intent(in) :: delta_loc(N_states, N_det_non_ref, 2)
integer, intent(in) :: ind, task_id
integer :: rc, i
rc = f77_zmq_send( zmq_socket_push, ind, 4, ZMQ_SNDMORE)
if(rc /= 4) stop "push"
rc = f77_zmq_send( zmq_socket_push, dress_detail, 8*N_states, ZMQ_SNDMORE)
if(rc /= 8*N_states) stop "push"
rc = f77_zmq_send( zmq_socket_push, delta_loc(1,1,1), 8*N_states*N_det_non_ref, ZMQ_SNDMORE)
if(rc /= 8*N_states*N_det_non_ref) stop "push"
rc = f77_zmq_send( zmq_socket_push, delta_loc(1,1,2), 8*N_states*N_det_non_ref, ZMQ_SNDMORE)
if(rc /= 8*N_states*N_det_non_ref) stop "push"
rc = f77_zmq_send( zmq_socket_push, task_id, 4, 0)
if(rc /= 4) stop "push"
! Activate is zmq_socket_push is a REQ
IRP_IF ZMQ_PUSH
IRP_ELSE
character*(2) :: ok
rc = f77_zmq_recv( zmq_socket_push, ok, 2, 0)
IRP_ENDIF
end subroutine
subroutine pull_dress_results(zmq_socket_pull, ind, dress_detail, delta_loc, task_id)
use f77_zmq
implicit none
integer(ZMQ_PTR), intent(in) :: zmq_socket_pull
double precision, intent(inout) :: dress_detail(N_states)
double precision, intent(inout) :: delta_loc(N_states, N_det_non_ref, 2)
double precision, allocatable :: dress_dress_mwen(:,:)
integer, intent(out) :: ind
integer, intent(out) :: task_id
integer :: rc, rn, i
allocate(dress_dress_mwen(N_states, N_det_non_ref))
rc = f77_zmq_recv( zmq_socket_pull, ind, 4, 0)
if(rc /= 4) stop "pull"
rc = f77_zmq_recv( zmq_socket_pull, dress_detail, N_states*8, 0)
if(rc /= 8*N_states) stop "pull"
rc = f77_zmq_recv( zmq_socket_pull, delta_loc(1,1,1), N_states*8*N_det_non_ref, 0)
if(rc /= 8*N_states*N_det_non_ref) stop "pull"
rc = f77_zmq_recv( zmq_socket_pull, delta_loc(1,1,2), N_states*8*N_det_non_ref, 0)
if(rc /= 8*N_states*N_det_non_ref) stop "pull"
rc = f77_zmq_recv( zmq_socket_pull, task_id, 4, 0)
if(rc /= 4) stop "pull"
! Activate is zmq_socket_pull is a REP
IRP_IF ZMQ_PUSH
IRP_ELSE
rc = f77_zmq_send( zmq_socket_pull, 'ok', 2, 0)
IRP_ENDIF
end subroutine