Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
mma_io_hdf5.f90
1
34
36submodule(mma) mma_io_hdf5
37#ifdef HAVE_HDF5
38 use hdf5
39#endif
40 use mpi_f08, only: mpi_info_null, mpi_integer8
41
42contains
43
44#ifdef HAVE_HDF5
45
51 module subroutine mma_write_hdf5(this, filename)
52 class(mma_t), intent(inout) :: this
53 character(len=*), intent(in) :: filename
54 integer(hid_t) :: fapl_id, xf_id, file_id, dset_id, filespace, memspace, &
55 attr_id, grp_id, mma_grp_id, str_type
56 integer(hid_t) :: H5T_NEKO_REAL
57 integer(hsize_t), dimension(1) :: ddim, dcount, doffset
58 integer :: ierr, info, drank
59 integer(kind=8) :: local_n8, prefix8, total_n8
60
61
62 ! Ensure device state is on host
63 call this%copy_from(device_to_host, sync = .true.)
64
65 call h5open_f(ierr)
66
67 ! Assign the correct HDF5 data type based on the neko real kind
68 select case (rp)
69 case (dp)
70 h5t_neko_real = h5t_native_double
71 case (sp)
72 h5t_neko_real = h5t_native_real
73 case default
74 call neko_error('mma: unsupported real kind for HDF5')
75 end select
76
77 ! Create file with MPIO access
78 call h5pcreate_f(h5p_file_access_f, fapl_id, ierr)
79 info = mpi_info_null%mpi_val
80 call h5pset_fapl_mpio_f(fapl_id, neko_comm%mpi_val, info, ierr)
81
82 call h5fcreate_f(trim(filename), h5f_acc_trunc_f, file_id, ierr, &
83 access_prp = fapl_id)
84 call h5gcreate_f(file_id, "MMA", mma_grp_id, ierr, &
85 lcpl_id = h5p_default_f, gcpl_id = h5p_default_f, &
86 gapl_id = h5p_default_f)
87
88 ! ------------------------------------------------------------------------ !
89 ! Write basic Parameters attributes
90
91 call h5gcreate_f(mma_grp_id, "Parameters", grp_id, ierr, &
92 lcpl_id = h5p_default_f, gcpl_id = h5p_default_f, &
93 gapl_id = h5p_default_f)
94
95 call h5screate_f(h5s_scalar_f, filespace, ierr)
96 ddim = 1
97
98 ! Integer-valued attributes
99
100 call h5acreate_f(grp_id, 'n', h5t_native_integer, filespace, attr_id, &
101 ierr, h5p_default_f, h5p_default_f)
102 call h5awrite_f(attr_id, h5t_native_integer, this%n, ddim, ierr)
103 call h5aclose_f(attr_id, ierr)
104
105 call h5acreate_f(grp_id, 'm', h5t_native_integer, filespace, attr_id, &
106 ierr, h5p_default_f, h5p_default_f)
107 call h5awrite_f(attr_id, h5t_native_integer, this%m, ddim, ierr)
108 call h5aclose_f(attr_id, ierr)
109
110 call h5acreate_f(grp_id, 'max_iter', h5t_native_integer, filespace, &
111 attr_id, ierr, h5p_default_f, h5p_default_f)
112 call h5awrite_f(attr_id, h5t_native_integer, this%max_iter, ddim, ierr)
113 call h5aclose_f(attr_id, ierr)
114
115 ! Real-valued attributes
116
117 call h5acreate_f(grp_id, 'asyinit', h5t_neko_real, filespace, attr_id, &
118 ierr, h5p_default_f, h5p_default_f)
119 call h5awrite_f(attr_id, h5t_neko_real, this%asyinit, ddim, ierr)
120 call h5aclose_f(attr_id, ierr)
121
122 call h5acreate_f(grp_id, 'asyincr', h5t_neko_real, filespace, attr_id, &
123 ierr, h5p_default_f, h5p_default_f)
124 call h5awrite_f(attr_id, h5t_neko_real, this%asyincr, ddim, ierr)
125 call h5aclose_f(attr_id, ierr)
126
127 call h5acreate_f(grp_id, 'asydecr', h5t_neko_real, filespace, attr_id, &
128 ierr, h5p_default_f, h5p_default_f)
129 call h5awrite_f(attr_id, h5t_neko_real, this%asydecr, ddim, ierr)
130 call h5aclose_f(attr_id, ierr)
131
132 call h5acreate_f(grp_id, 'epsimin', h5t_neko_real, filespace, attr_id, &
133 ierr, h5p_default_f, h5p_default_f)
134 call h5awrite_f(attr_id, h5t_neko_real, this%epsimin, ddim, ierr)
135 call h5aclose_f(attr_id, ierr)
136
137 ! String-valued attributes
138
139 ! Create the string type
140 call h5tcopy_f(h5t_fortran_s1, str_type, ierr)
141 call h5tset_strpad_f(str_type, h5t_str_spacepad_f, ierr)
142
143 ddim(1) = len_trim(this%bcknd)
144 call h5tset_size_f(str_type, ddim(1), ierr)
145 call h5acreate_f(grp_id, 'bcknd', str_type, filespace, attr_id, &
146 ierr, h5p_default_f, h5p_default_f)
147 call h5awrite_f(attr_id, str_type, this%bcknd, ddim, ierr)
148 call h5aclose_f(attr_id, ierr)
149
150 ddim(1) = len_trim(this%subsolver)
151 call h5tset_size_f(str_type, ddim(1), ierr)
152 call h5acreate_f(grp_id, 'subsolver', str_type, filespace, attr_id, &
153 ierr, h5p_default_f, h5p_default_f)
154 call h5awrite_f(attr_id, str_type, this%subsolver, ddim, ierr)
155 call h5aclose_f(attr_id, ierr)
156
157 ! Close the string type
158 call h5tclose_f(str_type, ierr)
159
160 ! Close the filespace and Parameters group
161 call h5sclose_f(filespace, ierr)
162 call h5gclose_f(grp_id, ierr)
163
164 ! ------------------------------------------------------------------------ !
165 ! Global arrays datasets
166
167 ddim(1) = 1
168 call h5screate_f(h5s_scalar_f, filespace, ierr)
169
170 call h5dcreate_f(mma_grp_id, 'a0', h5t_neko_real, filespace, dset_id, ierr)
171 call h5dwrite_f(dset_id, h5t_neko_real, this%a0, ddim, ierr)
172 call h5dclose_f(dset_id, ierr)
173
174 ! The next batch are vectors of size m
175 ddim(1) = this%m
176 drank = 1
177
178 call h5screate_simple_f(drank, ddim, filespace, ierr)
179
180 call h5dcreate_f(mma_grp_id, 'a', h5t_neko_real, filespace, dset_id, ierr)
181 call h5dwrite_f(dset_id, h5t_neko_real, this%a%x(1), ddim, ierr)
182 call h5dclose_f(dset_id, ierr)
183
184 call h5dcreate_f(mma_grp_id, 'c', h5t_neko_real, filespace, dset_id, ierr)
185 call h5dwrite_f(dset_id, h5t_neko_real, this%c%x(1), ddim, ierr)
186 call h5dclose_f(dset_id, ierr)
187
188 call h5dcreate_f(mma_grp_id, 'd', h5t_neko_real, filespace, dset_id, ierr)
189 call h5dwrite_f(dset_id, h5t_neko_real, this%d%x(1), ddim, ierr)
190 call h5dclose_f(dset_id, ierr)
191
192 ! Close the filespace and group
193 call h5sclose_f(filespace, ierr)
194
195 ! ------------------------------------------------------------------------ !
196 ! Write per-rank datasets
197
198 ! Define the sizes and offsets
199 local_n8 = int(this%n, 8)
200 total_n8 = int(this%n_global, 8)
201 call mpi_scan(local_n8, prefix8, 1, mpi_integer8, mpi_sum, neko_comm, ierr)
202
203 drank = 1
204 dcount(1) = local_n8
205 doffset(1) = prefix8 - local_n8
206 ddim(1) = total_n8
207
208 ! Create a file and memory space
209 call h5screate_simple_f(drank, ddim, filespace, ierr)
210 call h5screate_simple_f(drank, dcount, memspace, ierr)
211
212 ! Dataset transfer property list: collective
213 call h5pcreate_f(h5p_dataset_xfer_f, xf_id, ierr)
214 call h5pset_dxpl_mpio_f(xf_id, h5fd_mpio_collective_f, ierr)
215
216 ! Select hyperslab in the file space
217 call h5sselect_hyperslab_f(filespace, h5s_select_set_f, doffset, dcount, &
218 ierr)
219
220 ! Write xmin
221 call h5dcreate_f(mma_grp_id, 'xmin', h5t_neko_real, &
222 filespace, dset_id, ierr)
223 call h5dwrite_f(dset_id, h5t_neko_real, this%xmin%x(1), ddim, ierr, &
224 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
225 call h5dclose_f(dset_id, ierr)
226
227 ! Write xmax
228 call h5dcreate_f(mma_grp_id, 'xmax', h5t_neko_real, &
229 filespace, dset_id, ierr)
230 call h5dwrite_f(dset_id, h5t_neko_real, this%xmax%x(1), ddim, ierr, &
231 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
232 call h5dclose_f(dset_id, ierr)
233
234 ! Write xold1
235 call h5dcreate_f(mma_grp_id, 'xold1', h5t_neko_real, &
236 filespace, dset_id, ierr)
237 call h5dwrite_f(dset_id, h5t_neko_real, this%xold1%x(1), ddim, ierr, &
238 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
239 call h5dclose_f(dset_id, ierr)
240
241 ! Write xold2
242 call h5dcreate_f(mma_grp_id, 'xold2', h5t_neko_real, &
243 filespace, dset_id, ierr)
244 call h5dwrite_f(dset_id, h5t_neko_real, this%xold2%x(1), ddim, ierr, &
245 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
246 call h5dclose_f(dset_id, ierr)
247
248 ! Write low
249 call h5dcreate_f(mma_grp_id, 'low', h5t_neko_real, &
250 filespace, dset_id, ierr)
251 call h5dwrite_f(dset_id, h5t_neko_real, this%low%x(1), ddim, ierr, &
252 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
253 call h5dclose_f(dset_id, ierr)
254
255 ! Write upp
256 call h5dcreate_f(mma_grp_id, 'upp', h5t_neko_real, &
257 filespace, dset_id, ierr)
258 call h5dwrite_f(dset_id, h5t_neko_real, this%upp%x(1), ddim, ierr, &
259 file_space_id = filespace, mem_space_id = memspace, xfer_prp = xf_id)
260 call h5dclose_f(dset_id, ierr)
261
262 ! Close the spaces used
263 call h5sclose_f(filespace, ierr)
264 call h5sclose_f(memspace, ierr)
265 call h5pclose_f(xf_id, ierr)
266
267 ! ------------------------------------------------------------------------ !
268 ! Close the group and file
269
270 call h5gclose_f(mma_grp_id, ierr)
271 call h5fclose_f(file_id, ierr)
272 call h5pclose_f(fapl_id, ierr)
273 call h5close_f(ierr)
274
275 end subroutine mma_write_hdf5
276
280 module subroutine mma_read_hdf5(this, filename)
281 class(mma_t), intent(inout) :: this
282 character(len=*), intent(in) :: filename
283 integer(hid_t) :: fapl_id, file_id, dset_id, &
284 attr_id, grp_id, str_type
285 integer(hid_t) :: H5T_NEKO_REAL
286 integer(hsize_t), dimension(1) :: ddim
287 integer :: ierr
288 integer :: n, m, max_iter
289 real(kind=rp) :: asyinit, asyincr, asydecr, epsimin
290
291 character(len=12) :: bcknd, subsolver
292
293 ! Initialize reader values
294 n = -1
295 m = -1
296 max_iter = -1
297 asyinit = -1.0_rp
298 asyincr = -1.0_rp
299 asydecr = -1.0_rp
300 epsimin = -1.0_rp
301
302 bcknd = ''
303 subsolver = ''
304
305 ! Open file and prepare reading
306 call h5open_f(ierr)
307 call h5pcreate_f(h5p_file_access_f, fapl_id, ierr)
308 call h5pset_fapl_mpio_f(fapl_id, neko_comm%mpi_val, &
309 mpi_info_null%mpi_val, ierr)
310 call h5fopen_f(trim(filename), h5f_acc_rdonly_f, file_id, ierr, &
311 access_prp = fapl_id)
312
313 ! Assign the correct HDF5 data type based on the neko real kind
314 select case (rp)
315 case (dp)
316 h5t_neko_real = h5t_native_double
317 case (sp)
318 h5t_neko_real = h5t_native_real
319 case default
320 call neko_error('mma: unsupported real kind for HDF5')
321 end select
322
323 ! ------------------------------------------------------------------------ !
324 ! Read basic Parameters attributes
325
326 call h5gopen_f(file_id, 'MMA/Parameters', grp_id, ierr)
327
328 ddim(1) = 1
329 call h5aopen_f(grp_id, 'n', attr_id, ierr)
330 call h5aread_f(attr_id, h5t_native_integer, n, ddim, ierr)
331 call h5aclose_f(attr_id, ierr)
332
333 call h5aopen_f(grp_id, 'm', attr_id, ierr)
334 call h5aread_f(attr_id, h5t_native_integer, m, ddim, ierr)
335 call h5aclose_f(attr_id, ierr)
336
337 call h5aopen_f(grp_id, 'max_iter', attr_id, ierr)
338 call h5aread_f(attr_id, h5t_native_integer, max_iter, ddim, ierr)
339 call h5aclose_f(attr_id, ierr)
340
341 call h5aopen_f(grp_id, 'asyinit', attr_id, ierr)
342 call h5aread_f(attr_id, h5t_neko_real, asyinit, ddim, ierr)
343 call h5aclose_f(attr_id, ierr)
344
345 call h5aopen_f(grp_id, 'asyincr', attr_id, ierr)
346 call h5aread_f(attr_id, h5t_neko_real, asyincr, ddim, ierr)
347 call h5aclose_f(attr_id, ierr)
348
349 call h5aopen_f(grp_id, 'asydecr', attr_id, ierr)
350 call h5aread_f(attr_id, h5t_neko_real, asydecr, ddim, ierr)
351 call h5aclose_f(attr_id, ierr)
352
353 call h5aopen_f(grp_id, 'epsimin', attr_id, ierr)
354 call h5aread_f(attr_id, h5t_neko_real, epsimin, ddim, ierr)
355 call h5aclose_f(attr_id, ierr)
356
357 ! Read strings
358 call h5aopen_f(grp_id, 'bcknd', attr_id, ierr)
359 call h5aget_type_f(attr_id, str_type, ierr)
360 call h5aread_f(attr_id, str_type, bcknd, ddim, ierr)
361 call h5aclose_f(attr_id, ierr)
362
363 call h5aopen_f(grp_id, 'subsolver', attr_id, ierr)
364 call h5aget_type_f(attr_id, str_type, ierr)
365 call h5aread_f(attr_id, str_type, subsolver, ddim, ierr)
366 call h5aclose_f(attr_id, ierr)
367
368 call h5tclose_f(str_type, ierr)
369 call h5gclose_f(grp_id, ierr)
370
371 ! Ensure the MMA object is allocated with the same configuration as the file
372 if (n .ne. this%n) then
373 call neko_error('mma: mismatch in n during HDF5 read')
374 end if
375 if (m .ne. this%m) then
376 call neko_error('mma: mismatch in m during HDF5 read')
377 end if
378 if (max_iter .ne. this%max_iter) then
379 call neko_error('mma: mismatch in max_iter during HDF5 read')
380 end if
381 if (abs(asyinit - this%asyinit) .gt. 1.0e-12_rp) then
382 call neko_error('mma: mismatch in asyinit during HDF5 read')
383 end if
384 if (abs(asyincr - this%asyincr) .gt. 1.0e-12_rp) then
385 call neko_error('mma: mismatch in asyincr during HDF5 read')
386 end if
387 if (abs(asydecr - this%asydecr) .gt. 1.0e-12_rp) then
388 call neko_error('mma: mismatch in asydecr during HDF5 read')
389 end if
390 if (abs(epsimin - this%epsimin) .gt. 1.0e-12_rp) then
391 call neko_error('mma: mismatch in epsimin during HDF5 read')
392 end if
393 if (trim(bcknd) .ne. trim(this%bcknd)) then
394 call neko_error('mma: mismatch in bcknd during HDF5 read')
395 end if
396 if (trim(subsolver) .ne. trim(this%subsolver)) then
397 call neko_error('mma: mismatch in subsolver during HDF5 read')
398 end if
399
400 ! ------------------------------------------------------------------------ !
401 ! Global arrays datasets
402
403 call h5gopen_f(file_id, 'MMA', grp_id, ierr)
404
405 ! Read penalty parameters
406 ddim(1) = 1
407 call h5dopen_f(file_id, '/MMA/a0', dset_id, ierr)
408 call h5dread_f(dset_id, h5t_neko_real, this%a0, ddim, ierr)
409 call h5dclose_f(dset_id, ierr)
410
411 ddim(1) = m
412 call h5dopen_f(file_id, '/MMA/a', dset_id, ierr)
413 call h5dread_f(dset_id, h5t_neko_real, this%a%x(1), ddim, ierr)
414 call h5dclose_f(dset_id, ierr)
415
416 call h5dopen_f(file_id, '/MMA/c', dset_id, ierr)
417 call h5dread_f(dset_id, h5t_neko_real, this%c%x(1), ddim, ierr)
418 call h5dclose_f(dset_id, ierr)
419
420 call h5dopen_f(file_id, '/MMA/d', dset_id, ierr)
421 call h5dread_f(dset_id, h5t_neko_real, this%d%x(1), ddim, ierr)
422 call h5dclose_f(dset_id, ierr)
423
424 ! ------------------------------------------------------------------------ !
425 ! Read per-rank datasets
426
427 ddim(1) = this%n
428
429 call h5dopen_f(file_id, '/MMA/xmin', dset_id, ierr)
430 call h5dread_f(dset_id, h5t_neko_real, this%xmin%x(1), ddim, ierr)
431 call h5dclose_f(dset_id, ierr)
432
433 call h5dopen_f(file_id, '/MMA/xmax', dset_id, ierr)
434 call h5dread_f(dset_id, h5t_neko_real, this%xmax%x(1), ddim, ierr)
435 call h5dclose_f(dset_id, ierr)
436
437 call h5dopen_f(file_id, '/MMA/xold1', dset_id, ierr)
438 call h5dread_f(dset_id, h5t_neko_real, this%xold1%x(1), ddim, ierr)
439 call h5dclose_f(dset_id, ierr)
440
441 call h5dopen_f(file_id, '/MMA/xold2', dset_id, ierr)
442 call h5dread_f(dset_id, h5t_neko_real, this%xold2%x(1), ddim, ierr)
443 call h5dclose_f(dset_id, ierr)
444
445 call h5dopen_f(file_id, '/MMA/low', dset_id, ierr)
446 call h5dread_f(dset_id, h5t_neko_real, this%low%x(1), ddim, ierr)
447 call h5dclose_f(dset_id, ierr)
448
449 call h5dopen_f(file_id, '/MMA/upp', dset_id, ierr)
450 call h5dread_f(dset_id, h5t_neko_real, this%upp%x(1), ddim, ierr)
451 call h5dclose_f(dset_id, ierr)
452
453 ! ------------------------------------------------------------------------ !
454 ! Close the group and file
455
456 call h5gclose_f(grp_id, ierr)
457 call h5fclose_f(file_id, ierr)
458 call h5pclose_f(fapl_id, ierr)
459 call h5close_f(ierr)
460
461 ! Ensure device state is updated
462 call this%copy_from(host_to_device, sync = .true.)
463
464 end subroutine mma_read_hdf5
465
466#else
467
468 module subroutine mma_write_hdf5(this, filename)
469 class(mma_t), intent(inout) :: this
470 character(len=*), intent(in) :: filename
471 call neko_error('mma: HDF5 support not enabled rebuild with HAVE_HDF5')
472 end subroutine mma_write_hdf5
473
474 module subroutine mma_read_hdf5(this, filename)
475 class(mma_t), intent(inout) :: this
476 character(len=*), intent(in) :: filename
477 call neko_error('mma: HDF5 support not enabled rebuild with HAVE_HDF5')
478 end subroutine mma_read_hdf5
479
480#endif
481
482end submodule mma_io_hdf5
MMA module.
Definition mma.f90:69