Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleHDF5.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLEHDF5_H
2#define AMREX_PARTICLEHDF5_H
3#include <AMReX.H>
4#include <AMReX_Config.H>
5
7
8#ifdef AMREX_USE_HDF5
9#include "hdf5.h"
10
11#ifdef AMREX_USE_HDF5_ZFP
12#include "H5Zzfp_lib.h"
13#include "H5Zzfp_props.h"
14#endif
15
16#ifdef AMREX_USE_HDF5_SZ
17#include "H5Z_SZ.h"
18#endif
19
20namespace amrex {
21
22template <typename ParticleType, int NArrayReal, int NArrayInt,
23 template<class> class Allocator, class CellAssignor>
24void
25ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
26::CheckpointHDF5 (const std::string& dir,
27 const std::string& name, bool /*is_checkpoint*/,
28 const Vector<std::string>& real_comp_names,
29 const Vector<std::string>& int_comp_names,
30 const std::string& compression) const
31{
32 Vector<int> write_real_comp;
33 Vector<std::string> tmp_real_comp_names;
34
35 // skip the first spacedim comps for soa particles
36 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
37 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
38 {
39 write_real_comp.push_back(1);
40 if (real_comp_names.empty())
41 {
42 tmp_real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
43 }
44 else
45 {
46 tmp_real_comp_names.push_back(real_comp_names[i-first_rcomp]);
47 }
48 }
49
50 Vector<int> write_int_comp;
51 Vector<std::string> tmp_int_comp_names;
52 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
53 {
54 write_int_comp.push_back(1);
55 if (int_comp_names.empty())
56 {
57 tmp_int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
58 }
59 else
60 {
61 tmp_int_comp_names.push_back(int_comp_names[i]);
62 }
63 }
64
65 WriteHDF5ParticleData(dir, name, write_real_comp, write_int_comp,
66 tmp_real_comp_names, tmp_int_comp_names,
67 compression,
68 [=] AMREX_GPU_HOST_DEVICE (const ConstPTDType& ptd, int i) -> int
69 {
70 return ptd.id(i).is_valid();
71 }, true);
72}
73
74template <typename ParticleType, int NArrayReal, int NArrayInt,
75 template<class> class Allocator, class CellAssignor>
76void
77ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
78::CheckpointHDF5 (const std::string& dir, const std::string& name,
79 const std::string& compression) const
80{
81 Vector<int> write_real_comp;
82 Vector<std::string> real_comp_names;
83 // skip the first spacedim comps for soa particles
84 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
85 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
86 {
87 write_real_comp.push_back(1);
88 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
89 }
90
91 Vector<int> write_int_comp;
92 Vector<std::string> int_comp_names;
93 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
94 {
95 write_int_comp.push_back(1);
96 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
97 }
98
99 WriteHDF5ParticleData(dir, name, write_real_comp, write_int_comp,
100 real_comp_names, int_comp_names,
101 compression,
102 [=] AMREX_GPU_HOST_DEVICE (const ConstPTDType& ptd, int i) -> int
103 {
104 return ptd.id(i).is_valid();
105 });
106}
107
108template <typename ParticleType, int NArrayReal, int NArrayInt,
109 template<class> class Allocator, class CellAssignor>
110void
111ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
112::WritePlotFileHDF5 (const std::string& dir, const std::string& name,
113 const std::string& compression) const
114{
115 Vector<int> write_real_comp;
116 Vector<std::string> real_comp_names;
117 // skip the first spacedim comps for soa particles
118 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
119 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
120 {
121 write_real_comp.push_back(1);
122 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
123 }
124
125 Vector<int> write_int_comp;
126 Vector<std::string> int_comp_names;
127 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
128 {
129 write_int_comp.push_back(1);
130 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
131 }
132
133 WriteHDF5ParticleData(dir, name, write_real_comp, write_int_comp,
134 real_comp_names, int_comp_names,
135 compression,
136 [=] AMREX_GPU_HOST_DEVICE (const ConstPTDType& ptd, int i)
137 {
138 return ptd.id(i).is_valid();
139 });
140}
141
142template <typename ParticleType, int NArrayReal, int NArrayInt,
143 template<class> class Allocator, class CellAssignor>
144void
145ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
146::WritePlotFileHDF5 (const std::string& dir, const std::string& name,
147 const Vector<std::string>& real_comp_names,
148 const Vector<std::string>& int_comp_names,
149 const std::string& compression) const
150{
151 // skip the first spacedim comps for soa particles
152 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
153 AMREX_ALWAYS_ASSERT(static_cast<int>(real_comp_names.size()) == NStructReal + NumRealComps() - first_rcomp);
154 AMREX_ALWAYS_ASSERT(static_cast<int>(int_comp_names.size()) == NStructInt + NumIntComps() );
155
156 Vector<int> write_real_comp;
157 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
158
159 Vector<int> write_int_comp;
160 for (int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
161
162 WriteHDF5ParticleData(dir, name,
163 write_real_comp, write_int_comp,
164 real_comp_names, int_comp_names,
165 compression,
166 [=] AMREX_GPU_HOST_DEVICE (const ConstPTDType& ptd, int i)
167 {
168 return ptd.id(i).is_valid();
169 });
170}
171
172template <typename ParticleType, int NArrayReal, int NArrayInt,
173 template<class> class Allocator, class CellAssignor>
174void
175ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
176::WritePlotFileHDF5 (const std::string& dir, const std::string& name,
177 const Vector<std::string>& real_comp_names,
178 const std::string& compression) const
179{
180 // skip the first spacedim comps for soa particles
181 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
182 AMREX_ALWAYS_ASSERT(static_cast<int>(real_comp_names.size()) == NStructReal + NumRealComps() - first_rcomp);
183
184 Vector<int> write_real_comp;
185 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) write_real_comp.push_back(1);
186
187 Vector<int> write_int_comp;
188 for (int i = 0; i < NStructInt + NumIntComps(); ++i) write_int_comp.push_back(1);
189
190 Vector<std::string> int_comp_names;
191 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
192 {
193 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
194 }
195
196 WriteHDF5ParticleData(dir, name,
197 write_real_comp, write_int_comp,
198 real_comp_names, int_comp_names,
199 compression,
200 [=] AMREX_GPU_HOST_DEVICE (const ConstPTDType& ptd, int i)
201 {
202 return ptd.id(i).is_valid();
203 });
204}
205
206template <typename ParticleType, int NArrayReal, int NArrayInt,
207 template<class> class Allocator, class CellAssignor>
208void
209ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
210::WritePlotFileHDF5 (const std::string& dir,
211 const std::string& name,
212 const Vector<int>& write_real_comp,
213 const Vector<int>& write_int_comp,
214 const std::string& compression) const
215{
216 // skip the first spacedim comps for soa particles
217 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
218 AMREX_ALWAYS_ASSERT(static_cast<int>(write_real_comp.size()) == NStructReal + NumRealComps() - first_rcomp);
219 AMREX_ALWAYS_ASSERT(static_cast<int>(write_int_comp.size()) == NStructInt + NumIntComps() );
220
221 Vector<std::string> real_comp_names;
222 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
223 {
224 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
225 }
226
227 Vector<std::string> int_comp_names;
228 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
229 {
230 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
231 }
232
233 WriteHDF5ParticleData(dir, name, write_real_comp, write_int_comp,
234 real_comp_names, int_comp_names,
235 compression,
236 [=] AMREX_GPU_HOST_DEVICE (const ConstPTDType& ptd, int i)
237 {
238 return ptd.id(i).is_valid();
239 });
240}
241
242template <typename ParticleType, int NArrayReal, int NArrayInt,
243 template<class> class Allocator, class CellAssignor>
244void
245ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>::
246WritePlotFileHDF5 (const std::string& dir, const std::string& name,
247 const Vector<int>& write_real_comp,
248 const Vector<int>& write_int_comp,
249 const Vector<std::string>& real_comp_names,
250 const Vector<std::string>& int_comp_names,
251 const std::string& compression) const
252{
253 BL_PROFILE("ParticleContainer::WritePlotFile()");
254
255 WriteHDF5ParticleData(dir, name,
256 write_real_comp, write_int_comp,
257 real_comp_names, int_comp_names,
258 compression,
259 [=] AMREX_GPU_HOST_DEVICE (const ConstPTDType& ptd, int i)
260 {
261 return ptd.id(i).is_valid();
262 });
263}
264
265template <typename ParticleType, int NArrayReal, int NArrayInt,
266 template<class> class Allocator, class CellAssignor>
267template <class F, std::enable_if_t<!std::is_same_v<F, Vector<std::string>>>*>
268void
269ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
270::WritePlotFileHDF5 (const std::string& dir, const std::string& name,
271 const std::string& compression, F&& f) const
272{
273 Vector<int> write_real_comp;
274 Vector<std::string> real_comp_names;
275 // skip the first spacedim comps for soa particles
276 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
277 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
278 {
279 write_real_comp.push_back(1);
280 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
281 }
282
283 Vector<int> write_int_comp;
284 Vector<std::string> int_comp_names;
285 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
286 {
287 write_int_comp.push_back(1);
288 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
289 }
290
291 WriteHDF5ParticleData(dir, name, write_real_comp, write_int_comp,
292 real_comp_names, int_comp_names, compression,
293 std::forward<F>(f));
294}
295
296template <typename ParticleType, int NArrayReal, int NArrayInt,
297 template<class> class Allocator, class CellAssignor>
298template <class F>
299void
300ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
301::WritePlotFileHDF5 (const std::string& dir, const std::string& name,
302 const Vector<std::string>& real_comp_names,
303 const Vector<std::string>& int_comp_names,
304 const std::string& compression, F&& f) const
305{
306 // skip the first spacedim comps for soa particles
307 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
308 AMREX_ALWAYS_ASSERT(static_cast<int>(real_comp_names.size()) == NStructReal + NumRealComps() - first_rcomp);
309 AMREX_ALWAYS_ASSERT(static_cast<int>(int_comp_names.size()) == NStructInt + NumIntComps() );
310
311 Vector<int> write_real_comp;
312 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
313
314 Vector<int> write_int_comp;
315 for (int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
316
317 WriteHDF5ParticleData(dir, name,
318 write_real_comp, write_int_comp,
319 real_comp_names, int_comp_names,
320 compression, std::forward<F>(f));
321}
322
323template <typename ParticleType, int NArrayReal, int NArrayInt,
324 template<class> class Allocator, class CellAssignor>
325template <class F, std::enable_if_t<!std::is_same_v<F, Vector<std::string>>>*>
326void
327ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
328::WritePlotFileHDF5 (const std::string& dir, const std::string& name,
329 const Vector<std::string>& real_comp_names,
330 const std::string& compression, F&& f) const
331{
332 // skip the first spacedim comps for soa particles
333 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
334 AMREX_ALWAYS_ASSERT(static_cast<int>(real_comp_names.size()) == NStructReal + NumRealComps() - first_rcomp);
335
336 Vector<int> write_real_comp;
337 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
338
339 Vector<int> write_int_comp;
340 for (int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
341
342 Vector<std::string> int_comp_names;
343 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
344 {
345 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
346 }
347
348 WriteHDF5ParticleData(dir, name,
349 write_real_comp, write_int_comp,
350 real_comp_names, int_comp_names,
351 compression, std::forward<F>(f));
352}
353
354template <typename ParticleType, int NArrayReal, int NArrayInt,
355 template<class> class Allocator, class CellAssignor>
356template <class F>
357void
358ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
359::WritePlotFileHDF5 (const std::string& dir,
360 const std::string& name,
361 const Vector<int>& write_real_comp,
362 const Vector<int>& write_int_comp,
363 const std::string& compression, F&& f) const
364{
365 // skip the first spacedim comps for soa particles
366 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
367 AMREX_ALWAYS_ASSERT(static_cast<int>(write_real_comp.size()) == NStructReal + NumRealComps() - first_rcomp);
368 AMREX_ALWAYS_ASSERT(static_cast<int>(write_int_comp.size()) == NStructInt + NumIntComps() );
369
370 Vector<std::string> real_comp_names;
371 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
372 {
373 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
374 }
375
376 Vector<std::string> int_comp_names;
377 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
378 {
379 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
380 }
381
382 WriteHDF5ParticleData(dir, name, write_real_comp, write_int_comp,
383 real_comp_names, int_comp_names,
384 compression, std::forward<F>(f));
385}
386
387template <typename ParticleType, int NArrayReal, int NArrayInt,
388 template<class> class Allocator, class CellAssignor>
389template <class F>
390void
391ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>::
392WritePlotFileHDF5 (const std::string& dir, const std::string& name,
393 const Vector<int>& write_real_comp,
394 const Vector<int>& write_int_comp,
395 const Vector<std::string>& real_comp_names,
396 const Vector<std::string>& int_comp_names,
397 const std::string& compression, F&& f) const
398{
399 BL_PROFILE("ParticleContainer::WritePlotFile()");
400
401 WriteHDF5ParticleData(dir, name,
402 write_real_comp, write_int_comp,
403 real_comp_names, int_comp_names,
404 compression, std::forward<F>(f));
405}
406
407template <typename ParticleType, int NArrayReal, int NArrayInt,
408 template<class> class Allocator, class CellAssignor>
409template <class F>
410void
411ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
412::WriteHDF5ParticleData (const std::string& dir, const std::string& name,
413 const Vector<int>& write_real_comp,
414 const Vector<int>& write_int_comp,
415 const Vector<std::string>& real_comp_names,
416 const Vector<std::string>& int_comp_names,
417 const std::string& compression,
418 F&& f, bool is_checkpoint) const
419{
420 /* HDF5 async is implemented in WriteHDF5ParticleDataSync, enabled by compile flag */
421 /* if (AsyncOut::UseAsyncOut()) { */
422 /* WriteHDF5ParticleDataAsync(*this, dir, name, */
423 /* write_real_comp, write_int_comp, */
424 /* real_comp_names, int_comp_names); */
425 /* } else */
426 /* { */
427 WriteHDF5ParticleDataSync(*this, dir, name,
428 write_real_comp, write_int_comp,
429 real_comp_names, int_comp_names,
430 compression,
431 std::forward<F>(f), is_checkpoint);
432 /* } */
433}
434
435template <typename ParticleType, int NArrayReal, int NArrayInt,
436 template<class> class Allocator, class CellAssignor>
437void
438ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
439::CheckpointPreHDF5 ()
440{
441 if( ! usePrePost) {
442 return;
443 }
444
445 BL_PROFILE("ParticleContainer::CheckpointPre()");
446
447 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
448 Long nparticles = 0;
449 Long maxnextid = ParticleType::NextID();
450
451 for (int lev = 0; lev < m_particles.size(); lev++) {
452 const auto& pmap = m_particles[lev];
453 for (const auto& kv : pmap) {
454 const auto& aos = kv.second.GetArrayOfStructs();
455 for (int k = 0; k < aos.numParticles(); ++k) {
456 const ParticleType& p = aos[k];
457 if (p.id().is_valid()) {
458 //
459 // Only count (and checkpoint) valid particles.
460 //
461 nparticles++;
462 }
463 }
464 }
465 }
466 ParallelDescriptor::ReduceLongSum(nparticles, IOProcNumber);
467
468 ParticleType::NextID(maxnextid);
469 ParallelDescriptor::ReduceLongMax(maxnextid, IOProcNumber);
470
471 nparticlesPrePost = nparticles;
472 maxnextidPrePost = maxnextid;
473
474 nParticlesAtLevelPrePost.clear();
475 nParticlesAtLevelPrePost.resize(finestLevel() + 1, 0);
476 for(int lev(0); lev <= finestLevel(); ++lev) {
477 nParticlesAtLevelPrePost[lev] = NumberOfParticlesAtLevel(lev);
478 }
479
480 whichPrePost.clear();
481 whichPrePost.resize(finestLevel() + 1);
482 countPrePost.clear();
483 countPrePost.resize(finestLevel() + 1);
484 wherePrePost.clear();
485 wherePrePost.resize(finestLevel() + 1);
486
487 filePrefixPrePost.clear();
488 filePrefixPrePost.resize(finestLevel() + 1);
489}
490
491
492template <typename ParticleType, int NArrayReal, int NArrayInt,
493 template<class> class Allocator, class CellAssignor>
494void
495ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
496::CheckpointPostHDF5 ()
497{
498 if( ! usePrePost) {
499 return;
500 }
501
502 BL_PROFILE("ParticleContainer::CheckpointPostHDF5()");
503
504 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
505 std::ofstream HdrFile;
506 HdrFile.open(HdrFileNamePrePost.c_str(), std::ios::out | std::ios::app);
507
508 for(int lev(0); lev <= finestLevel(); ++lev) {
509 ParallelDescriptor::ReduceIntSum (whichPrePost[lev].dataPtr(), whichPrePost[lev].size(), IOProcNumber);
510 ParallelDescriptor::ReduceIntSum (countPrePost[lev].dataPtr(), countPrePost[lev].size(), IOProcNumber);
511 ParallelDescriptor::ReduceLongSum(wherePrePost[lev].dataPtr(), wherePrePost[lev].size(), IOProcNumber);
512
513
515 for(int j(0); j < whichPrePost[lev].size(); ++j) {
516 HdrFile << whichPrePost[lev][j] << ' ' << countPrePost[lev][j] << ' ' << wherePrePost[lev][j] << '\n';
517 }
518
519 const bool gotsome = (nParticlesAtLevelPrePost[lev] > 0);
520 if(gotsome && doUnlink) {
521// BL_PROFILE_VAR("PC<NNNN>::Checkpoint:unlink", unlink_post);
522 // Unlink any zero-length data files.
523 Vector<Long> cnt(nOutFilesPrePost,0);
524
525 for(int i(0), N = countPrePost[lev].size(); i < N; ++i) {
526 cnt[whichPrePost[lev][i]] += countPrePost[lev][i];
527 }
528
529 for(int i(0), N = cnt.size(); i < N; ++i) {
530 if(cnt[i] == 0) {
531 std::string FullFileName = NFilesIter::FileName(i, filePrefixPrePost[lev]);
532 FileSystem::Remove(FullFileName);
533 }
534 }
535 }
536 }
537 }
538
540 HdrFile.flush();
541 HdrFile.close();
542 if( ! HdrFile.good()) {
543 amrex::Abort("ParticleContainer::CheckpointPostHDF5(): problem writing HdrFile");
544 }
545 }
546}
547
548template <typename ParticleType, int NArrayReal, int NArrayInt,
549 template<class> class Allocator, class CellAssignor>
550void
551ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
552::WritePlotFilePreHDF5 ()
553{
555}
556
557
558template <typename ParticleType, int NArrayReal, int NArrayInt,
559 template<class> class Allocator, class CellAssignor>
560void
561ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
562::WritePlotFilePostHDF5 ()
563{
565}
566
567
568template <typename ParticleType, int NArrayReal, int NArrayInt,
569 template<class> class Allocator, class CellAssignor>
570void
571ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
572::WriteParticlesHDF5 (int lev, hid_t grp,
573 Vector<int>& /*which*/, Vector<int>& count, Vector<Long>& /*where*/,
574 const Vector<int>& write_real_comp,
575 const Vector<int>& write_int_comp,
576 const std::string& compression,
577 const Vector<std::map<std::pair<int, int>, IntVector>>& particle_io_flags,
578 bool is_checkpoint) const
579{
580 BL_PROFILE("ParticleContainer::WriteParticlesHDF5()");
581
582 // For a each grid, the tiles it contains
583 std::map<int, Vector<int> > tile_map;
584
585 int ret;
586 hid_t dxpl_col, dxpl_ind, dcpl_int, dcpl_real;
587 dxpl_col = H5Pcreate(H5P_DATASET_XFER);
588 dxpl_ind = H5Pcreate(H5P_DATASET_XFER);
589#ifdef AMREX_USE_MPI
590 H5Pset_dxpl_mpio(dxpl_col, H5FD_MPIO_COLLECTIVE);
591#endif
592 dcpl_int = H5Pcreate(H5P_DATASET_CREATE);
593 dcpl_real = H5Pcreate(H5P_DATASET_CREATE);
594
595 std::string mode_env, value_env;
596 double comp_value = -1.0;
597 std::string::size_type pos = compression.find('@');
598 if (pos != std::string::npos) {
599 mode_env = compression.substr(0, pos);
600 value_env = compression.substr(pos+1);
601 if (!value_env.empty()) {
602 comp_value = atof(value_env.c_str());
603 }
604 }
605
606 H5Pset_alloc_time(dcpl_int, H5D_ALLOC_TIME_LATE);
607 H5Pset_alloc_time(dcpl_real, H5D_ALLOC_TIME_LATE);
608
609 if (!mode_env.empty() && mode_env != "None") {
610 const char *chunk_env = NULL;
611 hsize_t chunk_dim = 1024;
612 chunk_env = getenv("HDF5_CHUNK_SIZE");
613 if (chunk_env != NULL) {
614 chunk_dim = atoi(chunk_env);
615 }
616
617 H5Pset_chunk(dcpl_int, 1, &chunk_dim);
618 H5Pset_chunk(dcpl_real, 1, &chunk_dim);
619
620#ifdef AMREX_USE_HDF5_ZFP
621 pos = mode_env.find("ZFP");
622 if (pos != std::string::npos) {
623 ret = H5Z_zfp_initialize();
624 if (ret < 0) { amrex::Abort("ZFP initialize failed!"); }
625 }
626#endif
627
628 if (mode_env == "ZLIB") {
629 H5Pset_shuffle(dcpl_int);
630 H5Pset_shuffle(dcpl_real);
631 H5Pset_deflate(dcpl_int, (int)comp_value);
632 H5Pset_deflate(dcpl_real, (int)comp_value);
633 }
634#ifdef AMREX_USE_HDF5_SZ
635 else if (mode_env == "SZ") {
636 ret = H5Z_SZ_Init((char*)value_env.c_str());
637 if (ret < 0) {
638 std::cout << "SZ config file:" << value_env.c_str() << std::endl;
639 amrex::Abort("SZ initialize failed, check SZ config file!");
640 }
641 }
642#endif
643#ifdef AMREX_USE_HDF5_ZFP
644 else if (mode_env == "ZFP_RATE") {
645 H5Pset_zfp_rate(dcpl_int, comp_value);
646 H5Pset_zfp_rate(dcpl_real, comp_value);
647 }
648 else if (mode_env == "ZFP_PRECISION") {
649 H5Pset_zfp_precision(dcpl_int, (unsigned int)comp_value);
650 H5Pset_zfp_precision(dcpl_real, (unsigned int)comp_value);
651 }
652 else if (mode_env == "ZFP_ACCURACY") {
653 H5Pset_zfp_accuracy(dcpl_int, comp_value);
654 H5Pset_zfp_accuracy(dcpl_real, comp_value);
655 }
656 else if (mode_env == "ZFP_REVERSIBLE") {
657 H5Pset_zfp_reversible(dcpl_int);
658 H5Pset_zfp_reversible(dcpl_real);
659 }
660#endif
661 if (ParallelDescriptor::MyProc() == 0) {
662 std::cout << "\nHDF5 particle checkpoint using " << mode_env << ", " <<
663 value_env << ", " << chunk_dim << std::endl;
664 }
665 }
666
667 for (const auto& kv : m_particles[lev])
668 {
669 const int grid = kv.first.first;
670 const int tile = kv.first.second;
671 tile_map[grid].push_back(tile);
672 const auto& pflags = particle_io_flags[lev].at(kv.first);
673
674 // Only write out valid particles.
675 count[grid] += particle_detail::countFlags(pflags);
676 }
677
678 MFInfo info;
679 info.SetAlloc(false);
680 MultiFab state(ParticleBoxArray(lev), ParticleDistributionMap(lev), 1,0,info);
681
682 int my_mfi_cnt = 0;
683 ULong my_mfi_int_total_size = 0, my_mfi_real_total_size = 0, int_size, real_size;
684 Vector<int> all_mfi_cnt(ParallelDescriptor::NProcs());
685 Vector<int> my_mfi_real_size;
686 Vector<int> my_mfi_int_size;
687 Vector<int> my_nparticles;
688 Vector<ULong> all_mfi_real_total_size(ParallelDescriptor::NProcs());
689 Vector<ULong> all_mfi_int_total_size(ParallelDescriptor::NProcs());
690 hid_t offset_id, offset_space, real_mem_space, real_dset_space, real_dset_id;
691 hid_t int_mem_space, int_dset_id, int_dset_space;
692 hsize_t total_mfi = 0, total_real_size = 0, total_int_size = 0, real_file_offset = 0, int_file_offset = 0;
693 hsize_t my_int_offset, my_int_count, my_real_offset, my_real_count;
694
695 // Count total number of components written (positions always included;
696 // write_real_comp covers the non-position components for pure SoA)
697 int real_comp_count = AMREX_SPACEDIM; // position values
698
699 for (int i = 0; i < static_cast<int>(write_real_comp.size()); ++i ) {
700 if (write_real_comp[i]) { ++real_comp_count; }
701 }
702
703 int int_comp_count = 2; // cpu and id values
704
705 for (int i = 0; i < NStructInt + NumIntComps(); ++i ) {
706 if (write_int_comp[i]) { ++int_comp_count; }
707 }
708
709 // Get the size for each mf so we know the amount of data from each rank
710 for (MFIter mfi(state); mfi.isValid(); ++mfi) {
711 const int grid = mfi.index();
712 if (count[grid] == 0)
713 continue;
714
715 int_size = count[grid] * int_comp_count;
716 my_mfi_int_size.push_back(int_size);
717 my_nparticles.push_back(count[grid]);
718 my_mfi_int_total_size += int_size;
719
720
721 real_size = count[grid] * real_comp_count;
722 my_mfi_real_size.push_back(real_size);
723 my_mfi_real_total_size += real_size;
724 my_mfi_cnt++;
725 }
726
727 #ifdef AMREX_USE_MPI
728 // Collect the number of mf and total size of mf from each rank
729 MPI_Allgather(&my_mfi_cnt, 1, ParallelDescriptor::Mpi_typemap<int>::type(), &(all_mfi_cnt[0]), 1,
731
732 for (int i = 0; i < ParallelDescriptor::NProcs(); i++)
733 total_mfi += all_mfi_cnt[i];
734
735 // Create the int data
736 MPI_Allgather(&my_mfi_int_total_size, 1, ParallelDescriptor::Mpi_typemap<ULong>::type(),
738 #else
739 all_mfi_cnt[0] = my_mfi_cnt;
740 all_mfi_int_total_size[0] = my_mfi_int_total_size;
741 #endif
742
743 int_file_offset = 0;
744 for (int i = 0; i < ParallelDescriptor::MyProc(); i++)
745 int_file_offset += all_mfi_int_total_size[i];
746 my_int_offset = int_file_offset;
747 my_int_count = 0;
748
749 for (int i = 0; i < ParallelDescriptor::NProcs(); i++)
750 total_int_size += all_mfi_int_total_size[i];
751
752 // SZ int compression seems to have issues at the moment
753/* #ifdef AMREX_USE_HDF5_SZ */
754/* if (mode_env == "SZ") { */
755/* size_t cd_nelmts; */
756/* unsigned int* cd_values = NULL; */
757/* unsigned filter_config; */
758/* SZ_metaDataToCdArray(&cd_nelmts, &cd_values, SZ_INT32, 0, 0, 0, 0, total_int_size); */
759/* H5Pset_filter(dcpl_int, H5Z_FILTER_SZ, H5Z_FLAG_MANDATORY, cd_nelmts, cd_values); */
760/* } */
761/* #endif */
762
763 hsize_t chunk_size;
764 if (H5Pget_layout(dcpl_int) == H5D_CHUNKED) {
765 if (H5Pget_chunk(dcpl_int, 1, &chunk_size) > -1) {
766 if (chunk_size > total_int_size) {
767 H5Pset_chunk(dcpl_int, 1, &total_int_size);
768 }
769 }
770 }
771
772 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
773#ifdef AMREX_USE_HDF5_ASYNC
774 int_dset_id = H5Dcreate_async(grp, "data:datatype=0", H5T_NATIVE_INT, int_dset_space, H5P_DEFAULT, dcpl_int, H5P_DEFAULT, es_par_g);
775#else
776 int_dset_id = H5Dcreate(grp, "data:datatype=0", H5T_NATIVE_INT, int_dset_space, H5P_DEFAULT, dcpl_int, H5P_DEFAULT);
777#endif
778
779 H5Sclose(int_dset_space);
780
781 // Create the real data
782 #ifdef AMREX_USE_MPI
783 MPI_Allgather(&my_mfi_real_total_size, 1, ParallelDescriptor::Mpi_typemap<ULong>::type(),
785 #else
786 all_mfi_real_total_size[0] = my_mfi_real_total_size;
787 #endif
788
789 for (int i = 0; i < ParallelDescriptor::NProcs(); i++)
790 total_real_size += all_mfi_real_total_size[i];
791
792#ifdef AMREX_USE_HDF5_SZ
793 if (mode_env == "SZ") {
794 size_t cd_nelmts;
795 unsigned int* cd_values = NULL;
796 unsigned filter_config;
797 if (sizeof(typename ParticleType::RealType) == 4) {
798 SZ_metaDataToCdArray(&cd_nelmts, &cd_values, SZ_FLOAT, 0, 0, 0, 0, total_real_size);
799 H5Pset_filter(dcpl_real, H5Z_FILTER_SZ, H5Z_FLAG_MANDATORY, cd_nelmts, cd_values);
800 }
801 else {
802 SZ_metaDataToCdArray(&cd_nelmts, &cd_values, SZ_DOUBLE, 0, 0, 0, 0, total_real_size);
803 H5Pset_filter(dcpl_real, H5Z_FILTER_SZ, H5Z_FLAG_MANDATORY, cd_nelmts, cd_values);
804 }
805 }
806#endif
807
808 if (H5Pget_layout(dcpl_real) == H5D_CHUNKED) {
809 if (H5Pget_chunk(dcpl_real, 1, &chunk_size) > -1) {
810 if (chunk_size > total_real_size) {
811 H5Pset_chunk(dcpl_real, 1, &total_real_size);
812 }
813 }
814 }
815
816 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
817 if (sizeof(typename ParticleType::RealType) == 4) {
818#ifdef AMREX_USE_HDF5_ASYNC
819 real_dset_id = H5Dcreate_async(grp, "data:datatype=1", H5T_NATIVE_FLOAT, real_dset_space,
820 H5P_DEFAULT, dcpl_real, H5P_DEFAULT, es_par_g);
821#else
822 real_dset_id = H5Dcreate(grp, "data:datatype=1", H5T_NATIVE_FLOAT, real_dset_space,
823 H5P_DEFAULT, dcpl_real, H5P_DEFAULT);
824#endif
825 }
826 else {
827#ifdef AMREX_USE_HDF5_ASYNC
828 real_dset_id = H5Dcreate_async(grp, "data:datatype=1", H5T_NATIVE_DOUBLE, real_dset_space,
829 H5P_DEFAULT, dcpl_real, H5P_DEFAULT, es_par_g);
830#else
831 real_dset_id = H5Dcreate(grp, "data:datatype=1", H5T_NATIVE_DOUBLE, real_dset_space,
832 H5P_DEFAULT, dcpl_real, H5P_DEFAULT);
833#endif
834 }
835 H5Sclose(real_dset_space);
836
837 real_file_offset = 0;
838 for (int i = 0; i < ParallelDescriptor::MyProc(); i++) {
839 real_file_offset += all_mfi_real_total_size[i];
840 }
841 my_real_offset = real_file_offset;
842 my_real_count = 0;
843
844 int max_mfi_count = 0, write_count = 0;
845 for (int i = 0; i < ParallelDescriptor::NProcs(); i++) {
846 if (max_mfi_count < all_mfi_cnt[i]) {
847 max_mfi_count = all_mfi_cnt[i];
848 }
849 }
850
851
852 for (MFIter mfi(state); mfi.isValid(); ++mfi)
853 {
854 const int grid = mfi.index();
855
856 if (count[grid] == 0) { continue; }
857
858 Vector<int> istuff;
859 Vector<ParticleReal> rstuff;
860 particle_detail::packIOData(istuff, rstuff, *this, lev, grid,
861 write_real_comp, write_int_comp,
862 particle_io_flags, tile_map[grid], count[grid],
863 is_checkpoint);
864
865 my_int_count = istuff.size();
866 int_mem_space = H5Screate_simple(1, &my_int_count, NULL);
867 /* std::cout << "Rank " << ParallelDescriptor::MyProc() << ": my_int_offset = " << */
868 /* my_int_offset << ", my_int_count = " << my_int_count << ", total_int_size = " << total_int_size << std::endl; */
869 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
870 H5Sselect_hyperslab (int_dset_space, H5S_SELECT_SET, &my_int_offset, NULL, &my_int_count, NULL);
871
872#ifdef AMREX_USE_HDF5_ASYNC
873 ret = H5Dwrite_async(int_dset_id, H5T_NATIVE_INT, int_mem_space, int_dset_space, dxpl_col, istuff.dataPtr(), es_par_g);
874#else
875 ret = H5Dwrite(int_dset_id, H5T_NATIVE_INT, int_mem_space, int_dset_space, dxpl_col, istuff.dataPtr());
876#endif
877 if (ret < 0) { amrex::Abort("H5Dwrite int_dset failed!"); }
878
879 H5Sclose(int_dset_space);
880 H5Sclose(int_mem_space);
881
882 my_int_offset += my_int_count;
883
884 my_real_count = rstuff.size();
885 real_mem_space = H5Screate_simple(1, &my_real_count, NULL);
886 /* std::cout << "Rank " << ParallelDescriptor::MyProc() << ": my_real_offset = " << */
887 /* my_real_offset << ", my_real_count = " << my_real_count << ", total_real_size = " << total_real_size << '\n'; */
888 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
889 H5Sselect_hyperslab (real_dset_space, H5S_SELECT_SET, &my_real_offset, NULL, &my_real_count, NULL);
890 if (sizeof(typename ParticleType::RealType) == 4) {
891#ifdef AMREX_USE_HDF5_ASYNC
892 ret = H5Dwrite_async(real_dset_id, H5T_NATIVE_FLOAT, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr(), es_par_g);
893#else
894 ret = H5Dwrite(real_dset_id, H5T_NATIVE_FLOAT, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr());
895#endif
896 } else {
897#ifdef AMREX_USE_HDF5_ASYNC
898 ret = H5Dwrite_async(real_dset_id, H5T_NATIVE_DOUBLE, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr(), es_par_g);
899#else
900 ret = H5Dwrite(real_dset_id, H5T_NATIVE_DOUBLE, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr());
901#endif
902 }
903
904 if (ret < 0) { amrex::Abort("H5Dwrite real_dset failed!"); }
905
906 H5Sclose(real_mem_space);
907 H5Sclose(real_dset_space);
908
909 my_real_offset += my_real_count;
910 write_count++;
911
912 } // end for (mfi)
913
914 // Dummy writes so that every rank participates to every possible H5Dwrite (collective)
915 while (write_count < max_mfi_count) {
916 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
917 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
918 H5Sselect_none(int_dset_space);
919 H5Sselect_none(real_dset_space);
920
921#ifdef AMREX_USE_HDF5_ASYNC
922 H5Dwrite_async(int_dset_id, H5T_NATIVE_INT, int_dset_space, int_dset_space, dxpl_col, NULL, es_par_g);
923 if (sizeof(typename ParticleType::RealType) == 4) {
924 H5Dwrite_async(real_dset_id, H5T_NATIVE_FLOAT, real_dset_space, real_dset_space, dxpl_col, NULL, es_par_g);
925 } else {
926 H5Dwrite_async(real_dset_id, H5T_NATIVE_DOUBLE, real_dset_space, real_dset_space, dxpl_col, NULL, es_par_g);
927 }
928#else
929 H5Dwrite(int_dset_id, H5T_NATIVE_INT, int_dset_space, int_dset_space, dxpl_col, NULL);
930 if (sizeof(typename ParticleType::RealType) == 4) {
931 H5Dwrite(real_dset_id, H5T_NATIVE_FLOAT, real_dset_space, real_dset_space, dxpl_col, NULL);
932 } else {
933 H5Dwrite(real_dset_id, H5T_NATIVE_DOUBLE, real_dset_space, real_dset_space, dxpl_col, NULL);
934 }
935#endif
936
937 H5Sclose(int_dset_space);
938 H5Sclose(real_dset_space);
939
940 write_count++;
941 }
942
943#ifdef AMREX_USE_HDF5_ASYNC
944 H5Dclose_async(real_dset_id, es_par_g);
945 H5Dclose_async(int_dset_id, es_par_g);
946#else
947 H5Dclose(real_dset_id);
948 H5Dclose(int_dset_id);
949#endif
950
951 // Create and write the size dataset
952 offset_space = H5Screate_simple(1, &total_mfi, NULL);
953#ifdef AMREX_USE_HDF5_ASYNC
954 offset_id = H5Dcreate_async(grp, "nparticles_grid", H5T_NATIVE_INT, offset_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT, es_par_g);
955#else
956 offset_id = H5Dcreate(grp, "nparticles_grid", H5T_NATIVE_INT, offset_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
957#endif
958
959 my_int_offset = 0;
960 for (int i = 0; i < ParallelDescriptor::MyProc(); i++) {
961 my_int_offset += all_mfi_cnt[i];
962 }
963 my_int_count = my_mfi_cnt;
964 int_mem_space = H5Screate_simple(1, &my_int_count, NULL);
965 /* std::cout << "Rank " << ParallelDescriptor::MyProc() << ": my_int_offset = " << */
966 /* my_int_offset << ", my_int_count = " << my_int_count << ", total_mfi = " << total_mfi << '\n'; */
967 H5Sselect_hyperslab (offset_space, H5S_SELECT_SET, &my_int_offset, NULL, &my_int_count, NULL);
968
969#ifdef AMREX_USE_HDF5_ASYNC
970 ret = H5Dwrite_async(offset_id, H5T_NATIVE_INT, int_mem_space, offset_space, dxpl_col, my_nparticles.dataPtr(), es_par_g);
971#else
972 ret = H5Dwrite(offset_id, H5T_NATIVE_INT, int_mem_space, offset_space, dxpl_col, my_nparticles.dataPtr());
973#endif
974 if (ret < 0) { amrex::Abort("H5Dwrite offset failed!"); }
975
976 H5Pclose(dcpl_int);
977 H5Pclose(dcpl_real);
978 H5Pclose(dxpl_col);
979 H5Pclose(dxpl_ind);
980
981 H5Sclose(int_mem_space);
982 H5Sclose(offset_space);
983
984#ifdef AMREX_USE_HDF5_ASYNC
985 H5Dclose_async(offset_id, es_par_g);
986#else
987 H5Dclose(offset_id);
988#endif
989
990 /* std::cout << "Rank " << ParallelDescriptor::MyProc() << ": done WriteParticlesHDF5" << std::endl; */
991 return;
992} // End WriteParticlesHDF5
993
994template <typename ParticleType, int NArrayReal, int NArrayInt,
995 template<class> class Allocator, class CellAssignor>
996void
997ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
998::RestartHDF5 (const std::string& dir, const std::string& file, bool /*is_checkpoint*/)
999{
1000 RestartHDF5(dir, file);
1001}
1002
1003template <typename ParticleType, int NArrayReal, int NArrayInt,
1004 template<class> class Allocator, class CellAssignor>
1005void
1006ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1007::RestartHDF5 (const std::string& dir, const std::string& file)
1008{
1009 BL_PROFILE("ParticleContainer::RestartHDF5()");
1010 AMREX_ASSERT(!dir.empty());
1011 AMREX_ASSERT(!file.empty());
1012
1013 const auto strttime = amrex::second();
1014
1015 std::string fullname = dir;
1016 if (!fullname.empty() && fullname[fullname.size()-1] != '/') {
1017 fullname += '/';
1018 }
1019
1020 fullname += file;
1021 fullname += ".h5";
1022
1023 hid_t fid, dset, grp, fapl, attr, atype, dspace, int_dset, real_dset;
1024 int ret;
1025
1026 fapl = H5Pcreate (H5P_FILE_ACCESS);
1027#ifdef AMREX_USE_MPI
1029#else
1030 SetHDF5fapl(fapl);
1031#endif
1032
1033 fid = H5Fopen(fullname.c_str(), H5F_ACC_RDONLY, fapl);
1034 if (fid < 0) {
1035 std::string msg("ParticleContainer::RestartHDF5(): unable to open file: ");
1036 msg += fullname;
1037 amrex::Abort(msg.c_str());
1038 }
1039
1040 attr = H5Aopen(fid, "version_name", H5P_DEFAULT);
1041 if (attr < 0) {
1042 std::string msg("ParticleContainer::RestartHDF5(): unable to open version attribute ");
1043 amrex::Abort(msg.c_str());
1044 }
1045
1046 atype = H5Aget_type(attr);
1047 if (atype < 0) {
1048 std::string msg("ParticleContainer::RestartHDF5(): unable to get type of attribute ");
1049 amrex::Abort(msg.c_str());
1050 }
1051
1052 size_t attr_len = H5Tget_size(atype);
1053 if (attr_len == 0) {
1054 std::string msg("ParticleContainer::RestartHDF5(): unable to get size of attribute ");
1055 amrex::Abort(msg.c_str());
1056 }
1057
1058 std::string version;
1059 version.resize(attr_len+1);
1060
1061 ret = H5Aread(attr, atype, &version[0]);
1062 H5Tclose(atype);
1063
1064 H5Aclose(attr);
1065
1066 AMREX_ASSERT(!version.empty());
1067
1068 // What do our version strings mean?
1069 // "Version_One_Dot_Zero" -- hard-wired to write out in double precision.
1070 // "Version_One_Dot_One" -- can write out either as either single or double precision.
1071 // Appended to the latter version string are either "_single" or "_double" to
1072 // indicate how the particles were written.
1073 // "Version_Two_Dot_Zero" -- this is the AMReX particle file format
1074 // "Version_Two_Dot_One" -- expanded particle ids to allow for 2**39-1 per proc
1075 std::string how;
1076 bool convert_ids = false;
1077 if (version.find("Version_Two_Dot_One") != std::string::npos) {
1078 convert_ids = true;
1079 }
1080 if (version.find("_single") != std::string::npos) {
1081 how = "single";
1082 }
1083 else if (version.find("_double") != std::string::npos) {
1084 how = "double";
1085 }
1086 else {
1087 std::string msg("ParticleContainer::Restart(): bad version string: ");
1088 msg += version;
1089 amrex::Error(version.c_str());
1090 }
1091
1092 std::string gname = "Chombo_global";
1093 std::string aname = "SpaceDim";
1094 int dm;
1095 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1096 if (grp < 0) {
1097 std::string msg("ParticleContainer::RestartHDF5(): unable to open group ");
1098 msg += gname;
1099 amrex::Abort(msg.c_str());
1100 }
1101 ret = ReadHDF5Attr(grp, aname.c_str(), &dm, H5T_NATIVE_INT);
1102 if (ret < 0) {
1103 std::string msg("ParticleContainer::RestartHDF5(): unable to read attribute ");
1104 msg += aname;
1105 amrex::Abort(msg.c_str());
1106 }
1107 H5Gclose(grp);
1108
1109 aname = "num_component_real";
1110 int nr;
1111 ret = ReadHDF5Attr(fid, aname.c_str(), &nr, H5T_NATIVE_INT);
1112 if (ret < 0) {
1113 std::string msg("ParticleContainer::RestartHDF5(): unable to read attribute ");
1114 msg += aname;
1115 amrex::Abort(msg.c_str());
1116 }
1117 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
1118 if (nr != nrc) {
1119 amrex::Abort("ParticleContainer::RestartHDF5(): nr not the expected value");
1120 }
1121
1122 aname = "num_component_int";
1123 int ni;
1124 ret = ReadHDF5Attr(fid, aname.c_str(), &ni, H5T_NATIVE_INT);
1125 if (ret < 0) {
1126 std::string msg("ParticleContainer::RestartHDF5(): unable to read attribute ");
1127 msg += aname;
1128 amrex::Abort(msg.c_str());
1129 }
1130 if (ni != NStructInt + NumIntComps()) {
1131 amrex::Abort("ParticleContainer::Restart(): ni != NStructInt");
1132 }
1133
1134 aname = "nparticles";
1135 Long nparticles;
1136 ret = ReadHDF5Attr(fid, aname.c_str(), &nparticles, H5T_NATIVE_LLONG);
1137 if (ret < 0) {
1138 std::string msg("ParticleContainer::RestartHDF5(): unable to read attribute ");
1139 msg += aname;
1140 amrex::Abort(msg.c_str());
1141 }
1142 AMREX_ASSERT(nparticles >= 0);
1143
1144 aname = "maxnextid";
1145 Long maxnextid;
1146 ret = ReadHDF5Attr(fid, aname.c_str(), &maxnextid, H5T_NATIVE_LLONG);
1147 if (ret < 0) {
1148 std::string msg("ParticleContainer::RestartHDF5(): unable to read attribute ");
1149 msg += aname;
1150 amrex::Abort(msg.c_str());
1151 }
1152 AMREX_ASSERT(maxnextid > 0);
1153 ParticleType::NextID(maxnextid);
1154
1155 aname = "finest_level";
1156 int finest_level_in_file;
1157 ret = ReadHDF5Attr(fid, aname.c_str(), &finest_level_in_file, H5T_NATIVE_INT);
1158 if (ret < 0) {
1159 std::string msg("ParticleContainer::RestartHDF5(): unable to read attribute ");
1160 msg += aname;
1161 amrex::Abort(msg.c_str());
1162 }
1163 AMREX_ASSERT(finest_level_in_file >= 0);
1164
1165 hid_t comp_dtype = H5Tcreate (H5T_COMPOUND, 2 * AMREX_SPACEDIM * sizeof(int));
1166 if (1 == AMREX_SPACEDIM) {
1167 H5Tinsert (comp_dtype, "lo_i", 0 * sizeof(int), H5T_NATIVE_INT);
1168 H5Tinsert (comp_dtype, "hi_i", 1 * sizeof(int), H5T_NATIVE_INT);
1169 }
1170 else if (2 == AMREX_SPACEDIM) {
1171 H5Tinsert (comp_dtype, "lo_i", 0 * sizeof(int), H5T_NATIVE_INT);
1172 H5Tinsert (comp_dtype, "lo_j", 1 * sizeof(int), H5T_NATIVE_INT);
1173 H5Tinsert (comp_dtype, "hi_i", 2 * sizeof(int), H5T_NATIVE_INT);
1174 H5Tinsert (comp_dtype, "hi_j", 3 * sizeof(int), H5T_NATIVE_INT);
1175 }
1176 else if (3 == AMREX_SPACEDIM) {
1177 H5Tinsert (comp_dtype, "lo_i", 0 * sizeof(int), H5T_NATIVE_INT);
1178 H5Tinsert (comp_dtype, "lo_j", 1 * sizeof(int), H5T_NATIVE_INT);
1179 H5Tinsert (comp_dtype, "lo_k", 2 * sizeof(int), H5T_NATIVE_INT);
1180 H5Tinsert (comp_dtype, "hi_i", 3 * sizeof(int), H5T_NATIVE_INT);
1181 H5Tinsert (comp_dtype, "hi_j", 4 * sizeof(int), H5T_NATIVE_INT);
1182 H5Tinsert (comp_dtype, "hi_k", 5 * sizeof(int), H5T_NATIVE_INT);
1183 }
1184
1185 // Determine whether this is a dual-grid restart or not.
1186 Vector<BoxArray> particle_box_arrays(finest_level_in_file + 1);
1187 bool dual_grid = false;
1188
1189 for (int lev = 0; lev <= finest_level_in_file; lev++)
1190 {
1191 if (lev > finestLevel())
1192 {
1193 dual_grid = true;
1194 break;
1195 }
1196
1197 gname = "level_" + std::to_string(lev);
1198 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1199 if (grp < 0) {
1200 std::string msg("ParticleContainer::RestartHDF5(): unable to open group ");
1201 msg += gname;
1202 amrex::Abort(msg.c_str());
1203 }
1204
1205 dset = H5Dopen(grp, "boxes", H5P_DEFAULT);
1206 dspace = H5Dget_space(dset);
1207 hsize_t ngrid;
1208 H5Sget_simple_extent_dims(dspace, &ngrid, NULL);
1209
1210 Vector<int> boxes(ngrid*AMREX_SPACEDIM*2);
1211 ret = H5Dread(dset, comp_dtype, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(boxes[0]));
1212 if (ret < 0) {
1213 std::string msg("ParticleContainer::RestartHDF5(): unable to read nparticles_grid dataset");
1214 amrex::Abort(msg.c_str());
1215 }
1216 H5Sclose(dspace);
1217 H5Dclose(dset);
1218 H5Gclose(grp);
1219
1220 // Build the BoxArray from the boxes stored in the file.
1221 // Note: the compound dtype lays out each box as
1222 // (lo_0, lo_1, lo_2, hi_0, hi_1, hi_2) stored consecutively,
1223 // so lo[d] is at index i*2*AMREX_SPACEDIM+d and
1224 // hi[d] is at index i*2*AMREX_SPACEDIM+AMREX_SPACEDIM+d.
1225 {
1226 BoxList bl;
1227 for (int i = 0; i < (int)ngrid; i++) {
1228 const int base = i * 2 * AMREX_SPACEDIM;
1229 Box tmp(IntVect{AMREX_D_DECL(boxes[base+0],
1230 boxes[base+1],
1231 boxes[base+2])},
1232 IntVect{AMREX_D_DECL(boxes[base+AMREX_SPACEDIM+0],
1233 boxes[base+AMREX_SPACEDIM+1],
1234 boxes[base+AMREX_SPACEDIM+2])});
1235 bl.push_back(tmp);
1236 }
1237 particle_box_arrays[lev] = BoxArray(std::move(bl));
1238 }
1239
1240 if (! particle_box_arrays[lev].CellEqual(ParticleBoxArray(lev)))
1241 dual_grid = true;
1242 }
1243
1244 if (dual_grid) {
1245 for (int lev = 0; lev <= finestLevel(); lev++) {
1246 SetParticleBoxArray(lev, particle_box_arrays[lev]);
1247 DistributionMapping pdm(particle_box_arrays[lev]);
1248 SetParticleDistributionMap(lev, pdm);
1249 }
1250 }
1251
1252 Vector<int> ngrids(finest_level_in_file+1);
1253 for (int lev = 0; lev <= finest_level_in_file; lev++) {
1254 gname = "level_" + std::to_string(lev);
1255 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1256 if (grp < 0) {
1257 std::string msg("ParticleContainer::RestartHDF5(): unable to open group ");
1258 msg += gname;
1259 amrex::Abort(msg.c_str());
1260 }
1261
1262 aname = "ngrids";
1263 ret = ReadHDF5Attr(grp, aname.c_str(), &ngrids[lev], H5T_NATIVE_INT);
1264 if (ret < 0) {
1265 std::string msg("ParticleContainer::RestartHDF5(): unable to read attribute ");
1266 msg += aname;
1267 amrex::Abort(msg.c_str());
1268 }
1269
1270 AMREX_ASSERT(ngrids[lev] > 0);
1271 if (lev <= finestLevel()) {
1272 AMREX_ASSERT(ngrids[lev] == int(ParticleBoxArray(lev).size()));
1273 }
1274
1275 H5Gclose(grp);
1276 }
1277
1278 resizeData();
1279
1280 if (finest_level_in_file > finestLevel()) {
1281 m_particles.resize(finest_level_in_file+1);
1282 }
1283
1284 for (int lev = 0; lev <= finest_level_in_file; lev++) {
1285
1286 gname = "level_" + std::to_string(lev);
1287 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1288 if (grp < 0) {
1289 std::string msg("ParticleContainer::RestartHDF5(): unable to open group ");
1290 msg += gname;
1291 amrex::Abort(msg.c_str());
1292 }
1293
1294 dset = H5Dopen(grp, "nparticles_grid", H5P_DEFAULT);
1295 Vector<int> count(ngrids[lev]);
1296 ret = H5Dread(dset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(count[0]));
1297 if (ret < 0) {
1298 std::string msg("ParticleContainer::RestartHDF5(): unable to read nparticles_grid dataset");
1299 amrex::Abort(msg.c_str());
1300 }
1301 H5Dclose(dset);
1302
1303 Vector<hsize_t> offset(ngrids[lev]);
1304 offset[0] = 0;
1305 for (int i = 1; i < ngrids[lev]; i++) {
1306 offset[i] = offset[i-1] + count[i-1];
1307 }
1308
1309 int_dset = H5Dopen(grp, "data:datatype=0", H5P_DEFAULT);
1310 if (int_dset < 0) {
1311 std::string msg("ParticleContainer::RestartHDF5(): unable to open int dataset");
1312 amrex::Abort(msg.c_str());
1313 }
1314 real_dset = H5Dopen(grp, "data:datatype=1", H5P_DEFAULT);
1315 if (real_dset < 0) {
1316 std::string msg("ParticleContainer::RestartHDF5(): unable to open int dataset");
1317 amrex::Abort(msg.c_str());
1318 }
1319
1320 Vector<int> grids_to_read;
1321 if (lev <= finestLevel()) {
1322 for (MFIter mfi(*m_dummy_mf[lev]); mfi.isValid(); ++mfi) {
1323 grids_to_read.push_back(mfi.index());
1324 }
1325 } else {
1326
1327 // we lost a level on restart. we still need to read in particles
1328 // on finer levels, and put them in the right place via Redistribute()
1329
1330 const int rank = ParallelDescriptor::MyProc();
1331 const int NReaders = MaxReaders();
1332 if (rank >= NReaders) {
1333 H5Dclose(int_dset);
1334 H5Dclose(real_dset);
1335 H5Gclose(grp);
1336 continue;
1337 }
1338
1339 const int Navg = ngrids[lev] / NReaders;
1340 const int Nleft = ngrids[lev] - Navg * NReaders;
1341
1342 int lo, hi;
1343 if (rank < Nleft) {
1344 lo = rank*(Navg + 1);
1345 hi = lo + Navg + 1;
1346 }
1347 else {
1348 lo = rank * Navg + Nleft;
1349 hi = lo + Navg;
1350 }
1351
1352 for (int i = lo; i < hi; ++i) {
1353 grids_to_read.push_back(i);
1354 }
1355 }
1356
1357 for(int igrid = 0; igrid < static_cast<int>(grids_to_read.size()); ++igrid) {
1358 const int grid = grids_to_read[igrid];
1359
1360 if (how == "single") {
1361 ReadParticlesHDF5<float>(offset[grid], count[grid], grid, lev, int_dset, real_dset, finest_level_in_file, convert_ids);
1362 }
1363 else if (how == "double") {
1364 ReadParticlesHDF5<double>(offset[grid], count[grid], grid, lev, int_dset, real_dset, finest_level_in_file, convert_ids);
1365 }
1366 else {
1367 std::string msg("ParticleContainer::Restart(): bad parameter: ");
1368 msg += how;
1369 amrex::Error(msg.c_str());
1370 }
1371 }
1372
1373 H5Dclose(int_dset);
1374 H5Dclose(real_dset);
1375 H5Gclose(grp);
1376 } // end for level
1377
1378 H5Tclose(comp_dtype);
1379 H5Pclose(fapl);
1380 H5Fclose(fid);
1381
1382 Redistribute();
1383
1384 AMREX_ASSERT(OK());
1385
1386 if (m_verbose > 1) {
1387 auto stoptime = amrex::second() - strttime;
1389 amrex::Print() << "ParticleContainer::Restart() time: " << stoptime << '\n';
1390 }
1391}
1392
1393// Read a batch of particles from the checkpoint file
1394template <typename ParticleType, int NArrayReal, int NArrayInt,
1395 template<class> class Allocator, class CellAssignor>
1396template <class RTYPE>
1397void
1398ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1399::ReadParticlesHDF5 (hsize_t offset, hsize_t cnt, int grd, int lev,
1400 hid_t int_dset, hid_t real_dset, int finest_level_in_file,
1401 bool convert_ids)
1402{
1403 BL_PROFILE("ParticleContainer::ReadParticlesHDF5()");
1404 AMREX_ASSERT(cnt > 0);
1405 AMREX_ASSERT(lev < int(m_particles.size()));
1406
1407 hid_t int_dspace, int_fspace, real_dspace, real_fspace;
1408
1409 // First read in the integer data in binary. We do not store
1410 // the m_lev and m_grid data on disk. We can easily recreate
1411 // that given the structure of the checkpoint file.
1412 const int iChunkSize = 2 + NStructInt + NumIntComps();
1413 Vector<int> istuff(cnt*iChunkSize);
1414 int_fspace = H5Dget_space(int_dset);
1415 hsize_t int_cnt = cnt*iChunkSize;
1416 hsize_t int_offset = offset*iChunkSize;
1417 int_dspace = H5Screate_simple(1, &int_cnt, NULL);
1418 H5Sselect_hyperslab (int_fspace, H5S_SELECT_SET, &int_offset, NULL, &int_cnt, NULL);
1419 H5Dread(int_dset, H5T_NATIVE_INT, int_dspace, int_fspace, H5P_DEFAULT, istuff.dataPtr());
1420
1421 H5Sclose(int_fspace);
1422 H5Sclose(int_dspace);
1423
1424 // Then the real data in binary.
1425 // For pure SoA, positions are stored in the first AMREX_SPACEDIM SoA real slots;
1426 // there are no struct reals beyond positions, so total = NStructReal + NumRealComps().
1427 const int rChunkSize = ParticleType::is_soa_particle ?
1428 NStructReal + NumRealComps() : AMREX_SPACEDIM + NStructReal + NumRealComps();
1429 Vector<RTYPE> rstuff(cnt*rChunkSize);
1430 real_fspace = H5Dget_space(real_dset);
1431 hsize_t real_cnt = cnt*rChunkSize;
1432 hsize_t real_offset = offset*rChunkSize;
1433 real_dspace = H5Screate_simple(1, &real_cnt, NULL);
1434 H5Sselect_hyperslab (real_fspace, H5S_SELECT_SET, &real_offset, NULL, &real_cnt, NULL);
1435 if (sizeof(RTYPE) == 4) {
1436 H5Dread(real_dset, H5T_NATIVE_FLOAT, real_dspace, real_fspace, H5P_DEFAULT, rstuff.dataPtr());
1437 } else {
1438 H5Dread(real_dset, H5T_NATIVE_DOUBLE, real_dspace, real_fspace, H5P_DEFAULT, rstuff.dataPtr());
1439 }
1440
1441 H5Sclose(real_fspace);
1442 H5Sclose(real_dspace);
1443
1444 // Now reassemble the particles.
1445 int* iptr = istuff.dataPtr();
1446 RTYPE* rptr = rstuff.dataPtr();
1447
1448 // Use a plain AoS-style particle struct for reading id/cpu/pos/struct data,
1449 // regardless of whether the container is pure SoA.
1450 Particle<NStructReal, NStructInt> ptemp;
1451 ParticleLocData pld;
1452
1453 Vector<std::map<std::pair<int, int>, Gpu::HostVector<ParticleType> > > host_particles;
1454 host_particles.reserve(15);
1455 host_particles.resize(finest_level_in_file+1);
1456
1457 Vector<std::map<std::pair<int, int>,
1458 std::vector<Gpu::HostVector<ParticleReal> > > > host_real_attribs;
1459 host_real_attribs.reserve(15);
1460 host_real_attribs.resize(finest_level_in_file+1);
1461
1462 Vector<std::map<std::pair<int, int>,
1463 std::vector<Gpu::HostVector<int> > > > host_int_attribs;
1464 host_int_attribs.reserve(15);
1465 host_int_attribs.resize(finest_level_in_file+1);
1466
1467 Vector<std::map<std::pair<int, int>, Gpu::HostVector<uint64_t> > > host_idcpu;
1468 host_idcpu.reserve(15);
1469 host_idcpu.resize(finest_level_in_file+1);
1470
1471 for (hsize_t i = 0; i < cnt; i++) {
1472 // note: for pure SoA particle layouts, we do write the id, cpu and positions as a struct
1473 // for backwards compatibility with readers
1474 if (convert_ids) {
1475 std::int32_t xi, yi;
1476 std::uint32_t xu, yu;
1477 xi = iptr[0];
1478 yi = iptr[1];
1479 std::memcpy(&xu, &xi, sizeof(xi));
1480 std::memcpy(&yu, &yi, sizeof(yi));
1481 ptemp.m_idcpu = ((std::uint64_t)xu) << 32 | yu;
1482 } else {
1483 ptemp.id() = iptr[0];
1484 ptemp.cpu() = iptr[1];
1485 }
1486 iptr += 2;
1487
1488 for (int j = 0; j < NStructInt; j++)
1489 {
1490 ptemp.idata(j) = *iptr;
1491 ++iptr;
1492 }
1493
1494 AMREX_ASSERT(ptemp.id().is_valid());
1495
1496 AMREX_D_TERM(ptemp.pos(0) = ParticleReal(rptr[0]);,
1497 ptemp.pos(1) = ParticleReal(rptr[1]);,
1498 ptemp.pos(2) = ParticleReal(rptr[2]););
1499
1500 rptr += AMREX_SPACEDIM;
1501
1502 for (int j = 0; j < NStructReal; j++)
1503 {
1504 ptemp.rdata(j) = ParticleReal(*rptr);
1505 ++rptr;
1506 }
1507
1508 locateParticle(ptemp, pld, 0, finestLevel(), 0);
1509
1510 std::pair<int, int> ind(grd, pld.m_tile);
1511
1512 host_real_attribs[lev][ind].resize(NumRealComps());
1513 host_int_attribs[lev][ind].resize(NumIntComps());
1514
1515 if constexpr (!ParticleType::is_soa_particle)
1516 {
1517 // add the struct
1518 host_particles[lev][ind].push_back(ptemp);
1519
1520 // add the real...
1521 for (int icomp = 0; icomp < NumRealComps(); icomp++) {
1522 host_real_attribs[lev][ind][icomp].push_back(ParticleReal(*rptr));
1523 ++rptr;
1524 }
1525
1526 // ... and int array data
1527 for (int icomp = 0; icomp < NumIntComps(); icomp++) {
1528 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1529 ++iptr;
1530 }
1531 } else {
1532 // we need to create this map entry for the loop below
1533 host_particles[lev][ind];
1534
1535 for (int j = 0; j < AMREX_SPACEDIM; j++) {
1536 host_real_attribs[lev][ind][j].push_back(ptemp.pos(j));
1537 }
1538
1539 host_idcpu[lev][ind].push_back(ptemp.m_idcpu);
1540
1541 // read remaining SoA real data (starting after positions)
1542 for (int icomp = AMREX_SPACEDIM; icomp < NumRealComps(); icomp++) {
1543 host_real_attribs[lev][ind][icomp].push_back(ParticleReal(*rptr));
1544 ++rptr;
1545 }
1546
1547 // ... and int array data
1548 for (int icomp = 0; icomp < NumIntComps(); icomp++) {
1549 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1550 ++iptr;
1551 }
1552 }
1553 }
1554
1555 for (int host_lev = 0; host_lev < static_cast<int>(host_particles.size()); ++host_lev)
1556 {
1557 for (auto& kv : host_particles[host_lev]) {
1558 auto grid = kv.first.first;
1559 auto tile = kv.first.second;
1560 const auto& src_tile = kv.second;
1561
1562 auto& dst_tile = DefineAndReturnParticleTile(host_lev, grid, tile);
1563 auto old_size = dst_tile.size();
1564 auto new_size = old_size;
1565 if constexpr (!ParticleType::is_soa_particle)
1566 {
1567 new_size += src_tile.size();
1568 } else {
1569 amrex::ignore_unused(src_tile);
1570 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].size();
1571 }
1572 dst_tile.resize(new_size);
1573
1574 if constexpr (!ParticleType::is_soa_particle)
1575 {
1576 Gpu::copyAsync(Gpu::hostToDevice, src_tile.begin(), src_tile.end(),
1577 dst_tile.GetArrayOfStructs().begin() + old_size);
1578 } else {
1580 host_idcpu[host_lev][std::make_pair(grid,tile)].begin(),
1581 host_idcpu[host_lev][std::make_pair(grid,tile)].end(),
1582 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1583 }
1584
1585 for (int i = 0; i < NumRealComps(); ++i) {
1587 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1588 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1589 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1590 }
1591
1592 for (int i = 0; i < NumIntComps(); ++i) {
1594 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1595 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1596 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
1597 }
1598 }
1599 }
1600
1602}
1603
1604}
1605
1606#endif /*AMREX_USE_HDF5*/
1607#endif /*AMREX_PARTICLEHDF5_H*/
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
void CheckpointPreHDF5()
Prepare metadata used by deferred HDF5 checkpoint writes.
void RestartHDF5(const std::string &dir, const std::string &file)
Restart from checkpoint.
void WriteHDF5ParticleData(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=false) const
Writes particle data to disk in the HDF5 format.
void CheckpointPostHDF5()
Finalize a checkpoint header handed off by CheckpointPreHDF5().
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
const std::string & FileName() const
Definition AMReX_NFiles.H:160
This class provides the user with a few print options.
Definition AMReX_Print.H:35
amrex_ulong ULong
Unsigned integer type guaranteed to be wider than unsigned int.
Definition AMReX_INT.H:32
amrex_particle_real ParticleReal
Floating Point Type for Particles.
Definition AMReX_REAL.H:90
amrex_long Long
Definition AMReX_INT.H:30
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
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
bool Remove(std::string const &filename)
Definition AMReX_FileSystem.cpp:190
void copyAsync(HostToDevice, InIter begin, InIter end, OutIter result) noexcept
A host-to-device copy routine. Note this is just a wrapper around memcpy, so it assumes contiguous st...
Definition AMReX_GpuContainers.H:228
static constexpr HostToDevice hostToDevice
Definition AMReX_GpuContainers.H:105
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
MPI_Comm Communicator() noexcept
Definition AMReX_ParallelDescriptor.H:223
void ReduceRealMax(Vector< std::reference_wrapper< Real > > const &)
Definition AMReX_ParallelDescriptor.cpp:1228
__host__ __device__ bool is_valid(const uint64_t idcpu) noexcept
Definition AMReX_Particle.H:146
Definition AMReX_Amr.cpp:49
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2006
double second() noexcept
Definition AMReX_Utility.cpp:940
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:234
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
const int[]
Definition AMReX_BLProfiler.cpp:1664
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
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2015
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