Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
redist.f90
Go to the documentation of this file.
1! Copyright (c) 2021-2025, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34module redist
35 use mesh_field, only : mesh_fld_t
38 use mpi_f08, only : mpi_status, mpi_allreduce, mpi_sendrecv, &
39 mpi_get_count, mpi_integer, mpi_max, mpi_in_place
40 use htable, only : htable_i4_t
41 use point, only : point_t
43 use curve, only : curve_t
44 use comm, only : pe_size, pe_rank, neko_comm
45 use mesh, only : mesh_t, neko_msh_max_zlbls
48 use element, only : element_t
49 use utils, only : neko_error
50 implicit none
51 private
52
53 public :: redist_mesh
54
55contains
56
58 subroutine redist_mesh(msh, parts)
59 type(mesh_t), intent(inout), target :: msh
60 type(mesh_fld_t), intent(in) :: parts
61 type(stack_nh_t), allocatable :: new_mesh_dist(:)
62 type(stack_nz_t), allocatable :: new_zone_dist(:)
63 type(stack_nc_t), allocatable :: new_curve_dist(:)
64 type(nmsh_hex_t) :: el
65 type(nmsh_hex_t), allocatable :: recv_buf_msh(:)
66 type(nmsh_zone_t), allocatable :: recv_buf_zone(:)
67 type(nmsh_curve_el_t), allocatable :: recv_buf_curve(:)
68 class(element_t), pointer :: ep
69 integer, allocatable :: recv_buf_idx(:), send_buf_idx(:)
70 type(mpi_status) :: status
71 integer :: i, j, k, ierr, max_recv_idx, label
72 integer :: src, dst, recv_size, gdim, tmp, new_el_idx, new_pel_idx
73 integer :: max_recv(3)
74 type(point_t) :: p(8)
75 type(htable_i4_t) :: el_map, glb_map
76 type(stack_i4_t) :: pe_lst
77
78
79 !
80 ! Reset possible periodic ids
81 !
82 call msh%reset_periodic_ids()
83
84 !
85 ! Extract new zone distributions
86 !
87
88 allocate(new_zone_dist(0:pe_size - 1))
89 do i = 0, pe_size - 1
90 call new_zone_dist(i)%init()
91 end do
92
93 call redist_zone(msh, msh%periodic, 5, parts, new_zone_dist)
94
95 do j = 1, neko_msh_max_zlbls
96 label = j
97 call redist_zone(msh, msh%labeled_zones(j), 7, parts, &
98 new_zone_dist, label)
99 end do
100
101 !
102 ! Extract new curve info. distributions
103 !
104 allocate(new_curve_dist(0:pe_size - 1))
105 do i = 0, pe_size - 1
106 call new_curve_dist(i)%init()
107 end do
108
109 call redist_curve(msh, msh%curve, parts, new_curve_dist)
110
111 !
112 ! Extract new mesh distribution
113 !
114
115 allocate(new_mesh_dist(0:(pe_size - 1)))
116 do i = 0, pe_size - 1
117 call new_mesh_dist(i)%init()
118 end do
119
120 do i = 1, msh%nelv
121 ep => msh%elements(i)%e
122 el%el_idx = ep%id()
123 do j = 1, 8
124 el%v(j)%v_idx = ep%pts(j)%p%id()
125 el%v(j)%v_xyz = ep%pts(j)%p%x
126 end do
127 call new_mesh_dist(parts%data(i))%push(el)
128 end do
129 nullify(ep)
130
131 gdim = msh%gdim
132 call msh%free()
133
134 max_recv = 0
135 do i = 0, pe_size - 1
136 max_recv(1) = max(max_recv(1), new_mesh_dist(i)%size())
137 max_recv(2) = max(max_recv(2), new_zone_dist(i)%size())
138 max_recv(3) = max(max_recv(3), new_curve_dist(i)%size())
139 end do
140
141 call mpi_allreduce(mpi_in_place, max_recv, 3, mpi_integer, &
142 mpi_max, neko_comm, ierr)
143
144 allocate(recv_buf_msh(max_recv(1)))
145 allocate(recv_buf_zone(max_recv(2)))
146 allocate(recv_buf_curve(max_recv(3)))
147
148 do i = 1, pe_size - 1
149 src = modulo(pe_rank - i + pe_size, pe_size)
150 dst = modulo(pe_rank + i, pe_size)
151
152 ! We should use the %array() procedure, which works great for
153 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
154 ! certain data types
155 select type (nmd_array => new_mesh_dist(dst)%data)
156 type is (nmsh_hex_t)
157 call mpi_sendrecv(nmd_array, &
158 new_mesh_dist(dst)%size(), mpi_nmsh_hex, dst, 0, recv_buf_msh, &
159 max_recv(1), mpi_nmsh_hex, src, 0, neko_comm, status, ierr)
160 end select
161 call mpi_get_count(status, mpi_nmsh_hex, recv_size, ierr)
162
163 do j = 1, recv_size
164 call new_mesh_dist(pe_rank)%push(recv_buf_msh(j))
165 end do
166
167 ! We should use the %array() procedure, which works great for
168 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
169 ! certain data types
170 select type (nzd_array => new_zone_dist(dst)%data)
171 type is (nmsh_zone_t)
172 call mpi_sendrecv(nzd_array, &
173 new_zone_dist(dst)%size(), mpi_nmsh_zone, dst, 1, recv_buf_zone,&
174 max_recv(2), mpi_nmsh_zone, src, 1, neko_comm, status, ierr)
175 end select
176 call mpi_get_count(status, mpi_nmsh_zone, recv_size, ierr)
177
178 do j = 1, recv_size
179 call new_zone_dist(pe_rank)%push(recv_buf_zone(j))
180 end do
181
182 call mpi_sendrecv(new_curve_dist(dst)%array(), &
183 new_curve_dist(dst)%size(), mpi_nmsh_curve, dst, 2, recv_buf_curve,&
184 max_recv(3), mpi_nmsh_curve, src, 2, neko_comm, status, ierr)
185 call mpi_get_count(status, mpi_nmsh_curve, recv_size, ierr)
186
187 do j = 1, recv_size
188 call new_curve_dist(pe_rank)%push(recv_buf_curve(j))
189 end do
190 end do
191
192 deallocate(recv_buf_msh)
193 deallocate(recv_buf_zone)
194 deallocate(recv_buf_curve)
195
196 do i = 0, pe_size - 1
197 if (i .ne. pe_rank) then
198 call new_mesh_dist(i)%free
199 call new_zone_dist(i)%free
200 call new_curve_dist(i)%free
201 end if
202 end do
203
204 !
205 ! Create a new mesh based on the given distribution
206 !
207 call msh%init(gdim, new_mesh_dist(pe_rank)%size())
208
209 call el_map%init(new_mesh_dist(pe_rank)%size())
210 call glb_map%init(new_mesh_dist(pe_rank)%size())
211
212 ! We should use the %array() procedure, which works great for
213 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
214 ! certain data types
215 select type (np => new_mesh_dist(pe_rank)%data)
216 type is (nmsh_hex_t)
217 do i = 1, new_mesh_dist(pe_rank)%size()
218 do j = 1, 8
219 call p(j)%init(np(i)%v(j)%v_xyz, np(i)%v(j)%v_idx)
220 end do
221 call msh%add_element(i, np(i)%el_idx, &
222 p(1), p(2), p(3), p(4), p(5), p(6), p(7), p(8))
223
224 if (el_map%get(np(i)%el_idx, tmp) .gt. 0) then
225 ! Old glb to new local
226 tmp = i
227 call el_map%set(np(i)%el_idx, tmp)
228
229 ! Old glb to new glb
230 tmp = msh%elements(i)%e%id()
231 call glb_map%set(np(i)%el_idx, tmp)
232 else
233 call neko_error('Global element id already defined')
234 end if
235 end do
236 end select
237 call new_mesh_dist(pe_rank)%free()
238 if (allocated(new_mesh_dist)) then
239 deallocate(new_mesh_dist)
240 end if
241
242
243 !
244 ! Figure out new mesh distribution (necessary for periodic zones)
245 !
246 call pe_lst%init()
247
248 ! We should use the %array() procedure, which works great for
249 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
250 ! certain data types
251 select type(zp => new_zone_dist(pe_rank)%data)
252 type is (nmsh_zone_t)
253 do i = 1, new_zone_dist(pe_rank)%size()
254 if (zp(i)%type .eq. 5) then
255 call pe_lst%push(zp(i)%p_e)
256 end if
257 end do
258 end select
259
260 max_recv_idx = 2 * pe_lst%size()
261 call mpi_allreduce(mpi_in_place, max_recv_idx, 1, mpi_integer, &
262 mpi_max, neko_comm, ierr)
263 allocate(recv_buf_idx(max_recv_idx))
264 allocate(send_buf_idx(max_recv_idx))
265
266 do i = 1, pe_size - 1
267 src = modulo(pe_rank - i + pe_size, pe_size)
268 dst = modulo(pe_rank + i, pe_size)
269
270 ! We should use the %array() procedure, which works great for
271 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
272 ! certain data types
273 select type (pe_lst_array => pe_lst%data)
274 type is (integer)
275 call mpi_sendrecv(pe_lst_array, &
276 pe_lst%size(), mpi_integer, dst, 0, recv_buf_idx, &
277 max_recv_idx, mpi_integer, src, 0, neko_comm, status, ierr)
278 end select
279 call mpi_get_count(status, mpi_integer, recv_size, ierr)
280
281 k = 0
282 do j = 1, recv_size
283 if (glb_map%get(recv_buf_idx(j), tmp) .eq. 0) then
284 k = k + 1
285 send_buf_idx(k) = recv_buf_idx(j)
286 k = k + 1
287 send_buf_idx(k) = tmp
288 end if
289 end do
290
291 call mpi_sendrecv(send_buf_idx, k, mpi_integer, src, 1, &
292 recv_buf_idx, max_recv_idx, mpi_integer, dst, 1, &
293 neko_comm, status, ierr)
294 call mpi_get_count(status, mpi_integer, recv_size, ierr)
295
296 do j = 1, recv_size, 2
297 call glb_map%set(recv_buf_idx(j), recv_buf_idx(j+1))
298 end do
299 end do
300 deallocate(recv_buf_idx)
301 deallocate(send_buf_idx)
302 call pe_lst%free()
303
304 !
305 ! Add zone data for new mesh distribution
306 !
307 select type (zp => new_zone_dist(pe_rank)%data)
308 type is (nmsh_zone_t)
309 do i = 1, new_zone_dist(pe_rank)%size()
310 if (el_map%get(zp(i)%e, new_el_idx) .gt. 0) then
311 call neko_error('Missing element after redistribution')
312 end if
313 select case(zp(i)%type)
314 case(5)
315 if (glb_map%get(zp(i)%p_e, new_pel_idx) .gt. 0) then
316 call neko_error('Missing periodic element after redistribution')
317 end if
318
319 call msh%mark_periodic_facet(zp(i)%f, new_el_idx, &
320 zp(i)%p_f, zp(i)%p_e, zp(i)%glb_pt_ids)
321 case(7)
322 call msh%mark_labeled_facet(zp(i)%f, new_el_idx, zp(i)%p_f)
323 end select
324 end do
325 do i = 1, new_zone_dist(pe_rank)%size()
326 if (el_map%get(zp(i)%e, new_el_idx) .gt. 0) then
327 call neko_error('Missing element after redistribution')
328 end if
329 select case(zp(i)%type)
330 case(5)
331 if (glb_map%get(zp(i)%p_e, new_pel_idx) .gt. 0) then
332 call neko_error('Missing periodic element after redistribution')
333 end if
334
335 call msh%apply_periodic_facet(zp(i)%f, new_el_idx, &
336 zp(i)%p_f, zp(i)%p_e, zp(i)%glb_pt_ids)
337 end select
338 end do
339 end select
340 call new_zone_dist(pe_rank)%free()
341 if (allocated(new_zone_dist)) then
342 deallocate(new_zone_dist)
343 end if
344
345
346 !
347 ! Add curve element information for new mesh distribution
348 !
349 select type (cp => new_curve_dist(pe_rank)%data)
350 type is (nmsh_curve_el_t)
351 do i = 1, new_curve_dist(pe_rank)%size()
352 if (el_map%get(cp(i)%e, new_el_idx) .gt. 0) then
353 call neko_error('Missing element after redistribution')
354 end if
355 call msh%mark_curve_element(new_el_idx, cp(i)%curve_data, cp(i)%type)
356 end do
357 end select
358 call new_curve_dist(pe_rank)%free()
359 if (allocated(new_curve_dist)) then
360 deallocate(new_curve_dist)
361 end if
362
363
364 call msh%finalize()
365
366 end subroutine redist_mesh
367
369 subroutine redist_zone(msh, z, type, parts, new_dist, label)
370 type(mesh_t), intent(inout) :: msh
371 class(facet_zone_t), intent(in) :: z
372 integer, intent(in) :: type
373 type(mesh_fld_t), intent(in) :: parts
374 type(stack_nz_t), intent(inout), allocatable :: new_dist(:)
375 integer, intent(in), optional :: label
376 type(nmsh_zone_t) :: nmsh_zone
377 integer :: i, j, zone_el, lbl
378
379 if (present(label)) then
380 lbl = label
381 else
382 lbl = 0
383 end if
384
385 select type(zp => z)
386 type is (facet_zone_periodic_t)
387 do i = 1, zp%size
388 zone_el = zp%facet_el(i)%x(2)
389 nmsh_zone%e = msh%elements(zp%facet_el(i)%x(2))%e%id()
390 nmsh_zone%f = zp%facet_el(i)%x(1)
391 nmsh_zone%p_e = zp%p_facet_el(i)%x(2)
392 nmsh_zone%p_f = zp%p_facet_el(i)%x(1)
393 nmsh_zone%glb_pt_ids = zp%p_ids(i)%x
394 nmsh_zone%type = type
395 call new_dist(parts%data(zone_el))%push(nmsh_zone)
396 end do
397 type is (facet_zone_t)
398 do i = 1, zp%size
399 zone_el = zp%facet_el(i)%x(2)
400 nmsh_zone%e = msh%elements(zp%facet_el(i)%x(2))%e%id()
401 nmsh_zone%f = zp%facet_el(i)%x(1)
402 nmsh_zone%p_f = lbl ! Labels are encoded in the periodic facet...
403 nmsh_zone%type = type
404 call new_dist(parts%data(zone_el))%push(nmsh_zone)
405 end do
406 end select
407 end subroutine redist_zone
408
410 subroutine redist_curve(msh, c, parts, new_dist)
411 type(mesh_t), intent(inout) :: msh
412 type(curve_t), intent(in) :: c
413 type(mesh_fld_t), intent(in) :: parts
414 type(stack_nc_t), intent(inout), allocatable :: new_dist(:)
415 type(nmsh_curve_el_t) :: nmsh_curve
416 integer :: i, j, curve_el
417
418 do i = 1, c%size
419 curve_el = c%curve_el(i)%el_idx
420 nmsh_curve%e = msh%elements(curve_el)%e%id()
421 nmsh_curve%curve_data = c%curve_el(i)%curve_data
422 nmsh_curve%type = c%curve_el(i)%curve_type
423 call new_dist(parts%data(curve_el))%push(nmsh_curve)
424 end do
425
426 end subroutine redist_curve
427
428end module redist
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Defines a domain as a subset of facets in a mesh.
Definition curve.f90:2
Defines a zone as a subset of facets in a mesh.
Implements a hash table ADT.
Definition htable.f90:36
Defines a mesh field.
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:63
MPI derived types.
Definition mpi_types.f90:34
type(mpi_datatype), public mpi_nmsh_hex
MPI derived type for 3D Neko nmsh data.
Definition mpi_types.f90:47
type(mpi_datatype), public mpi_nmsh_curve
MPI derived type for Neko nmsh curved elements.
Definition mpi_types.f90:50
type(mpi_datatype), public mpi_nmsh_zone
MPI derived type for Neko nmsh zone data.
Definition mpi_types.f90:49
Neko binary mesh format.
Definition nmsh.f90:2
Implements a point.
Definition point.f90:35
Redistribution routines.
Definition redist.f90:34
subroutine, public redist_mesh(msh, parts)
Redistribute a mesh msh according to new partitions.
Definition redist.f90:59
subroutine redist_zone(msh, z, type, parts, new_dist, label)
Fill redistribution list for zone data.
Definition redist.f90:370
subroutine redist_curve(msh, c, parts, new_dist)
Fill redistribution list for zone data.
Definition redist.f90:411
Implements a dynamic stack ADT.
Definition stack.f90:35
Utilities.
Definition utils.f90:35
Base type for an element.
Definition element.f90:44
Integer based hash table.
Definition htable.f90:82
Neko curve data.
Definition nmsh.f90:39
Neko hex element data.
Definition nmsh.f90:24
Neko zone data.
Definition nmsh.f90:29
A point in with coordinates .
Definition point.f90:43
Integer based stack.
Definition stack.f90:63
Neko curve info based stack.
Definition stack.f90:140
Neko hex element based stack.
Definition stack.f90:126
Neko zone based stack.
Definition stack.f90:133
#define max(a, b)
Definition tensor.cu:40