Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleIO.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLEIO_H
2#define AMREX_PARTICLEIO_H
3#include <AMReX_Config.H>
4
6#include <AMReX_VisMF.H>
7
8#include <iterator>
9
10namespace amrex {
11
12template <typename ParticleType, int NArrayReal, int NArrayInt,
13 template<class> class Allocator, class CellAssignor>
14void
15ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
16::WriteParticleRealData (void* data, size_t size, std::ostream& os) const
17{
18 if (sizeof(typename ParticleType::RealType) == 4) {
19 writeFloatData((float*) data, size, os, ParticleRealDescriptor);
20 }
21 else if (sizeof(typename ParticleType::RealType) == 8) {
22 writeDoubleData((double*) data, size, os, ParticleRealDescriptor);
23 }
24}
25
26template <typename ParticleType, int NArrayReal, int NArrayInt,
27 template<class> class Allocator, class CellAssignor>
28void
30::ReadParticleRealData (void* data, size_t size, std::istream& is)
31{
32 if (sizeof(typename ParticleType::RealType) == 4) {
33 readFloatData((float*) data, size, is, ParticleRealDescriptor);
34 }
35 else if (sizeof(typename ParticleType::RealType) == 8) {
36 readDoubleData((double*) data, size, is, ParticleRealDescriptor);
37 }
38}
39
41 template <typename P>
43 int operator() (const P& p) const {
44 return p.id().is_valid();
45 }
46};
47
49template <typename ParticleType, int NArrayReal, int NArrayInt,
50 template<class> class Allocator, class CellAssignor>
51void
52ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
53::Checkpoint (const std::string& dir,
54 const std::string& name, bool /*is_checkpoint*/,
55 const Vector<std::string>& real_comp_names,
56 const Vector<std::string>& int_comp_names) const
57{
58 Vector<int> write_real_comp;
59 Vector<std::string> tmp_real_comp_names;
60
61 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
62 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
63 {
64 write_real_comp.push_back(1);
65 if (real_comp_names.empty())
66 {
67 tmp_real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
68 }
69 else
70 {
71 tmp_real_comp_names.push_back(real_comp_names[i-first_rcomp]);
72 }
73 }
74
75 Vector<int> write_int_comp;
76 Vector<std::string> tmp_int_comp_names;
77 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
78 {
79 write_int_comp.push_back(1);
80 if (int_comp_names.empty())
81 {
82 tmp_int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
83 }
84 else
85 {
86 tmp_int_comp_names.push_back(int_comp_names[i]);
87 }
88 }
89
90 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
91 tmp_real_comp_names, tmp_int_comp_names,
92 FilterPositiveID{}, true);
93}
94
95template <typename ParticleType, int NArrayReal, int NArrayInt,
96 template<class> class Allocator, class CellAssignor>
97void
99::Checkpoint (const std::string& dir, const std::string& name,
100 const Vector<int>& write_real_comp,
101 const Vector<int>& write_int_comp,
102 const Vector<std::string>& real_comp_names,
103 const Vector<std::string>& int_comp_names) const
104{
105 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
106 real_comp_names, int_comp_names,
107 FilterPositiveID{}, true);
108}
109
111template <typename ParticleType, int NArrayReal, int NArrayInt,
112 template<class> class Allocator, class CellAssignor>
113void
115::WritePlotFile (const std::string& dir, const std::string& name) const
116{
117 Vector<int> write_real_comp;
118 Vector<std::string> real_comp_names;
119
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 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
136 real_comp_names, int_comp_names,
138}
139
141template <typename ParticleType, int NArrayReal, int NArrayInt,
142 template<class> class Allocator, class CellAssignor>
143void
145::WritePlotFile (const std::string& dir, const std::string& name,
146 const Vector<std::string>& real_comp_names,
147 const Vector<std::string>& int_comp_names) const
148{
149 if constexpr(ParticleType::is_soa_particle) {
150 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
151 } else {
152 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal);
153 }
154 AMREX_ASSERT( int_comp_names.size() == NStructInt + NumIntComps() );
155
156 Vector<int> write_real_comp;
157 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
158 for (int i = 0; i < nrc; ++i) {
159 write_real_comp.push_back(1);
160 }
161
162 Vector<int> write_int_comp;
163 for (int i = 0; i < NStructInt + NumIntComps(); ++i) {
164 write_int_comp.push_back(1);
165 }
166
167 WriteBinaryParticleData(dir, name,
168 write_real_comp, write_int_comp,
169 real_comp_names, int_comp_names,
171}
172
174template <typename ParticleType, int NArrayReal, int NArrayInt,
175 template<class> class Allocator, class CellAssignor>
176void
178::WritePlotFile (const std::string& dir, const std::string& name,
179 const Vector<std::string>& real_comp_names) const
180{
181 if constexpr(ParticleType::is_soa_particle) {
182 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
183 } else {
184 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal);
185 }
186
187 Vector<int> write_real_comp;
188 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
189 for (int i = 0; i < nrc; ++i) {
190 write_real_comp.push_back(1);
191 }
192
193 Vector<int> write_int_comp;
194 for (int i = 0; i < NStructInt + NumIntComps(); ++i) {
195 write_int_comp.push_back(1);
196 }
197
198 Vector<std::string> int_comp_names;
199 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
200 {
201 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
202 }
203
204 WriteBinaryParticleData(dir, name,
205 write_real_comp, write_int_comp,
206 real_comp_names, int_comp_names,
208}
209
211template <typename ParticleType, int NArrayReal, int NArrayInt,
212 template<class> class Allocator, class CellAssignor>
213void
215::WritePlotFile (const std::string& dir,
216 const std::string& name,
217 const Vector<int>& write_real_comp,
218 const Vector<int>& write_int_comp) const
219{
220
221 if constexpr(ParticleType::is_soa_particle) {
222 AMREX_ALWAYS_ASSERT(write_real_comp.size() == NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
223 } else {
224 AMREX_ALWAYS_ASSERT(write_real_comp.size() == NumRealComps() + NStructReal);
225 }
226 AMREX_ASSERT(write_int_comp.size() == NStructInt + NumIntComps() );
227
228 Vector<std::string> real_comp_names;
229 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
230 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
231 {
232 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
233 }
234
235 Vector<std::string> int_comp_names;
236 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
237 {
238 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
239 }
240
241 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
242 real_comp_names, int_comp_names,
244}
245
247template <typename ParticleType, int NArrayReal, int NArrayInt,
248 template<class> class Allocator, class CellAssignor>
249void
251WritePlotFile (const std::string& dir, const std::string& name,
252 const Vector<int>& write_real_comp,
253 const Vector<int>& write_int_comp,
254 const Vector<std::string>& real_comp_names,
255 const Vector<std::string>& int_comp_names) const
256{
257 BL_PROFILE("ParticleContainer::WritePlotFile()");
258
259 WriteBinaryParticleData(dir, name,
260 write_real_comp, write_int_comp,
261 real_comp_names, int_comp_names,
263}
264
266template <typename ParticleType, int NArrayReal, int NArrayInt,
267 template<class> class Allocator, class CellAssignor>
268template <class F>
269void
271::WritePlotFile (const std::string& dir, const std::string& name, F const& f) const
272{
273 Vector<int> write_real_comp;
274 Vector<std::string> real_comp_names;
275
276 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
277 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
278 {
279 write_real_comp.push_back(1);
280 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
281 }
282
283 Vector<int> write_int_comp;
284 Vector<std::string> int_comp_names;
285 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
286 {
287 write_int_comp.push_back(1);
288 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
289 }
290
291 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
292 real_comp_names, int_comp_names, f);
293}
294
296template <typename ParticleType, int NArrayReal, int NArrayInt,
297 template<class> class Allocator, class CellAssignor>
298template <class F>
299void
301::WritePlotFile (const std::string& dir, const std::string& name,
302 const Vector<std::string>& real_comp_names,
303 const Vector<std::string>& int_comp_names, F const& f) const
304{
305 if constexpr(ParticleType::is_soa_particle) {
306 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
307 } else {
308 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal);
309 }
310 AMREX_ASSERT( int_comp_names.size() == NStructInt + NumIntComps() );
311
312 Vector<int> write_real_comp;
313 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
314 for (int i = 0; i < nrc; ++i) {
315 write_real_comp.push_back(1);
316 }
317
318 Vector<int> write_int_comp;
319 for (int i = 0; i < NStructInt + NumIntComps(); ++i) {
320 write_int_comp.push_back(1);
321 }
322
323 WriteBinaryParticleData(dir, name,
324 write_real_comp, write_int_comp,
325 real_comp_names, int_comp_names, f);
326}
327
329template <typename ParticleType, int NArrayReal, int NArrayInt,
330 template<class> class Allocator, class CellAssignor>
331template <class F>
332void
334::WritePlotFile (const std::string& dir, const std::string& name,
335 const Vector<std::string>& real_comp_names, F const& f) const
336{
337 if constexpr(ParticleType::is_soa_particle) {
338 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
339 } else {
340 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal);
341 }
342
343 Vector<int> write_real_comp;
344 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
345 for (int i = 0; i < nrc; ++i) {
346 write_real_comp.push_back(1);
347 }
348
349 Vector<int> write_int_comp;
350 for (int i = 0; i < NStructInt + NumIntComps(); ++i) {
351 write_int_comp.push_back(1);
352 }
353
354 Vector<std::string> int_comp_names;
355 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
356 {
357 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
358 }
359
360 WriteBinaryParticleData(dir, name,
361 write_real_comp, write_int_comp,
362 real_comp_names, int_comp_names, f);
363}
364
366template <typename ParticleType, int NArrayReal, int NArrayInt,
367 template<class> class Allocator, class CellAssignor>
368template <class F>
369void
371::WritePlotFile (const std::string& dir,
372 const std::string& name,
373 const Vector<int>& write_real_comp,
374 const Vector<int>& write_int_comp, F const& f) const
375{
376 if constexpr(ParticleType::is_soa_particle) {
377 AMREX_ALWAYS_ASSERT(write_real_comp.size() == NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
378 } else {
379 AMREX_ALWAYS_ASSERT(write_real_comp.size() == NumRealComps() + NStructReal);
380 }
381 AMREX_ASSERT(write_int_comp.size() == NStructInt + NumIntComps() );
382
383 Vector<std::string> real_comp_names;
384 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
385 for (int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
386 {
387 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
388 }
389
390 Vector<std::string> int_comp_names;
391 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
392 {
393 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
394 }
395
396 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
397 real_comp_names, int_comp_names, f);
398}
399
401template <typename ParticleType, int NArrayReal, int NArrayInt,
402 template<class> class Allocator, class CellAssignor>
403template <class F>
404void
406WritePlotFile (const std::string& dir, const std::string& name,
407 const Vector<int>& write_real_comp,
408 const Vector<int>& write_int_comp,
409 const Vector<std::string>& real_comp_names,
410 const Vector<std::string>& int_comp_names,
411 F const& f) const
412{
413 BL_PROFILE("ParticleContainer::WritePlotFile()");
414
415 WriteBinaryParticleData(dir, name,
416 write_real_comp, write_int_comp,
417 real_comp_names, int_comp_names, f);
418}
419
420template <typename ParticleType, int NArrayReal, int NArrayInt,
421 template<class> class Allocator, class CellAssignor>
422template <class F>
423void
425::WriteBinaryParticleData (const std::string& dir, const std::string& name,
426 const Vector<int>& write_real_comp,
427 const Vector<int>& write_int_comp,
428 const Vector<std::string>& real_comp_names,
429 const Vector<std::string>& int_comp_names,
430 F const& f, bool is_checkpoint) const
431{
432 if (AsyncOut::UseAsyncOut()) {
433 WriteBinaryParticleDataAsync(*this, dir, name,
434 write_real_comp, write_int_comp,
435 real_comp_names, int_comp_names, is_checkpoint);
436 } else
437 {
438 WriteBinaryParticleDataSync(*this, dir, name,
439 write_real_comp, write_int_comp,
440 real_comp_names, int_comp_names,
441 f, is_checkpoint);
442 }
443}
444
445template <typename ParticleType, int NArrayReal, int NArrayInt,
446 template<class> class Allocator, class CellAssignor>
447void
450{
451 if( ! usePrePost) {
452 return;
453 }
454
455 BL_PROFILE("ParticleContainer::CheckpointPre()");
456
457 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
458 Long nparticles = 0;
459 Long maxnextid = ParticleType::NextID();
460
461 for (int lev = 0; lev < std::ssize(m_particles); lev++) {
462 const auto& pmap = m_particles[lev];
463 for (const auto& kv : pmap) {
464 const auto& aos = kv.second.GetArrayOfStructs();
465 for (int k = 0; k < aos.numParticles(); ++k) {
466 const ParticleType& p = aos[k];
467 if (p.id().is_valid()) {
468 //
469 // Only count (and checkpoint) valid particles.
470 //
471 nparticles++;
472 }
473 }
474 }
475 }
476 ParallelDescriptor::ReduceLongSum(nparticles, IOProcNumber);
477
478 ParticleType::NextID(maxnextid);
479 ParallelDescriptor::ReduceLongMax(maxnextid, IOProcNumber);
480
481 nparticlesPrePost = nparticles;
482 maxnextidPrePost = maxnextid;
483
484 nParticlesAtLevelPrePost.clear();
485 nParticlesAtLevelPrePost.resize(finestLevel() + 1, 0);
486 for(int lev(0); lev <= finestLevel(); ++lev) {
487 nParticlesAtLevelPrePost[lev] = NumberOfParticlesAtLevel(lev);
488 }
489
490 whichPrePost.clear();
491 whichPrePost.resize(finestLevel() + 1);
492 countPrePost.clear();
493 countPrePost.resize(finestLevel() + 1);
494 wherePrePost.clear();
495 wherePrePost.resize(finestLevel() + 1);
496
497 filePrefixPrePost.clear();
498 filePrefixPrePost.resize(finestLevel() + 1);
499}
500
501
502template <typename ParticleType, int NArrayReal, int NArrayInt,
503 template<class> class Allocator, class CellAssignor>
504void
507{
508 if( ! usePrePost) {
509 return;
510 }
511
512 BL_PROFILE("ParticleContainer::CheckpointPost()");
513
514 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
515 std::ofstream HdrFile;
516 HdrFile.open(HdrFileNamePrePost.c_str(), std::ios::out | std::ios::app);
517
518 for(int lev(0); lev <= finestLevel(); ++lev) {
519 ParallelDescriptor::ReduceIntSum (whichPrePost[lev].dataPtr(), whichPrePost[lev].size(), IOProcNumber);
520 ParallelDescriptor::ReduceIntSum (countPrePost[lev].dataPtr(), countPrePost[lev].size(), IOProcNumber);
521 ParallelDescriptor::ReduceLongSum(wherePrePost[lev].dataPtr(), wherePrePost[lev].size(), IOProcNumber);
522
523
525 for(int j(0); j < whichPrePost[lev].size(); ++j) {
526 HdrFile << whichPrePost[lev][j] << ' ' << countPrePost[lev][j] << ' ' << wherePrePost[lev][j] << '\n';
527 }
528
529 const bool gotsome = (nParticlesAtLevelPrePost[lev] > 0);
530 if(gotsome && doUnlink) {
531// BL_PROFILE_VAR("PC<NNNN>::Checkpoint:unlink", unlink_post);
532 // Unlink any zero-length data files.
533 Vector<Long> cnt(nOutFilesPrePost,0);
534
535 for(int i(0), N = countPrePost[lev].size(); i < N; ++i) {
536 cnt[whichPrePost[lev][i]] += countPrePost[lev][i];
537 }
538
539 for(int i = 0; i < std::ssize(cnt); ++i) {
540 if(cnt[i] == 0) {
541 std::string FullFileName = NFilesIter::FileName(i, filePrefixPrePost[lev]);
542 FileSystem::Remove(FullFileName);
543 }
544 }
545 }
546 }
547 }
548
550 HdrFile.flush();
551 HdrFile.close();
552 if( ! HdrFile.good()) {
553 amrex::Abort("ParticleContainer::CheckpointPost(): problem writing HdrFile");
554 }
555 }
556}
557
558template <typename ParticleType, int NArrayReal, int NArrayInt,
559 template<class> class Allocator, class CellAssignor>
560void
563{
564 CheckpointPre();
565}
566
567
568template <typename ParticleType, int NArrayReal, int NArrayInt,
569 template<class> class Allocator, class CellAssignor>
570void
573{
574 CheckpointPost();
575}
576
577
578template <typename ParticleType, int NArrayReal, int NArrayInt,
579 template<class> class Allocator, class CellAssignor>
580void
582::WriteParticles (int lev, std::ofstream& ofs, int fnum,
583 Vector<int>& which, Vector<int>& count, Vector<Long>& where,
584 const Vector<int>& write_real_comp,
585 const Vector<int>& write_int_comp,
586 const Vector<std::map<std::pair<int, int>, FlagsVector>>& particle_io_flags,
587 bool is_checkpoint) const
588{
589 BL_PROFILE("ParticleContainer::WriteParticles()");
590
591 // For a each grid, the tiles it contains
592 std::map<int, Vector<int> > tile_map;
593
594 for (const auto& kv : m_particles[lev])
595 {
596 const int grid = kv.first.first;
597 const int tile = kv.first.second;
598 tile_map[grid].push_back(tile);
599 const auto& pflags = particle_io_flags[lev].at(kv.first);
600
601 // Only write out valid particles.
602 count[grid] += particle_detail::countFlags(pflags);
603 }
604
605 MFInfo info;
606 info.SetAlloc(false);
607 MultiFab state(ParticleBoxArray(lev), ParticleDistributionMap(lev), 1,0,info);
608
609 for (MFIter mfi(state); mfi.isValid(); ++mfi)
610 {
611 const int grid = mfi.index();
612
613 which[grid] = fnum;
614 where[grid] = VisMF::FileOffset(ofs);
615
616 if (count[grid] <= 0) { continue; }
617
618 Vector<int> istuff;
620 particle_detail::packIOData(istuff, rstuff, *this, lev, grid,
621 write_real_comp, write_int_comp,
622 particle_io_flags, tile_map[grid], count[grid], is_checkpoint);
623
624 writeIntData(istuff.dataPtr(), istuff.size(), ofs);
626 ofs.flush(); // Some systems require this flush() (probably due to a bug)
627 }
628
629 WriteParticleRealData(rstuff.dataPtr(), rstuff.size(), ofs);
631 ofs.flush(); // Some systems require this flush() (probably due to a bug)
632 }
633 }
634}
635
636
637template <typename ParticleType, int NArrayReal, int NArrayInt,
638 template<class> class Allocator, class CellAssignor>
639void
641::Restart (const std::string& dir, const std::string& file, bool /*is_checkpoint*/)
642{
643 Restart(dir, file);
644}
645
646template <typename ParticleType, int NArrayReal, int NArrayInt,
647 template<class> class Allocator, class CellAssignor>
648void
650::Restart (const std::string& dir, const std::string& file)
651{
652 BL_PROFILE("ParticleContainer::Restart()");
653 AMREX_ASSERT(!dir.empty());
654 AMREX_ASSERT(!file.empty());
655
656 const auto strttime = amrex::second();
657
658 int DATA_Digits_Read(5);
659 ParmParse pp("particles");
660 pp.query("datadigits_read",DATA_Digits_Read);
661
662 std::string fullname = dir;
663 if (!fullname.empty() && fullname[fullname.size()-1] != '/') {
664 fullname += '/';
665 }
666 fullname += file;
667 std::string HdrFileName = fullname;
668 if (!HdrFileName.empty() && HdrFileName[HdrFileName.size()-1] != '/') {
669 HdrFileName += '/';
670 }
671 HdrFileName += "Header";
672
673 Vector<char> fileCharPtr;
674 ParallelDescriptor::ReadAndBcastFile(HdrFileName, fileCharPtr);
675 std::string fileCharPtrString(fileCharPtr.dataPtr());
676 std::istringstream HdrFile(fileCharPtrString, std::istringstream::in);
677
678 std::string version;
679 HdrFile >> version;
680 AMREX_ASSERT(!version.empty());
681
682 // What do our version strings mean?
683 // "Version_One_Dot_Zero" -- hard-wired to write out in double precision.
684 // "Version_One_Dot_One" -- can write out either as either single or double precision.
685 // Appended to the latter version string are either "_single" or "_double" to
686 // indicate how the particles were written.
687 // "Version_Two_Dot_Zero" -- this is the AMReX particle file format
688 // "Version_Two_Dot_One" -- expanded particle ids to allow for 2**39-1 per proc
689 std::string how;
690 bool convert_ids = false;
691 if (version.find("Version_Two_Dot_One") != std::string::npos) {
692 convert_ids = true;
693 }
694 if (version.find("Version_One_Dot_Zero") != std::string::npos) {
695 how = "double";
696 }
697 else if (version.find("Version_One_Dot_One") != std::string::npos ||
698 version.find("Version_Two_Dot_Zero") != std::string::npos ||
699 version.find("Version_Two_Dot_One") != std::string::npos) {
700 if (version.find("_single") != std::string::npos) {
701 how = "single";
702 }
703 else if (version.find("_double") != std::string::npos) {
704 how = "double";
705 }
706 else {
707 std::string msg("ParticleContainer::Restart(): bad version string: ");
708 msg += version;
709 amrex::Error(version.c_str());
710 }
711 }
712 else {
713 std::string msg("ParticleContainer::Restart(): unknown version string: ");
714 msg += version;
715 amrex::Abort(msg.c_str());
716 }
717
718 int dm;
719 HdrFile >> dm;
720 if (dm != AMREX_SPACEDIM) {
721 amrex::Abort("ParticleContainer::Restart(): dm != AMREX_SPACEDIM");
722 }
723
724 int nr;
725 HdrFile >> nr;
726 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
727 if (nr != nrc) {
728 amrex::Abort("ParticleContainer::Restart(): nr not the expected value");
729 }
730
731 std::string comp_name;
732 for (int i = 0; i < nr; ++i) {
733 HdrFile >> comp_name;
734 }
735
736 int ni;
737 HdrFile >> ni;
738 if (ni != NStructInt + NumIntComps()) {
739 amrex::Abort("ParticleContainer::Restart(): ni != NStructInt");
740 }
741
742 for (int i = 0; i < ni; ++i) {
743 HdrFile >> comp_name;
744 }
745
746 bool checkpoint;
747 HdrFile >> checkpoint;
748
749 Long nparticles;
750 HdrFile >> nparticles;
751 AMREX_ASSERT(nparticles >= 0);
752
753 Long maxnextid;
754 HdrFile >> maxnextid;
755 AMREX_ASSERT(maxnextid > 0);
756 ParticleType::NextID(maxnextid);
757
758 int finest_level_in_file;
759 HdrFile >> finest_level_in_file;
760 AMREX_ASSERT(finest_level_in_file >= 0);
761
762 // Determine whether this is a dual-grid restart or not.
763 // Size the vectors to cover both the levels present in the checkpoint
764 // file and the levels present in the current container, since the
765 // loops below iterate up to finestLevel()
766 int const max_level = std::max(finest_level_in_file, finestLevel());
767 Vector<DistributionMapping> old_dms(max_level + 1);
768 Vector<BoxArray> old_bas(max_level + 1);
769 Vector<BoxArray> particle_box_arrays(max_level + 1);
770 bool dual_grid = false;
771
772 bool have_pheaders = false;
773 for (int lev = 0; lev <= finest_level_in_file; lev++)
774 {
775 std::string phdr_name = fullname;
776 phdr_name += "/Level_";
777 phdr_name = amrex::Concatenate(phdr_name, lev, 1);
778 phdr_name += "/Particle_H";
779
780 if (amrex::FileExists(phdr_name)) {
781 have_pheaders = true;
782 break;
783 }
784 }
785
786 if (have_pheaders)
787 {
788 for (int lev = 0; lev <= finestLevel(); lev++)
789 {
790 old_dms[lev] = ParticleDistributionMap(lev);
791 old_bas[lev] = ParticleBoxArray(lev);
792 std::string phdr_name = fullname;
793 phdr_name += "/Level_";
794 phdr_name = amrex::Concatenate(phdr_name, lev, 1);
795 phdr_name += "/Particle_H";
796
797 if (! amrex::FileExists(phdr_name)) { continue; }
798
799 Vector<char> phdr_chars;
800 ParallelDescriptor::ReadAndBcastFile(phdr_name, phdr_chars);
801 std::string phdr_string(phdr_chars.dataPtr());
802 std::istringstream phdr_file(phdr_string, std::istringstream::in);
803
804 particle_box_arrays[lev].readFrom(phdr_file);
805 if (! particle_box_arrays[lev].CellEqual(ParticleBoxArray(lev))) { dual_grid = true; }
806 }
807 } else // if no particle box array information exists in the file, we assume a single grid restart
808 {
809 dual_grid = false;
810 }
811
812 if (dual_grid) {
813 for (int lev = 0; lev <= finestLevel(); lev++) {
814 // this can happen if there are no particles at a given level in the checkpoint
815 if (particle_box_arrays[lev].empty()) {
816 particle_box_arrays[lev] = BoxArray(Geom(lev).Domain());
817 }
818 SetParticleBoxArray(lev, particle_box_arrays[lev]);
819 DistributionMapping pdm(particle_box_arrays[lev]);
820 SetParticleDistributionMap(lev, pdm);
821 }
822 }
823
824 Vector<int> ngrids(finest_level_in_file+1);
825 for (int lev = 0; lev <= finest_level_in_file; lev++) {
826 HdrFile >> ngrids[lev];
827 AMREX_ASSERT(ngrids[lev] > 0);
829
830 resizeData();
831
832 if (finest_level_in_file > finestLevel()) {
833 m_particles.resize(finest_level_in_file+1);
834 }
835
836 for (int lev = 0; lev <= finest_level_in_file; lev++) {
837 Vector<int> which(ngrids[lev]);
838 Vector<int> count(ngrids[lev]);
839 Vector<Long> where(ngrids[lev]);
840 for (int i = 0; i < ngrids[lev]; i++) {
841 HdrFile >> which[i] >> count[i] >> where[i];
842 }
843
844 Vector<int> grids_to_read;
845 if (lev <= finestLevel()) {
846 for (MFIter mfi(*m_dummy_mf[lev]); mfi.isValid(); ++mfi) {
847 grids_to_read.push_back(mfi.index());
848 }
849 } else {
850
851 // we lost a level on restart. we still need to read in particles
852 // on finer levels, and put them in the right place via Redistribute()
854 const int rank = ParallelDescriptor::MyProc();
855 const int NReaders = MaxReaders();
856 if (rank >= NReaders) { continue; }
857
858 const int Navg = ngrids[lev] / NReaders;
859 const int Nleft = ngrids[lev] - Navg * NReaders;
860
861 int lo, hi;
862 if (rank < Nleft) {
863 lo = rank*(Navg + 1);
864 hi = lo + Navg + 1;
865 }
866 else {
867 lo = rank * Navg + Nleft;
868 hi = lo + Navg;
869 }
870
871 for (int i = lo; i < hi; ++i) {
872 grids_to_read.push_back(i);
873 }
874 }
875
876 for(int grid : grids_to_read) {
877 if (count[grid] <= 0) { continue; }
878
879 // The file names in the header file are relative.
880 std::string name = fullname;
881
882 if (!name.empty() && name[name.size()-1] != '/') {
883 name += '/';
884 }
885
886 name += "Level_";
887 name += amrex::Concatenate("", lev, 1);
888 name += '/';
889 name += DataPrefix();
890 name += amrex::Concatenate("", which[grid], DATA_Digits_Read);
891
892 std::ifstream ParticleFile;
893
894 ParticleFile.open(name.c_str(), std::ios::in | std::ios::binary);
896 if (!ParticleFile.good()) {
898 }
899
900 ParticleFile.seekg(where[grid], std::ios::beg);
901
902 // Use if constexpr to avoid instantiating the mis-matched
903 // type case and triggering the static_assert on the
904 // underlying copy calls
905 if (how == "single") {
906 if constexpr (std::is_same_v<ParticleReal, float>) {
907 ReadParticles<float>(count[grid], grid, lev, ParticleFile, finest_level_in_file, convert_ids);
908 } else {
909 amrex::Error("File contains single-precision data, while AMReX is compiled with ParticleReal==double");
910 }
911 }
912 else if (how == "double") {
913 if constexpr (std::is_same_v<ParticleReal, double>) {
914 ReadParticles<double>(count[grid], grid, lev, ParticleFile, finest_level_in_file, convert_ids);
915 } else {
916 amrex::Error("File contains double-precision data, while AMReX is compiled with ParticleReal==float");
917 }
918 }
919 else {
920 std::string msg("ParticleContainer::Restart(): bad parameter: ");
921 msg += how;
922 amrex::Error(msg.c_str());
923 }
925 ParticleFile.close();
926
927 if (!ParticleFile.good()) {
928 amrex::Abort("ParticleContainer::Restart(): problem reading particles");
929 }
930 }
931 }
932
933 if (dual_grid) {
934 for (int lev = 0; lev <= finestLevel(); lev++) {
935 SetParticleBoxArray(lev, old_bas[lev]);
936 SetParticleDistributionMap(lev, old_dms[lev]);
937 }
938 }
939
940 Redistribute();
941
942 AMREX_ASSERT(OK());
943
944 if (m_verbose > 1) {
945 auto stoptime = amrex::second() - strttime;
947 amrex::Print() << "ParticleContainer::Restart() time: " << stoptime << '\n';
948 }
949}
950
951// Read a batch of particles from the checkpoint file
952template <typename ParticleType, int NArrayReal, int NArrayInt,
953 template<class> class Allocator, class CellAssignor>
954template <class RTYPE>
955void
957::ReadParticles (int cnt, int grd, int lev, std::ifstream& ifs,
958 int finest_level_in_file, bool convert_ids)
959{
960 BL_PROFILE("ParticleContainer::ReadParticles()");
961 AMREX_ASSERT(cnt > 0);
962 AMREX_ASSERT(lev < std::ssize(m_particles));
963
964 // First read in the integer data in binary. We do not store
965 // the m_lev and m_grid data on disk. We can easily recreate
966 // that given the structure of the checkpoint file.
967 const int iChunkSize = 2 + NStructInt + NumIntComps();
968 Vector<int> istuff(std::size_t(cnt)*iChunkSize);
969 readIntData(istuff.dataPtr(), istuff.size(), ifs, FPC::NativeIntDescriptor());
970
971 // Then the real data in binary.
972 const int rChunkSize = ParticleType::is_soa_particle ? NStructReal + NumRealComps() : AMREX_SPACEDIM + NStructReal + NumRealComps();
973 Vector<RTYPE> rstuff(std::size_t(cnt)*rChunkSize);
974 ReadParticleRealData(rstuff.dataPtr(), rstuff.size(), ifs);
975
976 // Now reassemble the particles.
977 int* iptr = istuff.dataPtr();
978 RTYPE* rptr = rstuff.dataPtr();
979
981 ParticleLocData pld;
982
984 host_particles.reserve(15);
985 host_particles.resize(finest_level_in_file+1);
986
988 std::vector<Gpu::HostVector<RTYPE> > > > host_real_attribs;
989 host_real_attribs.reserve(15);
990 host_real_attribs.resize(finest_level_in_file+1);
991
993 std::vector<Gpu::HostVector<int> > > > host_int_attribs;
994 host_int_attribs.reserve(15);
995 host_int_attribs.resize(finest_level_in_file+1);
996
998 host_idcpu.reserve(15);
999 host_idcpu.resize(finest_level_in_file+1);
1001 for (int i = 0; i < cnt; i++) {
1002 // note: for pure SoA particle layouts, we do write the id, cpu and positions as a struct
1003 // for backwards compatibility with readers
1004 if (convert_ids) {
1005 std::int32_t xi, yi;
1006 std::uint32_t xu, yu;
1007 xi = iptr[0];
1008 yi = iptr[1];
1009 std::memcpy(&xu, &xi, sizeof(xi));
1010 std::memcpy(&yu, &yi, sizeof(yi));
1011 ptemp.m_idcpu = ((std::uint64_t)xu) << 32 | yu;
1012 } else {
1013 ptemp.id() = iptr[0];
1014 ptemp.cpu() = iptr[1];
1015 }
1016 iptr += 2;
1017
1018 for (int j = 0; j < NStructInt; j++)
1019 {
1020 ptemp.idata(j) = *iptr;
1021 ++iptr;
1022 }
1023
1024 AMREX_ASSERT(ptemp.id().is_valid());
1025
1026 AMREX_D_TERM(ptemp.pos(0) = ParticleReal(rptr[0]);,
1027 ptemp.pos(1) = ParticleReal(rptr[1]);,
1028 ptemp.pos(2) = ParticleReal(rptr[2]););
1029
1030 rptr += AMREX_SPACEDIM;
1031
1032 for (int j = 0; j < NStructReal; j++)
1033 {
1034 ptemp.rdata(j) = ParticleReal(*rptr);
1035 ++rptr;
1036 }
1037
1038 locateParticle(ptemp, pld, 0, finestLevel(), 0);
1039
1040 std::pair<int, int> ind(grd, pld.m_tile);
1041
1042 host_real_attribs[lev][ind].resize(NumRealComps());
1043 host_int_attribs[lev][ind].resize(NumIntComps());
1044
1045 // add the struct
1046 if constexpr(!ParticleType::is_soa_particle)
1047 {
1048 host_particles[lev][ind].push_back(ptemp);
1049
1050 // add the real...
1051 for (int icomp = 0; icomp < NumRealComps(); icomp++) {
1052 host_real_attribs[lev][ind][icomp].push_back(*rptr);
1053 ++rptr;
1054 }
1055
1056 // ... and int array data
1057 for (int icomp = 0; icomp < NumIntComps(); icomp++) {
1058 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1059 ++iptr;
1060 }
1061 } else {
1062 host_particles[lev][ind];
1063
1064 for (int j = 0; j < AMREX_SPACEDIM; j++) {
1065 host_real_attribs[lev][ind][j].push_back(ptemp.pos(j));
1066 }
1067
1068 host_idcpu[lev][ind].push_back(ptemp.m_idcpu);
1069
1070 // read all other SoA
1071 // add the real...
1072 for (int icomp = AMREX_SPACEDIM; icomp < NumRealComps(); icomp++) {
1073 host_real_attribs[lev][ind][icomp].push_back(*rptr);
1074 ++rptr;
1075 }
1076
1077 // ... and int array data
1078 for (int icomp = 0; icomp < NumIntComps(); icomp++) {
1079 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1080 ++iptr;
1081 }
1082 }
1083 }
1084
1085 for (int host_lev = 0; host_lev < std::ssize(host_particles); ++host_lev)
1086 {
1087 for (auto& kv : host_particles[host_lev]) {
1088 auto grid = kv.first.first;
1089 auto tile = kv.first.second;
1090 const auto& src_tile = kv.second;
1091
1092 auto& dst_tile = DefineAndReturnParticleTile(host_lev, grid, tile);
1093 auto old_size = dst_tile.size();
1094 auto new_size = old_size;
1095 if constexpr(!ParticleType::is_soa_particle)
1096 {
1097 new_size += src_tile.size();
1098 } else {
1099 amrex::ignore_unused(src_tile);
1100 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].size();
1101 }
1102 dst_tile.resize(new_size);
1103
1104 if constexpr(!ParticleType::is_soa_particle)
1105 {
1106 Gpu::copyAsync(Gpu::hostToDevice, src_tile.begin(), src_tile.end(),
1107 dst_tile.GetArrayOfStructs().begin() + old_size);
1108 } else {
1110 host_idcpu[host_lev][std::make_pair(grid,tile)].begin(),
1111 host_idcpu[host_lev][std::make_pair(grid,tile)].end(),
1112 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1113 }
1114
1115 for (int i = 0; i < NumRealComps(); ++i) { // NOLINT(readability-misleading-indentation)
1117 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1118 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1119 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1120 }
1121
1122 for (int i = 0; i < NumIntComps(); ++i) {
1124 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1125 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1126 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
1127 }
1128 }
1129 }
1130
1132}
1133
1134template <typename ParticleType, int NArrayReal, int NArrayInt,
1135 template<class> class Allocator, class CellAssignor>
1136void
1137ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1138::WriteAsciiFile (const std::string& filename)
1139{
1140 BL_PROFILE("ParticleContainer::WriteAsciiFile()");
1141 AMREX_ASSERT(!filename.empty());
1142
1143 const auto strttime = amrex::second();
1144 //
1145 // Count # of valid particles.
1146 //
1147 Long nparticles = 0;
1148
1149 for (int lev = 0; lev < std::ssize(m_particles); lev++) {
1150 auto& pmap = m_particles[lev];
1151 for (const auto& kv : pmap) {
1152 const auto& aos = kv.second.GetArrayOfStructs();
1153 auto np = aos.numParticles();
1154 Gpu::HostVector<ParticleType> host_aos(np);
1155 Gpu::copy(Gpu::deviceToHost, aos.begin(), aos.begin() + np, host_aos.begin());
1156 for (int k = 0; k < np; ++k) {
1157 const ParticleType& p = host_aos[k];
1158 if (p.id().is_valid()) {
1159 //
1160 // Only count (and checkpoint) valid particles.
1161 //
1162 nparticles++;
1163 }
1164 }
1165 }
1166 }
1167
1168 //
1169 // And send count to I/O processor.
1170 //
1172
1174 {
1175 //
1176 // Have I/O processor open file and write out particle metadata.
1177 //
1178 std::ofstream File;
1179
1180 File.open(filename.c_str(), std::ios::out|std::ios::trunc);
1181
1182 if (!File.good()) {
1183 amrex::FileOpenFailed(filename);
1184 }
1185
1186 File << nparticles << '\n';
1187 File << NStructReal << '\n';
1188 File << NStructInt << '\n';
1189 File << NumRealComps() << '\n';
1190 File << NumIntComps() << '\n';
1191
1192 File.flush();
1193
1194 File.close();
1195
1196 if (!File.good()) {
1197 amrex::Abort("ParticleContainer::WriteAsciiFile(): problem writing file");
1198 }
1199 }
1200
1202
1203 const int MyProc = ParallelDescriptor::MyProc();
1204
1205 for (int proc = 0; proc < ParallelDescriptor::NProcs(); proc++)
1206 {
1207 if (MyProc == proc)
1208 {
1209 //
1210 // Each CPU opens the file for appending and adds its particles.
1211 //
1213
1214 std::ofstream File;
1215
1216 File.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
1217
1218 File.open(filename.c_str(), std::ios::out|std::ios::app);
1219
1220 File.precision(15);
1221
1222 if (!File.good()) {
1223 amrex::FileOpenFailed(filename);
1224 }
1225
1226 for (int lev = 0; lev < std::ssize(m_particles); lev++) {
1227 auto& pmap = m_particles[lev];
1228 for (const auto& kv : pmap) {
1229 ParticleTile<ParticleType, NArrayReal, NArrayInt,
1231 pinned_ptile.define(NumRuntimeRealComps(), NumRuntimeIntComps(),
1232 nullptr, nullptr, The_Pinned_Arena());
1233 pinned_ptile.resize(kv.second.numParticles());
1234 amrex::copyParticles(pinned_ptile, kv.second);
1235 const auto& host_aos = pinned_ptile.GetArrayOfStructs();
1236 const auto& host_soa = pinned_ptile.GetStructOfArrays();
1237
1238 auto np = host_aos.numParticles();
1239 for (int index = 0; index < np; ++index) {
1240 const ParticleType* it = &host_aos[index];
1241 if (it->id().is_valid()) {
1242
1243 // write out the particle struct first...
1244 AMREX_D_TERM(File << it->pos(0) << ' ',
1245 << it->pos(1) << ' ',
1246 << it->pos(2) << ' ');
1247
1248 for (int i = 0; i < NStructReal; i++) {
1249 File << it->rdata(i) << ' ';
1250 }
1251
1252 File << it->id() << ' ';
1253 File << it->cpu() << ' ';
1254
1255 for (int i = 0; i < NStructInt; i++) {
1256 File << it->idata(i) << ' ';
1257 }
1258
1259 // then the particle attributes.
1260 for (int i = 0; i < NumRealComps(); i++) {
1261 File << host_soa.GetRealData(i)[index] << ' ';
1262 }
1263
1264 for (int i = 0; i < NumIntComps(); i++) {
1265 File << host_soa.GetIntData(i)[index] << ' ';
1266 }
1267
1268 File << '\n';
1269 }
1270 }
1271 }
1272 }
1273
1274 File.flush();
1275
1276 File.close();
1277
1278 if (!File.good()) {
1279 amrex::Abort("ParticleContainer::WriteAsciiFile(): problem writing file");
1280 }
1281
1282 }
1283
1285 }
1286
1287 if (m_verbose > 1)
1288 {
1289 auto stoptime = amrex::second() - strttime;
1290
1292
1293 amrex::Print() << "ParticleContainer::WriteAsciiFile() time: " << stoptime << '\n';
1294 }
1295}
1296
1297}
1298
1299#endif /*AMREX_PARTICLEIO_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
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
static const IntDescriptor & NativeIntDescriptor()
Returns a constant reference to an IntDescriptor describing the native "int" under which AMReX was co...
Definition AMReX_FPC.cpp:76
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
const std::string & FileName() const
Definition AMReX_NFiles.H:160
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
iterator begin() noexcept
Definition AMReX_PODVector.H:674
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:349
int query(std::string_view name, bool &ref, int ival=FIRST) const
Same as querykth() but searches for the last occurrence of name.
Definition AMReX_ParmParse.cpp:1947
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:149
T_ParticleType ParticleType
Definition AMReX_ParticleContainer.H:151
Definition AMReX_GpuAllocators.H:156
This class provides the user with a few print options.
Definition AMReX_Print.H:35
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
T * dataPtr() noexcept
get access to the underlying data pointer
Definition AMReX_Vector.H:49
Long size() const noexcept
Definition AMReX_Vector.H:53
static constexpr int IO_Buffer_Size
We try to do I/O with buffers of this size.
Definition AMReX_VisMFBuffer.H:16
static Long FileOffset(std::ostream &os)
The file offset of the passed ostream.
Definition AMReX_VisMF.cpp:580
static bool GetNoFlushAfterWrite()
Definition AMReX_VisMF.H:284
amrex_particle_real ParticleReal
Floating Point Type for Particles.
Definition AMReX_REAL.H:90
amrex_long Long
Definition AMReX_INT.H:30
void WritePlotFile(const std::string &dir, const std::string &name) const
This version of WritePlotFile writes all components and assigns component names.
Definition AMReX_ParticleIO.H:115
Arena * The_Pinned_Arena()
Definition AMReX_Arena.cpp:860
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
void Barrier(const std::string &)
Definition AMReX_ParallelDescriptor.cpp:1215
void ReduceIntSum(int &)
Definition AMReX_ParallelDescriptor.cpp:1265
void ReadAndBcastFile(const std::string &filename, Vector< char > &charBuf, bool bExitOnError, const MPI_Comm &comm)
Definition AMReX_ParallelDescriptor.cpp:1495
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 UseAsyncOut()
Definition AMReX_AsyncOut.cpp:70
bool Remove(std::string const &filename)
Remove a file, symbolic link, or empty directory.
Definition AMReX_FileSystem.cpp:41
void copy(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:128
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 DeviceToHost deviceToHost
Definition AMReX_GpuContainers.H:106
static constexpr HostToDevice hostToDevice
Definition AMReX_GpuContainers.H:105
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
void ReduceRealMax(Vector< std::reference_wrapper< Real > > const &)
Definition AMReX_ParallelDescriptor.cpp:1228
Definition AMReX_Amr.cpp:50
void writeIntData(const From *data, std::size_t size, std::ostream &os, const amrex::IntDescriptor &id)
Definition AMReX_IntConv.H:23
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
void WriteBinaryParticleDataSync(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, F const &f, bool is_checkpoint)
Definition AMReX_WriteBinaryParticleData.H:485
void FileOpenFailed(const std::string &file)
Output a message and abort when couldn't open the file.
Definition AMReX_Utility.cpp:137
void copyParticles(DstTile &dst, const SrcTile &src) noexcept
Copy particles from src to dst. This version copies all the particles, writing them to the beginning ...
Definition AMReX_ParticleTransformation.H:222
void readFloatData(float *data, std::size_t size, std::istream &is, const RealDescriptor &rd)
Definition AMReX_VectorIO.cpp:136
void writeFloatData(const float *data, std::size_t size, std::ostream &os, const RealDescriptor &rd=FPC::Native32RealDescriptor())
Definition AMReX_VectorIO.cpp:130
bool FileExists(const std::string &filename)
Check if a file already exists. Return true if the filename is an existing file, directory,...
Definition AMReX_Utility.cpp:145
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2006
std::string Concatenate(const std::string &root, int num, int mindigits)
Returns rootNNNN where NNNN == num.
Definition AMReX_String.cpp:35
double second() noexcept
Definition AMReX_Utility.cpp:940
void writeDoubleData(const double *data, std::size_t size, std::ostream &os, const RealDescriptor &rd=FPC::Native64RealDescriptor())
Definition AMReX_VectorIO.cpp:142
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:235
void readDoubleData(double *data, std::size_t size, std::istream &is, const RealDescriptor &rd)
Definition AMReX_VectorIO.cpp:148
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
void WriteBinaryParticleDataAsync(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, bool is_checkpoint)
Definition AMReX_WriteBinaryParticleData.H:784
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2015
void readIntData(To *data, std::size_t size, std::istream &is, const amrex::IntDescriptor &id)
Definition AMReX_IntConv.H:38
Definition AMReX_ParticleIO.H:40
__host__ __device__ int operator()(const P &p) const
Definition AMReX_ParticleIO.H:43
FabArray memory allocation information.
Definition AMReX_FabArray.H:67
MFInfo & SetAlloc(bool a) noexcept
Definition AMReX_FabArray.H:74
uint64_t m_idcpu
Definition AMReX_Particle.H:359
__host__ __device__ bool is_valid() const noexcept
Definition AMReX_Particle.H:252
A struct used for storing a particle's position in the AMR hierarchy.
Definition AMReX_ParticleContainer.H:93
int m_tile
Definition AMReX_ParticleContainer.H:96
Definition AMReX_ParticleTile.H:764
The struct used to store particles.
Definition AMReX_Particle.H:405
__host__ __device__ RealVect pos() const &
Definition AMReX_Particle.H:456
__host__ __device__ int & idata(int index) &
Definition AMReX_Particle.H:545
__host__ __device__ ParticleCPUWrapper cpu() &
Definition AMReX_Particle.H:424
__host__ __device__ ParticleIDWrapper id() &
Definition AMReX_Particle.H:427
__host__ __device__ RealType & rdata(int index) &
Definition AMReX_Particle.H:474