Block-Structured AMR Software Framework
AMReX_ParticleContainer.H
Go to the documentation of this file.
1 #ifndef AMREX_PARTICLE_CONTAINER_H_
2 #define AMREX_PARTICLE_CONTAINER_H_
3 #include <AMReX_Config.H>
4 
6 #include <AMReX_ParmParse.H>
7 #include <AMReX_ParGDB.H>
8 #include <AMReX_REAL.H>
9 #include <AMReX_IntVect.H>
10 #include <AMReX_Array.H>
11 #include <AMReX_Vector.H>
12 #include <AMReX_Utility.H>
13 #include <AMReX_Geometry.H>
14 #include <AMReX_VisMF.H>
15 #include <AMReX_RealBox.H>
16 #include <AMReX_Print.H>
17 #include <AMReX_MultiFabUtil.H>
18 #include <AMReX_NFiles.H>
19 #include <AMReX_VectorIO.H>
20 #include <AMReX_Particle_mod_K.H>
21 #include <AMReX_ParticleMPIUtil.H>
22 #include <AMReX_StructOfArrays.H>
23 #include <AMReX_ArrayOfStructs.H>
24 #include <AMReX_Particle.H>
25 #include <AMReX_ParticleTile.H>
26 #include <AMReX_TypeTraits.H>
27 #include <AMReX_GpuContainers.H>
28 #include <AMReX_ParticleUtil.H>
29 #include <AMReX_ParticleReduce.H>
32 #include <AMReX_ParticleLocator.H>
33 #include <AMReX_Scan.H>
34 #include <AMReX_DenseBins.H>
35 #include <AMReX_SparseBins.H>
37 #include <AMReX_ParticleMesh.H>
38 #include <AMReX_OpenMP.H>
39 #include <AMReX_ParIter.H>
40 
41 #ifdef AMREX_LAZY
42 #include <AMReX_Lazy.H>
43 #endif
44 
45 #ifdef AMREX_USE_OMP
46 #include <omp.h>
47 #endif
48 
49 #ifdef AMREX_USE_HDF5
50 #include <hdf5.h>
51 #endif
52 
53 #include <algorithm>
54 #include <array>
55 #include <cstring>
56 #include <fstream>
57 #include <iostream>
58 #include <limits>
59 #include <map>
60 #include <memory>
61 #include <numeric>
62 #include <random>
63 #include <string>
64 #include <tuple>
65 #include <type_traits>
66 #include <utility>
67 #include <vector>
68 
69 namespace amrex {
70 
71 #ifdef AMREX_USE_HDF5_ASYNC
72  extern hid_t es_par_g;
73 #endif
79 {
80  using RealType = ParticleReal;
81  int m_lev;
82  int m_grid;
84  RealType m_data[1 + AMREX_SPACEDIM];
85 };
86 
91 {
92  int m_lev = -1;
93  int m_grid = -1;
94  int m_tile = -1;
99 };
100 
114 template<int NStructReal, int NStructInt, int NArrayReal, int NArrayInt>
116 {
117  std::array<double, NStructReal> real_struct_data;
118  std::array<int, NStructInt > int_struct_data;
119  std::array<double, NArrayReal > real_array_data;
120  std::array<int, NArrayInt > int_array_data;
121 };
122 
123 template <bool is_const, typename T_ParticleType, int NArrayReal, int NArrayInt,
124  template<class> class Allocator, class CellAssignor>
125 class ParIterBase_impl;
126 
140 template <typename T_ParticleType, int T_NArrayReal, int T_NArrayInt,
141  template<class> class Allocator=DefaultAllocator,
142  class T_CellAssignor=DefaultAssignor>
143 
145 {
146 public:
147  using ParticleType = T_ParticleType;
148  using ConstParticleType = typename ParticleType::ConstType;
149  using CellAssignor = T_CellAssignor;
150 
152  static constexpr int NStructReal = ParticleType::NReal;
154  static constexpr int NStructInt = ParticleType::NInt;
156  static constexpr int NArrayReal = T_NArrayReal;
158  static constexpr int NArrayInt = T_NArrayInt;
160 
161 private:
162  friend class ParIterBase_impl<true,ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>;
163  friend class ParIterBase_impl<false,ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>;
164 
165 public:
167  template <typename T>
168  using AllocatorType = Allocator<T>;
173 
174 #ifdef AMREX_SINGLE_PRECISION_PARTICLES
176 #else
178 #endif
179 
183 
186  using ParticleLevel = std::map<std::pair<int, int>, ParticleTileType>;
187  using AoS = typename ParticleTileType::AoS;
188  using SoA = typename ParticleTileType::SoA;
189 
190  using RealVector = typename SoA::RealVector;
191  using IntVector = typename SoA::IntVector;
192  using ParticleVector = typename AoS::ParticleVector;
196 
200  :
202  h_redistribute_real_comp(AMREX_SPACEDIM + NStructReal + NArrayReal, true),
204  {
205  Initialize ();
206  }
207 
216  :
218  h_redistribute_real_comp(AMREX_SPACEDIM + NStructReal + NArrayReal, true),
220  {
221  Initialize ();
224  }
225 
234  const DistributionMapping & dmap,
235  const BoxArray & ba)
236  :
237  ParticleContainerBase(geom, dmap, ba),
238  h_redistribute_real_comp(AMREX_SPACEDIM + NStructReal + NArrayReal, true),
240  {
241  Initialize ();
244  }
245 
256  const Vector<DistributionMapping> & dmap,
257  const Vector<BoxArray> & ba,
258  const Vector<int> & rr)
259  :
260  ParticleContainerBase(geom, dmap, ba, rr),
261  h_redistribute_real_comp(AMREX_SPACEDIM + NStructReal + NArrayReal, true),
263  {
264  Initialize ();
267  }
268 
278  const Vector<DistributionMapping> & dmap,
279  const Vector<BoxArray> & ba,
280  const Vector<IntVect> & rr)
281  :
282  ParticleContainerBase(geom, dmap, ba, rr),
283  h_redistribute_real_comp(AMREX_SPACEDIM + NStructReal + NArrayReal, true),
285  {
286  Initialize ();
289  }
290 
291  ~ParticleContainer_impl () override = default;
292 
295 
297  ParticleContainer_impl& operator= ( ParticleContainer_impl && ) noexcept = default;
298 
299 
307  void Define (ParGDBBase* gdb)
308  {
310  reserveData();
311  resizeData();
312  }
313 
321  void Define (const Geometry & geom,
322  const DistributionMapping & dmap,
323  const BoxArray & ba)
324  {
325  this->ParticleContainerBase::Define(geom, dmap, ba);
326  reserveData();
327  resizeData();
328  }
329 
339  void Define (const Vector<Geometry> & geom,
340  const Vector<DistributionMapping> & dmap,
341  const Vector<BoxArray> & ba,
342  const Vector<int> & rr)
343  {
344  this->ParticleContainerBase::Define(geom, dmap, ba, rr);
345  reserveData();
346  resizeData();
347  }
348 
358  void Define (const Vector<Geometry> & geom,
359  const Vector<DistributionMapping> & dmap,
360  const Vector<BoxArray> & ba,
361  const Vector<IntVect> & rr)
362  {
363  this->ParticleContainerBase::Define(geom, dmap, ba, rr);
364  reserveData();
365  resizeData();
366  }
367 
369  int numLocalTilesAtLevel (int lev) const {
370  return (lev < m_particles.size()) ? m_particles[lev].size() : 0;
371  }
372 
377  void reserveData () override;
378 
385  void resizeData () override;
386 
387  void InitFromAsciiFile (const std::string& file, int extradata,
388  const IntVect* Nrep = nullptr);
389 
390  void InitFromBinaryFile (const std::string& file, int extradata);
391 
392  void InitFromBinaryMetaFile (const std::string& file, int extradata);
393 
410  void InitRandom (Long icount, ULong iseed,
411  const ParticleInitData& pdata,
412  bool serialize = false, RealBox bx = RealBox());
413 
414 
430  void InitRandomPerBox (Long icount, ULong iseed, const ParticleInitData& pdata);
431 
432 
447  void InitOnePerCell (Real x_off, Real y_off, Real z_off,
448  const ParticleInitData& pdata);
449 
450 
462  void InitNRandomPerCell (int n_per_cell, const ParticleInitData& pdata);
463 
464  void Increment (MultiFab& mf, int level);
465 
466  Long IncrementWithTotal (MultiFab& mf, int level, bool local = false);
467 
504  void Redistribute (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
505  bool remove_negative=true);
506 
507 
519  template <class index_type>
520  void ReorderParticles (int lev, const MFIter& mfi, const index_type* permutations);
521 
536  void SortParticlesForDeposition (IntVect idx_type);
537 
541  void SortParticlesByCell ();
542 
549  void SortParticlesByBin (IntVect bin_size);
550 
561  bool OK (int lev_min = 0, int lev_max = -1, int nGrow = 0) const;
562 
563  std::array<Long, 3> ByteSpread () const;
564 
565  std::array<Long, 3> PrintCapacity () const;
566 
567  void ShrinkToFit ();
568 
579  Long NumberOfParticlesAtLevel (int level, bool only_valid = true, bool only_local = false) const;
580 
581  Vector<Long> NumberOfParticlesInGrid (int level, bool only_valid = true, bool only_local = false) const;
582 
591  Long TotalNumberOfParticles (bool only_valid=true, bool only_local=false) const;
592 
593 
601  void RemoveParticlesAtLevel (int level);
602 
604 
612  void CreateVirtualParticles (int level, AoS& virts) const;
613 
622  void CreateGhostParticles (int level, int ngrow, AoS& ghosts) const;
623 
631  void AddParticlesAtLevel (AoS& particles, int level, int nGrow=0);
632 
640  void CreateVirtualParticles (int level, ParticleTileType& virts) const;
641 
650  void CreateGhostParticles (int level, int ngrow, ParticleTileType& ghosts) const;
651 
659  void AddParticlesAtLevel (ParticleTileType& particles, int level, int nGrow=0);
660 
661 
665  void clearParticles ();
666 
667 
676  template <class PCType,
677  std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0>
678  void copyParticles (const PCType& other, bool local=false);
679 
687  template <class PCType,
688  std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0>
689  void addParticles (const PCType& other, bool local=false);
690 
705  template <class F, class PCType,
706  std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0,
707  std::enable_if_t<! std::is_integral_v<F>, int> bar = 0>
708  void copyParticles (const PCType& other, F&&f, bool local=false);
709 
723  template <class F, class PCType,
724  std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0,
725  std::enable_if_t<! std::is_integral_v<F>, int> bar = 0>
726  void addParticles (const PCType& other, F const& f, bool local=false);
727 
735  void WriteParticleRealData (void* data, size_t size, std::ostream& os) const;
736 
745  void ReadParticleRealData (void* data, size_t size, std::istream& is);
746 
755  void Checkpoint (const std::string& dir, const std::string& name,
756  const Vector<std::string>& real_comp_names = Vector<std::string>(),
757  const Vector<std::string>& int_comp_names = Vector<std::string>()) const
758  {
759  Checkpoint(dir, name, true, real_comp_names, int_comp_names);
760  }
761 
767  void Checkpoint (const std::string& dir, const std::string& name, bool is_checkpoint,
768  const Vector<std::string>& real_comp_names = Vector<std::string>(),
769  const Vector<std::string>& int_comp_names = Vector<std::string>()) const;
770 
785  template <class F>
786  void WriteBinaryParticleData (const std::string& dir,
787  const std::string& name,
788  const Vector<int>& write_real_comp,
789  const Vector<int>& write_int_comp,
790  const Vector<std::string>& real_comp_names,
791  const Vector<std::string>& int_comp_names,
792  F&& f, bool is_checkpoint=false) const;
793 
794  void CheckpointPre ();
795 
796  void CheckpointPost ();
797 
804  void Restart (const std::string& dir, const std::string& file);
805 
813  void Restart (const std::string& dir, const std::string& file, bool is_checkpoint);
814 
821  void WritePlotFile (const std::string& dir, const std::string& name) const;
822 
834  template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>* = nullptr>
835  void WritePlotFile (const std::string& dir, const std::string& name, F&& f) const;
836 
845  void WritePlotFile (const std::string& dir, const std::string& name,
846  const Vector<std::string>& real_comp_names,
847  const Vector<std::string>& int_comp_names) const;
848 
862  template <class F>
863  void WritePlotFile (const std::string& dir, const std::string& name,
864  const Vector<std::string>& real_comp_names,
865  const Vector<std::string>& int_comp_names, F&& f) const;
866 
875  void WritePlotFile (const std::string& dir, const std::string& name,
876  const Vector<std::string>& real_comp_names) const;
877 
891  template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>* = nullptr>
892  void WritePlotFile (const std::string& dir, const std::string& name,
893  const Vector<std::string>& real_comp_names, F&& f) const;
894 
904  void WritePlotFile (const std::string& dir,
905  const std::string& name,
906  const Vector<int>& write_real_comp,
907  const Vector<int>& write_int_comp) const;
908 
923  template <class F>
924  void WritePlotFile (const std::string& dir,
925  const std::string& name,
926  const Vector<int>& write_real_comp,
927  const Vector<int>& write_int_comp, F&& f) const;
928 
942  void WritePlotFile (const std::string& dir,
943  const std::string& name,
944  const Vector<int>& write_real_comp,
945  const Vector<int>& write_int_comp,
946  const Vector<std::string>& real_comp_names,
947  const Vector<std::string>& int_comp_names) const;
948 
967  template <class F>
968  void WritePlotFile (const std::string& dir,
969  const std::string& name,
970  const Vector<int>& write_real_comp,
971  const Vector<int>& write_int_comp,
972  const Vector<std::string>& real_comp_names,
973  const Vector<std::string>& int_comp_names,
974  F&& f) const;
975 
976  void WritePlotFilePre ();
977 
978  void WritePlotFilePost ();
979 
980  void WriteAsciiFile (const std::string& file);
981 
986  const Vector<ParticleLevel>& GetParticles () const { return m_particles; }
987 
993 
1007  const ParticleLevel& GetParticles (int lev) const { return m_particles[lev]; }
1008 
1022  ParticleLevel & GetParticles (int lev) { return m_particles[lev]; }
1023 
1046  const ParticleTileType& ParticlesAt (int lev, int grid, int tile) const
1047  { return m_particles[lev].at(std::make_pair(grid, tile)); }
1048 
1071  ParticleTileType& ParticlesAt (int lev, int grid, int tile)
1072  { return m_particles[lev].at(std::make_pair(grid, tile)); }
1073 
1095  template <class Iterator>
1096  const ParticleTileType& ParticlesAt (int lev, const Iterator& iter) const
1097  { return ParticlesAt(lev, iter.index(), iter.LocalTileIndex()); }
1098 
1120  template <class Iterator>
1121  ParticleTileType& ParticlesAt (int lev, const Iterator& iter)
1122  { return ParticlesAt(lev, iter.index(), iter.LocalTileIndex()); }
1123 
1146  ParticleTileType& DefineAndReturnParticleTile (int lev, int grid, int tile)
1147  {
1148  m_particles[lev][std::make_pair(grid, tile)].define(NumRuntimeRealComps(), NumRuntimeIntComps(), &m_soa_rdata_names, &m_soa_idata_names);
1149 
1150  return ParticlesAt(lev, grid, tile);
1151  }
1152 
1174  template <class Iterator>
1175  ParticleTileType& DefineAndReturnParticleTile (int lev, const Iterator& iter)
1176  {
1177  auto index = std::make_pair(iter.index(), iter.LocalTileIndex());
1178  m_particles[lev][index].define(NumRuntimeRealComps(), NumRuntimeIntComps());
1179  return ParticlesAt(lev, iter);
1180  }
1181 
1192  void AssignDensity (int rho_index,
1193  Vector<std::unique_ptr<MultiFab> >& mf_to_be_filled,
1194  int lev_min, int ncomp, int finest_level, int ngrow=2) const;
1195 
1196  void AssignCellDensitySingleLevel (int rho_index, MultiFab& mf, int level,
1197  int ncomp=1, int particle_lvl_offset = 0) const;
1198 
1199  template <typename P, typename Assignor=CellAssignor>
1200  IntVect Index (const P& p, int lev) const;
1201 
1211  ParticleLocData Reset (ParticleType& prt, bool update, bool verbose=true,
1212  ParticleLocData pld = ParticleLocData()) const;
1213 
1219  template <typename P>
1220  bool PeriodicShift (P& p) const;
1221 
1223 
1225 
1226  void SetUsePrePost (bool tf) const {
1227  usePrePost = tf;
1228  }
1229  bool GetUsePrePost () const {
1230  return usePrePost;
1231  }
1232 
1233  int GetMaxNextIDPrePost () const { return maxnextidPrePost; }
1234  Long GetNParticlesPrePost () const { return nparticlesPrePost; }
1235 
1236  void SetUseUnlink (bool tf) const {
1237  doUnlink = tf;
1238  }
1239 
1240  bool GetUseUnlink () const {
1241  return doUnlink;
1242  }
1243 
1244  void RedistributeCPU (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
1245  bool remove_negative=true);
1246 
1247  void RedistributeGPU (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
1248  bool remove_negative=true);
1249 
1250  Long superParticleSize() const { return superparticle_size; }
1251 
1252  void AddRealComp (std::string const & name, int communicate=1)
1253  {
1254  m_soa_rdata_names.push_back(name);
1255 
1256  m_runtime_comps_defined = true;
1258  h_redistribute_real_comp.push_back(communicate);
1259  SetParticleSize();
1260  this->resizeData();
1261 
1262  // resize runtime SoA
1263  for (int lev = 0; lev < numLevels(); ++lev) {
1264  for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1265  auto& tile = DefineAndReturnParticleTile(lev, pti);
1266  auto np = tile.numParticles();
1267  if (np > 0) {
1268  auto& soa = tile.GetStructOfArrays();
1269  soa.resize(np);
1270  }
1271  }
1272  }
1273  }
1274 
1275  void AddRealComp (int communicate=1)
1276  {
1277  AddRealComp(getDefaultCompNameReal<ParticleType>(NArrayReal+m_num_runtime_real), communicate);
1278  }
1279 
1280  void AddIntComp (std::string const & name, int communicate=1)
1281  {
1282  m_soa_idata_names.push_back(name);
1283 
1284  m_runtime_comps_defined = true;
1286  h_redistribute_int_comp.push_back(communicate);
1287  SetParticleSize();
1288  this->resizeData();
1289 
1290  // resize runtime SoA
1291  for (int lev = 0; lev < numLevels(); ++lev) {
1292  for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1293  auto& tile = DefineAndReturnParticleTile(lev, pti);
1294  auto np = tile.numParticles();
1295  if (np > 0) {
1296  auto& soa = tile.GetStructOfArrays();
1297  soa.resize(np);
1298  }
1299  }
1300  }
1301  }
1302 
1303  void AddIntComp (int communicate=1)
1304  {
1305  AddIntComp(getDefaultCompNameInt<ParticleType>(NArrayInt+m_num_runtime_int), communicate);
1306  }
1307 
1308  int NumRuntimeRealComps () const { return m_num_runtime_real; }
1309  int NumRuntimeIntComps () const { return m_num_runtime_int; }
1310 
1311  int NumRealComps () const { return NArrayReal + NumRuntimeRealComps(); }
1312  int NumIntComps () const { return NArrayInt + NumRuntimeIntComps() ; }
1313 
1319  void ResizeRuntimeRealComp (int new_size, bool communicate);
1320 
1326  void ResizeRuntimeIntComp (int new_size, bool communicate);
1327 
1329  template <template<class> class NewAllocator=amrex::DefaultAllocator>
1331 
1342  template <template<class> class NewAllocator=amrex::DefaultAllocator>
1344  make_alike () const
1345  {
1347 
1348  // add runtime real comps to tmp
1349  for (int ic = 0; ic < this->NumRuntimeRealComps(); ++ic) { tmp.AddRealComp(false); }
1350 
1351  // add runtime int comps to tmp
1352  for (int ic = 0; ic < this->NumRuntimeIntComps(); ++ic) { tmp.AddIntComp(false); }
1353 
1354  return tmp;
1355  }
1356 
1359 
1362  mutable bool usePrePost;
1363  mutable bool doUnlink;
1365  mutable int nOutFilesPrePost;
1371  mutable std::string HdrFileNamePrePost;
1373 
1374 protected:
1375 
1389  template <typename P>
1390  bool Where (const P& prt, ParticleLocData& pld,
1391  int lev_min = 0, int lev_max = -1, int nGrow=0, int local_grid=-1) const;
1392 
1393 
1404  template <typename P>
1406  int lev_min = 0, int lev_max = -1, int local_grid=-1) const;
1407 
1408 public:
1409  void
1410  WriteParticles (int level, std::ofstream& ofs, int fnum,
1411  Vector<int>& which, Vector<int>& count, Vector<Long>& where,
1412  const Vector<int>& write_real_comp, const Vector<int>& write_int_comp,
1413  const Vector<std::map<std::pair<int, int>,IntVector>>& particle_io_flags, bool is_checkpoint) const;
1414 #ifdef AMREX_USE_HDF5
1415 #include "AMReX_ParticlesHDF5.H"
1416 #endif
1417 
1419  void SetSoACompileTimeNames (std::vector<std::string> const & rdata_name, std::vector<std::string> const & idata_name);
1420 
1422  std::vector<std::string> GetRealSoANames () const {return m_soa_rdata_names;}
1423 
1425  std::vector<std::string> GetIntSoANames () const {return m_soa_idata_names;}
1426 
1427 protected:
1428 
1429  template <class RTYPE>
1430  void ReadParticles (int cnt, int grd, int lev, std::ifstream& ifs, int finest_level_in_file, bool convert_ids);
1431 
1432  void SetParticleSize ();
1433 
1435 
1436 private:
1437  virtual void particlePostLocate (ParticleType& /*p*/, const ParticleLocData& /*pld*/,
1438  const int /*lev*/) {}
1439 
1440  virtual void correctCellVectors (int /*old_index*/, int /*new_index*/,
1441  int /*grid*/, const ParticleType& /*p*/) {}
1442 
1443  void RedistributeMPI (std::map<int, Vector<char> >& not_ours,
1444  int lev_min = 0, int lev_max = 0, int nGrow = 0, int local=0);
1445 
1446  template <typename P>
1448  int lev_min, int lev_max, int nGrow, int local_grid=-1) const;
1449 
1450  void Initialize ();
1451 
1455 
1459 
1460  // names of both compile-time and runtime Real and Int SoA data
1461  std::vector<std::string> m_soa_rdata_names;
1462  std::vector<std::string> m_soa_idata_names;
1463 };
1464 
1465 template <int T_NStructReal, int T_NStructInt, int T_NArrayReal, int T_NArrayInt, template<class> class Allocator, class CellAssignor>
1466 using ParticleContainer = ParticleContainer_impl<Particle<T_NStructReal, T_NStructInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
1467 
1468 template <int T_NArrayReal, int T_NArrayInt, template<class> class Allocator=DefaultAllocator, class CellAssignor=DefaultAssignor>
1469 using ParticleContainerPureSoA = ParticleContainer_impl<SoAParticle<T_NArrayReal, T_NArrayInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
1470 
1471 
1472 #include "AMReX_ParticleInit.H"
1473 #include "AMReX_ParticleContainerI.H"
1474 #include "AMReX_ParticleIO.H"
1475 
1476 #ifdef AMREX_USE_HDF5
1477 #include "AMReX_ParticleHDF5.H"
1478 #endif
1479 
1480 }
1481 
1482 #endif /*_PARTICLES_H_*/
#define AMREX_D_DECL(a, b, c)
Definition: AMReX_SPACE.H:104
Definition: AMReX_GpuAllocators.H:114
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:549
A container for storing items in a set of bins.
Definition: AMReX_DenseBins.H:77
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
static const RealDescriptor & Native64RealDescriptor()
Definition: AMReX_FPC.cpp:137
static const RealDescriptor & Native32RealDescriptor()
NativeRealDescriptor is equivalent to Native32RealDescriptor if BL_FLOAT is used. Otherwise,...
Definition: AMReX_FPC.cpp:124
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
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
Definition: AMReX_ParIter.H:142
Definition: AMReX_ParGDB.H:13
Definition: AMReX_ParIter.H:32
Definition: AMReX_ParIter.H:113
Definition: AMReX_ParticleContainerBase.H:23
ParGDBBase * m_gdb
Definition: AMReX_ParticleContainerBase.H:269
void Define(ParGDBBase *gdb)
Definition: AMReX_ParticleContainerBase.H:83
int numLevels() const
the number of defined levels in the ParticleContainer
Definition: AMReX_ParticleContainerBase.H:222
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition: AMReX_ParticleContainer.H:145
void AssignDensity(int rho_index, Vector< std::unique_ptr< MultiFab > > &mf_to_be_filled, int lev_min, int ncomp, int finest_level, int ngrow=2) const
Functions depending the layout of the data. Use with caution.
Definition: AMReX_AmrParticles.H:17
void SetUseUnlink(bool tf) const
Definition: AMReX_ParticleContainer.H:1236
Vector< ParticleLevel > m_particles
Definition: AMReX_ParticleContainer.H:1458
Vector< Long > nParticlesAtLevelPrePost
Definition: AMReX_ParticleContainer.H:1367
void WriteBinaryParticleData(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 &&f, bool is_checkpoint=false) const
Writes particle data to disk in the AMReX native format.
void InitFromBinaryFile(const std::string &file, int extradata)
Definition: AMReX_ParticleInit.H:483
bool GetLevelDirectoriesCreated() const
Definition: AMReX_ParticleContainer.H:1224
size_t superparticle_size
Definition: AMReX_ParticleContainer.H:1456
void AddIntComp(std::string const &name, int communicate=1)
Definition: AMReX_ParticleContainer.H:1280
std::array< Long, 3 > ByteSpread() const
Definition: AMReX_ParticleContainerI.H:501
void clearParticles()
Clear all the particles in this container. This does not free memory.
Definition: AMReX_ParticleContainerI.H:1010
void copyParticles(const PCType &other, bool local=false)
Copy particles from other to this ParticleContainer. Will clear all the particles from this container...
bool doUnlink
Definition: AMReX_ParticleContainer.H:1363
Vector< Vector< int > > whichPrePost
Definition: AMReX_ParticleContainer.H:1368
void WriteParticleRealData(void *data, size_t size, std::ostream &os) const
Write a contiguous chunk of real particle data to an ostream.
Definition: AMReX_ParticleIO.H:11
const ParticleTileType & ParticlesAt(int lev, const Iterator &iter) const
Return the ParticleTile for level "lev" and Iterator "iter". Const version.
Definition: AMReX_ParticleContainer.H:1096
ParticleContainer_impl(const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< IntVect > &rr)
Same as the above, but accepts different refinement ratios in each direction.
Definition: AMReX_ParticleContainer.H:277
std::string HdrFileNamePrePost
Definition: AMReX_ParticleContainer.H:1371
Vector< Vector< Long > > wherePrePost
Definition: AMReX_ParticleContainer.H:1370
void WritePlotFile(const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names, F &&f) const
This version of WritePlotFile writes all components and allows the user to specify the names of the c...
RealDescriptor ParticleRealDescriptor
Definition: AMReX_ParticleContainer.H:177
ParticleTileType & DefineAndReturnParticleTile(int lev, int grid, int tile)
Define and return the ParticleTile for level "lev", grid "grid" and tile "tile.".
Definition: AMReX_ParticleContainer.H:1146
int nOutFilesPrePost
Definition: AMReX_ParticleContainer.H:1365
int maxnextidPrePost
Definition: AMReX_ParticleContainer.H:1364
void WritePlotFile(const std::string &dir, const std::string &name, F &&f) const
This version of WritePlotFile writes all components and assigns component names.
ParticleContainer_impl(ParticleContainer_impl &&) noexcept=default
void Checkpoint(const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names=Vector< std::string >(), const Vector< std::string > &int_comp_names=Vector< std::string >()) const
Writes a particle checkpoint to file, suitable for restarting.
Definition: AMReX_ParticleContainer.H:755
void addParticles(const PCType &other, bool local=false)
Add particles from other to this ParticleContainer. local controls whether or not to call Redistribut...
std::map< std::pair< int, int >, ParticleTileType > ParticleLevel
Definition: AMReX_ParticleContainer.H:186
bool levelDirectoriesCreated
Variables for i/o optimization saved for pre and post checkpoint.
Definition: AMReX_ParticleContainer.H:1361
ContainerLike< NewAllocator > make_alike() const
Definition: AMReX_ParticleContainer.H:1344
void WriteAsciiFile(const std::string &file)
Definition: AMReX_ParticleIO.H:1104
std::vector< std::string > m_soa_idata_names
Definition: AMReX_ParticleContainer.H:1462
static constexpr int NArrayInt
Number of extra integer components stored in struct-of-array form.
Definition: AMReX_ParticleContainer.H:158
void Initialize()
Definition: AMReX_ParticleContainerI.H:37
static constexpr int NArrayReal
Number of extra Real components stored in struct-of-array form.
Definition: AMReX_ParticleContainer.H:156
int GetMaxNextIDPrePost() const
Definition: AMReX_ParticleContainer.H:1233
std::array< Long, 3 > PrintCapacity() const
Definition: AMReX_ParticleContainerI.H:545
bool GetUsePrePost() const
Definition: AMReX_ParticleContainer.H:1229
ParticleContainer_impl & operator=(const ParticleContainer_impl &)=delete
bool PeriodicShift(P &p) const
Returns true if the particle was shifted.
typename ParticleTileType::AoS AoS
Definition: AMReX_ParticleContainer.H:187
ParticleInitType< NStructReal, NStructInt, NArrayReal, NArrayInt > ParticleInitData
Definition: AMReX_ParticleContainer.H:182
typename SoA::RealVector RealVector
Definition: AMReX_ParticleContainer.H:190
typename ParticleType::ConstType ConstParticleType
Definition: AMReX_ParticleContainer.H:148
ParticleContainer_impl(const ParticleContainer_impl &)=delete
void WriteParticles(int level, std::ofstream &ofs, int fnum, Vector< int > &which, Vector< int > &count, Vector< Long > &where, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::map< std::pair< int, int >, IntVector >> &particle_io_flags, bool is_checkpoint) const
Definition: AMReX_ParticleIO.H:556
ParticleContainer_impl(const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
Construct a particle container using a given Geometry, DistributionMapping, and BoxArray....
Definition: AMReX_ParticleContainer.H:233
Long NumberOfParticlesAtLevel(int level, bool only_valid=true, bool only_local=false) const
Returns # of particles at specified the level.
Definition: AMReX_ParticleContainerI.H:455
void AssignCellDensitySingleLevel(int rho_index, MultiFab &mf, int level, int ncomp=1, int particle_lvl_offset=0) const
Definition: AMReX_ParticleContainerI.H:2401
void Increment(MultiFab &mf, int level)
Definition: AMReX_ParticleContainerI.H:605
void SetLevelDirectoriesCreated(bool tf)
Definition: AMReX_ParticleContainer.H:1222
void ReorderParticles(int lev, const MFIter &mfi, const index_type *permutations)
Reorder particles on the tile given by lev and mfi using a the permutations array.
typename AoS::ParticleVector ParticleVector
Definition: AMReX_ParticleContainer.H:192
const ParticleLevel & GetParticles(int lev) const
Return the ParticleLevel for level "lev". Const version.
Definition: AMReX_ParticleContainer.H:1007
void Define(const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
Define a default-constructed ParticleContainer using a ParGDB object. Single-level version.
Definition: AMReX_ParticleContainer.H:321
static constexpr int NStructInt
Number of extra integer components in the particle struct.
Definition: AMReX_ParticleContainer.H:154
ParticleContainer_impl()
Default constructor - construct an empty particle container that has no concept of a level hierarchy....
Definition: AMReX_ParticleContainer.H:199
void SetUsePrePost(bool tf) const
Definition: AMReX_ParticleContainer.H:1226
int NumRealComps() const
Definition: AMReX_ParticleContainer.H:1311
void WritePlotFilePre()
Definition: AMReX_ParticleIO.H:536
bool EnforcePeriodicWhere(P &prt, ParticleLocData &pld, int lev_min=0, int lev_max=-1, int local_grid=-1) const
Checks whether the particle has crossed a periodic boundary in such a way that it is on levels lev_mi...
void RemoveParticlesAtLevel(int level)
The Following methods are for managing Virtual and Ghost Particles.
Definition: AMReX_ParticleContainerI.H:642
void InitFromAsciiFile(const std::string &file, int extradata, const IntVect *Nrep=nullptr)
Definition: AMReX_ParticleInit.H:37
ParticleTileType & DefineAndReturnParticleTile(int lev, const Iterator &iter)
Define and return the ParticleTile for level "lev", and Iterator "iter".
Definition: AMReX_ParticleContainer.H:1175
Long TotalNumberOfParticles(bool only_valid=true, bool only_local=false) const
Returns # of particles at all levels.
Definition: AMReX_ParticleContainerI.H:385
Vector< int > h_redistribute_int_comp
Definition: AMReX_ParticleContainer.H:1358
std::vector< std::string > m_soa_rdata_names
Definition: AMReX_ParticleContainer.H:1461
typename ParticleTileType::SoA SoA
Definition: AMReX_ParticleContainer.H:188
virtual void correctCellVectors(int, int, int, const ParticleType &)
Definition: AMReX_ParticleContainer.H:1440
void WritePlotFile(const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, F &&f) const
This version of WritePlotFile writes all components and allows the user to specify the names of the c...
bool OK(int lev_min=0, int lev_max=-1, int nGrow=0) const
OK checks that all particles are in the right places (for some value of right)
Definition: AMReX_ParticleContainerI.H:2340
void locateParticle(P &p, ParticleLocData &pld, int lev_min, int lev_max, int nGrow, int local_grid=-1) const
int NumRuntimeIntComps() const
Definition: AMReX_ParticleContainer.H:1309
Long GetNParticlesPrePost() const
Definition: AMReX_ParticleContainer.H:1234
void InitRandom(Long icount, ULong iseed, const ParticleInitData &pdata, bool serialize=false, RealBox bx=RealBox())
This initializes the particle container with icount randomly distributed particles....
Definition: AMReX_ParticleInit.H:968
void ResizeRuntimeRealComp(int new_size, bool communicate)
Definition: AMReX_ParticleContainerI.H:2541
T_CellAssignor CellAssignor
Definition: AMReX_ParticleContainer.H:149
void RemoveParticlesNotAtFinestLevel()
Definition: AMReX_ParticleContainerI.H:656
DenseBins< typename ParticleTileType::ParticleTileDataType > m_bins
Definition: AMReX_ParticleContainer.H:1434
void ReadParticleRealData(void *data, size_t size, std::istream &is)
Read a contiguous chunk of real particle data from an istream.
Definition: AMReX_ParticleIO.H:25
void RedistributeCPU(int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
Definition: AMReX_ParticleContainerI.H:1479
int m_num_runtime_int
Definition: AMReX_ParticleContainer.H:1454
ParticleContainer_impl(ParGDBBase *gdb)
Construct a particle container using a ParGDB object. The container will track changes in the grid st...
Definition: AMReX_ParticleContainer.H:215
void Restart(const std::string &dir, const std::string &file)
Restart from checkpoint.
Definition: AMReX_ParticleIO.H:620
void WritePlotFilePost()
Definition: AMReX_ParticleIO.H:546
void Define(const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< IntVect > &rr)
Define a default-constructed ParticleContainer using a ParGDB object. Multi-level version.
Definition: AMReX_ParticleContainer.H:358
void ReadParticles(int cnt, int grd, int lev, std::ifstream &ifs, int finest_level_in_file, bool convert_ids)
void CreateVirtualParticles(int level, AoS &virts) const
Creates virtual particles for a given level that represent in some capacity all particles at finer le...
Definition: AMReX_ParticleContainerI.H:723
int NumRuntimeRealComps() const
Definition: AMReX_ParticleContainer.H:1308
void Checkpoint(const std::string &dir, const std::string &name, bool is_checkpoint, const Vector< std::string > &real_comp_names=Vector< std::string >(), const Vector< std::string > &int_comp_names=Vector< std::string >()) const
Writes a particle checkpoint to file, suitable for restarting. This version allows the particle compo...
void AddIntComp(int communicate=1)
Definition: AMReX_ParticleContainer.H:1303
Long nparticlesPrePost
Definition: AMReX_ParticleContainer.H:1366
void InitNRandomPerCell(int n_per_cell, const ParticleInitData &pdata)
This initializes the particle container with n_per_cell randomly distributed particles per cell,...
Definition: AMReX_ParticleInit.H:1547
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:93
void WritePlotFile(const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, F &&f) const
This version of WritePlotFile assigns component names, but allows the user to pass in a vector of int...
int m_num_runtime_real
Definition: AMReX_ParticleContainer.H:1453
void copyParticles(const PCType &other, F &&f, bool local=false)
Copy particles from other to this ParticleContainer. Will clear all the particles from this container...
void CheckpointPre()
Definition: AMReX_ParticleIO.H:423
ParticleLocData Reset(ParticleType &prt, bool update, bool verbose=true, ParticleLocData pld=ParticleLocData()) const
Updates a particle's location (Where), tries to periodic shift any particles that have left the domai...
Definition: AMReX_ParticleContainerI.H:296
Vector< std::string > filePrefixPrePost
Definition: AMReX_ParticleContainer.H:1372
void RedistributeGPU(int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
Definition: AMReX_ParticleContainerI.H:1283
int numLocalTilesAtLevel(int lev) const
The total number of tiles on this rank on this level.
Definition: AMReX_ParticleContainer.H:369
void reserveData() override
This reserves data in the vector of dummy MultiFabs used by the ParticleContainer for the maximum num...
Definition: AMReX_ParticleContainerI.H:330
void WritePlotFile(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 &&f) const
This is the most general version of WritePlotFile, which takes component names and flags for whether ...
std::vector< std::string > GetRealSoANames() const
Definition: AMReX_ParticleContainer.H:1422
static constexpr int NStructReal
Number of extra Real components in the particle struct.
Definition: AMReX_ParticleContainer.H:152
Vector< int > h_redistribute_real_comp
Definition: AMReX_ParticleContainer.H:1357
std::vector< std::string > GetIntSoANames() const
Definition: AMReX_ParticleContainer.H:1425
typename Particle< NStructReal, NStructInt >::RealType RealType
The type of the Real data.
Definition: AMReX_ParticleContainer.H:172
bool GetUseUnlink() const
Definition: AMReX_ParticleContainer.H:1240
Vector< Long > NumberOfParticlesInGrid(int level, bool only_valid=true, bool only_local=false) const
Definition: AMReX_ParticleContainerI.H:400
int NumIntComps() const
Definition: AMReX_ParticleContainer.H:1312
Allocator< T > AllocatorType
The memory allocator in use.
Definition: AMReX_ParticleContainer.H:168
void RedistributeMPI(std::map< int, Vector< char > > &not_ours, int lev_min=0, int lev_max=0, int nGrow=0, int local=0)
Definition: AMReX_ParticleContainerI.H:1993
ParticleLevel & GetParticles(int lev)
Return the ParticleLevel for level "lev". Non-const version.
Definition: AMReX_ParticleContainer.H:1022
void CreateGhostParticles(int level, int ngrow, ParticleTileType &ghosts) const
Create ghost particles for a given level that are copies of particles near coarse->fine boundaries in...
void AddParticlesAtLevel(ParticleTileType &particles, int level, int nGrow=0)
Add particles from a pbox to the grid at this level.
void ShrinkToFit()
Definition: AMReX_ParticleContainerI.H:585
void addParticles(const PCType &other, F const &f, bool local=false)
Add particles from other to this ParticleContainer. local controls whether or not to call Redistribut...
Vector< ParticleLevel > & GetParticles()
Return the underlying Vector (over AMR levels) of ParticleLevels. Non-const version.
Definition: AMReX_ParticleContainer.H:992
bool Where(const P &prt, ParticleLocData &pld, int lev_min=0, int lev_max=-1, int nGrow=0, int local_grid=-1) const
Checks a particle's location on levels lev_min and higher. Returns false if the particle does not exi...
void Restart(const std::string &dir, const std::string &file, bool is_checkpoint)
Older version, for backwards compatibility.
void Define(ParGDBBase *gdb)
Define a default-constructed ParticleContainer using a ParGDB object. The container will track change...
Definition: AMReX_ParticleContainer.H:307
ParticleContainer_impl(const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
Construct a particle container using a given Geometry, DistributionMapping, BoxArray and Vector of re...
Definition: AMReX_ParticleContainer.H:255
Long superParticleSize() const
Definition: AMReX_ParticleContainer.H:1250
typename SoA::IntVector IntVector
Definition: AMReX_ParticleContainer.H:191
int num_int_comm_comps
Definition: AMReX_ParticleContainer.H:1457
size_t particle_size
Definition: AMReX_ParticleContainer.H:1456
void SortParticlesByBin(IntVect bin_size)
Sort the particles on each tile by groups of cells, given an IntVect bin_size.
Definition: AMReX_ParticleContainerI.H:1220
virtual void particlePostLocate(ParticleType &, const ParticleLocData &, const int)
Definition: AMReX_ParticleContainer.H:1437
void CreateGhostParticles(int level, int ngrow, AoS &ghosts) const
Create ghost particles for a given level that are copies of particles near coarse->fine boundaries in...
void Define(const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
Define a default-constructed ParticleContainer using a ParGDB object. Multi-level version.
Definition: AMReX_ParticleContainer.H:339
void SetSoACompileTimeNames(std::vector< std::string > const &rdata_name, std::vector< std::string > const &idata_name)
Definition: AMReX_ParticleContainerI.H:84
void AddParticlesAtLevel(AoS &particles, int level, int nGrow=0)
Add particles from a pbox to the grid at this level.
void SortParticlesByCell()
Sort the particles on each tile by cell, using Fortran ordering.
Definition: AMReX_ParticleContainerI.H:1211
ParticleTileType & ParticlesAt(int lev, int grid, int tile)
Return the ParticleTile for level "lev", grid "grid" and tile "tile." Non-const version.
Definition: AMReX_ParticleContainer.H:1071
void ResizeRuntimeIntComp(int new_size, bool communicate)
Definition: AMReX_ParticleContainerI.H:2567
void InitFromBinaryMetaFile(const std::string &file, int extradata)
Definition: AMReX_ParticleInit.H:929
void resizeData() override
This resizes the vector of dummy MultiFabs used by the ParticleContainer for the current number of le...
Definition: AMReX_ParticleContainerI.H:339
void AddRealComp(std::string const &name, int communicate=1)
Definition: AMReX_ParticleContainer.H:1252
int num_real_comm_comps
Definition: AMReX_ParticleContainer.H:1457
ParticleTileType & ParticlesAt(int lev, const Iterator &iter)
Return the ParticleTile for level "lev" and Iterator "iter". Non-const version.
Definition: AMReX_ParticleContainer.H:1121
IntVect Index(const P &p, int lev) const
void Redistribute(int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
Redistribute puts all the particles back in the right places (for some value of right)
Definition: AMReX_ParticleContainerI.H:1100
Long IncrementWithTotal(MultiFab &mf, int level, bool local=false)
Definition: AMReX_ParticleContainerI.H:632
const Vector< ParticleLevel > & GetParticles() const
Return the underlying Vector (over AMR levels) of ParticleLevels. Const version.
Definition: AMReX_ParticleContainer.H:986
void InitOnePerCell(Real x_off, Real y_off, Real z_off, const ParticleInitData &pdata)
This initializes the particle container with one particle per cell, where the other particle data and...
Definition: AMReX_ParticleInit.H:1463
const ParticleTileType & ParticlesAt(int lev, int grid, int tile) const
Return the ParticleTile for level "lev", grid "grid" and tile "tile." Const version.
Definition: AMReX_ParticleContainer.H:1046
void SortParticlesForDeposition(IntVect idx_type)
Sort particles on each tile such that particles adjacent in memory are likely to map to adjacent cell...
Definition: AMReX_ParticleContainerI.H:1253
void AddRealComp(int communicate=1)
Definition: AMReX_ParticleContainer.H:1275
~ParticleContainer_impl() override=default
void CheckpointPost()
Definition: AMReX_ParticleIO.H:480
T_ParticleType ParticleType
Definition: AMReX_ParticleContainer.H:147
bool m_runtime_comps_defined
Definition: AMReX_ParticleContainer.H:1452
Vector< Vector< int > > countPrePost
Definition: AMReX_ParticleContainer.H:1369
void InitRandomPerBox(Long icount, ULong iseed, const ParticleInitData &pdata)
This initializes the container with icount randomly distributed particles per box,...
Definition: AMReX_ParticleInit.H:1362
bool usePrePost
Definition: AMReX_ParticleContainer.H:1362
void SetParticleSize()
Definition: AMReX_ParticleContainerI.H:11
A Box with real dimensions. A RealBox is OK iff volume >= 0.
Definition: AMReX_RealBox.H:21
A Descriptor of the Real Type.
Definition: AMReX_FabConv.H:105
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition: AMReX_Vector.H:27
Long size() const noexcept
Definition: AMReX_Vector.H:50
AMREX_GPU_HOST_DEVICE Long size(T const &b) noexcept
integer version
Definition: AMReX_GpuRange.H:26
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
Definition: AMReX_Amr.cpp:49
amrex::ArenaAllocator< T > DefaultAllocator
Definition: AMReX_GpuAllocators.H:194
int verbose
Definition: AMReX_DistributionMapping.cpp:36
Definition: AMReX_ParticleUtil.H:432
A struct used for communicating particle data across processes during multi-level operations.
Definition: AMReX_ParticleContainer.H:79
ParticleReal RealType
Definition: AMReX_ParticleContainer.H:80
int m_grid
Definition: AMReX_ParticleContainer.H:82
RealType m_data[1+AMREX_SPACEDIM]
Definition: AMReX_ParticleContainer.H:84
IntVect m_cell
Definition: AMReX_ParticleContainer.H:83
int m_lev
Definition: AMReX_ParticleContainer.H:81
A struct used to pass initial data into the various Init methods. This struct is used to pass initial...
Definition: AMReX_ParticleContainer.H:116
std::array< int, NStructInt > int_struct_data
Definition: AMReX_ParticleContainer.H:118
std::array< int, NArrayInt > int_array_data
Definition: AMReX_ParticleContainer.H:120
std::array< double, NArrayReal > real_array_data
Definition: AMReX_ParticleContainer.H:119
std::array< double, NStructReal > real_struct_data
Definition: AMReX_ParticleContainer.H:117
A struct used for storing a particle's position in the AMR hierarchy.
Definition: AMReX_ParticleContainer.H:91
Box m_grown_gridbox
Definition: AMReX_ParticleContainer.H:98
IntVect m_cell
Definition: AMReX_ParticleContainer.H:95
int m_grid
Definition: AMReX_ParticleContainer.H:93
int m_tile
Definition: AMReX_ParticleContainer.H:94
int m_lev
Definition: AMReX_ParticleContainer.H:92
Box m_tilebox
Definition: AMReX_ParticleContainer.H:97
Box m_gridbox
Definition: AMReX_ParticleContainer.H:96
Definition: AMReX_ParticleTile.H:693
std::conditional_t< ParticleType::is_soa_particle, StructOfArrays< NArrayReal, NArrayInt, Allocator, true >, StructOfArrays< NArrayReal, NArrayInt, Allocator, false > > SoA
Definition: AMReX_ParticleTile.H:716
std::conditional_t< ParticleType::is_soa_particle, ThisParticleTileHasNoAoS, ArrayOfStructs< ParticleType, Allocator > > AoS
Definition: AMReX_ParticleTile.H:710
The struct used to store particles.
Definition: AMReX_Particle.H:295
ParticleReal RealType
The floating point type used for the particles.
Definition: AMReX_Particle.H:307