10
0
mirror of https://github.com/LCPQ/quantum_package synced 2025-01-09 20:48:47 +01:00

corrected ZMQ bug - removed stack array

This commit is contained in:
Yann Garniron 2016-10-06 10:54:40 +02:00
parent 1afee06817
commit e2154fb892
3 changed files with 5 additions and 28 deletions

View File

@ -18,9 +18,9 @@ subroutine davidson_process(blockb, blocke, N, idx, vt, st)
integer :: i, j, sh, sh2, exa, ext, org_i, org_j, istate, ni, endi integer :: i, j, sh, sh2, exa, ext, org_i, org_j, istate, ni, endi
integer(bit_kind) :: sorted_i(N_int) integer(bit_kind) :: sorted_i(N_int)
double precision :: s2, hij double precision :: s2, hij
logical :: wrotten(dav_size) logical, allocatable :: wrotten(:)
allocate(wrotten(dav_size))
wrotten = .false. wrotten = .false.
do sh = blockb, blocke do sh = blockb, blocke
@ -51,7 +51,7 @@ subroutine davidson_process(blockb, blocke, N, idx, vt, st)
ext = ext + popcnt(xor(sorted_i(ni), sorted_(ni,j,1))) ext = ext + popcnt(xor(sorted_i(ni), sorted_(ni,j,1)))
end do end do
if(ext <= 4) then if(ext <= 4) then
call i_h_j (dav_det(1,1,org_j),dav_det(1,1,org_i),n_int,hij) ! psi_det call i_h_j (dav_det(1,1,org_j),dav_det(1,1,org_i),n_int,hij)
call get_s2(dav_det(1,1,org_j),dav_det(1,1,org_i),n_int,s2) call get_s2(dav_det(1,1,org_j),dav_det(1,1,org_i),n_int,s2)
if(.not. wrotten(org_i)) then if(.not. wrotten(org_i)) then
wrotten(org_i) = .true. wrotten(org_i) = .true.
@ -116,7 +116,7 @@ subroutine davidson_init(zmq_to_qp_run_socket)
use f77_zmq use f77_zmq
implicit none implicit none
integer(ZMQ_PTR), intent(out) :: zmq_to_qp_run_socket ! zmq_to_qp_run_socket integer(ZMQ_PTR), intent(out) :: zmq_to_qp_run_socket
touch dav_size dav_det dav_ut touch dav_size dav_det dav_ut
call new_parallel_job(zmq_to_qp_run_socket,"davidson") call new_parallel_job(zmq_to_qp_run_socket,"davidson")
@ -184,8 +184,6 @@ subroutine davidson_run_slave(thread,iproc)
end if end if
call davidson_slave_work(zmq_to_qp_run_socket, zmq_socket_push, worker_id) call davidson_slave_work(zmq_to_qp_run_socket, zmq_socket_push, worker_id)
! print *, "done slavin'"
!call sleep(1)
call disconnect_from_taskserver(zmq_to_qp_run_socket,zmq_socket_push,worker_id) 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_to_qp_run_socket(zmq_to_qp_run_socket)
call end_zmq_push_socket(zmq_socket_push,thread) call end_zmq_push_socket(zmq_socket_push,thread)
@ -444,7 +442,6 @@ subroutine davidson_miniserver_end()
rc = f77_zmq_recv(requester, buf, 3, 0) rc = f77_zmq_recv(requester, buf, 3, 0)
rc = f77_zmq_close(requester) rc = f77_zmq_close(requester)
rc = f77_zmq_ctx_destroy(context) rc = f77_zmq_ctx_destroy(context)
! print *, "closed miniserver"
end subroutine end subroutine

View File

@ -249,7 +249,7 @@ subroutine H_S2_u_0_nstates(v_0,s_0,u_0,H_jj,S2_jj,n,keys_tmp,Nint,N_st,sze_8)
call davidson_init(handler) call davidson_init(handler)
do sh=shortcut(0,1),1,-1 do sh=shortcut(0,1),1,-1
workload += (shortcut(sh+1,1) - shortcut(sh,1))**2 workload += (shortcut(sh+1,1) - shortcut(sh,1))**2
if(workload > 1000) then if(workload > 10000) then
blocke = sh blocke = sh
call davidson_add_task(handler, blocke, blockb) call davidson_add_task(handler, blocke, blockb)
blockb = sh-1 blockb = sh-1

View File

@ -473,26 +473,6 @@ subroutine end_zmq_push_socket(zmq_socket_push,thread)
integer :: rc integer :: rc
character*(8), external :: zmq_port character*(8), external :: zmq_port
if (thread == 1) then
rc = f77_zmq_disconnect(zmq_socket_push,zmq_socket_push_inproc_address)
! if (rc /= 0) then
! print *, 'f77_zmq_disconnect(zmq_socket_push,zmq_socket_push_inproc_address)'
! stop 'error'
! endif
else
rc = f77_zmq_disconnect(zmq_socket_push,zmq_socket_push_tcp_address)
if (rc /= 0) then
print *, 'f77_zmq_disconnect(zmq_socket_push,zmq_socket_push_tcp_address)'
stop 'error'
endif
endif
! rc = f77_zmq_setsockopt(zmq_socket_push,ZMQ_LINGER,20000,4)
! if (rc /= 0) then
! stop 'Unable to set ZMQ_LINGER on push socket'
! endif
rc = f77_zmq_close(zmq_socket_push) rc = f77_zmq_close(zmq_socket_push)
if (rc /= 0) then if (rc /= 0) then
print *, 'f77_zmq_close(zmq_socket_push)' print *, 'f77_zmq_close(zmq_socket_push)'