Neko  0.8.99
A portable framework for high-order spectral element flow simulations
nmsh_file.f90
Go to the documentation of this file.
1 ! Copyright (c) 2019-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 nmsh_file
35  use generic_file, only: generic_file_t
36  use comm, only: neko_comm, pe_rank, pe_size
37  use num_types, only: rp, dp, i4, i8
38  use mesh, only: mesh_t, neko_msh_max_zlbls
39  use utils, only: neko_error
40  use point, only: point_t
41  use tuple, only: tuple4_i4_t
43  use element, only: element_t
44  use datadist, only: linear_dist_t
47  use mpi_f08, only: mpi_wtime, mpi_status, mpi_file, mpi_offset_kind, &
48  mpi_mode_wronly, mpi_mode_create, mpi_mode_rdonly, mpi_info_null, &
49  mpi_file_open, mpi_file_close, mpi_file_read_all, mpi_file_write_all, &
50  mpi_file_write_at_all, mpi_file_read_at_all, mpi_integer, mpi_sum, &
51  mpi_exscan, mpi_barrier, mpi_type_size
52  use logger, only: neko_log, log_size
53  implicit none
54 
55  private
58  integer, parameter :: max_write_nel = 8000000
60  type, public, extends(generic_file_t) :: nmsh_file_t
61  contains
62  procedure :: read => nmsh_file_read
63  procedure :: write => nmsh_file_write
64  end type nmsh_file_t
65 
66 contains
67 
69  subroutine nmsh_file_read(this, data)
70  class(nmsh_file_t) :: this
71  class(*), target, intent(inout) :: data
72  type(nmsh_hex_t), allocatable :: nmsh_hex(:)
73  type(nmsh_quad_t), allocatable :: nmsh_quad(:)
74  type(nmsh_zone_t), allocatable :: nmsh_zone(:)
75  type(nmsh_curve_el_t), allocatable :: nmsh_curve(:)
76  type(mesh_t), pointer :: msh
77  type(mpi_status) :: status
78  type(mpi_file) :: fh
79  integer (kind=MPI_OFFSET_KIND) :: mpi_offset, mpi_el_offset
80  integer :: i, j, ierr, element_offset
81  integer :: nmsh_quad_size, nmsh_hex_size, nmsh_zone_size
82  integer :: nelv, gdim, nzones, ncurves
83  integer :: el_idx
84  type(point_t), target :: p(8)
85  type(linear_dist_t) :: dist
86  character(len=LOG_SIZE) :: log_buf
87  real(kind=rp) :: t_start, t_end
88 
89  call this%check_exists()
90 
91  select type (data)
92  type is (mesh_t)
93  msh => data
94  class default
95  call neko_error('Invalid output data')
96  end select
97 
98 
99  call neko_log%section("Mesh")
100  call neko_log%message('Reading a binary Neko file ' // this%fname)
101 
102  call mpi_type_size(mpi_nmsh_hex, nmsh_hex_size, ierr)
103  call mpi_type_size(mpi_nmsh_quad, nmsh_quad_size, ierr)
104  call mpi_type_size(mpi_nmsh_zone, nmsh_zone_size, ierr)
105 
106  call mpi_file_open(neko_comm, trim(this%fname), &
107  mpi_mode_rdonly, mpi_info_null, fh, ierr)
108 
109  if (ierr > 0) then
110  call neko_error('Could not open the mesh file ' // this%fname // &
111  'for reading!')
112  end if
113  call mpi_file_read_all(fh, nelv, 1, mpi_integer, status, ierr)
114  call mpi_file_read_all(fh, gdim, 1, mpi_integer, status, ierr)
115 
116  write(log_buf,1) gdim, nelv
117 1 format('gdim = ', i1, ', nelements = ', i9)
118  call neko_log%message(log_buf)
119 
120  if (gdim .eq. 2) then
121  call mpi_file_close(fh, ierr)
122  call nmsh_file_read_2d(this, msh)
123  else
124 
125  dist = linear_dist_t(nelv, pe_rank, pe_size, neko_comm)
126  nelv = dist%num_local()
127  element_offset = dist%start_idx()
128 
129  call msh%init(gdim, nelv)
130 
131  call neko_log%message('Reading elements')
132 
133  if (msh%gdim .eq. 2) then
134  allocate(nmsh_quad(msh%nelv))
135  mpi_offset = int(2 * mpi_integer_size, i8) + &
136  int(element_offset, i8) * int(nmsh_quad_size, i8)
137  call mpi_file_read_at_all(fh, mpi_offset, &
138  nmsh_quad, msh%nelv, mpi_nmsh_quad, status, ierr)
139  do i = 1, nelv
140  do j = 1, 4
141  p(j) = point_t(nmsh_quad(i)%v(j)%v_xyz, nmsh_quad(i)%v(j)%v_idx)
142  end do
143  ! swap vertices to keep symmetric vertex numbering in neko
144  call msh%add_element(i, p(1), p(2), p(4), p(3))
145  end do
146  deallocate(nmsh_quad)
147  mpi_el_offset = int(2 * mpi_integer_size, i8) + &
148  int(dist%num_global(), i8) * int(nmsh_quad_size, i8)
149  else if (msh%gdim .eq. 3) then
150  allocate(nmsh_hex(msh%nelv))
151  mpi_offset = int(2 * mpi_integer_size, i8) + &
152  int(element_offset, i8) * int(nmsh_hex_size, i8)
153  call mpi_file_read_at_all(fh, mpi_offset, &
154  nmsh_hex, msh%nelv, mpi_nmsh_hex, status, ierr)
155  do i = 1, nelv
156  do j = 1, 8
157  p(j) = point_t(nmsh_hex(i)%v(j)%v_xyz, nmsh_hex(i)%v(j)%v_idx)
158  end do
159  ! swap vertices to keep symmetric vertex numbering in neko
160  call msh%add_element(i, &
161  p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
162  end do
163  deallocate(nmsh_hex)
164  mpi_el_offset = int(2 * mpi_integer_size, i8) + &
165  int(dist%num_global(), i8) * int(nmsh_hex_size, i8)
166  else
167  if (pe_rank .eq. 0) call neko_error('Invalid dimension of mesh')
168  end if
169  call neko_log%message('Reading BC/zone data')
170 
171  mpi_offset = mpi_el_offset
172  call mpi_file_read_at_all(fh, mpi_offset, &
173  nzones, 1, mpi_integer, status, ierr)
174  if (nzones .gt. 0) then
175  allocate(nmsh_zone(nzones))
176 
182  mpi_offset = mpi_el_offset + int(mpi_integer_size, i8)
183  call mpi_file_read_at_all(fh, mpi_offset, &
184  nmsh_zone, nzones, mpi_nmsh_zone, status, ierr)
185 
186  do i = 1, nzones
187  el_idx = nmsh_zone(i)%e
188  if (el_idx .gt. msh%offset_el .and. &
189  el_idx .le. msh%offset_el + msh%nelv) then
190  el_idx = el_idx - msh%offset_el
191  select case (nmsh_zone(i)%type)
192  case (1)
193  call msh%mark_wall_facet(nmsh_zone(i)%f, el_idx)
194  case (2)
195  call msh%mark_inlet_facet(nmsh_zone(i)%f, el_idx)
196  case (3)
197  call msh%mark_outlet_facet(nmsh_zone(i)%f, el_idx)
198  case (4)
199  call msh%mark_sympln_facet(nmsh_zone(i)%f, el_idx)
200  case (5)
201  call msh%mark_periodic_facet(nmsh_zone(i)%f, el_idx, &
202  nmsh_zone(i)%p_f, nmsh_zone(i)%p_e, &
203  nmsh_zone(i)%glb_pt_ids)
204  case (6)
205  call msh%mark_outlet_normal_facet(nmsh_zone(i)%f, el_idx)
206  case (7)
207  call msh%mark_labeled_facet(nmsh_zone(i)%f, el_idx, &
208  nmsh_zone(i)%p_f)
209  end select
210  end if
211  end do
212  !Apply facets, important that marking is finished
213  do i = 1, nzones
214  el_idx = nmsh_zone(i)%e
215  if (el_idx .gt. msh%offset_el .and. &
216  el_idx .le. msh%offset_el + msh%nelv) then
217  el_idx = el_idx - msh%offset_el
218  select case (nmsh_zone(i)%type)
219  case (5)
220  call msh%apply_periodic_facet(nmsh_zone(i)%f, el_idx, &
221  nmsh_zone(i)%p_f, nmsh_zone(i)%p_e, &
222  nmsh_zone(i)%glb_pt_ids)
223  end select
224  end if
225  end do
226 
227  deallocate(nmsh_zone)
228  end if
229  call neko_log%message('Reading deformation data')
230 
231  mpi_offset = mpi_el_offset + int(mpi_integer_size, i8) + &
232  int(nzones, i8)*int(nmsh_zone_size, i8)
233  call mpi_file_read_at_all(fh, mpi_offset, &
234  ncurves, 1, mpi_integer, status, ierr)
235 
236  if (ncurves .gt. 0) then
237 
238  allocate(nmsh_curve(ncurves))
239  mpi_offset = mpi_el_offset + int(2 * mpi_integer_size, i8) + &
240  int(nzones, i8)*int(nmsh_zone_size, i8)
241  call mpi_file_read_at_all(fh, mpi_offset, &
242  nmsh_curve, ncurves, mpi_nmsh_curve, status, ierr)
243 
244  do i = 1, ncurves
245  el_idx = nmsh_curve(i)%e - msh%offset_el
246  if (el_idx .gt. 0 .and. &
247  el_idx .le. msh%nelv) then
248  call msh%mark_curve_element(el_idx, &
249  nmsh_curve(i)%curve_data, nmsh_curve(i)%type)
250  end if
251 
252  end do
253 
254  deallocate(nmsh_curve)
255  end if
256 
257  call mpi_file_close(fh, ierr)
258  call neko_log%message('Mesh read, setting up connectivity')
259 
260  t_start = mpi_wtime()
261  call msh%finalize()
262  call mpi_barrier(neko_comm, ierr)
263  t_end = mpi_wtime()
264  write(log_buf, '(A)') 'Done setting up mesh and connectivity'
265  call neko_log%message(log_buf)
266  write(log_buf, '(A,F9.6)') &
267  'Mesh and connectivity setup (excluding read) time (s): ', &
268  t_end - t_start
269  call neko_log%message(log_buf)
270 
271  call neko_log%end_section()
272  end if
273 
274  end subroutine nmsh_file_read
275 
277  subroutine nmsh_file_read_2d(this, msh)
278  class(nmsh_file_t) :: this
279  type(mesh_t), pointer, intent(inout) :: msh
280  type(nmsh_hex_t), allocatable :: nmsh_hex(:)
281  type(nmsh_quad_t), allocatable :: nmsh_quad(:)
282  type(nmsh_zone_t), allocatable :: nmsh_zone(:)
283  type(nmsh_curve_el_t), allocatable :: nmsh_curve(:)
284  type(mpi_status) :: status
285  type(mpi_file) :: fh
286  integer (kind=MPI_OFFSET_KIND) :: mpi_offset, mpi_el_offset
287  integer :: i, j, ierr, element_offset, id
288  integer :: nmsh_quad_size, nmsh_hex_size, nmsh_zone_size
289  integer :: nelv, gdim, nzones, ncurves, ids(4)
290  integer :: el_idx
291  type(point_t) :: p(8)
292  type(linear_dist_t) :: dist
293  character(len=LOG_SIZE) :: log_buf
294  real(kind=rp) :: depth = 1d0
295  real(kind=dp) :: coord(3)
296  type(tuple4_i4_t) :: glb_pt_ids
297 
298 
299  call mpi_type_size(mpi_nmsh_hex, nmsh_hex_size, ierr)
300  call mpi_type_size(mpi_nmsh_quad, nmsh_quad_size, ierr)
301  call mpi_type_size(mpi_nmsh_zone, nmsh_zone_size, ierr)
302 
303  call mpi_file_open(neko_comm, trim(this%fname), &
304  mpi_mode_rdonly, mpi_info_null, fh, ierr)
305  call mpi_file_read_all(fh, nelv, 1, mpi_integer, status, ierr)
306  call mpi_file_read_all(fh, gdim, 1, mpi_integer, status, ierr)
307 
308  write(log_buf,2) gdim
309 2 format('gdim = ', i1, ', no full 2d support, creating thin slab')
310  call neko_log%message(log_buf)
311  gdim = 3
312 
313  dist = linear_dist_t(nelv, pe_rank, pe_size, neko_comm)
314  nelv = dist%num_local()
315  element_offset = dist%start_idx()
316 
317  call msh%init(gdim, nelv)
318 
319  allocate(nmsh_quad(msh%nelv))
320  mpi_offset = int(2 * mpi_integer_size, i8) + &
321  int(element_offset, i8) * int(nmsh_quad_size, i8)
322  call mpi_file_read_at_all(fh, mpi_offset, &
323  nmsh_quad, msh%nelv, mpi_nmsh_quad, status, ierr)
324  do i = 1, nelv
325  do j = 1, 4
326  coord = nmsh_quad(i)%v(j)%v_xyz
327  coord(3) = 0_rp
328  p(j) = point_t(coord, nmsh_quad(i)%v(j)%v_idx)
329  end do
330  do j = 1, 4
331  coord = nmsh_quad(i)%v(j)%v_xyz
332  coord(3) = depth
333  id = nmsh_quad(i)%v(j)%v_idx+msh%glb_nelv*8
334  p(j+4) = point_t(coord, id)
335  end do
336  ! swap vertices to keep symmetric vertex numbering in neko
337  call msh%add_element(i, &
338  p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
339  end do
340  deallocate(nmsh_quad)
341  mpi_el_offset = int(2 * mpi_integer_size, i8) + &
342  int(dist%num_global(), i8) * int(nmsh_quad_size, i8)
343 
344  mpi_offset = mpi_el_offset
345  call mpi_file_read_at_all(fh, mpi_offset, &
346  nzones, 1, mpi_integer, status, ierr)
347  if (nzones .gt. 0) then
348  allocate(nmsh_zone(nzones))
349 
355  mpi_offset = mpi_el_offset + int(mpi_integer_size, i8)
356  call mpi_file_read_at_all(fh, mpi_offset, &
357  nmsh_zone, nzones, mpi_nmsh_zone, status, ierr)
358 
359  do i = 1, nzones
360  el_idx = nmsh_zone(i)%e
361  if (el_idx .gt. msh%offset_el .and. &
362  el_idx .le. msh%offset_el + msh%nelv) then
363  el_idx = el_idx - msh%offset_el
364  select case (nmsh_zone(i)%type)
365  case (1)
366  call msh%mark_wall_facet(nmsh_zone(i)%f, el_idx)
367  case (2)
368  call msh%mark_inlet_facet(nmsh_zone(i)%f, el_idx)
369  case (3)
370  call msh%mark_outlet_facet(nmsh_zone(i)%f, el_idx)
371  case (4)
372  call msh%mark_sympln_facet(nmsh_zone(i)%f, el_idx)
373  case (5)
374  nmsh_zone(i)%glb_pt_ids(3) = nmsh_zone(i)%glb_pt_ids(1) + &
375  msh%glb_nelv * 8
376  nmsh_zone(i)%glb_pt_ids(4) = nmsh_zone(i)%glb_pt_ids(2) + &
377  msh%glb_nelv * 8
378  if (nmsh_zone(i)%f .eq. 1 .or. nmsh_zone(i)%f .eq. 2) then
379  ids(1) = nmsh_zone(i)%glb_pt_ids(1)
380  ids(2) = nmsh_zone(i)%glb_pt_ids(3)
381  ids(3) = nmsh_zone(i)%glb_pt_ids(4)
382  ids(4) = nmsh_zone(i)%glb_pt_ids(2)
383  else
384  ids(1) = nmsh_zone(i)%glb_pt_ids(1)
385  ids(2) = nmsh_zone(i)%glb_pt_ids(2)
386  ids(3) = nmsh_zone(i)%glb_pt_ids(4)
387  ids(4) = nmsh_zone(i)%glb_pt_ids(3)
388  end if
389  nmsh_zone(i)%glb_pt_ids = ids
390  call msh%mark_periodic_facet(nmsh_zone(i)%f, el_idx, &
391  nmsh_zone(i)%p_f, nmsh_zone(i)%p_e, ids)
392  case (6)
393  call msh%mark_outlet_normal_facet(nmsh_zone(i)%f, el_idx)
394  case (7)
395  call msh%mark_labeled_facet(nmsh_zone(i)%f, el_idx, &
396  nmsh_zone(i)%p_f)
397  end select
398  end if
399  end do
400  !Apply facets, important that marking is finished
401  do i = 1, nzones
402  el_idx = nmsh_zone(i)%e
403  if (el_idx .gt. msh%offset_el .and. &
404  el_idx .le. msh%offset_el + msh%nelv) then
405  el_idx = el_idx - msh%offset_el
406  select case (nmsh_zone(i)%type)
407  case (5)
408  call msh%apply_periodic_facet(nmsh_zone(i)%f, el_idx, &
409  nmsh_zone(i)%p_f, nmsh_zone(i)%p_e, &
410  nmsh_zone(i)%glb_pt_ids)
411  end select
412  end if
413  end do
414  !Do the same for extruded 3d points
415  do el_idx = 1, nelv
416  call msh%elements(el_idx)%e%facet_order(glb_pt_ids,5)
417  call msh%mark_periodic_facet(6, el_idx, &
418  5, el_idx, glb_pt_ids%x)
419  call msh%elements(el_idx)%e%facet_order(glb_pt_ids,5)
420  call msh%mark_periodic_facet(5, el_idx, &
421  6, el_idx, glb_pt_ids%x)
422  end do
423  do el_idx = 1, nelv
424  call msh%elements(el_idx)%e%facet_order(glb_pt_ids,5)
425  call msh%apply_periodic_facet(6, el_idx, &
426  5, el_idx, glb_pt_ids%x)
427  call msh%elements(el_idx)%e%facet_order(glb_pt_ids,5)
428  call msh%apply_periodic_facet(5, el_idx, &
429  6, el_idx, glb_pt_ids%x)
430  end do
431 
432  deallocate(nmsh_zone)
433  end if
434 
435  mpi_offset = mpi_el_offset + &
436  int(mpi_integer_size, i8) + int(nzones, i8)*int(nmsh_zone_size, i8)
437  call mpi_file_read_at_all(fh, mpi_offset, &
438  ncurves, 1, mpi_integer, status, ierr)
439 
440  if (ncurves .gt. 0) then
441 
442  allocate(nmsh_curve(ncurves))
443  mpi_offset = mpi_el_offset + &
444  int(2*mpi_integer_size, i8) + &
445  int(nzones, i8)*int(nmsh_zone_size, i8)
446  call mpi_file_read_at_all(fh, mpi_offset, &
447  nmsh_curve, ncurves, mpi_nmsh_curve, status, ierr)
448 
449  do i = 1, ncurves
450  el_idx = nmsh_curve(i)%e - msh%offset_el
451  if (el_idx .gt. 0 .and. &
452  el_idx .le. msh%nelv) then
453  call msh%mark_curve_element(el_idx, &
454  nmsh_curve(i)%curve_data, nmsh_curve(i)%type)
455  end if
456 
457  end do
458 
459  deallocate(nmsh_curve)
460  end if
461 
462  call mpi_file_close(fh, ierr)
463 
464  call msh%finalize()
465 
466  call neko_log%end_section()
467 
468  end subroutine nmsh_file_read_2d
469 
470 
472  subroutine nmsh_file_write(this, data, t)
473  class(nmsh_file_t), intent(inout) :: this
474  class(*), target, intent(in) :: data
475  real(kind=rp), intent(in), optional :: t
476  type(nmsh_quad_t), allocatable :: nmsh_quad(:)
477  type(nmsh_hex_t), allocatable :: nmsh_hex(:)
478  type(nmsh_zone_t), allocatable :: nmsh_zone(:)
479  type(nmsh_curve_el_t), allocatable :: nmsh_curve(:)
480  type(mesh_t), pointer :: msh
481  type(mpi_status) :: status
482  type(mpi_file) :: fh
483  integer (kind=MPI_OFFSET_KIND) :: mpi_offset, mpi_el_offset
484  integer :: i, j, ierr, nelgv, element_offset, k
485  integer :: nmsh_quad_size, nmsh_hex_size, nmsh_zone_size, nmsh_curve_size
486  integer :: nzones, ncurves
487  class(element_t), pointer :: ep
488  integer(i4), dimension(8), parameter :: vcyc_to_sym = [1, 2, 4, 3, 5, &
489  & 6, 8, 7] ! cyclic to symmetric vertex mapping
490 
491  select type (data)
492  type is (mesh_t)
493  msh => data
494  class default
495  call neko_error('Invalid output data')
496  end select
497 
498  call mpi_type_size(mpi_nmsh_quad, nmsh_quad_size, ierr)
499  call mpi_type_size(mpi_nmsh_hex, nmsh_hex_size, ierr)
500  call mpi_type_size(mpi_nmsh_zone, nmsh_zone_size, ierr)
501  call mpi_type_size(mpi_nmsh_curve, nmsh_curve_size, ierr)
502 
503  call mpi_reduce(msh%nelv, nelgv, 1, mpi_integer, &
504  mpi_sum, 0, neko_comm, ierr)
505  element_offset = 0
506  call mpi_exscan(msh%nelv, element_offset, 1, &
507  mpi_integer, mpi_sum, neko_comm, ierr)
508 
509  call neko_log%message('Writing data as a binary Neko file ' // this%fname)
510 
511  call mpi_file_open(neko_comm, trim(this%fname), &
512  mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
513 
514  call mpi_file_write_all(fh, nelgv, 1, mpi_integer, status, ierr)
515  call mpi_file_write_all(fh, msh%gdim, 1, mpi_integer, status, ierr)
516 
517  call msh%reset_periodic_ids()
518 
519  if (msh%gdim .eq. 2) then
520  allocate(nmsh_quad(msh%nelv))
521  do i = 1, msh%nelv
522  ep => msh%elements(i)%e
523  nmsh_quad(i)%el_idx = ep%id()
524  do j = 1, 4
525  nmsh_quad(i)%v(j)%v_idx = ep%pts(vcyc_to_sym(j))%p%id()
526  nmsh_quad(i)%v(j)%v_xyz = ep%pts(vcyc_to_sym(j))%p%x
527  end do
528  end do
529  mpi_offset = int(2 * mpi_integer_size, i8) + &
530  int(element_offset, i8) * int(nmsh_quad_size, i8)
531  call mpi_file_write_at_all(fh, mpi_offset, &
532  nmsh_quad, msh%nelv, mpi_nmsh_quad, status, ierr)
533  deallocate(nmsh_quad)
534  mpi_el_offset = int(2 * mpi_integer_size, i8) + &
535  int(nelgv, i8) * int(nmsh_quad_size, i8)
536  else if (msh%gdim .eq. 3) then
537  allocate(nmsh_hex(msh%nelv))
538  do i = 1, msh%nelv
539  ep => msh%elements(i)%e
540  nmsh_hex(i)%el_idx = ep%id()
541  do j = 1, 8
542  nmsh_hex(i)%v(j)%v_idx = ep%pts(vcyc_to_sym(j))%p%id()
543  nmsh_hex(i)%v(j)%v_xyz = ep%pts(vcyc_to_sym(j))%p%x
544  end do
545  end do
546  mpi_offset = int(2 * mpi_integer_size, i8) + &
547  int(element_offset, i8) * int(nmsh_hex_size, i8)
548  call mpi_file_write_at_all(fh, mpi_offset, &
549  nmsh_hex, min(msh%nelv, max_write_nel), mpi_nmsh_hex, status, ierr)
550  do i = 1, msh%nelv/max_write_nel
551  mpi_offset = int(2 * mpi_integer_size, i8) + &
552  int(element_offset+i*max_write_nel, i8) * int(nmsh_hex_size, i8)
553  call mpi_file_write_at_all(fh, mpi_offset, &
554  nmsh_hex(i*max_write_nel+1), &
555  min(msh%nelv-i*max_write_nel, max_write_nel), &
556  mpi_nmsh_hex, status, ierr)
557  end do
558  deallocate(nmsh_hex)
559  mpi_el_offset = int(2 * mpi_integer_size, i8) + &
560  int(nelgv, i8) * int(nmsh_hex_size, i8)
561  else
562  call neko_error('Invalid dimension of mesh')
563  end if
564 
565  nzones = msh%wall%size + msh%inlet%size + msh%outlet%size + &
566  msh%sympln%size + msh%periodic%size + msh%outlet_normal%size
567 
568  do i = 1, neko_msh_max_zlbls
569  nzones = nzones + msh%labeled_zones(i)%size
570  end do
571  mpi_offset = mpi_el_offset
572  call mpi_file_write_at_all(fh, mpi_offset, &
573  nzones, 1, mpi_integer, status, ierr)
574 
575  if (nzones .gt. 0) then
576  allocate(nmsh_zone(nzones))
577 
578  nmsh_zone(:)%type = 0
579 
580  j = 1
581  do i = 1, msh%wall%size
582  nmsh_zone(j)%e = msh%wall%facet_el(i)%x(2) + msh%offset_el
583  nmsh_zone(j)%f = msh%wall%facet_el(i)%x(1)
584  nmsh_zone(j)%type = 1
585  j = j + 1
586  end do
587 
588  do i = 1, msh%inlet%size
589  nmsh_zone(j)%e = msh%inlet%facet_el(i)%x(2) + msh%offset_el
590  nmsh_zone(j)%f = msh%inlet%facet_el(i)%x(1)
591  nmsh_zone(j)%type = 2
592  j = j + 1
593  end do
594 
595  do i = 1, msh%outlet%size
596  nmsh_zone(j)%e = msh%outlet%facet_el(i)%x(2) + msh%offset_el
597  nmsh_zone(j)%f = msh%outlet%facet_el(i)%x(1)
598  nmsh_zone(j)%type = 3
599  j = j + 1
600  end do
601 
602  do i = 1, msh%sympln%size
603  nmsh_zone(j)%e = msh%sympln%facet_el(i)%x(2) + msh%offset_el
604  nmsh_zone(j)%f = msh%sympln%facet_el(i)%x(1)
605  nmsh_zone(j)%type = 4
606  j = j + 1
607  end do
608 
609  do i = 1, msh%periodic%size
610  nmsh_zone(j)%e = msh%periodic%facet_el(i)%x(2) + msh%offset_el
611  nmsh_zone(j)%f = msh%periodic%facet_el(i)%x(1)
612  nmsh_zone(j)%p_e = msh%periodic%p_facet_el(i)%x(2)
613  nmsh_zone(j)%p_f = msh%periodic%p_facet_el(i)%x(1)
614  nmsh_zone(j)%glb_pt_ids = msh%periodic%p_ids(i)%x
615  nmsh_zone(j)%type = 5
616  j = j + 1
617  end do
618  do i = 1, msh%outlet_normal%size
619  nmsh_zone(j)%e = msh%outlet_normal%facet_el(i)%x(2) + msh%offset_el
620  nmsh_zone(j)%f = msh%outlet_normal%facet_el(i)%x(1)
621  nmsh_zone(j)%type = 6
622  j = j + 1
623  end do
624  do k = 1, neko_msh_max_zlbls
625  do i = 1, msh%labeled_zones(k)%size
626  nmsh_zone(j)%e = msh%labeled_zones(k)%facet_el(i)%x(2) + &
627  msh%offset_el
628  nmsh_zone(j)%f = msh%labeled_zones(k)%facet_el(i)%x(1)
629  nmsh_zone(j)%p_f = k
630  nmsh_zone(j)%type = 7
631  j = j + 1
632  end do
633  end do
634 
635 
636  mpi_offset = mpi_el_offset + int(mpi_integer_size, i8)
637  call mpi_file_write_at_all(fh, mpi_offset, &
638  nmsh_zone, nzones, mpi_nmsh_zone, status, ierr)
639 
640  deallocate(nmsh_zone)
641  end if
642 
643  ncurves = msh%curve%size
644  mpi_offset = mpi_el_offset + int(mpi_integer_size, i8) + &
645  int(nzones, i8)*int(nmsh_zone_size, i8)
646 
647  call mpi_file_write_at_all(fh, mpi_offset, &
648  ncurves, 1, mpi_integer, status, ierr)
649 
650  if (ncurves .gt. 0) then
651  allocate(nmsh_curve(ncurves))
652  do i = 1, ncurves
653  nmsh_curve(i)%type = 0
654  end do
655 
656  do i = 1, ncurves
657  nmsh_curve(i)%e = msh%curve%curve_el(i)%el_idx + msh%offset_el
658  nmsh_curve(i)%curve_data = msh%curve%curve_el(i)%curve_data
659  nmsh_curve(i)%type = msh%curve%curve_el(i)%curve_type
660  end do
661  mpi_offset = mpi_el_offset + int(2*mpi_integer_size, i8) + &
662  int(nzones, i8)*int(nmsh_zone_size, i8)
663  call mpi_file_write_at_all(fh, mpi_offset, &
664  nmsh_curve, ncurves, mpi_nmsh_curve, status, ierr)
665 
666  deallocate(nmsh_curve)
667  end if
668 
669  call mpi_file_sync(fh, ierr)
670  call mpi_file_close(fh, ierr)
671  call neko_log%message('Done')
672 
673  end subroutine nmsh_file_write
674 
675 end module nmsh_file
676 
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:26
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
integer pe_size
MPI size of communicator.
Definition: comm.F90:29
Defines practical data distributions.
Definition: datadist.f90:34
Logging routines.
Definition: log.f90:34
type(log_t), public neko_log
Global log stream.
Definition: log.f90:61
integer, parameter, public log_size
Definition: log.f90:40
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_quad
MPI derived type for 2D Neko nmsh data.
Definition: mpi_types.f90:43
type(mpi_datatype), public mpi_nmsh_curve
MPI derived type for Neko nmsh curved elements.
Definition: mpi_types.f90:45
integer, public mpi_integer_size
Size of MPI type integer.
Definition: mpi_types.f90:63
type(mpi_datatype), public mpi_nmsh_zone
MPI derived type for Neko nmsh zone data.
Definition: mpi_types.f90:44
Neko binary mesh data.
Definition: nmsh_file.f90:34
subroutine nmsh_file_write(this, data, t)
Write a mesh from to a binary Neko nmsh file.
Definition: nmsh_file.f90:473
subroutine nmsh_file_read(this, data)
Load a mesh from a binary Neko nmsh file.
Definition: nmsh_file.f90:70
integer, parameter max_write_nel
Specifices the maximum number of elements any rank is allowed to write (for nmsh)....
Definition: nmsh_file.f90:58
subroutine nmsh_file_read_2d(this, msh)
Load a mesh from a binary Neko nmsh file.
Definition: nmsh_file.f90:278
Neko binary mesh format.
Definition: nmsh.f90:2
integer, parameter, public i8
Definition: num_types.f90:7
integer, parameter, public i4
Definition: num_types.f90:6
integer, parameter, public dp
Definition: num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Implements a point.
Definition: point.f90:35
Implements a n-tuple.
Definition: tuple.f90:34
Utilities.
Definition: utils.f90:35
Load-balanced linear distribution .
Definition: datadist.f90:50
Base type for an element.
Definition: element.f90:44
A generic file handler.
Neko curve data.
Definition: nmsh.f90:39
Neko hex element data.
Definition: nmsh.f90:24
Neko quad element data.
Definition: nmsh.f90:19
Neko zone data.
Definition: nmsh.f90:29
Interface for Neko nmsh files.
Definition: nmsh_file.f90:60
A point in with coordinates .
Definition: point.f90:43
Integer based 4-tuple.
Definition: tuple.f90:69