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