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