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