1subroutine global_gather_&
3&(obj, z, n, mpi_comm_rows, mpi_comm_cols, npc_n, np_prev, np_next, &
13#ifdef WITH_OPENMP_TRADITIONAL
18 integer(kind=ik),
intent(in) :: mpi_comm_cols, mpi_comm_rows
19 integer(kind=ik),
intent(in) :: npc_n, np_prev, np_next
21 integer(kind=MPI_KIND) :: mpierr, np_rowsMPI, np_colsMPI
23 integer(kind=ik) :: n, np_rows, np_cols
24 real(kind=real_datatype) :: z(n)
25 real(kind=real_datatype) :: tmp(n)
26 integer(kind=ik) :: np
27 integer(kind=MPI_KIND) :: allreduce_request1, allreduce_request2
28 logical :: useNonBlockingCollectivesCols
29 logical :: useNonBlockingCollectivesRows
30 integer(kind=c_int) :: non_blocking_collectives_rows, error, &
31 non_blocking_collectives_cols
36 call obj%get(
"nbc_row_global_gather", non_blocking_collectives_rows, error)
37 if (error .ne. elpa_ok)
then
38 write(error_unit,*)
"Problem getting option for non blocking collectives for rows in global_gather. Aborting..."
43 call obj%get(
"nbc_col_global_gather", non_blocking_collectives_cols, error)
44 if (error .ne. elpa_ok)
then
45 write(error_unit,*)
"Problem getting option for non blocking collectives for cols in global_gather. Aborting..."
50 if (non_blocking_collectives_rows .eq. 1)
then
51 usenonblockingcollectivesrows = .true.
53 usenonblockingcollectivesrows = .false.
56 if (non_blocking_collectives_cols .eq. 1)
then
57 usenonblockingcollectivescols = .true.
59 usenonblockingcollectivescols = .false.
63 call obj%timer%start(
"mpi_communication")
64 call mpi_comm_size(int(mpi_comm_rows,kind=mpi_kind) ,np_rowsmpi, mpierr)
65 np_rows = int(np_rowsmpi,kind=c_int)
68 call mpi_comm_size(int(mpi_comm_cols,kind=mpi_kind) ,np_colsmpi, mpierr)
69 np_cols = int(np_colsmpi,kind=c_int)
71 call obj%timer%stop(
"mpi_communication")
75 if (npc_n==1 .and. np_rows==1)
return
79 if (usenonblockingcollectivesrows)
then
80 call obj%timer%start(
"mpi_nbc_communication")
81 call mpi_iallreduce(z, tmp, int(n,kind=mpi_kind), mpi_real_precision, mpi_sum, int(mpi_comm_rows,kind=mpi_kind), &
82 allreduce_request1, mpierr)
83 call mpi_wait(allreduce_request1, mpi_status_ignore, mpierr)
84 call obj%timer%stop(
"mpi_nbc_communication")
86 call obj%timer%start(
"mpi_communication")
87 call mpi_allreduce(z, tmp, int(n,kind=mpi_kind), mpi_real_precision, mpi_sum, int(mpi_comm_rows,kind=mpi_kind), &
89 call obj%timer%stop(
"mpi_communication")
101 if (npc_n==np_cols)
then
103 if (usenonblockingcollectivescols)
then
104 call obj%timer%start(
"mpi_nbc_communication")
105 call mpi_iallreduce(tmp, z, int(n,kind=mpi_kind), mpi_real_precision, mpi_sum, int(mpi_comm_cols,kind=mpi_kind), &
106 allreduce_request2, mpierr)
107 call mpi_wait(allreduce_request2, mpi_status_ignore, mpierr)
108 call obj%timer%stop(
"mpi_nbc_communication")
110 call obj%timer%start(
"mpi_communication")
111 call mpi_allreduce(tmp, z, int(n,kind=mpi_kind), mpi_real_precision, mpi_sum, int(mpi_comm_cols,kind=mpi_kind), &
113 call obj%timer%stop(
"mpi_communication")
127 call obj%timer%start(
"mpi_communication")
128 call mpi_sendrecv_replace(z, int(n,kind=mpi_kind), mpi_real_precision, int(np_next,kind=mpi_kind), &
129 1112_mpi_kind+int(np,kind=mpi_kind), &
130 int(np_prev,kind=mpi_kind), 1112_mpi_kind+int(np,kind=mpi_kind), &
131 int(mpi_comm_cols,kind=mpi_kind), mpi_status_ignore, mpierr)
132 call obj%timer%stop(
"mpi_communication")