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