71 class(*),
target,
intent(inout) :: data
76 type(
mesh_t),
pointer :: msh
77 type(mpi_status) :: status
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
86 character(len=LOG_SIZE) :: log_buf
87 real(kind=
rp) :: t_start, t_end
89 call this%check_exists()
100 call neko_log%message(
'Reading a binary Neko file ' // this%fname)
106 call mpi_file_open(
neko_comm, trim(this%fname), &
107 mpi_mode_rdonly, mpi_info_null, fh, ierr)
110 call neko_error(
'Could not open the mesh file ' // this%fname // &
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)
116 write(log_buf,1) gdim, nelv
1171
format(
'gdim = ', i1,
', nelements = ', i9)
120 if (gdim .eq. 2)
then
121 call mpi_file_close(fh, ierr)
126 nelv = dist%num_local()
127 element_offset = dist%start_idx()
129 call msh%init(gdim, nelv)
131 call neko_log%message(
'Reading elements')
133 if (msh%gdim .eq. 2)
then
134 allocate(nmsh_quad(msh%nelv))
136 int(element_offset,
i8) * int(nmsh_quad_size,
i8)
137 call mpi_file_read_at_all(fh, mpi_offset, &
141 p(j) =
point_t(nmsh_quad(i)%v(j)%v_xyz, nmsh_quad(i)%v(j)%v_idx)
144 call msh%add_element(i, nmsh_quad(i)%el_idx, &
145 p(1), p(2), p(4), p(3))
147 deallocate(nmsh_quad)
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))
153 int(element_offset,
i8) * int(nmsh_hex_size,
i8)
154 call mpi_file_read_at_all(fh, mpi_offset, &
158 p(j) =
point_t(nmsh_hex(i)%v(j)%v_xyz, nmsh_hex(i)%v(j)%v_idx)
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))
166 int(dist%num_global(),
i8) * int(nmsh_hex_size,
i8)
170 call neko_log%message(
'Reading BC/zone data')
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))
184 call mpi_file_read_at_all(fh, mpi_offset, &
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)
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)
196 call msh%mark_labeled_facet(nmsh_zone(i)%f, el_idx, &
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)
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)
214 deallocate(nmsh_zone)
216 call neko_log%message(
'Reading deformation data')
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)
223 if (ncurves .gt. 0)
then
225 allocate(nmsh_curve(ncurves))
227 int(nzones,
i8)*int(nmsh_zone_size,
i8)
228 call mpi_file_read_at_all(fh, mpi_offset, &
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)
240 deallocate(nmsh_curve)
243 call mpi_file_close(fh, ierr)
244 call neko_log%message(
'Mesh read, setting up connectivity')
246 t_start = mpi_wtime()
250 write(log_buf,
'(A)')
'Done setting up mesh and connectivity'
252 write(log_buf,
'(A,F9.6)') &
253 'Mesh and connectivity setup (excluding read) time (s): ', &
265 type(
mesh_t),
pointer,
intent(inout) :: msh
270 type(mpi_status) :: status
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
279 character(len=LOG_SIZE) :: log_buf
280 real(kind=
rp) :: depth = 1d0
281 real(kind=
dp) :: coord(3)
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)
294 write(log_buf,2) gdim
2952
format(
'gdim = ', i1,
', no full 2d support, creating thin slab')
300 nelv = dist%num_local()
301 element_offset = dist%start_idx()
303 call msh%init(gdim, nelv)
305 allocate(nmsh_quad(msh%nelv))
307 int(element_offset,
i8) * int(nmsh_quad_size,
i8)
308 call mpi_file_read_at_all(fh, mpi_offset, &
312 coord = nmsh_quad(i)%v(j)%v_xyz
314 p(j) =
point_t(coord, nmsh_quad(i)%v(j)%v_idx)
317 coord = nmsh_quad(i)%v(j)%v_xyz
319 id = nmsh_quad(i)%v(j)%v_idx+msh%glb_nelv*8
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))
326 deallocate(nmsh_quad)
328 int(dist%num_global(),
i8) * int(nmsh_quad_size,
i8)
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))
342 call mpi_file_read_at_all(fh, mpi_offset, &
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)
350 nmsh_zone(i)%glb_pt_ids(3) = nmsh_zone(i)%glb_pt_ids(1) + &
352 nmsh_zone(i)%glb_pt_ids(4) = nmsh_zone(i)%glb_pt_ids(2) + &
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)
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)
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)
369 call msh%mark_labeled_facet(nmsh_zone(i)%f, el_idx, &
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)
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)
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)
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)
404 deallocate(nmsh_zone)
407 mpi_offset = mpi_el_offset + &
409 call mpi_file_read_at_all(fh, mpi_offset, &
410 ncurves, 1, mpi_integer, status, ierr)
412 if (ncurves .gt. 0)
then
414 allocate(nmsh_curve(ncurves))
415 mpi_offset = mpi_el_offset + &
417 int(nzones,
i8)*int(nmsh_zone_size,
i8)
418 call mpi_file_read_at_all(fh, mpi_offset, &
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)
429 deallocate(nmsh_curve)
432 call mpi_file_close(fh, ierr)
444 class(*),
target,
intent(in) :: data
445 real(kind=
rp),
intent(in),
optional :: t
450 type(
mesh_t),
pointer :: msh
451 type(mpi_status) :: status
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
460 integer(i4),
dimension(8),
parameter :: vcyc_to_sym = &
461 [1, 2, 4, 3, 5, 6, 8, 7]
475 call neko_log%message(
'Writing data as a binary Neko file ' // this%fname)
477 call mpi_file_open(
neko_comm, trim(this%fname), &
478 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
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)
483 call msh%reset_periodic_ids()
485 if (msh%gdim .eq. 2)
then
486 allocate(nmsh_quad(msh%nelv))
488 ep => msh%elements(i)%e
489 nmsh_quad(i)%el_idx = ep%id()
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
496 int(msh%offset_el,
i8) * int(nmsh_quad_size,
i8)
497 call mpi_file_write_at_all(fh, mpi_offset, &
499 deallocate(nmsh_quad)
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))
505 ep => msh%elements(i)%e
506 nmsh_hex(i)%el_idx = ep%id()
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
513 int(msh%offset_el,
i8) * int(nmsh_hex_size,
i8)
514 call mpi_file_write_at_all(fh, mpi_offset, &
519 call mpi_file_write_at_all(fh, mpi_offset, &
526 int(msh%glb_nelv,
i8) * int(nmsh_hex_size,
i8)
531 nzones = msh%periodic%size
533 nzones = nzones + msh%labeled_zones(i)%size
536 call mpi_allreduce(nzones, nzones_glb, 1, &
540 call mpi_exscan(nzones, nzones_offset, 1, &
543 mpi_offset = mpi_el_offset
544 call mpi_file_write_at_all(fh, mpi_offset, &
545 nzones_glb, 1, mpi_integer, status, ierr)
547 if (nzones_glb .gt. 0)
then
548 allocate(nmsh_zone(nzones))
550 if (nzones .gt. 0)
then
551 nmsh_zone(:)%type = 0
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
565 do i = 1, msh%labeled_zones(k)%size
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)
570 nmsh_zone(j)%type = 7
577 int(nzones_offset,
i8) * int(nmsh_zone_size,
i8)
578 call mpi_file_write_at_all(fh, mpi_offset, &
581 deallocate(nmsh_zone)
584 ncurves = msh%curve%size
587 call mpi_allreduce(ncurves, ncurves_glb, 1, &
591 call mpi_exscan(ncurves, ncurves_offset, 1, &
595 int(nzones_glb,
i8)*int(nmsh_zone_size,
i8)
597 call mpi_file_write_at_all(fh, mpi_offset, &
598 ncurves_glb, 1, mpi_integer, status, ierr)
600 if (ncurves_glb .gt. 0)
then
602 allocate(nmsh_curve(ncurves))
605 nmsh_curve(i)%type = 0
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
615 int(nzones_glb,
i8) * int(nmsh_zone_size,
i8) + &
616 int(ncurves_offset,
i8) * int(nmsh_curve_size,
i8)
618 call mpi_file_write_at_all(fh, mpi_offset, &
620 deallocate(nmsh_curve)
623 call mpi_file_sync(fh, ierr)
624 call mpi_file_close(fh, ierr)
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)