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