Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
hdf5_file.F90
Go to the documentation of this file.
1! Copyright (c) 2024-2025, 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, dp, sp
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
46 use vector, only : vector_t
47 use matrix, only : matrix_t
48 use datadist, only : linear_dist_t
49 use comm, only : pe_rank, pe_size, neko_comm
50 use mpi_f08, only : mpi_info_null, mpi_allreduce, mpi_allgather, &
51 mpi_in_place, mpi_integer, mpi_sum, mpi_max, mpi_comm_size, mpi_exscan, &
52 mpi_barrier, mpi_integer8, mpi_scan
53#ifdef HAVE_HDF5
54 use hdf5
55#endif
56 implicit none
57 private
58
60 type, public, extends(generic_file_t) :: hdf5_file_t
61
62 ! HDF5 members
63#ifdef HAVE_HDF5
64 integer(hid_t) :: file_id = -1_hid_t
65 integer(hid_t) :: active_group_id = -1_hid_t
66 integer(hid_t) :: plist_id = -1_hid_t
67#endif
68 character(len=1) :: mode
69 integer :: precision = -1
70 integer :: offset = 0
71 integer :: count = 0
72
73 contains
74 ! General methods for reading/writing HDF5 files
75 procedure :: read => hdf5_file_read
76 procedure :: write => hdf5_file_write
77 procedure :: set_overwrite => hdf5_file_set_overwrite
78 ! Granular methods for dealing with HDF5 files
79 procedure :: open => hdf5_file_open
80 procedure :: close => hdf5_file_close
81 procedure :: set_active_group => hdf5_file_set_group
82 procedure :: set_precision => hdf5_file_set_precision
83 procedure, pass(this) :: write_vector => hdf5_file_write_vector
84 procedure, pass(this) :: write_matrix => hdf5_file_write_matrix
85 procedure, pass(this) :: write_field => hdf5_file_write_field
86 procedure, pass(this) :: write_int_attribute => &
88 procedure, pass(this) :: write_rp_attribute => hdf5_file_write_rp_attribute
89 procedure, pass(this) :: read_vector => hdf5_file_read_vector
90 procedure, pass(this) :: read_matrix => hdf5_file_read_matrix
91 procedure, pass(this) :: read_int_attribute => hdf5_file_read_int_attribute
92 procedure, pass(this) :: read_rp_attribute => hdf5_file_read_rp_attribute
93 procedure :: write_dataset => hdf5_file_write_dataset
94 procedure :: read_dataset => hdf5_file_read_dataset
95 procedure :: write_attribute => hdf5_file_write_attribute
96 procedure :: read_attribute => hdf5_file_read_attribute
97 end type hdf5_file_t
98
99contains
100
102 subroutine hdf5_file_set_overwrite(this, overwrite)
103 class(hdf5_file_t), intent(inout) :: this
104 logical, intent(in) :: overwrite
105 this%overwrite = overwrite
106 end subroutine hdf5_file_set_overwrite
107
109 function file_get_fname(this) result(base_fname)
110 class(hdf5_file_t), intent(in) :: this
111 character(len=1024) :: base_fname
112 character(len=1024) :: fname
113 character(len=1024) :: path, name, suffix
114
115 fname = trim(this%get_base_fname())
116 call filename_split(fname, path, name, suffix)
117
118 ! Append a counter
119 !write(base_fname, '(A,A,"_",I0,A)') &
120 ! trim(path), trim(name), this%get_start_counter(), trim(suffix)
121
122 ! Do not append anything
123 base_fname = trim(fname)
124
125 end function file_get_fname
126
128 subroutine hdf5_file_set_precision(this, precision)
129 class(hdf5_file_t), intent(inout) :: this
130 integer, intent(in) :: precision
131 this%precision = precision
132 end subroutine hdf5_file_set_precision
133
134
135#ifdef HAVE_HDF5
136
137 ! ===============
138 ! General methods
139 ! ===============
140
142 subroutine hdf5_file_write(this, data, t)
143 class(hdf5_file_t), intent(inout) :: this
144 class(*), target, intent(in) :: data
145 real(kind=rp), intent(in), optional :: t
146 type(mesh_t), pointer :: msh
147 type(dofmap_t), pointer :: dof
148 type(field_ptr_t), allocatable :: fp(:)
149 type(field_series_ptr_t), allocatable :: fsp(:)
150 real(kind=rp), pointer :: dtlag(:)
151 real(kind=rp), pointer :: tlag(:)
152 integer :: ierr, info, drank, i, j
153 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
154 integer(hid_t) :: filespace, memspace
155 integer(hid_t) :: H5T_NEKO_REAL
156 integer(hsize_t), dimension(1) :: ddim, dcount, doffset
157 integer :: suffix_pos
158 character(len=5) :: id_str
159 character(len=1024) :: fname
160
161 call hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
162
163 if (.not. this%overwrite) call this%increment_counter()
164 fname = trim(this%get_fname())
165
166 call h5open_f(ierr)
167
168 call hdf5_file_determine_real(h5t_neko_real)
169
170 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
171 info = mpi_info_null%mpi_val
172 call h5pset_fapl_mpio_f(plist_id, neko_comm%mpi_val, info, ierr)
173
174 call h5fcreate_f(fname, h5f_acc_trunc_f, &
175 file_id, ierr, access_prp = plist_id)
176
177 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
178 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
179
180 call h5screate_f(h5s_scalar_f, filespace, ierr)
181 ddim = 1
182
183 if (present(t)) then
184 call h5acreate_f(file_id, "Time", h5t_neko_real, filespace, attr_id, &
185 ierr, h5p_default_f, h5p_default_f)
186 call h5awrite_f(attr_id, h5t_neko_real, t, ddim, ierr)
187 call h5aclose_f(attr_id, ierr)
188 end if
189
190 if (associated(dof)) then
191 call h5acreate_f(file_id, "Lx", h5t_native_integer, filespace, attr_id, &
192 ierr, h5p_default_f, h5p_default_f)
193 call h5awrite_f(attr_id, h5t_native_integer, dof%Xh%lx, ddim, ierr)
194 call h5aclose_f(attr_id, ierr)
195 end if
196
197 if (associated(msh)) then
198 call h5gcreate_f(file_id, "Mesh", grp_id, ierr, &
199 lcpl_id = h5p_default_f, gcpl_id = h5p_default_f, &
200 gapl_id = h5p_default_f)
201
202 call h5acreate_f(grp_id, "Elements", h5t_native_integer, filespace, &
203 attr_id, ierr, h5p_default_f, h5p_default_f)
204 call h5awrite_f(attr_id, h5t_native_integer, msh%glb_nelv, ddim, ierr)
205 call h5aclose_f(attr_id, ierr)
206
207 call h5acreate_f(grp_id, "Dimension", h5t_native_integer, filespace, &
208 attr_id, ierr, h5p_default_f, h5p_default_f)
209 call h5awrite_f(attr_id, h5t_native_integer, msh%gdim, ddim, ierr)
210 call h5aclose_f(attr_id, ierr)
211
212 call h5gclose_f(grp_id, ierr)
213 end if
214
215
216 call h5sclose_f(filespace, ierr)
217
218 !
219 ! Write restart group (tlag, dtlag)
220 !
221 if (associated(tlag) .and. associated(dtlag)) then
222 call h5gcreate_f(file_id, "Restart", grp_id, ierr, &
223 lcpl_id = h5p_default_f, gcpl_id = h5p_default_f, &
224 gapl_id = h5p_default_f)
225
226 drank = 1
227 ddim = size(tlag)
228 doffset(1) = 0
229 if (pe_rank .eq. 0) then
230 dcount = size(tlag)
231 else
232 dcount = 0
233 end if
234
235 call h5screate_simple_f(drank, ddim, filespace, ierr)
236
237 call h5dcreate_f(grp_id, 'tlag', h5t_neko_real, &
238 filespace, dset_id, ierr)
239 call h5dget_space_f(dset_id, filespace, ierr)
240 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
241 doffset, dcount, ierr)
242 call h5dwrite_f(dset_id, h5t_neko_real, tlag, &
243 ddim, ierr, xfer_prp = plist_id)
244 call h5dclose_f(dset_id, ierr)
245
246 call h5dcreate_f(grp_id, 'dtlag', h5t_neko_real, &
247 filespace, dset_id, ierr)
248 call h5dget_space_f(dset_id, filespace, ierr)
249 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
250 doffset, dcount, ierr)
251 call h5dwrite_f(dset_id, h5t_neko_real, dtlag, &
252 ddim, ierr, xfer_prp = plist_id)
253 call h5dclose_f(dset_id, ierr)
254
255 call h5sclose_f(filespace, ierr)
256 call h5gclose_f(grp_id, ierr)
257
258 end if
259
260
261 !
262 ! Write fields group
263 !
264 if (allocated(fp) .or. allocated(fsp)) then
265 call h5gcreate_f(file_id, "Fields", grp_id, ierr, &
266 lcpl_id = h5p_default_f, gcpl_id = h5p_default_f, &
267 gapl_id = h5p_default_f)
268
269 dcount(1) = int(dof%size(), 8)
270 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
271 ddim = int(dof%size(), 8)
272 drank = 1
273 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
274 mpi_integer8, mpi_sum, neko_comm, ierr)
275
276 call h5screate_simple_f(drank, ddim, filespace, ierr)
277 call h5screate_simple_f(drank, dcount, memspace, ierr)
278
279
280 if (allocated(fp)) then
281 do i = 1, size(fp)
282 call h5dcreate_f(grp_id, fp(i)%ptr%name, h5t_neko_real, &
283 filespace, dset_id, ierr)
284 call h5dget_space_f(dset_id, filespace, ierr)
285 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
286 doffset, dcount, ierr)
287 call h5dwrite_f(dset_id, h5t_neko_real, &
288 fp(i)%ptr%x(1,1,1,1), &
289 ddim, ierr, file_space_id = filespace, &
290 mem_space_id = memspace, xfer_prp = plist_id)
291 call h5dclose_f(dset_id, ierr)
292 end do
293 deallocate(fp)
294 end if
295
296 if (allocated(fsp)) then
297 do i = 1, size(fsp)
298 do j = 1, fsp(i)%ptr%size()
299 call h5dcreate_f(grp_id, fsp(i)%ptr%lf(j)%name, &
300 h5t_neko_real, filespace, dset_id, ierr)
301 call h5dget_space_f(dset_id, filespace, ierr)
302 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
303 doffset, dcount, ierr)
304 call h5dwrite_f(dset_id, h5t_neko_real, &
305 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
306 ddim, ierr, file_space_id = filespace, &
307 mem_space_id = memspace, xfer_prp = plist_id)
308 call h5dclose_f(dset_id, ierr)
309 end do
310 end do
311 deallocate(fsp)
312 end if
313
314 call h5gclose_f(grp_id, ierr)
315 call h5sclose_f(filespace, ierr)
316 call h5sclose_f(memspace, ierr)
317 end if
318
319 call h5pclose_f(plist_id, ierr)
320 call h5fclose_f(file_id, ierr)
321
322 call h5close_f(ierr)
323
324 end subroutine hdf5_file_write
325
327 subroutine hdf5_file_read(this, data)
328 class(hdf5_file_t) :: this
329 class(*), target, intent(inout) :: data
330 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
331 integer(hid_t) :: filespace, memspace
332 integer(hid_t) :: H5T_NEKO_REAL
333 integer(hsize_t), dimension(1) :: ddim, dcount, doffset
334 integer :: i,j, ierr, info, glb_nelv, gdim, lx, drank
335 type(mesh_t), pointer :: msh
336 type(dofmap_t), pointer :: dof
337 type(field_ptr_t), allocatable :: fp(:)
338 type(field_series_ptr_t), allocatable :: fsp(:)
339 real(kind=rp), pointer :: dtlag(:)
340 real(kind=rp), pointer :: tlag(:)
341 real(kind=rp) :: t
342 character(len=1024) :: fname
343
344 fname = trim(this%get_fname())
345
346 call hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
347
348 call h5open_f(ierr)
349
350 call hdf5_file_determine_real(h5t_neko_real)
351
352 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
353 info = mpi_info_null%mpi_val
354 call h5pset_fapl_mpio_f(plist_id, neko_comm%mpi_val, info, ierr)
355
356 call h5fopen_f(fname, h5f_acc_rdonly_f, &
357 file_id, ierr, access_prp = plist_id)
358
359 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
360 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
361
362 ddim = 1
363 call h5aopen_name_f(file_id, 'Time', attr_id, ierr)
364 call h5aread_f(attr_id, h5t_neko_real, t, ddim, ierr)
365 call h5aclose_f(attr_id, ierr)
366
367 select type (data)
368 type is (chkp_t)
369 data%t = t
370 end select
371
372 call h5aopen_name_f(file_id, 'Lx', attr_id, ierr)
373 call h5aread_f(attr_id, h5t_native_integer, lx, ddim, ierr)
374 call h5aclose_f(attr_id, ierr)
375
376 call h5gopen_f(file_id, 'Mesh', grp_id, ierr, gapl_id = h5p_default_f)
377
378 call h5aopen_name_f(grp_id, 'Elements', attr_id, ierr)
379 call h5aread_f(attr_id, h5t_native_integer, glb_nelv, ddim, ierr)
380 call h5aclose_f(attr_id, ierr)
381
382 call h5aopen_name_f(grp_id, 'Dimension', attr_id, ierr)
383 call h5aread_f(attr_id, h5t_native_integer, gdim, ddim, ierr)
384 call h5aclose_f(attr_id, ierr)
385 call h5gclose_f(grp_id, ierr)
386
387
388 if (associated(tlag) .and. associated(dtlag)) then
389 drank = 1
390 ddim = size(tlag)
391 doffset(1) = 0
392 if (pe_rank .eq. 0) then
393 dcount = size(tlag)
394 else
395 dcount = 0
396 end if
397
398 call h5gopen_f(file_id, 'Restart', grp_id, ierr, &
399 gapl_id = h5p_default_f)
400 call h5dopen_f(grp_id, 'tlag', dset_id, ierr)
401 call h5dget_space_f(dset_id, filespace, ierr)
402 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
403 doffset, dcount, ierr)
404 call h5dread_f(dset_id, h5t_neko_real, tlag, ddim, ierr, &
405 xfer_prp = plist_id)
406 call h5dclose_f(dset_id, ierr)
407 call h5sclose_f(filespace, ierr)
408
409 call h5dopen_f(grp_id, 'dtlag', dset_id, ierr)
410 call h5dget_space_f(dset_id, filespace, ierr)
411 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
412 doffset, dcount, ierr)
413 call h5dread_f(dset_id, h5t_neko_real, dtlag, ddim, ierr, &
414 xfer_prp = plist_id)
415 call h5dclose_f(dset_id, ierr)
416 call h5sclose_f(filespace, ierr)
417
418 call h5gclose_f(grp_id, ierr)
419 end if
420
421 if (allocated(fp) .or. allocated(fsp)) then
422 call h5gopen_f(file_id, 'Fields', grp_id, ierr, gapl_id = h5p_default_f)
423
424 dcount(1) = int(dof%size(), 8)
425 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
426 ddim = int(dof%size(), 8)
427 drank = 1
428
429 dcount(1) = int(dof%size(), 8)
430 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
431 ddim = int(dof%size(), 8)
432 drank = 1
433 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
434 mpi_integer8, mpi_sum, neko_comm, ierr)
435
436 call h5screate_simple_f(drank, dcount, memspace, ierr)
437
438 if (allocated(fp)) then
439 do i = 1, size(fp)
440 call h5dopen_f(grp_id, fp(i)%ptr%name, dset_id, ierr)
441 call h5dget_space_f(dset_id, filespace, ierr)
442 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
443 doffset, dcount, ierr)
444 call h5dread_f(dset_id, h5t_neko_real, &
445 fp(i)%ptr%x(1,1,1,1), &
446 ddim, ierr, file_space_id = filespace, &
447 mem_space_id = memspace, xfer_prp = plist_id)
448 call h5dclose_f(dset_id, ierr)
449 call h5sclose_f(filespace, ierr)
450 end do
451 end if
452
453 if (allocated(fsp)) then
454 do i = 1, size(fsp)
455 do j = 1, fsp(i)%ptr%size()
456 call h5dopen_f(grp_id, fsp(i)%ptr%lf(j)%name, dset_id, ierr)
457 call h5dget_space_f(dset_id, filespace, ierr)
458 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
459 doffset, dcount, ierr)
460 call h5dread_f(dset_id, h5t_neko_real, &
461 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
462 ddim, ierr, file_space_id = filespace, &
463 mem_space_id = memspace, xfer_prp = plist_id)
464 call h5dclose_f(dset_id, ierr)
465 call h5sclose_f(filespace, ierr)
466 end do
467 end do
468 end if
469 call h5sclose_f(memspace, ierr)
470 call h5gclose_f(grp_id, ierr)
471 end if
472
473 call h5pclose_f(plist_id, ierr)
474 call h5fclose_f(file_id, ierr)
475
476 call h5close_f(ierr)
477
478 end subroutine hdf5_file_read
479
480
481 subroutine hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
482 class(*), target, intent(in) :: data
483 type(mesh_t), pointer, intent(inout) :: msh
484 type(dofmap_t), pointer, intent(inout) :: dof
485 type(field_ptr_t), allocatable, intent(inout) :: fp(:)
486 type(field_series_ptr_t), allocatable, intent(inout) :: fsp(:)
487 real(kind=rp), pointer, intent(inout) :: dtlag(:)
488 real(kind=rp), pointer, intent(inout) :: tlag(:)
489 integer :: i, j, fp_size, fp_cur, fsp_size, fsp_cur, scalar_count, ab_count
490 character(len=32) :: scalar_name
491
492 select type (data)
493 type is (field_t)
494 dof => data%dof
495 msh => data%msh
496 fp_size = 1
497 allocate(fp(fp_size))
498 fp(1)%ptr => data
499
500 nullify(dtlag)
501 nullify(tlag)
502
503 type is (field_list_t)
504
505 if (data%size() .gt. 0) then
506 allocate(fp(data%size()))
507
508 dof => data%dof(1)
509 msh => data%msh(1)
510
511 do i = 1, data%size()
512 fp(i)%ptr => data%items(i)%ptr
513 end do
514 else
515 call neko_error('Empty field list')
516 end if
517
518 nullify(dtlag)
519 nullify(tlag)
520
521 type is (chkp_t)
522
523 if ( .not. associated(data%u) .or. &
524 .not. associated(data%v) .or. &
525 .not. associated(data%w) .or. &
526 .not. associated(data%p) ) then
527 call neko_error('Checkpoint not initialized')
528 end if
529
530 fp_size = 4
531
532 if (allocated(data%scalar_lags%items) .and. &
533 data%scalar_lags%size() > 0) then
534 scalar_count = data%scalar_lags%size()
535 else if (associated(data%s)) then
536 scalar_count = 1
537 else
538 scalar_count = 0
539 end if
540
541 if (scalar_count .gt. 1) then
542 fp_size = fp_size + scalar_count
543
544 ! Add abx1 and abx2 fields for each scalar
545 fp_size = fp_size + (scalar_count * 2)
546 else if (associated(data%s)) then
547 ! Single scalar support
548 fp_size = fp_size + 1
549 if (associated(data%abs1)) then
550 fp_size = fp_size + 2
551 end if
552 end if
553
554 if (associated(data%abx1)) then
555 fp_size = fp_size + 6
556 end if
557
558 allocate(fp(fp_size))
559
560 fsp_size = 0
561 if (associated(data%ulag)) then
562 fsp_size = fsp_size + 3
563 end if
564
565 if (scalar_count .gt. 1) then
566 if (allocated(data%scalar_lags%items)) then
567 fsp_size = fsp_size + data%scalar_lags%size()
568 end if
569 else if (associated(data%slag)) then
570 fsp_size = fsp_size + 1
571 end if
572
573 if (fsp_size .gt. 0) then
574 allocate(fsp(fsp_size))
575 fsp_cur = 1
576 end if
577
578 dof => data%u%dof
579 msh => data%u%msh
580
581 fp(1)%ptr => data%u
582 fp(2)%ptr => data%v
583 fp(3)%ptr => data%w
584 fp(4)%ptr => data%p
585
586 fp_cur = 5
587
588 if (scalar_count .gt. 1) then
589 block
590 type(field_series_t), pointer :: slag
591 do i = 1, scalar_count
592 slag => data%scalar_lags%get(i)
593 fp(fp_cur)%ptr => slag%f
594 fp_cur = fp_cur + 1
595 end do
596 end block
597
598 do i = 1, scalar_count
599 fp(fp_cur)%ptr => data%scalar_abx1(i)%ptr
600 fp_cur = fp_cur + 1
601 fp(fp_cur)%ptr => data%scalar_abx2(i)%ptr
602 fp_cur = fp_cur + 1
603 end do
604 else if (associated(data%s)) then
605 ! Single scalar support
606 fp(fp_cur)%ptr => data%s
607 fp_cur = fp_cur + 1
608
609 if (associated(data%abs1)) then
610 fp(fp_cur)%ptr => data%abs1
611 fp(fp_cur+1)%ptr => data%abs2
612 fp_cur = fp_cur + 2
613 end if
614 end if
615
616 if (associated(data%abx1)) then
617 fp(fp_cur)%ptr => data%abx1
618 fp(fp_cur+1)%ptr => data%abx2
619 fp(fp_cur+2)%ptr => data%aby1
620 fp(fp_cur+3)%ptr => data%aby2
621 fp(fp_cur+4)%ptr => data%abz1
622 fp(fp_cur+5)%ptr => data%abz2
623 fp_cur = fp_cur + 6
624 end if
625
626 if (associated(data%ulag)) then
627 fsp(fsp_cur)%ptr => data%ulag
628 fsp(fsp_cur+1)%ptr => data%vlag
629 fsp(fsp_cur+2)%ptr => data%wlag
630 fsp_cur = fsp_cur + 3
631 end if
632
633
634 if (scalar_count .gt. 1) then
635 if (allocated(data%scalar_lags%items)) then
636 do j = 1, data%scalar_lags%size()
637 fsp(fsp_cur)%ptr => data%scalar_lags%get(j)
638 fsp_cur = fsp_cur + 1
639 end do
640 end if
641 else if (associated(data%slag)) then
642 fsp(fsp_cur)%ptr => data%slag
643 fsp_cur = fsp_cur + 1
644 end if
645
646 if (associated(data%tlag)) then
647 tlag => data%tlag
648 dtlag => data%dtlag
649 end if
650
651 class default
652 call neko_log%error('Invalid data')
653 end select
654
655 end subroutine hdf5_file_determine_data
656
660 subroutine hdf5_file_determine_real(H5T_NEKO_REAL)
661 integer(hid_t), intent(inout) :: H5T_NEKO_REAL
662 select case (rp)
663 case (dp)
664 h5t_neko_real = h5t_native_double
665 case (sp)
666 h5t_neko_real = h5t_native_real
667 case default
668 call neko_error("Unsupported real type")
669 end select
670 end subroutine hdf5_file_determine_real
671
672 ! ================
673 ! Granular methods
674 ! ================
675
677 subroutine hdf5_file_open(this, mode)
678 class(hdf5_file_t), intent(inout) :: this
679 character(len=1), intent(in) :: mode
680 integer :: ierr, mpi_info, mpi_comm, i, n_fields, counter
681 logical :: file_exists
682 character(len=1024) :: fname
683 character(len=LOG_SIZE) :: log_buf
684
685 ! Set the mode for the file
686 this%mode = mode
687
688 ! Ensure precision is set and are valid.
689 if (this%precision .gt. rp) then
690 this%precision = rp
691 call neko_warning('Requested precision is higher than working precision')
692 else if (this%precision .eq. -1) then
693 this%precision = rp
694 end if
695
696 fname = trim(file_get_fname(this))
697 counter = this%get_counter() - this%get_start_counter()
698
699 ! Set the configuration for MPI IO
700 mpi_info = mpi_info_null%mpi_val
701 mpi_comm = neko_comm%mpi_val
702 call h5open_f(ierr)
703 call h5pcreate_f(h5p_file_access_f, this%plist_id, ierr)
704 call h5pset_fapl_mpio_f(this%plist_id, mpi_comm, mpi_info, ierr)
705
706 ! Open the file
707 inquire(file = fname, exist = file_exists)
708 if (file_exists) then
709 call h5fopen_f(fname, h5f_acc_rdwr_f, this%file_id, ierr, &
710 access_prp = this%plist_id)
711 else
712 call h5fcreate_f(fname, h5f_acc_trunc_f, &
713 this%file_id, ierr, access_prp = this%plist_id)
714 end if
715
716 ! Set the active group to the root of the file
717 call this%set_active_group()
718
719 write (log_buf, *) "Opened HDF5 file: ", trim(fname), " with counter: ", &
720 counter
721 call neko_log%message(log_buf, lvl = neko_log_debug)
722
723 end subroutine hdf5_file_open
724
726 subroutine hdf5_file_close(this)
727 class(hdf5_file_t), intent(inout) :: this
728 integer :: ierr
729
730 if (this%active_group_id .ne. -1_hid_t .and. &
731 this%active_group_id .ne. this%file_id) then
732 call h5gclose_f(this%active_group_id, ierr)
733 end if
734 this%active_group_id = -1_hid_t
735
736 call h5pclose_f(this%plist_id, ierr)
737 this%plist_id = -1_hid_t
738 call h5fclose_f(this%file_id, ierr)
739 this%file_id = -1_hid_t
740 call h5close_f(ierr)
741
742 call neko_log%message("Closed HDF5 file: " // trim(this%get_fname()), &
743 lvl = neko_log_debug)
744
745 end subroutine hdf5_file_close
746
747
751 subroutine hdf5_file_set_group(this, group_name_path)
752 class(hdf5_file_t), intent(inout) :: this
753 character(len=*), intent(in), optional :: group_name_path
754 character(len=1000), allocatable :: group_name(:)
755
756 integer(hid_t) :: current_id, group_id
757 integer :: ierr, i, j, num_groups, name_len, group_loc
758 logical :: group_exists
759
760
761 ! Close previous active group if one is open
762 if (this%active_group_id .ne. -1_hid_t .and. this%active_group_id .ne. &
763 this%file_id) then
764 call h5gclose_f(this%active_group_id, ierr)
765 end if
766 this%active_group_id = -1_hid_t
767
768 ! Start from root location = file
769 current_id = this%file_id
770 ! Return the root directory if no group name is given
771 if (.not. present(group_name_path)) then
772 this%active_group_id = current_id
773 return
774 end if
775
776 ! Split the input string into group names using "/" as a delimiter
777 name_len = len(trim(group_name_path))
778 ! Count how many groups
779 num_groups = 1 ! There is at least one group if this was passed
780 do i = 1, name_len
781 if (group_name_path .eq. "/") then
782 num_groups = num_groups + 1
783 end if
784 end do
785
786 ! Allocate the group array and populate it
787 allocate(group_name(num_groups))
788 j = 1
789 group_loc = 1
790 do i = 1, name_len
791 if (group_name_path .eq. "/") then
792 group_name(group_loc) = group_name_path(j:i-1)
793 group_loc = group_loc + 1
794 j = i + 1
795 end if
796 end do
797 if (j .ne. name_len) then
798 group_name(group_loc) = group_name_path(j:name_len)
799 end if
800
801 ! Iterate over the groups in the path
802 do i = 1, num_groups
803 call h5lexists_f(current_id, trim(group_name(i)), group_exists, ierr)
804
805 ! Only create groups if they dont exist and we are in write mode "w"
806 if (group_exists) then
807 call h5gopen_f(current_id, trim(group_name(i)), group_id, ierr)
808 else
809 if (this%mode == "r") then
810 call neko_error("Group " // trim(group_name(i)) // &
811 " does not exist in file " // trim(file_get_fname(this)))
812 end if
813 call h5gcreate_f(current_id, trim(group_name(i)), group_id, ierr)
814 end if
815
816 ! Close previous location only if it was an opened group, not the file
817 if (i > 1) then
818 call h5gclose_f(current_id, ierr)
819 end if
820
821 current_id = group_id
822 end do
823
824 this%active_group_id = current_id
825 end subroutine hdf5_file_set_group
826
827
828 subroutine hdf5_file_write_dataset(this, data)
829 class(hdf5_file_t), intent(inout) :: this
830 class(*), intent(inout) :: data
831
832 select type (d => data)
833 type is (vector_t)
834 call this%write_vector(d)
835 type is (matrix_t)
836 call this%write_matrix(d)
837 type is (field_t)
838 call this%write_field(d)
839 class default
840 call neko_error("write_dataset not implemented for this data type")
841 end select
842 end subroutine hdf5_file_write_dataset
843
844 subroutine hdf5_file_read_dataset(this, data_name, data, strategy)
845 class(hdf5_file_t), intent(inout) :: this
846 character(len=*), intent(in) :: data_name
847 class(*), intent(inout) :: data
848 character(len=*), intent(in), optional :: strategy
849
850 select type (d => data)
851 type is (vector_t)
852 call this%read_vector(data_name, d, strategy)
853 type is (matrix_t)
854 call this%read_matrix(data_name, d, strategy)
855 type is (field_t)
856 call neko_error("Reading a field_t is not supported yet")
857 class default
858 call neko_error("read_dataset not implemented for this data type")
859 end select
860 end subroutine hdf5_file_read_dataset
861
862 subroutine hdf5_file_write_attribute(this, data_name, data)
863 class(hdf5_file_t), intent(inout) :: this
864 character(len=*), intent(in) :: data_name
865 class(*), intent(inout) :: data
866
867 select type (d => data)
868 type is (integer)
869 call this%write_int_attribute(data_name, d)
870 type is (real(kind=rp))
871 call this%write_rp_attribute(data_name, d)
872 class default
873 call neko_error("write_attribute not implemented for this data type")
874 end select
875 end subroutine hdf5_file_write_attribute
876
877 subroutine hdf5_file_read_attribute(this, data_name, data, exist)
878 class(hdf5_file_t), intent(inout) :: this
879 character(len=*), intent(in) :: data_name
880 class(*), intent(inout) :: data
881 logical, intent(inout) :: exist
882
883 select type (d => data)
884 type is (integer)
885 call this%read_int_attribute(data_name, d, exist)
886 type is (real(kind=rp))
887 call this%read_rp_attribute(data_name, d, exist)
888 class default
889 call neko_error("read_attribute not implemented for this data type")
890 end select
891 end subroutine hdf5_file_read_attribute
892
893
894 subroutine hdf5_file_write_vector(this, vec)
895 class(hdf5_file_t), intent(inout) :: this
896 type(vector_t), intent(inout) :: vec
897 integer :: ierr, counts, offset, total_count, dset_rank, max_count
898 integer(hsize_t) :: append_offset
899 integer(hid_t) :: precision_hdf
900 integer(hid_t) :: xf_id, filespace, dset_id, memspace, dcpl_id
901 integer(hsize_t), dimension(1) :: dcount, doffset
902 integer(hsize_t), dimension(1) :: ddims, ddims_max, chunkdims
903 integer(hsize_t), dimension(1) :: tempddims, tempmaxddims
904 logical :: dset_exists
905 real(kind=sp), allocatable :: write_buffer_sp(:) ! Write buffer single
906 real(kind=dp), allocatable :: write_buffer_dp(:) ! Write buffer double
907
908 ! ===============
909 ! Get vector info
910 ! ===============
911 counts = vec%size()
912 append_offset = 0_hsize_t
913 offset = 0
914 total_count = 0
915 max_count = 0
916 call mpi_scan(counts, offset, 1, mpi_integer, &
917 mpi_sum, neko_comm, ierr)
918 offset = offset - counts ! Not using exclusive scan
919 call mpi_allreduce(counts, total_count, 1, mpi_integer, &
920 mpi_sum, neko_comm, ierr)
921 call mpi_allreduce(counts, max_count, 1, mpi_integer, &
922 mpi_max, neko_comm, ierr)
923
924 ! ===============
925 ! Configure MPIIO
926 ! ===============
927 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
928 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
929 precision_hdf = h5kind_to_type(this%precision, h5_real_kind)
930
931 ! ===================
932 ! Create the data set
933 ! ===================
934 dset_rank = 1 ! rank 1 array, i.e. a vector
935 ddims = [int(total_count, hsize_t)] ! global size of the vector
936
937 ! Enable chunking to be able to append
938 chunkdims = [max(int(max_count, hsize_t), 1_hsize_t)]
939 ddims_max = [h5s_unlimited_f] ! allow unlimited size for appending
940 call h5lexists_f(this%active_group_id, trim(vec%name), dset_exists, ierr)
941 if (dset_exists) then
942 if (this%overwrite) then
943 ! retrieve the dset id for the existing data set
944 call h5dopen_f(this%active_group_id, trim(vec%name), dset_id, ierr)
945 else
946 ! Retreive the existing data set
947 call h5dopen_f(this%active_group_id, trim(vec%name), dset_id, ierr)
948 ! Retrieve the current filespace (shape space)
949 call h5dget_space_f(dset_id, filespace, ierr)
950 ! Get the current shape
951 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
952 ierr)
953 ! Clean up the opened file space
954 call h5sclose_f(filespace, ierr)
955 ! Overwrite the new full shape
956 ddims(1) = ddims(1) + tempddims(1) ! New size
957 append_offset = tempddims(1) ! current size which is the offset
958 ! Extend the data set to the new shape
959 call h5dset_extent_f(dset_id, ddims, ierr)
960 end if
961 else
962 ! create file space of this shape
963 call h5screate_simple_f(dset_rank, ddims, filespace, ierr, ddims_max)
964 ! Create chunk property list (needed to be able to append)
965 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
966 call h5pset_chunk_f(dcpl_id, dset_rank, chunkdims, ierr)
967 ! create the data set with the given shape
968 call h5dcreate_f(this%active_group_id, trim(vec%name), precision_hdf, &
969 filespace, dset_id, ierr, dcpl_id = dcpl_id)
970 ! clean opened ids
971 call h5sclose_f(filespace, ierr)
972 call h5pclose_f(dcpl_id, ierr)
973 end if
974
975 ! ===========================
976 ! Set up writing the data set
977 ! ===========================
978 dcount = [int(counts, hsize_t)] ! local size of the vector
979
980 ! offset for this rank in the global vector
981 doffset = [int(offset, hsize_t) + append_offset]
982 ! Get the total file space (shape) of the data set
983 call h5dget_space_f(dset_id, filespace, ierr)
984 ! Get only the slice where my rank writes
985 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
986 ierr)
987 ! Create the corresponding memory space (buffer) for my local data
988 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
989
990
991 ! =======================
992 ! Cast and write the data
993 ! =======================
994 if (this%precision == sp) then
995 allocate(write_buffer_sp(vec%size()))
996 if (vec%size() > 0) write_buffer_sp = real(vec%x, kind=sp)
997 ! Write the data
998 call h5dwrite_f(dset_id, precision_hdf, write_buffer_sp, dcount, ierr, &
999 file_space_id = filespace, mem_space_id = memspace, &
1000 xfer_prp = xf_id)
1001 deallocate(write_buffer_sp)
1002 else if (this%precision == dp) then
1003 allocate(write_buffer_dp(vec%size()))
1004 if (vec%size() > 0) write_buffer_dp = real(vec%x, kind=dp)
1005 ! Write the data
1006 call h5dwrite_f(dset_id, precision_hdf, write_buffer_dp, dcount, ierr, &
1007 file_space_id = filespace, mem_space_id = memspace, &
1008 xfer_prp = xf_id)
1009 deallocate(write_buffer_dp)
1010 else
1011 call neko_error("Unsupported precision")
1012 end if
1013
1014 ! =======================
1015 ! Clean up
1016 ! =======================
1017 call h5pclose_f(xf_id, ierr)
1018 call h5sclose_f(memspace, ierr)
1019 call h5sclose_f(filespace, ierr)
1020 call h5dclose_f(dset_id, ierr)
1021
1022 end subroutine hdf5_file_write_vector
1023
1024 subroutine hdf5_file_write_matrix(this, mat)
1025 class(hdf5_file_t), intent(inout) :: this
1026 type(matrix_t), intent(inout) :: mat
1027 integer :: ierr, counts, offset, total_count, dset_rank, strides, max_count
1028 integer(hsize_t) :: append_offset
1029 integer(hid_t) :: precision_hdf
1030 integer(hid_t) :: xf_id, filespace, dset_id, memspace, dcpl_id
1031 integer(hsize_t), dimension(2) :: dcount, doffset
1032 integer(hsize_t), dimension(2) :: ddims, ddims_max, chunkdims
1033 integer(hsize_t), dimension(2) :: tempddims, tempmaxddims
1034 logical :: dset_exists
1035 real(kind=sp), allocatable :: write_buffer_sp(:,:) ! Write buffer single
1036 real(kind=dp), allocatable :: write_buffer_dp(:,:) ! Write buffer double
1037
1038 ! ===============
1039 ! Get Matrix info
1040 ! ===============
1041 strides = mat%get_nrows()
1042 counts = mat%get_ncols()
1043 append_offset = 0_hsize_t
1044 total_count = 0
1045 max_count = 0
1046 offset = 0
1047 call mpi_scan(counts, offset, 1, mpi_integer, &
1048 mpi_sum, neko_comm, ierr)
1049 offset = offset - counts ! Not using exclusive scan
1050 call mpi_allreduce(counts, total_count, 1, mpi_integer, &
1051 mpi_sum, neko_comm, ierr)
1052 call mpi_allreduce(counts, max_count, 1, mpi_integer, &
1053 mpi_max, neko_comm, ierr)
1054
1055 ! ===============
1056 ! Configure MPIIO
1057 ! ===============
1058 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1059 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1060 precision_hdf = h5kind_to_type(this%precision, h5_real_kind)
1061
1062 ! ===================
1063 ! Create the data set
1064 ! ===================
1065 dset_rank = 2 ! rank 2 array, i.e. a matrix
1066 ! global size of the matrix
1067 ddims = [int(strides, hsize_t), int(total_count, hsize_t)]
1068 chunkdims = [int(strides, hsize_t), max(int(max_count, hsize_t), 1_hsize_t)]
1069 ddims_max = [int(strides, hsize_t), h5s_unlimited_f]
1070 call h5lexists_f(this%active_group_id, trim(mat%name), dset_exists, ierr)
1071 if (dset_exists) then
1072 if (this%overwrite) then
1073
1074 if (pe_rank .eq. 0) then
1075 call neko_warning("Dataset " // trim(mat%name) // &
1076 " already exists and wil be overwritten")
1077 end if
1078 ! retrieve the dset id for the existing data set
1079 call h5dopen_f(this%active_group_id, trim(mat%name), dset_id, ierr)
1080 else
1081 call h5dopen_f(this%active_group_id, trim(mat%name), dset_id, ierr)
1082 call h5dget_space_f(dset_id, filespace, ierr)
1083 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
1084 ierr)
1085 call h5sclose_f(filespace, ierr)
1086 ddims(2) = ddims(2) + tempddims(2)
1087 append_offset = tempddims(2)
1088 call h5dset_extent_f(dset_id, ddims, ierr)
1089 end if
1090 else
1091 ! create file space of this shape
1092 call h5screate_simple_f(dset_rank, ddims, filespace, ierr, ddims_max)
1093 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
1094 call h5pset_chunk_f(dcpl_id, dset_rank, chunkdims, ierr)
1095 ! create the data set with the given shape
1096 call h5dcreate_f(this%active_group_id, trim(mat%name), precision_hdf, &
1097 filespace, dset_id, ierr, dcpl_id = dcpl_id)
1098 call h5sclose_f(filespace, ierr)
1099 call h5pclose_f(dcpl_id, ierr)
1100 end if
1101
1102 ! ===========================
1103 ! Set up writing the data set
1104 ! ===========================
1105 ! local size of the matrix
1106 dcount = [int(strides, hsize_t), int(counts, hsize_t)]
1107 ! offset for this rank in the global matrix
1108 doffset = [0_hsize_t, int(offset, hsize_t) + append_offset]
1109 ! Get the total file space (shape) of the data set
1110 call h5dget_space_f(dset_id, filespace, ierr)
1111 ! Get only the slice where my rank writes
1112 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
1113 ierr)
1114 ! Create the corresponding memory space (buffer) for my local data
1115 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
1116
1117 ! =======================
1118 ! Cast and write the data
1119 ! =======================
1120 if (this%precision == sp) then
1121 allocate(write_buffer_sp(mat%get_nrows(), mat%get_ncols()))
1122 if (mat%size() > 0) write_buffer_sp = real(mat%x, kind=sp)
1123 ! Write the data
1124 call h5dwrite_f(dset_id, precision_hdf, write_buffer_sp, dcount, ierr, &
1125 file_space_id = filespace, mem_space_id = memspace, &
1126 xfer_prp = xf_id)
1127 deallocate(write_buffer_sp)
1128 else if (this%precision == dp) then
1129 allocate(write_buffer_dp(mat%get_nrows(), mat%get_ncols()))
1130 if (mat%size() > 0) write_buffer_dp = real(mat%x, kind=dp)
1131 ! Write the data
1132 call h5dwrite_f(dset_id, precision_hdf, write_buffer_dp, dcount, ierr, &
1133 file_space_id = filespace, mem_space_id = memspace, &
1134 xfer_prp = xf_id)
1135 deallocate(write_buffer_dp)
1136 else
1137 call neko_error("Unsupported precision")
1138 end if
1139
1140 ! =======================
1141 ! Clean up
1142 ! =======================
1143 call h5pclose_f(xf_id, ierr)
1144 call h5sclose_f(memspace, ierr)
1145 call h5sclose_f(filespace, ierr)
1146 call h5dclose_f(dset_id, ierr)
1147
1148 end subroutine hdf5_file_write_matrix
1149
1150 subroutine hdf5_file_write_field(this, field)
1151 class(hdf5_file_t), intent(inout) :: this
1152 type(field_t), intent(inout) :: field
1153 integer :: ierr, counts, offset, total_count, dset_rank, max_count
1154 integer :: stride_ax_1, stride_ax_2, stride_ax_3
1155 integer(hsize_t) :: append_offset
1156 integer(hid_t) :: precision_hdf
1157 integer(hid_t) :: xf_id, filespace, dset_id, memspace, dcpl_id
1158 integer(hsize_t), dimension(4) :: dcount, doffset
1159 integer(hsize_t), dimension(4) :: ddims, ddims_max, chunkdims
1160 integer(hsize_t), dimension(4) :: tempddims, tempmaxddims
1161 logical :: dset_exists
1162 real(kind=sp), allocatable :: write_buffer_sp(:,:,:,:) ! Write buffer single
1163 real(kind=dp), allocatable :: write_buffer_dp(:,:,:,:) ! Write buffer double
1164
1165 ! ==============
1166 ! Get Field info
1167 ! ==============
1168 stride_ax_1 = field%Xh%lx
1169 stride_ax_2 = field%Xh%ly
1170 stride_ax_3 = field%Xh%lz
1171 counts = field%msh%nelv
1172 append_offset = 0_hsize_t
1173 total_count = field%msh%glb_nelv
1174 max_count = 0
1175 offset = field%msh%offset_el
1176 call mpi_allreduce(counts, max_count, 1, mpi_integer, &
1177 mpi_max, neko_comm, ierr)
1178
1179 ! ===============
1180 ! Configure MPIIO
1181 ! ===============
1182 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1183 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1184 precision_hdf = h5kind_to_type(this%precision, h5_real_kind)
1185
1186 ! ===================
1187 ! Create the data set
1188 ! ===================
1189 dset_rank = 4 ! rank 4 array, i.e. a 4D tensor
1190 ddims = [int(stride_ax_1, hsize_t), &
1191 int(stride_ax_2, hsize_t), &
1192 int(stride_ax_3, hsize_t), &
1193 int(total_count, hsize_t)] ! global size of the tensor
1194 chunkdims = [int(stride_ax_1, hsize_t), &
1195 int(stride_ax_2, hsize_t), &
1196 int(stride_ax_3, hsize_t), &
1197 max(int(max_count, hsize_t), 1_hsize_t)]
1198 ddims_max = [int(stride_ax_1, hsize_t), &
1199 int(stride_ax_2, hsize_t), &
1200 int(stride_ax_3, hsize_t), &
1201 h5s_unlimited_f]
1202 call h5lexists_f(this%active_group_id, trim(field%name), dset_exists, ierr)
1203 if (dset_exists) then
1204 if (this%overwrite) then
1205 ! retrieve the dset id for the existing data set
1206 if (pe_rank .eq. 0) then
1207 call neko_warning("Overwriting dataset: " // trim(field%name))
1208 end if
1209 call h5dopen_f(this%active_group_id, trim(field%name), dset_id, ierr)
1210 else
1211 call h5dopen_f(this%active_group_id, trim(field%name), dset_id, ierr)
1212 call h5dget_space_f(dset_id, filespace, ierr)
1213 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
1214 ierr)
1215 call h5sclose_f(filespace, ierr)
1216 ddims(4) = ddims(4) + tempddims(4)
1217 append_offset = tempddims(4)
1218 call h5dset_extent_f(dset_id, ddims, ierr)
1219 end if
1220 else
1221 ! create file space of this shape
1222 call h5screate_simple_f(dset_rank, ddims, filespace, ierr, ddims_max)
1223 call h5pcreate_f(h5p_dataset_create_f, dcpl_id, ierr)
1224 call h5pset_chunk_f(dcpl_id, dset_rank, chunkdims, ierr)
1225 ! create the data set with the given shape
1226 call h5dcreate_f(this%active_group_id, trim(field%name), precision_hdf, &
1227 filespace, dset_id, ierr, dcpl_id = dcpl_id)
1228 call h5sclose_f(filespace, ierr)
1229 call h5pclose_f(dcpl_id, ierr)
1230 end if
1231
1232 ! ===========================
1233 ! Set up writing the data set
1234 ! ===========================
1235 dcount = [int(stride_ax_1, hsize_t), &
1236 int(stride_ax_2, hsize_t), &
1237 int(stride_ax_3, hsize_t), &
1238 int(counts, hsize_t)] ! local size of the tensor
1239 doffset = [0_hsize_t, 0_hsize_t, 0_hsize_t, &
1240 int(offset, hsize_t) + append_offset] ! offset in the global tensor
1241 ! Get the total file space (shape) of the data set
1242 call h5dget_space_f(dset_id, filespace, ierr)
1243 ! Get only the slice where my rank writes
1244 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
1245 ierr)
1246 ! Create the corresponding memory space (buffer) for my local data
1247 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
1248
1249 ! =======================
1250 ! Cast and write the data
1251 ! =======================
1252 if (this%precision == sp) then
1253 allocate(write_buffer_sp(field%Xh%lx, field%Xh%ly, field%Xh%lz, &
1254 field%msh%nelv))
1255 if (field%msh%nelv > 0) write_buffer_sp = real(field%x, kind=sp)
1256 ! Write the data
1257 call h5dwrite_f(dset_id, precision_hdf, write_buffer_sp, dcount, ierr, &
1258 file_space_id = filespace, mem_space_id = memspace, &
1259 xfer_prp = xf_id)
1260 deallocate(write_buffer_sp)
1261 else if (this%precision == dp) then
1262 allocate(write_buffer_dp(field%Xh%lx, field%Xh%ly, field%Xh%lz, &
1263 field%msh%nelv))
1264 if (field%msh%nelv > 0) write_buffer_dp = real(field%x, kind=dp)
1265 ! Write the data
1266 call h5dwrite_f(dset_id, precision_hdf, write_buffer_dp, dcount, ierr, &
1267 file_space_id = filespace, mem_space_id = memspace, &
1268 xfer_prp = xf_id)
1269 deallocate(write_buffer_dp)
1270 else
1271 call neko_error("Unsupported precision")
1272 end if
1273
1274 ! =======================
1275 ! Clean up
1276 ! =======================
1277 call h5pclose_f(xf_id, ierr)
1278 call h5sclose_f(memspace, ierr)
1279 call h5sclose_f(filespace, ierr)
1280 call h5dclose_f(dset_id, ierr)
1281
1282 end subroutine hdf5_file_write_field
1283
1285 subroutine hdf5_file_read_vector(this, data_name, vec, strategy)
1286 class(hdf5_file_t) :: this
1287 character(len=*), intent(in) :: data_name
1288 type(vector_t), intent(inout) :: vec
1289 character(len=*), intent(in), optional :: strategy
1290 character(len=1000) :: strategy_
1291 integer :: ierr, counts, offset, total_count, dset_rank
1292 integer(hid_t) :: precision_hdf
1293 integer(hid_t) :: xf_id, filespace, dset_id, memspace
1294 integer(hsize_t), dimension(1) :: dcount, doffset
1295 integer(hsize_t), dimension(1) :: tempddims, tempmaxddims
1296 integer :: temprank
1297 logical :: dset_exists
1298 type(linear_dist_t) :: dist
1299
1300 ! Set up strategy
1301 if (present(strategy)) then
1302 if (trim(strategy) .eq. "linear" .or. &
1303 trim(strategy) .eq. "rank_0") then
1304 strategy_ = strategy
1305 else
1306 call neko_error("Unsupported strategy: " // trim(strategy))
1307 end if
1308 else
1309 strategy_ = "linear"
1310 end if
1311
1312 ! Free the input
1313 call vec%free()
1314
1315 ! ===============
1316 ! Configure MPIIO
1317 ! ===============
1318 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1319 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1320 precision_hdf = h5kind_to_type(rp, h5_real_kind)
1321
1322 ! ===================
1323 ! Get the data set info
1324 ! ===================
1325 call h5lexists_f(this%active_group_id, trim(data_name), dset_exists, ierr)
1326 if (dset_exists) then
1327 ! Open the data set
1328 call h5dopen_f(this%active_group_id, trim(data_name), dset_id, ierr)
1329 ! Get the current rank of the dataset
1330 call h5dget_space_f(dset_id, filespace, ierr)
1331 call h5sget_simple_extent_ndims_f(filespace, temprank, ierr)
1332 if (temprank .ne. 1) then
1333 call neko_error("Dataset " // trim(data_name) // &
1334 " is not a rank 1 vector in file " // trim(file_get_fname(this)))
1335 end if
1336 ! Get the current shape and close the filespace
1337 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
1338 ierr)
1339 call h5sclose_f(filespace, ierr)
1340 else
1341 call neko_error("Dataset " // trim(data_name) // &
1342 " does not exist in current group " // trim(file_get_fname(this)))
1343 end if
1344
1345 ! =============================
1346 ! Perform the data distribution
1347 ! =============================
1348 total_count = int(tempddims(1))
1349 if (strategy_ .eq. "linear") then
1350 dist = linear_dist_t(total_count, pe_rank, pe_size, neko_comm)
1351 counts = dist%num_local()
1352 offset = 0
1353 else if (strategy_ .eq. "rank_0") then
1354 if (pe_rank .eq. 0) then
1355 counts = total_count
1356 else
1357 counts = 0
1358 end if
1359 offset = 0
1360 end if
1361 call mpi_exscan(counts, offset, 1, mpi_integer, &
1362 mpi_sum, neko_comm, ierr)
1363
1364 ! ===========================
1365 ! Set up reading the data set
1366 ! ===========================
1367 dset_rank = 1 ! rank 1 array, i.e. a vector
1368 dcount = [int(counts, hsize_t)] ! local size of the vector
1369 doffset = [int(offset, hsize_t)] ! offset for this rank in the global vector
1370 ! Get the total file space (shape) of the data set
1371 call h5dget_space_f(dset_id, filespace, ierr)
1372 ! Get only the slice where my rank reads
1373 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
1374 ierr)
1375 ! Create the corresponding memory space (buffer) for my local data
1376 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
1377
1378 ! =============================
1379 ! Allocate data. HDF5 will cast
1380 ! =============================
1381 call vec%init(counts, trim(data_name)) ! this is rp
1382 call h5dread_f(dset_id, precision_hdf, vec%x, dcount, ierr, &
1383 file_space_id = filespace, mem_space_id = memspace, &
1384 xfer_prp = xf_id)
1385
1386 ! =======================
1387 ! Clean up
1388 ! =======================
1389 call h5pclose_f(xf_id, ierr)
1390 call h5sclose_f(memspace, ierr)
1391 call h5sclose_f(filespace, ierr)
1392 call h5dclose_f(dset_id, ierr)
1393
1394 end subroutine hdf5_file_read_vector
1395
1397 subroutine hdf5_file_read_matrix(this, data_name, mat, strategy)
1398 class(hdf5_file_t) :: this
1399 character(len=*), intent(in) :: data_name
1400 type(matrix_t), intent(inout) :: mat
1401 character(len=*), intent(in), optional :: strategy
1402 character(len=1000) :: strategy_
1403 integer :: ierr, counts, offset, total_count, dset_rank
1404 integer(hid_t) :: precision_hdf
1405 integer(hid_t) :: xf_id, filespace, dset_id, memspace
1406 integer(hsize_t), dimension(2) :: dcount, doffset
1407 integer(hsize_t), dimension(2) :: tempddims, tempmaxddims
1408 integer :: temprank
1409 logical :: dset_exists
1410 type(linear_dist_t) :: dist
1411
1412 ! Set up strategy
1413 if (present(strategy)) then
1414 if (trim(strategy) .eq. "linear" .or. &
1415 trim(strategy) .eq. "rank_0") then
1416 strategy_ = strategy
1417 else
1418 call neko_error("Unsupported strategy: " // trim(strategy))
1419 end if
1420 else
1421 strategy_ = "linear"
1422 end if
1423
1424 ! Free the input
1425 call mat%free()
1426
1427 ! ===============
1428 ! Configure MPIIO
1429 ! ===============
1430 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
1431 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
1432 precision_hdf = h5kind_to_type(rp, h5_real_kind)
1433
1434 ! ===================
1435 ! Get the data set info
1436 ! ===================
1437 call h5lexists_f(this%active_group_id, trim(data_name), dset_exists, ierr)
1438 if (dset_exists) then
1439 ! Openr the data set
1440 call h5dopen_f(this%active_group_id, trim(data_name), dset_id, ierr)
1441 ! Get the current rank of the of the dataset
1442 call h5dget_space_f(dset_id, filespace, ierr)
1443 call h5sget_simple_extent_ndims_f(filespace, temprank, ierr)
1444 if (temprank .ne. 2) then
1445 call neko_error("Dataset " // trim(data_name) // &
1446 " is not a rank 2 matrix in file " // trim(file_get_fname(this)))
1447 end if
1448 ! Get the current shape and close the filespace
1449 call h5sget_simple_extent_dims_f(filespace, tempddims, tempmaxddims, &
1450 ierr)
1451 call h5sclose_f(filespace, ierr)
1452 else
1453 call neko_error("Dataset " // trim(data_name) &
1454 // " does not exist in current group " // &
1455 trim(file_get_fname(this)))
1456 end if
1457
1458 ! =============================
1459 ! Perform the data distribution
1460 ! =============================
1461 total_count = int(tempddims(2))
1462 if (strategy_ .eq. "linear") then
1463 dist = linear_dist_t(total_count, pe_rank, pe_size, neko_comm)
1464 counts = dist%num_local()
1465 offset = 0
1466 else if (strategy_ .eq. "rank_0") then
1467 if (pe_rank .eq. 0) then
1468 counts = total_count
1469 else
1470 counts = 0
1471 end if
1472 offset = 0
1473 end if
1474 call mpi_scan(counts, offset, 1, mpi_integer, &
1475 mpi_sum, neko_comm, ierr)
1476 offset = offset - counts ! Not using exclusive scan
1477
1478 ! ===========================
1479 ! Set up reading the data set
1480 ! ===========================
1481 dset_rank = 2 ! rank 2 array, i.e. a matrix
1482 dcount = [int(tempddims(1), hsize_t), int(counts, hsize_t)] ! local size
1483 doffset = [0_hsize_t, int(offset, hsize_t)] ! offset in the global matrix
1484 ! Get the total file space (shape) of the data set
1485 call h5dget_space_f(dset_id, filespace, ierr)
1486 ! Get only the slice where my rank reads
1487 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
1488 ierr)
1489 ! Create the corresponding memory space (buffer) for my local data
1490 call h5screate_simple_f(dset_rank, dcount, memspace, ierr)
1491
1492 ! =============================
1493 ! Allocate data. HDF5 will cast
1494 ! =============================
1495 call mat%init(int(tempddims(1)), counts, trim(data_name)) ! this is rp
1496 call h5dread_f(dset_id, precision_hdf, mat%x, dcount, ierr, &
1497 file_space_id = filespace, mem_space_id = memspace, &
1498 xfer_prp = xf_id)
1499
1500 ! =======================
1501 ! Clean up
1502 ! =======================
1503 call h5pclose_f(xf_id, ierr)
1504 call h5sclose_f(memspace, ierr)
1505 call h5sclose_f(filespace, ierr)
1506 call h5dclose_f(dset_id, ierr)
1507
1508
1509 end subroutine hdf5_file_read_matrix
1510
1512 subroutine hdf5_file_write_int_attribute(this, attr_name, attr)
1513 class(hdf5_file_t), intent(inout) :: this
1514 character(len=*), intent(in) :: attr_name
1515 integer, intent(in) :: attr
1516 integer :: ierr
1517 integer(hid_t) :: filespace, attr_id
1518 integer(hsize_t), dimension(1) :: dcount
1519 logical :: attr_exists
1520
1521 ! ====================
1522 ! Create the attribute
1523 ! ====================
1524 dcount = [int(1, hsize_t)]
1525 call h5aexists_f(this%active_group_id, trim(attr_name), attr_exists, ierr)
1526 if (attr_exists) then
1527 ! retrieve the attr id for the existing attribute
1528 call h5aopen_f(this%active_group_id, trim(attr_name), attr_id, ierr)
1529 else
1530 ! create file space of this shape
1531 call h5screate_f(h5s_scalar_f, filespace, ierr)
1532 ! create the data set with the given shape
1533 call h5acreate_f(this%active_group_id, trim(attr_name), &
1534 h5t_native_integer, &
1535 filespace, attr_id, ierr, h5p_default_f, h5p_default_f)
1536 call h5sclose_f(filespace, ierr)
1537 end if
1538
1539 ! ===========================
1540 ! Set up writing the data set
1541 ! ===========================
1542 call h5awrite_f(attr_id, h5t_native_integer, attr, dcount, ierr)
1543
1544 ! =======================
1545 ! Clean up
1546 ! =======================
1547 call h5aclose_f(attr_id, ierr)
1548
1549 end subroutine hdf5_file_write_int_attribute
1550
1552 subroutine hdf5_file_write_rp_attribute(this, attr_name, attr)
1553 class(hdf5_file_t), intent(inout) :: this
1554 character(len=*), intent(in) :: attr_name
1555 real(kind=rp), intent(in) :: attr
1556 integer :: ierr
1557 integer(hid_t) :: precision_hdf
1558 integer(hid_t) :: filespace, attr_id
1559 integer(hsize_t), dimension(1) :: dcount
1560 logical :: attr_exists
1561
1562 ! Get the precision
1563 precision_hdf = h5kind_to_type(rp, h5_real_kind)
1564
1565 ! ====================
1566 ! Create the attribute
1567 ! ====================
1568 dcount = [int(1, hsize_t)]
1569 call h5aexists_f(this%active_group_id, trim(attr_name), attr_exists, ierr)
1570 if (attr_exists) then
1571 ! retrieve the attr id for the existing attribute
1572 call h5aopen_f(this%active_group_id, trim(attr_name), attr_id, ierr)
1573 else
1574 ! create file space of this shape
1575 call h5screate_f(h5s_scalar_f, filespace, ierr)
1576 ! create the data set with the given shape
1577 call h5acreate_f(this%active_group_id, trim(attr_name), precision_hdf, &
1578 filespace, attr_id, ierr, h5p_default_f, h5p_default_f)
1579 call h5sclose_f(filespace, ierr)
1580 end if
1581
1582 ! ===========================
1583 ! Set up writing the data set
1584 ! ===========================
1585 call h5awrite_f(attr_id, precision_hdf, attr, dcount, ierr)
1586
1587 ! =======================
1588 ! Clean up
1589 ! =======================
1590 call h5aclose_f(attr_id, ierr)
1591
1592 end subroutine hdf5_file_write_rp_attribute
1593
1595 subroutine hdf5_file_read_int_attribute(this, attr_name, attr, attr_exists)
1596 class(hdf5_file_t), intent(inout) :: this
1597 character(len=*), intent(in) :: attr_name
1598 integer, intent(inout) :: attr
1599 logical, intent(inout) :: attr_exists
1600 integer :: ierr
1601 integer(hid_t) :: filespace, attr_id
1602 integer(hsize_t), dimension(1) :: dcount
1603
1604 ! ====================
1605 ! Create the attribute
1606 ! ====================
1607 dcount = [int(1, hsize_t)]
1608 call h5aexists_f(this%active_group_id, trim(attr_name), attr_exists, ierr)
1609 if (attr_exists) then
1610 ! retrieve the attr id for the existing attribute
1611 call h5aopen_f(this%active_group_id, trim(attr_name), attr_id, ierr)
1612 else
1613 return
1614 end if
1615
1616 ! ===========================
1617 ! Set up writing the data set
1618 ! ===========================
1619 call h5aread_f(attr_id, h5t_native_integer, attr, dcount, ierr)
1620
1621 ! =======================
1622 ! Clean up
1623 ! =======================
1624 call h5aclose_f(attr_id, ierr)
1625
1626 end subroutine hdf5_file_read_int_attribute
1627
1629 subroutine hdf5_file_read_rp_attribute(this, attr_name, attr, attr_exists)
1630 class(hdf5_file_t), intent(inout) :: this
1631 character(len=*), intent(in) :: attr_name
1632 real(kind=rp), intent(inout) :: attr
1633 logical, intent(inout) :: attr_exists
1634 integer :: ierr
1635 integer(hid_t) :: precision_hdf
1636 integer(hid_t) :: filespace, attr_id
1637 integer(hsize_t), dimension(1) :: dcount
1638
1639 ! Get the precision
1640 precision_hdf = h5kind_to_type(rp, h5_real_kind)
1641
1642 ! ====================
1643 ! Create the attribute
1644 ! ====================
1645 dcount = [int(1, hsize_t)]
1646 call h5aexists_f(this%active_group_id, trim(attr_name), attr_exists, ierr)
1647 if (attr_exists) then
1648 ! retrieve the attr id for the existing attribute
1649 call h5aopen_f(this%active_group_id, trim(attr_name), attr_id, ierr)
1650 else
1651 return
1652 end if
1653
1654 ! ===========================
1655 ! Set up writing the data set
1656 ! ===========================
1657 call h5aread_f(attr_id, precision_hdf, attr, dcount, ierr)
1658
1659 ! =======================
1660 ! Clean up
1661 ! =======================
1662 call h5aclose_f(attr_id, ierr)
1663
1664 end subroutine hdf5_file_read_rp_attribute
1665
1666#else
1667
1669 subroutine hdf5_file_open(this, mode)
1670 class(hdf5_file_t), intent(inout) :: this
1671 character(len=1), intent(in) :: mode
1672 call neko_error('Neko needs to be built with HDF5 support')
1673 end subroutine hdf5_file_open
1674
1676 subroutine hdf5_file_close(this)
1677 class(hdf5_file_t), intent(inout) :: this
1678 call neko_error('Neko needs to be built with HDF5 support')
1679 end subroutine hdf5_file_close
1680
1682 subroutine hdf5_file_set_group(this, group_name_path)
1683 class(hdf5_file_t), intent(inout) :: this
1684 character(len=*), intent(in), optional :: group_name_path
1685 call neko_error('Neko needs to be built with HDF5 support')
1686 end subroutine hdf5_file_set_group
1687
1689 subroutine hdf5_file_write(this, data, t)
1690 class(hdf5_file_t), intent(inout) :: this
1691 class(*), target, intent(in) :: data
1692 real(kind=rp), intent(in), optional :: t
1693 call neko_error('Neko needs to be built with HDF5 support')
1694 end subroutine hdf5_file_write
1695
1697 subroutine hdf5_file_read(this, data)
1698 class(hdf5_file_t) :: this
1699 class(*), target, intent(inout) :: data
1700 call neko_error('Neko needs to be built with HDF5 support')
1701 end subroutine hdf5_file_read
1702
1703 subroutine hdf5_file_write_dataset(this, data)
1704 class(hdf5_file_t), intent(inout) :: this
1705 class(*), intent(inout) :: data
1706 call neko_error('Neko needs to be built with HDF5 support')
1707 end subroutine hdf5_file_write_dataset
1708
1709 subroutine hdf5_file_read_dataset(this, data_name, data, strategy)
1710 class(hdf5_file_t), intent(inout) :: this
1711 character(len=*), intent(in) :: data_name
1712 class(*), intent(inout) :: data
1713 character(len=*), intent(in), optional :: strategy
1714 call neko_error('Neko needs to be built with HDF5 support')
1715 end subroutine hdf5_file_read_dataset
1716
1717 subroutine hdf5_file_write_attribute(this, data_name, data)
1718 class(hdf5_file_t), intent(inout) :: this
1719 character(len=*), intent(in) :: data_name
1720 class(*), intent(inout) :: data
1721 call neko_error('Neko needs to be built with HDF5 support')
1722 end subroutine hdf5_file_write_attribute
1723
1724 subroutine hdf5_file_read_attribute(this, data_name, data, exist)
1725 class(hdf5_file_t), intent(inout) :: this
1726 character(len=*), intent(in) :: data_name
1727 class(*), intent(inout) :: data
1728 logical, intent(inout) :: exist
1729 call neko_error('Neko needs to be built with HDF5 support')
1730 end subroutine hdf5_file_read_attribute
1731
1732 subroutine hdf5_file_write_vector(this, vec)
1733 class(hdf5_file_t), intent(inout) :: this
1734 type(vector_t), intent(inout) :: vec
1735 call neko_error('Neko needs to be built with HDF5 support')
1736 end subroutine hdf5_file_write_vector
1737
1738 subroutine hdf5_file_write_matrix(this, mat)
1739 class(hdf5_file_t), intent(inout) :: this
1740 type(matrix_t), intent(inout) :: mat
1741 call neko_error('Neko needs to be built with HDF5 support')
1742 end subroutine hdf5_file_write_matrix
1743
1744 subroutine hdf5_file_write_field(this, fld)
1745 class(hdf5_file_t), intent(inout) :: this
1746 type(field_t), intent(inout) :: fld
1747 call neko_error('Neko needs to be built with HDF5 support')
1748 end subroutine hdf5_file_write_field
1749
1750 subroutine hdf5_file_read_vector(this, data_name, vec, strategy)
1751 class(hdf5_file_t) :: this
1752 character(len=*), intent(in) :: data_name
1753 type(vector_t), intent(inout) :: vec
1754 character(len=*), intent(in), optional :: strategy
1755 call neko_error('Neko needs to be built with HDF5 support')
1756 end subroutine hdf5_file_read_vector
1757
1758 subroutine hdf5_file_read_matrix(this, data_name, mat, strategy)
1759 class(hdf5_file_t) :: this
1760 character(len=*), intent(in) :: data_name
1761 type(matrix_t), intent(inout) :: mat
1762 character(len=*), intent(in), optional :: strategy
1763 call neko_error('Neko needs to be built with HDF5 support')
1764 end subroutine hdf5_file_read_matrix
1765
1766 subroutine hdf5_file_write_int_attribute(this, attr_name, attr)
1767 class(hdf5_file_t), intent(inout) :: this
1768 character(len=*), intent(in) :: attr_name
1769 integer, intent(in) :: attr
1770 call neko_error('Neko needs to be built with HDF5 support')
1771 end subroutine hdf5_file_write_int_attribute
1772
1773 subroutine hdf5_file_write_rp_attribute(this, attr_name, attr)
1774 class(hdf5_file_t), intent(inout) :: this
1775 character(len=*), intent(in) :: attr_name
1776 real(kind=rp), intent(in) :: attr
1777 call neko_error('Neko needs to be built with HDF5 support')
1778 end subroutine hdf5_file_write_rp_attribute
1779
1780 subroutine hdf5_file_read_int_attribute(this, attr_name, attr, attr_exists)
1781 class(hdf5_file_t), intent(inout) :: this
1782 character(len=*), intent(in) :: attr_name
1783 integer, intent(inout) :: attr
1784 logical, intent(inout) :: attr_exists
1785 call neko_error('Neko needs to be built with HDF5 support')
1786 end subroutine hdf5_file_read_int_attribute
1787
1788 subroutine hdf5_file_read_rp_attribute(this, attr_name, attr, attr_exists)
1789 class(hdf5_file_t), intent(inout) :: this
1790 character(len=*), intent(in) :: attr_name
1791 real(kind=rp), intent(inout) :: attr
1792 logical, intent(inout) :: attr_exists
1793 call neko_error('Neko needs to be built with HDF5 support')
1794 end subroutine hdf5_file_read_rp_attribute
1795
1796#endif
1797
1798end module hdf5_file
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
Defines practical data distributions.
Definition datadist.f90:34
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
HDF5 file format.
Definition hdf5_file.F90:34
subroutine hdf5_file_open(this, mode)
Open a HDF5 file in a given mode.
character(len=1024) function file_get_fname(this)
Return the file name with the start counter.
subroutine hdf5_file_write_field(this, fld)
subroutine hdf5_file_read(this, data)
Read data in HDF5 format.
subroutine hdf5_file_write(this, data, t)
Write data in HDF5 format.
subroutine hdf5_file_read_dataset(this, data_name, data, strategy)
subroutine hdf5_file_write_int_attribute(this, attr_name, attr)
subroutine hdf5_file_write_attribute(this, data_name, data)
subroutine hdf5_file_read_vector(this, data_name, vec, strategy)
subroutine hdf5_file_write_vector(this, vec)
subroutine hdf5_file_read_matrix(this, data_name, mat, strategy)
subroutine hdf5_file_read_int_attribute(this, attr_name, attr, attr_exists)
subroutine hdf5_file_write_rp_attribute(this, attr_name, attr)
subroutine hdf5_file_set_precision(this, precision)
Set the precision for the output (single or double)
subroutine hdf5_file_read_rp_attribute(this, attr_name, attr, attr_exists)
subroutine hdf5_file_read_attribute(this, data_name, data, exist)
subroutine hdf5_file_close(this)
Close the file.
subroutine hdf5_file_set_group(this, group_name_path)
Set the active group for HDF5 files.
subroutine hdf5_file_write_matrix(this, mat)
subroutine hdf5_file_write_dataset(this, data)
subroutine hdf5_file_set_overwrite(this, overwrite)
Set the overwrite flag for HDF5 files.
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_debug
Debug log level.
Definition log.f90:89
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Defines a matrix.
Definition matrix.f90:34
Defines a mesh.
Definition mesh.f90:34
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
subroutine, public filename_split(fname, path, name, suffix)
Extract file name components.
Definition utils.f90:125
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:392
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition utils.f90:68
Defines a vector.
Definition vector.f90:34
Load-balanced linear distribution .
Definition datadist.f90:50
field_ptr_t, To easily obtain a pointer to a field
Definition field.f90:83
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.
Definition hdf5_file.F90:60
#define max(a, b)
Definition tensor.cu:40