10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-06-26 07:02:14 +02:00

miniserver when needed - untested

This commit is contained in:
Yann Garniron 2016-10-03 10:42:47 +02:00
parent c431da3830
commit 272482e8fa
2 changed files with 22 additions and 13 deletions

View File

@ -22,7 +22,7 @@ subroutine davidson_process(block, N, idx, vt, st)
vt = 0d0
st = 0d0
N = N_det
N = dav_size ! N_det
do i=1,N
idx(i) = i
end do
@ -56,8 +56,8 @@ subroutine davidson_process(block, N, idx, vt, st)
ext = ext + popcnt(xor(sorted_i(ni), sorted_(ni,j,1)))
end do
if(ext <= 4) then
call i_h_j (psi_det(1,1,org_j),psi_det(1,1,org_i),n_int,hij)
call get_s2(psi_det(1,1,org_j),psi_det(1,1,org_i),n_int,s2)
call i_h_j (dav_det(1,1,org_j),dav_det(1,1,org_i),n_int,hij) ! psi_det
call get_s2(dav_det(1,1,org_j),dav_det(1,1,org_i),n_int,s2)
do istate=1,N_states
vt (istate,org_i) = vt (istate,org_i) + hij*dav_ut(istate,org_j)
vt (istate,org_j) = vt (istate,org_j) + hij*dav_ut(istate,org_i)
@ -100,6 +100,7 @@ subroutine davidson_init(zmq_to_qp_run_socket)
integer(ZMQ_PTR), intent(out) :: zmq_to_qp_run_socket ! zmq_to_qp_run_socket
touch dav_size dav_det dav_ut
call new_parallel_job(zmq_to_qp_run_socket,'davidson')
end subroutine
@ -195,7 +196,7 @@ subroutine davidson_slave_work(zmq_to_qp_run_socket, zmq_socket_push, worker_id)
allocate(vt(N_states, N_det))
allocate(st(N_states, N_det))
call hobo_get()
!call hobo_get()
do
call get_task_from_taskserver(zmq_to_qp_run_socket,worker_id, task_id, task)
@ -347,9 +348,9 @@ subroutine davidson_run(zmq_to_qp_run_socket , v0, s0)
call davidson_collector(zmq_collector, zmq_socket_pull , v0, s0)
call end_zmq_to_qp_run_socket(zmq_collector)
call end_zmq_pull_socket(zmq_socket_pull)
call hobo_server_end()
call davidson_miniserver_end()
else if(i==1) then
call hobo_server()
call davidson_miniserver_run()
else
call davidson_slave_inproc(i)
endif
@ -359,7 +360,7 @@ end subroutine
subroutine hobo_server()
subroutine davidson_miniserver_run()
use f77_zmq
implicit none
integer(ZMQ_PTR) context
@ -378,9 +379,9 @@ subroutine hobo_server()
do
rc = f77_zmq_recv(responder, buffer, 5, 0)
if (buffer(1:rc) /= 'end') then
rc = f77_zmq_send (responder, N_det, 4, ZMQ_SNDMORE)
rc = f77_zmq_send (responder, psi_det, 16*N_int*N_det, ZMQ_SNDMORE)
rc = f77_zmq_send (responder, ut, 8*N_det*N_states, 0)
rc = f77_zmq_send (responder, dav_size, 4, ZMQ_SNDMORE)
rc = f77_zmq_send (responder, dav_det, 16*N_int*N_det, ZMQ_SNDMORE)
rc = f77_zmq_send (responder, dav_ut, 8*N_det*N_states, 0)
else
rc = f77_zmq_send (responder, "end", 3, 0)
exit
@ -392,7 +393,7 @@ subroutine hobo_server()
end subroutine
subroutine hobo_server_end()
subroutine davidson_miniserver_end()
implicit none
use f77_zmq
@ -414,7 +415,7 @@ subroutine hobo_server_end()
end subroutine
subroutine hobo_get()
subroutine davidson_miniserver_get()
implicit none
use f77_zmq
@ -440,6 +441,8 @@ subroutine hobo_get()
rc = f77_zmq_recv(requester, dav_ut, 8*dav_size*N_states, 0)
rc = f77_zmq_close(requester)
rc = f77_zmq_ctx_destroy(context)
touch dav_det dav_ut
end subroutine

View File

@ -227,6 +227,8 @@ subroutine H_S2_u_0_nstates(v_0,s_0,u_0,H_jj,S2_jj,n,keys_tmp,Nint,N_st,sze_8)
v_0 = 0.d0
s_0 = 0.d0
if(n /= N_det) stop "n /= N_det"
do i=1,n
do istate=1,N_st
ut(istate,i) = u_0(i,istate)
@ -234,7 +236,11 @@ subroutine H_S2_u_0_nstates(v_0,s_0,u_0,H_jj,S2_jj,n,keys_tmp,Nint,N_st,sze_8)
enddo
call sort_dets_ab_v(keys_tmp, sorted(1,1,1), sort_idx(1,1), shortcut(0,1), version(1,1,1), n, Nint)
call sort_dets_ba_v(keys_tmp, sorted(1,1,2), sort_idx(1,2), shortcut(0,2), version(1,1,2), n, Nint)
dav_size = n
dav_det = psi_det
dav_ut = ut
call davidson_init(handler)
do sh=shortcut(0,1),1,-1