76 class(*),
target,
intent(inout) :: data
77 type(
mesh_t),
pointer :: msh
79 character(len=5) :: hdr_ver
80 character(len=54) :: hdr_str
81 character(len=80) :: hdr_full
82 integer :: nel, ndim, nelv, ierr, nBCre2
83 type(mpi_status) :: status
85 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
88 integer :: ncurv, nbcs
92 character(len=1024) :: map_fname
94 integer :: re2_data_xy_size
95 integer :: re2_data_xyz_size
96 integer :: re2_data_cv_size
97 integer :: re2_data_bc_size
99 character(len=LOG_SIZE) :: log_buf
101 call this%check_exists()
111 open(newunit=file_unit,
file=trim(this%fname), status=
'old', iostat=ierr)
112 call neko_log%message(
'Reading binary NEKTON file ' // this%fname)
114 read(file_unit,
'(a80)') hdr_full
115 read(hdr_full,
'(a5)') hdr_ver
117 if (hdr_ver .eq.
'#v004')
then
118 read(hdr_full,
'(a5,i16,i3,i16,i4,a36)') hdr_ver, nel, ndim, nelv, nbcre2, hdr_str
120 else if (hdr_ver .eq.
'#v002' .or. hdr_ver .eq.
'#v003')
then
121 read(hdr_full,
'(a5,i9,i3,i9,a54)') hdr_ver, nel, ndim, nelv, hdr_str
123 else if (hdr_ver .eq.
'#v001')
then
124 read(hdr_full,
'(a5,i9,i3,i9,a54)') hdr_ver, nel, ndim, nelv, hdr_str
139 write(log_buf,1) ndim, nelv
1401
format(
'gdim = ', i1,
', nelements =', i7)
146 inquire(
file=map_fname, exist=read_map)
148 call nm%init(nelv, 2**ndim)
152 call neko_log%warning(
'No NEKTON map file found')
155 call mpi_file_open(
neko_comm, trim(this%fname), &
156 mpi_mode_rdonly, mpi_info_null, fh, ierr)
158 if (ierr .ne. 0)
then
159 call neko_log%error(
"Can't open binary NEKTON file ")
164 call msh%init(ndim, dist)
169 call mpi_file_read_at_all(fh, mpi_offset, test, 1, mpi_real, status, ierr)
173 call neko_error(
'Invalid endian of re2 file, byte swap not implemented yet')
177 mpi_offset, re2_data_xy_size, re2_data_xyz_size, v2_format)
182 if (ndim .eq. 2)
then
183 mpi_offset = mpi_offset + int(dist%num_global(),
i8) * int(re2_data_xy_size,
i8)
185 mpi_offset = mpi_offset + int(dist%num_global(),
i8) * int(re2_data_xyz_size,
i8)
191 call mpi_file_read_at_all(fh, mpi_offset, t2, 1, mpi_double_precision, status, ierr)
195 mpi_offset = mpi_offset + int(ncurv,
i8) * int(re2_data_cv_size,
i8)
196 call mpi_file_read_at_all(fh, mpi_offset, t2, 1, mpi_double_precision, status, ierr)
202 call mpi_file_read_at_all(fh, mpi_offset, ncurv, 1, mpi_integer, status, ierr)
205 mpi_offset = mpi_offset + int(ncurv,
i8) * int(re2_data_cv_size,
i8)
206 call mpi_file_read_at_all(fh, mpi_offset, nbcs, 1, mpi_integer, status, ierr)
213 call mpi_file_close(fh, ierr)
224 class(*),
target,
intent(in) :: data
225 real(kind=
rp),
intent(in),
optional :: t
226 type(
re2v1_xy_t),
allocatable :: re2_data_xy(:)
228 type(
mesh_t),
pointer :: msh
229 character(len=5),
parameter :: RE2_HDR_VER =
'#v001'
230 character(len=54),
parameter :: RE2_HDR_STR =
'RE2 exported by NEKO'
231 integer :: i, j, ierr, nelgv
232 type(mpi_status) :: status
234 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
235 integer :: element_offset
236 integer :: re2_data_xy_size
237 integer :: re2_data_xyz_size
248 call mpi_reduce(msh%nelv, nelgv, 1, &
249 mpi_integer, mpi_sum, 0,
neko_comm, ierr)
251 call mpi_exscan(msh%nelv, element_offset, 1, &
254 call neko_log%message(
'Writing data as a binary NEKTON file ' // this%fname)
257 open(unit=9,
file=trim(this%fname), status=
'new', iostat=ierr)
258 write(9,
'(a5,i9,i3,i9,a54)') re2_hdr_ver, nelgv, msh%gdim,&
264 call mpi_file_open(
neko_comm, trim(this%fname), &
265 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
269 mpi_real, status, ierr)
272 if (msh%gdim .eq. 2)
then
273 allocate(re2_data_xy(msh%nelv))
275 re2_data_xy(i)%rgroup = 1.0
277 re2_data_xy(i)%x(j) =
real(msh%elements(i)%e%pts(j)%p%x(1))
278 re2_data_xy(i)%y(j) =
real(msh%elements(i)%e%pts(j)%p%x(2))
281 mpi_offset = mpi_offset + element_offset * re2_data_xy_size
282 call mpi_file_write_at(fh, mpi_offset, &
285 deallocate(re2_data_xy)
287 else if (msh%gdim .eq. 3)
then
288 allocate(re2_data_xyz(msh%nelv))
290 re2_data_xyz(i)%rgroup = 1.0
292 re2_data_xyz(i)%x(j) =
real(msh%elements(i)%e%pts(j)%p%x(1))
293 re2_data_xyz(i)%y(j) =
real(msh%elements(i)%e%pts(j)%p%x(2))
294 re2_data_xyz(i)%z(j) =
real(msh%elements(i)%e%pts(j)%p%x(3))
297 mpi_offset = mpi_offset + element_offset * re2_data_xyz_size
298 call mpi_file_write_at(fh, mpi_offset, &
301 deallocate(re2_data_xyz)
306 call mpi_file_close(fh, ierr)
314 mpi_offset, re2_data_xy_size, re2_data_xyz_size, v2_format)
315 type(
mesh_t),
intent(inout) :: msh
316 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
317 integer,
intent(inout) :: ndim
318 integer,
intent(inout) :: nel
319 type(mpi_file),
intent(inout) :: fh
320 integer,
intent(in) :: re2_data_xy_size
321 integer,
intent(in) :: re2_data_xyz_size
322 logical,
intent(in) :: v2_format
324 integer :: element_offset
325 type(
re2v1_xy_t),
allocatable :: re2v1_data_xy(:)
326 type(
re2v1_xyz_t),
allocatable :: re2v1_data_xyz(:)
327 type(
re2v2_xy_t),
allocatable :: re2v2_data_xy(:)
328 type(
re2v2_xyz_t),
allocatable :: re2v2_data_xyz(:)
329 type(mpi_status) :: status
332 integer :: pt_idx, nelv
333 integer :: i, j, ierr
335 nelv = dist%num_local()
336 element_offset = dist%start_idx()
338 call htp%init(2*nel, ndim)
340 if (ndim .eq. 2)
then
341 mpi_offset = mpi_offset + element_offset * re2_data_xy_size
342 if (.not. v2_format)
then
343 allocate(re2v1_data_xy(nelv))
344 call mpi_file_read_at_all(fh, mpi_offset, &
349 real(re2v1_data_xy(i)%y(j),dp), 0.0d0)
353 if(mod(i,nelv/10) .eq. 0)
write(*,*) i,
'elements read'
356 call msh%add_element(i, i, p(1), p(2), p(4), p(3))
358 deallocate(re2v1_data_xy)
360 allocate(re2v2_data_xy(nelv))
361 call mpi_file_read_at_all(fh, mpi_offset, &
365 p(j) =
point_t(re2v2_data_xy(i)%x(j), &
366 re2v2_data_xy(i)%y(j), 0.0d0)
370 if(mod(i,nelv/10) .eq. 0)
write(*,*) i,
'elements read'
373 call msh%add_element(i, i, p(1), p(2), p(4), p(3))
375 deallocate(re2v2_data_xy)
377 else if (ndim .eq. 3)
then
378 mpi_offset = mpi_offset + element_offset * re2_data_xyz_size
379 if (.not. v2_format)
then
380 allocate(re2v1_data_xyz(nelv))
381 call mpi_file_read_at_all(fh, mpi_offset, &
386 real(re2v1_data_xyz(i)%y(j),dp),&
387 real(re2v1_data_xyz(i)%z(j),dp))
391 if(mod(i,nelv/100) .eq. 0)
write(*,*) i,
'elements read'
394 call msh%add_element(i, i, &
395 p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
397 deallocate(re2v1_data_xyz)
399 allocate(re2v2_data_xyz(nelv))
400 call mpi_file_read_at_all(fh, mpi_offset, &
404 p(j) =
point_t(re2v2_data_xyz(i)%x(j), &
405 re2v2_data_xyz(i)%y(j),&
406 re2v2_data_xyz(i)%z(j))
410 if(mod(i,nelv/100) .eq. 0)
write(*,*) i,
'elements read'
413 call msh%add_element(i, i, &
414 p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
416 deallocate(re2v2_data_xyz)
424 type(
mesh_t),
intent(inout) :: msh
425 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
426 integer,
intent(inout) :: ncurve
428 type(mpi_file),
intent(inout) :: fh
429 logical,
intent(in) :: v2_format
430 type(mpi_status) :: status
431 integer :: i, j, l, ierr, el_idx, id
434 real(kind=
dp),
allocatable :: curve_data(:,:,:)
435 integer,
allocatable :: curve_type(:,:)
436 logical,
allocatable :: curve_element(:)
437 character(len=1) :: chtemp
438 logical :: curve_skip = .false.
440 allocate(curve_data(5,12,msh%nelv))
441 allocate(curve_element(msh%nelv))
442 allocate(curve_type(12,msh%nelv))
444 curve_element(i) = .false.
448 curve_data(l,j,i) = 0d0
453 call neko_log%message(
'Reading curved data')
454 if (.not. v2_format)
then
455 allocate(re2v1_data_curve(ncurve))
456 call mpi_file_read_at_all(fh, mpi_offset, re2v1_data_curve, ncurve, &
459 allocate(re2v2_data_curve(ncurve))
460 call mpi_file_read_at_all(fh, mpi_offset, re2v2_data_curve, ncurve, &
466 el_idx = re2v2_data_curve(i)%elem - dist%start_idx()
467 id = re2v2_data_curve(i)%zone
468 chtemp = re2v2_data_curve(i)%type
470 curve_data(j,id, el_idx) = re2v2_data_curve(i)%point(j)
473 el_idx = re2v1_data_curve(i)%elem - dist%start_idx()
474 id = re2v1_data_curve(i)%zone
475 chtemp = re2v1_data_curve(i)%type
477 curve_data(j,id, el_idx) =
real(re2v1_data_curve(i)%point(j),
dp)
481 curve_element(el_idx) = .true.
483 select case(trim(chtemp))
485 curve_type(id,el_idx) = 1
486 call neko_log%warning(
'curve type s not supported, treating mesh as non-curved')
490 curve_type(id,el_idx) = 2
491 call neko_log%warning(
'curve type e not supported, treating mesh as non-curved')
495 curve_type(id,el_idx) = 3
497 curve_type(id,el_idx) = 4
499 write(*,*) chtemp,
'curve type not supported yet, treating mesh as non-curved',id, el_idx
506 deallocate(re2v2_data_curve)
508 deallocate(re2v1_data_curve)
510 if (.not. curve_skip)
then
511 do el_idx = 1, msh%nelv
512 if (curve_element(el_idx))
then
513 call msh%mark_curve_element(el_idx, &
514 curve_data(1,1,el_idx), curve_type(1,el_idx))
519 deallocate(curve_data)
520 deallocate(curve_element)
521 deallocate(curve_type)
526 type(
mesh_t),
intent(inout) :: msh
527 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
528 integer,
intent(inout) :: nbcs
530 type(mpi_file),
intent(inout) :: fh
531 logical,
intent(in) :: v2_format
532 type(mpi_status) :: status
534 integer :: sym_facet, label
535 integer :: p_el_idx, p_facet
536 integer :: i, j, ierr, el_idx
537 integer,
parameter,
dimension(6) :: facet_map = (/3, 2, 4, 1, 5, 6/)
539 type(
re2v1_bc_t),
allocatable :: re2v1_data_bc(:)
540 type(
re2v2_bc_t),
allocatable :: re2v2_data_bc(:)
541 character(len=LOG_SIZE) :: log_buf
544 integer :: user_labeled_zones(NEKO_MSH_MAX_ZLBLS)
545 integer :: labeled_zone_offsets(NEKO_MSH_MAX_ZLBLS)
546 integer :: total_labeled_zone_offset, current_internal_zone
547 total_labeled_zone_offset = 0
548 labeled_zone_offsets = 0
549 user_labeled_zones = 0
550 current_internal_zone = 1
553 call neko_log%message(
"Reading boundary conditions")
554 if (.not. v2_format)
then
555 allocate(re2v1_data_bc(nbcs))
556 call mpi_file_read_at_all(fh, mpi_offset, re2v1_data_bc, nbcs, &
559 allocate(re2v2_data_bc(nbcs))
560 call mpi_file_read_at_all(fh, mpi_offset, re2v2_data_bc, nbcs, &
574 el_idx = int(re2v2_data_bc(i)%elem) - dist%start_idx()
575 sym_facet = facet_map(int(re2v2_data_bc(i)%face))
577 select case(trim(re2v2_data_bc(i)%type))
578 case (
'MSH',
'msh',
'EXO',
'exo')
580 label = int(re2v2_data_bc(i)%bc_data(5))
582 if (label .lt. 1 .or. label .gt. neko_msh_max_zlbls)
then
583 call neko_error(
'Invalid label id (valid range [1,...,20])')
586 if (user_labeled_zones(label) .eq. 0)
then
587 write (log_buf,
"(A,I2,A)")
"Labeled zone ", label, &
590 user_labeled_zones(label) = 1
593 call msh%mark_labeled_facet(sym_facet, el_idx, label)
597 total_labeled_zone_offset =
max(total_labeled_zone_offset, label)
604 el_idx = int(re2v2_data_bc(i)%elem) - dist%start_idx()
605 sym_facet = facet_map(int(re2v2_data_bc(i)%face))
606 select case(trim(re2v2_data_bc(i)%type))
607 case (
'MSH',
'msh',
'EXO',
'exo')
612 current_internal_zone = current_internal_zone + 1
616 re2v2_data_bc(i)%type, &
625 current_internal_zone = current_internal_zone + 1
629 re2v2_data_bc(i)%type, &
637 current_internal_zone = current_internal_zone + 1
641 re2v2_data_bc(i)%type, &
649 current_internal_zone = current_internal_zone + 1
653 re2v2_data_bc(i)%type, &
661 current_internal_zone = current_internal_zone + 1
665 re2v2_data_bc(i)%type, &
670 case (
's',
'sl',
'sh',
'shl',
'S',
'SL',
'SH',
'SHL')
673 current_internal_zone = current_internal_zone + 1
677 re2v2_data_bc(i)%type, &
684 p_el_idx = int(re2v2_data_bc(i)%bc_data(1))
685 p_facet = facet_map(int(re2v2_data_bc(i)%bc_data(2)))
686 call msh%get_facet_ids(sym_facet, el_idx, pids)
687 call msh%mark_periodic_facet(sym_facet, el_idx, &
688 p_facet, p_el_idx, pids)
690 write (*,*) trim(re2v2_data_bc(i)%type),
'bc type not supported yet'
691 write (*,*) re2v2_data_bc(i)%bc_data
701 el_idx = re2v2_data_bc(i)%elem - dist%start_idx()
702 sym_facet = facet_map(int(re2v2_data_bc(i)%face))
703 select case(trim(re2v2_data_bc(i)%type))
705 p_el_idx = int(re2v2_data_bc(i)%bc_data(1))
706 p_facet = facet_map(int(re2v2_data_bc(i)%bc_data(2)))
707 call msh%create_periodic_ids(sym_facet, el_idx, &
713 deallocate(re2v2_data_bc)
722 el_idx = re2v1_data_bc(i)%elem - dist%start_idx()
723 sym_facet = facet_map(re2v1_data_bc(i)%face)
725 select case(trim(re2v1_data_bc(i)%type))
726 case (
'MSH',
'msh',
'EXO',
'exo')
728 label = int(re2v1_data_bc(i)%bc_data(5))
730 if (label .lt. 1 .or. label .gt. neko_msh_max_zlbls)
then
731 call neko_error(
'Invalid label id (valid range [1,...,20])')
734 if (user_labeled_zones(label) .eq. 0)
then
735 write (log_buf,
"(A,I2,A,I3)")
"Labeled zone ", label, &
738 user_labeled_zones(label) = 1
743 total_labeled_zone_offset =
max(total_labeled_zone_offset, label)
745 call msh%mark_labeled_facet(sym_facet, el_idx, label)
752 el_idx = re2v1_data_bc(i)%elem - dist%start_idx()
753 sym_facet = facet_map(re2v1_data_bc(i)%face)
754 select case(trim(re2v1_data_bc(i)%type))
755 case (
'MSH',
'msh',
'EXO',
'exo')
760 current_internal_zone = current_internal_zone + 1
764 re2v1_data_bc(i)%type, &
773 current_internal_zone = current_internal_zone + 1
777 re2v1_data_bc(i)%type, &
785 current_internal_zone = current_internal_zone + 1
789 re2v1_data_bc(i)%type, &
797 current_internal_zone = current_internal_zone + 1
801 re2v1_data_bc(i)%type, &
809 current_internal_zone = current_internal_zone + 1
813 re2v1_data_bc(i)%type, &
818 case (
's',
'sl',
'sh',
'shl',
'S',
'SL',
'SH',
'SHL')
821 current_internal_zone = current_internal_zone + 1
825 re2v1_data_bc(i)%type, &
832 p_el_idx = int(re2v1_data_bc(i)%bc_data(1))
833 p_facet = facet_map(int(re2v1_data_bc(i)%bc_data(2)))
834 call msh%get_facet_ids(sym_facet, el_idx, pids)
835 call msh%mark_periodic_facet(sym_facet, el_idx, &
836 p_facet, p_el_idx, pids)
838 write (*,*) re2v1_data_bc(i)%type,
'bc type not supported yet'
839 write (*,*) re2v1_data_bc(i)%bc_data
849 el_idx = re2v1_data_bc(i)%elem - dist%start_idx()
850 sym_facet = facet_map(re2v1_data_bc(i)%face)
851 select case(trim(re2v1_data_bc(i)%type))
853 p_el_idx = int(re2v1_data_bc(i)%bc_data(1))
854 p_facet = facet_map(int(re2v1_data_bc(i)%bc_data(2)))
855 call msh%create_periodic_ids(sym_facet, el_idx, &
862 deallocate(re2v1_data_bc)