10
1
mirror of https://gitlab.com/scemama/qmcchem.git synced 2024-06-02 11:25:18 +02:00

Fixed zmq_ezfio

This commit is contained in:
Anthony Scemama 2018-09-12 09:34:51 +02:00
parent 3c0d404661
commit cc0620ee32

View File

@ -27,7 +27,7 @@ subroutine zmq_ezfio_has(cmd_in,exists)
character*(4) :: buffer_size_char
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, &
len(buffer_size_char), ZMQ_SNDMORE)
len(buffer_size_char), 0)
read(buffer_size_char(1:rc),*) buffer_size
logical :: w
@ -59,7 +59,7 @@ subroutine zmq_ezfio_get_logical(cmd_in,w,d)
rc = f77_zmq_send(zmq_to_dataserver_socket, 'Ezfio', 5, ZMQ_SNDMORE)
rc = f77_zmq_send(zmq_to_dataserver_socket, cmd, len_trim(cmd), 0)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), ZMQ_SNDMORE)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), 0)
read(buffer_size_char(1:rc),*) buffer_size
allocate (buffer(buffer_size+1))
buffer = ' '
@ -107,7 +107,7 @@ subroutine zmq_ezfio_get_double_precision(cmd_in,w,d)
rc = f77_zmq_send(zmq_to_dataserver_socket, 'Ezfio', 5, ZMQ_SNDMORE)
rc = f77_zmq_send(zmq_to_dataserver_socket, cmd, len_trim(cmd), 0)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), ZMQ_SNDMORE)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), 0)
read(buffer_size_char(1:rc),*) buffer_size
allocate (buffer(buffer_size+1))
buffer = ' '
@ -153,7 +153,7 @@ subroutine zmq_ezfio_get_integer(cmd_in,w,d)
rc = f77_zmq_send(zmq_to_dataserver_socket, 'Ezfio', 5, ZMQ_SNDMORE)
rc = f77_zmq_send(zmq_to_dataserver_socket, cmd, len_trim(cmd), 0)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), ZMQ_SNDMORE)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), 0)
read(buffer_size_char(1:rc),*) buffer_size
allocate (buffer(buffer_size+1))
buffer = ' '
@ -199,7 +199,7 @@ subroutine zmq_ezfio_get_integer8(cmd_in,w,d)
rc = f77_zmq_send(zmq_to_dataserver_socket, 'Ezfio', 5, ZMQ_SNDMORE)
rc = f77_zmq_send(zmq_to_dataserver_socket, cmd, len_trim(cmd), 0)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), ZMQ_SNDMORE)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), 0)
read(buffer_size_char(1:rc),*) buffer_size
allocate (buffer(buffer_size+1))
buffer = ' '
@ -247,7 +247,7 @@ subroutine zmq_ezfio_get_real(cmd_in,w,d)
rc = f77_zmq_send(zmq_to_dataserver_socket, 'Ezfio', 5, ZMQ_SNDMORE)
rc = f77_zmq_send(zmq_to_dataserver_socket, cmd, len_trim(cmd), 0)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), ZMQ_SNDMORE)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), 0)
read(buffer_size_char(1:rc),*) buffer_size
allocate (buffer(buffer_size+1))
buffer = ' '
@ -296,7 +296,7 @@ subroutine zmq_ezfio_get_character(cmd_in,w,d)
rc = f77_zmq_send(zmq_to_dataserver_socket, 'Ezfio', 5, ZMQ_SNDMORE)
rc = f77_zmq_send(zmq_to_dataserver_socket, cmd, len_trim(cmd), 0)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), ZMQ_SNDMORE)
rc = f77_zmq_recv(zmq_to_dataserver_socket, buffer_size_char, len(buffer_size_char), 0)
read(buffer_size_char(1:rc),*) buffer_size
allocate (buffer(buffer_size+1))
buffer = ' '