10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-12-22 20:35:19 +01:00
quantum_package/plugins/Full_CI_ZMQ/run_pt2_slave.irp.f

202 lines
5.1 KiB
Fortran
Raw Normal View History

2017-01-12 15:46:10 +01:00
subroutine run_pt2_slave(thread,iproc,energy)
use f77_zmq
use selection_types
implicit none
double precision, intent(in) :: energy(N_states_diag)
integer, intent(in) :: thread, iproc
integer :: rc, i
2017-11-27 19:44:29 +01:00
integer :: worker_id, ctask, ltask
character*(512), allocatable :: task(:)
integer, allocatable :: task_id(:)
2017-01-12 15:46:10 +01:00
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
2017-05-05 15:08:51 +02:00
type(selection_buffer) :: buf
2017-01-12 15:46:10 +01:00
logical :: done
2017-11-27 19:44:29 +01:00
double precision,allocatable :: pt2(:,:)
integer :: n_tasks, k, n_tasks_max
integer, allocatable :: i_generator(:), subset(:)
2017-01-12 15:46:10 +01:00
2017-11-27 19:44:29 +01:00
n_tasks_max = N_det_generators/100+1
allocate(task_id(n_tasks_max), task(n_tasks_max))
allocate(pt2(N_states,n_tasks_max), i_generator(n_tasks_max), subset(n_tasks_max))
2017-11-29 15:15:10 +01:00
2017-01-12 15:46:10 +01:00
zmq_to_qp_run_socket = new_zmq_to_qp_run_socket()
2017-11-29 15:15:10 +01:00
integer, external :: connect_to_taskserver
if (connect_to_taskserver(zmq_to_qp_run_socket,worker_id,thread) == -1) then
2017-01-12 15:46:10 +01:00
call end_zmq_to_qp_run_socket(zmq_to_qp_run_socket)
return
2017-11-29 15:15:10 +01:00
endif
zmq_socket_push = new_zmq_push_socket(thread)
2017-01-12 15:46:10 +01:00
buf%N = 0
2017-11-27 23:18:18 +01:00
n_tasks = 0
2017-11-27 19:44:29 +01:00
call create_selection_buffer(1, 2, buf)
done = .False.
do while (.not.done)
n_tasks = min(n_tasks+1,n_tasks_max)
2017-01-12 15:46:10 +01:00
2017-11-29 15:15:10 +01:00
integer, external :: get_tasks_from_taskserver
if (get_tasks_from_taskserver(zmq_to_qp_run_socket,worker_id, task_id, task, n_tasks) == -1) then
exit
endif
2017-11-27 19:44:29 +01:00
done = task_id(n_tasks) == 0
if (done) n_tasks = n_tasks-1
2017-11-27 23:18:18 +01:00
if (n_tasks == 0) exit
2017-11-27 19:44:29 +01:00
do k=1,n_tasks
read (task(k),*) subset(k), i_generator(k)
enddo
do k=1,n_tasks
pt2(:,k) = 0.d0
buf%cur = 0
call select_connected(i_generator(k),energy,pt2(1,k),buf,subset(k))
enddo
2017-11-29 15:15:10 +01:00
integer, external :: tasks_done_to_taskserver
if (tasks_done_to_taskserver(zmq_to_qp_run_socket,worker_id,task_id,n_tasks) == -1) then
done = .true.
endif
2017-11-27 19:44:29 +01:00
call push_pt2_results(zmq_socket_push, i_generator, pt2, task_id, n_tasks)
2017-01-12 15:46:10 +01:00
end do
2017-11-29 15:15:10 +01:00
integer, external :: disconnect_from_taskserver
2017-12-05 18:54:10 +01:00
if (disconnect_from_taskserver(zmq_to_qp_run_socket,worker_id) == -1) then
2017-11-29 15:15:10 +01:00
continue
endif
2017-01-12 15:46:10 +01:00
call end_zmq_to_qp_run_socket(zmq_to_qp_run_socket)
call end_zmq_push_socket(zmq_socket_push,thread)
2017-11-27 19:44:29 +01:00
call delete_selection_buffer(buf)
2017-01-12 15:46:10 +01:00
end subroutine
2017-11-27 19:44:29 +01:00
subroutine push_pt2_results(zmq_socket_push, index, pt2, task_id, n_tasks)
2017-01-12 15:46:10 +01:00
use f77_zmq
use selection_types
implicit none
integer(ZMQ_PTR), intent(in) :: zmq_socket_push
2017-11-27 19:44:29 +01:00
double precision, intent(in) :: pt2(N_states,n_tasks)
integer, intent(in) :: n_tasks, index(n_tasks), task_id(n_tasks)
2017-01-12 15:46:10 +01:00
integer :: rc
2017-11-27 19:44:29 +01:00
rc = f77_zmq_send( zmq_socket_push, n_tasks, 4, ZMQ_SNDMORE)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
return
endif
2017-11-24 01:02:53 +01:00
if(rc /= 4) stop 'push'
2017-01-12 15:46:10 +01:00
2017-11-27 19:44:29 +01:00
rc = f77_zmq_send( zmq_socket_push, index, 4*n_tasks, ZMQ_SNDMORE)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
return
endif
2017-11-27 19:44:29 +01:00
if(rc /= 4*n_tasks) stop 'push'
2017-01-12 15:46:10 +01:00
2017-11-27 19:44:29 +01:00
rc = f77_zmq_send( zmq_socket_push, pt2, 8*N_states*n_tasks, ZMQ_SNDMORE)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
return
endif
2017-11-27 19:44:29 +01:00
if(rc /= 8*N_states*n_tasks) stop 'push'
2017-01-12 15:46:10 +01:00
2017-11-27 19:44:29 +01:00
rc = f77_zmq_send( zmq_socket_push, task_id, n_tasks*4, 0)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
return
endif
2017-11-27 19:44:29 +01:00
if(rc /= 4*n_tasks) stop 'push'
2017-01-12 15:46:10 +01:00
! Activate is zmq_socket_push is a REQ
2017-05-16 16:31:35 +02:00
IRP_IF ZMQ_PUSH
IRP_ELSE
character*(2) :: ok
rc = f77_zmq_recv( zmq_socket_push, ok, 2, 0)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
return
endif
2017-11-24 01:02:53 +01:00
if ((rc /= 2).and.(ok(1:2) /= 'ok')) then
print *, irp_here//': error in receiving ok'
stop -1
endif
2017-05-16 16:31:35 +02:00
IRP_ENDIF
2017-01-12 15:46:10 +01:00
end subroutine
2017-11-27 19:44:29 +01:00
subroutine pull_pt2_results(zmq_socket_pull, index, pt2, task_id, n_tasks)
2017-01-12 15:46:10 +01:00
use f77_zmq
use selection_types
implicit none
integer(ZMQ_PTR), intent(in) :: zmq_socket_pull
2017-11-27 19:44:29 +01:00
double precision, intent(inout) :: pt2(N_states,*)
integer, intent(out) :: index(*)
integer, intent(out) :: n_tasks, task_id(*)
2017-01-12 15:46:10 +01:00
integer :: rc, rn, i
2017-11-27 19:44:29 +01:00
rc = f77_zmq_recv( zmq_socket_pull, n_tasks, 4, 0)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
n_tasks = 1
task_id(1) = 0
endif
2017-11-24 01:02:53 +01:00
if(rc /= 4) stop 'pull'
2017-01-12 15:46:10 +01:00
2017-11-27 19:44:29 +01:00
rc = f77_zmq_recv( zmq_socket_pull, index, 4*n_tasks, 0)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
n_tasks = 1
task_id(1) = 0
endif
2017-11-27 19:44:29 +01:00
if(rc /= 4*n_tasks) stop 'pull'
2017-01-12 15:46:10 +01:00
2017-11-27 19:44:29 +01:00
rc = f77_zmq_recv( zmq_socket_pull, pt2, N_states*8*n_tasks, 0)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
n_tasks = 1
task_id(1) = 0
endif
2017-11-27 19:44:29 +01:00
if(rc /= 8*N_states*n_tasks) stop 'pull'
2017-01-12 15:46:10 +01:00
2017-11-27 19:44:29 +01:00
rc = f77_zmq_recv( zmq_socket_pull, task_id, n_tasks*4, 0)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
n_tasks = 1
task_id(1) = 0
endif
2017-11-27 19:44:29 +01:00
if(rc /= 4*n_tasks) stop 'pull'
2017-01-12 15:46:10 +01:00
! Activate is zmq_socket_pull is a REP
2017-05-16 16:31:35 +02:00
IRP_IF ZMQ_PUSH
IRP_ELSE
rc = f77_zmq_send( zmq_socket_pull, 'ok', 2, 0)
2017-11-27 23:18:18 +01:00
if (rc == -1) then
n_tasks = 1
task_id(1) = 0
endif
2017-11-24 01:02:53 +01:00
if (rc /= 2) then
print *, irp_here//': error in sending ok'
stop -1
endif
2017-05-16 16:31:35 +02:00
IRP_ENDIF
2017-01-12 15:46:10 +01:00
end subroutine
2017-04-12 18:26:57 +02:00
BEGIN_PROVIDER [ double precision, pt2_workload, (N_det_generators) ]
2017-01-12 15:46:10 +01:00
integer :: i
2017-04-12 18:26:57 +02:00
do i=1,N_det_generators
pt2_workload(i) = dfloat(N_det_generators - i + 1)**2
2017-01-12 15:46:10 +01:00
end do
pt2_workload = pt2_workload / sum(pt2_workload)
END_PROVIDER