10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-12-23 21:03:56 +01:00

Can diagonalize much larger spaces (ZMQ)

This commit is contained in:
Anthony Scemama 2018-09-19 15:52:13 +02:00
parent a968003838
commit 4308fa9a1f
3 changed files with 128 additions and 7 deletions

View File

@ -36,6 +36,7 @@ subroutine run_wf
double precision :: t0, t1
integer, external :: zmq_get_dvector, zmq_get_N_det_generators
integer, external :: zmq_get8_dvector
integer, external :: zmq_get_ivector
integer, external :: zmq_get_psi, zmq_get_N_det_selectors
integer, external :: zmq_get_N_states_diag

View File

@ -82,11 +82,14 @@ subroutine davidson_slave_work(zmq_to_qp_run_socket, zmq_socket_push, N_st, sze,
double precision, allocatable :: energy(:)
integer, external :: zmq_get_dvector
integer, external :: zmq_get_dmatrix
allocate(u_t(N_st,N_det))
allocate (energy(N_st))
if (zmq_get_dvector(zmq_to_qp_run_socket, worker_id, 'u_t', u_t, size(u_t)) == -1) then
! Warning : dimensions are permuted for performance considerations, It is OK
! since we get the full matrix
if (zmq_get_dmatrix(zmq_to_qp_run_socket, worker_id, 'u_t', u_t, size(u_t,2), size(u_t,1) ) == -1) then
print *, irp_here, ': Unable to get u_t'
deallocate(u_t,energy)
return
@ -313,6 +316,7 @@ subroutine H_S2_u_0_nstates_zmq(v_0,s_0,u_0,N_st,sze)
double precision :: energy(N_st)
integer, external :: zmq_put_dvector, zmq_put_psi, zmq_put_N_states_diag
integer, external :: zmq_put_dmatrix
energy = 0.d0
@ -325,7 +329,9 @@ subroutine H_S2_u_0_nstates_zmq(v_0,s_0,u_0,N_st,sze)
if (zmq_put_dvector(zmq_to_qp_run_socket,1,'energy',energy,size(energy)) == -1) then
stop 'Unable to put energy on ZMQ server'
endif
if (zmq_put_dvector(zmq_to_qp_run_socket, 1, 'u_t', u_t, size(u_t)) == -1) then
! Warning : dimensions are permuted for performance considerations, It is OK
! since we get the full matrix
if (zmq_put_dmatrix(zmq_to_qp_run_socket, 1, 'u_t', u_t, size(u_t,2),size(u_t,1) ) == -1) then
stop 'Unable to put u_t on ZMQ server'
endif

View File

@ -218,17 +218,20 @@ integer function zmq_put8_dvector(zmq_to_qp_run_socket, worker_id, name, x, size
rc = f77_zmq_send(zmq_to_qp_run_socket,trim(msg),len(trim(msg)),ZMQ_SNDMORE)
if (rc /= len(trim(msg))) then
zmq_put8_dvector = -1
print *, 'Failed in put_data'
return
endif
rc = f77_zmq_send(zmq_to_qp_run_socket,x,size_x*8,0)
if (rc /= size_x*8) then
rc = f77_zmq_send8(zmq_to_qp_run_socket,x,size_x*8_8,0)
if (rc /= size_x*8_8) then
print *, 'Failed in send ', rc, size_x*8, size_x, N_det
zmq_put8_dvector = -1
return
endif
rc = f77_zmq_recv(zmq_to_qp_run_socket,msg,len(msg),0)
if (msg(1:rc) /= 'put_data_reply ok') then
print *, 'Failed in recv ', rc
zmq_put8_dvector = -1
return
endif
@ -270,7 +273,7 @@ integer function zmq_get8_dvector(zmq_to_qp_run_socket, worker_id, name, x, size
go to 10
endif
rc = f77_zmq_recv(zmq_to_qp_run_socket,x,size_x*8,0)
rc = f77_zmq_recv8(zmq_to_qp_run_socket,x,size_x*8,0)
if (rc /= size_x*8) then
print *, irp_here, 'rc /= size_x*8', rc, size_x*8
zmq_get8_dvector = -1
@ -300,6 +303,117 @@ end
integer function zmq_put_dmatrix(zmq_to_qp_run_socket, worker_id, name, x, size_x1, size_x2)
use f77_zmq
implicit none
BEGIN_DOC
! Put a float vector on the qp_run scheduler
END_DOC
integer(ZMQ_PTR), intent(in) :: zmq_to_qp_run_socket
integer, intent(in) :: worker_id
character*(*) :: name
integer, intent(in) :: size_x1, size_x2
double precision, intent(in) :: x(size_x1, size_x2)
integer*8 :: rc
integer :: j
character*(256) :: msg
zmq_put_dmatrix = 0
do j=1,size_x2
write(msg,'(A,1X,I8,1X,A,I8.8)') 'put_data '//trim(zmq_state), worker_id, trim(name), j
rc = f77_zmq_send(zmq_to_qp_run_socket,trim(msg),len(trim(msg)),ZMQ_SNDMORE)
if (rc /= len(trim(msg))) then
zmq_put_dmatrix = -1
print *, 'Failed in put_data', rc, j
return
endif
rc = f77_zmq_send8(zmq_to_qp_run_socket,x(1,j),size_x1*8_8,0)
if (rc /= size_x1*8_8) then
print *, 'Failed in send ', rc, j
zmq_put_dmatrix = -1
return
endif
rc = f77_zmq_recv(zmq_to_qp_run_socket,msg,len(msg),0)
if (msg(1:rc) /= 'put_data_reply ok') then
print *, 'Failed in recv ', rc, j
zmq_put_dmatrix = -1
return
endif
enddo
end
integer function zmq_get_dmatrix(zmq_to_qp_run_socket, worker_id, name, x, size_x1, size_x2)
use f77_zmq
implicit none
BEGIN_DOC
! Get a float vector from the qp_run scheduler
END_DOC
integer(ZMQ_PTR), intent(in) :: zmq_to_qp_run_socket
integer, intent(in) :: worker_id
integer, intent(in) :: size_x1, size_x2
character*(*), intent(in) :: name
double precision, intent(out) :: x(size_x1,size_x2)
integer*8 :: rc
integer :: j
character*(256) :: msg
PROVIDE zmq_state
! Success
zmq_get_dmatrix = 0
if (mpi_master) then
do j=1, size_x2
write(msg,'(A,1X,I8,1X,A,I8.8)') 'get_data '//trim(zmq_state), worker_id, trim(name),j
rc = f77_zmq_send(zmq_to_qp_run_socket,trim(msg),len(trim(msg)),0)
if (rc /= len(trim(msg))) then
zmq_get_dmatrix = -1
print *, irp_here, 'rc /= len(trim(msg))', rc, len(trim(msg))
go to 10
endif
rc = f77_zmq_recv(zmq_to_qp_run_socket,msg,len(msg),0)
if (msg(1:14) /= 'get_data_reply') then
print *, irp_here, 'msg(1:14) /= get_data_reply', msg(1:14)
zmq_get_dmatrix = -1
go to 10
endif
rc = f77_zmq_recv8(zmq_to_qp_run_socket,x(1,j),size_x1*8,0)
if (rc /= size_x1*8) then
print *, irp_here, 'rc /= size_x1*8', rc, size_x1*8
zmq_get_dmatrix = -1
go to 10
endif
enddo
endif
10 continue
IRP_IF MPI_DEBUG
print *, irp_here, mpi_rank
call MPI_BARRIER(MPI_COMM_WORLD, ierr)
IRP_ENDIF
IRP_IF MPI
integer :: ierr
include 'mpif.h'
call MPI_BCAST (zmq_get_dmatrix, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
if (ierr /= MPI_SUCCESS) then
print *, irp_here//': Unable to broadcast zmq_get_dmatrix'
stop -1
endif
call MPI_BARRIER(MPI_COMM_WORLD,ierr)
call broadcast_chunks_double(x, int(size_x1,8)*int(size_x2,8))
IRP_ENDIF
end
integer function zmq_put8_ivector(zmq_to_qp_run_socket, worker_id, name, x, size_x)
use f77_zmq
implicit none
@ -323,7 +437,7 @@ integer function zmq_put8_ivector(zmq_to_qp_run_socket, worker_id, name, x, size
return
endif
rc = f77_zmq_send(zmq_to_qp_run_socket,x,size_x*4,0)
rc = f77_zmq_send8(zmq_to_qp_run_socket,x,size_x*4,0)
if (rc /= size_x*4) then
zmq_put8_ivector = -1
return
@ -364,7 +478,7 @@ integer function zmq_get8_ivector(zmq_to_qp_run_socket, worker_id, name, x, size
go to 10
endif
rc = f77_zmq_recv(zmq_to_qp_run_socket,msg,len(msg),0)
rc = f77_zmq_recv8(zmq_to_qp_run_socket,msg,len(msg),0)
if (msg(1:14) /= 'get_data_reply') then
zmq_get8_ivector = -1
go to 10