48 use mpi_f08,
only : mpi_info_null, mpi_allreduce, mpi_allgather, &
49 mpi_in_place, mpi_integer, mpi_sum, mpi_max, mpi_comm_size, mpi_exscan, &
55 h5open_f, h5close_f, &
56 h5fcreate_f, h5fopen_f, h5fclose_f, &
57 h5gcreate_f, h5gopen_f, h5gclose_f, &
58 h5acreate_f, h5aopen_f, h5awrite_f, h5aclose_f, h5aexists_f, &
59 h5dcreate_f, h5dopen_f, h5dwrite_f, h5dclose_f, &
60 h5tcopy_f, h5tclose_f, h5tset_strpad_f, h5tset_size_f, &
61 h5screate_f, h5screate_simple_f, h5sclose_f, &
62 h5sselect_hyperslab_f, h5sselect_all_f, h5sget_simple_extent_dims_f, &
63 h5dget_space_f, h5dset_extent_f, &
64 h5pcreate_f, h5pclose_f, h5pset_fapl_mpio_f, h5pset_dxpl_mpio_f, &
65 h5lexists_f, h5ldelete_f, &
66 h5p_file_access_f, h5p_dataset_xfer_f, h5p_dataset_create_f, &
67 h5f_acc_trunc_f, h5f_acc_rdwr_f, h5pset_chunk_f, &
68 h5t_std_u8le, h5t_native_integer, h5t_fortran_s1, h5t_str_nullterm_f, &
69 h5kind_to_type, h5_real_kind, h5_integer_kind, h5pset_virtual_f, &
70 h5s_scalar_f, h5s_select_set_f, h5fd_mpio_collective_f, h5p_default_f, &
78 logical :: amr_enabled = .false.
79 logical :: subdivide = .false.
80 integer :: precision = -1
101 logical,
intent(in) :: overwrite
102 this%overwrite = overwrite
108 this%amr_enabled = .false.
114 integer,
intent(in) :: precision
115 this%precision = precision
121 character(len=1024) :: base_fname
122 character(len=1024) :: fname
123 character(len=1024) :: path, name, suffix
125 fname = trim(this%get_base_fname())
128 write(base_fname,
'(A,A,"_",I0,A)') &
129 trim(path), trim(name), this%get_start_counter(), trim(suffix)
140 logical,
intent(in) :: subdivide
141 this%subdivide = subdivide
152 class(*),
target,
intent(in) :: data
153 real(kind=
rp),
intent(in),
optional :: t
154 type(
mesh_t),
pointer :: msh
157 integer :: ierr, mpi_info, mpi_comm, i, n_fields
158 integer(hid_t) :: plist_id, file_id, attr_id, vtkhdf_grp
159 integer(hid_t) :: filespace, memspace
160 integer(hsize_t),
dimension(1) :: vdims
161 integer :: lx, ly, lz
162 integer :: local_points, local_cells, local_conn
163 integer :: total_points, total_cells, total_conn
164 integer :: point_offset
165 integer :: max_local_points
166 integer,
allocatable :: part_points(:), part_cells(:), part_conns(:)
167 character(len=1024) :: fname
168 character(len=16) :: type_str
183 n_fields = data%size()
184 allocate(fp(n_fields))
186 fp(i)%ptr => data%items(i)%ptr
189 call neko_error(
'Invalid data type for vtkhdf_file_write')
193 if (.not.
associated(msh))
then
194 call neko_error(
'Mesh must be associated for vtkhdf_file_write')
196 if (dof%Xh%lx < 2 .or. dof%Xh%ly < 2)
then
197 call neko_error(
'VTKHDF linear output requires lx, ly >= 2')
199 if (msh%gdim .eq. 3 .and. dof%Xh%lz < 2)
then
200 call neko_error(
'VTKHDF linear output requires lz >= 2 in 3D')
202 if (msh%gdim .lt. 2 .or. msh%gdim .gt. 3)
then
203 call neko_error(
'VTKHDF output only supports 2D and 3D meshes')
207 if (this%precision .gt.
rp)
then
209 call neko_warning(
'Requested precision is higher than working precision')
210 else if (this%precision .eq. -1)
then
214 call this%increment_counter()
215 fname = trim(this%get_vtkhdf_fname())
216 counter = this%get_counter() - this%get_start_counter()
218 mpi_info = mpi_info_null%mpi_val
222 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
223 call h5pset_fapl_mpio_f(plist_id, mpi_comm, mpi_info, ierr)
225 if (counter .eq. 0)
then
227 call h5fcreate_f(fname, h5f_acc_trunc_f, &
228 file_id, ierr, access_prp = plist_id)
230 call h5fopen_f(fname, h5f_acc_rdwr_f, file_id, ierr, &
231 access_prp = plist_id)
235 call h5lexists_f(file_id,
"VTKHDF", exists, ierr)
237 call h5gopen_f(file_id,
"VTKHDF", vtkhdf_grp, ierr)
239 call h5gcreate_f(file_id,
"VTKHDF", vtkhdf_grp, ierr)
243 call h5screate_simple_f(1, vdims, filespace, ierr)
244 call h5acreate_f(vtkhdf_grp,
"Version", h5t_native_integer, filespace, &
246 call h5awrite_f(attr_id, h5t_native_integer,
vtkhdf_version, vdims, ierr)
247 call h5aclose_f(attr_id, ierr)
248 call h5sclose_f(filespace, ierr)
251 type_str =
"UnstructuredGrid"
253 call h5screate_f(h5s_scalar_f, filespace, ierr)
255 call h5tcopy_f(h5t_fortran_s1, memspace, ierr)
256 call h5tset_size_f(memspace, int(len_trim(type_str), kind=hsize_t), ierr)
257 call h5tset_strpad_f(memspace, h5t_str_nullterm_f, ierr)
259 call h5acreate_f(vtkhdf_grp,
"Type", memspace, filespace, attr_id, ierr)
260 call h5awrite_f(attr_id, memspace, [type_str], vdims, ierr)
261 call h5aclose_f(attr_id, ierr)
263 call h5tclose_f(memspace, ierr)
264 call h5sclose_f(filespace, ierr)
268 call vtkhdf_write_steps(vtkhdf_grp, counter, t)
271 if (
associated(msh))
then
272 call vtkhdf_write_mesh(vtkhdf_grp, dof, msh, &
273 this%amr_enabled, counter, this%subdivide, t)
277 if (n_fields > 0)
then
278 call vtkhdf_write_pointdata(vtkhdf_grp, fp, this%precision, counter, &
282 call h5gclose_f(vtkhdf_grp, ierr)
283 call h5pclose_f(plist_id, ierr)
284 call h5fclose_f(file_id, ierr)
287 if (
allocated(fp))
deallocate(fp)
300 subroutine vtkhdf_write_mesh(vtkhdf_grp, dof, msh, amr, &
301 counter, subdivide, t)
303 type(
mesh_t),
intent(in) :: msh
304 integer(hid_t),
intent(in) :: vtkhdf_grp
305 logical,
intent(in) :: amr
306 integer,
intent(in) :: counter
307 logical,
intent(in) :: subdivide
308 real(kind=
rp),
intent(in),
optional :: t
310 integer(kind=1) :: VTK_cell_type
311 integer :: ierr, i, ii, jj, kk, el, local_idx
312 integer :: lx, ly, lz, npts_per_cell, nodes_per_cell, cells_per_element
313 integer :: local_points, local_cells, local_conn
314 integer :: total_points, total_cells, total_conn
315 integer :: point_offset, max_local_points
316 integer :: total_offsets, cell_offset, conn_offset, offsets_offset
317 integer :: max_local_cells, max_local_conn
318 integer(hid_t) :: xf_id, dset_id, dcpl_id, grp_id, attr_id
319 integer(hid_t) :: filespace, memspace, H5T_NEKO_DOUBLE
320 integer(hsize_t),
dimension(1) :: dcount, vdims, maxdims, doffset, chunkdims
321 integer(hsize_t),
dimension(2) :: dcount2, vdims2, maxdims2, doffset2
322 integer(kind=i8) :: i8_value
324 integer,
dimension(3) :: component_sizes
325 integer,
dimension(3) :: component_offsets
326 integer,
dimension(3) :: component_max_sizes
332 if (subdivide .and. msh%gdim .eq. 3)
then
334 cells_per_element = (lx - 1) * (ly - 1) * (lz - 1)
336 else if (subdivide .and. msh%gdim .eq. 2)
then
338 cells_per_element = (lx - 1) * (ly - 1)
340 else if (msh%gdim .eq. 3)
then
342 cells_per_element = 1
343 nodes_per_cell = lx * ly * lz
344 else if (msh%gdim .eq. 2)
then
346 cells_per_element = 1
347 nodes_per_cell = lx * ly
351 local_points = dof%size()
352 local_cells = msh%nelv * cells_per_element
353 local_conn = local_cells * nodes_per_cell
355 total_points = dof%global_size()
356 total_cells = msh%glb_nelv * cells_per_element
357 total_conn = total_cells * nodes_per_cell
359 component_sizes = [local_points, local_cells, local_conn]
360 component_offsets = 0
361 component_max_sizes = 0
363 call mpi_exscan(component_sizes, component_offsets, 3, mpi_integer, &
365 call mpi_allreduce(component_sizes, component_max_sizes, 3, mpi_integer, &
368 point_offset = component_offsets(1)
369 cell_offset = component_offsets(2)
370 conn_offset = component_offsets(3)
371 max_local_points = component_max_sizes(1)
372 max_local_cells = component_max_sizes(2)
373 max_local_conn = component_max_sizes(3)
375 offsets_offset = cell_offset +
pe_rank
376 total_offsets = total_cells +
pe_size
379 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
380 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
388 integer(hsize_t),
dimension(1) :: nof_dims, nof_maxdims
389 integer(hsize_t),
dimension(1) :: nof_count, nof_offset, nof_chunk
390 integer(hid_t) :: nof_filespace, nof_memspace, nof_dcpl
391 integer(kind=i8) :: nof_values(3)
393 nof_count(1) = 1_hsize_t
394 nof_offset(1) = int(counter, hsize_t) * int(
pe_size, hsize_t) &
396 nof_chunk(1) =
max(1_hsize_t, int(
pe_size, hsize_t))
397 nof_values = [int(local_points, kind=
i8), int(local_cells, kind=
i8), &
398 int(local_conn, kind=
i8)]
400 call h5pcreate_f(h5p_dataset_create_f, nof_dcpl, ierr)
401 call h5pset_chunk_f(nof_dcpl, 1, nof_chunk, ierr)
403 call vtkhdf_write_numberof(vtkhdf_grp,
"NumberOfPoints", &
404 nof_values(1), nof_offset, nof_count, nof_dcpl, &
405 counter, xf_id, ierr)
406 call vtkhdf_write_numberof(vtkhdf_grp,
"NumberOfCells", &
407 nof_values(2), nof_offset, nof_count, nof_dcpl, &
408 counter, xf_id, ierr)
409 call vtkhdf_write_numberof(vtkhdf_grp,
"NumberOfConnectivityIds", &
410 nof_values(3), nof_offset, nof_count, nof_dcpl, &
411 counter, xf_id, ierr)
413 call h5pclose_f(nof_dcpl, ierr)
417 call h5lexists_f(vtkhdf_grp,
"Points", exists, ierr)
418 if (.not. exists)
then
420 vdims2 = [3_hsize_t, int(total_points, hsize_t)]
421 maxdims2 = [3_hsize_t, h5s_unlimited_f]
422 chunkdims(1) = int(
max(1, min(max_local_points, total_points)), hsize_t)
423 dcount2 = [3_hsize_t, int(local_points, hsize_t)]
424 doffset2 = [0_hsize_t, int(point_offset, hsize_t)]
425 h5t_neko_double = h5kind_to_type(
dp, h5_real_kind)
427 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
428 call h5screate_simple_f(2, dcount2, memspace, ierr)
429 call h5screate_simple_f(2, vdims2, filespace, ierr, maxdims2)
431 call h5pset_chunk_f(dcpl_id, 2, [3_hsize_t, chunkdims(1)], ierr)
432 call h5dcreate_f(vtkhdf_grp,
"Points", h5t_neko_double, &
433 filespace, dset_id, ierr, dcpl_id = dcpl_id)
434 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
435 doffset2, dcount2, ierr)
438 real(kind=
dp),
allocatable :: coords(:,:)
440 allocate(coords(3, local_points))
441 do concurrent(local_idx = 1:local_points)
446 coords(1, local_idx) = dof%x(idx(1), idx(2), idx(3), idx(4))
447 coords(2, local_idx) = dof%y(idx(1), idx(2), idx(3), idx(4))
448 coords(3, local_idx) = dof%z(idx(1), idx(2), idx(3), idx(4))
451 call h5dwrite_f(dset_id, h5t_neko_double, coords, dcount2, ierr, &
452 file_space_id = filespace, mem_space_id = memspace, &
457 call h5dclose_f(dset_id, ierr)
458 call h5sclose_f(filespace, ierr)
459 call h5sclose_f(memspace, ierr)
460 call h5pclose_f(dcpl_id, ierr)
464 call h5lexists_f(vtkhdf_grp,
"Connectivity", exists, ierr)
465 if (exists)
call h5ldelete_f(vtkhdf_grp,
"Connectivity", ierr)
467 vdims = int(total_conn, hsize_t)
468 maxdims = h5s_unlimited_f
469 chunkdims = int(
max(1, min(max_local_conn, total_conn)), hsize_t)
470 dcount = int(local_conn, hsize_t)
471 doffset = int(conn_offset, hsize_t)
473 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
474 call h5screate_simple_f(1, dcount, memspace, ierr)
475 call h5screate_simple_f(1, vdims, filespace, ierr, maxdims)
477 call h5pset_chunk_f(dcpl_id, 1, chunkdims, ierr)
478 call h5dcreate_f(vtkhdf_grp,
"Connectivity", h5t_native_integer, &
479 filespace, dset_id, ierr, dcpl_id = dcpl_id)
480 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
481 doffset, dcount, ierr)
484 integer,
allocatable :: connectivity(:)
486 allocate(connectivity(local_conn))
487 call vtkhdf_build_connectivity(connectivity, vtk_cell_type, msh, dof, &
489 call h5dwrite_f(dset_id, h5t_native_integer, connectivity, dcount, &
490 ierr, file_space_id = filespace, mem_space_id = memspace, &
492 deallocate(connectivity)
495 call h5dclose_f(dset_id, ierr)
496 call h5sclose_f(filespace, ierr)
497 call h5sclose_f(memspace, ierr)
498 call h5pclose_f(dcpl_id, ierr)
501 call h5lexists_f(vtkhdf_grp,
"Offsets", exists, ierr)
502 if (exists)
call h5ldelete_f(vtkhdf_grp,
"Offsets", ierr)
504 vdims = int(total_offsets, hsize_t)
505 maxdims = h5s_unlimited_f
506 chunkdims = int(
max(1, min(max_local_cells + 1, total_offsets)), hsize_t)
507 dcount = int(local_cells + 1, hsize_t)
508 doffset = int(offsets_offset, hsize_t)
510 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
511 call h5screate_simple_f(1, dcount, memspace, ierr)
512 call h5screate_simple_f(1, vdims, filespace, ierr, maxdims)
514 call h5pset_chunk_f(dcpl_id, 1, chunkdims, ierr)
515 call h5dcreate_f(vtkhdf_grp,
"Offsets", h5t_native_integer, &
516 filespace, dset_id, ierr, dcpl_id = dcpl_id)
517 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
518 doffset, dcount, ierr)
521 integer,
allocatable :: offsets(:)
523 allocate(offsets(local_cells + 1))
524 do concurrent(i = 1:local_cells)
525 offsets(i) = (i - 1) * nodes_per_cell
527 offsets(local_cells + 1) = local_conn
528 call h5dwrite_f(dset_id, h5t_native_integer, offsets, dcount, ierr, &
529 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
533 call h5dclose_f(dset_id, ierr)
534 call h5sclose_f(filespace, ierr)
535 call h5sclose_f(memspace, ierr)
536 call h5pclose_f(dcpl_id, ierr)
539 call h5lexists_f(vtkhdf_grp,
"Types", exists, ierr)
540 if (exists)
call h5ldelete_f(vtkhdf_grp,
"Types", ierr)
542 vdims = int(total_cells, hsize_t)
543 maxdims = h5s_unlimited_f
544 chunkdims = int(
max(1, min(max_local_cells, total_cells)), hsize_t)
545 dcount = int(local_cells, hsize_t)
546 doffset = int(cell_offset, hsize_t)
548 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
549 call h5screate_simple_f(1, dcount, memspace, ierr)
550 call h5screate_simple_f(1, vdims, filespace, ierr, maxdims)
552 call h5pset_chunk_f(dcpl_id, 1, chunkdims, ierr)
553 call h5dcreate_f(vtkhdf_grp,
"Types", h5t_std_u8le, &
554 filespace, dset_id, ierr, dcpl_id = dcpl_id)
555 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
556 doffset, dcount, ierr)
559 integer(kind=1),
allocatable :: cell_types(:)
560 allocate(cell_types(local_cells), source=vtk_cell_type)
561 call h5dwrite_f(dset_id, h5t_std_u8le, cell_types, dcount, ierr, &
562 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
563 deallocate(cell_types)
566 call h5dclose_f(dset_id, ierr)
567 call h5sclose_f(filespace, ierr)
568 call h5sclose_f(memspace, ierr)
569 call h5pclose_f(dcpl_id, ierr)
573 call h5gopen_f(vtkhdf_grp,
"Steps", grp_id, ierr)
576 call vtkhdf_write_i8_at(grp_id,
"NumberOfParts", int(pe_size, kind=i8), &
582 i8_value = int(counter - 1, kind=i8) * int(pe_size, kind=i8)
583 call vtkhdf_write_i8_at(grp_id,
"PartOffsets", i8_value, counter)
585 i8_value = int(counter - 1, kind=i8) * int(total_points, kind=i8)
586 call vtkhdf_write_i8_at(grp_id,
"PointOffsets", i8_value, counter)
588 i8_value = int(counter - 1, kind=i8) * int(total_cells, kind=i8)
589 call vtkhdf_write_i8_at(grp_id,
"CellOffsets", i8_value, counter)
591 i8_value = int(counter - 1, kind=i8) * int(total_conn, kind=i8)
592 call vtkhdf_write_i8_at(grp_id,
"ConnectivityIdOffsets", i8_value, &
597 call vtkhdf_write_i8_at(grp_id,
"PartOffsets", i8_value, counter)
598 call vtkhdf_write_i8_at(grp_id,
"PointOffsets", i8_value, counter)
599 call vtkhdf_write_i8_at(grp_id,
"CellOffsets", i8_value, counter)
600 call vtkhdf_write_i8_at(grp_id,
"ConnectivityIdOffsets", i8_value, &
604 call h5gclose_f(grp_id, ierr)
607 call h5pclose_f(xf_id, ierr)
609 end subroutine vtkhdf_write_mesh
617 subroutine vtkhdf_write_steps(vtkhdf_grp, counter, t)
618 integer(hid_t),
intent(in) :: vtkhdf_grp
619 integer,
intent(in) :: counter
620 real(kind=rp),
intent(in) :: t
622 integer(hid_t) :: xf_id, H5T_NEKO_DOUBLE
624 integer(hid_t) :: grp_id, dset_id, dcpl_id, filespace, memspace, attr_id
625 integer(hsize_t),
dimension(1) :: step_dims, step_maxdims
626 integer(hsize_t),
dimension(1) :: step_count, step_offset, chunkdims, ddim
627 real(kind=dp),
dimension(1) :: time_value
628 integer(kind=i8) :: i8_value
629 logical :: exists, attr_exists
632 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
633 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
634 h5t_neko_double = h5kind_to_type(dp, h5_real_kind)
637 call h5lexists_f(vtkhdf_grp,
"Steps", exists, ierr)
639 call h5gopen_f(vtkhdf_grp,
"Steps", grp_id, ierr)
641 call h5gcreate_f(vtkhdf_grp,
"Steps", grp_id, ierr)
645 call h5lexists_f(grp_id,
"Values", exists, ierr)
647 call h5dopen_f(grp_id,
"Values", dset_id, ierr)
648 call h5dget_space_f(dset_id, filespace, ierr)
649 call h5sget_simple_extent_dims_f(filespace, step_dims, step_maxdims, &
651 call h5sclose_f(filespace, ierr)
654 if (step_dims(1) .eq. counter)
then
655 step_dims(1) = int(counter + 1, hsize_t)
656 call h5dset_extent_f(dset_id, step_dims, ierr)
657 else if (step_dims(1) .lt. counter)
then
658 call neko_error(
"VTKHDF: Time steps written out of order.")
661 step_dims(1) = 1_hsize_t
662 step_maxdims(1) = h5s_unlimited_f
663 chunkdims(1) = 1_hsize_t
665 call h5screate_simple_f(1, step_dims, filespace, ierr, step_maxdims)
666 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
667 call h5pset_chunk_f(dcpl_id, 1, chunkdims, ierr)
668 call h5dcreate_f(grp_id,
"Values", h5t_neko_double, &
669 filespace, dset_id, ierr, dcpl_id = dcpl_id)
670 call h5sclose_f(filespace, ierr)
671 call h5pclose_f(dcpl_id, ierr)
674 step_count(1) = 1_hsize_t
675 step_offset(1) = int(counter, hsize_t)
677 call h5dget_space_f(dset_id, filespace, ierr)
678 call h5screate_simple_f(1, step_count, memspace, ierr)
679 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
680 step_offset, step_count, ierr)
682 time_value(1) =
real(t, kind=dp)
683 call h5dwrite_f(dset_id, h5t_neko_double, time_value, step_count, ierr, &
684 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
686 call h5sclose_f(memspace, ierr)
687 call h5sclose_f(filespace, ierr)
688 call h5dclose_f(dset_id, ierr)
692 call h5aexists_f(grp_id,
"NSteps", attr_exists, ierr)
693 if (attr_exists)
then
694 call h5aopen_f(grp_id,
"NSteps", attr_id, ierr)
696 call h5screate_f(h5s_scalar_f, filespace, ierr)
697 call h5acreate_f(grp_id,
"NSteps", h5t_native_integer, filespace, &
698 attr_id, ierr, h5p_default_f, h5p_default_f)
699 call h5sclose_f(filespace, ierr)
702 call h5awrite_f(attr_id, h5t_native_integer, counter + 1, ddim, ierr)
704 call h5aclose_f(attr_id, ierr)
705 call h5gclose_f(grp_id, ierr)
706 call h5pclose_f(xf_id, ierr)
708 end subroutine vtkhdf_write_steps
724 subroutine vtkhdf_write_pointdata(vtkhdf_grp, fp, precision, counter, &
726 integer(hid_t),
intent(in) :: vtkhdf_grp
727 type(field_ptr_t),
intent(in) :: fp(:)
728 integer,
intent(in) :: precision
729 integer,
intent(in) :: counter
730 character(len=*),
intent(in) :: fname
731 real(kind=rp),
intent(in),
optional :: t
733 integer(kind=i8) :: time_offset
734 integer :: local_points, point_offset, total_points
735 integer(hid_t) :: precision_hdf
736 integer :: ierr, i, j
738 integer(hid_t) :: pointdata_grp, grp_id, step_grp_id
739 integer(hid_t) :: dset_id, dcpl_id, filespace
740 integer(hsize_t),
dimension(1) :: pd_dims1, pd_maxdims1
741 integer(hsize_t),
dimension(2) :: pd_dims2, pd_maxdims2
742 type(field_t),
pointer :: fld, u, v, w
743 character(len=128) :: field_name
744 logical :: exists, is_vector
747 character(len=1024) :: ext_fname, ext_path, src_pattern
748 character(len=1024) :: main_path, main_name, main_suffix
749 integer(hid_t) :: ext_file_id, ext_plist_id, vds_src_space
750 integer(hid_t) :: write_target
751 integer :: mpi_info, mpi_comm
754 integer :: fields_written
755 character(len=128),
allocatable :: name_list(:)
756 logical,
allocatable :: vector_list(:)
758 mpi_info = mpi_info_null%mpi_val
759 mpi_comm = neko_comm%mpi_val
764 local_points = fp(1)%ptr%dof%size()
765 total_points = fp(1)%ptr%dof%global_size()
767 call mpi_exscan(local_points, point_offset, 1, mpi_integer, &
768 mpi_sum, neko_comm, ierr)
772 if (
associated(fp(i)%ptr))
then
773 call fp(i)%ptr%copy_from(device_to_host, sync = i .eq. n_fields)
778 allocate(name_list(n_fields))
782 call h5lexists_f(vtkhdf_grp,
"PointData", exists, ierr)
783 if (.not. exists)
then
784 call h5gcreate_f(vtkhdf_grp,
"PointData", pointdata_grp, ierr)
785 call h5gclose_f(pointdata_grp, ierr)
794 call filename_split(fname, main_path, main_name, main_suffix)
795 write(ext_path,
'(A,A,".data/")') trim(main_path), trim(main_name)
796 write(ext_fname,
'(A,I0,".h5")') trim(ext_path), counter
797 write(src_pattern,
'(A,".data/%b.h5")') trim(main_name)
799 if (pe_rank == 0)
then
800 inquire(
file = trim(ext_path), exist = exists)
801 if (.not. exists)
then
802 call execute_command_line(
"mkdir -p '" // trim(ext_path) //
"'")
805 call mpi_barrier(neko_comm, ierr)
807 call h5pcreate_f(h5p_file_access_f, ext_plist_id, ierr)
808 call h5pset_fapl_mpio_f(ext_plist_id, mpi_comm, mpi_info, ierr)
809 call h5fcreate_f(trim(ext_fname), h5f_acc_trunc_f, write_target, ierr, &
810 access_prp = ext_plist_id)
811 call h5pclose_f(ext_plist_id, ierr)
815 call h5gopen_f(vtkhdf_grp,
"PointData", write_target, ierr)
823 field_name = fld%name
824 if (field_name .eq.
'p') field_name =
'Pressure'
828 if (field_name .eq.
'u' .or. &
829 field_name .eq.
'v' .or. &
830 field_name .eq.
'w')
then
835 select case (fp(j)%ptr%name)
845 if (
associated(u) .and.
associated(v) .and.
associated(w))
then
847 field_name =
'Velocity'
854 do j = 1, fields_written
855 if (trim(name_list(j)) .eq. trim(field_name))
then
865 fields_written = fields_written + 1
866 name_list(fields_written) = field_name
871 call write_vector_field(write_target, field_name, u%x, v%x, w%x, &
872 local_points, precision, total_points, point_offset)
874 call write_scalar_field(write_target, field_name, fld%x, &
875 local_points, precision, total_points, point_offset)
881 call h5fclose_f(write_target, ierr)
883 call h5gclose_f(write_target, ierr)
892 call h5gopen_f(vtkhdf_grp,
"Steps", step_grp_id, ierr)
893 time_offset = int(counter, kind=i8) * int(total_points, kind=i8)
896 call h5lexists_f(step_grp_id,
"PointDataOffsets", exists, ierr)
898 call h5gopen_f(step_grp_id,
"PointDataOffsets", grp_id, ierr)
900 call h5gcreate_f(step_grp_id,
"PointDataOffsets", grp_id, ierr)
902 do i = 1, fields_written
903 call vtkhdf_write_i8_at(grp_id, trim(name_list(i)), &
904 time_offset, counter)
906 call h5gclose_f(grp_id, ierr)
907 call h5gclose_f(step_grp_id, ierr)
910 call h5gopen_f(vtkhdf_grp,
"PointData", pointdata_grp, ierr)
912 do i = 1, fields_written
913 field_name = name_list(i)
916 if (counter .eq. 0)
then
918 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
919 precision_hdf = h5kind_to_type(precision, h5_real_kind)
922 pd_dims2 = [3_hsize_t, int(total_points, hsize_t)]
923 call h5screate_simple_f(2, pd_dims2, vds_src_space, ierr)
924 call h5sselect_all_f(vds_src_space, ierr)
926 pd_maxdims2 = [3_hsize_t, h5s_unlimited_f]
927 call h5screate_simple_f(2, pd_dims2, filespace, ierr, &
930 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
931 [0_hsize_t, 0_hsize_t], &
932 [1_hsize_t, h5s_unlimited_f], &
934 stride = [3_hsize_t, int(total_points, hsize_t)], &
935 block = [3_hsize_t, int(total_points, hsize_t)])
937 call h5pset_virtual_f(dcpl_id, filespace, trim(src_pattern), &
938 trim(field_name), vds_src_space, ierr)
939 call h5sclose_f(vds_src_space, ierr)
941 call h5dcreate_f(pointdata_grp, trim(field_name), &
942 precision_hdf, filespace, dset_id, ierr, &
944 call h5sclose_f(filespace, ierr)
946 pd_dims1 = int(total_points, hsize_t)
947 call h5screate_simple_f(1, pd_dims1, vds_src_space, ierr)
948 call h5sselect_all_f(vds_src_space, ierr)
950 pd_maxdims1(1) = h5s_unlimited_f
951 call h5screate_simple_f(1, pd_dims1, filespace, ierr, &
954 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
958 stride = [int(total_points, hsize_t)], &
959 block = [int(total_points, hsize_t)])
961 call h5pset_virtual_f(dcpl_id, filespace, trim(src_pattern), &
962 trim(field_name), vds_src_space, ierr)
963 call h5sclose_f(vds_src_space, ierr)
965 call h5dcreate_f(pointdata_grp, trim(field_name), &
966 precision_hdf, filespace, dset_id, ierr, &
968 call h5sclose_f(filespace, ierr)
971 call h5pclose_f(dcpl_id, ierr)
972 call h5dclose_f(dset_id, ierr)
976 call h5dopen_f(pointdata_grp, trim(field_name), dset_id, ierr)
979 pd_dims2 = [3_hsize_t, &
980 int((counter + 1) * total_points, hsize_t)]
981 call h5dset_extent_f(dset_id, pd_dims2, ierr)
983 pd_dims1 = int((counter + 1) * total_points, hsize_t)
984 call h5dset_extent_f(dset_id, pd_dims1, ierr)
987 call h5dclose_f(dset_id, ierr)
991 call h5gclose_f(pointdata_grp, ierr)
997 deallocate(name_list)
1000 end subroutine vtkhdf_write_pointdata
1017 subroutine vtkhdf_build_connectivity(conn, vtk_type, msh, dof, subdivide)
1018 integer,
intent(inout) :: conn(:)
1019 integer(kind=1),
intent(in) :: vtk_type
1020 type(mesh_t),
intent(in) :: msh
1021 type(dofmap_t),
intent(in) :: dof
1022 logical,
intent(in) :: subdivide
1023 integer :: lx, ly, lz, nelv
1024 integer :: ie, ii, n_pts_per_elem, n_conn_per_elem
1025 integer,
allocatable :: node_order(:)
1031 n_pts_per_elem = lx * ly * lz
1033 if (subdivide .and. vtk_type .eq. 12)
then
1035 else if (subdivide .and. vtk_type .eq. 9)
then
1038 node_order = vtk_ordering(vtk_type, lx, ly, lz)
1041 n_conn_per_elem =
size(node_order)
1043 do concurrent(ie = 1:nelv, ii = 1:n_conn_per_elem)
1045 integer :: idx, base
1046 idx = (ie - 1) * n_conn_per_elem
1047 base = (ie - 1) * n_pts_per_elem
1048 conn(idx + ii) = base + node_order(ii)
1052 deallocate(node_order)
1054 end subroutine vtkhdf_build_connectivity
1068 subroutine vtkhdf_write_numberof(grp, dset_name, value, offset, cnt, &
1069 dcpl, index, xf_id, ierr)
1070 integer(hid_t),
intent(in) :: grp, dcpl, xf_id
1071 character(len=*),
intent(in) :: dset_name
1072 integer(kind=i8),
intent(in) :: value
1073 integer,
intent(in):: index
1074 integer(hsize_t),
dimension(1),
intent(in) :: offset, cnt
1075 integer,
intent(out) :: ierr
1077 integer(hid_t) :: dset_id, fspace, mspace
1078 integer(hsize_t),
dimension(1) :: dims, maxdims
1079 integer(hid_t) :: H5T_NEKO_INTEGER
1082 h5t_neko_integer = h5kind_to_type(i8, h5_integer_kind)
1084 call h5lexists_f(grp, dset_name, exists, ierr)
1086 call h5dopen_f(grp, dset_name, dset_id, ierr)
1087 dims(1) = int(index + 1, hsize_t) * int(pe_size, hsize_t)
1088 call h5dset_extent_f(dset_id, dims, ierr)
1090 dims(1) = int(pe_size, hsize_t)
1091 maxdims(1) = h5s_unlimited_f
1092 call h5screate_simple_f(1, dims, fspace, ierr, maxdims)
1093 call h5dcreate_f(grp, dset_name, h5t_neko_integer, &
1094 fspace, dset_id, ierr, dcpl_id = dcpl)
1095 call h5sclose_f(fspace, ierr)
1098 call h5dget_space_f(dset_id, fspace, ierr)
1099 call h5screate_simple_f(1, cnt, mspace, ierr)
1100 call h5sselect_hyperslab_f(fspace, h5s_select_set_f, offset, cnt, ierr)
1102 call h5dwrite_f(dset_id, h5t_neko_integer,
value, cnt, ierr, &
1103 file_space_id = fspace, mem_space_id = mspace, xfer_prp = xf_id)
1105 call h5sclose_f(mspace, ierr)
1106 call h5sclose_f(fspace, ierr)
1107 call h5dclose_f(dset_id, ierr)
1108 end subroutine vtkhdf_write_numberof
1117 subroutine vtkhdf_write_i8_at(grp_id, name, value, index)
1118 integer(hid_t),
intent(in) :: grp_id
1119 character(len=*),
intent(in) :: name
1120 integer(kind=i8),
intent(in) :: value
1121 integer,
intent(in) :: index
1124 integer(hid_t) :: dset_id, dcpl_id, xf_id, filespace, memspace
1125 integer(hsize_t),
dimension(1) :: dims, maxdims, count, offset, chunkdims
1126 integer(hid_t) :: H5T_NEKO_INTEGER
1129 h5t_neko_integer = h5kind_to_type(i8, h5_integer_kind)
1132 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1133 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1135 call h5lexists_f(grp_id, name, exists, ierr)
1137 call h5dopen_f(grp_id, name, dset_id, ierr)
1138 call h5dget_space_f(dset_id, filespace, ierr)
1139 call h5sget_simple_extent_dims_f(filespace, dims, maxdims, ierr)
1140 call h5sclose_f(filespace, ierr)
1142 if (index .eq. dims(1))
then
1143 dims(1) = int(index + 1, hsize_t)
1144 call h5dset_extent_f(dset_id, dims, ierr)
1145 else if (index .gt. dims(1))
then
1146 call neko_error(
"VTKHDF: Values written out of order.")
1150 maxdims = h5s_unlimited_f
1151 chunkdims = 1_hsize_t
1153 call h5screate_simple_f(1, dims, filespace, ierr, maxdims)
1154 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
1155 call h5pset_chunk_f(dcpl_id, 1, chunkdims, ierr)
1156 call h5dcreate_f(grp_id, name, h5t_neko_integer, &
1157 filespace, dset_id, ierr, dcpl_id = dcpl_id)
1158 call h5sclose_f(filespace, ierr)
1159 call h5pclose_f(dcpl_id, ierr)
1163 offset = int(index, hsize_t)
1165 call h5dget_space_f(dset_id, filespace, ierr)
1166 call h5screate_simple_f(1, count, memspace, ierr)
1167 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, offset, count, ierr)
1168 call h5dwrite_f(dset_id, h5t_neko_integer,
value, count, ierr, &
1169 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
1171 call h5sclose_f(memspace, ierr)
1172 call h5sclose_f(filespace, ierr)
1173 call h5dclose_f(dset_id, ierr)
1174 call h5pclose_f(xf_id, ierr)
1176 end subroutine vtkhdf_write_i8_at
1189 subroutine write_scalar_field(hdf_root, name, x, n_local, &
1190 precision, n_total, offset)
1191 integer,
intent(in) :: n_local
1192 integer(hid_t),
intent(in) :: hdf_root
1193 character(len=*),
intent(in) :: name
1194 real(kind=rp),
dimension(n_local),
intent(in) :: x
1195 integer,
intent(in),
optional :: precision
1196 integer,
intent(in),
optional :: n_total, offset
1198 integer(hsize_t),
dimension(1) :: dims, dcount, doffset
1199 integer(hid_t) :: xf_id, dset_id, filespace, memspace, precision_hdf
1200 integer :: i, ierr, precision_local, n_tot, off
1203 dcount = int(n_local, hsize_t)
1205 if (
present(n_total))
then
1206 dims = int(n_total, hsize_t)
1208 call mpi_allreduce(n_local, n_tot, 1, mpi_integer, mpi_sum, neko_comm, &
1210 dims = int(n_tot, hsize_t)
1213 if (
present(offset))
then
1214 doffset = int(offset, hsize_t)
1216 call mpi_exscan(n_local, off, 1, mpi_integer, mpi_sum, neko_comm, &
1218 doffset = int(off, hsize_t)
1221 if (
present(precision))
then
1222 precision_local = precision
1224 precision_local = rp
1226 precision_hdf = h5kind_to_type(precision_local, h5_real_kind)
1229 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1230 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1232 call h5screate_simple_f(1, dims, filespace, ierr)
1233 call h5screate_simple_f(1, dcount, memspace, ierr)
1234 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
1235 doffset, dcount, ierr)
1237 call h5dcreate_f(hdf_root, trim(name), precision_hdf, &
1238 filespace, dset_id, ierr)
1240 if (precision_local .eq. rp)
then
1241 call h5dwrite_f(dset_id, precision_hdf, x, dcount, ierr, &
1242 file_space_id = filespace, mem_space_id = memspace, &
1245 else if (precision_local .eq. sp)
then
1247 real(kind=sp),
allocatable :: x_sp(:)
1249 allocate(x_sp(n_local))
1250 do concurrent(i = 1:n_local)
1251 x_sp(i) =
real(x(i), sp)
1254 call h5dwrite_f(dset_id, precision_hdf, x_sp, dcount, ierr, &
1255 file_space_id = filespace, mem_space_id = memspace, &
1260 else if (precision_local .eq. dp)
then
1262 real(kind=dp),
allocatable :: x_dp(:)
1264 allocate(x_dp(n_local))
1265 do concurrent(i = 1:n_local)
1266 x_dp(i) =
real(x(i), dp)
1269 call h5dwrite_f(dset_id, precision_hdf, x_dp, dcount, ierr, &
1270 file_space_id = filespace, mem_space_id = memspace, &
1275 else if (precision_local .eq. qp)
then
1277 real(kind=qp),
allocatable :: x_qp(:)
1279 allocate(x_qp(n_local))
1280 do concurrent(i = 1:n_local)
1281 x_qp(i) =
real(x(i), qp)
1284 call h5dwrite_f(dset_id, precision_hdf, x_qp, dcount, ierr, &
1285 file_space_id = filespace, mem_space_id = memspace, &
1291 call neko_error(
"Unsupported precision in HDF5 write_scalar_field")
1294 call h5sclose_f(filespace, ierr)
1295 call h5sclose_f(memspace, ierr)
1296 call h5dclose_f(dset_id, ierr)
1297 call h5pclose_f(xf_id, ierr)
1298 end subroutine write_scalar_field
1313 subroutine write_vector_field(hdf_root, name, u, v, w, n_local, &
1314 precision, n_total, offset)
1315 integer,
intent(in) :: n_local
1316 integer(hid_t),
intent(in) :: hdf_root
1317 character(len=*),
intent(in) :: name
1318 real(kind=rp),
dimension(n_local),
intent(in) :: u, v, w
1319 integer,
intent(in),
optional :: precision
1320 integer,
intent(in),
optional :: n_total, offset
1322 integer(hsize_t),
dimension(2) :: dims, dcount, doffset
1323 integer(hid_t) :: xf_id, dset_id, filespace, memspace, precision_hdf
1324 integer :: i, ierr, precision_local, n_tot, off
1327 dcount = [3_hsize_t, int(n_local, hsize_t)]
1329 if (
present(n_total))
then
1330 dims = [3_hsize_t, int(n_total, hsize_t)]
1332 call mpi_allreduce(n_local, n_tot, 1, mpi_integer, mpi_sum, neko_comm, &
1334 dims = [3_hsize_t, int(n_tot, hsize_t)]
1337 if (
present(offset))
then
1338 doffset = [0_hsize_t, int(offset, hsize_t)]
1340 call mpi_exscan(n_local, off, 1, mpi_integer, mpi_sum, neko_comm, &
1342 doffset = [0_hsize_t, int(off, hsize_t)]
1345 if (
present(precision))
then
1346 precision_local = precision
1348 precision_local = rp
1350 precision_hdf = h5kind_to_type(precision_local, h5_real_kind)
1353 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1354 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1356 call h5screate_simple_f(2, dims, filespace, ierr)
1357 call h5screate_simple_f(2, dcount, memspace, ierr)
1358 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
1359 doffset, dcount, ierr)
1361 call h5dcreate_f(hdf_root, trim(name), precision_hdf, &
1362 filespace, dset_id, ierr)
1364 if (precision_local .eq. sp)
then
1366 real(kind=sp),
allocatable :: f(:,:)
1368 allocate(f(3, n_local))
1369 do concurrent(i = 1:n_local)
1370 f(1, i) =
real(u(i), sp)
1371 f(2, i) =
real(v(i), sp)
1372 f(3, i) =
real(w(i), sp)
1375 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1376 file_space_id = filespace, mem_space_id = memspace, &
1381 else if (precision_local .eq. dp)
then
1383 real(kind=dp),
allocatable :: f(:,:)
1385 allocate(f(3, n_local))
1386 do concurrent(i = 1:n_local)
1387 f(1, i) =
real(u(i), dp)
1388 f(2, i) =
real(v(i), dp)
1389 f(3, i) =
real(w(i), dp)
1392 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1393 file_space_id = filespace, mem_space_id = memspace, &
1398 else if (precision_local .eq. qp)
then
1400 real(kind=qp),
allocatable :: f(:,:)
1402 allocate(f(3, n_local))
1403 do concurrent(i = 1:n_local)
1404 f(1, i) =
real(u(i), qp)
1405 f(2, i) =
real(v(i), qp)
1406 f(3, i) =
real(w(i), qp)
1409 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1410 file_space_id = filespace, mem_space_id = memspace, &
1416 call neko_error(
"Unsupported precision in HDF5 write_vector_field")
1419 call h5sclose_f(filespace, ierr)
1420 call h5sclose_f(memspace, ierr)
1421 call h5dclose_f(dset_id, ierr)
1422 call h5pclose_f(xf_id, ierr)
1423 end subroutine write_vector_field
1428 class(*),
target,
intent(inout) :: data
1430 call neko_error(
'VTKHDF file reading is not yet implemented')
1441 class(*),
target,
intent(in) :: data
1442 real(kind=rp),
intent(in),
optional :: t
1443 call neko_error(
'Neko needs to be built with HDF5 support')
1449 class(*),
target,
intent(inout) :: data
1450 call neko_error(
'Neko needs to be built with HDF5 support')
1467 integer,
intent(in) :: lx, ly, lz
1468 integer :: node_order(8 * (lx - 1) * (ly - 1) * (lz - 1))
1469 integer :: ii, jj, kk, idx
1476 node_order(idx + 1) = (kk - 1) * lx * ly + (jj - 1) * lx + ii - 1
1477 node_order(idx + 2) = (kk - 1) * lx * ly + (jj - 1) * lx + ii
1478 node_order(idx + 3) = (kk - 1) * lx * ly + jj * lx + ii
1479 node_order(idx + 4) = (kk - 1) * lx * ly + jj * lx + ii - 1
1480 node_order(idx + 5) = kk * lx * ly + (jj - 1) * lx + ii - 1
1481 node_order(idx + 6) = kk * lx * ly + (jj - 1) * lx + ii
1482 node_order(idx + 7) = kk * lx * ly + jj * lx + ii
1483 node_order(idx + 8) = kk * lx * ly + jj * lx + ii - 1
1499 integer,
intent(in) :: lx, ly
1500 integer :: node_order(4 * (lx - 1) * (ly - 1))
1501 integer :: ii, jj, idx
1507 node_order(idx + 1) = (jj - 1) * lx + ii - 1
1508 node_order(idx + 2) = (jj - 1) * lx + ii
1509 node_order(idx + 3) = jj * lx + ii
1510 node_order(idx + 4) = jj * lx + ii - 1
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
integer, public pe_size
MPI size of communicator.
integer, public pe_rank
MPI rank.
type(mpi_comm), public neko_comm
MPI communicator.
Device abstraction, common interface for various accelerators.
integer, parameter, public device_to_host
Defines a mapping of the degrees of freedom.
Contains the field_serties_t type.
Module for file I/O operations.
type(log_t), public neko_log
Global log stream.
integer, parameter, public i8
integer, parameter, public qp
integer, parameter, public dp
integer, parameter, public sp
integer, parameter, public rp
Global precision used in computations.
pure integer function, public linear_index(i, j, k, l, lx, ly, lz)
Compute the address of a (i,j,k,l) array with sizes (1:lx, 1:ly, 1:lz, :)
subroutine, public filename_split(fname, path, name, suffix)
Extract file name components.
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
VTK Module containing utilities for VTK file handling.
integer function, dimension(:), allocatable, public vtk_ordering(cell_type, lx, ly, lz)
Get the VTK node ordering for a given cell type. For Lagrange cells, returns an array mapping VTK nod...
subroutine vtkhdf_file_set_overwrite(this, overwrite)
Set the overwrite flag for HDF5 files.
subroutine vtkhdf_file_enable_amr(this)
Enable support for Adaptive Mesh Refinement.
subroutine vtkhdf_file_set_precision(this, precision)
Set the precision for VTKHDF output (single or double)
subroutine vtkhdf_file_write(this, data, t)
Write data in HDF5 format (no HDF5 support)
pure integer function, dimension(8 *(lx - 1) *(ly - 1) *(lz - 1)) subdivide_to_hex_ordering(lx, ly, lz)
Build linear hexahedron sub-cell node ordering for a spectral element. Returns an array of 0-based te...
integer, dimension(2), parameter vtkhdf_version
subroutine vtkhdf_file_set_subdivide(this, subdivide)
Enable or disable subdivision of spectral elements into linear sub-cells. When subdivision is enabled...
subroutine vtkhdf_file_read(this, data)
Read data in HDF5 format (no HDF5 support)
character(len=1024) function vtkhdf_file_get_fname(this)
Return the file name with the start counter.
pure integer function, dimension(4 *(lx - 1) *(ly - 1)) subdivide_to_quad_ordering(lx, ly)
Build linear quadrilateral sub-cell node ordering for a spectral element. Returns an array of 0-based...
field_ptr_t, To easily obtain a pointer to a field
field_list_t, To be able to group fields together
A wrapper for a pointer to a field_series_t.
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
Interface for HDF5 files.