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