mirror of
https://github.com/LCPQ/quantum_package
synced 2024-06-02 11:25:26 +02:00
274 lines
8.6 KiB
Plaintext
274 lines
8.6 KiB
Plaintext
subroutine four_index_transform_zmq(map_a,map_c,matrix_B,LDB, &
|
|
i_start, j_start, k_start, l_start, &
|
|
i_end , j_end , k_end , l_end , &
|
|
a_start, b_start, c_start, d_start, &
|
|
a_end , b_end , c_end , d_end )
|
|
implicit none
|
|
use f77_zmq
|
|
use map_module
|
|
BEGIN_DOC
|
|
! Performs a four-index transformation of map_a(N^4) into map_c(M^4) using b(NxM)
|
|
! C_{abcd} = \sum_{ijkl} A_{ijkl}.B_{ia}.B_{jb}.B_{kc}.B_{ld}
|
|
! Loops run over *_start->*_end
|
|
END_DOC
|
|
type(map_type), intent(in) :: map_a
|
|
type(map_type), intent(inout) :: map_c
|
|
integer, intent(in) :: LDB
|
|
double precision, intent(in) :: matrix_B(LDB,*)
|
|
integer, intent(in) :: i_start, j_start, k_start, l_start
|
|
integer, intent(in) :: i_end , j_end , k_end , l_end
|
|
integer, intent(in) :: a_start, b_start, c_start, d_start
|
|
integer, intent(in) :: a_end , b_end , c_end , d_end
|
|
|
|
double precision, allocatable :: T(:,:), U(:,:,:), V(:,:)
|
|
double precision, allocatable :: T2d(:,:), V2d(:,:)
|
|
integer :: i_max, j_max, k_max, l_max
|
|
integer :: i_min, j_min, k_min, l_min
|
|
integer :: i, j, k, l, ik, ll
|
|
integer :: l_start_block, l_end_block, l_block
|
|
integer :: a, b, c, d
|
|
double precision, external :: get_ao_bielec_integral
|
|
integer*8 :: ii
|
|
integer(key_kind) :: idx
|
|
real(integral_kind) :: tmp
|
|
integer(key_kind), allocatable :: key(:)
|
|
real(integral_kind), allocatable :: value(:)
|
|
integer*8, allocatable :: l_pointer(:)
|
|
|
|
ASSERT (k_start == i_start)
|
|
ASSERT (l_start == j_start)
|
|
ASSERT (a_start == c_start)
|
|
ASSERT (b_start == d_start)
|
|
|
|
i_min = min(i_start,a_start)
|
|
i_max = max(i_end ,a_end )
|
|
j_min = min(j_start,b_start)
|
|
j_max = max(j_end ,b_end )
|
|
k_min = min(k_start,c_start)
|
|
k_max = max(k_end ,c_end )
|
|
l_min = min(l_start,d_start)
|
|
l_max = max(l_end ,d_end )
|
|
|
|
ASSERT (0 < i_max)
|
|
ASSERT (0 < j_max)
|
|
ASSERT (0 < k_max)
|
|
ASSERT (0 < l_max)
|
|
ASSERT (LDB >= i_max)
|
|
ASSERT (LDB >= j_max)
|
|
ASSERT (LDB >= k_max)
|
|
ASSERT (LDB >= l_max)
|
|
|
|
|
|
integer(ZMQ_PTR) :: zmq_to_qp_run_socket
|
|
call new_parallel_job(zmq_to_qp_run_socket,'four_idx')
|
|
|
|
integer*8 :: new_size
|
|
new_size = max(1024_8, 5_8 * map_a % n_elements )
|
|
|
|
integer :: npass
|
|
integer*8 :: tempspace
|
|
|
|
tempspace = (new_size * 14_8) / (1024_8 * 1024_8)
|
|
npass = min(l_end-l_start,1 + tempspace / 2048) ! 2 GiB of scratch space
|
|
l_block = (l_end-l_start)/npass
|
|
|
|
! Create tasks
|
|
! ============
|
|
|
|
character(len=64), allocatable :: task
|
|
|
|
do l_start_block = l_start, l_end, l_block
|
|
l_end_block = min(l_end, l_start_block+l_block-1)
|
|
write(task,'I10,X,I10') l_start_block, l_end_block
|
|
call add_task_to_taskserver(zmq_to_qp_run_socket,trim(task))
|
|
enddo
|
|
|
|
call zmq_set_running(zmq_to_qp_run_socket)
|
|
|
|
PROVIDE nproc
|
|
|
|
call omp_set_nested(.True.)
|
|
integer :: ithread
|
|
!$OMP PARALLEL NUM_THREADS(2) PRIVATE(ithread)
|
|
ithread = omp_get_thread_num()
|
|
if (ithread==0) then
|
|
call four_idx_collector(zmq_to_qp_run_socket,map_c)
|
|
else
|
|
!TODO : Put strings of map_a and matrix_b on server and broadcast
|
|
call four_index_transform_slave_inproc(map_a,map_c,matrix_B,LDB, &
|
|
i_start, j_start, k_start, l_start_block, &
|
|
i_end , j_end , k_end , l_end_block , &
|
|
a_start, b_start, c_start, d_start, &
|
|
a_end , b_end , c_end , d_end, 1 )
|
|
endif
|
|
!$OMP END PARALLEL
|
|
|
|
call end_parallel_job(zmq_to_qp_run_socket, 'four_idx')
|
|
|
|
|
|
end
|
|
|
|
|
|
subroutine four_idx_slave_work(zmq_to_qp_run_socket, worker_id)
|
|
use f77_zmq
|
|
implicit none
|
|
|
|
integer(ZMQ_PTR),intent(in) :: zmq_to_qp_run_socket
|
|
integer,intent(in) :: worker_id
|
|
integer :: task_id
|
|
character*(512) :: msg
|
|
|
|
integer :: i_start, j_start, k_start, l_start_block
|
|
integer :: i_end , j_end , k_end , l_end_block
|
|
integer :: a_start, b_start, c_start, d_start
|
|
integer :: a_end , b_end , c_end , d_end
|
|
|
|
!TODO : get map_a and matrix_B from server
|
|
do
|
|
call get_task_from_taskserver(zmq_to_qp_run_socket,worker_id, task_id, msg)
|
|
if(task_id == 0) exit
|
|
read (msg,*) LDB, &
|
|
i_start, j_start, k_start, l_start_block, &
|
|
i_end , j_end , k_end , l_end_block , &
|
|
a_start, b_start, c_start, d_start, &
|
|
a_end , b_end , c_end , d_end
|
|
|
|
call four_index_transform_slave(map_a,map_c,matrix_B,LDB, &
|
|
i_start, j_start, k_start, l_start_block, &
|
|
i_end , j_end , k_end , l_end_block , &
|
|
a_start, b_start, c_start, d_start, &
|
|
a_end , b_end , c_end , d_end, zmq_to_qp_run_socket, &
|
|
task_id)
|
|
call task_done_to_taskserver(zmq_to_qp_run_socket,worker_id,task_id)
|
|
|
|
enddo
|
|
end
|
|
|
|
|
|
BEGIN_PROVIDER [ integer, nthreads_four_idx ]
|
|
implicit none
|
|
BEGIN_DOC
|
|
! Number of threads for 4-index transformation
|
|
END_DOC
|
|
nthreads_four_idx = nproc
|
|
character*(32) :: env
|
|
call getenv('NTHREADS_FOUR_IDX',env)
|
|
if (trim(env) /= '') then
|
|
read(env,*) nthreads_four_idx
|
|
endif
|
|
call write_int(6,nthreads_davidson,'Number of threads for 4-index transformation')
|
|
END_PROVIDER
|
|
|
|
|
|
|
|
subroutine four_idx_collector(zmq_to_qp_run_socket,map_c)
|
|
use f77_zmq
|
|
use map_module
|
|
implicit none
|
|
type(map_type), intent(inout) :: map_c
|
|
|
|
integer :: more
|
|
integer(ZMQ_PTR), external :: new_zmq_pull_socket
|
|
integer(ZMQ_PTR) :: zmq_socket_pull
|
|
|
|
|
|
more = 1
|
|
zmq_socket_pull = new_zmq_pull_socket()
|
|
|
|
do while (more == 1)
|
|
call four_idx_pull_results(zmq_socket_pull, map_c, task_id)
|
|
call zmq_delete_task(zmq_to_qp_run_socket,zmq_socket_pull,task_id,more)
|
|
enddo
|
|
|
|
call end_zmq_pull_socket(zmq_socket_pull)
|
|
|
|
end
|
|
|
|
|
|
subroutine four_idx_pull_results(zmq_socket_pull, map_c, task_id)
|
|
use f77_zmq
|
|
use map_module
|
|
implicit none
|
|
type(map_type), intent(inout) :: map_c
|
|
integer(ZMQ_PTR), intent(inout) :: zmq_socket_pull
|
|
|
|
integer, intent(out) :: task_id
|
|
|
|
integer :: rc, sze
|
|
integer*8 :: rc8
|
|
|
|
|
|
rc = f77_zmq_recv( zmq_socket_pull, task_id, 4, 0)
|
|
if(rc /= 4) stop "four_idx_pull_results failed to pull task_id"
|
|
|
|
rc = f77_zmq_recv( zmq_socket_pull, sze, 4, 0)
|
|
if(rc /= 4) stop "four_idx_pull_results failed to pull sze"
|
|
|
|
integer(key_kind), allocatable :: key(:)
|
|
real(integral_kind), allocatable :: value(:)
|
|
|
|
allocate(key(sze), value(sze))
|
|
|
|
rc8 = f77_zmq_recv8( zmq_socket_pull, key, key_kind*sze, 0)
|
|
if(rc8 /= key_kind*sze) stop "four_idx_pull_results failed to pull key"
|
|
|
|
rc8 = f77_zmq_recv8( zmq_socket_pull, value, integral_kind*sze, 0)
|
|
if(rc8 /= integral_kind*sze) stop "four_idx_pull_results failed to pull value"
|
|
|
|
! Activate if zmq_socket_pull is a REP
|
|
IRP_IF ZMQ_PUSH
|
|
IRP_ELSE
|
|
rc = f77_zmq_send( zmq_socket_pull, 0, 4, 0)
|
|
if (rc /= 4) then
|
|
print *, irp_here, ' : f77_zmq_send (zmq_socket_pull,...'
|
|
stop 'error'
|
|
endif
|
|
IRP_ENDIF
|
|
|
|
call map_update(map_c, key, value, sze, 1.d-15) ! TODO : threshold
|
|
|
|
deallocate(key, value)
|
|
end
|
|
|
|
|
|
|
|
subroutine four_idx_push_results(zmq_socket_push, key, value, sze, task_id)
|
|
use f77_zmq
|
|
use map_module
|
|
implicit none
|
|
integer, intent(in) :: sze
|
|
integer(key_kind), intent(in) :: key(sze)
|
|
real(integral_kind), intent(in) :: value(sze)
|
|
integer(ZMQ_PTR), intent(in) :: zmq_socket_push
|
|
integer, intent(in) :: task_id
|
|
|
|
integer :: rc, sze
|
|
integer*8 :: rc8
|
|
|
|
|
|
rc = f77_zmq_send( zmq_socket_push, task_id, 4, ZMQ_SNDMORE)
|
|
if(rc /= 4) stop "four_idx_push_results failed to push task_id"
|
|
|
|
rc = f77_zmq_send( zmq_socket_push, sze, 4, ZMQ_SNDMORE)
|
|
if(rc /= 4) stop "four_idx_push_results failed to push sze"
|
|
|
|
rc8 = f77_zmq_send8( zmq_socket_push, key, key_kind*sze, ZMQ_SNDMORE)
|
|
if(rc8 /= key_kind*sze) stop "four_idx_push_results failed to push key"
|
|
|
|
rc8 = f77_zmq_send8( zmq_socket_push, value, integral_kind*sze, 0)
|
|
if(rc8 /= integral_kind*sze) stop "four_idx_push_results failed to push value"
|
|
|
|
! Activate if zmq_socket_push is a REP
|
|
IRP_IF ZMQ_PUSH
|
|
IRP_ELSE
|
|
rc = f77_zmq_send( zmq_socket_push, 0, 4, 0)
|
|
if (rc /= 4) then
|
|
print *, irp_here, ' : f77_zmq_send (zmq_socket_push,...'
|
|
stop 'error'
|
|
endif
|
|
IRP_ENDIF
|
|
|
|
end
|
|
|
|
|