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