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