Neko 1.99.1
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
39 use mesh, only : mesh_t
40 use field, only : field_t, field_ptr_t
41 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, neko_comm
47 use mpi_f08, only : mpi_info_null, mpi_allreduce, mpi_in_place, &
48 mpi_integer8, mpi_sum
49#ifdef HAVE_HDF5
50 use hdf5
51#endif
52 implicit none
53 private
54
56 type, public, extends(generic_file_t) :: hdf5_file_t
57 contains
58 procedure :: read => hdf5_file_read
59 procedure :: write => hdf5_file_write
60 procedure :: set_overwrite => hdf5_file_set_overwrite
61 end type hdf5_file_t
62
63contains
64
65#ifdef HAVE_HDF5
66
68 subroutine hdf5_file_write(this, data, t)
69 class(hdf5_file_t), intent(inout) :: this
70 class(*), target, intent(in) :: data
71 real(kind=rp), intent(in), optional :: t
72 type(mesh_t), pointer :: msh
73 type(dofmap_t), pointer :: dof
74 type(field_ptr_t), allocatable :: fp(:)
75 type(field_series_ptr_t), allocatable :: fsp(:)
76 real(kind=rp), pointer :: dtlag(:)
77 real(kind=rp), pointer :: tlag(:)
78 integer :: ierr, info, drank, i, j
79 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
80 integer(hid_t) :: filespace, memspace
81 integer(hid_t) :: h5t_neko_real
82 integer(hsize_t), dimension(1) :: ddim, dcount, doffset
83 integer :: suffix_pos
84 character(len=5) :: id_str
85 character(len=1024) :: fname
86
87 call hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
88
89 if (.not. this%overwrite) call this%increment_counter()
90 fname = trim(this%get_fname())
91
92 call h5open_f(ierr)
93
94 call hdf5_file_determine_real(h5t_neko_real)
95
96 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
97 info = mpi_info_null%mpi_val
98 call h5pset_fapl_mpio_f(plist_id, neko_comm%mpi_val, info, ierr)
99
100 call h5fcreate_f(fname, h5f_acc_trunc_f, &
101 file_id, ierr, access_prp = plist_id)
102
103 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
104 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
105
106 call h5screate_f(h5s_scalar_f, filespace, ierr)
107 ddim = 1
108
109 if (present(t)) then
110 call h5acreate_f(file_id, "Time", h5t_neko_real, filespace, attr_id, &
111 ierr, h5p_default_f, h5p_default_f)
112 call h5awrite_f(attr_id, h5t_neko_real, t, ddim, ierr)
113 call h5aclose_f(attr_id, ierr)
114 end if
115
116 if (associated(dof)) then
117 call h5acreate_f(file_id, "Lx", h5t_native_integer, filespace, attr_id, &
118 ierr, h5p_default_f, h5p_default_f)
119 call h5awrite_f(attr_id, h5t_native_integer, dof%Xh%lx, ddim, ierr)
120 call h5aclose_f(attr_id, ierr)
121 end if
122
123 if (associated(msh)) then
124 call h5gcreate_f(file_id, "Mesh", grp_id, ierr, lcpl_id=h5p_default_f, &
125 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
126
127 call h5acreate_f(grp_id, "Elements", h5t_native_integer, filespace, attr_id, &
128 ierr, h5p_default_f, h5p_default_f)
129 call h5awrite_f(attr_id, h5t_native_integer, msh%glb_nelv, ddim, ierr)
130 call h5aclose_f(attr_id, ierr)
131
132 call h5acreate_f(grp_id, "Dimension", h5t_native_integer, filespace, attr_id, &
133 ierr, h5p_default_f, h5p_default_f)
134 call h5awrite_f(attr_id, h5t_native_integer, msh%gdim, ddim, ierr)
135 call h5aclose_f(attr_id, ierr)
136
137 call h5gclose_f(grp_id, ierr)
138 end if
139
140
141 call h5sclose_f(filespace, ierr)
142
143 !
144 ! Write restart group (tlag, dtlag)
145 !
146 if (associated(tlag) .and. associated(dtlag)) then
147 call h5gcreate_f(file_id, "Restart", grp_id, ierr, lcpl_id=h5p_default_f, &
148 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
149
150 drank = 1
151 ddim = size(tlag)
152 doffset(1) = 0
153 if (pe_rank .eq. 0) then
154 dcount = size(tlag)
155 else
156 dcount = 0
157 end if
158
159 call h5screate_simple_f(drank, ddim, filespace, ierr)
160
161 call h5dcreate_f(grp_id,'tlag', h5t_neko_real, &
162 filespace, dset_id, ierr)
163 call h5dget_space_f(dset_id, filespace, ierr)
164 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
165 doffset, dcount, ierr)
166 call h5dwrite_f(dset_id, h5t_neko_real, tlag, &
167 ddim, ierr, xfer_prp = plist_id)
168 call h5dclose_f(dset_id, ierr)
169
170 call h5dcreate_f(grp_id,'dtlag', h5t_neko_real, &
171 filespace, dset_id, ierr)
172 call h5dget_space_f(dset_id, filespace, ierr)
173 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
174 doffset, dcount, ierr)
175 call h5dwrite_f(dset_id, h5t_neko_real, dtlag, &
176 ddim, ierr, xfer_prp = plist_id)
177 call h5dclose_f(dset_id, ierr)
178
179 call h5sclose_f(filespace, ierr)
180 call h5gclose_f(grp_id, ierr)
181
182 end if
183
184
185 !
186 ! Write fields group
187 !
188 if (allocated(fp) .or. allocated(fsp)) then
189 call h5gcreate_f(file_id, "Fields", grp_id, ierr, lcpl_id=h5p_default_f, &
190 gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
191
192 dcount(1) = int(dof%size(), 8)
193 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
194 ddim = int(dof%size(), 8)
195 drank = 1
196 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
197 mpi_integer8, mpi_sum, neko_comm, ierr)
198
199 call h5screate_simple_f(drank, ddim, filespace, ierr)
200 call h5screate_simple_f(drank, dcount, memspace, ierr)
201
202
203 if (allocated(fp)) then
204 do i = 1, size(fp)
205 call h5dcreate_f(grp_id, fp(i)%ptr%name, h5t_neko_real, &
206 filespace, dset_id, ierr)
207 call h5dget_space_f(dset_id, filespace, ierr)
208 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
209 doffset, dcount, ierr)
210 call h5dwrite_f(dset_id, h5t_neko_real, &
211 fp(i)%ptr%x(1,1,1,1), &
212 ddim, ierr, file_space_id = filespace, &
213 mem_space_id = memspace, xfer_prp = plist_id)
214 call h5dclose_f(dset_id, ierr)
215 end do
216 deallocate(fp)
217 end if
218
219 if (allocated(fsp)) then
220 do i = 1, size(fsp)
221 do j = 1, fsp(i)%ptr%size()
222 call h5dcreate_f(grp_id, fsp(i)%ptr%lf(j)%name, &
223 h5t_neko_real, filespace, dset_id, ierr)
224 call h5dget_space_f(dset_id, filespace, ierr)
225 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
226 doffset, dcount, ierr)
227 call h5dwrite_f(dset_id, h5t_neko_real, &
228 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
229 ddim, ierr, file_space_id = filespace, &
230 mem_space_id = memspace, xfer_prp = plist_id)
231 call h5dclose_f(dset_id, ierr)
232 end do
233 end do
234 deallocate(fsp)
235 end if
236
237 call h5gclose_f(grp_id, ierr)
238 call h5sclose_f(filespace, ierr)
239 call h5sclose_f(memspace, ierr)
240 end if
241
242 call h5pclose_f(plist_id, ierr)
243 call h5fclose_f(file_id, ierr)
244
245 call h5close_f(ierr)
246
247 end subroutine hdf5_file_write
248
250 subroutine hdf5_file_read(this, data)
251 class(hdf5_file_t) :: this
252 class(*), target, intent(inout) :: data
253 integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
254 integer(hid_t) :: filespace, memspace
255 integer(hid_t) :: h5t_neko_real
256 integer(hsize_t), dimension(1) :: ddim, dcount, doffset
257 integer :: i,j, ierr, info, glb_nelv, gdim, lx, drank
258 type(mesh_t), pointer :: msh
259 type(dofmap_t), pointer :: dof
260 type(field_ptr_t), allocatable :: fp(:)
261 type(field_series_ptr_t), allocatable :: fsp(:)
262 real(kind=rp), pointer :: dtlag(:)
263 real(kind=rp), pointer :: tlag(:)
264 real(kind=rp) :: t
265 character(len=1024) :: fname
266
267 fname = trim(this%get_fname())
268
269 call hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
270
271 call h5open_f(ierr)
272
273 call hdf5_file_determine_real(h5t_neko_real)
274
275 call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
276 info = mpi_info_null%mpi_val
277 call h5pset_fapl_mpio_f(plist_id, neko_comm%mpi_val, info, ierr)
278
279 call h5fopen_f(fname, h5f_acc_rdonly_f, &
280 file_id, ierr, access_prp = plist_id)
281
282 call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
283 call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
284
285 ddim = 1
286 call h5aopen_name_f(file_id, 'Time', attr_id, ierr)
287 call h5aread_f(attr_id, h5t_neko_real, t, ddim, ierr)
288 call h5aclose_f(attr_id, ierr)
289
290 select type(data)
291 type is(chkp_t)
292 data%t = t
293 end select
294
295 call h5aopen_name_f(file_id, 'Lx', attr_id, ierr)
296 call h5aread_f(attr_id, h5t_native_integer, lx, ddim, ierr)
297 call h5aclose_f(attr_id, ierr)
298
299 call h5gopen_f(file_id, 'Mesh', grp_id, ierr, gapl_id=h5p_default_f)
300
301 call h5aopen_name_f(grp_id, 'Elements', attr_id, ierr)
302 call h5aread_f(attr_id, h5t_native_integer, glb_nelv, ddim, ierr)
303 call h5aclose_f(attr_id, ierr)
304
305 call h5aopen_name_f(grp_id, 'Dimension', attr_id, ierr)
306 call h5aread_f(attr_id, h5t_native_integer, gdim, ddim, ierr)
307 call h5aclose_f(attr_id, ierr)
308 call h5gclose_f(grp_id, ierr)
309
310
311 if (associated(tlag) .and. associated(dtlag)) then
312 drank = 1
313 ddim = size(tlag)
314 doffset(1) = 0
315 if (pe_rank .eq. 0) then
316 dcount = size(tlag)
317 else
318 dcount = 0
319 end if
320
321 call h5gopen_f(file_id, 'Restart', grp_id, ierr, gapl_id=h5p_default_f)
322 call h5dopen_f(grp_id, 'tlag', dset_id, ierr)
323 call h5dget_space_f(dset_id, filespace, ierr)
324 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
325 doffset, dcount, ierr)
326 call h5dread_f(dset_id, h5t_neko_real, tlag, ddim, ierr, xfer_prp=plist_id)
327 call h5dclose_f(dset_id, ierr)
328 call h5sclose_f(filespace, ierr)
329
330 call h5dopen_f(grp_id, 'dtlag', dset_id, ierr)
331 call h5dget_space_f(dset_id, filespace, ierr)
332 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
333 doffset, dcount, ierr)
334 call h5dread_f(dset_id, h5t_neko_real, dtlag, ddim, ierr, xfer_prp=plist_id)
335 call h5dclose_f(dset_id, ierr)
336 call h5sclose_f(filespace, ierr)
337
338 call h5gclose_f(grp_id, ierr)
339 end if
340
341 if (allocated(fp) .or. allocated(fsp)) then
342 call h5gopen_f(file_id, 'Fields', grp_id, ierr, gapl_id=h5p_default_f)
343
344 dcount(1) = int(dof%size(), 8)
345 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
346 ddim = int(dof%size(), 8)
347 drank = 1
348
349 dcount(1) = int(dof%size(), 8)
350 doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
351 ddim = int(dof%size(), 8)
352 drank = 1
353 call mpi_allreduce(mpi_in_place, ddim(1), 1, &
354 mpi_integer8, mpi_sum, neko_comm, ierr)
355
356 call h5screate_simple_f(drank, dcount, memspace, ierr)
357
358 if (allocated(fp)) then
359 do i = 1, size(fp)
360 call h5dopen_f(grp_id, fp(i)%ptr%name, dset_id, ierr)
361 call h5dget_space_f(dset_id, filespace, ierr)
362 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
363 doffset, dcount, ierr)
364 call h5dread_f(dset_id, h5t_neko_real, &
365 fp(i)%ptr%x(1,1,1,1), &
366 ddim, ierr, file_space_id = filespace, &
367 mem_space_id = memspace, xfer_prp=plist_id)
368 call h5dclose_f(dset_id, ierr)
369 call h5sclose_f(filespace, ierr)
370 end do
371 end if
372
373 if (allocated(fsp)) then
374 do i = 1, size(fsp)
375 do j = 1, fsp(i)%ptr%size()
376 call h5dopen_f(grp_id, fsp(i)%ptr%lf(j)%name, dset_id, ierr)
377 call h5dget_space_f(dset_id, filespace, ierr)
378 call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
379 doffset, dcount, ierr)
380 call h5dread_f(dset_id, h5t_neko_real, &
381 fsp(i)%ptr%lf(j)%x(1,1,1,1), &
382 ddim, ierr, file_space_id = filespace, &
383 mem_space_id = memspace, xfer_prp=plist_id)
384 call h5dclose_f(dset_id, ierr)
385 call h5sclose_f(filespace, ierr)
386 end do
387 end do
388 end if
389 call h5sclose_f(memspace, ierr)
390 call h5gclose_f(grp_id, ierr)
391 end if
392
393 call h5pclose_f(plist_id, ierr)
394 call h5fclose_f(file_id, ierr)
395
396 call h5close_f(ierr)
397
398 end subroutine hdf5_file_read
399
400
401 subroutine hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
402 class(*), target, intent(in) :: data
403 type(mesh_t), pointer, intent(inout) :: msh
404 type(dofmap_t), pointer, intent(inout) :: dof
405 type(field_ptr_t), allocatable, intent(inout) :: fp(:)
406 type(field_series_ptr_t), allocatable, intent(inout) :: fsp(:)
407 real(kind=rp), pointer, intent(inout) :: dtlag(:)
408 real(kind=rp), pointer, intent(inout) :: tlag(:)
409 integer :: i, j, fp_size, fp_cur, fsp_size, fsp_cur, scalar_count, ab_count
410 character(len=32) :: scalar_name
411
412 select type(data)
413 type is (field_t)
414 dof => data%dof
415 msh => data%msh
416 fp_size = 1
417 allocate(fp(fp_size))
418 fp(1)%ptr => data
419
420 nullify(dtlag)
421 nullify(tlag)
422
423 type is (field_list_t)
424
425 if (data%size() .gt. 0) then
426 allocate(fp(data%size()))
427
428 dof => data%dof(1)
429 msh => data%msh(1)
430
431 do i = 1, data%size()
432 fp(i)%ptr => data%items(i)%ptr
433 end do
434 else
435 call neko_error('Empty field list')
436 end if
437
438 nullify(dtlag)
439 nullify(tlag)
440
441 type is(chkp_t)
442
443 if ( .not. associated(data%u) .or. &
444 .not. associated(data%v) .or. &
445 .not. associated(data%w) .or. &
446 .not. associated(data%p) ) then
447 call neko_error('Checkpoint not initialized')
448 end if
449
450 fp_size = 4
451
452 if (allocated(data%scalar_lags%items) .and. data%scalar_lags%size() > 0) then
453 scalar_count = data%scalar_lags%size()
454 else if (associated(data%s)) then
455 scalar_count = 1
456 else
457 scalar_count = 0
458 end if
459
460 if (scalar_count .gt. 1) then
461 fp_size = fp_size + scalar_count
462
463 ! Add abx1 and abx2 fields for each scalar
464 fp_size = fp_size + (scalar_count * 2)
465 else if (associated(data%s)) then
466 ! Single scalar support
467 fp_size = fp_size + 1
468 if (associated(data%abs1)) then
469 fp_size = fp_size + 2
470 end if
471 end if
472
473 if (associated(data%abx1)) then
474 fp_size = fp_size + 6
475 end if
476
477 allocate(fp(fp_size))
478
479 fsp_size = 0
480 if (associated(data%ulag)) then
481 fsp_size = fsp_size + 3
482 end if
483
484 if (scalar_count .gt. 1) then
485 if (allocated(data%scalar_lags%items)) then
486 fsp_size = fsp_size + data%scalar_lags%size()
487 end if
488 else if (associated(data%slag)) then
489 fsp_size = fsp_size + 1
490 end if
491
492 if (fsp_size .gt. 0) then
493 allocate(fsp(fsp_size))
494 fsp_cur = 1
495 end if
496
497 dof => data%u%dof
498 msh => data%u%msh
499
500 fp(1)%ptr => data%u
501 fp(2)%ptr => data%v
502 fp(3)%ptr => data%w
503 fp(4)%ptr => data%p
504
505 fp_cur = 5
506
507 if (scalar_count .gt. 1) then
508 block
509 type(field_series_t), pointer :: slag
510 do i = 1, scalar_count
511 slag => data%scalar_lags%get(i)
512 fp(fp_cur)%ptr => slag%f
513 fp_cur = fp_cur + 1
514 end do
515 end block
516
517 do i = 1, scalar_count
518 fp(fp_cur)%ptr => data%scalar_abx1(i)%ptr
519 fp_cur = fp_cur + 1
520 fp(fp_cur)%ptr => data%scalar_abx2(i)%ptr
521 fp_cur = fp_cur + 1
522 end do
523 else if (associated(data%s)) then
524 ! Single scalar support
525 fp(fp_cur)%ptr => data%s
526 fp_cur = fp_cur + 1
527
528 if (associated(data%abs1)) then
529 fp(fp_cur)%ptr => data%abs1
530 fp(fp_cur+1)%ptr => data%abs2
531 fp_cur = fp_cur + 2
532 end if
533 end if
534
535 if (associated(data%abx1)) then
536 fp(fp_cur)%ptr => data%abx1
537 fp(fp_cur+1)%ptr => data%abx2
538 fp(fp_cur+2)%ptr => data%aby1
539 fp(fp_cur+3)%ptr => data%aby2
540 fp(fp_cur+4)%ptr => data%abz1
541 fp(fp_cur+5)%ptr => data%abz2
542 fp_cur = fp_cur + 6
543 end if
544
545 if (associated(data%ulag)) then
546 fsp(fsp_cur)%ptr => data%ulag
547 fsp(fsp_cur+1)%ptr => data%vlag
548 fsp(fsp_cur+2)%ptr => data%wlag
549 fsp_cur = fsp_cur + 3
550 end if
551
552
553 if (scalar_count .gt. 1) then
554 if (allocated(data%scalar_lags%items)) then
555 do j = 1, data%scalar_lags%size()
556 fsp(fsp_cur)%ptr => data%scalar_lags%get(j)
557 fsp_cur = fsp_cur + 1
558 end do
559 end if
560 else if (associated(data%slag)) then
561 fsp(fsp_cur)%ptr => data%slag
562 fsp_cur = fsp_cur + 1
563 end if
564
565 if (associated(data%tlag)) then
566 tlag => data%tlag
567 dtlag => data%dtlag
568 end if
569
570 class default
571 call neko_log%error('Invalid data')
572 end select
573
574 end subroutine hdf5_file_determine_data
575
579 subroutine hdf5_file_determine_real(H5T_NEKO_REAL)
580 integer(hid_t), intent(inout) :: h5t_neko_real
581 select case(rp)
582 case(dp)
583 h5t_neko_real = h5t_native_double
584 case(sp)
585 h5t_neko_real = h5t_native_real
586 case default
587 call neko_error("Unsupported real type")
588 end select
589 end subroutine hdf5_file_determine_real
590
592 subroutine hdf5_file_set_overwrite(this, overwrite)
593 class(hdf5_file_t), intent(inout) :: this
594 logical, intent(in) :: overwrite
595 this%overwrite = overwrite
596 end subroutine hdf5_file_set_overwrite
597
598#else
599
601 subroutine hdf5_file_write(this, data, t)
602 class(hdf5_file_t), intent(inout) :: this
603 class(*), target, intent(in) :: data
604 real(kind=rp), intent(in), optional :: t
605 call neko_error('Neko needs to be built with HDF5 support')
606 end subroutine hdf5_file_write
607
609 subroutine hdf5_file_read(this, data)
610 class(hdf5_file_t) :: this
611 class(*), target, intent(inout) :: data
612 call neko_error('Neko needs to be built with HDF5 support')
613 end subroutine hdf5_file_read
614
616 subroutine hdf5_file_set_overwrite(this, overwrite)
617 class(hdf5_file_t), intent(inout) :: this
618 logical, intent(in) :: overwrite
619 call neko_error('Neko needs to be built with HDF5 support')
620 end subroutine hdf5_file_set_overwrite
621
622#endif
623
624end module hdf5_file
Defines a checkpoint.
Definition comm.F90:1
integer, public pe_rank
MPI rank.
Definition comm.F90:55
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
HDF5 file format.
Definition hdf5_file.F90:34
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_set_overwrite(this, overwrite)
Set the overwrite flag for HDF5 files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
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
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition utils.f90:58
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.
Definition hdf5_file.F90:56