Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
vtkhdf_file.F90
Go to the documentation of this file.
1! Copyright (c) 2026, 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!
35 use num_types, only : rp, sp, dp, qp, i8
37 use checkpoint, only : chkp_t
40 use mesh, only : mesh_t
41 use field, only : field_t, field_ptr_t
42 use field_list, only : field_list_t
44 use dofmap, only : dofmap_t
45 use logger, only : neko_log
46 use comm, only : pe_rank, pe_size, neko_comm
47 use device, only : device_to_host
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, &
50 mpi_barrier, mpi_logical
51 use vtk, only : vtk_ordering
52#ifdef HAVE_HDF5
53 use hdf5, only : &
54 hid_t, hsize_t, &
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, &
71 h5s_unlimited_f
72#endif
73 implicit none
74 private
75
77 type, public, extends(generic_file_t) :: vtkhdf_file_t
78 logical :: amr_enabled = .false.
79 logical :: subdivide = .false.
80 integer :: precision = -1
81 contains
82 procedure :: get_vtkhdf_fname => vtkhdf_file_get_fname
83 procedure :: read => vtkhdf_file_read
84 procedure :: write => vtkhdf_file_write
85 procedure :: set_overwrite => vtkhdf_file_set_overwrite
86 procedure :: enable_amr => vtkhdf_file_enable_amr
87 procedure :: set_precision => vtkhdf_file_set_precision
88 procedure :: set_subdivide => vtkhdf_file_set_subdivide
89 end type vtkhdf_file_t
90
91 integer, dimension(2), parameter :: vtkhdf_version = [2, 6]
92
93contains
94
95 ! -------------------------------------------------------------------------- !
96 ! Well defined subroutines
97
99 subroutine vtkhdf_file_set_overwrite(this, overwrite)
100 class(vtkhdf_file_t), intent(inout) :: this
101 logical, intent(in) :: overwrite
102 this%overwrite = overwrite
103 end subroutine vtkhdf_file_set_overwrite
104
106 subroutine vtkhdf_file_enable_amr(this)
107 class(vtkhdf_file_t), intent(inout) :: this
108 this%amr_enabled = .false.
109 end subroutine vtkhdf_file_enable_amr
110
112 subroutine vtkhdf_file_set_precision(this, precision)
113 class(vtkhdf_file_t), intent(inout) :: this
114 integer, intent(in) :: precision
115 this%precision = precision
116 end subroutine vtkhdf_file_set_precision
117
119 function vtkhdf_file_get_fname(this) result(base_fname)
120 class(vtkhdf_file_t), intent(in) :: this
121 character(len=1024) :: base_fname
122 character(len=1024) :: fname
123 character(len=1024) :: path, name, suffix
124
125 fname = trim(this%get_base_fname())
126 call filename_split(fname, path, name, suffix)
127
128 write(base_fname, '(A,A,"_",I0,A)') &
129 trim(path), trim(name), this%get_start_counter(), trim(suffix)
130
131 end function vtkhdf_file_get_fname
132
138 subroutine vtkhdf_file_set_subdivide(this, subdivide)
139 class(vtkhdf_file_t), intent(inout) :: this
140 logical, intent(in) :: subdivide
141 this%subdivide = subdivide
142 end subroutine vtkhdf_file_set_subdivide
143
144#ifdef HAVE_HDF5
145 ! -------------------------------------------------------------------------- !
146 ! HDF5 Required subroutines
147
150 subroutine vtkhdf_file_write(this, data, t)
151 class(vtkhdf_file_t), intent(inout) :: this
152 class(*), target, intent(in) :: data
153 real(kind=rp), intent(in), optional :: t
154 type(mesh_t), pointer :: msh
155 type(dofmap_t), pointer :: dof
156 type(field_ptr_t), allocatable :: fp(:)
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, H5T_NEKO_STRING
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
169 logical :: exists
170 integer :: counter
171
172 ! Determine mesh and field data
173 select type(data)
174 type is (field_t)
175 msh => data%msh
176 dof => data%dof
177 n_fields = 1
178 allocate(fp(1))
179 fp(1)%ptr => data
180 type is (field_list_t)
181 msh => data%msh(1)
182 dof => data%dof(1)
183 n_fields = data%size()
184 allocate(fp(n_fields))
185 do i = 1, n_fields
186 fp(i)%ptr => data%items(i)%ptr
187 end do
188 class default
189 call neko_error('Invalid data type for vtkhdf_file_write')
190 end select
191
192 ! Check conditions to ensure the input data is supported.
193 if (.not. associated(msh)) then
194 call neko_error('Mesh must be associated for vtkhdf_file_write')
195 end if
196 if (dof%Xh%lx < 2 .or. dof%Xh%ly < 2) then
197 call neko_error('VTKHDF linear output requires lx, ly >= 2')
198 end if
199 if (msh%gdim .eq. 3 .and. dof%Xh%lz < 2) then
200 call neko_error('VTKHDF linear output requires lz >= 2 in 3D')
201 end if
202 if (msh%gdim .lt. 2 .or. msh%gdim .gt. 3) then
203 call neko_error('VTKHDF output only supports 2D and 3D meshes')
204 end if
205
206 ! Ensure precision is set and are valid.
207 if (this%precision .gt. rp) then
208 this%precision = rp
209 call neko_warning('Requested precision is higher than working precision')
210 else if (this%precision .eq. -1) then
211 this%precision = rp
212 end if
213
214 call this%increment_counter()
215 fname = trim(this%get_vtkhdf_fname())
216 counter = this%get_counter() - this%get_start_counter()
217
218 mpi_info = mpi_info_null%mpi_val
219 mpi_comm = neko_comm%mpi_val
220
221 call h5open_f(ierr)
222 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
223 call h5pset_fapl_mpio_f(plist_id, mpi_comm, mpi_info, ierr)
224
225 if (counter .eq. 0) then
226 ! First write: always create a fresh file to avoid stale data
227 call h5fcreate_f(fname, h5f_acc_trunc_f, &
228 file_id, ierr, access_prp = plist_id)
229 else
230 call h5fopen_f(fname, h5f_acc_rdwr_f, file_id, ierr, &
231 access_prp = plist_id)
232 end if
233
234 ! Create/open VTKHDF root group with vtkhdf_version and type attributes
235 call h5lexists_f(file_id, "VTKHDF", exists, ierr)
236 if (exists) then
237 call h5gopen_f(file_id, "VTKHDF", vtkhdf_grp, ierr)
238 else
239 call h5gcreate_f(file_id, "VTKHDF", vtkhdf_grp, ierr)
240
241 ! Write Version attribute
242 vdims = 2
243 call h5screate_simple_f(1, vdims, filespace, ierr)
244 call h5acreate_f(vtkhdf_grp, "Version", h5t_native_integer, filespace, &
245 attr_id, ierr)
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)
249
250 ! Write Type attribute "UnstructuredGrid" as a fixed-length string
251 type_str = "UnstructuredGrid"
252 vdims = 1
253 call h5screate_f(h5s_scalar_f, filespace, ierr)
254
255 call h5tcopy_f(h5t_fortran_s1, h5t_neko_string, ierr)
256 call h5tset_size_f(h5t_neko_string, int(len_trim(type_str), kind=hsize_t), ierr)
257 call h5tset_strpad_f(h5t_neko_string, h5t_str_nullterm_f, ierr)
258
259 call h5acreate_f(vtkhdf_grp, "Type", h5t_neko_string, filespace, attr_id, ierr)
260 call h5awrite_f(attr_id, h5t_neko_string, [type_str], vdims, ierr)
261 call h5aclose_f(attr_id, ierr)
262
263 call h5tclose_f(h5t_neko_string, ierr)
264 call h5sclose_f(filespace, ierr)
265 end if
266
267 if (present(t)) then
268 call vtkhdf_write_steps(vtkhdf_grp, counter, t)
269 end if
270
271 if (associated(msh)) then
272 call vtkhdf_write_mesh(vtkhdf_grp, dof, msh, &
273 this%amr_enabled, counter, this%subdivide, t)
274 end if
275
276 ! Write field data in PointData group
277 if (n_fields > 0) then
278 call vtkhdf_write_pointdata(vtkhdf_grp, fp, this%precision, counter, &
279 fname, t)
280 end if
281
282 call h5gclose_f(vtkhdf_grp, ierr)
283 call h5pclose_f(plist_id, ierr)
284 call h5fclose_f(file_id, ierr)
285 call h5close_f(ierr)
286
287 if (allocated(fp)) deallocate(fp)
288
289 end subroutine vtkhdf_file_write
290
300 subroutine vtkhdf_write_mesh(vtkhdf_grp, dof, msh, amr, &
301 counter, subdivide, t)
302 type(dofmap_t), intent(in) :: dof
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
309
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
323 logical :: exists
324 integer, dimension(3) :: component_sizes
325 integer, dimension(3) :: component_offsets
326 integer, dimension(3) :: component_max_sizes
327
328 lx = dof%Xh%lx
329 ly = dof%Xh%ly
330 lz = dof%Xh%lz
331
332 if (subdivide .and. msh%gdim .eq. 3) then
333 vtk_cell_type = 12 ! VTK_HEXAHEDRON
334 cells_per_element = (lx - 1) * (ly - 1) * (lz - 1)
335 nodes_per_cell = 8
336 else if (subdivide .and. msh%gdim .eq. 2) then
337 vtk_cell_type = 9 ! VTK_QUAD
338 cells_per_element = (lx - 1) * (ly - 1)
339 nodes_per_cell = 4
340 else if (msh%gdim .eq. 3) then
341 vtk_cell_type = 72 ! VTK_LAGRANGE_HEXAHEDRON
342 cells_per_element = 1
343 nodes_per_cell = lx * ly * lz
344 else if (msh%gdim .eq. 2) then
345 vtk_cell_type = 70 ! VTK_LAGRANGE_QUADRILATERAL
346 cells_per_element = 1
347 nodes_per_cell = lx * ly
348 end if
349
350 ! --- Build the number of cells and the connectivity
351 local_points = dof%size()
352 local_cells = msh%nelv * cells_per_element
353 local_conn = local_cells * nodes_per_cell
354
355 total_points = dof%global_size()
356 total_cells = msh%glb_nelv * cells_per_element
357 total_conn = total_cells * nodes_per_cell
358
359 component_sizes = [local_points, local_cells, local_conn]
360 component_offsets = 0
361 component_max_sizes = 0
362
363 call mpi_exscan(component_sizes, component_offsets, 3, mpi_integer, &
364 mpi_sum, neko_comm, ierr)
365 call mpi_allreduce(component_sizes, component_max_sizes, 3, mpi_integer, &
366 mpi_max, neko_comm, ierr)
367
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)
374
375 offsets_offset = cell_offset + pe_rank
376 total_offsets = total_cells + pe_size
377
378 ! Create collective transfer property list
379 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
380 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
381
382 ! --- NumberOfPoints, NumberOfCells, NumberOfConnectivityIds ---
383 ! These datasets must accumulate nPieces entries per timestep,
384 ! giving a total size of nSteps * nPieces. VTK's reader computes
385 ! numberOfPieces = dims[0] / nSteps, so missing entries cause
386 ! garbage reads.
387 block
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)
392
393 nof_count(1) = 1_hsize_t
394 nof_offset(1) = int(counter, hsize_t) * int(pe_size, hsize_t) &
395 + int(pe_rank, 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)]
399
400 call h5pcreate_f(h5p_dataset_create_f, nof_dcpl, ierr)
401 call h5pset_chunk_f(nof_dcpl, 1, nof_chunk, ierr)
402
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)
412
413 call h5pclose_f(nof_dcpl, ierr)
414 end block
415
416 ! --- Points dataset (global coordinates) ---
417 call h5lexists_f(vtkhdf_grp, "Points", exists, ierr)
418 if (.not. exists) then
419
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)
426
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)
430
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)
436
437 block
438 real(kind=dp), allocatable :: coords(:,:)
439
440 allocate(coords(3, local_points))
441 do concurrent(local_idx = 1:local_points)
442 block
443 integer :: idx(4)
444 idx = nonlinear_index(local_idx, lx, ly, lz)
445
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))
449 end block
450 end do
451 call h5dwrite_f(dset_id, h5t_neko_double, coords, dcount2, ierr, &
452 file_space_id = filespace, mem_space_id = memspace, &
453 xfer_prp = xf_id)
454 deallocate(coords)
455 end block
456
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)
461 end if
462
463 ! --- Connectivity dataset ---
464 call h5lexists_f(vtkhdf_grp, "Connectivity", exists, ierr)
465 if (exists) call h5ldelete_f(vtkhdf_grp, "Connectivity", ierr)
466
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)
472
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)
476
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)
482
483 block
484 integer, allocatable :: connectivity(:)
485
486 allocate(connectivity(local_conn))
487 call vtkhdf_build_connectivity(connectivity, vtk_cell_type, msh, dof, &
488 subdivide)
489 call h5dwrite_f(dset_id, h5t_native_integer, connectivity, dcount, &
490 ierr, file_space_id = filespace, mem_space_id = memspace, &
491 xfer_prp = xf_id)
492 deallocate(connectivity)
493 end block
494
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)
499
500 ! --- Offsets dataset ---
501 call h5lexists_f(vtkhdf_grp, "Offsets", exists, ierr)
502 if (exists) call h5ldelete_f(vtkhdf_grp, "Offsets", ierr)
503
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)
509
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)
513
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)
519
520 block
521 integer, allocatable :: offsets(:)
522
523 allocate(offsets(local_cells + 1))
524 do concurrent(i = 1:local_cells)
525 offsets(i) = (i - 1) * nodes_per_cell
526 end do
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)
530 deallocate(offsets)
531 end block
532
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)
537
538 ! --- Types dataset (VTK cell types) ---
539 call h5lexists_f(vtkhdf_grp, "Types", exists, ierr)
540 if (exists) call h5ldelete_f(vtkhdf_grp, "Types", ierr)
541
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)
547
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)
551
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)
557
558 block
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)
564 end block
565
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)
570
571 if (present(t)) then
572 ! Open Steps group
573 call h5gopen_f(vtkhdf_grp, "Steps", grp_id, ierr)
574
575 ! --- NumberOfParts ---
576 call vtkhdf_write_i8_at(grp_id, "NumberOfParts", int(pe_size, kind=i8), &
577 counter)
578
579 ! --- PartOffsets ---
580 i8_value = 0
581 if (amr) then
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)
584
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)
587
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)
590
591 i8_value = int(counter - 1, kind=i8) * int(total_conn, kind=i8)
592 call vtkhdf_write_i8_at(grp_id, "ConnectivityIdOffsets", i8_value, &
593 counter)
594
595 else
596 i8_value = 0_i8
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, &
601 counter)
602 end if
603
604 call h5gclose_f(grp_id, ierr)
605 end if
606
607 call h5pclose_f(xf_id, ierr)
608
609 end subroutine vtkhdf_write_mesh
610
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
621
622 integer(hid_t) :: xf_id, H5T_NEKO_DOUBLE
623 integer :: ierr
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
630
631 ! Create collective transfer property list
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)
635
636 ! Create or open Steps group
637 call h5lexists_f(vtkhdf_grp, "Steps", exists, ierr)
638 if (exists) then
639 call h5gopen_f(vtkhdf_grp, "Steps", grp_id, ierr)
640 else
641 call h5gcreate_f(vtkhdf_grp, "Steps", grp_id, ierr)
642 end if
643
644 ! --- Values dataset (time values, real type) ---
645 call h5lexists_f(grp_id, "Values", exists, ierr)
646 if (exists) then
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, &
650 ierr)
651 call h5sclose_f(filespace, ierr)
652
653 ! We have not written this timestep yet, expand the array
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.")
659 end if
660 else
661 step_dims(1) = 1_hsize_t
662 step_maxdims(1) = h5s_unlimited_f
663 chunkdims(1) = 1_hsize_t
664
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)
672 end if
673
674 step_count(1) = 1_hsize_t
675 step_offset(1) = int(counter, hsize_t)
676
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)
681
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)
685
686 call h5sclose_f(memspace, ierr)
687 call h5sclose_f(filespace, ierr)
688 call h5dclose_f(dset_id, ierr)
689
690 ! --- NSteps attribute ---
691 ddim(1) = 1_hsize_t
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)
695 else
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)
700 end if
701
702 call h5awrite_f(attr_id, h5t_native_integer, counter + 1, ddim, ierr)
703
704 call h5aclose_f(attr_id, ierr)
705 call h5gclose_f(grp_id, ierr)
706 call h5pclose_f(xf_id, ierr)
707
708 end subroutine vtkhdf_write_steps
709
724 subroutine vtkhdf_write_pointdata(vtkhdf_grp, fp, precision, counter, &
725 fname, t)
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
732
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
737 integer :: n_fields
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
745
746 ! VDS and per-timestep external file variables
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, attr_id, H5T_NEKO_STRING
751 integer :: mpi_info, mpi_comm
752
753 ! Collected field info for VDS phase
754 integer :: fields_written
755 character(len=128), allocatable :: name_list(:)
756 logical, allocatable :: vector_list(:)
757
758 mpi_info = mpi_info_null%mpi_val
759 mpi_comm = neko_comm%mpi_val
760
761 n_fields = size(fp)
762
763 ! Compute local/global point counts and MPI offsets
764 local_points = fp(1)%ptr%dof%size()
765 total_points = fp(1)%ptr%dof%global_size()
766 point_offset = 0
767 call mpi_exscan(local_points, point_offset, 1, mpi_integer, &
768 mpi_sum, neko_comm, ierr)
769
770 ! Sync all the fields
771 do i = 1, n_fields
772 if (associated(fp(i)%ptr)) then
773 call fp(i)%ptr%copy_from(device_to_host, sync = i .eq. n_fields)
774 end if
775 end do
776
777 fields_written = 0
778 allocate(name_list(n_fields))
779 allocate(vector_list(n_fields))
780
781 ! Create PointData group if missing
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)
786 end if
787
788 ! ------------------------------------------------------------------------ !
789 ! Construct the target where data is written
790
791 if (present(t)) then
792
793 ! Derive base path from main filename for external files
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)
798
799 if (pe_rank .eq. 0) inquire(file = trim(ext_path), exist = exists)
800 call mpi_bcast(exists, 1, mpi_logical, 0, neko_comm, ierr)
801 if (.not. exists) then
802 if (pe_rank .eq. 0) then
803 call execute_command_line("mkdir -p '" // trim(ext_path) // "'")
804 end if
805 call mpi_barrier(neko_comm, ierr)
806 end if
807
808 call h5pcreate_f(h5p_file_access_f, ext_plist_id, ierr)
809 call h5pset_fapl_mpio_f(ext_plist_id, mpi_comm, mpi_info, ierr)
810 call h5fcreate_f(trim(ext_fname), h5f_acc_trunc_f, write_target, ierr, &
811 access_prp = ext_plist_id)
812 call h5pclose_f(ext_plist_id, ierr)
813
814 else
815 ! Non-temporal: write directly into the main file's PointData group
816 call h5gopen_f(vtkhdf_grp, "PointData", write_target, ierr)
817 end if
818
819 ! ------------------------------------------------------------------------ !
820 ! Write field data
821
822 do i = 1, n_fields
823 fld => fp(i)%ptr
824 field_name = fld%name
825 if (field_name .eq. 'p') field_name = 'Pressure'
826
827 ! Determine if this is a velocity component to group as a vector
828 is_vector = .false.
829 if (field_name .eq. 'u' .or. &
830 field_name .eq. 'v' .or. &
831 field_name .eq. 'w') then
832 u => null()
833 v => null()
834 w => null()
835 do j = 1, n_fields
836 select case (fp(j)%ptr%name)
837 case ('u')
838 u => fp(j)%ptr
839 case ('v')
840 v => fp(j)%ptr
841 case ('w')
842 w => fp(j)%ptr
843 end select
844 end do
845
846 if (associated(u) .and. associated(v) .and. associated(w)) then
847 is_vector = .true.
848 field_name = 'Velocity'
849 end if
850 end if
851
852 ! Skip duplicate fields (e.g. fluid_rho added by both fluid output
853 ! and field_writer when sharing the same output file)
854 exists = .false.
855 do j = 1, fields_written
856 if (trim(name_list(j)) .eq. trim(field_name)) then
857 exists = .true.
858 exit
859 end if
860 end do
861
862 ! Skip duplicate fields
863 if (exists) cycle
864
865 ! Track unique field names for VDS phase
866 fields_written = fields_written + 1
867 name_list(fields_written) = field_name
868 vector_list(fields_written) = is_vector
869
870 ! Write field data to the target
871 if (is_vector) then
872 call write_vector_field(write_target, field_name, u%x, v%x, w%x, &
873 local_points, precision, total_points, point_offset)
874 else
875 call write_scalar_field(write_target, field_name, fld%x, &
876 local_points, precision, total_points, point_offset)
877 end if
878 end do
879
880 ! Close write target
881 if (present(t)) then
882 call h5fclose_f(write_target, ierr)
883 else
884 call h5gclose_f(write_target, ierr)
885 end if
886
887 ! ------------------------------------------------------------------------ !
888 ! Manage temporal datasets through VDS
889
890 if (present(t)) then
891
892 ! Temporal: set up per-timestep external file as write target
893 call h5gopen_f(vtkhdf_grp, "Steps", step_grp_id, ierr)
894 time_offset = int(counter, kind=i8) * int(total_points, kind=i8)
895
896 ! Write PointDataOffsets under Steps for each unique field
897 call h5lexists_f(step_grp_id, "PointDataOffsets", exists, ierr)
898 if (exists) then
899 call h5gopen_f(step_grp_id, "PointDataOffsets", grp_id, ierr)
900 else
901 call h5gcreate_f(step_grp_id, "PointDataOffsets", grp_id, ierr)
902 end if
903 do i = 1, fields_written
904 call vtkhdf_write_i8_at(grp_id, trim(name_list(i)), &
905 time_offset, counter)
906 end do
907 call h5gclose_f(grp_id, ierr)
908 call h5gclose_f(step_grp_id, ierr)
909
910 ! ===== Create or extend VDS in main file =====
911 call h5gopen_f(vtkhdf_grp, "PointData", pointdata_grp, ierr)
912
913 do i = 1, fields_written
914 field_name = name_list(i)
915 is_vector = vector_list(i)
916
917 if (counter .eq. 0) then
918 ! First write: create VDS with pattern-based mapping
919 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
920 precision_hdf = h5kind_to_type(precision, h5_real_kind)
921
922 if (is_vector) then
923 pd_dims2 = [3_hsize_t, int(total_points, hsize_t)]
924 call h5screate_simple_f(2, pd_dims2, vds_src_space, ierr)
925 call h5sselect_all_f(vds_src_space, ierr)
926
927 pd_maxdims2 = [3_hsize_t, h5s_unlimited_f]
928 call h5screate_simple_f(2, pd_dims2, filespace, ierr, &
929 pd_maxdims2)
930
931 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
932 [0_hsize_t, 0_hsize_t], &
933 [1_hsize_t, h5s_unlimited_f], &
934 ierr, &
935 stride = [3_hsize_t, int(total_points, hsize_t)], &
936 block = [3_hsize_t, int(total_points, hsize_t)])
937
938 call h5pset_virtual_f(dcpl_id, filespace, trim(src_pattern), &
939 trim(field_name), vds_src_space, ierr)
940 call h5sclose_f(vds_src_space, ierr)
941
942 call h5dcreate_f(pointdata_grp, trim(field_name), &
943 precision_hdf, filespace, dset_id, ierr, &
944 dcpl_id = dcpl_id)
945 call h5sclose_f(filespace, ierr)
946 else
947 pd_dims1 = int(total_points, hsize_t)
948 call h5screate_simple_f(1, pd_dims1, vds_src_space, ierr)
949 call h5sselect_all_f(vds_src_space, ierr)
950
951 pd_maxdims1(1) = h5s_unlimited_f
952 call h5screate_simple_f(1, pd_dims1, filespace, ierr, &
953 pd_maxdims1)
954
955 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
956 [0_hsize_t], &
957 [h5s_unlimited_f], &
958 ierr, &
959 stride = [int(total_points, hsize_t)], &
960 block = [int(total_points, hsize_t)])
961
962 call h5pset_virtual_f(dcpl_id, filespace, trim(src_pattern), &
963 trim(field_name), vds_src_space, ierr)
964 call h5sclose_f(vds_src_space, ierr)
965
966 call h5dcreate_f(pointdata_grp, trim(field_name), &
967 precision_hdf, filespace, dset_id, ierr, &
968 dcpl_id = dcpl_id)
969 call h5sclose_f(filespace, ierr)
970 end if
971
972 call h5pclose_f(dcpl_id, ierr)
973 call h5dclose_f(dset_id, ierr)
974
975 else
976 ! Subsequent write: extend VDS to include new timestep
977 call h5dopen_f(pointdata_grp, trim(field_name), dset_id, ierr)
978
979 if (is_vector) then
980 pd_dims2 = [3_hsize_t, &
981 int(counter + 1, hsize_t) * int(total_points, hsize_t)]
982 call h5dset_extent_f(dset_id, pd_dims2, ierr)
983 else
984 pd_dims1 = int(counter + 1, hsize_t) * int(total_points, hsize_t)
985 call h5dset_extent_f(dset_id, pd_dims1, ierr)
986 end if
987
988 call h5dclose_f(dset_id, ierr)
989 end if
990 end do
991
992 call h5gclose_f(pointdata_grp, ierr)
993 end if
994
995 ! ------------------------------------------------------------------------ !
996 ! Add attributes to the PointData
997
998 do i = 1, fields_written
999 field_name = name_list(i)
1000 is_vector = vector_list(i)
1001 pd_dims1 = 1_hsize_t
1002
1003 call h5gopen_f(vtkhdf_grp, "PointData", pointdata_grp, ierr)
1004 call h5dopen_f(pointdata_grp, trim(field_name), dset_id, ierr)
1005
1006 call h5aexists_f(dset_id, "Attribute", exists, ierr)
1007 if (exists) then
1008 call h5dclose_f(dset_id, ierr)
1009 call h5gclose_f(pointdata_grp, ierr)
1010 cycle
1011 end if
1012
1013 ! Write the attribute what this PointData represents
1014 call h5screate_f(h5s_scalar_f, filespace, ierr)
1015
1016 call h5tcopy_f(h5t_fortran_s1, h5t_neko_string, ierr)
1017 call h5tset_size_f(h5t_neko_string, 6_hsize_t, ierr)
1018 call h5tset_strpad_f(h5t_neko_string, h5t_str_nullterm_f, ierr)
1019
1020 call h5acreate_f(dset_id, "Attribute", h5t_neko_string, filespace, &
1021 attr_id, ierr)
1022 if (is_vector) then
1023 call h5awrite_f(attr_id, h5t_neko_string, ["Vector"], pd_dims1, ierr)
1024 else
1025 call h5awrite_f(attr_id, h5t_neko_string, ["Scalar"], pd_dims1, ierr)
1026 end if
1027
1028 call h5aclose_f(attr_id, ierr)
1029 call h5tclose_f(h5t_neko_string, ierr)
1030 call h5sclose_f(filespace, ierr)
1031 call h5dclose_f(dset_id, ierr)
1032 call h5gclose_f(pointdata_grp, ierr)
1033
1034 end do
1035
1036 ! ------------------------------------------------------------------------ !
1037 ! Cleanup before returning
1038
1039 deallocate(name_list)
1040 deallocate(vector_list)
1041
1042 end subroutine vtkhdf_write_pointdata
1043
1044 ! -------------------------------------------------------------------------- !
1045 ! Helper functions and routines
1046
1059 subroutine vtkhdf_build_connectivity(conn, vtk_type, msh, dof, subdivide)
1060 integer, intent(inout) :: conn(:)
1061 integer(kind=1), intent(in) :: vtk_type
1062 type(mesh_t), intent(in) :: msh
1063 type(dofmap_t), intent(in) :: dof
1064 logical, intent(in) :: subdivide
1065 integer :: lx, ly, lz, nelv
1066 integer :: ie, ii, n_pts_per_elem, n_conn_per_elem
1067 integer, allocatable :: node_order(:)
1068
1069 nelv = msh%nelv
1070 lx = dof%Xh%lx
1071 ly = dof%Xh%ly
1072 lz = dof%Xh%lz
1073 n_pts_per_elem = lx * ly * lz
1074
1075 if (subdivide .and. vtk_type .eq. 12) then
1076 node_order = subdivide_to_hex_ordering(lx, ly, lz)
1077 else if (subdivide .and. vtk_type .eq. 9) then
1078 node_order = subdivide_to_quad_ordering(lx, ly)
1079 else
1080 node_order = vtk_ordering(vtk_type, lx, ly, lz)
1081 end if
1082
1083 n_conn_per_elem = size(node_order)
1084
1085 do concurrent(ie = 1:nelv, ii = 1:n_conn_per_elem)
1086 block
1087 integer :: idx, base
1088 idx = (ie - 1) * n_conn_per_elem
1089 base = (ie - 1) * n_pts_per_elem
1090 conn(idx + ii) = base + node_order(ii)
1091 end block
1092 end do
1093
1094 deallocate(node_order)
1095
1096 end subroutine vtkhdf_build_connectivity
1097
1110 subroutine vtkhdf_write_numberof(grp, dset_name, value, offset, cnt, &
1111 dcpl, index, xf_id, ierr)
1112 integer(hid_t), intent(in) :: grp, dcpl, xf_id
1113 character(len=*), intent(in) :: dset_name
1114 integer(kind=i8), intent(in) :: value
1115 integer, intent(in):: index
1116 integer(hsize_t), dimension(1), intent(in) :: offset, cnt
1117 integer, intent(out) :: ierr
1118
1119 integer(hid_t) :: dset_id, fspace, mspace
1120 integer(hsize_t), dimension(1) :: dims, maxdims
1121 integer(hid_t) :: H5T_NEKO_INTEGER
1122 logical :: exists
1123
1124 h5t_neko_integer = h5kind_to_type(i8, h5_integer_kind)
1125
1126 call h5lexists_f(grp, dset_name, exists, ierr)
1127 if (exists) then
1128 call h5dopen_f(grp, dset_name, dset_id, ierr)
1129 dims(1) = int(index + 1, hsize_t) * int(pe_size, hsize_t)
1130 call h5dset_extent_f(dset_id, dims, ierr)
1131 else
1132 dims(1) = int(pe_size, hsize_t)
1133 maxdims(1) = h5s_unlimited_f
1134 call h5screate_simple_f(1, dims, fspace, ierr, maxdims)
1135 call h5dcreate_f(grp, dset_name, h5t_neko_integer, &
1136 fspace, dset_id, ierr, dcpl_id = dcpl)
1137 call h5sclose_f(fspace, ierr)
1138 end if
1139
1140 call h5dget_space_f(dset_id, fspace, ierr)
1141 call h5screate_simple_f(1, cnt, mspace, ierr)
1142 call h5sselect_hyperslab_f(fspace, h5s_select_set_f, offset, cnt, ierr)
1143
1144 call h5dwrite_f(dset_id, h5t_neko_integer, value, cnt, ierr, &
1145 file_space_id = fspace, mem_space_id = mspace, xfer_prp = xf_id)
1146
1147 call h5sclose_f(mspace, ierr)
1148 call h5sclose_f(fspace, ierr)
1149 call h5dclose_f(dset_id, ierr)
1150 end subroutine vtkhdf_write_numberof
1151
1159 subroutine vtkhdf_write_i8_at(grp_id, name, value, index)
1160 integer(hid_t), intent(in) :: grp_id
1161 character(len=*), intent(in) :: name
1162 integer(kind=i8), intent(in) :: value
1163 integer, intent(in) :: index
1164
1165 integer :: ierr
1166 integer(hid_t) :: dset_id, dcpl_id, xf_id, filespace, memspace
1167 integer(hsize_t), dimension(1) :: dims, maxdims, count, offset, chunkdims
1168 integer(hid_t) :: H5T_NEKO_INTEGER
1169 logical :: exists
1170
1171 h5t_neko_integer = h5kind_to_type(i8, h5_integer_kind)
1172
1173 ! Create collective transfer property list
1174 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1175 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1176
1177 call h5lexists_f(grp_id, name, exists, ierr)
1178 if (exists) then
1179 call h5dopen_f(grp_id, name, dset_id, ierr)
1180 call h5dget_space_f(dset_id, filespace, ierr)
1181 call h5sget_simple_extent_dims_f(filespace, dims, maxdims, ierr)
1182 call h5sclose_f(filespace, ierr)
1183
1184 if (index .eq. dims(1)) then
1185 dims(1) = int(index + 1, hsize_t)
1186 call h5dset_extent_f(dset_id, dims, ierr)
1187 else if (index .gt. dims(1)) then
1188 call neko_error("VTKHDF: Values written out of order.")
1189 end if
1190 else
1191 dims = 1_hsize_t
1192 maxdims = h5s_unlimited_f
1193 chunkdims = 1_hsize_t
1194
1195 call h5screate_simple_f(1, dims, filespace, ierr, maxdims)
1196 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
1197 call h5pset_chunk_f(dcpl_id, 1, chunkdims, ierr)
1198 call h5dcreate_f(grp_id, name, h5t_neko_integer, &
1199 filespace, dset_id, ierr, dcpl_id = dcpl_id)
1200 call h5sclose_f(filespace, ierr)
1201 call h5pclose_f(dcpl_id, ierr)
1202 end if
1203
1204 count = 1_hsize_t
1205 offset = int(index, hsize_t)
1206
1207 call h5dget_space_f(dset_id, filespace, ierr)
1208 call h5screate_simple_f(1, count, memspace, ierr)
1209 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, offset, count, ierr)
1210 call h5dwrite_f(dset_id, h5t_neko_integer, value, count, ierr, &
1211 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
1212
1213 call h5sclose_f(memspace, ierr)
1214 call h5sclose_f(filespace, ierr)
1215 call h5dclose_f(dset_id, ierr)
1216 call h5pclose_f(xf_id, ierr)
1217
1218 end subroutine vtkhdf_write_i8_at
1219
1231 subroutine write_scalar_field(hdf_root, name, x, n_local, &
1232 precision, n_total, offset)
1233 integer, intent(in) :: n_local
1234 integer(hid_t), intent(in) :: hdf_root
1235 character(len=*), intent(in) :: name
1236 real(kind=rp), dimension(n_local), intent(in) :: x
1237 integer, intent(in), optional :: precision
1238 integer, intent(in), optional :: n_total, offset
1239
1240 integer(hsize_t), dimension(1) :: dims, dcount, doffset
1241 integer(hid_t) :: xf_id, dset_id, filespace, memspace, precision_hdf
1242 integer :: i, ierr, precision_local, n_tot, off
1243
1244 ! Setup data sizes, offsets and precision
1245 dcount = int(n_local, hsize_t)
1246
1247 if (present(n_total)) then
1248 dims = int(n_total, hsize_t)
1249 else
1250 call mpi_allreduce(n_local, n_tot, 1, mpi_integer, mpi_sum, neko_comm, &
1251 ierr)
1252 dims = int(n_tot, hsize_t)
1253 end if
1254
1255 if (present(offset)) then
1256 doffset = int(offset, hsize_t)
1257 else
1258 call mpi_exscan(n_local, off, 1, mpi_integer, mpi_sum, neko_comm, &
1259 ierr)
1260 doffset = int(off, hsize_t)
1261 end if
1262
1263 if (present(precision)) then
1264 precision_local = precision
1265 else
1266 precision_local = rp
1267 end if
1268 precision_hdf = h5kind_to_type(precision_local, h5_real_kind)
1269
1270 ! Prepare memory and filespaces
1271 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1272 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1273
1274 call h5screate_simple_f(1, dims, filespace, ierr)
1275 call h5screate_simple_f(1, dcount, memspace, ierr)
1276 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
1277 doffset, dcount, ierr)
1278
1279 call h5dcreate_f(hdf_root, trim(name), precision_hdf, &
1280 filespace, dset_id, ierr)
1281
1282 if (precision_local .eq. rp) then
1283 call h5dwrite_f(dset_id, precision_hdf, x, dcount, ierr, &
1284 file_space_id = filespace, mem_space_id = memspace, &
1285 xfer_prp = xf_id)
1286
1287 else if (precision_local .eq. sp) then
1288 block
1289 real(kind=sp), allocatable :: x_sp(:)
1290
1291 allocate(x_sp(n_local))
1292 do concurrent(i = 1:n_local)
1293 x_sp(i) = real(x(i), sp)
1294 end do
1295
1296 call h5dwrite_f(dset_id, precision_hdf, x_sp, dcount, ierr, &
1297 file_space_id = filespace, mem_space_id = memspace, &
1298 xfer_prp = xf_id)
1299 deallocate(x_sp)
1300 end block
1301
1302 else if (precision_local .eq. dp) then
1303 block
1304 real(kind=dp), allocatable :: x_dp(:)
1305
1306 allocate(x_dp(n_local))
1307 do concurrent(i = 1:n_local)
1308 x_dp(i) = real(x(i), dp)
1309 end do
1310
1311 call h5dwrite_f(dset_id, precision_hdf, x_dp, dcount, ierr, &
1312 file_space_id = filespace, mem_space_id = memspace, &
1313 xfer_prp = xf_id)
1314 deallocate(x_dp)
1315 end block
1316
1317 else if (precision_local .eq. qp) then
1318 block
1319 real(kind=qp), allocatable :: x_qp(:)
1320
1321 allocate(x_qp(n_local))
1322 do concurrent(i = 1:n_local)
1323 x_qp(i) = real(x(i), qp)
1324 end do
1325
1326 call h5dwrite_f(dset_id, precision_hdf, x_qp, dcount, ierr, &
1327 file_space_id = filespace, mem_space_id = memspace, &
1328 xfer_prp = xf_id)
1329 deallocate(x_qp)
1330 end block
1331
1332 else
1333 call neko_error("Unsupported precision in HDF5 write_scalar_field")
1334 end if
1335
1336 call h5sclose_f(filespace, ierr)
1337 call h5sclose_f(memspace, ierr)
1338 call h5dclose_f(dset_id, ierr)
1339 call h5pclose_f(xf_id, ierr)
1340 end subroutine write_scalar_field
1341
1355 subroutine write_vector_field(hdf_root, name, u, v, w, n_local, &
1356 precision, n_total, offset)
1357 integer, intent(in) :: n_local
1358 integer(hid_t), intent(in) :: hdf_root
1359 character(len=*), intent(in) :: name
1360 real(kind=rp), dimension(n_local), intent(in) :: u, v, w
1361 integer, intent(in), optional :: precision
1362 integer, intent(in), optional :: n_total, offset
1363
1364 integer(hsize_t), dimension(2) :: dims, dcount, doffset
1365 integer(hid_t) :: xf_id, dset_id, filespace, memspace, precision_hdf
1366 integer :: i, ierr, precision_local, n_tot, off
1367
1368 ! Setup data sizes, offsets and precision
1369 dcount = [3_hsize_t, int(n_local, hsize_t)]
1370
1371 if (present(n_total)) then
1372 dims = [3_hsize_t, int(n_total, hsize_t)]
1373 else
1374 call mpi_allreduce(n_local, n_tot, 1, mpi_integer, mpi_sum, neko_comm, &
1375 ierr)
1376 dims = [3_hsize_t, int(n_tot, hsize_t)]
1377 end if
1378
1379 if (present(offset)) then
1380 doffset = [0_hsize_t, int(offset, hsize_t)]
1381 else
1382 call mpi_exscan(n_local, off, 1, mpi_integer, mpi_sum, neko_comm, &
1383 ierr)
1384 doffset = [0_hsize_t, int(off, hsize_t)]
1385 end if
1386
1387 if (present(precision)) then
1388 precision_local = precision
1389 else
1390 precision_local = rp
1391 end if
1392 precision_hdf = h5kind_to_type(precision_local, h5_real_kind)
1393
1394 ! Prepare memory and filespaces
1395 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1396 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1397
1398 call h5screate_simple_f(2, dims, filespace, ierr)
1399 call h5screate_simple_f(2, dcount, memspace, ierr)
1400 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
1401 doffset, dcount, ierr)
1402
1403 call h5dcreate_f(hdf_root, trim(name), precision_hdf, filespace, dset_id, &
1404 ierr)
1405
1406 if (precision_local .eq. sp) then
1407 block
1408 real(kind=sp), allocatable :: f(:,:)
1409
1410 allocate(f(3, n_local))
1411 do concurrent(i = 1:n_local)
1412 f(1, i) = real(u(i), sp)
1413 f(2, i) = real(v(i), sp)
1414 f(3, i) = real(w(i), sp)
1415 end do
1416
1417 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1418 file_space_id = filespace, mem_space_id = memspace, &
1419 xfer_prp = xf_id)
1420 deallocate(f)
1421 end block
1422
1423 else if (precision_local .eq. dp) then
1424 block
1425 real(kind=dp), allocatable :: f(:,:)
1426
1427 allocate(f(3, n_local))
1428 do concurrent(i = 1:n_local)
1429 f(1, i) = real(u(i), dp)
1430 f(2, i) = real(v(i), dp)
1431 f(3, i) = real(w(i), dp)
1432 end do
1433
1434 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1435 file_space_id = filespace, mem_space_id = memspace, &
1436 xfer_prp = xf_id)
1437 deallocate(f)
1438 end block
1439
1440 else if (precision_local .eq. qp) then
1441 block
1442 real(kind=qp), allocatable :: f(:,:)
1443
1444 allocate(f(3, n_local))
1445 do concurrent(i = 1:n_local)
1446 f(1, i) = real(u(i), qp)
1447 f(2, i) = real(v(i), qp)
1448 f(3, i) = real(w(i), qp)
1449 end do
1450
1451 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1452 file_space_id = filespace, mem_space_id = memspace, &
1453 xfer_prp = xf_id)
1454 deallocate(f)
1455 end block
1456
1457 else
1458 call neko_error("Unsupported precision in HDF5 write_vector_field")
1459 end if
1460
1461 call h5sclose_f(filespace, ierr)
1462 call h5sclose_f(memspace, ierr)
1463 call h5dclose_f(dset_id, ierr)
1464 call h5pclose_f(xf_id, ierr)
1465 end subroutine write_vector_field
1466
1468 subroutine vtkhdf_file_read(this, data)
1469 class(vtkhdf_file_t) :: this
1470 class(*), target, intent(inout) :: data
1471
1472 call neko_error('VTKHDF file reading is not yet implemented')
1473
1474 end subroutine vtkhdf_file_read
1475
1476#else
1477 ! -------------------------------------------------------------------------- !
1478 ! Dummy functions and subroutines
1479
1481 subroutine vtkhdf_file_write(this, data, t)
1482 class(vtkhdf_file_t), intent(inout) :: this
1483 class(*), target, intent(in) :: data
1484 real(kind=rp), intent(in), optional :: t
1485 call neko_error('Neko needs to be built with HDF5 support')
1486 end subroutine vtkhdf_file_write
1487
1489 subroutine vtkhdf_file_read(this, data)
1490 class(vtkhdf_file_t) :: this
1491 class(*), target, intent(inout) :: data
1492 call neko_error('Neko needs to be built with HDF5 support')
1493 end subroutine vtkhdf_file_read
1494
1495#endif
1496
1497 ! -------------------------------------------------------------------------- !
1498 ! Sub-cell node ordering functions for VTK compatibility
1499
1508 pure function subdivide_to_hex_ordering(lx, ly, lz) result(node_order)
1509 integer, intent(in) :: lx, ly, lz
1510 integer :: node_order(8 * (lx - 1) * (ly - 1) * (lz - 1))
1511 integer :: ii, jj, kk, idx
1512
1513 idx = 0
1514
1515 do ii = 1, lx - 1
1516 do jj = 1, ly - 1
1517 do kk = 1, lz - 1
1518 node_order(idx + 1) = (kk - 1) * lx * ly + (jj - 1) * lx + ii - 1
1519 node_order(idx + 2) = (kk - 1) * lx * ly + (jj - 1) * lx + ii
1520 node_order(idx + 3) = (kk - 1) * lx * ly + jj * lx + ii
1521 node_order(idx + 4) = (kk - 1) * lx * ly + jj * lx + ii - 1
1522 node_order(idx + 5) = kk * lx * ly + (jj - 1) * lx + ii - 1
1523 node_order(idx + 6) = kk * lx * ly + (jj - 1) * lx + ii
1524 node_order(idx + 7) = kk * lx * ly + jj * lx + ii
1525 node_order(idx + 8) = kk * lx * ly + jj * lx + ii - 1
1526 idx = idx + 8
1527 end do
1528 end do
1529 end do
1530
1531 end function subdivide_to_hex_ordering
1532
1540 pure function subdivide_to_quad_ordering(lx, ly) result(node_order)
1541 integer, intent(in) :: lx, ly
1542 integer :: node_order(4 * (lx - 1) * (ly - 1))
1543 integer :: ii, jj, idx
1544
1545 idx = 0
1546
1547 do jj = 1, ly - 1
1548 do ii = 1, lx - 1
1549 node_order(idx + 1) = (jj - 1) * lx + ii - 1
1550 node_order(idx + 2) = (jj - 1) * lx + ii
1551 node_order(idx + 3) = jj * lx + ii
1552 node_order(idx + 4) = jj * lx + ii - 1
1553 idx = idx + 4
1554 end do
1555 end do
1556
1557 end function subdivide_to_quad_ordering
1558
1559end module vtkhdf_file
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
double real
Defines a checkpoint.
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public device_to_host
Definition device.F90:47
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public i8
Definition num_types.f90:7
integer, parameter, public qp
Definition num_types.f90:10
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Utilities.
Definition utils.f90:35
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, :)
Definition utils.f90:237
subroutine, public filename_split(fname, path, name, suffix)
Extract file name components.
Definition utils.f90:115
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:346
VTK Module containing utilities for VTK file handling.
Definition vtk.f90:42
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...
Definition vtk.f90:60
VTKHDF file format.
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
Definition field.f90:82
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 ...
A generic file handler.
Interface for HDF5 files.
#define max(a, b)
Definition tensor.cu:40