Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
redist.f90
Go to the documentation of this file.
1! Copyright (c) 2021, 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
37 use mpi_f08
38 use htable, only : htable_i4_t
39 use point, only : point_t
41 use curve, only : curve_t
42 use comm
43 use mesh, only : mesh_t, neko_msh_max_zlbls
46 use element, only : element_t
47 implicit none
48 private
49
50 public :: redist_mesh
51
52contains
53
55 subroutine redist_mesh(msh, parts)
56 type(mesh_t), intent(inout), target :: msh
57 type(mesh_fld_t), intent(in) :: parts
58 type(stack_nh_t), allocatable :: new_mesh_dist(:)
59 type(stack_nz_t), allocatable :: new_zone_dist(:)
60 type(stack_nc_t), allocatable :: new_curve_dist(:)
61 type(nmsh_hex_t) :: el
62 type(nmsh_hex_t), allocatable :: recv_buf_msh(:)
63 type(nmsh_zone_t), allocatable :: recv_buf_zone(:)
64 type(nmsh_curve_el_t), allocatable :: recv_buf_curve(:)
65 class(element_t), pointer :: ep
66 integer, allocatable :: recv_buf_idx(:), send_buf_idx(:)
67 type(mpi_status) :: status
68 integer :: i, j, k, ierr, max_recv_idx, label
69 integer :: src, dst, recv_size, gdim, tmp, new_el_idx, new_pel_idx
70 integer :: max_recv(3)
71 type(point_t) :: p(8)
72 type(htable_i4_t) :: el_map, glb_map
73 type(stack_i4_t) :: pe_lst
74
75
76 !
77 ! Reset possible periodic ids
78 !
79 call msh%reset_periodic_ids()
80
81 !
82 ! Extract new zone distributions
83 !
84
85 allocate(new_zone_dist(0:pe_size - 1))
86 do i = 0, pe_size - 1
87 call new_zone_dist(i)%init()
88 end do
89
90 call redist_zone(msh, msh%wall, 1, parts, new_zone_dist)
91 call redist_zone(msh, msh%inlet, 2, parts, new_zone_dist)
92 call redist_zone(msh, msh%outlet, 3, parts, new_zone_dist)
93 call redist_zone(msh, msh%sympln, 4, parts, new_zone_dist)
94 call redist_zone(msh, msh%periodic, 5, parts, new_zone_dist)
95 call redist_zone(msh, msh%outlet_normal, 6, parts, new_zone_dist)
96
97 do j = 1, neko_msh_max_zlbls
98 label = j
99 call redist_zone(msh, msh%labeled_zones(j), 7, parts, &
100 new_zone_dist, label)
101 end do
102
103 !
104 ! Extract new curve info. distributions
105 !
106 allocate(new_curve_dist(0:pe_size - 1))
107 do i = 0, pe_size - 1
108 call new_curve_dist(i)%init()
109 end do
110
111 call redist_curve(msh, msh%curve, parts, new_curve_dist)
112
113 !
114 ! Extract new mesh distribution
115 !
116
117 allocate(new_mesh_dist(0:(pe_size - 1)))
118 do i = 0, pe_size - 1
119 call new_mesh_dist(i)%init()
120 end do
121
122 do i = 1, msh%nelv
123 ep => msh%elements(i)%e
124 el%el_idx = ep%id()
125 do j = 1, 8
126 el%v(j)%v_idx = ep%pts(j)%p%id()
127 el%v(j)%v_xyz = ep%pts(j)%p%x
128 end do
129 call new_mesh_dist(parts%data(i))%push(el)
130 end do
131
132
133 gdim = msh%gdim
134 call msh%free()
135
136 max_recv = 0
137 do i = 0, pe_size - 1
138 max_recv(1) = max(max_recv(1), new_mesh_dist(i)%size())
139 max_recv(2) = max(max_recv(2), new_zone_dist(i)%size())
140 max_recv(3) = max(max_recv(3), new_curve_dist(i)%size())
141 end do
142
143 call mpi_allreduce(mpi_in_place, max_recv, 3, mpi_integer, &
144 mpi_max, neko_comm, ierr)
145
146 allocate(recv_buf_msh(max_recv(1)))
147 allocate(recv_buf_zone(max_recv(2)))
148 allocate(recv_buf_curve(max_recv(3)))
149
150 do i = 1, pe_size - 1
151 src = modulo(pe_rank - i + pe_size, pe_size)
152 dst = modulo(pe_rank + i, pe_size)
153
154 ! We should use the %array() procedure, which works great for
155 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
156 ! certain data types
157 select type (nmd_array => new_mesh_dist(dst)%data)
158 type is (nmsh_hex_t)
159 call mpi_sendrecv(nmd_array, &
160 new_mesh_dist(dst)%size(), mpi_nmsh_hex, dst, 0, recv_buf_msh, &
161 max_recv(1), mpi_nmsh_hex, src, 0, neko_comm, status, ierr)
162 end select
163 call mpi_get_count(status, mpi_nmsh_hex, recv_size, ierr)
164
165 do j = 1, recv_size
166 call new_mesh_dist(pe_rank)%push(recv_buf_msh(j))
167 end do
168
169 ! We should use the %array() procedure, which works great for
170 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
171 ! certain data types
172 select type (nzd_array => new_zone_dist(dst)%data)
173 type is (nmsh_zone_t)
174 call mpi_sendrecv(nzd_array, &
175 new_zone_dist(dst)%size(), mpi_nmsh_zone, dst, 1, recv_buf_zone,&
176 max_recv(2), mpi_nmsh_zone, src, 1, neko_comm, status, ierr)
177 end select
178 call mpi_get_count(status, mpi_nmsh_zone, recv_size, ierr)
179
180 do j = 1, recv_size
181 call new_zone_dist(pe_rank)%push(recv_buf_zone(j))
182 end do
183
184 call mpi_sendrecv(new_curve_dist(dst)%array(), &
185 new_curve_dist(dst)%size(), mpi_nmsh_curve, dst, 2, recv_buf_curve,&
186 max_recv(3), mpi_nmsh_curve, src, 2, neko_comm, status, ierr)
187 call mpi_get_count(status, mpi_nmsh_curve, recv_size, ierr)
188
189 do j = 1, recv_size
190 call new_curve_dist(pe_rank)%push(recv_buf_curve(j))
191 end do
192 end do
193
194 deallocate(recv_buf_msh)
195 deallocate(recv_buf_zone)
196 deallocate(recv_buf_curve)
197
198 do i = 0, pe_size - 1
199 if (i .ne. pe_rank) then
200 call new_mesh_dist(i)%free
201 call new_zone_dist(i)%free
202 call new_curve_dist(i)%free
203 end if
204 end do
205
206 !
207 ! Create a new mesh based on the given distribution
208 !
209 call msh%init(gdim, new_mesh_dist(pe_rank)%size())
210
211 call el_map%init(new_mesh_dist(pe_rank)%size())
212 call glb_map%init(new_mesh_dist(pe_rank)%size())
213
214 ! We should use the %array() procedure, which works great for
215 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
216 ! certain data types
217 select type (np => new_mesh_dist(pe_rank)%data)
218 type is (nmsh_hex_t)
219 do i = 1, new_mesh_dist(pe_rank)%size()
220 do j = 1, 8
221 p(j) = point_t(np(i)%v(j)%v_xyz, np(i)%v(j)%v_idx)
222 end do
223 call msh%add_element(i, np(i)%el_idx, &
224 p(1), p(2), p(3), p(4), p(5), p(6), p(7), p(8))
225
226 if (el_map%get(np(i)%el_idx, tmp) .gt. 0) then
227 ! Old glb to new local
228 tmp = i
229 call el_map%set(np(i)%el_idx, tmp)
230
231 ! Old glb to new glb
232 tmp = msh%elements(i)%e%id()
233 call glb_map%set(np(i)%el_idx, tmp)
234 else
235 call neko_error('Global element id already defined')
236 end if
237 end do
238 end select
239 call new_mesh_dist(pe_rank)%free()
240
241
242 !
243 ! Figure out new mesh distribution (necessary for periodic zones)
244 !
245 call pe_lst%init()
246
247 ! We should use the %array() procedure, which works great for
248 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
249 ! certain data types
250 select type(zp => new_zone_dist(pe_rank)%data)
251 type is (nmsh_zone_t)
252 do i = 1, new_zone_dist(pe_rank)%size()
253 if (zp(i)%type .eq. 5) then
254 call pe_lst%push(zp(i)%p_e)
255 end if
256 end do
257 end select
258
259 max_recv_idx = 2 * pe_lst%size()
260 call mpi_allreduce(mpi_in_place, max_recv_idx, 1, mpi_integer, &
261 mpi_max, neko_comm, ierr)
262 allocate(recv_buf_idx(max_recv_idx))
263 allocate(send_buf_idx(max_recv_idx))
264
265 do i = 1, pe_size - 1
266 src = modulo(pe_rank - i + pe_size, pe_size)
267 dst = modulo(pe_rank + i, pe_size)
268
269 ! We should use the %array() procedure, which works great for
270 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
271 ! certain data types
272 select type (pe_lst_array => pe_lst%data)
273 type is (integer)
274 call mpi_sendrecv(pe_lst_array, &
275 pe_lst%size(), mpi_integer, dst, 0, recv_buf_idx, &
276 max_recv_idx, mpi_integer, src, 0, neko_comm, status, ierr)
277 end select
278 call mpi_get_count(status, mpi_integer, recv_size, ierr)
279
280 k = 0
281 do j = 1, recv_size
282 if (glb_map%get(recv_buf_idx(j), tmp) .eq. 0) then
283 k = k + 1
284 send_buf_idx(k) = recv_buf_idx(j)
285 k = k + 1
286 send_buf_idx(k) = tmp
287 end if
288 end do
289
290 call mpi_sendrecv(send_buf_idx, k, mpi_integer, src, 1, &
291 recv_buf_idx, max_recv_idx, mpi_integer, dst, 1, &
292 neko_comm, status, ierr)
293 call mpi_get_count(status, mpi_integer, recv_size, ierr)
294
295 do j = 1, recv_size, 2
296 call glb_map%set(recv_buf_idx(j), recv_buf_idx(j+1))
297 end do
298 end do
299 deallocate(recv_buf_idx)
300 deallocate(send_buf_idx)
301 call pe_lst%free()
302
303 !
304 ! Add zone data for new mesh distribution
305 !
306 select type (zp => new_zone_dist(pe_rank)%data)
307 type is (nmsh_zone_t)
308 do i = 1, new_zone_dist(pe_rank)%size()
309 if (el_map%get(zp(i)%e, new_el_idx) .gt. 0) then
310 call neko_error('Missing element after redistribution')
311 end if
312 select case(zp(i)%type)
313 case(1)
314 call msh%mark_wall_facet(zp(i)%f, new_el_idx)
315 case(2)
316 call msh%mark_inlet_facet(zp(i)%f, new_el_idx)
317 case(3)
318 call msh%mark_outlet_facet(zp(i)%f, new_el_idx)
319 case(4)
320 call msh%mark_sympln_facet(zp(i)%f, new_el_idx)
321 case(5)
322 if (glb_map%get(zp(i)%p_e, new_pel_idx) .gt. 0) then
323 call neko_error('Missing periodic element after redistribution')
324 end if
325
326 call msh%mark_periodic_facet(zp(i)%f, new_el_idx, &
327 zp(i)%p_f, zp(i)%p_e, zp(i)%glb_pt_ids)
328 case(6)
329 call msh%mark_outlet_normal_facet(zp(i)%f, new_el_idx)
330 case(7)
331 call msh%mark_labeled_facet(zp(i)%f, new_el_idx, zp(i)%p_f)
332 end select
333 end do
334 do i = 1, new_zone_dist(pe_rank)%size()
335 if (el_map%get(zp(i)%e, new_el_idx) .gt. 0) then
336 call neko_error('Missing element after redistribution')
337 end if
338 select case(zp(i)%type)
339 case(5)
340 if (glb_map%get(zp(i)%p_e, new_pel_idx) .gt. 0) then
341 call neko_error('Missing periodic element after redistribution')
342 end if
343
344 call msh%apply_periodic_facet(zp(i)%f, new_el_idx, &
345 zp(i)%p_f, zp(i)%p_e, zp(i)%glb_pt_ids)
346 end select
347 end do
348 end select
349 call new_zone_dist(pe_rank)%free()
350
351
352 !
353 ! Add curve element information for new mesh distribution
354 !
355 select type (cp => new_curve_dist(pe_rank)%data)
356 type is (nmsh_curve_el_t)
357 do i = 1, new_curve_dist(pe_rank)%size()
358 if (el_map%get(cp(i)%e, new_el_idx) .gt. 0) then
359 call neko_error('Missing element after redistribution')
360 end if
361 call msh%mark_curve_element(new_el_idx, cp(i)%curve_data, cp(i)%type)
362 end do
363 end select
364 call new_curve_dist(pe_rank)%free()
365
366
367 call msh%finalize()
368
369 end subroutine redist_mesh
370
372 subroutine redist_zone(msh, z, type, parts, new_dist, label)
373 type(mesh_t), intent(inout) :: msh
374 class(facet_zone_t), intent(in) :: z
375 integer, intent(in) :: type
376 type(mesh_fld_t), intent(in) :: parts
377 type(stack_nz_t), intent(inout), allocatable :: new_dist(:)
378 integer, intent(in), optional :: label
379 type(nmsh_zone_t) :: nmsh_zone
380 integer :: i, j, zone_el, lbl
381
382 if (present(label)) then
383 lbl = label
384 else
385 lbl = 0
386 end if
387
388 select type(zp => z)
389 type is (facet_zone_periodic_t)
390 do i = 1, zp%size
391 zone_el = zp%facet_el(i)%x(2)
392 nmsh_zone%e = msh%elements(zp%facet_el(i)%x(2))%e%id()
393 nmsh_zone%f = zp%facet_el(i)%x(1)
394 nmsh_zone%p_e = zp%p_facet_el(i)%x(2)
395 nmsh_zone%p_f = zp%p_facet_el(i)%x(1)
396 nmsh_zone%glb_pt_ids = zp%p_ids(i)%x
397 nmsh_zone%type = type
398 call new_dist(parts%data(zone_el))%push(nmsh_zone)
399 end do
400 type is (facet_zone_t)
401 do i = 1, zp%size
402 zone_el = zp%facet_el(i)%x(2)
403 nmsh_zone%e = msh%elements(zp%facet_el(i)%x(2))%e%id()
404 nmsh_zone%f = zp%facet_el(i)%x(1)
405 nmsh_zone%p_f = lbl ! Labels are encoded in the periodic facet...
406 nmsh_zone%type = type
407 call new_dist(parts%data(zone_el))%push(nmsh_zone)
408 end do
409 end select
410 end subroutine redist_zone
411
413 subroutine redist_curve(msh, c, parts, new_dist)
414 type(mesh_t), intent(inout) :: msh
415 type(curve_t), intent(in) :: c
416 type(mesh_fld_t), intent(in) :: parts
417 type(stack_nc_t), intent(inout), allocatable :: new_dist(:)
418 type(nmsh_curve_el_t) :: nmsh_curve
419 integer :: i, j, curve_el
420
421 do i = 1, c%size
422 curve_el = c%curve_el(i)%el_idx
423 nmsh_curve%e = msh%elements(curve_el)%e%id()
424 nmsh_curve%curve_data = c%curve_el(i)%curve_data
425 nmsh_curve%type = c%curve_el(i)%curve_type
426 call new_dist(parts%data(curve_el))%push(nmsh_curve)
427 end do
428
429 end subroutine redist_curve
430
431end module redist
Definition comm.F90:1
integer pe_rank
MPI rank.
Definition comm.F90:50
type(mpi_comm) neko_comm
MPI communicator.
Definition comm.F90:38
integer pe_size
MPI size of communicator.
Definition comm.F90:53
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.
subroutine, public mesh_field_free(fld)
subroutine, public mesh_field_init(fld, msh, fld_name)
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:56
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:42
type(mpi_datatype), public mpi_nmsh_curve
MPI derived type for Neko nmsh curved elements.
Definition mpi_types.f90:45
type(mpi_datatype), public mpi_nmsh_zone
MPI derived type for Neko nmsh zone data.
Definition mpi_types.f90:44
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:56
subroutine redist_zone(msh, z, type, parts, new_dist, label)
Fill redistribution list for zone data.
Definition redist.f90:373
subroutine redist_curve(msh, c, parts, new_dist)
Fill redistribution list for zone data.
Definition redist.f90:414
Implements a dynamic stack ADT.
Definition stack.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