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
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!
36 use comm, only: neko_comm, pe_rank, pe_size
37 use num_types, only: rp, dp, i4, i8
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, mpi_allreduce, mpi_file_sync
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
66contains
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, el_idx_glb
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
1171 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, nmsh_quad(i)%el_idx, &
145 p(1), p(2), p(4), p(3))
146 end do
147 deallocate(nmsh_quad)
148 mpi_el_offset = int(2 * mpi_integer_size, i8) + &
149 int(dist%num_global(), i8) * int(nmsh_quad_size, i8)
150 else if (msh%gdim .eq. 3) then
151 allocate(nmsh_hex(msh%nelv))
152 mpi_offset = int(2 * mpi_integer_size, i8) + &
153 int(element_offset, i8) * int(nmsh_hex_size, i8)
154 call mpi_file_read_at_all(fh, mpi_offset, &
155 nmsh_hex, msh%nelv, mpi_nmsh_hex, status, ierr)
156 do i = 1, nelv
157 do j = 1, 8
158 p(j) = point_t(nmsh_hex(i)%v(j)%v_xyz, nmsh_hex(i)%v(j)%v_idx)
159 end do
160 ! swap vertices to keep symmetric vertex numbering in neko
161 call msh%add_element(i, nmsh_hex(i)%el_idx, &
162 p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
163 end do
164 deallocate(nmsh_hex)
165 mpi_el_offset = int(2 * mpi_integer_size, i8) + &
166 int(dist%num_global(), i8) * int(nmsh_hex_size, i8)
167 else
168 if (pe_rank .eq. 0) call neko_error('Invalid dimension of mesh')
169 end if
170 call neko_log%message('Reading BC/zone data')
171
172 mpi_offset = mpi_el_offset
173 call mpi_file_read_at_all(fh, mpi_offset, &
174 nzones, 1, mpi_integer, status, ierr)
175 if (nzones .gt. 0) then
176 allocate(nmsh_zone(nzones))
177
183 mpi_offset = mpi_el_offset + int(mpi_integer_size, i8)
184 call mpi_file_read_at_all(fh, mpi_offset, &
185 nmsh_zone, nzones, mpi_nmsh_zone, status, ierr)
186
187 do i = 1, nzones
188 el_idx_glb = nmsh_zone(i)%e
189 if (msh%htel%get(el_idx_glb, el_idx) .eq. 0) then
190 select case (nmsh_zone(i)%type)
191 case (5)
192 call msh%mark_periodic_facet(nmsh_zone(i)%f, el_idx, &
193 nmsh_zone(i)%p_f, nmsh_zone(i)%p_e, &
194 nmsh_zone(i)%glb_pt_ids)
195 case (7)
196 call msh%mark_labeled_facet(nmsh_zone(i)%f, el_idx, &
197 nmsh_zone(i)%p_f)
198 end select
199 end if
200 end do
201 !Apply facets, important that marking is finished
202 do i = 1, nzones
203 el_idx_glb = nmsh_zone(i)%e
204 if (msh%htel%get(el_idx_glb, el_idx) .eq. 0) then
205 select case (nmsh_zone(i)%type)
206 case (5)
207 call msh%apply_periodic_facet(nmsh_zone(i)%f, el_idx, &
208 nmsh_zone(i)%p_f, nmsh_zone(i)%p_e, &
209 nmsh_zone(i)%glb_pt_ids)
210 end select
211 end if
212 end do
213
214 deallocate(nmsh_zone)
215 end if
216 call neko_log%message('Reading deformation data')
217
218 mpi_offset = mpi_el_offset + int(mpi_integer_size, i8) + &
219 int(nzones, i8)*int(nmsh_zone_size, i8)
220 call mpi_file_read_at_all(fh, mpi_offset, &
221 ncurves, 1, mpi_integer, status, ierr)
222
223 if (ncurves .gt. 0) then
224
225 allocate(nmsh_curve(ncurves))
226 mpi_offset = mpi_el_offset + int(2 * mpi_integer_size, i8) + &
227 int(nzones, i8)*int(nmsh_zone_size, i8)
228 call mpi_file_read_at_all(fh, mpi_offset, &
229 nmsh_curve, ncurves, mpi_nmsh_curve, status, ierr)
230
231 do i = 1, ncurves
232 el_idx_glb = nmsh_curve(i)%e
233 if (msh%htel%get(el_idx_glb, el_idx) .eq. 0) then
234 call msh%mark_curve_element(el_idx, &
235 nmsh_curve(i)%curve_data, nmsh_curve(i)%type)
236 end if
237
238 end do
239
240 deallocate(nmsh_curve)
241 end if
242
243 call mpi_file_close(fh, ierr)
244 call neko_log%message('Mesh read, setting up connectivity')
245
246 t_start = mpi_wtime()
247 call msh%finalize()
248 call mpi_barrier(neko_comm, ierr)
249 t_end = mpi_wtime()
250 write(log_buf, '(A)') 'Done setting up mesh and connectivity'
251 call neko_log%message(log_buf)
252 write(log_buf, '(A,F9.6)') &
253 'Mesh and connectivity setup (excluding read) time (s): ', &
254 t_end - t_start
255 call neko_log%message(log_buf)
256
257 call neko_log%end_section()
258 end if
259
260 end subroutine nmsh_file_read
261
263 subroutine nmsh_file_read_2d(this, msh)
264 class(nmsh_file_t) :: this
265 type(mesh_t), pointer, intent(inout) :: msh
266 type(nmsh_hex_t), allocatable :: nmsh_hex(:)
267 type(nmsh_quad_t), allocatable :: nmsh_quad(:)
268 type(nmsh_zone_t), allocatable :: nmsh_zone(:)
269 type(nmsh_curve_el_t), allocatable :: nmsh_curve(:)
270 type(mpi_status) :: status
271 type(mpi_file) :: fh
272 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, mpi_el_offset
273 integer :: i, j, ierr, element_offset, id
274 integer :: nmsh_quad_size, nmsh_hex_size, nmsh_zone_size
275 integer :: nelv, gdim, nzones, ncurves, ids(4)
276 integer :: el_idx_glb, el_idx
277 type(point_t) :: p(8)
278 type(linear_dist_t) :: dist
279 character(len=LOG_SIZE) :: log_buf
280 real(kind=rp) :: depth = 1d0
281 real(kind=dp) :: coord(3)
282 type(tuple4_i4_t) :: glb_pt_ids
283
284
285 call mpi_type_size(mpi_nmsh_hex, nmsh_hex_size, ierr)
286 call mpi_type_size(mpi_nmsh_quad, nmsh_quad_size, ierr)
287 call mpi_type_size(mpi_nmsh_zone, nmsh_zone_size, ierr)
288
289 call mpi_file_open(neko_comm, trim(this%fname), &
290 mpi_mode_rdonly, mpi_info_null, fh, ierr)
291 call mpi_file_read_all(fh, nelv, 1, mpi_integer, status, ierr)
292 call mpi_file_read_all(fh, gdim, 1, mpi_integer, status, ierr)
293
294 write(log_buf,2) gdim
2952 format('gdim = ', i1, ', no full 2d support, creating thin slab')
296 call neko_log%message(log_buf)
297 gdim = 3
298
299 dist = linear_dist_t(nelv, pe_rank, pe_size, neko_comm)
300 nelv = dist%num_local()
301 element_offset = dist%start_idx()
302
303 call msh%init(gdim, nelv)
304
305 allocate(nmsh_quad(msh%nelv))
306 mpi_offset = int(2 * mpi_integer_size, i8) + &
307 int(element_offset, i8) * int(nmsh_quad_size, i8)
308 call mpi_file_read_at_all(fh, mpi_offset, &
309 nmsh_quad, msh%nelv, mpi_nmsh_quad, status, ierr)
310 do i = 1, nelv
311 do j = 1, 4
312 coord = nmsh_quad(i)%v(j)%v_xyz
313 coord(3) = 0_rp
314 p(j) = point_t(coord, nmsh_quad(i)%v(j)%v_idx)
315 end do
316 do j = 1, 4
317 coord = nmsh_quad(i)%v(j)%v_xyz
318 coord(3) = depth
319 id = nmsh_quad(i)%v(j)%v_idx+msh%glb_nelv*8
320 p(j+4) = point_t(coord, id)
321 end do
322 ! swap vertices to keep symmetric vertex numbering in neko
323 call msh%add_element(i, nmsh_quad(i)%el_idx, &
324 p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
325 end do
326 deallocate(nmsh_quad)
327 mpi_el_offset = int(2 * mpi_integer_size, i8) + &
328 int(dist%num_global(), i8) * int(nmsh_quad_size, i8)
329
330 mpi_offset = mpi_el_offset
331 call mpi_file_read_at_all(fh, mpi_offset, &
332 nzones, 1, mpi_integer, status, ierr)
333 if (nzones .gt. 0) then
334 allocate(nmsh_zone(nzones))
335
341 mpi_offset = mpi_el_offset + int(mpi_integer_size, i8)
342 call mpi_file_read_at_all(fh, mpi_offset, &
343 nmsh_zone, nzones, mpi_nmsh_zone, status, ierr)
344
345 do i = 1, nzones
346 el_idx_glb = nmsh_zone(i)%e
347 if (msh%htel%get(el_idx_glb, el_idx) .eq. 0) then
348 select case (nmsh_zone(i)%type)
349 case (5)
350 nmsh_zone(i)%glb_pt_ids(3) = nmsh_zone(i)%glb_pt_ids(1) + &
351 msh%glb_nelv * 8
352 nmsh_zone(i)%glb_pt_ids(4) = nmsh_zone(i)%glb_pt_ids(2) + &
353 msh%glb_nelv * 8
354 if (nmsh_zone(i)%f .eq. 1 .or. nmsh_zone(i)%f .eq. 2) then
355 ids(1) = nmsh_zone(i)%glb_pt_ids(1)
356 ids(2) = nmsh_zone(i)%glb_pt_ids(3)
357 ids(3) = nmsh_zone(i)%glb_pt_ids(4)
358 ids(4) = nmsh_zone(i)%glb_pt_ids(2)
359 else
360 ids(1) = nmsh_zone(i)%glb_pt_ids(1)
361 ids(2) = nmsh_zone(i)%glb_pt_ids(2)
362 ids(3) = nmsh_zone(i)%glb_pt_ids(4)
363 ids(4) = nmsh_zone(i)%glb_pt_ids(3)
364 end if
365 nmsh_zone(i)%glb_pt_ids = ids
366 call msh%mark_periodic_facet(nmsh_zone(i)%f, el_idx, &
367 nmsh_zone(i)%p_f, nmsh_zone(i)%p_e, ids)
368 case (7)
369 call msh%mark_labeled_facet(nmsh_zone(i)%f, el_idx, &
370 nmsh_zone(i)%p_f)
371 end select
372 end if
373 end do
374 !Apply facets, important that marking is finished
375 do i = 1, nzones
376 el_idx_glb = nmsh_zone(i)%e
377 if (msh%htel%get(el_idx_glb, el_idx) .eq. 0) then
378 select case (nmsh_zone(i)%type)
379 case (5)
380 call msh%apply_periodic_facet(nmsh_zone(i)%f, el_idx, &
381 nmsh_zone(i)%p_f, nmsh_zone(i)%p_e, &
382 nmsh_zone(i)%glb_pt_ids)
383 end select
384 end if
385 end do
386 !Do the same for extruded 3d points
387 do el_idx = 1, nelv
388 call msh%elements(el_idx)%e%facet_order(glb_pt_ids,5)
389 call msh%mark_periodic_facet(6, el_idx, &
390 5, el_idx, glb_pt_ids%x)
391 call msh%elements(el_idx)%e%facet_order(glb_pt_ids,5)
392 call msh%mark_periodic_facet(5, el_idx, &
393 6, el_idx, glb_pt_ids%x)
394 end do
395 do el_idx = 1, nelv
396 call msh%elements(el_idx)%e%facet_order(glb_pt_ids,5)
397 call msh%apply_periodic_facet(6, el_idx, &
398 5, el_idx, glb_pt_ids%x)
399 call msh%elements(el_idx)%e%facet_order(glb_pt_ids,5)
400 call msh%apply_periodic_facet(5, el_idx, &
401 6, el_idx, glb_pt_ids%x)
402 end do
403
404 deallocate(nmsh_zone)
405 end if
406
407 mpi_offset = mpi_el_offset + &
408 int(mpi_integer_size, i8) + int(nzones, i8)*int(nmsh_zone_size, i8)
409 call mpi_file_read_at_all(fh, mpi_offset, &
410 ncurves, 1, mpi_integer, status, ierr)
411
412 if (ncurves .gt. 0) then
413
414 allocate(nmsh_curve(ncurves))
415 mpi_offset = mpi_el_offset + &
416 int(2*mpi_integer_size, i8) + &
417 int(nzones, i8)*int(nmsh_zone_size, i8)
418 call mpi_file_read_at_all(fh, mpi_offset, &
419 nmsh_curve, ncurves, mpi_nmsh_curve, status, ierr)
420
421 do i = 1, ncurves
422 el_idx_glb = nmsh_curve(i)%e
423 if (msh%htel%get(el_idx_glb, el_idx) .eq. 0) then
424 call msh%mark_curve_element(el_idx, &
425 nmsh_curve(i)%curve_data, nmsh_curve(i)%type)
426 end if
427 end do
428
429 deallocate(nmsh_curve)
430 end if
431
432 call mpi_file_close(fh, ierr)
433
434 call msh%finalize()
435
436 call neko_log%end_section()
437
438 end subroutine nmsh_file_read_2d
439
440
442 subroutine nmsh_file_write(this, data, t)
443 class(nmsh_file_t), intent(inout) :: this
444 class(*), target, intent(in) :: data
445 real(kind=rp), intent(in), optional :: t
446 type(nmsh_quad_t), allocatable :: nmsh_quad(:)
447 type(nmsh_hex_t), allocatable :: nmsh_hex(:)
448 type(nmsh_zone_t), allocatable :: nmsh_zone(:)
449 type(nmsh_curve_el_t), allocatable :: nmsh_curve(:)
450 type(mesh_t), pointer :: msh
451 type(mpi_status) :: status
452 type(mpi_file) :: fh
453 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, mpi_el_offset
454 integer :: i, j, ierr, k
455 integer :: nmsh_quad_size, nmsh_hex_size, nmsh_zone_size, nmsh_curve_size
456 integer :: nzones, nzones_glb, nzones_offset
457 integer :: ncurves, ncurves_glb, ncurves_offset
458 integer :: el_idx, el_idx_glb
459 class(element_t), pointer :: ep
460 integer(i4), dimension(8), parameter :: vcyc_to_sym = &
461 [1, 2, 4, 3, 5, 6, 8, 7] ! cyclic to symmetric vertex mapping
462
463 select type (data)
464 type is (mesh_t)
465 msh => data
466 class default
467 call neko_error('Invalid output data')
468 end select
469
470 call mpi_type_size(mpi_nmsh_quad, nmsh_quad_size, ierr)
471 call mpi_type_size(mpi_nmsh_hex, nmsh_hex_size, ierr)
472 call mpi_type_size(mpi_nmsh_zone, nmsh_zone_size, ierr)
473 call mpi_type_size(mpi_nmsh_curve, nmsh_curve_size, ierr)
474
475 call neko_log%message('Writing data as a binary Neko file ' // this%fname)
476
477 call mpi_file_open(neko_comm, trim(this%fname), &
478 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
479
480 call mpi_file_write_all(fh, msh%glb_nelv, 1, mpi_integer, status, ierr)
481 call mpi_file_write_all(fh, msh%gdim, 1, mpi_integer, status, ierr)
482
483 call msh%reset_periodic_ids()
484
485 if (msh%gdim .eq. 2) then
486 allocate(nmsh_quad(msh%nelv))
487 do i = 1, msh%nelv
488 ep => msh%elements(i)%e
489 nmsh_quad(i)%el_idx = ep%id()
490 do j = 1, 4
491 nmsh_quad(i)%v(j)%v_idx = ep%pts(vcyc_to_sym(j))%p%id()
492 nmsh_quad(i)%v(j)%v_xyz = ep%pts(vcyc_to_sym(j))%p%x
493 end do
494 end do
495 mpi_offset = int(2 * mpi_integer_size, i8) + &
496 int(msh%offset_el, i8) * int(nmsh_quad_size, i8)
497 call mpi_file_write_at_all(fh, mpi_offset, &
498 nmsh_quad, msh%nelv, mpi_nmsh_quad, status, ierr)
499 deallocate(nmsh_quad)
500 mpi_el_offset = int(2 * mpi_integer_size, i8) + &
501 int(msh%glb_nelv, i8) * int(nmsh_quad_size, i8)
502 else if (msh%gdim .eq. 3) then
503 allocate(nmsh_hex(msh%nelv))
504 do i = 1, msh%nelv
505 ep => msh%elements(i)%e
506 nmsh_hex(i)%el_idx = ep%id()
507 do j = 1, 8
508 nmsh_hex(i)%v(j)%v_idx = ep%pts(vcyc_to_sym(j))%p%id()
509 nmsh_hex(i)%v(j)%v_xyz = ep%pts(vcyc_to_sym(j))%p%x
510 end do
511 end do
512 mpi_offset = int(2 * mpi_integer_size, i8) + &
513 int(msh%offset_el, i8) * int(nmsh_hex_size, i8)
514 call mpi_file_write_at_all(fh, mpi_offset, &
515 nmsh_hex, min(msh%nelv, max_write_nel), mpi_nmsh_hex, status, ierr)
516 do i = 1, msh%nelv/max_write_nel
517 mpi_offset = int(2 * mpi_integer_size, i8) + &
518 int(msh%offset_el+i*max_write_nel, i8) * int(nmsh_hex_size, i8)
519 call mpi_file_write_at_all(fh, mpi_offset, &
520 nmsh_hex(i*max_write_nel+1), &
521 min(msh%nelv-i*max_write_nel, max_write_nel), &
522 mpi_nmsh_hex, status, ierr)
523 end do
524 deallocate(nmsh_hex)
525 mpi_el_offset = int(2 * mpi_integer_size, i8) + &
526 int(msh%glb_nelv, i8) * int(nmsh_hex_size, i8)
527 else
528 call neko_error('Invalid dimension of mesh')
529 end if
530
531 nzones = msh%periodic%size
532 do i = 1, neko_msh_max_zlbls
533 nzones = nzones + msh%labeled_zones(i)%size
534 end do
535
536 call mpi_allreduce(nzones, nzones_glb, 1, &
537 mpi_integer, mpi_sum, neko_comm, ierr)
538
539 nzones_offset = 0
540 call mpi_exscan(nzones, nzones_offset, 1, &
541 mpi_integer, mpi_sum, neko_comm, ierr)
542
543 mpi_offset = mpi_el_offset
544 call mpi_file_write_at_all(fh, mpi_offset, &
545 nzones_glb, 1, mpi_integer, status, ierr)
546
547 if (nzones_glb .gt. 0) then
548 allocate(nmsh_zone(nzones))
549
550 if (nzones .gt. 0) then
551 nmsh_zone(:)%type = 0
552
553 j = 1
554 do i = 1, msh%periodic%size
555 nmsh_zone(j)%e = msh%elements(msh%periodic%facet_el(i)%x(2))%e%id()
556 nmsh_zone(j)%f = msh%periodic%facet_el(i)%x(1)
557 nmsh_zone(j)%p_e = msh%periodic%p_facet_el(i)%x(2)
558 nmsh_zone(j)%p_f = msh%periodic%p_facet_el(i)%x(1)
559 nmsh_zone(j)%glb_pt_ids = msh%periodic%p_ids(i)%x
560 nmsh_zone(j)%type = 5
561 j = j + 1
562 end do
563
564 do k = 1, neko_msh_max_zlbls
565 do i = 1, msh%labeled_zones(k)%size
566 nmsh_zone(j)%e = &
567 msh%elements(msh%labeled_zones(k)%facet_el(i)%x(2))%e%id()
568 nmsh_zone(j)%f = msh%labeled_zones(k)%facet_el(i)%x(1)
569 nmsh_zone(j)%p_f = k
570 nmsh_zone(j)%type = 7
571 j = j + 1
572 end do
573 end do
574 end if
575
576 mpi_offset = mpi_el_offset + int(mpi_integer_size, i8) + &
577 int(nzones_offset, i8) * int(nmsh_zone_size, i8)
578 call mpi_file_write_at_all(fh, mpi_offset, &
579 nmsh_zone, nzones, mpi_nmsh_zone, status, ierr)
580
581 deallocate(nmsh_zone)
582 end if
583
584 ncurves = msh%curve%size
585
586
587 call mpi_allreduce(ncurves, ncurves_glb, 1, &
588 mpi_integer, mpi_sum, neko_comm, ierr)
589
590 ncurves_offset = 0
591 call mpi_exscan(ncurves, ncurves_offset, 1, &
592 mpi_integer, mpi_sum, neko_comm, ierr)
593
594 mpi_offset = mpi_el_offset + int(mpi_integer_size, i8) + &
595 int(nzones_glb, i8)*int(nmsh_zone_size, i8)
596
597 call mpi_file_write_at_all(fh, mpi_offset, &
598 ncurves_glb, 1, mpi_integer, status, ierr)
599
600 if (ncurves_glb .gt. 0) then
601
602 allocate(nmsh_curve(ncurves))
603
604 do i = 1, ncurves
605 nmsh_curve(i)%type = 0
606 end do
607
608 do i = 1, ncurves
609 nmsh_curve(i)%e = msh%elements(msh%curve%curve_el(i)%el_idx)%e%id()
610 nmsh_curve(i)%curve_data = msh%curve%curve_el(i)%curve_data
611 nmsh_curve(i)%type = msh%curve%curve_el(i)%curve_type
612 end do
613
614 mpi_offset = mpi_el_offset + int(2*mpi_integer_size, i8) + &
615 int(nzones_glb, i8) * int(nmsh_zone_size, i8) + &
616 int(ncurves_offset, i8) * int(nmsh_curve_size, i8)
617
618 call mpi_file_write_at_all(fh, mpi_offset, &
619 nmsh_curve, ncurves, mpi_nmsh_curve, status, ierr)
620 deallocate(nmsh_curve)
621 end if
622
623 call mpi_file_sync(fh, ierr)
624 call mpi_file_close(fh, ierr)
625 call neko_log%message('Done')
626
627 !
628 ! Re-apply periodic facets
629 ! (necessary if the mesh is going to be used after I/O)
630 !
631 do i = 1, msh%periodic%size
632 el_idx_glb = msh%elements(msh%periodic%facet_el(i)%x(2))%e%id()
633 if (msh%htel%get(el_idx_glb, el_idx) .eq. 0) then
634 call msh%apply_periodic_facet(msh%periodic%facet_el(i)%x(1), el_idx, &
635 msh%periodic%p_facet_el(i)%x(1), msh%periodic%p_facet_el(i)%x(2), &
636 msh%periodic%p_ids(i)%x)
637 end if
638 end do
639
640 end subroutine nmsh_file_write
641
642end module nmsh_file
643
Definition comm.F90:1
integer pe_rank
MPI rank.
Definition comm.F90:51
type(mpi_comm) neko_comm
MPI communicator.
Definition comm.F90:38
integer pe_size
MPI size of communicator.
Definition comm.F90:54
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:65
integer, parameter, public log_size
Definition log.f90:42
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.
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.
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