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 + NArrayInt );
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, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>*>
267void
269::WritePlotFile (const std::string& dir, const std::string& name, F&& 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,
291 std::forward<F>(f));
292}
293
295template <typename ParticleType, int NArrayReal, int NArrayInt,
296 template<class> class Allocator, class CellAssignor>
297template <class F>
298void
300::WritePlotFile (const std::string& dir, const std::string& name,
301 const Vector<std::string>& real_comp_names,
302 const Vector<std::string>& int_comp_names, F&& f) const
303{
304 if constexpr(ParticleType::is_soa_particle) {
305 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
306 } else {
307 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal);
308 }
309 AMREX_ASSERT( int_comp_names.size() == NStructInt + NArrayInt );
310
311 Vector<int> write_real_comp;
312 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
313 for (int i = 0; i < nrc; ++i) {
314 write_real_comp.push_back(1);
315 }
316
317 Vector<int> write_int_comp;
318 for (int i = 0; i < NStructInt + NumIntComps(); ++i) {
319 write_int_comp.push_back(1);
320 }
321
322 WriteBinaryParticleData(dir, name,
323 write_real_comp, write_int_comp,
324 real_comp_names, int_comp_names,
325 std::forward<F>(f));
326}
327
329template <typename ParticleType, int NArrayReal, int NArrayInt,
330 template<class> class Allocator, class CellAssignor>
331template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>*>
332void
334::WritePlotFile (const std::string& dir, const std::string& name,
335 const Vector<std::string>& real_comp_names, F&& f) const
336{
337 if constexpr(ParticleType::is_soa_particle) {
338 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
339 } else {
340 AMREX_ALWAYS_ASSERT(real_comp_names.size() == NumRealComps() + NStructReal);
341 }
342
343 Vector<int> write_real_comp;
344 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
345 for (int i = 0; i < nrc; ++i) {
346 write_real_comp.push_back(1);
347 }
348
349 Vector<int> write_int_comp;
350 for (int i = 0; i < NStructInt + NumIntComps(); ++i) {
351 write_int_comp.push_back(1);
352 }
353
354 Vector<std::string> int_comp_names;
355 for (int i = 0; i < NStructInt + NumIntComps(); ++i )
356 {
357 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
358 }
359
360 WriteBinaryParticleData(dir, name,
361 write_real_comp, write_int_comp,
362 real_comp_names, int_comp_names,
363 std::forward<F>(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&& 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,
399 std::forward<F>(f));
400}
401
403template <typename ParticleType, int NArrayReal, int NArrayInt,
404 template<class> class Allocator, class CellAssignor>
405template <class F>
406void
408WritePlotFile (const std::string& dir, const std::string& name,
409 const Vector<int>& write_real_comp,
410 const Vector<int>& write_int_comp,
411 const Vector<std::string>& real_comp_names,
412 const Vector<std::string>& int_comp_names,
413 F&& f) const
414{
415 BL_PROFILE("ParticleContainer::WritePlotFile()");
416
417 WriteBinaryParticleData(dir, name,
418 write_real_comp, write_int_comp,
419 real_comp_names, int_comp_names,
420 std::forward<F>(f));
421}
422
423template <typename ParticleType, int NArrayReal, int NArrayInt,
424 template<class> class Allocator, class CellAssignor>
425template <class F>
426void
428::WriteBinaryParticleData (const std::string& dir, const std::string& name,
429 const Vector<int>& write_real_comp,
430 const Vector<int>& write_int_comp,
431 const Vector<std::string>& real_comp_names,
432 const Vector<std::string>& int_comp_names,
433 F&& f, bool is_checkpoint) const
434{
435 if (AsyncOut::UseAsyncOut()) {
436 WriteBinaryParticleDataAsync(*this, dir, name,
437 write_real_comp, write_int_comp,
438 real_comp_names, int_comp_names, is_checkpoint);
439 } else
440 {
441 WriteBinaryParticleDataSync(*this, dir, name,
442 write_real_comp, write_int_comp,
443 real_comp_names, int_comp_names,
444 std::forward<F>(f), is_checkpoint);
445 }
446}
447
448template <typename ParticleType, int NArrayReal, int NArrayInt,
449 template<class> class Allocator, class CellAssignor>
450void
453{
454 if( ! usePrePost) {
455 return;
456 }
457
458 BL_PROFILE("ParticleContainer::CheckpointPre()");
459
460 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
461 Long nparticles = 0;
462 Long maxnextid = ParticleType::NextID();
463
464 for (int lev = 0; lev < m_particles.size(); lev++) {
465 const auto& pmap = m_particles[lev];
466 for (const auto& kv : pmap) {
467 const auto& aos = kv.second.GetArrayOfStructs();
468 for (int k = 0; k < aos.numParticles(); ++k) {
469 const ParticleType& p = aos[k];
470 if (p.id().is_valid()) {
471 //
472 // Only count (and checkpoint) valid particles.
473 //
474 nparticles++;
475 }
476 }
477 }
478 }
479 ParallelDescriptor::ReduceLongSum(nparticles, IOProcNumber);
480
481 ParticleType::NextID(maxnextid);
482 ParallelDescriptor::ReduceLongMax(maxnextid, IOProcNumber);
483
484 nparticlesPrePost = nparticles;
485 maxnextidPrePost = int(maxnextid);
486
487 nParticlesAtLevelPrePost.clear();
488 nParticlesAtLevelPrePost.resize(finestLevel() + 1, 0);
489 for(int lev(0); lev <= finestLevel(); ++lev) {
490 nParticlesAtLevelPrePost[lev] = NumberOfParticlesAtLevel(lev);
491 }
492
493 whichPrePost.clear();
494 whichPrePost.resize(finestLevel() + 1);
495 countPrePost.clear();
496 countPrePost.resize(finestLevel() + 1);
497 wherePrePost.clear();
498 wherePrePost.resize(finestLevel() + 1);
499
500 filePrefixPrePost.clear();
501 filePrefixPrePost.resize(finestLevel() + 1);
502}
503
504
505template <typename ParticleType, int NArrayReal, int NArrayInt,
506 template<class> class Allocator, class CellAssignor>
507void
510{
511 if( ! usePrePost) {
512 return;
513 }
514
515 BL_PROFILE("ParticleContainer::CheckpointPost()");
516
517 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
518 std::ofstream HdrFile;
519 HdrFile.open(HdrFileNamePrePost.c_str(), std::ios::out | std::ios::app);
520
521 for(int lev(0); lev <= finestLevel(); ++lev) {
522 ParallelDescriptor::ReduceIntSum (whichPrePost[lev].dataPtr(), whichPrePost[lev].size(), IOProcNumber);
523 ParallelDescriptor::ReduceIntSum (countPrePost[lev].dataPtr(), countPrePost[lev].size(), IOProcNumber);
524 ParallelDescriptor::ReduceLongSum(wherePrePost[lev].dataPtr(), wherePrePost[lev].size(), IOProcNumber);
525
526
528 for(int j(0); j < whichPrePost[lev].size(); ++j) {
529 HdrFile << whichPrePost[lev][j] << ' ' << countPrePost[lev][j] << ' ' << wherePrePost[lev][j] << '\n';
530 }
531
532 const bool gotsome = (nParticlesAtLevelPrePost[lev] > 0);
533 if(gotsome && doUnlink) {
534// BL_PROFILE_VAR("PC<NNNN>::Checkpoint:unlink", unlink_post);
535 // Unlink any zero-length data files.
536 Vector<Long> cnt(nOutFilesPrePost,0);
537
538 for(int i(0), N = countPrePost[lev].size(); i < N; ++i) {
539 cnt[whichPrePost[lev][i]] += countPrePost[lev][i];
540 }
541
542 for(int i(0), N = int(cnt.size()); i < N; ++i) {
543 if(cnt[i] == 0) {
544 std::string FullFileName = NFilesIter::FileName(i, filePrefixPrePost[lev]);
545 FileSystem::Remove(FullFileName);
546 }
547 }
548 }
549 }
550 }
551
553 HdrFile.flush();
554 HdrFile.close();
555 if( ! HdrFile.good()) {
556 amrex::Abort("ParticleContainer::CheckpointPost(): problem writing HdrFile");
557 }
558 }
559}
560
561template <typename ParticleType, int NArrayReal, int NArrayInt,
562 template<class> class Allocator, class CellAssignor>
563void
566{
567 CheckpointPre();
568}
569
570
571template <typename ParticleType, int NArrayReal, int NArrayInt,
572 template<class> class Allocator, class CellAssignor>
573void
576{
577 CheckpointPost();
578}
579
580
581template <typename ParticleType, int NArrayReal, int NArrayInt,
582 template<class> class Allocator, class CellAssignor>
583void
585::WriteParticles (int lev, std::ofstream& ofs, int fnum,
586 Vector<int>& which, Vector<int>& count, Vector<Long>& where,
587 const Vector<int>& write_real_comp,
588 const Vector<int>& write_int_comp,
589 const Vector<std::map<std::pair<int, int>, IntVector>>& particle_io_flags,
590 bool is_checkpoint) const
591{
592 BL_PROFILE("ParticleContainer::WriteParticles()");
593
594 // For a each grid, the tiles it contains
595 std::map<int, Vector<int> > tile_map;
596
597 for (const auto& kv : m_particles[lev])
598 {
599 const int grid = kv.first.first;
600 const int tile = kv.first.second;
601 tile_map[grid].push_back(tile);
602 const auto& pflags = particle_io_flags[lev].at(kv.first);
603
604 // Only write out valid particles.
605 count[grid] += particle_detail::countFlags(pflags);
606 }
607
608 MFInfo info;
609 info.SetAlloc(false);
610 MultiFab state(ParticleBoxArray(lev), ParticleDistributionMap(lev), 1,0,info);
611
612 for (MFIter mfi(state); mfi.isValid(); ++mfi)
613 {
614 const int grid = mfi.index();
615
616 which[grid] = fnum;
617 where[grid] = VisMF::FileOffset(ofs);
618
619 if (count[grid] == 0) { continue; }
620
621 Vector<int> istuff;
623 particle_detail::packIOData(istuff, rstuff, *this, lev, grid,
624 write_real_comp, write_int_comp,
625 particle_io_flags, tile_map[grid], count[grid], is_checkpoint);
626
627 writeIntData(istuff.dataPtr(), istuff.size(), ofs);
629 ofs.flush(); // Some systems require this flush() (probably due to a bug)
630 }
631
632 WriteParticleRealData(rstuff.dataPtr(), rstuff.size(), ofs);
634 ofs.flush(); // Some systems require this flush() (probably due to a bug)
635 }
636 }
637}
638
639
640template <typename ParticleType, int NArrayReal, int NArrayInt,
641 template<class> class Allocator, class CellAssignor>
642void
644::Restart (const std::string& dir, const std::string& file, bool /*is_checkpoint*/)
645{
646 Restart(dir, file);
647}
648
649template <typename ParticleType, int NArrayReal, int NArrayInt,
650 template<class> class Allocator, class CellAssignor>
651void
653::Restart (const std::string& dir, const std::string& file)
654{
655 BL_PROFILE("ParticleContainer::Restart()");
656 AMREX_ASSERT(!dir.empty());
657 AMREX_ASSERT(!file.empty());
658
659 const auto strttime = amrex::second();
660
661 int DATA_Digits_Read(5);
662 ParmParse pp("particles");
663 pp.query("datadigits_read",DATA_Digits_Read);
664
665 std::string fullname = dir;
666 if (!fullname.empty() && fullname[fullname.size()-1] != '/') {
667 fullname += '/';
668 }
669 fullname += file;
670 std::string HdrFileName = fullname;
671 if (!HdrFileName.empty() && HdrFileName[HdrFileName.size()-1] != '/') {
672 HdrFileName += '/';
673 }
674 HdrFileName += "Header";
675
676 Vector<char> fileCharPtr;
677 ParallelDescriptor::ReadAndBcastFile(HdrFileName, fileCharPtr);
678 std::string fileCharPtrString(fileCharPtr.dataPtr());
679 std::istringstream HdrFile(fileCharPtrString, std::istringstream::in);
680
681 std::string version;
682 HdrFile >> version;
683 AMREX_ASSERT(!version.empty());
684
685 // What do our version strings mean?
686 // "Version_One_Dot_Zero" -- hard-wired to write out in double precision.
687 // "Version_One_Dot_One" -- can write out either as either single or double precision.
688 // Appended to the latter version string are either "_single" or "_double" to
689 // indicate how the particles were written.
690 // "Version_Two_Dot_Zero" -- this is the AMReX particle file format
691 // "Version_Two_Dot_One" -- expanded particle ids to allow for 2**39-1 per proc
692 std::string how;
693 bool convert_ids = false;
694 if (version.find("Version_Two_Dot_One") != std::string::npos) {
695 convert_ids = true;
696 }
697 if (version.find("Version_One_Dot_Zero") != std::string::npos) {
698 how = "double";
699 }
700 else if (version.find("Version_One_Dot_One") != std::string::npos ||
701 version.find("Version_Two_Dot_Zero") != std::string::npos ||
702 version.find("Version_Two_Dot_One") != std::string::npos) {
703 if (version.find("_single") != std::string::npos) {
704 how = "single";
705 }
706 else if (version.find("_double") != std::string::npos) {
707 how = "double";
708 }
709 else {
710 std::string msg("ParticleContainer::Restart(): bad version string: ");
711 msg += version;
712 amrex::Error(version.c_str());
713 }
714 }
715 else {
716 std::string msg("ParticleContainer::Restart(): unknown version string: ");
717 msg += version;
718 amrex::Abort(msg.c_str());
719 }
720
721 int dm;
722 HdrFile >> dm;
723 if (dm != AMREX_SPACEDIM) {
724 amrex::Abort("ParticleContainer::Restart(): dm != AMREX_SPACEDIM");
725 }
726
727 int nr;
728 HdrFile >> nr;
729 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
730 if (nr != nrc) {
731 amrex::Abort("ParticleContainer::Restart(): nr not the expected value");
732 }
733
734 std::string comp_name;
735 for (int i = 0; i < nr; ++i) {
736 HdrFile >> comp_name;
737 }
738
739 int ni;
740 HdrFile >> ni;
741 if (ni != NStructInt + NumIntComps()) {
742 amrex::Abort("ParticleContainer::Restart(): ni != NStructInt");
743 }
744
745 for (int i = 0; i < ni; ++i) {
746 HdrFile >> comp_name;
747 }
748
749 bool checkpoint;
750 HdrFile >> checkpoint;
751
752 Long nparticles;
753 HdrFile >> nparticles;
754 AMREX_ASSERT(nparticles >= 0);
755
756 Long maxnextid;
757 HdrFile >> maxnextid;
758 AMREX_ASSERT(maxnextid > 0);
759 ParticleType::NextID(maxnextid);
760
761 int finest_level_in_file;
762 HdrFile >> finest_level_in_file;
763 AMREX_ASSERT(finest_level_in_file >= 0);
764
765 // Determine whether this is a dual-grid restart or not.
766 Vector<DistributionMapping> old_dms(finest_level_in_file + 1);
767 Vector<BoxArray> old_bas(finest_level_in_file + 1);
768 Vector<BoxArray> particle_box_arrays(finest_level_in_file + 1);
769 bool dual_grid = false;
770
771 bool have_pheaders = false;
772 for (int lev = 0; lev <= finest_level_in_file; lev++)
773 {
774 std::string phdr_name = fullname;
775 phdr_name += "/Level_";
776 phdr_name = amrex::Concatenate(phdr_name, lev, 1);
777 phdr_name += "/Particle_H";
778
779 if (amrex::FileExists(phdr_name)) {
780 have_pheaders = true;
781 break;
782 }
783 }
784
785 if (have_pheaders)
786 {
787 for (int lev = 0; lev <= finestLevel(); lev++)
788 {
789 old_dms[lev] = ParticleDistributionMap(lev);
790 old_bas[lev] = ParticleBoxArray(lev);
791 std::string phdr_name = fullname;
792 phdr_name += "/Level_";
793 phdr_name = amrex::Concatenate(phdr_name, lev, 1);
794 phdr_name += "/Particle_H";
795
796 if (! amrex::FileExists(phdr_name)) { continue; }
797
798 Vector<char> phdr_chars;
799 ParallelDescriptor::ReadAndBcastFile(phdr_name, phdr_chars);
800 std::string phdr_string(phdr_chars.dataPtr());
801 std::istringstream phdr_file(phdr_string, std::istringstream::in);
802
803 particle_box_arrays[lev].readFrom(phdr_file);
804 if (! particle_box_arrays[lev].CellEqual(ParticleBoxArray(lev))) { dual_grid = true; }
805 }
806 } else // if no particle box array information exists in the file, we assume a single grid restart
807 {
808 dual_grid = false;
809 }
810
811 if (dual_grid) {
812 for (int lev = 0; lev <= finestLevel(); lev++) {
813 // this can happen if there are no particles at a given level in the checkpoint
814 if (particle_box_arrays[lev].empty()) {
815 particle_box_arrays[lev] = BoxArray(Geom(lev).Domain());
817 SetParticleBoxArray(lev, particle_box_arrays[lev]);
818 DistributionMapping pdm(particle_box_arrays[lev]);
819 SetParticleDistributionMap(lev, pdm);
820 }
821 }
822
823 Vector<int> ngrids(finest_level_in_file+1);
824 for (int lev = 0; lev <= finest_level_in_file; lev++) {
825 HdrFile >> ngrids[lev];
826 AMREX_ASSERT(ngrids[lev] > 0);
827 }
828
829 resizeData();
830
831 if (finest_level_in_file > finestLevel()) {
832 m_particles.resize(finest_level_in_file+1);
833 }
834
835 for (int lev = 0; lev <= finest_level_in_file; lev++) {
836 Vector<int> which(ngrids[lev]);
837 Vector<int> count(ngrids[lev]);
838 Vector<Long> where(ngrids[lev]);
839 for (int i = 0; i < ngrids[lev]; i++) {
840 HdrFile >> which[i] >> count[i] >> where[i];
841 }
842
843 Vector<int> grids_to_read;
844 if (lev <= finestLevel()) {
845 for (MFIter mfi(*m_dummy_mf[lev]); mfi.isValid(); ++mfi) {
846 grids_to_read.push_back(mfi.index());
847 }
848 } else {
849
850 // we lost a level on restart. we still need to read in particles
851 // on finer levels, and put them in the right place via Redistribute()
852
853 const int rank = ParallelDescriptor::MyProc();
854 const int NReaders = MaxReaders();
855 if (rank >= NReaders) { return; }
856
857 const int Navg = ngrids[lev] / NReaders;
858 const int Nleft = ngrids[lev] - Navg * NReaders;
859
860 int lo, hi;
861 if (rank < Nleft) {
862 lo = rank*(Navg + 1);
863 hi = lo + Navg + 1;
864 }
865 else {
866 lo = rank * Navg + Nleft;
867 hi = lo + Navg;
868 }
869
870 for (int i = lo; i < hi; ++i) {
871 grids_to_read.push_back(i);
872 }
873 }
874
875 for(int grid : grids_to_read) {
876 if (count[grid] <= 0) { continue; }
877
878 // The file names in the header file are relative.
879 std::string name = fullname;
880
881 if (!name.empty() && name[name.size()-1] != '/') {
882 name += '/';
883 }
884
885 name += "Level_";
886 name += amrex::Concatenate("", lev, 1);
887 name += '/';
888 name += DataPrefix();
889 name += amrex::Concatenate("", which[grid], DATA_Digits_Read);
890
891 std::ifstream ParticleFile;
892
893 ParticleFile.open(name.c_str(), std::ios::in | std::ios::binary);
894
895 if (!ParticleFile.good()) {
898
899 ParticleFile.seekg(where[grid], std::ios::beg);
900
901 // Use if constexpr to avoid instantiating the mis-matched
902 // type case and triggering the static_assert on the
903 // underlying copy calls
904 if (how == "single") {
905 if constexpr (std::is_same_v<ParticleReal, float>) {
906 ReadParticles<float>(count[grid], grid, lev, ParticleFile, finest_level_in_file, convert_ids);
907 } else {
908 amrex::Error("File contains single-precision data, while AMReX is compiled with ParticleReal==double");
909 }
910 }
911 else if (how == "double") {
912 if constexpr (std::is_same_v<ParticleReal, double>) {
913 ReadParticles<double>(count[grid], grid, lev, ParticleFile, finest_level_in_file, convert_ids);
914 } else {
915 amrex::Error("File contains double-precision data, while AMReX is compiled with ParticleReal==float");
916 }
917 }
918 else {
919 std::string msg("ParticleContainer::Restart(): bad parameter: ");
920 msg += how;
921 amrex::Error(msg.c_str());
922 }
923
924 ParticleFile.close();
925
926 if (!ParticleFile.good()) {
927 amrex::Abort("ParticleContainer::Restart(): problem reading particles");
928 }
929 }
930 }
931
932 if (dual_grid) {
933 for (int lev = 0; lev <= finestLevel(); lev++) {
934 SetParticleBoxArray(lev, old_bas[lev]);
935 SetParticleDistributionMap(lev, old_dms[lev]);
936 }
937 }
938
939 Redistribute();
940
941 AMREX_ASSERT(OK());
942
943 if (m_verbose > 1) {
944 auto stoptime = amrex::second() - strttime;
946 amrex::Print() << "ParticleContainer::Restart() time: " << stoptime << '\n';
947 }
948}
949
950// Read a batch of particles from the checkpoint file
951template <typename ParticleType, int NArrayReal, int NArrayInt,
952 template<class> class Allocator, class CellAssignor>
953template <class RTYPE>
954void
955ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
956::ReadParticles (int cnt, int grd, int lev, std::ifstream& ifs,
957 int finest_level_in_file, bool convert_ids)
958{
959 BL_PROFILE("ParticleContainer::ReadParticles()");
960 AMREX_ASSERT(cnt > 0);
961 AMREX_ASSERT(lev < int(m_particles.size()));
962
963 // First read in the integer data in binary. We do not store
964 // the m_lev and m_grid data on disk. We can easily recreate
965 // that given the structure of the checkpoint file.
966 const int iChunkSize = 2 + NStructInt + NumIntComps();
967 Vector<int> istuff(std::size_t(cnt)*iChunkSize);
968 readIntData(istuff.dataPtr(), istuff.size(), ifs, FPC::NativeIntDescriptor());
969
970 // Then the real data in binary.
971 const int rChunkSize = ParticleType::is_soa_particle ? NStructReal + NumRealComps() : AMREX_SPACEDIM + NStructReal + NumRealComps();
972 Vector<RTYPE> rstuff(std::size_t(cnt)*rChunkSize);
973 ReadParticleRealData(rstuff.dataPtr(), rstuff.size(), ifs);
974
975 // Now reassemble the particles.
976 int* iptr = istuff.dataPtr();
977 RTYPE* rptr = rstuff.dataPtr();
978
980 ParticleLocData pld;
981
983 host_particles.reserve(15);
984 host_particles.resize(finest_level_in_file+1);
985
987 std::vector<Gpu::HostVector<RTYPE> > > > host_real_attribs;
988 host_real_attribs.reserve(15);
989 host_real_attribs.resize(finest_level_in_file+1);
992 std::vector<Gpu::HostVector<int> > > > host_int_attribs;
993 host_int_attribs.reserve(15);
994 host_int_attribs.resize(finest_level_in_file+1);
995
997 host_idcpu.reserve(15);
998 host_idcpu.resize(finestLevel()+1);
999
1000 for (int i = 0; i < cnt; i++) {
1001 // note: for pure SoA particle layouts, we do write the id, cpu and positions as a struct
1002 // for backwards compatibility with readers
1003 if (convert_ids) {
1004 std::int32_t xi, yi;
1005 std::uint32_t xu, yu;
1006 xi = iptr[0];
1007 yi = iptr[1];
1008 std::memcpy(&xu, &xi, sizeof(xi));
1009 std::memcpy(&yu, &yi, sizeof(yi));
1010 ptemp.m_idcpu = ((std::uint64_t)xu) << 32 | yu;
1011 } else {
1012 ptemp.id() = iptr[0];
1013 ptemp.cpu() = iptr[1];
1014 }
1015 iptr += 2;
1016
1017 for (int j = 0; j < NStructInt; j++)
1018 {
1019 ptemp.idata(j) = *iptr;
1020 ++iptr;
1021 }
1022
1023 AMREX_ASSERT(ptemp.id().is_valid());
1024
1025 AMREX_D_TERM(ptemp.pos(0) = ParticleReal(rptr[0]);,
1026 ptemp.pos(1) = ParticleReal(rptr[1]);,
1027 ptemp.pos(2) = ParticleReal(rptr[2]););
1028
1029 rptr += AMREX_SPACEDIM;
1030
1031 for (int j = 0; j < NStructReal; j++)
1032 {
1033 ptemp.rdata(j) = ParticleReal(*rptr);
1034 ++rptr;
1035 }
1036
1037 locateParticle(ptemp, pld, 0, finestLevel(), 0);
1038
1039 std::pair<int, int> ind(grd, pld.m_tile);
1040
1041 host_real_attribs[lev][ind].resize(NumRealComps());
1042 host_int_attribs[lev][ind].resize(NumIntComps());
1043
1044 // add the struct
1045 if constexpr(!ParticleType::is_soa_particle)
1046 {
1047 host_particles[lev][ind].push_back(ptemp);
1048
1049 // add the real...
1050 for (int icomp = 0; icomp < NumRealComps(); icomp++) {
1051 host_real_attribs[lev][ind][icomp].push_back(*rptr);
1052 ++rptr;
1053 }
1054
1055 // ... and int array data
1056 for (int icomp = 0; icomp < NumIntComps(); icomp++) {
1057 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1058 ++iptr;
1059 }
1060 } else {
1061 host_particles[lev][ind];
1062
1063 for (int j = 0; j < AMREX_SPACEDIM; j++) {
1064 host_real_attribs[lev][ind][j].push_back(ptemp.pos(j));
1065 }
1066
1067 host_idcpu[lev][ind].push_back(ptemp.m_idcpu);
1068
1069 // read all other SoA
1070 // add the real...
1071 for (int icomp = AMREX_SPACEDIM; icomp < NumRealComps(); icomp++) {
1072 host_real_attribs[lev][ind][icomp].push_back(*rptr);
1073 ++rptr;
1074 }
1075
1076 // ... and int array data
1077 for (int icomp = 0; icomp < NumIntComps(); icomp++) {
1078 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1079 ++iptr;
1080 }
1081 }
1082 }
1083
1084 for (int host_lev = 0; host_lev < static_cast<int>(host_particles.size()); ++host_lev)
1085 {
1086 for (auto& kv : host_particles[host_lev]) {
1087 auto grid = kv.first.first;
1088 auto tile = kv.first.second;
1089 const auto& src_tile = kv.second;
1090
1091 auto& dst_tile = DefineAndReturnParticleTile(host_lev, grid, tile);
1092 auto old_size = dst_tile.size();
1093 auto new_size = old_size;
1094 if constexpr(!ParticleType::is_soa_particle)
1095 {
1096 new_size += src_tile.size();
1097 } else {
1098 amrex::ignore_unused(src_tile);
1099 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].size();
1100 }
1101 dst_tile.resize(new_size);
1102
1103 if constexpr(!ParticleType::is_soa_particle)
1104 {
1105 Gpu::copyAsync(Gpu::hostToDevice, src_tile.begin(), src_tile.end(),
1106 dst_tile.GetArrayOfStructs().begin() + old_size);
1107 } else {
1109 host_idcpu[host_lev][std::make_pair(grid,tile)].begin(),
1110 host_idcpu[host_lev][std::make_pair(grid,tile)].end(),
1111 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1112 }
1113
1114 for (int i = 0; i < NumRealComps(); ++i) { // NOLINT(readability-misleading-indentation)
1116 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1117 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1118 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1119 }
1120
1121 for (int i = 0; i < NumIntComps(); ++i) {
1123 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1124 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1125 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
1126 }
1127 }
1128 }
1129
1131}
1132
1133template <typename ParticleType, int NArrayReal, int NArrayInt,
1134 template<class> class Allocator, class CellAssignor>
1135void
1136ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1137::WriteAsciiFile (const std::string& filename)
1138{
1139 BL_PROFILE("ParticleContainer::WriteAsciiFile()");
1140 AMREX_ASSERT(!filename.empty());
1141
1142 const auto strttime = amrex::second();
1143 //
1144 // Count # of valid particles.
1145 //
1146 Long nparticles = 0;
1147
1148 for (int lev = 0; lev < m_particles.size(); lev++) {
1149 auto& pmap = m_particles[lev];
1150 for (const auto& kv : pmap) {
1151 const auto& aos = kv.second.GetArrayOfStructs();
1152 auto np = aos.numParticles();
1153 Gpu::HostVector<ParticleType> host_aos(np);
1154 Gpu::copy(Gpu::deviceToHost, aos.begin(), aos.begin() + np, host_aos.begin());
1155 for (int k = 0; k < np; ++k) {
1156 const ParticleType& p = host_aos[k];
1157 if (p.id().is_valid()) {
1158 //
1159 // Only count (and checkpoint) valid particles.
1160 //
1161 nparticles++;
1162 }
1163 }
1164 }
1165 }
1166
1167 //
1168 // And send count to I/O processor.
1169 //
1171
1173 {
1174 //
1175 // Have I/O processor open file and write out particle metadata.
1176 //
1177 std::ofstream File;
1178
1179 File.open(filename.c_str(), std::ios::out|std::ios::trunc);
1180
1181 if (!File.good()) {
1182 amrex::FileOpenFailed(filename);
1183 }
1184
1185 File << nparticles << '\n';
1186 File << NStructReal << '\n';
1187 File << NStructInt << '\n';
1188 File << NumRealComps() << '\n';
1189 File << NumIntComps() << '\n';
1190
1191 File.flush();
1192
1193 File.close();
1194
1195 if (!File.good()) {
1196 amrex::Abort("ParticleContainer::WriteAsciiFile(): problem writing file");
1197 }
1198 }
1199
1201
1202 const int MyProc = ParallelDescriptor::MyProc();
1203
1204 for (int proc = 0; proc < ParallelDescriptor::NProcs(); proc++)
1205 {
1206 if (MyProc == proc)
1207 {
1208 //
1209 // Each CPU opens the file for appending and adds its particles.
1210 //
1212
1213 std::ofstream File;
1214
1215 File.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
1216
1217 File.open(filename.c_str(), std::ios::out|std::ios::app);
1218
1219 File.precision(15);
1220
1221 if (!File.good()) {
1222 amrex::FileOpenFailed(filename);
1223 }
1224
1225 for (int lev = 0; lev < m_particles.size(); lev++) {
1226 auto& pmap = m_particles[lev];
1227 for (const auto& kv : pmap) {
1228 ParticleTile<ParticleType, NArrayReal, NArrayInt,
1230 pinned_ptile.define(NumRuntimeRealComps(), NumRuntimeIntComps(),
1231 nullptr, nullptr, The_Pinned_Arena());
1232 pinned_ptile.resize(kv.second.numParticles());
1233 amrex::copyParticles(pinned_ptile, kv.second);
1234 const auto& host_aos = pinned_ptile.GetArrayOfStructs();
1235 const auto& host_soa = pinned_ptile.GetStructOfArrays();
1236
1237 auto np = host_aos.numParticles();
1238 for (int index = 0; index < np; ++index) {
1239 const ParticleType* it = &host_aos[index];
1240 if (it->id().is_valid()) {
1241
1242 // write out the particle struct first...
1243 AMREX_D_TERM(File << it->pos(0) << ' ',
1244 << it->pos(1) << ' ',
1245 << it->pos(2) << ' ');
1246
1247 for (int i = 0; i < NStructReal; i++) {
1248 File << it->rdata(i) << ' ';
1249 }
1250
1251 File << it->id() << ' ';
1252 File << it->cpu() << ' ';
1253
1254 for (int i = 0; i < NStructInt; i++) {
1255 File << it->idata(i) << ' ';
1256 }
1257
1258 // then the particle attributes.
1259 for (int i = 0; i < NumRealComps(); i++) {
1260 File << host_soa.GetRealData(i)[index] << ' ';
1261 }
1262
1263 for (int i = 0; i < NumIntComps(); i++) {
1264 File << host_soa.GetIntData(i)[index] << ' ';
1265 }
1266
1267 File << '\n';
1268 }
1269 }
1270 }
1271 }
1272
1273 File.flush();
1274
1275 File.close();
1276
1277 if (!File.good()) {
1278 amrex::Abort("ParticleContainer::WriteAsciiFile(): problem writing file");
1279 }
1280
1281 }
1282
1284 }
1285
1286 if (m_verbose > 1)
1287 {
1288 auto stoptime = amrex::second() - strttime;
1289
1291
1292 amrex::Print() << "ParticleContainer::WriteAsciiFile() time: " << stoptime << '\n';
1293 }
1294}
1295
1296}
1297
1298#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:567
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
static const IntDescriptor & NativeIntDescriptor()
Returns a constant reference to an IntDescriptor describing the native "int" under which AMReX was co...
Definition AMReX_FPC.cpp:76
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:63
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:147
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:345
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:1395
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:148
typename SoA::IntVector IntVector
Definition AMReX_ParticleContainer.H:196
T_ParticleType ParticleType
Definition AMReX_ParticleContainer.H:150
Definition AMReX_GpuAllocators.H:144
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
@ IO_Buffer_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
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:126
void Barrier(const std::string &)
Definition AMReX_ParallelDescriptor.cpp:1211
void ReduceIntSum(int &)
Definition AMReX_ParallelDescriptor.cpp:1261
void ReadAndBcastFile(const std::string &filename, Vector< char > &charBuf, bool bExitOnError, const MPI_Comm &comm)
Definition AMReX_ParallelDescriptor.cpp:1491
int NProcs() noexcept
Definition AMReX_ParallelDescriptor.H:246
void ReduceLongSum(Long &)
Definition AMReX_ParallelDescriptor.cpp:1232
void ReduceLongMax(Long &)
Definition AMReX_ParallelDescriptor.cpp:1233
bool IOProcessor() noexcept
Is this CPU the I/O Processor? To get the rank number, call IOProcessorNumber()
Definition AMReX_ParallelDescriptor.H:281
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:263
std::string const & name()
Definition AMReX_Machine.cpp:46
int IOProcessorNumber() noexcept
Definition AMReX_ParallelDescriptor.H:271
void ReduceRealMax(Vector< std::reference_wrapper< Real > > const &)
Definition AMReX_ParallelDescriptor.cpp:1224
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:138
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:467
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:161
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:224
Arena * The_Pinned_Arena()
Definition AMReX_Arena.cpp:823
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:230
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:765
__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:92
int m_tile
Definition AMReX_ParticleContainer.H:95
Definition AMReX_ParticleTile.H:721
The struct used to store particles.
Definition AMReX_Particle.H:404
__host__ __device__ RealVect pos() const &
Definition AMReX_Particle.H:454
__host__ __device__ int & idata(int index) &
Definition AMReX_Particle.H:543
__host__ __device__ ParticleCPUWrapper cpu() &
Definition AMReX_Particle.H:422
__host__ __device__ ParticleIDWrapper id() &
Definition AMReX_Particle.H:425
__host__ __device__ RealType & rdata(int index) &
Definition AMReX_Particle.H:472