Neko  0.9.0
A portable framework for high-order spectral element flow simulations
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 !
34 module redist
36  use neko_mpi_types
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 
52 contains
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  class(nmsh_hex_t), pointer :: np(:)
63  type(nmsh_hex_t), allocatable :: recv_buf_msh(:)
64  class(nmsh_zone_t), pointer :: zp(:)
65  type(nmsh_zone_t), allocatable :: recv_buf_zone(:)
66  type(nmsh_curve_el_t), pointer :: cp(:)
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%wall, 1, parts, new_zone_dist)
94  call redist_zone(msh, msh%inlet, 2, parts, new_zone_dist)
95  call redist_zone(msh, msh%outlet, 3, parts, new_zone_dist)
96  call redist_zone(msh, msh%sympln, 4, parts, new_zone_dist)
97  call redist_zone(msh, msh%periodic, 5, parts, new_zone_dist)
98 
99  do j = 1, neko_msh_max_zlbls
100  label = j
101  call redist_zone(msh, msh%labeled_zones(j), 7, parts, &
102  new_zone_dist, label)
103  end do
104 
105  !
106  ! Extract new curve info. distributions
107  !
108  allocate(new_curve_dist(0:pe_size - 1))
109  do i = 0, pe_size - 1
110  call new_curve_dist(i)%init()
111  end do
112 
113  call redist_curve(msh, msh%curve, parts, new_curve_dist)
114 
115  !
116  ! Extract new mesh distribution
117  !
118 
119  allocate(new_mesh_dist(0:(pe_size - 1)))
120  do i = 0, pe_size - 1
121  call new_mesh_dist(i)%init()
122  end do
123 
124  do i = 1, msh%nelv
125  ep => msh%elements(i)%e
126  el%el_idx = ep%id()
127  do j = 1, 8
128  el%v(j)%v_idx = ep%pts(j)%p%id()
129  el%v(j)%v_xyz = ep%pts(j)%p%x
130  end do
131  call new_mesh_dist(parts%data(i))%push(el)
132  end do
133 
134 
135  gdim = msh%gdim
136  call msh%free()
137 
138  max_recv = 0
139  do i = 0, pe_size - 1
140  max_recv(1) = max(max_recv(1), new_mesh_dist(i)%size())
141  max_recv(2) = max(max_recv(2), new_zone_dist(i)%size())
142  max_recv(3) = max(max_recv(3), new_curve_dist(i)%size())
143  end do
144 
145  call mpi_allreduce(mpi_in_place, max_recv, 3, mpi_integer, &
146  mpi_max, neko_comm, ierr)
147  allocate(recv_buf_msh(max_recv(1)))
148  allocate(recv_buf_zone(max_recv(2)))
149  allocate(recv_buf_curve(max_recv(3)))
150 
151  do i = 1, pe_size - 1
152  src = modulo(pe_rank - i + pe_size, pe_size)
153  dst = modulo(pe_rank + i, pe_size)
154 
155  ! We should use the %array() procedure, which works great for
156  ! GNU, Intel and NEC, but it breaks horribly on Cray when using
157  ! certain data types
158  select type (nmd_array => new_mesh_dist(dst)%data)
159  type is (nmsh_hex_t)
160  call mpi_sendrecv(nmd_array, &
161  new_mesh_dist(dst)%size(), mpi_nmsh_hex, dst, 0, recv_buf_msh, &
162  max_recv(1), mpi_nmsh_hex, src, 0, neko_comm, status, ierr)
163  end select
164  call mpi_get_count(status, mpi_nmsh_hex, recv_size, ierr)
165 
166  do j = 1, recv_size
167  call new_mesh_dist(pe_rank)%push(recv_buf_msh(j))
168  end do
169 
170  ! We should use the %array() procedure, which works great for
171  ! GNU, Intel and NEC, but it breaks horribly on Cray when using
172  ! certain data types
173  select type (nzd_array => new_zone_dist(dst)%data)
174  type is (nmsh_zone_t)
175  call mpi_sendrecv(nzd_array, &
176  new_zone_dist(dst)%size(), mpi_nmsh_zone, dst, 1, recv_buf_zone,&
177  max_recv(2), mpi_nmsh_zone, src, 1, neko_comm, status, ierr)
178  end select
179  call mpi_get_count(status, mpi_nmsh_zone, recv_size, ierr)
180 
181  do j = 1, recv_size
182  call new_zone_dist(pe_rank)%push(recv_buf_zone(j))
183  end do
184 
185  call mpi_sendrecv(new_curve_dist(dst)%array(), &
186  new_curve_dist(dst)%size(), mpi_nmsh_curve, dst, 2, recv_buf_curve,&
187  max_recv(3), mpi_nmsh_curve, src, 2, neko_comm, status, ierr)
188  call mpi_get_count(status, mpi_nmsh_curve, recv_size, ierr)
189 
190  do j = 1, recv_size
191  call new_curve_dist(pe_rank)%push(recv_buf_curve(j))
192  end do
193  end do
194 
195  deallocate(recv_buf_msh)
196  deallocate(recv_buf_zone)
197  deallocate(recv_buf_curve)
198 
199  do i = 0, pe_size - 1
200  if (i .ne. pe_rank) then
201  call new_mesh_dist(i)%free
202  call new_zone_dist(i)%free
203  call new_curve_dist(i)%free
204  end if
205  end do
206 
207  !
208  ! Create a new mesh based on the given distribution
209  !
210  call msh%init(gdim, new_mesh_dist(pe_rank)%size())
211 
212  call el_map%init(new_mesh_dist(pe_rank)%size())
213  call glb_map%init(new_mesh_dist(pe_rank)%size())
214 
215  ! We should use the %array() procedure, which works great for
216  ! GNU, Intel and NEC, but it breaks horribly on Cray when using
217  ! certain data types
218  select type (np => new_mesh_dist(pe_rank)%data)
219  type is (nmsh_hex_t)
220  do i = 1, new_mesh_dist(pe_rank)%size()
221  do j = 1, 8
222  p(j) = point_t(np(i)%v(j)%v_xyz, np(i)%v(j)%v_idx)
223  end do
224  call msh%add_element(i, &
225  p(1), p(2), p(3), p(4), p(5), p(6), p(7), p(8))
226 
227  if (el_map%get(np(i)%el_idx, tmp) .gt. 0) then
228  ! Old glb to new local
229  tmp = i
230  call el_map%set(np(i)%el_idx, tmp)
231 
232  ! Old glb to new glb
233  tmp = msh%elements(i)%e%id()
234  call glb_map%set(np(i)%el_idx, tmp)
235  else
236  call neko_error('Global element id already defined')
237  end if
238  end do
239  end select
240  call new_mesh_dist(pe_rank)%free()
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  zp => new_zone_dist(pe_rank)%array()
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, new_pel_idx, zp(i)%glb_pt_ids)
328  case(7)
329  call msh%mark_labeled_facet(zp(i)%f, new_el_idx, zp(i)%p_f)
330  end select
331  end do
332  do i = 1, new_zone_dist(pe_rank)%size()
333  if (el_map%get(zp(i)%e, new_el_idx) .gt. 0) then
334  call neko_error('Missing element after redistribution')
335  end if
336  select case(zp(i)%type)
337  case(5)
338  if (glb_map%get(zp(i)%p_e, new_pel_idx) .gt. 0) then
339  call neko_error('Missing periodic element after redistribution')
340  end if
341 
342  call msh%apply_periodic_facet(zp(i)%f, new_el_idx, &
343  zp(i)%p_f, new_pel_idx, zp(i)%glb_pt_ids)
344  end select
345  end do
346 
347  call new_zone_dist(pe_rank)%free()
348 
349  !
350  ! Add curve element information for new mesh distribution
351  !
352  cp => new_curve_dist(pe_rank)%array()
353  do i = 1, new_curve_dist(pe_rank)%size()
354  if (el_map%get(cp(i)%e, new_el_idx) .gt. 0) then
355  call neko_error('Missing element after redistribution')
356  end if
357  call msh%mark_curve_element(new_el_idx, cp(i)%curve_data, cp(i)%type)
358  end do
359  call new_curve_dist(pe_rank)%free()
360 
361  call msh%finalize()
362 
363  end subroutine redist_mesh
364 
366  subroutine redist_zone(msh, z, type, parts, new_dist, label)
367  type(mesh_t), intent(inout) :: msh
368  class(facet_zone_t), intent(in) :: z
369  integer, intent(in) :: type
370  type(mesh_fld_t), intent(in) :: parts
371  type(stack_nz_t), intent(inout), allocatable :: new_dist(:)
372  integer, intent(in), optional :: label
373  type(nmsh_zone_t) :: nmsh_zone
374  integer :: i, j, zone_el, lbl
375 
376  if (present(label)) then
377  lbl = label
378  else
379  lbl = 0
380  end if
381 
382  select type(zp => z)
383  type is (facet_zone_periodic_t)
384  do i = 1, zp%size
385  zone_el = zp%facet_el(i)%x(2)
386  nmsh_zone%e = zp%facet_el(i)%x(2) + msh%offset_el
387  nmsh_zone%f = zp%facet_el(i)%x(1)
388  nmsh_zone%p_e = zp%p_facet_el(i)%x(2)
389  nmsh_zone%p_f = zp%p_facet_el(i)%x(1)
390  nmsh_zone%glb_pt_ids = zp%p_ids(i)%x
391  nmsh_zone%type = type
392  call new_dist(parts%data(zone_el))%push(nmsh_zone)
393  end do
394  type is (facet_zone_t)
395  do i = 1, zp%size
396  zone_el = zp%facet_el(i)%x(2)
397  nmsh_zone%e = zp%facet_el(i)%x(2) + msh%offset_el
398  nmsh_zone%f = zp%facet_el(i)%x(1)
399  nmsh_zone%p_f = lbl ! Labels are encoded in the periodic facet...
400  nmsh_zone%type = type
401  call new_dist(parts%data(zone_el))%push(nmsh_zone)
402  end do
403  end select
404  end subroutine redist_zone
405 
407  subroutine redist_curve(msh, c, parts, new_dist)
408  type(mesh_t), intent(inout) :: msh
409  type(curve_t), intent(in) :: c
410  type(mesh_fld_t), intent(in) :: parts
411  type(stack_nc_t), intent(inout), allocatable :: new_dist(:)
412  type(nmsh_curve_el_t) :: nmsh_curve
413  integer :: i, j, curve_el
414 
415  do i = 1, c%size
416  curve_el = c%curve_el(i)%el_idx
417  nmsh_curve%e = curve_el + msh%offset_el
418  nmsh_curve%curve_data = c%curve_el(i)%curve_data
419  nmsh_curve%type = c%curve_el(i)%curve_type
420  call new_dist(parts%data(curve_el))%push(nmsh_curve)
421  end do
422 
423  end subroutine redist_curve
424 
425 end module redist
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:28
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
integer pe_size
MPI size of communicator.
Definition: comm.F90:31
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.
Definition: facet_zone.f90:34
Implements a hash table ADT.
Definition: htable.f90:36
Defines a mesh field.
Definition: mesh_field.f90:35
subroutine, public mesh_field_free(fld)
Definition: mesh_field.f90:73
subroutine, public mesh_field_init(fld, msh, fld_name)
Definition: mesh_field.f90:52
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:367
subroutine redist_curve(msh, c, parts, new_dist)
Fill redistribution list for zone data.
Definition: redist.f90:408
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