Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
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>
24#include <AMReX_Particle.H>
25#include <AMReX_ParticleTile.H>
27#include <AMReX_TypeTraits.H>
28#include <AMReX_GpuContainers.H>
29#include <AMReX_ParticleUtil.H>
34#include <AMReX_Scan.H>
35#include <AMReX_DenseBins.H>
36#include <AMReX_SparseBins.H>
38#include <AMReX_ParticleMesh.H>
39#include <AMReX_OpenMP.H>
40#include <AMReX_ParIter.H>
41
42#ifdef AMREX_LAZY
43#include <AMReX_Lazy.H>
44#endif
45
46#ifdef AMREX_USE_OMP
47#include <omp.h>
48#endif
49
50#ifdef AMREX_USE_HDF5
51#include <hdf5.h>
52#endif
53
54#include <algorithm>
55#include <array>
56#include <cstring>
57#include <fstream>
58#include <iostream>
59#include <limits>
60#include <map>
61#include <memory>
62#include <numeric>
63#include <random>
64#include <stdexcept>
65#include <string>
66#include <tuple>
67#include <type_traits>
68#include <utility>
69#include <vector>
70
71namespace amrex {
72
73#ifdef AMREX_USE_HDF5_ASYNC
74 extern hid_t es_par_g;
75#endif
81{
83 int m_lev;
84 int m_grid;
86 RealType m_data[1 + AMREX_SPACEDIM];
87};
88
102
117template<int NStructReal, int NStructInt, int NArrayReal, int NArrayInt>
119{
120 std::array<double, NStructReal> real_struct_data;
121 std::array<int, NStructInt > int_struct_data;
122 std::array<double, NArrayReal > real_array_data;
123 std::array<int, NArrayInt > int_array_data;
124};
125
126template <bool is_const, typename T_ParticleType, int NArrayReal, int NArrayInt,
127 template<class> class Allocator, class CellAssignor>
128class ParIterBase_impl;
129
144template <typename T_ParticleType, int T_NArrayReal, int T_NArrayInt,
145 template<class> class Allocator=DefaultAllocator,
146 class T_CellAssignor=DefaultAssignor>
147
149{
150public:
151 using ParticleType = T_ParticleType;
152 using ConstParticleType = typename ParticleType::ConstType;
153 using CellAssignor = T_CellAssignor;
154
156 static constexpr int NStructReal = ParticleType::NReal;
158 static constexpr int NStructInt = ParticleType::NInt;
160 static constexpr int NArrayReal = T_NArrayReal;
162 static constexpr int NArrayInt = T_NArrayInt;
164
165private:
166 friend class ParIterBase_impl<true,ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>;
167 friend class ParIterBase_impl<false,ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>;
168
169public:
171 template <typename T>
172 using AllocatorType = Allocator<T>;
177
178#ifdef AMREX_SINGLE_PRECISION_PARTICLES
180#else
182#endif
183
184 static constexpr bool is_rtsoa_pc = ParticleType::is_rtsoa_particle;
185
187 using ParticleTileType = std::conditional_t<
191 >;
193
196 using ParticleLevel = std::map<std::pair<int, int>, ParticleTileType>;
197 using PTDType = typename ParticleTileType::ParticleTileDataType;
198 using ConstPTDType = typename ParticleTileType::ConstParticleTileDataType;
199 using AoS = typename ParticleTileType::AoS;
200 using SoA = typename ParticleTileType::SoA;
201
202 using RealVector = typename SoA::RealVector;
203 using IntVector = typename SoA::IntVector;
204 using ParticleVector = typename AoS::ParticleVector;
208
209 static constexpr bool has_polymorphic_allocator =
211
212 static_assert(!ParticleType::is_soa_particle || T_NArrayReal >= AMREX_SPACEDIM || is_rtsoa_pc,
213 "Pure SoA particle containers require at least AMREX_SPACEDIM real components for positions");
214
218 :
220 {
221 Initialize ();
222 }
223
239
248 const DistributionMapping & dmap,
249 const BoxArray & ba)
250 :
251 ParticleContainerBase(geom, dmap, ba)
252 {
253 Initialize ();
256 }
257
268 const Vector<DistributionMapping> & dmap,
269 const Vector<BoxArray> & ba,
270 const Vector<int> & rr)
271 :
272 ParticleContainerBase(geom, dmap, ba, rr)
273 {
274 Initialize ();
277 }
278
288 const Vector<DistributionMapping> & dmap,
289 const Vector<BoxArray> & ba,
290 const Vector<IntVect> & rr)
291 :
292 ParticleContainerBase(geom, dmap, ba, rr)
293 {
294 Initialize ();
297 }
298
299 ~ParticleContainer_impl () override = default;
300
303
305 ParticleContainer_impl& operator= ( ParticleContainer_impl && ) noexcept = default;
306
307
315 void Define (ParGDBBase* gdb)
316 {
318 reserveData();
319 resizeData();
320 }
321
329 void Define (const Geometry & geom,
330 const DistributionMapping & dmap,
331 const BoxArray & ba)
332 {
333 this->ParticleContainerBase::Define(geom, dmap, ba);
334 reserveData();
335 resizeData();
336 }
337
347 void Define (const Vector<Geometry> & geom,
348 const Vector<DistributionMapping> & dmap,
349 const Vector<BoxArray> & ba,
350 const Vector<int> & rr)
351 {
352 this->ParticleContainerBase::Define(geom, dmap, ba, rr);
353 reserveData();
354 resizeData();
355 }
356
366 void Define (const Vector<Geometry> & geom,
367 const Vector<DistributionMapping> & dmap,
368 const Vector<BoxArray> & ba,
369 const Vector<IntVect> & rr)
370 {
371 this->ParticleContainerBase::Define(geom, dmap, ba, rr);
372 reserveData();
373 resizeData();
374 }
375
377 int numLocalTilesAtLevel (int lev) const {
378 return (lev < m_particles.size()) ? m_particles[lev].size() : 0;
379 }
380
385 void reserveData () override;
386
393 void resizeData () override;
394
395 void InitFromAsciiFile (const std::string& file, int extradata,
396 const IntVect* Nrep = nullptr);
397
398 void InitFromBinaryFile (const std::string& file, int extradata);
399
400 void InitFromBinaryMetaFile (const std::string& file, int extradata);
401
418 void InitRandom (Long icount, ULong iseed,
419 const ParticleInitData& pdata,
420 bool serialize = false, RealBox bx = RealBox());
421
422
438 void InitRandomPerBox (Long icount, ULong iseed, const ParticleInitData& pdata);
439
440
455 void InitOnePerCell (Real x_off, Real y_off, Real z_off,
456 const ParticleInitData& pdata);
457
458
470 void InitNRandomPerCell (int n_per_cell, const ParticleInitData& pdata);
471
472 void Increment (MultiFab& mf, int level);
473
474 Long IncrementWithTotal (MultiFab& mf, int level, bool local = false);
475
513 void Redistribute (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
514 bool remove_negative=true);
515
516
528 template <class index_type>
529 void ReorderParticles (int lev, const MFIter& mfi, const index_type* permutations);
530
545 void SortParticlesForDeposition (IntVect idx_type);
546
550 void SortParticlesByCell ();
551
558 void SortParticlesByBin (IntVect bin_size);
559
570 bool OK (int lev_min = 0, int lev_max = -1, int nGrow = 0) const;
571
572 std::array<Long, 3> ByteSpread () const;
573
574 std::array<Long, 3> PrintCapacity () const;
575
576 void ShrinkToFit ();
577
588 Long NumberOfParticlesAtLevel (int level, bool only_valid = true, bool only_local = false) const;
589
590 Vector<Long> NumberOfParticlesInGrid (int level, bool only_valid = true, bool only_local = false) const;
591
595 template <typename I, std::enable_if_t<std::is_integral_v<I> && (sizeof(I) >= sizeof(Long)), int> = 0>
596 void CapacityOfParticlesInGrid (LayoutData<I>& mem, int lev) const;
597
606 Long TotalNumberOfParticles (bool only_valid=true, bool only_local=false) const;
607
608
616 void RemoveParticlesAtLevel (int level);
617
619
627 void CreateVirtualParticles (int level, AoS& virts) const;
628
637 void CreateGhostParticles (int level, int ngrow, AoS& ghosts) const;
638
646 void AddParticlesAtLevel (AoS& particles, int level, int nGrow=0);
647
655 void CreateVirtualParticles (int level, ParticleTileType& virts) const;
656
665 void CreateGhostParticles (int level, int ngrow, ParticleTileType& ghosts) const;
666
674 void AddParticlesAtLevel (ParticleTileType& particles, int level, int nGrow=0);
675
676
680 void clearParticles ();
681
682
691 template <class PCType,
692 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0>
693 void copyParticles (const PCType& other, bool local=false);
694
702 template <class PCType,
703 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0>
704 void addParticles (const PCType& other, bool local=false);
705
720 template <class F, class PCType,
721 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0,
722 std::enable_if_t<! std::is_integral_v<F>, int> bar = 0>
723 void copyParticles (const PCType& other, F&&f, bool local=false);
724
738 template <class F, class PCType,
739 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0,
740 std::enable_if_t<! std::is_integral_v<F>, int> bar = 0>
741 void addParticles (const PCType& other, F const& f, bool local=false);
742
750 void WriteParticleRealData (void* data, size_t size, std::ostream& os) const;
751
759 void ReadParticleRealData (void* data, size_t size, std::istream& is);
760
769 void Checkpoint (const std::string& dir, const std::string& name,
770 const Vector<std::string>& real_comp_names = Vector<std::string>(),
771 const Vector<std::string>& int_comp_names = Vector<std::string>()) const
772 {
773 Checkpoint(dir, name, true, real_comp_names, int_comp_names);
774 }
775
781 void Checkpoint (const std::string& dir, const std::string& name, bool is_checkpoint,
782 const Vector<std::string>& real_comp_names = Vector<std::string>(),
783 const Vector<std::string>& int_comp_names = Vector<std::string>()) const;
784
797 void Checkpoint (const std::string& dir, const std::string& name,
798 const Vector<int>& write_real_comp,
799 const Vector<int>& write_int_comp,
800 const Vector<std::string>& real_comp_names,
801 const Vector<std::string>& int_comp_names) const;
802
817 template <class F>
818 void WriteBinaryParticleData (const std::string& dir,
819 const std::string& name,
820 const Vector<int>& write_real_comp,
821 const Vector<int>& write_int_comp,
822 const Vector<std::string>& real_comp_names,
823 const Vector<std::string>& int_comp_names,
824 F&& f, bool is_checkpoint=false) const;
825
826 void CheckpointPre ();
827
828 void CheckpointPost ();
829
836 void Restart (const std::string& dir, const std::string& file);
837
845 void Restart (const std::string& dir, const std::string& file, bool is_checkpoint);
846
853 void WritePlotFile (const std::string& dir, const std::string& name) const;
854
866 template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>* = nullptr>
867 void WritePlotFile (const std::string& dir, const std::string& name, F&& f) const;
868
877 void WritePlotFile (const std::string& dir, const std::string& name,
878 const Vector<std::string>& real_comp_names,
879 const Vector<std::string>& int_comp_names) const;
880
894 template <class F>
895 void WritePlotFile (const std::string& dir, const std::string& name,
896 const Vector<std::string>& real_comp_names,
897 const Vector<std::string>& int_comp_names, F&& f) const;
898
907 void WritePlotFile (const std::string& dir, const std::string& name,
908 const Vector<std::string>& real_comp_names) const;
909
923 template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>* = nullptr>
924 void WritePlotFile (const std::string& dir, const std::string& name,
925 const Vector<std::string>& real_comp_names, F&& f) const;
926
936 void WritePlotFile (const std::string& dir,
937 const std::string& name,
938 const Vector<int>& write_real_comp,
939 const Vector<int>& write_int_comp) const;
940
955 template <class F>
956 void WritePlotFile (const std::string& dir,
957 const std::string& name,
958 const Vector<int>& write_real_comp,
959 const Vector<int>& write_int_comp, F&& f) const;
960
974 void WritePlotFile (const std::string& dir,
975 const std::string& name,
976 const Vector<int>& write_real_comp,
977 const Vector<int>& write_int_comp,
978 const Vector<std::string>& real_comp_names,
979 const Vector<std::string>& int_comp_names) const;
980
999 template <class F>
1000 void WritePlotFile (const std::string& dir,
1001 const std::string& name,
1002 const Vector<int>& write_real_comp,
1003 const Vector<int>& write_int_comp,
1004 const Vector<std::string>& real_comp_names,
1005 const Vector<std::string>& int_comp_names,
1006 F&& f) const;
1007
1008 void WritePlotFilePre ();
1009
1010 void WritePlotFilePost ();
1011
1012 void WriteAsciiFile (const std::string& file);
1013
1018 const Vector<ParticleLevel>& GetParticles () const { return m_particles; }
1019
1024 Vector <ParticleLevel>& GetParticles () { return m_particles; }
1025
1039 const ParticleLevel& GetParticles (int lev) const { return m_particles[lev]; }
1040
1054 ParticleLevel & GetParticles (int lev) { return m_particles[lev]; }
1055
1078 const ParticleTileType& ParticlesAt (int lev, int grid, int tile) const
1079 { return m_particles[lev].at(std::make_pair(grid, tile)); }
1080
1103 ParticleTileType& ParticlesAt (int lev, int grid, int tile)
1104 { return m_particles[lev].at(std::make_pair(grid, tile)); }
1105
1127 template <class Iterator>
1128 const ParticleTileType& ParticlesAt (int lev, const Iterator& iter) const
1129 { return ParticlesAt(lev, iter.index(), iter.LocalTileIndex()); }
1130
1152 template <class Iterator>
1153 ParticleTileType& ParticlesAt (int lev, const Iterator& iter)
1154 { return ParticlesAt(lev, iter.index(), iter.LocalTileIndex()); }
1155
1178 ParticleTileType& DefineAndReturnParticleTile (int lev, int grid, int tile)
1179 {
1180 m_particles[lev][std::make_pair(grid, tile)].define(
1182 &m_soa_rdata_names, &m_soa_idata_names,
1183 arena()
1184 );
1185 return ParticlesAt(lev, grid, tile);
1186 }
1187
1209 template <class Iterator>
1210 ParticleTileType& DefineAndReturnParticleTile (int lev, const Iterator& iter)
1211 {
1212 auto index = std::make_pair(iter.index(), iter.LocalTileIndex());
1213 m_particles[lev][index].define(
1215 &m_soa_rdata_names, &m_soa_idata_names,
1216 arena()
1217 );
1218 return ParticlesAt(lev, iter);
1219 }
1220
1231 void AssignDensity (int rho_index,
1232 Vector<std::unique_ptr<MultiFab> >& mf_to_be_filled,
1233 int lev_min, int ncomp, int finest_level, int ngrow=2) const;
1234
1235 void AssignCellDensitySingleLevel (int rho_index, MultiFab& mf, int level,
1236 int ncomp=1, int particle_lvl_offset = 0) const;
1237
1238 template <typename P, typename Assignor=CellAssignor>
1239 IntVect Index (const P& p, int lev) const;
1240
1250 ParticleLocData Reset (ParticleType& prt, bool update, bool verbose=true,
1251 ParticleLocData pld = ParticleLocData()) const;
1252
1258 template <typename P>
1259 bool PeriodicShift (P& p) const;
1260
1262
1264
1265 void SetUsePrePost (bool tf) const {
1266 usePrePost = tf;
1267 }
1268 bool GetUsePrePost () const {
1269 return usePrePost;
1270 }
1271
1274
1275 void SetUseUnlink (bool tf) const {
1276 doUnlink = tf;
1277 }
1278
1279 bool GetUseUnlink () const {
1280 return doUnlink;
1281 }
1282
1283 void RedistributeCPU (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
1284 bool remove_negative=true);
1285
1286 void RedistributeGPU (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
1287 bool remove_negative=true);
1288
1289 void ReserveForRedistribute (ParticleCopyPlan const& plan);
1290
1291 Long superParticleSize() const { return superparticle_size; }
1292
1293 void AddRealComp (std::string const & name, int communicate=1)
1294 {
1295 // names must be unique
1296 auto const it = std::find(m_soa_rdata_names.begin(), m_soa_rdata_names.end(), name);
1297 if (it != m_soa_rdata_names.end()) {
1298 throw std::runtime_error("AddRealComp: name '" + name + "' is already present in the SoA.");
1299 }
1300 m_soa_rdata_names.push_back(name);
1301
1302 m_runtime_comps_defined = true;
1303 m_num_runtime_real++;
1304 h_redistribute_real_comp.push_back(communicate);
1306 this->resizeData();
1307
1308 // resize runtime SoA
1309 for (int lev = 0; lev < numLevels(); ++lev) {
1310 for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1311 auto& tile = DefineAndReturnParticleTile(lev, pti);
1312 auto np = tile.numParticles();
1313 if (np > 0) {
1314 auto& soa = tile.GetStructOfArrays();
1315 soa.resize(np);
1316 }
1317 }
1318 }
1319 }
1320
1321 void AddRealComp (int communicate=1)
1322 {
1323 AddRealComp(getDefaultCompNameReal<ParticleType>(NArrayReal+m_num_runtime_real), communicate);
1324 }
1325
1326 void AddIntComp (std::string const & name, int communicate=1)
1327 {
1328 // names must be unique
1329 auto const it = std::find(m_soa_idata_names.begin(), m_soa_idata_names.end(), name);
1330 if (it != m_soa_idata_names.end()) {
1331 throw std::runtime_error("AddIntComp: name '" + name + "' is already present in the SoA.");
1332 }
1333 m_soa_idata_names.push_back(name);
1334
1335 m_runtime_comps_defined = true;
1336 m_num_runtime_int++;
1337 h_redistribute_int_comp.push_back(communicate);
1339 this->resizeData();
1340
1341 // resize runtime SoA
1342 for (int lev = 0; lev < numLevels(); ++lev) {
1343 for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1344 auto& tile = DefineAndReturnParticleTile(lev, pti);
1345 auto np = tile.numParticles();
1346 if (np > 0) {
1347 auto& soa = tile.GetStructOfArrays();
1348 soa.resize(np);
1349 }
1350 }
1351 }
1352 }
1353
1354 void AddIntComp (int communicate=1)
1355 {
1356 AddIntComp(getDefaultCompNameInt<ParticleType>(NArrayInt+m_num_runtime_int), communicate);
1357 }
1358
1359 int NumRuntimeRealComps () const { return m_num_runtime_real; }
1360 int NumRuntimeIntComps () const { return m_num_runtime_int; }
1361
1362 int NumRealComps () const { return NArrayReal + NumRuntimeRealComps(); }
1363 int NumIntComps () const { return NArrayInt + NumRuntimeIntComps() ; }
1364
1370 void ResizeRuntimeRealComp (int new_size, bool communicate);
1371
1377 void ResizeRuntimeIntComp (int new_size, bool communicate);
1378
1380 template <template<class> class NewAllocator=amrex::DefaultAllocator>
1382
1393 template <template<class> class NewAllocator=Allocator>
1395 make_alike () const
1396 {
1398
1399 std::vector<std::string> const real_names = this->GetRealSoANames();
1400 std::vector<std::string> real_ct_names(NArrayReal);
1401 for (int ic = 0; ic < NArrayReal; ++ic) { real_ct_names.at(ic) = real_names[ic]; }
1402
1403 std::vector<std::string> const int_names = this->GetIntSoANames();
1404 std::vector<std::string> int_ct_names(NArrayInt);
1405 for (int ic = 0; ic < NArrayInt; ++ic) { int_ct_names.at(ic) = int_names[ic]; }
1406
1407 tmp.SetSoACompileTimeNames(real_ct_names, int_ct_names);
1408
1409 // add runtime real comps to tmp
1410 for (int ic = 0; ic < this->NumRuntimeRealComps(); ++ic) {
1411 tmp.AddRealComp(real_names.at(ic + NArrayReal));
1412 }
1413
1414 // add runtime int comps to tmp
1415 for (int ic = 0; ic < this->NumRuntimeIntComps(); ++ic) {
1416 tmp.AddIntComp(int_names.at(ic + NArrayInt));
1417 }
1418
1421
1422 return tmp;
1423 }
1424
1427
1430 mutable bool usePrePost;
1431 mutable bool doUnlink;
1433 mutable int nOutFilesPrePost;
1439 mutable std::string HdrFileNamePrePost;
1441
1442protected:
1443
1457 template <typename P>
1458 bool Where (const P& prt, ParticleLocData& pld,
1459 int lev_min = 0, int lev_max = -1, int nGrow=0, int local_grid=-1) const;
1460
1461
1472 template <typename P>
1473 bool EnforcePeriodicWhere (P& prt, ParticleLocData& pld,
1474 int lev_min = 0, int lev_max = -1, int local_grid=-1) const;
1475
1476public:
1478
1479 void
1480 WriteParticles (int level, std::ofstream& ofs, int fnum,
1481 Vector<int>& which, Vector<int>& count, Vector<Long>& where,
1482 const Vector<int>& write_real_comp, const Vector<int>& write_int_comp,
1483 const Vector<std::map<std::pair<int, int>,FlagsVector>>& particle_io_flags, bool is_checkpoint) const;
1484#ifdef AMREX_USE_HDF5
1485#include "AMReX_ParticlesHDF5.H"
1486#endif
1487
1488public: // NOLINT(readability-redundant-access-specifiers)
1490 void SetSoACompileTimeNames (std::vector<std::string> const & rdata_name, std::vector<std::string> const & idata_name);
1491
1493 std::vector<std::string> GetRealSoANames () const {return m_soa_rdata_names;}
1494
1496 std::vector<std::string> GetIntSoANames () const {return m_soa_idata_names;}
1497
1503 bool HasRealComp (std::string const & name);
1504
1510 bool HasIntComp (std::string const & name);
1511
1519 int GetRealCompIndex (std::string const & name);
1520
1528 int GetIntCompIndex (std::string const & name);
1529
1530protected:
1531
1532 template <class RTYPE>
1533 void ReadParticles (int cnt, int grd, int lev, std::ifstream& ifs, int finest_level_in_file, bool convert_ids);
1534
1535 void SetParticleSize ();
1536
1538
1539private:
1540 virtual void particlePostLocate (ParticleType& /*p*/, const ParticleLocData& /*pld*/,
1541 const int /*lev*/) {}
1542
1543 virtual void correctCellVectors (int /*old_index*/, int /*new_index*/,
1544 int /*grid*/, const ParticleType& /*p*/) {}
1545
1546 void RedistributeMPI (std::map<int, Vector<char> >& not_ours,
1547 int lev_min = 0, int lev_max = 0, int nGrow = 0, int local=0);
1548
1549 template <typename P>
1550 void locateParticle (P& p, ParticleLocData& pld,
1551 int lev_min, int lev_max, int nGrow, int local_grid=-1) const;
1552
1553 void Initialize ();
1554
1555 bool m_runtime_comps_defined{false};
1556 int m_num_runtime_real{0};
1557 int m_num_runtime_int{0};
1558
1559 size_t particle_size, superparticle_size;
1560 int num_real_comm_comps, num_int_comm_comps;
1561 Vector<ParticleLevel> m_particles;
1562
1563 // names of both compile-time and runtime Real and Int SoA data
1564 std::vector<std::string> m_soa_rdata_names;
1565 std::vector<std::string> m_soa_idata_names;
1566};
1567
1569template <int T_NStructReal, int T_NStructInt, int T_NArrayReal, int T_NArrayInt, template<class> class Allocator, class CellAssignor>
1570using ParticleContainer = ParticleContainer_impl<Particle<T_NStructReal, T_NStructInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
1571
1572template <int T_NArrayReal, int T_NArrayInt, template<class> class Allocator=DefaultAllocator, class CellAssignor=DefaultAssignor>
1573using ParticleContainerPureSoA = ParticleContainer_impl<SoAParticle<T_NArrayReal, T_NArrayInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
1574
1575template <class RType=ParticleReal, class IType=int, class CellAssignor=DefaultAssignor>
1577
1578}
1579
1580#include "AMReX_ParticleInit.H"
1582#include "AMReX_ParticleIO.H"
1583
1584#ifdef AMREX_USE_HDF5
1585#include "AMReX_ParticleHDF5.H"
1586#endif
1587
1588#endif /*_PARTICLES_H_*/
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Definition AMReX_GpuAllocators.H:121
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
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:43
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:74
a one-thingy-per-box distributed object
Definition AMReX_LayoutData.H:13
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
Definition AMReX_ParIter.H:147
Definition AMReX_ParGDB.H:13
Definition AMReX_ParIter.H:35
Definition AMReX_ParIter.H:118
Definition AMReX_ParticleContainerBase.H:23
ParGDBBase * m_gdb
Definition AMReX_ParticleContainerBase.H:278
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
Arena * arena() const
Definition AMReX_ParticleContainerBase.H:248
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:149
ParticleTileType & DefineAndReturnParticleTile(int lev, const Iterator &iter)
Define and return the ParticleTile for level "lev", and Iterator "iter".
Definition AMReX_ParticleContainer.H:1210
void WritePlotFilePre()
Definition AMReX_ParticleIO.H:565
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:36
typename ParticleTileType::ConstParticleTileDataType ConstPTDType
Definition AMReX_ParticleContainer.H:198
void addParticles(const PCType &other, bool local=false)
Add particles from other to this ParticleContainer. local controls whether or not to call Redistribut...
Definition AMReX_ParticleContainerI.H:1155
void SetUseUnlink(bool tf) const
Definition AMReX_ParticleContainer.H:1275
Vector< Long > nParticlesAtLevelPrePost
Definition AMReX_ParticleContainer.H:1435
bool GetLevelDirectoriesCreated() const
Definition AMReX_ParticleContainer.H:1263
IntVect Index(const P &p, int lev) const
Definition AMReX_ParticleContainerI.H:201
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:1237
void AddIntComp(std::string const &name, int communicate=1)
Definition AMReX_ParticleContainer.H:1326
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.
Definition AMReX_ParticleContainerI.H:1268
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...
Definition AMReX_ParticleContainerI.H:297
void ResizeRuntimeIntComp(int new_size, bool communicate)
Definition AMReX_ParticleContainerI.H:2486
bool doUnlink
Definition AMReX_ParticleContainer.H:1431
Vector< Vector< int > > whichPrePost
Definition AMReX_ParticleContainer.H:1436
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:287
std::string HdrFileNamePrePost
Definition AMReX_ParticleContainer.H:1439
Vector< Vector< Long > > wherePrePost
Definition AMReX_ParticleContainer.H:1438
RealDescriptor ParticleRealDescriptor
Definition AMReX_ParticleContainer.H:181
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:2259
ParticleTileType & ParticlesAt(int lev, const Iterator &iter)
Return the ParticleTile for level "lev" and Iterator "iter". Non-const version.
Definition AMReX_ParticleContainer.H:1153
int nOutFilesPrePost
Definition AMReX_ParticleContainer.H:1433
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:28
void InitRandomPerBox(Long icount, ULong iseed, const ParticleInitData &pdata)
This initializes the container with icount randomly distributed particles per box,...
Definition AMReX_ParticleInit.H:1364
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:769
bool PeriodicShift(P &p) const
Returns true if the particle was shifted.
Definition AMReX_ParticleContainerI.H:376
std::map< std::pair< int, int >, ParticleTileType > ParticleLevel
Definition AMReX_ParticleContainer.H:196
bool levelDirectoriesCreated
Variables for i/o optimization saved for pre and post checkpoint.
Definition AMReX_ParticleContainer.H:1429
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...
Definition AMReX_ParticleContainerI.H:216
void RedistributeCPU(int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
Definition AMReX_ParticleContainerI.H:1600
static constexpr int NArrayInt
Number of extra integer components stored in struct-of-array form.
Definition AMReX_ParticleContainer.H:162
void InitFromAsciiFile(const std::string &file, int extradata, const IntVect *Nrep=nullptr)
Definition AMReX_ParticleInit.H:39
static constexpr int NArrayReal
Number of extra Real components stored in struct-of-array form.
Definition AMReX_ParticleContainer.H:160
const ParticleTileType & ParticlesAt(int lev, const Iterator &iter) const
Return the ParticleTile for level "lev" and Iterator "iter". Const version.
Definition AMReX_ParticleContainer.H:1128
void WritePlotFilePost()
Definition AMReX_ParticleIO.H:575
bool GetUsePrePost() const
Definition AMReX_ParticleContainer.H:1268
void SetSoACompileTimeNames(std::vector< std::string > const &rdata_name, std::vector< std::string > const &idata_name)
Definition AMReX_ParticleContainerI.H:109
typename ParticleTileType::AoS AoS
Definition AMReX_ParticleContainer.H:199
ParticleInitType< NStructReal, NStructInt, NArrayReal, NArrayInt > ParticleInitData
Definition AMReX_ParticleContainer.H:192
typename SoA::RealVector RealVector
Definition AMReX_ParticleContainer.H:202
typename ParticleType::ConstType ConstParticleType
Definition AMReX_ParticleContainer.H:152
ParticleContainer_impl(const ParticleContainer_impl &)=delete
void WriteAsciiFile(const std::string &file)
Definition AMReX_ParticleIO.H:1141
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:247
void SetParticleSize()
Definition AMReX_ParticleContainerI.H:16
void CheckpointPost()
Definition AMReX_ParticleIO.H:509
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:1549
void SetLevelDirectoriesCreated(bool tf)
Definition AMReX_ParticleContainer.H:1261
void RemoveParticlesAtLevel(int level)
The Following methods are for managing Virtual and Ghost Particles.
Definition AMReX_ParticleContainerI.H:757
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:1078
typename AoS::ParticleVector ParticleVector
Definition AMReX_ParticleContainer.H:204
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:329
static constexpr int NStructInt
Number of extra integer components in the particle struct.
Definition AMReX_ParticleContainer.H:158
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:1334
ParticleContainer_impl()
Default constructor - construct an empty particle container that has no concept of a level hierarchy....
Definition AMReX_ParticleContainer.H:217
void SetUsePrePost(bool tf) const
Definition AMReX_ParticleContainer.H:1265
int NumRealComps() const
Definition AMReX_ParticleContainer.H:1362
Long GetMaxNextIDPrePost() const
Definition AMReX_ParticleContainer.H:1272
void clearParticles()
Clear all the particles in this container. This does not free memory.
Definition AMReX_ParticleContainerI.H:1128
void ReadParticles(int cnt, int grd, int lev, std::ifstream &ifs, int finest_level_in_file, bool convert_ids)
Definition AMReX_ParticleIO.H:960
void Increment(MultiFab &mf, int level)
Definition AMReX_ParticleContainerI.H:720
Vector< int > h_redistribute_int_comp
Definition AMReX_ParticleContainer.H:1426
typename ParticleTileType::SoA SoA
Definition AMReX_ParticleContainer.H:200
void AssignCellDensitySingleLevel(int rho_index, MultiFab &mf, int level, int ncomp=1, int particle_lvl_offset=0) const
Definition AMReX_ParticleContainerI.H:2320
int NumRuntimeIntComps() const
Definition AMReX_ParticleContainer.H:1360
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...
Definition AMReX_ParticleContainerI.H:1083
Long GetNParticlesPrePost() const
Definition AMReX_ParticleContainer.H:1273
void SortParticlesByBin(IntVect bin_size)
Sort the particles on each tile by groups of cells, given an IntVect bin_size.
Definition AMReX_ParticleContainerI.H:1301
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:1465
bool HasIntComp(std::string const &name)
Definition AMReX_ParticleContainerI.H:148
T_CellAssignor CellAssignor
Definition AMReX_ParticleContainer.H:153
const Vector< ParticleLevel > & GetParticles() const
Return the underlying Vector (over AMR levels) of ParticleLevels. Const version.
Definition AMReX_ParticleContainer.H:1018
const ParticleLevel & GetParticles(int lev) const
Return the ParticleLevel for level "lev". Const version.
Definition AMReX_ParticleContainer.H:1039
std::vector< std::string > GetRealSoANames() const
Definition AMReX_ParticleContainer.H:1493
DenseBins< typename ParticleTileType::ParticleTileDataType > m_bins
Definition AMReX_ParticleContainer.H:1537
void CheckpointPre()
Definition AMReX_ParticleIO.H:452
void ShrinkToFit()
Definition AMReX_ParticleContainerI.H:700
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:231
bool HasRealComp(std::string const &name)
Definition AMReX_ParticleContainerI.H:140
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:970
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:366
void resizeData() override
This resizes the vector of dummy MultiFabs used by the ParticleContainer for the current number of le...
Definition AMReX_ParticleContainerI.H:435
int NumRuntimeRealComps() const
Definition AMReX_ParticleContainer.H:1359
void ResizeRuntimeRealComp(int new_size, bool communicate)
Definition AMReX_ParticleContainerI.H:2460
void AddIntComp(int communicate=1)
Definition AMReX_ParticleContainer.H:1354
Long nparticlesPrePost
Definition AMReX_ParticleContainer.H:1434
Vector< std::string > filePrefixPrePost
Definition AMReX_ParticleContainer.H:1440
int numLocalTilesAtLevel(int lev) const
The total number of tiles on this rank on this level.
Definition AMReX_ParticleContainer.H:377
Long IncrementWithTotal(MultiFab &mf, int level, bool local=false)
Definition AMReX_ParticleContainerI.H:747
static constexpr int NStructReal
Number of extra Real components in the particle struct.
Definition AMReX_ParticleContainer.H:156
Long NumberOfParticlesAtLevel(int level, bool only_valid=true, bool only_local=false) const
Returns # of particles at specified the level.
Definition AMReX_ParticleContainerI.H:551
void SortParticlesByCell()
Sort the particles on each tile by cell, using Fortran ordering.
Definition AMReX_ParticleContainerI.H:1292
Vector< int > h_redistribute_real_comp
Definition AMReX_ParticleContainer.H:1425
typename Particle< NStructReal, NStructInt >::RealType RealType
The type of the Real data.
Definition AMReX_ParticleContainer.H:176
static constexpr bool is_rtsoa_pc
Definition AMReX_ParticleContainer.H:184
bool GetUseUnlink() const
Definition AMReX_ParticleContainer.H:1279
int GetRealCompIndex(std::string const &name)
Definition AMReX_ParticleContainerI.H:161
Long maxnextidPrePost
Definition AMReX_ParticleContainer.H:1432
void InitFromBinaryMetaFile(const std::string &file, int extradata)
Definition AMReX_ParticleInit.H:931
int NumIntComps() const
Definition AMReX_ParticleContainer.H:1363
void ReserveForRedistribute(ParticleCopyPlan const &plan)
Definition AMReX_ParticleContainerI.H:1560
typename ParticleTileType::ParticleTileDataType PTDType
Definition AMReX_ParticleContainer.H:197
void RedistributeGPU(int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
Definition AMReX_ParticleContainerI.H:1364
ParticleContainer_impl & operator=(const ParticleContainer_impl &)=delete
Allocator< T > AllocatorType
The memory allocator in use.
Definition AMReX_ParticleContainer.H:172
ParticleLevel & GetParticles(int lev)
Return the ParticleLevel for level "lev". Non-const version.
Definition AMReX_ParticleContainer.H:1054
int GetIntCompIndex(std::string const &name)
Definition AMReX_ParticleContainerI.H:183
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:1178
std::vector< std::string > GetIntSoANames() const
Definition AMReX_ParticleContainer.H:1496
void RemoveParticlesNotAtFinestLevel()
Definition AMReX_ParticleContainerI.H:771
Vector< Long > NumberOfParticlesInGrid(int level, bool only_valid=true, bool only_local=false) const
Definition AMReX_ParticleContainerI.H:496
void Define(ParGDBBase *gdb)
Define a default-constructed ParticleContainer using a ParGDB object. The container will track change...
Definition AMReX_ParticleContainer.H:315
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:267
void Restart(const std::string &dir, const std::string &file)
Restart from checkpoint.
Definition AMReX_ParticleIO.H:653
Long superParticleSize() const
Definition AMReX_ParticleContainer.H:1291
typename SoA::IntVector IntVector
Definition AMReX_ParticleContainer.H:203
void copyParticles(const PCType &other, bool local=false)
Copy particles from other to this ParticleContainer. Will clear all the particles from this container...
Definition AMReX_ParticleContainerI.H:1144
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:392
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:347
void CapacityOfParticlesInGrid(LayoutData< I > &mem, int lev) const
Return capacity of memory for particles at specific grid.
Definition AMReX_ParticleContainerI.H:593
std::array< Long, 3 > PrintCapacity() const
Definition AMReX_ParticleContainerI.H:660
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:838
static constexpr bool has_polymorphic_allocator
Definition AMReX_ParticleContainer.H:209
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 >, FlagsVector > > &particle_io_flags, bool is_checkpoint) const
Definition AMReX_ParticleIO.H:585
std::conditional_t< is_rtsoa_pc, ParticleTileRT< typename ParticleType::RealType, typename ParticleType::IntType >, ParticleTile< ParticleType, NArrayReal, NArrayInt, Allocator > > ParticleTileType
Definition AMReX_ParticleContainer.H:191
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:1103
void InitFromBinaryFile(const std::string &file, int extradata)
Definition AMReX_ParticleInit.H:485
void AddRealComp(std::string const &name, int communicate=1)
Definition AMReX_ParticleContainer.H:1293
std::array< Long, 3 > ByteSpread() const
Definition AMReX_ParticleContainerI.H:616
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:14
ContainerLike< NewAllocator > make_alike() const
Definition AMReX_ParticleContainer.H:1395
void AddRealComp(int communicate=1)
Definition AMReX_ParticleContainer.H:1321
void AddParticlesAtLevel(AoS &particles, int level, int nGrow=0)
Add particles from a pbox to the grid at this level.
Definition AMReX_ParticleContainerI.H:2274
~ParticleContainer_impl() override=default
Long TotalNumberOfParticles(bool only_valid=true, bool only_local=false) const
Returns # of particles at all levels.
Definition AMReX_ParticleContainerI.H:481
void reserveData() override
This reserves data in the vector of dummy MultiFabs used by the ParticleContainer for the maximum num...
Definition AMReX_ParticleContainerI.H:426
Vector< ParticleLevel > & GetParticles()
Return the underlying Vector (over AMR levels) of ParticleLevels. Non-const version.
Definition AMReX_ParticleContainer.H:1024
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.
Definition AMReX_ParticleIO.H:428
T_ParticleType ParticleType
Definition AMReX_ParticleContainer.H:151
Vector< Vector< int > > countPrePost
Definition AMReX_ParticleContainer.H:1437
bool usePrePost
Definition AMReX_ParticleContainer.H:1430
Definition AMReX_GpuAllocators.H:156
A Box with real dimensions.
Definition AMReX_RealBox.H:28
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:28
Long size() const noexcept
Definition AMReX_Vector.H:53
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_ulong ULong
Unsigned integer type guaranteed to be wider than unsigned int.
Definition AMReX_INT.H:32
amrex_particle_real ParticleReal
Floating Point Type for Particles.
Definition AMReX_REAL.H:90
amrex_long Long
Definition AMReX_INT.H:30
void WritePlotFile(const std::string &dir, const std::string &name) const
This version of WritePlotFile writes all components and assigns component names.
Definition AMReX_ParticleIO.H:113
ParticleContainer_impl< Particle< T_NStructReal, T_NStructInt >, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor > ParticleContainer
Definition AMReX_ParIter.H:29
Definition AMReX_Amr.cpp:49
amrex::ArenaAllocator< T > DefaultAllocator
Definition AMReX_GpuAllocators.H:205
Definition AMReX_ParticleUtil.H:390
Definition AMReX_GpuAllocators.H:184
A struct used for communicating particle data across processes during multi-level operations.
Definition AMReX_ParticleContainer.H:81
RealType m_data[1+3]
Definition AMReX_ParticleContainer.H:86
ParticleReal RealType
Definition AMReX_ParticleContainer.H:82
int m_grid
Definition AMReX_ParticleContainer.H:84
IntVect m_cell
Definition AMReX_ParticleContainer.H:85
int m_lev
Definition AMReX_ParticleContainer.H:83
Definition AMReX_ParticleCommunication.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:119
std::array< int, NStructInt > int_struct_data
Definition AMReX_ParticleContainer.H:121
std::array< int, NArrayInt > int_array_data
Definition AMReX_ParticleContainer.H:123
std::array< double, NArrayReal > real_array_data
Definition AMReX_ParticleContainer.H:122
std::array< double, NStructReal > real_struct_data
Definition AMReX_ParticleContainer.H:120
A struct used for storing a particle's position in the AMR hierarchy.
Definition AMReX_ParticleContainer.H:93
Box m_grown_gridbox
Definition AMReX_ParticleContainer.H:100
IntVect m_cell
Definition AMReX_ParticleContainer.H:97
int m_grid
Definition AMReX_ParticleContainer.H:95
int m_tile
Definition AMReX_ParticleContainer.H:96
int m_lev
Definition AMReX_ParticleContainer.H:94
Box m_tilebox
Definition AMReX_ParticleContainer.H:99
Box m_gridbox
Definition AMReX_ParticleContainer.H:98
Definition AMReX_ParticleTileRT.H:382
Definition AMReX_ParticleTile.H:723
The struct used to store particles.
Definition AMReX_Particle.H:405
ParticleReal RealType
The floating point type used for the particles.
Definition AMReX_Particle.H:418