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

ZMQ checks

This commit is contained in:
Anthony Scemama 2017-03-25 11:57:06 +01:00
parent 2cd4a513dc
commit 3d21999c7e
4 changed files with 29 additions and 6 deletions

View File

@ -23,7 +23,7 @@ subroutine run_wf
integer(ZMQ_PTR), external :: new_zmq_to_qp_run_socket
integer(ZMQ_PTR) :: zmq_to_qp_run_socket
double precision :: energy(N_states)
character*(64) :: states(2)
character*(64) :: states(4)
integer :: rc, i
logical :: force_update

View File

@ -1,6 +1,6 @@
! DO NOT MODIFY BY HAND
! Created by $QP_ROOT/scripts/ezfio_interface/ei_handler.py
! from file /panfs/panasas/cnt0024/cpq1738/scemama/workdir/quantum_package/src/mrcc_selected/EZFIO.cfg
! from file /ccc/work/cont003/gen1738/scemama/quantum_package/src/mrcc_selected/EZFIO.cfg
BEGIN_PROVIDER [ double precision, thresh_dressed_ci ]

View File

@ -145,13 +145,14 @@ subroutine davidson_collect(N, idx, vt, st , v0t, s0t)
end subroutine
subroutine davidson_init(zmq_to_qp_run_socket,u,n0,n,n_st,update_dets)
subroutine davidson_init(zmq_to_qp_run_socket,dets_in,u,n0,n,n_st,update_dets)
use f77_zmq
implicit none
integer(ZMQ_PTR), intent(out) :: zmq_to_qp_run_socket
integer, intent(in) :: n0,n, n_st, update_dets
double precision, intent(in) :: u(n0,n_st)
integer(bit_kind), intent(in) :: dets_in(N_int,2,n)
integer :: i,k
@ -160,8 +161,8 @@ subroutine davidson_init(zmq_to_qp_run_socket,u,n0,n,n_st,update_dets)
touch dav_size
do i=1,dav_size
do k=1,N_int
dav_det(k,1,i) = psi_det(k,1,i)
dav_det(k,2,i) = psi_det(k,2,i)
dav_det(k,1,i) = dets_in(k,1,i)
dav_det(k,2,i) = dets_in(k,2,i)
enddo
enddo
touch dav_det
@ -527,18 +528,40 @@ subroutine davidson_miniserver_get(force_update)
rc = f77_zmq_connect(requester,address)
rc = f77_zmq_send(requester, 'ut', 2, 0)
rc = f77_zmq_recv(requester, update_dets, 4, 0)
if (rc /= 4) then
print *, irp_here, ': f77_zmq_recv(requester, update_dets, 4, 0)'
print *, irp_here, ': rc = ', rc
endif
rc = f77_zmq_recv(requester, dav_size, 4, 0)
if (rc /= 4) then
print *, irp_here, ': f77_zmq_recv(requester, dav_size, 4, 0)'
print *, irp_here, ': rc = ', rc
endif
if (update_dets == 1 .or. force_update) then
TOUCH dav_size
endif
rc = f77_zmq_recv(requester, dav_ut, 8*dav_size*N_states_diag, 0)
if (rc /= 8*dav_size*N_states_diag) then
print *, irp_here, ': f77_zmq_recv(requester, dav_ut, 8*dav_size*N_states_diag, 0)'
print *, irp_here, ': rc = ', rc
endif
SOFT_TOUCH dav_ut
if (update_dets == 1 .or. force_update) then
rc = f77_zmq_send(requester, 'det', 3, 0)
rc = f77_zmq_recv(requester, dav_size, 4, 0)
if (rc /= 4) then
print *, irp_here, ': f77_zmq_recv(requester, dav_size, 4, 0)'
print *, irp_here, ': rc = ', rc
endif
rc = f77_zmq_recv(requester, dav_det, 16*N_int*dav_size, 0)
if (rc /= 16*N_int*dav_size) then
print *, irp_here, ': f77_zmq_recv(requester, dav_det, 16*N_int*dav_size, 0)'
print *, irp_here, ': rc = ', rc
endif
SOFT_TOUCH dav_det
endif

View File

@ -308,7 +308,7 @@ subroutine H_S2_u_0_nstates_zmq(v_0,s_0,u_0,H_jj,S2_jj,n,keys_tmp,Nint,N_st,sze_
v_0 = 0.d0
s_0 = 0.d0
call davidson_init(handler,u_0,size(u_0,1),n,N_st,update_dets)
call davidson_init(handler,keys_tmp,u_0,size(u_0,1),n,N_st,update_dets)
ave_workload = 0.d0
do sh=1,shortcut_(0,1)