Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_WriteBinaryParticleDataHDF5.H
Go to the documentation of this file.
1#ifndef AMREX_WRITE_BINARY_PARTICLE_DATA_HDF5_H
2#define AMREX_WRITE_BINARY_PARTICLE_DATA_HDF5_H
3
4#include <iterator>
5
13namespace amrex {
14
25static int CreateWriteHDF5Attr(hid_t loc, const char *name, hsize_t n, void *data, hid_t dtype)
26{
27 herr_t ret;
28 hid_t attr, attr_space;
29 hsize_t dims = n;
30
31 attr_space = H5Screate_simple(1, &dims, NULL);
32
33 attr = H5Acreate(loc, name, dtype, attr_space, H5P_DEFAULT, H5P_DEFAULT);
34 if (attr < 0) {
35 amrex::Print() << " Error with H5Acreate [" << name << "]\n";
36 H5Sclose(attr_space);
37 return -1;
38 }
39
40 ret = H5Awrite(attr, dtype, (void*)data);
41 if (ret < 0) {
42 amrex::Print() << " Error with H5Awrite [" << name << "]\n";
43 H5Sclose(attr_space);
44 H5Aclose(attr);
45 return -1;
46 }
47 H5Sclose(attr_space);
48 H5Aclose(attr);
49 return 1;
50}
51
60static int CreateWriteHDF5AttrString(hid_t loc, const char *name, const char* str)
61{
62 hid_t attr, atype, space;
63 herr_t ret;
64
65 AMREX_ASSERT(name);
66 AMREX_ASSERT(str);
67
68 space = H5Screate(H5S_SCALAR);
69 atype = H5Tcopy(H5T_C_S1);
70 H5Tset_size(atype, strlen(str)+1);
71 H5Tset_strpad(atype,H5T_STR_NULLTERM);
72 attr = H5Acreate(loc, name, atype, space, H5P_DEFAULT, H5P_DEFAULT);
73 if (attr < 0) {
74 amrex::Print() << " Error with H5Acreate [" << name << "]\n";
75 return -1;
76 }
77
78 ret = H5Awrite(attr, atype, str);
79 if (ret < 0) {
80 amrex::Print() << " Error with H5Awrite [" << name << "]\n";
81 return -1;
82 }
83
84 H5Tclose(atype);
85 H5Sclose(space);
86 H5Aclose(attr);
87
88 return 1;
89}
90
100static int ReadHDF5Attr(hid_t loc, const char *name, void *data, hid_t dtype)
101{
102 herr_t ret;
103 hid_t attr;
104
105 attr = H5Aopen(loc, name, H5P_DEFAULT);
106 if (attr < 0) {
107 amrex::Print() << " Error with H5Aopen [" << name << "]\n";
108 return -1;
109 }
110
111 ret = H5Aread(attr, dtype, data);
112 if (ret < 0) {
113 amrex::Print() << " Error with H5Aread [" << name << "]\n";
114 return -1;
115 }
116 H5Aclose(attr);
117 return 1;
118}
119
120#ifdef AMREX_USE_MPI
127static void SetHDF5fapl(hid_t fapl, MPI_Comm comm)
128#else
134static void SetHDF5fapl(hid_t fapl)
135#endif
136{
137#ifdef AMREX_USE_MPI
138 H5Pset_fapl_mpio(fapl, comm, MPI_INFO_NULL);
139
140 // Alignment and metadata block size
141 int alignment = 16 * 1024 * 1024;
142 int blocksize = 4 * 1024 * 1024;
143 H5Pset_alignment(fapl, alignment, alignment);
144 H5Pset_meta_block_size(fapl, blocksize);
145
146 // Collective metadata ops
147 H5Pset_coll_metadata_write(fapl, true);
148 H5Pset_all_coll_metadata_ops(fapl, true);
149
150 // Defer cache flush
151 H5AC_cache_config_t cache_config;
152 cache_config.version = H5AC__CURR_CACHE_CONFIG_VERSION;
153 H5Pget_mdc_config(fapl, &cache_config);
154 cache_config.set_initial_size = 1;
155 cache_config.initial_size = 16 * 1024 * 1024;
156 cache_config.evictions_enabled = 0;
157 cache_config.incr_mode = H5C_incr__off;
158 cache_config.flash_incr_mode = H5C_flash_incr__off;
159 cache_config.decr_mode = H5C_decr__off;
160 H5Pset_mdc_config (fapl, &cache_config);
161#else
162 H5Pset_fapl_sec2(fapl);
163#endif
164
165}
183template <class PC, class F, std::enable_if_t<IsParticleContainer<PC>::value, int> foo = 0>
184void WriteHDF5ParticleDataSync (PC const& pc,
185 const std::string& dir, const std::string& name,
186 const Vector<int>& write_real_comp,
187 const Vector<int>& write_int_comp,
188 const Vector<std::string>& real_comp_names,
189 const Vector<std::string>& int_comp_names,
190 const std::string& compression,
191 F&& f, bool is_checkpoint)
192{
193 BL_PROFILE("WriteHDF5ParticleDataSync()");
194 AMREX_ASSERT(pc.OK());
195
196 AMREX_ASSERT(sizeof(typename PC::ParticleType::RealType) == 4 ||
197 sizeof(typename PC::ParticleType::RealType) == 8);
198
199 constexpr int NStructReal = PC::NStructReal;
200 constexpr int NStructInt = PC::NStructInt;
201
202 const int NProcs = ParallelDescriptor::NProcs();
203 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
204
205 {
206 int nrc = PC::ParticleType::is_soa_particle ?
207 pc.NumRealComps() + NStructReal - AMREX_SPACEDIM : pc.NumRealComps() + NStructReal;
208 AMREX_ALWAYS_ASSERT(std::ssize(real_comp_names) == nrc);
209 }
210 AMREX_ALWAYS_ASSERT( int_comp_names.size() == pc.NumIntComps() + NStructInt);
211
212#ifdef AMREX_USE_HDF5_ASYNC
213 // For HDF5 async VOL, block and wait previous tasks have all completed
214 if (es_par_g != 0)
215 async_vol_es_wait_particle();
216 else {
217 ExecOnFinalize(async_vol_es_wait_close_particle);
218 es_par_g = H5EScreate();
219 }
220#endif
221
222 std::string pdir = dir;
223 if ( ! pdir.empty() && pdir[pdir.size()-1] != '/') { pdir += '/'; }
224 pdir += name;
225
226 if ( ! pc.GetLevelDirectoriesCreated()) {
228 {
229 if ( ! amrex::UtilCreateDirectory(pdir, 0755))
230 {
232 }
233 }
235 }
236
237 /* std::ofstream HdrFile; */
238 hid_t fid, grp, fapl, comp_dtype;
239 int status;
240
241 Long nparticles = 0;
242 Long maxnextid;
243
244 comp_dtype = H5Tcreate (H5T_COMPOUND, 2 * AMREX_SPACEDIM * sizeof(int));
245 if (1 == AMREX_SPACEDIM) {
246 H5Tinsert (comp_dtype, "lo_i", 0 * sizeof(int), H5T_NATIVE_INT);
247 H5Tinsert (comp_dtype, "hi_i", 1 * sizeof(int), H5T_NATIVE_INT);
248 }
249 else if (2 == AMREX_SPACEDIM) {
250 H5Tinsert (comp_dtype, "lo_i", 0 * sizeof(int), H5T_NATIVE_INT);
251 H5Tinsert (comp_dtype, "lo_j", 1 * sizeof(int), H5T_NATIVE_INT);
252 H5Tinsert (comp_dtype, "hi_i", 2 * sizeof(int), H5T_NATIVE_INT);
253 H5Tinsert (comp_dtype, "hi_j", 3 * sizeof(int), H5T_NATIVE_INT);
254 }
255 else if (3 == AMREX_SPACEDIM) {
256 H5Tinsert (comp_dtype, "lo_i", 0 * sizeof(int), H5T_NATIVE_INT);
257 H5Tinsert (comp_dtype, "lo_j", 1 * sizeof(int), H5T_NATIVE_INT);
258 H5Tinsert (comp_dtype, "lo_k", 2 * sizeof(int), H5T_NATIVE_INT);
259 H5Tinsert (comp_dtype, "hi_i", 3 * sizeof(int), H5T_NATIVE_INT);
260 H5Tinsert (comp_dtype, "hi_j", 4 * sizeof(int), H5T_NATIVE_INT);
261 H5Tinsert (comp_dtype, "hi_k", 5 * sizeof(int), H5T_NATIVE_INT);
262 }
263
264 // evaluate f for every particle to determine which ones to output
265 Vector<std::map<std::pair<int, int>, typename PC::IntVector > >
266 particle_io_flags(pc.GetParticles().size());
267 for (int lev = 0; lev < std::ssize(pc.GetParticles()); lev++)
268 {
269 const auto& pmap = pc.GetParticles(lev);
270 for (const auto& kv : pmap)
271 {
272 auto& flags = particle_io_flags[lev][kv.first];
273 particle_detail::fillFlags(flags, kv.second, f);
274 }
275 }
276
278
279 if(pc.GetUsePrePost())
280 {
281 nparticles = pc.GetNParticlesPrePost();
282 maxnextid = pc.GetMaxNextIDPrePost();
283 }
284 else
285 {
286 nparticles = particle_detail::countFlags(particle_io_flags, pc);
287 maxnextid = PC::ParticleType::NextID();
288 ParallelDescriptor::ReduceLongSum(nparticles, IOProcNumber);
289 PC::ParticleType::NextID(maxnextid);
290 ParallelDescriptor::ReduceLongMax(maxnextid, IOProcNumber);
291 }
292
293 std::string HDF5FileName = pdir;
294 if ( ! HDF5FileName.empty() && HDF5FileName[HDF5FileName.size()-1] != '/')
295 HDF5FileName += '/';
296
297 HDF5FileName += name;
298 HDF5FileName += ".h5";
299 pc.HdrFileNamePrePost = HDF5FileName;
300
301 BL_PROFILE_VAR("H5FileCreate", h5fg);
302
304 {
305 int set_stripe = 0;
306 char setstripe[1024];
307 int stripe_count = 128;
308 int stripe_size = 1;
309 char *stripe_count_str = getenv("HDF5_STRIPE_COUNT");
310 char *stripe_size_str = getenv("HDF5_STRIPE_SIZE");
311 if (stripe_count_str) {
312 stripe_count = atoi(stripe_count_str);
313 set_stripe = 1;
314 }
315 if (stripe_size_str) {
316 stripe_size = atoi(stripe_size_str);
317 set_stripe = 1;
318 }
319 if (set_stripe == 1) {
320 snprintf(setstripe, sizeof setstripe, "lfs setstripe -c %d -S %dm %s", stripe_count, stripe_size, pdir.c_str());
321 std::cout << "Setting stripe parameters for HDF5 output: " << setstripe << std::endl;
322 amrex::ignore_unused(std::system(setstripe));
323 }
324
325 fid = H5Fcreate(HDF5FileName.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
326 if (fid < 0) { amrex::FileOpenFailed(HDF5FileName.c_str()); }
327
328 //
329 // First thing written is our Checkpoint/Restart version string.
330 // We append "_single" or "_double" to the version string indicating
331 // whether we're using "float" or "double" floating point data in the
332 // particles so that we can Restart from the checkpoint files.
333 //
334 std::string versionName = is_checkpoint ? PC::CheckpointVersion() : PC::PlotfileVersion();
335 if (sizeof(typename PC::ParticleType::RealType) == 4) {
336 versionName += "_single";
337 } else {
338 versionName += "_double";
339 }
340
341 CreateWriteHDF5AttrString(fid, "version_name", versionName.c_str());
342
343 int num_output_real = 0;
344 for (int i = 0; i < std::ssize(write_real_comp); ++i) {
345 if (write_real_comp[i]) { ++num_output_real; }
346 }
347
348 int num_output_int = 0;
349 for (int i = 0; i < pc.NumIntComps() + NStructInt; ++i) {
350 if (write_int_comp[i]) { ++num_output_int; }
351 }
352
353 // AMREX_SPACEDIM and N for sanity checking.
354 int ndim = AMREX_SPACEDIM;
355 grp = H5Gcreate(fid, "Chombo_global", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
356 CreateWriteHDF5Attr(grp, "SpaceDim", 1, &ndim, H5T_NATIVE_INT);
357 H5Gclose(grp);
358
359 // The number of extra real parameters
360 int ncomp = num_output_real + num_output_int;
361 CreateWriteHDF5Attr(fid, "num_component", 1, &ncomp, H5T_NATIVE_INT);
362 CreateWriteHDF5Attr(fid, "num_component_real", 1, &num_output_real, H5T_NATIVE_INT);
363
364 char comp_name[128];
365
366 // Real component names
367 int real_comp_count = AMREX_SPACEDIM; // position values
368
369 for (int i = 0; i < std::ssize(write_real_comp); ++i ) {
370 if (write_real_comp[i]) {
371 /* HdrFile << real_comp_names[i] << '\n'; */
372 snprintf(comp_name, sizeof comp_name, "real_component_%d", real_comp_count);
373 CreateWriteHDF5AttrString(fid, comp_name, real_comp_names[i].c_str());
374 ++real_comp_count;
375 }
376 }
377
378 // The number of extra int parameters
379 CreateWriteHDF5Attr(fid, "num_component_int", 1, &num_output_int, H5T_NATIVE_INT);
380
381 // int component names
382 int int_comp_count = 2;
383
384 for (int i = 0; i < NStructInt + pc.NumIntComps(); ++i ) {
385 if (write_int_comp[i]) {
386 snprintf(comp_name, sizeof comp_name, "int_component_%d", int_comp_count);
387 CreateWriteHDF5AttrString(fid, comp_name, int_comp_names[i].c_str());
388 ++int_comp_count;
389 }
390 }
391
392 // The total number of particles.
393 CreateWriteHDF5Attr(fid, "nparticles", 1, &nparticles, H5T_NATIVE_LLONG);
394
395 // The value of nextid that we need to restore on restart.
396 CreateWriteHDF5Attr(fid, "maxnextid", 1, &maxnextid, H5T_NATIVE_LLONG);
397
398 // Then the finest level of the AMR hierarchy.
399 int finest_level = pc.finestLevel();
400 CreateWriteHDF5Attr(fid, "finest_level", 1, &finest_level, H5T_NATIVE_INT);
401
402 char level_name[128];
403 int ngrids;
404
405 hid_t dcpl_id = H5Pcreate(H5P_DATASET_CREATE);
406 H5Pset_fill_time(dcpl_id, H5D_FILL_TIME_NEVER);
407
408 for (int lev = 0; lev <= finest_level; ++lev) {
409 snprintf(level_name, sizeof level_name, "level_%d", lev);
410
411 grp = H5Gcreate(fid, level_name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
412 if (grp < 0) {
413 std::cout << "H5Gcreate [" << level_name << "] failed!" << std::endl;
414 continue;
415 }
416
417 // Then the number of grids at each level.
418 ngrids = pc.ParticleBoxArray(lev).size();
419 CreateWriteHDF5Attr(grp, "ngrids", 1, &ngrids, H5T_NATIVE_INT);
420
421 /* int ratio = 1; */
422 /* CreateWriteHDF5Attr(grp, "ref_ratio", 1, &ratio); */
423
424 int mfs_size = 2 * AMREX_SPACEDIM;
425 hsize_t mfs_dim = (hsize_t)ngrids;
426
427 hid_t mfs_dset_space = H5Screate_simple(1, &mfs_dim, NULL);
428 hid_t mfs_dset= H5Dcreate(grp, "boxes", comp_dtype, mfs_dset_space, H5P_DEFAULT, dcpl_id, H5P_DEFAULT);
429
430 Vector<int> vbox(ngrids * mfs_size);
431 for(int j = 0; j < pc.ParticleBoxArray(lev).size(); ++j) {
432 for(int i = 0; i < AMREX_SPACEDIM; ++i) {
433 vbox[(j * mfs_size) + i] = pc.ParticleBoxArray(lev)[j].smallEnd(i);
434 vbox[(j * mfs_size) + i + AMREX_SPACEDIM] = pc.ParticleBoxArray(lev)[j].bigEnd(i);
435 }
436 }
437
438 status = H5Dwrite(mfs_dset, comp_dtype, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(vbox[0]));
439 if (status < 0) {
440 std::string msg("ParticleContainer::WriteHDF5ParticleDataSync(): unable to write boxes dataset");
441 amrex::Abort(msg.c_str());
442 }
443
444 H5Sclose(mfs_dset_space);
445 H5Dclose(mfs_dset);
446
447 H5Gclose(grp);
448 }
449 H5Pclose(dcpl_id);
450 H5Fclose(fid);
451 }
452
455
456 // We want to write the data out in parallel.
457 // We'll allow up to nOutFiles active writers at a time.
458 int nOutFiles(-1);
459
460 ParmParse pp("particles");
461 pp.query("particles_nfiles",nOutFiles);
462 if(nOutFiles == -1) nOutFiles = NProcs;
463 /* nOutFiles = std::max(1, std::min(nOutFiles,NProcs)); */
464 pc.nOutFilesPrePost = nOutFiles;
465
466 fapl = H5Pcreate (H5P_FILE_ACCESS);
467#ifdef AMREX_USE_MPI
469#else
470 SetHDF5fapl(fapl);
471#endif
472
473#ifdef AMREX_USE_HDF5_ASYNC
474 fid = H5Fopen_async(HDF5FileName.c_str(), H5F_ACC_RDWR, fapl, es_par_g);
475#else
476 fid = H5Fopen(HDF5FileName.c_str(), H5F_ACC_RDWR, fapl);
477#endif
478 if (fid < 0) {
479 FileOpenFailed(HDF5FileName.c_str());
480 }
481
482 char level_name[64];
483 for (int lev = 0; lev <= pc.finestLevel(); lev++)
484 {
485 snprintf(level_name, sizeof level_name, "level_%d", lev);
486#ifdef AMREX_USE_HDF5_ASYNC
487 grp = H5Gopen_async(fid, level_name, H5P_DEFAULT, es_par_g);
488#else
489 grp = H5Gopen(fid, level_name, H5P_DEFAULT);
490#endif
491
492 bool gotsome;
493 if(pc.usePrePost) {
494 gotsome = (pc.nParticlesAtLevelPrePost[lev] > 0);
495 } else {
496 gotsome = (pc.NumberOfParticlesAtLevel(lev) > 0);
497 }
498
499 MFInfo info;
500 info.SetAlloc(false);
501 MultiFab state(pc.ParticleBoxArray(lev),
502 pc.ParticleDistributionMap(lev),
503 1,0,info);
504
505 // We eventually want to write out the file name and the offset
506 // into that file into which each grid of particles is written.
507 Vector<int> which(state.size(),0);
508 Vector<int > count(state.size(),0);
509 Vector<Long> where(state.size(),0);
510
511 if(pc.usePrePost) {
512 pc.filePrefixPrePost[lev] = HDF5FileName;
513 }
514
515 if (gotsome)
516 {
517 pc.WriteParticlesHDF5(lev, grp, which, count, where,
518 write_real_comp, write_int_comp,
519 compression,
520 particle_io_flags, is_checkpoint);
521
522 if(pc.usePrePost) {
523 pc.whichPrePost[lev] = which;
524 pc.countPrePost[lev] = count;
525 pc.wherePrePost[lev] = where;
526 } else {
527 ParallelDescriptor::ReduceIntSum (which.dataPtr(), which.size(), IOProcNumber);
528 ParallelDescriptor::ReduceIntSum (count.dataPtr(), count.size(), IOProcNumber);
529 ParallelDescriptor::ReduceLongSum(where.dataPtr(), where.size(), IOProcNumber);
530 }
531 }
532
533#ifdef AMREX_USE_HDF5_ASYNC
534 H5Gclose_async(grp, es_par_g);
535#else
536 H5Gclose(grp);
537#endif
538 } // end for (lev...)
539
540 H5Tclose(comp_dtype);
541 H5Pclose(fapl);
542
543#ifdef AMREX_USE_HDF5_ASYNC
544 H5Fclose_async(fid, es_par_g);
545#else
546 H5Fclose(fid);
547#endif
548
549 return;
550} // End WriteHDF5ParticleDataSync
551
552}
553
554#endif /*AMREX_WRITE_BINARY_PARTICLE_DATA_HDF5_H*/
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_PROFILE_VAR_STOP(vname)
Definition AMReX_BLProfiler.H:563
#define BL_PROFILE_VAR(fname, vname)
Definition AMReX_BLProfiler.H:560
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
int size() const noexcept
Return the number of FABs in the FabArray.
Definition AMReX_FabArrayBase.H:110
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:349
int query(std::string_view name, bool &ref, int ival=FIRST) const
Same as querykth() but searches for the last occurrence of name.
Definition AMReX_ParmParse.cpp:1947
This class provides the user with a few print options.
Definition AMReX_Print.H:35
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
T * dataPtr() noexcept
get access to the underlying data pointer
Definition AMReX_Vector.H:49
Long size() const noexcept
Definition AMReX_Vector.H:53
amrex_long Long
Definition AMReX_INT.H:30
void Barrier(const std::string &)
Definition AMReX_ParallelDescriptor.cpp:1215
void ReduceIntSum(int &)
Definition AMReX_ParallelDescriptor.cpp:1265
int NProcs() noexcept
Definition AMReX_ParallelDescriptor.H:255
void ReduceLongSum(Long &)
Definition AMReX_ParallelDescriptor.cpp:1236
int IOProcessorNumber() noexcept
The MPI rank number of the I/O Processor (probably rank 0). This rank is usually used to write to std...
Definition AMReX_ParallelDescriptor.H:279
void ReduceLongMax(Long &)
Definition AMReX_ParallelDescriptor.cpp:1237
bool IOProcessor() noexcept
Is this CPU the I/O Processor? To get the rank number, call IOProcessorNumber()
Definition AMReX_ParallelDescriptor.H:289
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
MPI_Comm Communicator() noexcept
Definition AMReX_ParallelDescriptor.H:223
int MPI_Comm
Definition AMReX_ccse-mpi.H:51
Definition AMReX_Amr.cpp:50
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
void FileOpenFailed(const std::string &file)
Output a message and abort when couldn't open the file.
Definition AMReX_Utility.cpp:137
void CreateDirectoryFailed(const std::string &dir)
Output a message and abort when couldn't create the directory.
Definition AMReX_Utility.cpp:129
bool UtilCreateDirectory(const std::string &path, mode_t mode, bool verbose=false)
Creates the specified directories. path may be either a full pathname or a relative pathname....
Definition AMReX_Utility.cpp:116
static int CreateWriteHDF5AttrString(hid_t loc, const char *name, const char *str)
Definition AMReX_PlotFileUtilHDF5.cpp:79
static int CreateWriteHDF5Attr(hid_t loc, const char *name, hsize_t n, void *data, hid_t dtype)
Create and write a numeric attribute on the given HDF5 object.
Definition AMReX_WriteBinaryParticleDataHDF5.H:25
static void SetHDF5fapl(hid_t fapl, MPI_Comm comm)
Definition AMReX_PlotFileUtilHDF5.cpp:111
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
void ExecOnFinalize(std::function< void()>)
We maintain a stack of functions that need to be called in Finalize(). The functions are called in LI...
Definition AMReX.cpp:331
void WriteHDF5ParticleDataSync(PC const &pc, const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, const std::string &compression, F &&f, bool is_checkpoint)
Synchronous HDF5 particle writer shared by checkpoint and plotfile routines.
Definition AMReX_WriteBinaryParticleDataHDF5.H:184
static int ReadHDF5Attr(hid_t loc, const char *name, void *data, hid_t dtype)
Read an attribute from an HDF5 object into the supplied buffer.
Definition AMReX_WriteBinaryParticleDataHDF5.H:100
FabArray memory allocation information.
Definition AMReX_FabArray.H:67
MFInfo & SetAlloc(bool a) noexcept
Definition AMReX_FabArray.H:74