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