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
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, 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
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, 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)
258
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)
262
263 call h5tclose_f(memspace, 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
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 == 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) // "'")
803 end if
804 end if
805 call mpi_barrier(neko_comm, ierr)
806
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)
812
813 else
814 ! Non-temporal: write directly into the main file's PointData group
815 call h5gopen_f(vtkhdf_grp, "PointData", write_target, ierr)
816 end if
817
818 ! ------------------------------------------------------------------------ !
819 ! Write field data
820
821 do i = 1, n_fields
822 fld => fp(i)%ptr
823 field_name = fld%name
824 if (field_name .eq. 'p') field_name = 'Pressure'
825
826 ! Determine if this is a velocity component to group as a vector
827 is_vector = .false.
828 if (field_name .eq. 'u' .or. &
829 field_name .eq. 'v' .or. &
830 field_name .eq. 'w') then
831 u => null()
832 v => null()
833 w => null()
834 do j = 1, n_fields
835 select case (fp(j)%ptr%name)
836 case ('u')
837 u => fp(j)%ptr
838 case ('v')
839 v => fp(j)%ptr
840 case ('w')
841 w => fp(j)%ptr
842 end select
843 end do
844
845 if (associated(u) .and. associated(v) .and. associated(w)) then
846 is_vector = .true.
847 field_name = 'Velocity'
848 end if
849 end if
850
851 ! Skip duplicate fields (e.g. fluid_rho added by both fluid output
852 ! and field_writer when sharing the same output file)
853 exists = .false.
854 do j = 1, fields_written
855 if (trim(name_list(j)) .eq. trim(field_name)) then
856 exists = .true.
857 exit
858 end if
859 end do
860
861 ! Skip duplicate fields
862 if (exists) cycle
863
864 ! Track unique field names for VDS phase
865 fields_written = fields_written + 1
866 name_list(fields_written) = field_name
867 vector_list(fields_written) = is_vector
868
869 ! Write field data to the target
870 if (is_vector) then
871 call write_vector_field(write_target, field_name, u%x, v%x, w%x, &
872 local_points, precision, total_points, point_offset)
873 else
874 call write_scalar_field(write_target, field_name, fld%x, &
875 local_points, precision, total_points, point_offset)
876 end if
877 end do
878
879 ! Close write target
880 if (present(t)) then
881 call h5fclose_f(write_target, ierr)
882 else
883 call h5gclose_f(write_target, ierr)
884 end if
885
886 ! ------------------------------------------------------------------------ !
887 ! Manage temporal datasets through VDS
888
889 if (present(t)) then
890
891 ! Temporal: set up per-timestep external file as write target
892 call h5gopen_f(vtkhdf_grp, "Steps", step_grp_id, ierr)
893 time_offset = int(counter, kind=i8) * int(total_points, kind=i8)
894
895 ! Write PointDataOffsets under Steps for each unique field
896 call h5lexists_f(step_grp_id, "PointDataOffsets", exists, ierr)
897 if (exists) then
898 call h5gopen_f(step_grp_id, "PointDataOffsets", grp_id, ierr)
899 else
900 call h5gcreate_f(step_grp_id, "PointDataOffsets", grp_id, ierr)
901 end if
902 do i = 1, fields_written
903 call vtkhdf_write_i8_at(grp_id, trim(name_list(i)), &
904 time_offset, counter)
905 end do
906 call h5gclose_f(grp_id, ierr)
907 call h5gclose_f(step_grp_id, ierr)
908
909 ! ===== Create or extend VDS in main file =====
910 call h5gopen_f(vtkhdf_grp, "PointData", pointdata_grp, ierr)
911
912 do i = 1, fields_written
913 field_name = name_list(i)
914 is_vector = vector_list(i)
915
916 if (counter .eq. 0) then
917 ! First write: create VDS with pattern-based mapping
918 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
919 precision_hdf = h5kind_to_type(precision, h5_real_kind)
920
921 if (is_vector) then
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)
925
926 pd_maxdims2 = [3_hsize_t, h5s_unlimited_f]
927 call h5screate_simple_f(2, pd_dims2, filespace, ierr, &
928 pd_maxdims2)
929
930 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
931 [0_hsize_t, 0_hsize_t], &
932 [1_hsize_t, h5s_unlimited_f], &
933 ierr, &
934 stride = [3_hsize_t, int(total_points, hsize_t)], &
935 block = [3_hsize_t, int(total_points, hsize_t)])
936
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)
940
941 call h5dcreate_f(pointdata_grp, trim(field_name), &
942 precision_hdf, filespace, dset_id, ierr, &
943 dcpl_id = dcpl_id)
944 call h5sclose_f(filespace, ierr)
945 else
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)
949
950 pd_maxdims1(1) = h5s_unlimited_f
951 call h5screate_simple_f(1, pd_dims1, filespace, ierr, &
952 pd_maxdims1)
953
954 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
955 [0_hsize_t], &
956 [h5s_unlimited_f], &
957 ierr, &
958 stride = [int(total_points, hsize_t)], &
959 block = [int(total_points, hsize_t)])
960
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)
964
965 call h5dcreate_f(pointdata_grp, trim(field_name), &
966 precision_hdf, filespace, dset_id, ierr, &
967 dcpl_id = dcpl_id)
968 call h5sclose_f(filespace, ierr)
969 end if
970
971 call h5pclose_f(dcpl_id, ierr)
972 call h5dclose_f(dset_id, ierr)
973
974 else
975 ! Subsequent write: extend VDS to include new timestep
976 call h5dopen_f(pointdata_grp, trim(field_name), dset_id, ierr)
977
978 if (is_vector) then
979 pd_dims2 = [3_hsize_t, &
980 int((counter + 1) * total_points, hsize_t)]
981 call h5dset_extent_f(dset_id, pd_dims2, ierr)
982 else
983 pd_dims1 = int((counter + 1) * total_points, hsize_t)
984 call h5dset_extent_f(dset_id, pd_dims1, ierr)
985 end if
986
987 call h5dclose_f(dset_id, ierr)
988 end if
989 end do
990
991 call h5gclose_f(pointdata_grp, ierr)
992 end if
993
994 ! ------------------------------------------------------------------------ !
995 ! Cleanup before returning
996
997 deallocate(name_list)
998 deallocate(vector_list)
999
1000 end subroutine vtkhdf_write_pointdata
1001
1002 ! -------------------------------------------------------------------------- !
1003 ! Helper functions and routines
1004
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(:)
1026
1027 nelv = msh%nelv
1028 lx = dof%Xh%lx
1029 ly = dof%Xh%ly
1030 lz = dof%Xh%lz
1031 n_pts_per_elem = lx * ly * lz
1032
1033 if (subdivide .and. vtk_type .eq. 12) then
1034 node_order = subdivide_to_hex_ordering(lx, ly, lz)
1035 else if (subdivide .and. vtk_type .eq. 9) then
1036 node_order = subdivide_to_quad_ordering(lx, ly)
1037 else
1038 node_order = vtk_ordering(vtk_type, lx, ly, lz)
1039 end if
1040
1041 n_conn_per_elem = size(node_order)
1042
1043 do concurrent(ie = 1:nelv, ii = 1:n_conn_per_elem)
1044 block
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)
1049 end block
1050 end do
1051
1052 deallocate(node_order)
1053
1054 end subroutine vtkhdf_build_connectivity
1055
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
1076
1077 integer(hid_t) :: dset_id, fspace, mspace
1078 integer(hsize_t), dimension(1) :: dims, maxdims
1079 integer(hid_t) :: H5T_NEKO_INTEGER
1080 logical :: exists
1081
1082 h5t_neko_integer = h5kind_to_type(i8, h5_integer_kind)
1083
1084 call h5lexists_f(grp, dset_name, exists, ierr)
1085 if (exists) then
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)
1089 else
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)
1096 end if
1097
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)
1101
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)
1104
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
1109
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
1122
1123 integer :: ierr
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
1127 logical :: exists
1128
1129 h5t_neko_integer = h5kind_to_type(i8, h5_integer_kind)
1130
1131 ! Create collective transfer property list
1132 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1133 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1134
1135 call h5lexists_f(grp_id, name, exists, ierr)
1136 if (exists) then
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)
1141
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.")
1147 end if
1148 else
1149 dims = 1_hsize_t
1150 maxdims = h5s_unlimited_f
1151 chunkdims = 1_hsize_t
1152
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)
1160 end if
1161
1162 count = 1_hsize_t
1163 offset = int(index, hsize_t)
1164
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)
1170
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)
1175
1176 end subroutine vtkhdf_write_i8_at
1177
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
1197
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
1201
1202 ! Setup data sizes, offsets and precision
1203 dcount = int(n_local, hsize_t)
1204
1205 if (present(n_total)) then
1206 dims = int(n_total, hsize_t)
1207 else
1208 call mpi_allreduce(n_local, n_tot, 1, mpi_integer, mpi_sum, neko_comm, &
1209 ierr)
1210 dims = int(n_tot, hsize_t)
1211 end if
1212
1213 if (present(offset)) then
1214 doffset = int(offset, hsize_t)
1215 else
1216 call mpi_exscan(n_local, off, 1, mpi_integer, mpi_sum, neko_comm, &
1217 ierr)
1218 doffset = int(off, hsize_t)
1219 end if
1220
1221 if (present(precision)) then
1222 precision_local = precision
1223 else
1224 precision_local = rp
1225 end if
1226 precision_hdf = h5kind_to_type(precision_local, h5_real_kind)
1227
1228 ! Prepare memory and filespaces
1229 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1230 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1231
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)
1236
1237 call h5dcreate_f(hdf_root, trim(name), precision_hdf, &
1238 filespace, dset_id, ierr)
1239
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, &
1243 xfer_prp = xf_id)
1244
1245 else if (precision_local .eq. sp) then
1246 block
1247 real(kind=sp), allocatable :: x_sp(:)
1248
1249 allocate(x_sp(n_local))
1250 do concurrent(i = 1:n_local)
1251 x_sp(i) = real(x(i), sp)
1252 end do
1253
1254 call h5dwrite_f(dset_id, precision_hdf, x_sp, dcount, ierr, &
1255 file_space_id = filespace, mem_space_id = memspace, &
1256 xfer_prp = xf_id)
1257 deallocate(x_sp)
1258 end block
1259
1260 else if (precision_local .eq. dp) then
1261 block
1262 real(kind=dp), allocatable :: x_dp(:)
1263
1264 allocate(x_dp(n_local))
1265 do concurrent(i = 1:n_local)
1266 x_dp(i) = real(x(i), dp)
1267 end do
1268
1269 call h5dwrite_f(dset_id, precision_hdf, x_dp, dcount, ierr, &
1270 file_space_id = filespace, mem_space_id = memspace, &
1271 xfer_prp = xf_id)
1272 deallocate(x_dp)
1273 end block
1274
1275 else if (precision_local .eq. qp) then
1276 block
1277 real(kind=qp), allocatable :: x_qp(:)
1278
1279 allocate(x_qp(n_local))
1280 do concurrent(i = 1:n_local)
1281 x_qp(i) = real(x(i), qp)
1282 end do
1283
1284 call h5dwrite_f(dset_id, precision_hdf, x_qp, dcount, ierr, &
1285 file_space_id = filespace, mem_space_id = memspace, &
1286 xfer_prp = xf_id)
1287 deallocate(x_qp)
1288 end block
1289
1290 else
1291 call neko_error("Unsupported precision in HDF5 write_scalar_field")
1292 end if
1293
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
1299
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
1321
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
1325
1326 ! Setup data sizes, offsets and precision
1327 dcount = [3_hsize_t, int(n_local, hsize_t)]
1328
1329 if (present(n_total)) then
1330 dims = [3_hsize_t, int(n_total, hsize_t)]
1331 else
1332 call mpi_allreduce(n_local, n_tot, 1, mpi_integer, mpi_sum, neko_comm, &
1333 ierr)
1334 dims = [3_hsize_t, int(n_tot, hsize_t)]
1335 end if
1336
1337 if (present(offset)) then
1338 doffset = [0_hsize_t, int(offset, hsize_t)]
1339 else
1340 call mpi_exscan(n_local, off, 1, mpi_integer, mpi_sum, neko_comm, &
1341 ierr)
1342 doffset = [0_hsize_t, int(off, hsize_t)]
1343 end if
1344
1345 if (present(precision)) then
1346 precision_local = precision
1347 else
1348 precision_local = rp
1349 end if
1350 precision_hdf = h5kind_to_type(precision_local, h5_real_kind)
1351
1352 ! Prepare memory and filespaces
1353 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1354 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1355
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)
1360
1361 call h5dcreate_f(hdf_root, trim(name), precision_hdf, &
1362 filespace, dset_id, ierr)
1363
1364 if (precision_local .eq. sp) then
1365 block
1366 real(kind=sp), allocatable :: f(:,:)
1367
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)
1373 end do
1374
1375 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1376 file_space_id = filespace, mem_space_id = memspace, &
1377 xfer_prp = xf_id)
1378 deallocate(f)
1379 end block
1380
1381 else if (precision_local .eq. dp) then
1382 block
1383 real(kind=dp), allocatable :: f(:,:)
1384
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)
1390 end do
1391
1392 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1393 file_space_id = filespace, mem_space_id = memspace, &
1394 xfer_prp = xf_id)
1395 deallocate(f)
1396 end block
1397
1398 else if (precision_local .eq. qp) then
1399 block
1400 real(kind=qp), allocatable :: f(:,:)
1401
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)
1407 end do
1408
1409 call h5dwrite_f(dset_id, precision_hdf, f, dcount, ierr, &
1410 file_space_id = filespace, mem_space_id = memspace, &
1411 xfer_prp = xf_id)
1412 deallocate(f)
1413 end block
1414
1415 else
1416 call neko_error("Unsupported precision in HDF5 write_vector_field")
1417 end if
1418
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
1424
1426 subroutine vtkhdf_file_read(this, data)
1427 class(vtkhdf_file_t) :: this
1428 class(*), target, intent(inout) :: data
1429
1430 call neko_error('VTKHDF file reading is not yet implemented')
1431
1432 end subroutine vtkhdf_file_read
1433
1434#else
1435 ! -------------------------------------------------------------------------- !
1436 ! Dummy functions and subroutines
1437
1439 subroutine vtkhdf_file_write(this, data, t)
1440 class(vtkhdf_file_t), intent(inout) :: this
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')
1444 end subroutine vtkhdf_file_write
1445
1447 subroutine vtkhdf_file_read(this, data)
1448 class(vtkhdf_file_t) :: this
1449 class(*), target, intent(inout) :: data
1450 call neko_error('Neko needs to be built with HDF5 support')
1451 end subroutine vtkhdf_file_read
1452
1453#endif
1454
1455 ! -------------------------------------------------------------------------- !
1456 ! Sub-cell node ordering functions for VTK compatibility
1457
1466 pure function subdivide_to_hex_ordering(lx, ly, lz) result(node_order)
1467 integer, intent(in) :: lx, ly, lz
1468 integer :: node_order(8 * (lx - 1) * (ly - 1) * (lz - 1))
1469 integer :: ii, jj, kk, idx
1470
1471 idx = 0
1472
1473 do ii = 1, lx - 1
1474 do jj = 1, ly - 1
1475 do kk = 1, lz - 1
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
1484 idx = idx + 8
1485 end do
1486 end do
1487 end do
1488
1489 end function subdivide_to_hex_ordering
1490
1498 pure function subdivide_to_quad_ordering(lx, ly) result(node_order)
1499 integer, intent(in) :: lx, ly
1500 integer :: node_order(4 * (lx - 1) * (ly - 1))
1501 integer :: ii, jj, idx
1502
1503 idx = 0
1504
1505 do jj = 1, ly - 1
1506 do ii = 1, lx - 1
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
1511 idx = idx + 4
1512 end do
1513 end do
1514
1515 end function subdivide_to_quad_ordering
1516
1517end 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