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