5 #include "sll_memory.h"
6 #include "sll_working_precision.h"
13 use sll_m_remapper,
only: &
14 sll_o_get_layout_collective, &
15 sll_o_get_layout_i_max, &
16 sll_o_get_layout_i_min, &
17 sll_o_get_layout_j_max, &
18 sll_o_get_layout_j_min, &
19 sll_o_get_layout_k_max, &
20 sll_o_get_layout_k_min, &
43 type(sll_t_layout_2d),
pointer :: layout
44 sll_int32,
intent(in) :: collective_size
45 sll_int32,
dimension(collective_size),
intent(out) :: disps
46 sll_int32 :: imin, imax
47 sll_int32 :: jmin, jmax
48 sll_int32 :: size_i, size_j
55 do rank = 1, collective_size - 1
56 imin = sll_o_get_layout_i_min(layout, rank - 1)
57 imax = sll_o_get_layout_i_max(layout, rank - 1)
58 jmin = sll_o_get_layout_j_min(layout, rank - 1)
59 jmax = sll_o_get_layout_j_max(layout, rank - 1)
60 size_i = imax - imin + 1
61 size_j = jmax - jmin + 1
62 counter = counter + size_i*size_j
63 disps(rank + 1) = counter
70 type(sll_t_layout_2d),
pointer :: layout
71 sll_real64,
dimension(:, :),
intent(in) ::
data
72 sll_real64,
dimension(:),
intent(out) :: buffer
74 sll_int32 :: data_size
77 sll_int32 :: imin, imax
78 sll_int32 :: jmin, jmax
79 sll_int32 :: size_i, size_j
83 col => sll_o_get_layout_collective(layout)
85 data_size =
size(
data, 1)*
size(
data, 2)
87 imin = sll_o_get_layout_i_min(layout, myrank)
88 imax = sll_o_get_layout_i_max(layout, myrank)
89 jmin = sll_o_get_layout_j_min(layout, myrank)
90 jmax = sll_o_get_layout_j_max(layout, myrank)
91 size_i = imax - imin + 1
92 size_j = jmax - jmin + 1
98 buffer(counter) =
data(i, j)
105 type(sll_t_layout_2d),
pointer :: layout
106 sll_real32,
dimension(:, :),
intent(in) ::
data
107 sll_real32,
dimension(:),
intent(out) :: buffer
109 sll_int32 :: data_size
112 sll_int32 :: imin, imax
113 sll_int32 :: jmin, jmax
114 sll_int32 :: size_i, size_j
118 col => sll_o_get_layout_collective(layout)
120 data_size =
size(
data, 1)*
size(
data, 2)
122 imin = sll_o_get_layout_i_min(layout, myrank)
123 imax = sll_o_get_layout_i_max(layout, myrank)
124 jmin = sll_o_get_layout_j_min(layout, myrank)
125 jmax = sll_o_get_layout_j_max(layout, myrank)
126 size_i = imax - imin + 1
127 size_j = jmax - jmin + 1
132 counter = counter + 1
133 buffer(counter) =
data(i, j)
142 type(sll_t_layout_3d),
pointer :: layout
143 sll_real64,
dimension(:, :, :),
intent(in) ::
data
144 sll_real64,
dimension(:),
intent(out) :: buffer
146 sll_int32 :: data_size
149 sll_int32 :: imin, imax
150 sll_int32 :: jmin, jmax
151 sll_int32 :: kmin, kmax
152 sll_int32 :: size_i, size_j, size_k
156 col => sll_o_get_layout_collective(layout)
158 data_size =
size(
data, 1)*
size(
data, 2)*
size(
data, 3)
159 imin = sll_o_get_layout_i_min(layout, myrank)
160 imax = sll_o_get_layout_i_max(layout, myrank)
161 jmin = sll_o_get_layout_j_min(layout, myrank)
162 jmax = sll_o_get_layout_j_max(layout, myrank)
163 kmin = sll_o_get_layout_k_min(layout, myrank)
164 kmax = sll_o_get_layout_k_max(layout, myrank)
165 size_i = imax - imin + 1
166 size_j = jmax - jmin + 1
167 size_k = kmax - kmin + 1
172 counter = counter + 1
173 buffer(counter) =
data(i, j, k)
182 type(sll_t_layout_2d),
pointer :: layout
183 sll_int32,
intent(in) :: n
184 sll_int32,
dimension(n) :: rc
186 sll_int32 :: imin, imax
187 sll_int32 :: jmin, jmax
188 sll_int32 :: size_i, size_j
191 imin = sll_o_get_layout_i_min(layout, i)
192 imax = sll_o_get_layout_i_max(layout, i)
193 jmin = sll_o_get_layout_j_min(layout, i)
194 jmax = sll_o_get_layout_j_max(layout, i)
195 size_i = imax - imin + 1
196 size_j = jmax - jmin + 1
197 rc(i + 1) = size_i*size_j
204 type(sll_t_layout_2d),
pointer :: layout
205 sll_real64,
dimension(:, :),
intent(out) ::
data
206 sll_real64,
dimension(:),
intent(in) :: buffer
211 sll_int32 :: imin, imax
212 sll_int32 :: jmin, jmax
215 col => sll_o_get_layout_collective(layout)
220 do box = 0, col_sz - 1
221 imin = sll_o_get_layout_i_min(layout, box)
222 imax = sll_o_get_layout_i_max(layout, box)
223 jmin = sll_o_get_layout_j_min(layout, box)
224 jmax = sll_o_get_layout_j_max(layout, box)
229 data(i, j) = buffer(pos)
239 type(sll_t_layout_3d),
pointer :: layout
240 sll_real64,
dimension(:, :, :),
intent(out) ::
data
241 sll_real64,
dimension(:),
intent(in) :: buffer
245 sll_int32 :: imin, imax
246 sll_int32 :: jmin, jmax
247 sll_int32 :: kmin, kmax
248 sll_int32 :: size_i, size_j, size_k
251 col => sll_o_get_layout_collective(layout)
257 imin = sll_o_get_layout_i_min(layout, myrank)
258 imax = sll_o_get_layout_i_max(layout, myrank)
259 jmin = sll_o_get_layout_j_min(layout, myrank)
260 jmax = sll_o_get_layout_j_max(layout, myrank)
261 kmin = sll_o_get_layout_k_min(layout, myrank)
262 kmax = sll_o_get_layout_k_max(layout, myrank)
263 size_i = imax - imin + 1
264 size_j = jmax - jmin + 1
265 size_k = kmax - kmin + 1
272 data(i, j, k) = buffer(pos)
subroutine, public sll_s_load_buffer_2d(layout, data, buffer)
this subroutine loads 2D array into a 1D buffer for use in collective communications
subroutine, public sll_s_unload_buffer_3d(layout, buffer, data)
this routine takes 1D array and stores it in a 3D array warning definition changes are done wrt sll_s...
integer(kind=i32) function, dimension(n), public sll_f_receive_counts_array_2d(layout, n)
this function is also a helper for collective routines for example gatherv should be change into subr...
subroutine, public sll_s_load_buffer_3d(layout, data, buffer)
this subroutine loads 3D array into a 1D buffer for use in collective communications
subroutine, public sll_s_load_buffer32_2d(layout, data, buffer)
subroutine, public sll_s_compute_displacements_array_2d(layout, collective_size, disps)
the objective of this subroutine is to help to prepare data for gatherv operations for example
subroutine, public sll_s_unload_buffer_2d(layout, buffer, data)
integer(kind=i32) function, public sll_f_get_collective_rank(col)
Determines the rank of the calling process in the communicator.
integer(kind=i32) function, public sll_f_get_collective_size(col)
Determines the size of the group associated with a communicator.
Wrapper around the communicator.