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