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>
26#include <AMReX_TypeTraits.H>
27#include <AMReX_GpuContainers.H>
28#include <AMReX_ParticleUtil.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 <stdexcept>
64#include <string>
65#include <tuple>
66#include <type_traits>
67#include <utility>
68#include <vector>
69
70namespace amrex {
71
72#ifdef AMREX_USE_HDF5_ASYNC
73 extern hid_t es_par_g;
74#endif
80{
82 int m_lev;
83 int m_grid;
85 RealType m_data[1 + AMREX_SPACEDIM];
86};
87
101
116template<int NStructReal, int NStructInt, int NArrayReal, int NArrayInt>
118{
119 std::array<double, NStructReal> real_struct_data;
120 std::array<int, NStructInt > int_struct_data;
121 std::array<double, NArrayReal > real_array_data;
122 std::array<int, NArrayInt > int_array_data;
123};
124
125template <bool is_const, typename T_ParticleType, int NArrayReal, int NArrayInt,
126 template<class> class Allocator, class CellAssignor>
127class ParIterBase_impl;
128
143template <typename T_ParticleType, int T_NArrayReal, int T_NArrayInt,
144 template<class> class Allocator=DefaultAllocator,
145 class T_CellAssignor=DefaultAssignor>
146
148{
149public:
150 using ParticleType = T_ParticleType;
151 using ConstParticleType = typename ParticleType::ConstType;
152 using CellAssignor = T_CellAssignor;
153
155 static constexpr int NStructReal = ParticleType::NReal;
157 static constexpr int NStructInt = ParticleType::NInt;
159 static constexpr int NArrayReal = T_NArrayReal;
161 static constexpr int NArrayInt = T_NArrayInt;
163
164private:
165 friend class ParIterBase_impl<true,ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>;
166 friend class ParIterBase_impl<false,ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>;
167
168public:
170 template <typename T>
171 using AllocatorType = Allocator<T>;
176
177#ifdef AMREX_SINGLE_PRECISION_PARTICLES
179#else
181#endif
182
186
189 using ParticleLevel = std::map<std::pair<int, int>, ParticleTileType>;
192 using AoS = typename ParticleTileType::AoS;
193 using SoA = typename ParticleTileType::SoA;
194
195 using RealVector = typename SoA::RealVector;
196 using IntVector = typename SoA::IntVector;
197 using ParticleVector = typename AoS::ParticleVector;
201
202 static constexpr bool has_polymorphic_allocator =
204
208 :
210 {
211 Initialize ();
212 }
213
229
238 const DistributionMapping & dmap,
239 const BoxArray & ba)
240 :
241 ParticleContainerBase(geom, dmap, ba)
242 {
243 Initialize ();
246 }
247
258 const Vector<DistributionMapping> & dmap,
259 const Vector<BoxArray> & ba,
260 const Vector<int> & rr)
261 :
262 ParticleContainerBase(geom, dmap, ba, rr)
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 {
284 Initialize ();
287 }
288
289 ~ParticleContainer_impl () override = default;
290
293
295 ParticleContainer_impl& operator= ( ParticleContainer_impl && ) noexcept = default;
296
297
305 void Define (ParGDBBase* gdb)
306 {
308 reserveData();
309 resizeData();
310 }
311
319 void Define (const Geometry & geom,
320 const DistributionMapping & dmap,
321 const BoxArray & ba)
322 {
323 this->ParticleContainerBase::Define(geom, dmap, ba);
324 reserveData();
325 resizeData();
326 }
327
337 void Define (const Vector<Geometry> & geom,
338 const Vector<DistributionMapping> & dmap,
339 const Vector<BoxArray> & ba,
340 const Vector<int> & rr)
341 {
342 this->ParticleContainerBase::Define(geom, dmap, ba, rr);
343 reserveData();
344 resizeData();
345 }
346
356 void Define (const Vector<Geometry> & geom,
357 const Vector<DistributionMapping> & dmap,
358 const Vector<BoxArray> & ba,
359 const Vector<IntVect> & rr)
360 {
361 this->ParticleContainerBase::Define(geom, dmap, ba, rr);
362 reserveData();
363 resizeData();
364 }
365
367 int numLocalTilesAtLevel (int lev) const {
368 return (lev < m_particles.size()) ? m_particles[lev].size() : 0;
369 }
370
375 void reserveData () override;
376
383 void resizeData () override;
384
385 void InitFromAsciiFile (const std::string& file, int extradata,
386 const IntVect* Nrep = nullptr);
387
388 void InitFromBinaryFile (const std::string& file, int extradata);
389
390 void InitFromBinaryMetaFile (const std::string& file, int extradata);
391
408 void InitRandom (Long icount, ULong iseed,
409 const ParticleInitData& pdata,
410 bool serialize = false, RealBox bx = RealBox());
411
412
428 void InitRandomPerBox (Long icount, ULong iseed, const ParticleInitData& pdata);
429
430
445 void InitOnePerCell (Real x_off, Real y_off, Real z_off,
446 const ParticleInitData& pdata);
447
448
460 void InitNRandomPerCell (int n_per_cell, const ParticleInitData& pdata);
461
462 void Increment (MultiFab& mf, int level);
463
464 Long IncrementWithTotal (MultiFab& mf, int level, bool local = false);
465
503 void Redistribute (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
504 bool remove_negative=true);
505
506
518 template <class index_type>
519 void ReorderParticles (int lev, const MFIter& mfi, const index_type* permutations);
520
535 void SortParticlesForDeposition (IntVect idx_type);
536
540 void SortParticlesByCell ();
541
548 void SortParticlesByBin (IntVect bin_size);
549
560 bool OK (int lev_min = 0, int lev_max = -1, int nGrow = 0) const;
561
562 std::array<Long, 3> ByteSpread () const;
563
564 std::array<Long, 3> PrintCapacity () const;
565
566 void ShrinkToFit ();
567
578 Long NumberOfParticlesAtLevel (int level, bool only_valid = true, bool only_local = false) const;
579
580 Vector<Long> NumberOfParticlesInGrid (int level, bool only_valid = true, bool only_local = false) const;
581
585 template <typename I, std::enable_if_t<std::is_integral_v<I> && (sizeof(I) >= sizeof(Long)), int> = 0>
586 void CapacityOfParticlesInGrid (LayoutData<I>& mem, int lev) const;
587
596 Long TotalNumberOfParticles (bool only_valid=true, bool only_local=false) const;
597
598
606 void RemoveParticlesAtLevel (int level);
607
609
617 void CreateVirtualParticles (int level, AoS& virts) const;
618
627 void CreateGhostParticles (int level, int ngrow, AoS& ghosts) const;
628
636 void AddParticlesAtLevel (AoS& particles, int level, int nGrow=0);
637
645 void CreateVirtualParticles (int level, ParticleTileType& virts) const;
646
655 void CreateGhostParticles (int level, int ngrow, ParticleTileType& ghosts) const;
656
664 void AddParticlesAtLevel (ParticleTileType& particles, int level, int nGrow=0);
665
666
670 void clearParticles ();
671
672
681 template <class PCType,
682 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0>
683 void copyParticles (const PCType& other, bool local=false);
684
692 template <class PCType,
693 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0>
694 void addParticles (const PCType& other, bool local=false);
695
710 template <class F, class PCType,
711 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0,
712 std::enable_if_t<! std::is_integral_v<F>, int> bar = 0>
713 void copyParticles (const PCType& other, F&&f, bool local=false);
714
728 template <class F, class PCType,
729 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0,
730 std::enable_if_t<! std::is_integral_v<F>, int> bar = 0>
731 void addParticles (const PCType& other, F const& f, bool local=false);
732
740 void WriteParticleRealData (void* data, size_t size, std::ostream& os) const;
741
749 void ReadParticleRealData (void* data, size_t size, std::istream& is);
750
759 void Checkpoint (const std::string& dir, const std::string& name,
760 const Vector<std::string>& real_comp_names = Vector<std::string>(),
761 const Vector<std::string>& int_comp_names = Vector<std::string>()) const
762 {
763 Checkpoint(dir, name, true, real_comp_names, int_comp_names);
764 }
765
771 void Checkpoint (const std::string& dir, const std::string& name, bool is_checkpoint,
772 const Vector<std::string>& real_comp_names = Vector<std::string>(),
773 const Vector<std::string>& int_comp_names = Vector<std::string>()) const;
774
787 void Checkpoint (const std::string& dir, 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) const;
792
807 template <class F>
808 void WriteBinaryParticleData (const std::string& dir,
809 const std::string& name,
810 const Vector<int>& write_real_comp,
811 const Vector<int>& write_int_comp,
812 const Vector<std::string>& real_comp_names,
813 const Vector<std::string>& int_comp_names,
814 F&& f, bool is_checkpoint=false) const;
815
816 void CheckpointPre ();
817
818 void CheckpointPost ();
819
826 void Restart (const std::string& dir, const std::string& file);
827
835 void Restart (const std::string& dir, const std::string& file, bool is_checkpoint);
836
843 void WritePlotFile (const std::string& dir, const std::string& name) const;
844
856 template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>* = nullptr>
857 void WritePlotFile (const std::string& dir, const std::string& name, F&& f) const;
858
867 void WritePlotFile (const std::string& dir, const std::string& name,
868 const Vector<std::string>& real_comp_names,
869 const Vector<std::string>& int_comp_names) const;
870
884 template <class F>
885 void WritePlotFile (const std::string& dir, const std::string& name,
886 const Vector<std::string>& real_comp_names,
887 const Vector<std::string>& int_comp_names, F&& f) const;
888
897 void WritePlotFile (const std::string& dir, const std::string& name,
898 const Vector<std::string>& real_comp_names) const;
899
913 template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>* = nullptr>
914 void WritePlotFile (const std::string& dir, const std::string& name,
915 const Vector<std::string>& real_comp_names, F&& f) const;
916
926 void WritePlotFile (const std::string& dir,
927 const std::string& name,
928 const Vector<int>& write_real_comp,
929 const Vector<int>& write_int_comp) const;
930
945 template <class F>
946 void WritePlotFile (const std::string& dir,
947 const std::string& name,
948 const Vector<int>& write_real_comp,
949 const Vector<int>& write_int_comp, F&& f) const;
950
964 void WritePlotFile (const std::string& dir,
965 const std::string& name,
966 const Vector<int>& write_real_comp,
967 const Vector<int>& write_int_comp,
968 const Vector<std::string>& real_comp_names,
969 const Vector<std::string>& int_comp_names) const;
970
989 template <class F>
990 void WritePlotFile (const std::string& dir,
991 const std::string& name,
992 const Vector<int>& write_real_comp,
993 const Vector<int>& write_int_comp,
994 const Vector<std::string>& real_comp_names,
995 const Vector<std::string>& int_comp_names,
996 F&& f) const;
997
998 void WritePlotFilePre ();
999
1000 void WritePlotFilePost ();
1001
1002 void WriteAsciiFile (const std::string& file);
1003
1008 const Vector<ParticleLevel>& GetParticles () const { return m_particles; }
1009
1014 Vector <ParticleLevel>& GetParticles () { return m_particles; }
1015
1029 const ParticleLevel& GetParticles (int lev) const { return m_particles[lev]; }
1030
1044 ParticleLevel & GetParticles (int lev) { return m_particles[lev]; }
1045
1068 const ParticleTileType& ParticlesAt (int lev, int grid, int tile) const
1069 { return m_particles[lev].at(std::make_pair(grid, tile)); }
1070
1093 ParticleTileType& ParticlesAt (int lev, int grid, int tile)
1094 { return m_particles[lev].at(std::make_pair(grid, tile)); }
1095
1117 template <class Iterator>
1118 const ParticleTileType& ParticlesAt (int lev, const Iterator& iter) const
1119 { return ParticlesAt(lev, iter.index(), iter.LocalTileIndex()); }
1120
1142 template <class Iterator>
1143 ParticleTileType& ParticlesAt (int lev, const Iterator& iter)
1144 { return ParticlesAt(lev, iter.index(), iter.LocalTileIndex()); }
1145
1168 ParticleTileType& DefineAndReturnParticleTile (int lev, int grid, int tile)
1169 {
1170 m_particles[lev][std::make_pair(grid, tile)].define(
1172 &m_soa_rdata_names, &m_soa_idata_names,
1173 arena()
1174 );
1175 return ParticlesAt(lev, grid, tile);
1176 }
1177
1199 template <class Iterator>
1200 ParticleTileType& DefineAndReturnParticleTile (int lev, const Iterator& iter)
1201 {
1202 auto index = std::make_pair(iter.index(), iter.LocalTileIndex());
1203 m_particles[lev][index].define(
1205 &m_soa_rdata_names, &m_soa_idata_names,
1206 arena()
1207 );
1208 return ParticlesAt(lev, iter);
1209 }
1210
1221 void AssignDensity (int rho_index,
1222 Vector<std::unique_ptr<MultiFab> >& mf_to_be_filled,
1223 int lev_min, int ncomp, int finest_level, int ngrow=2) const;
1224
1225 void AssignCellDensitySingleLevel (int rho_index, MultiFab& mf, int level,
1226 int ncomp=1, int particle_lvl_offset = 0) const;
1227
1228 template <typename P, typename Assignor=CellAssignor>
1229 IntVect Index (const P& p, int lev) const;
1230
1240 ParticleLocData Reset (ParticleType& prt, bool update, bool verbose=true,
1241 ParticleLocData pld = ParticleLocData()) const;
1242
1248 template <typename P>
1249 bool PeriodicShift (P& p) const;
1250
1252
1254
1255 void SetUsePrePost (bool tf) const {
1256 usePrePost = tf;
1257 }
1258 bool GetUsePrePost () const {
1259 return usePrePost;
1260 }
1261
1262 int GetMaxNextIDPrePost () const { return maxnextidPrePost; }
1264
1265 void SetUseUnlink (bool tf) const {
1266 doUnlink = tf;
1267 }
1268
1269 bool GetUseUnlink () const {
1270 return doUnlink;
1271 }
1272
1273 void RedistributeCPU (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
1274 bool remove_negative=true);
1275
1276 void RedistributeGPU (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
1277 bool remove_negative=true);
1278
1279 void ReserveForRedistribute (ParticleCopyPlan const& plan);
1280
1281 Long superParticleSize() const { return superparticle_size; }
1282
1283 void AddRealComp (std::string const & name, int communicate=1)
1284 {
1285 // names must be unique
1286 auto const it = std::find(m_soa_rdata_names.begin(), m_soa_rdata_names.end(), name);
1287 if (it != m_soa_rdata_names.end()) {
1288 throw std::runtime_error("AddRealComp: name '" + name + "' is already present in the SoA.");
1289 }
1290 m_soa_rdata_names.push_back(name);
1291
1292 m_runtime_comps_defined = true;
1293 m_num_runtime_real++;
1294 h_redistribute_real_comp.push_back(communicate);
1296 this->resizeData();
1297
1298 // resize runtime SoA
1299 for (int lev = 0; lev < numLevels(); ++lev) {
1300 for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1301 auto& tile = DefineAndReturnParticleTile(lev, pti);
1302 auto np = tile.numParticles();
1303 if (np > 0) {
1304 auto& soa = tile.GetStructOfArrays();
1305 soa.resize(np);
1306 }
1307 }
1308 }
1309 }
1310
1311 void AddRealComp (int communicate=1)
1312 {
1313 AddRealComp(getDefaultCompNameReal<ParticleType>(NArrayReal+m_num_runtime_real), communicate);
1314 }
1315
1316 void AddIntComp (std::string const & name, int communicate=1)
1317 {
1318 // names must be unique
1319 auto const it = std::find(m_soa_idata_names.begin(), m_soa_idata_names.end(), name);
1320 if (it != m_soa_idata_names.end()) {
1321 throw std::runtime_error("AddIntComp: name '" + name + "' is already present in the SoA.");
1322 }
1323 m_soa_idata_names.push_back(name);
1324
1325 m_runtime_comps_defined = true;
1326 m_num_runtime_int++;
1327 h_redistribute_int_comp.push_back(communicate);
1329 this->resizeData();
1330
1331 // resize runtime SoA
1332 for (int lev = 0; lev < numLevels(); ++lev) {
1333 for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1334 auto& tile = DefineAndReturnParticleTile(lev, pti);
1335 auto np = tile.numParticles();
1336 if (np > 0) {
1337 auto& soa = tile.GetStructOfArrays();
1338 soa.resize(np);
1339 }
1340 }
1341 }
1342 }
1343
1344 void AddIntComp (int communicate=1)
1345 {
1346 AddIntComp(getDefaultCompNameInt<ParticleType>(NArrayInt+m_num_runtime_int), communicate);
1347 }
1348
1349 int NumRuntimeRealComps () const { return m_num_runtime_real; }
1350 int NumRuntimeIntComps () const { return m_num_runtime_int; }
1351
1352 int NumRealComps () const { return NArrayReal + NumRuntimeRealComps(); }
1353 int NumIntComps () const { return NArrayInt + NumRuntimeIntComps() ; }
1354
1360 void ResizeRuntimeRealComp (int new_size, bool communicate);
1361
1367 void ResizeRuntimeIntComp (int new_size, bool communicate);
1368
1370 template <template<class> class NewAllocator=amrex::DefaultAllocator>
1372
1383 template <template<class> class NewAllocator=Allocator>
1385 make_alike () const
1386 {
1388
1389 std::vector<std::string> const real_names = this->GetRealSoANames();
1390 std::vector<std::string> real_ct_names(NArrayReal);
1391 for (int ic = 0; ic < NArrayReal; ++ic) { real_ct_names.at(ic) = real_names[ic]; }
1392
1393 std::vector<std::string> const int_names = this->GetIntSoANames();
1394 std::vector<std::string> int_ct_names(NArrayInt);
1395 for (int ic = 0; ic < NArrayInt; ++ic) { int_ct_names.at(ic) = int_names[ic]; }
1396
1397 tmp.SetSoACompileTimeNames(real_ct_names, int_ct_names);
1398
1399 // add runtime real comps to tmp
1400 for (int ic = 0; ic < this->NumRuntimeRealComps(); ++ic) {
1401 tmp.AddRealComp(real_names.at(ic + NArrayReal));
1402 }
1403
1404 // add runtime int comps to tmp
1405 for (int ic = 0; ic < this->NumRuntimeIntComps(); ++ic) {
1406 tmp.AddIntComp(int_names.at(ic + NArrayInt));
1407 }
1408
1411
1412 return tmp;
1413 }
1414
1417
1420 mutable bool usePrePost;
1421 mutable bool doUnlink;
1423 mutable int nOutFilesPrePost;
1429 mutable std::string HdrFileNamePrePost;
1431
1432protected:
1433
1447 template <typename P>
1448 bool Where (const P& prt, ParticleLocData& pld,
1449 int lev_min = 0, int lev_max = -1, int nGrow=0, int local_grid=-1) const;
1450
1451
1462 template <typename P>
1463 bool EnforcePeriodicWhere (P& prt, ParticleLocData& pld,
1464 int lev_min = 0, int lev_max = -1, int local_grid=-1) const;
1465
1466public:
1467 void
1468 WriteParticles (int level, std::ofstream& ofs, int fnum,
1469 Vector<int>& which, Vector<int>& count, Vector<Long>& where,
1470 const Vector<int>& write_real_comp, const Vector<int>& write_int_comp,
1471 const Vector<std::map<std::pair<int, int>,IntVector>>& particle_io_flags, bool is_checkpoint) const;
1472#ifdef AMREX_USE_HDF5
1473#include "AMReX_ParticlesHDF5.H"
1474#endif
1475
1477 void SetSoACompileTimeNames (std::vector<std::string> const & rdata_name, std::vector<std::string> const & idata_name);
1478
1480 std::vector<std::string> GetRealSoANames () const {return m_soa_rdata_names;}
1481
1483 std::vector<std::string> GetIntSoANames () const {return m_soa_idata_names;}
1484
1490 bool HasRealComp (std::string const & name);
1491
1497 bool HasIntComp (std::string const & name);
1498
1506 int GetRealCompIndex (std::string const & name);
1507
1515 int GetIntCompIndex (std::string const & name);
1516
1517protected:
1518
1519 template <class RTYPE>
1520 void ReadParticles (int cnt, int grd, int lev, std::ifstream& ifs, int finest_level_in_file, bool convert_ids);
1521
1522 void SetParticleSize ();
1523
1525
1526private:
1527 virtual void particlePostLocate (ParticleType& /*p*/, const ParticleLocData& /*pld*/,
1528 const int /*lev*/) {}
1529
1530 virtual void correctCellVectors (int /*old_index*/, int /*new_index*/,
1531 int /*grid*/, const ParticleType& /*p*/) {}
1532
1533 void RedistributeMPI (std::map<int, Vector<char> >& not_ours,
1534 int lev_min = 0, int lev_max = 0, int nGrow = 0, int local=0);
1535
1536 template <typename P>
1537 void locateParticle (P& p, ParticleLocData& pld,
1538 int lev_min, int lev_max, int nGrow, int local_grid=-1) const;
1539
1540 void Initialize ();
1541
1542 bool m_runtime_comps_defined{false};
1543 int m_num_runtime_real{0};
1544 int m_num_runtime_int{0};
1545
1546 size_t particle_size, superparticle_size;
1547 int num_real_comm_comps, num_int_comm_comps;
1548 Vector<ParticleLevel> m_particles;
1549
1550 // names of both compile-time and runtime Real and Int SoA data
1551 std::vector<std::string> m_soa_rdata_names;
1552 std::vector<std::string> m_soa_idata_names;
1553};
1554
1556template <int T_NStructReal, int T_NStructInt, int T_NArrayReal, int T_NArrayInt, template<class> class Allocator, class CellAssignor>
1557using ParticleContainer = ParticleContainer_impl<Particle<T_NStructReal, T_NStructInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
1558
1559template <int T_NArrayReal, int T_NArrayInt, template<class> class Allocator=DefaultAllocator, class CellAssignor=DefaultAssignor>
1560using ParticleContainerPureSoA = ParticleContainer_impl<SoAParticle<T_NArrayReal, T_NArrayInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
1561
1562}
1563
1564#include "AMReX_ParticleInit.H"
1566#include "AMReX_ParticleIO.H"
1567
1568#ifdef AMREX_USE_HDF5
1569#include "AMReX_ParticleHDF5.H"
1570#endif
1571
1572#endif /*_PARTICLES_H_*/
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Definition AMReX_GpuAllocators.H:115
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:567
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:144
Definition AMReX_ParGDB.H:13
Definition AMReX_ParIter.H:32
Definition AMReX_ParIter.H:115
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:148
ParticleTileType & DefineAndReturnParticleTile(int lev, const Iterator &iter)
Define and return the ParticleTile for level "lev", and Iterator "iter".
Definition AMReX_ParticleContainer.H:1200
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:17
typename ParticleTileType::ConstParticleTileDataType ConstPTDType
Definition AMReX_ParticleContainer.H:191
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:1265
Vector< Long > nParticlesAtLevelPrePost
Definition AMReX_ParticleContainer.H:1425
bool GetLevelDirectoriesCreated() const
Definition AMReX_ParticleContainer.H:1253
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:1236
void AddIntComp(std::string const &name, int communicate=1)
Definition AMReX_ParticleContainer.H:1316
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:1261
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:2785
bool doUnlink
Definition AMReX_ParticleContainer.H:1421
Vector< Vector< int > > whichPrePost
Definition AMReX_ParticleContainer.H:1426
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:1429
Vector< Vector< Long > > wherePrePost
Definition AMReX_ParticleContainer.H:1428
RealDescriptor ParticleRealDescriptor
Definition AMReX_ParticleContainer.H:180
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:2558
ParticleTileType & ParticlesAt(int lev, const Iterator &iter)
Return the ParticleTile for level "lev" and Iterator "iter". Non-const version.
Definition AMReX_ParticleContainer.H:1143
int nOutFilesPrePost
Definition AMReX_ParticleContainer.H:1423
int maxnextidPrePost
Definition AMReX_ParticleContainer.H:1422
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:759
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:189
bool levelDirectoriesCreated
Variables for i/o optimization saved for pre and post checkpoint.
Definition AMReX_ParticleContainer.H:1419
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:1671
static constexpr int NArrayInt
Number of extra integer components stored in struct-of-array form.
Definition AMReX_ParticleContainer.H:161
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:159
int GetMaxNextIDPrePost() const
Definition AMReX_ParticleContainer.H:1262
const ParticleTileType & ParticlesAt(int lev, const Iterator &iter) const
Return the ParticleTile for level "lev" and Iterator "iter". Const version.
Definition AMReX_ParticleContainer.H:1118
void WritePlotFilePost()
Definition AMReX_ParticleIO.H:575
bool GetUsePrePost() const
Definition AMReX_ParticleContainer.H:1258
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:192
ParticleInitType< NStructReal, NStructInt, NArrayReal, NArrayInt > ParticleInitData
Definition AMReX_ParticleContainer.H:185
typename SoA::RealVector RealVector
Definition AMReX_ParticleContainer.H:195
typename ParticleType::ConstType ConstParticleType
Definition AMReX_ParticleContainer.H:151
ParticleContainer_impl(const ParticleContainer_impl &)=delete
void WriteAsciiFile(const std::string &file)
Definition AMReX_ParticleIO.H:1137
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:237
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:1251
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:1068
typename AoS::ParticleVector ParticleVector
Definition AMReX_ParticleContainer.H:197
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:319
static constexpr int NStructInt
Number of extra integer components in the particle struct.
Definition AMReX_ParticleContainer.H:157
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:1405
ParticleContainer_impl()
Default constructor - construct an empty particle container that has no concept of a level hierarchy....
Definition AMReX_ParticleContainer.H:207
void SetUsePrePost(bool tf) const
Definition AMReX_ParticleContainer.H:1255
int NumRealComps() const
Definition AMReX_ParticleContainer.H:1352
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:956
void Increment(MultiFab &mf, int level)
Definition AMReX_ParticleContainerI.H:720
Vector< int > h_redistribute_int_comp
Definition AMReX_ParticleContainer.H:1416
typename ParticleTileType::SoA SoA
Definition AMReX_ParticleContainer.H:193
void AssignCellDensitySingleLevel(int rho_index, MultiFab &mf, int level, int ncomp=1, int particle_lvl_offset=0) const
Definition AMReX_ParticleContainerI.H:2619
int NumRuntimeIntComps() const
Definition AMReX_ParticleContainer.H:1350
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:1263
void SortParticlesByBin(IntVect bin_size)
Sort the particles on each tile by groups of cells, given an IntVect bin_size.
Definition AMReX_ParticleContainerI.H:1372
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:152
const Vector< ParticleLevel > & GetParticles() const
Return the underlying Vector (over AMR levels) of ParticleLevels. Const version.
Definition AMReX_ParticleContainer.H:1008
const ParticleLevel & GetParticles(int lev) const
Return the ParticleLevel for level "lev". Const version.
Definition AMReX_ParticleContainer.H:1029
std::vector< std::string > GetRealSoANames() const
Definition AMReX_ParticleContainer.H:1480
DenseBins< typename ParticleTileType::ParticleTileDataType > m_bins
Definition AMReX_ParticleContainer.H:1524
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:221
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:356
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:1349
void ResizeRuntimeRealComp(int new_size, bool communicate)
Definition AMReX_ParticleContainerI.H:2759
void AddIntComp(int communicate=1)
Definition AMReX_ParticleContainer.H:1344
Long nparticlesPrePost
Definition AMReX_ParticleContainer.H:1424
Vector< std::string > filePrefixPrePost
Definition AMReX_ParticleContainer.H:1430
int numLocalTilesAtLevel(int lev) const
The total number of tiles on this rank on this level.
Definition AMReX_ParticleContainer.H:367
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:155
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:1363
Vector< int > h_redistribute_real_comp
Definition AMReX_ParticleContainer.H:1415
typename Particle< NStructReal, NStructInt >::RealType RealType
The type of the Real data.
Definition AMReX_ParticleContainer.H:175
ParticleTile< ParticleType, NArrayReal, NArrayInt, Allocator > ParticleTileType
Definition AMReX_ParticleContainer.H:184
bool GetUseUnlink() const
Definition AMReX_ParticleContainer.H:1269
int GetRealCompIndex(std::string const &name)
Definition AMReX_ParticleContainerI.H:161
void InitFromBinaryMetaFile(const std::string &file, int extradata)
Definition AMReX_ParticleInit.H:931
int NumIntComps() const
Definition AMReX_ParticleContainer.H:1353
void ReserveForRedistribute(ParticleCopyPlan const &plan)
Definition AMReX_ParticleContainerI.H:1631
typename ParticleTileType::ParticleTileDataType PTDType
Definition AMReX_ParticleContainer.H:190
void RedistributeGPU(int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
Definition AMReX_ParticleContainerI.H:1435
ParticleContainer_impl & operator=(const ParticleContainer_impl &)=delete
Allocator< T > AllocatorType
The memory allocator in use.
Definition AMReX_ParticleContainer.H:171
ParticleLevel & GetParticles(int lev)
Return the ParticleLevel for level "lev". Non-const version.
Definition AMReX_ParticleContainer.H:1044
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:1168
std::vector< std::string > GetIntSoANames() const
Definition AMReX_ParticleContainer.H:1483
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:305
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:257
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:1281
typename SoA::IntVector IntVector
Definition AMReX_ParticleContainer.H:196
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
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:585
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:337
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:202
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:1093
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:1283
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:1385
void AddRealComp(int communicate=1)
Definition AMReX_ParticleContainer.H:1311
void AddParticlesAtLevel(AoS &particles, int level, int nGrow=0)
Add particles from a pbox to the grid at this level.
Definition AMReX_ParticleContainerI.H:2573
~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:1014
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:150
Vector< Vector< int > > countPrePost
Definition AMReX_ParticleContainer.H:1427
bool usePrePost
Definition AMReX_ParticleContainer.H:1420
A Box with real dimensions.
Definition AMReX_RealBox.H:26
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:26
Definition AMReX_Amr.cpp:49
amrex::ArenaAllocator< T > DefaultAllocator
Definition AMReX_GpuAllocators.H:199
Definition AMReX_ParticleUtil.H:394
Definition AMReX_GpuAllocators.H:178
A struct used for communicating particle data across processes during multi-level operations.
Definition AMReX_ParticleContainer.H:80
RealType m_data[1+3]
Definition AMReX_ParticleContainer.H:85
ParticleReal RealType
Definition AMReX_ParticleContainer.H:81
int m_grid
Definition AMReX_ParticleContainer.H:83
IntVect m_cell
Definition AMReX_ParticleContainer.H:84
int m_lev
Definition AMReX_ParticleContainer.H:82
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:118
std::array< int, NStructInt > int_struct_data
Definition AMReX_ParticleContainer.H:120
std::array< int, NArrayInt > int_array_data
Definition AMReX_ParticleContainer.H:122
std::array< double, NArrayReal > real_array_data
Definition AMReX_ParticleContainer.H:121
std::array< double, NStructReal > real_struct_data
Definition AMReX_ParticleContainer.H:119
A struct used for storing a particle's position in the AMR hierarchy.
Definition AMReX_ParticleContainer.H:92
Box m_grown_gridbox
Definition AMReX_ParticleContainer.H:99
IntVect m_cell
Definition AMReX_ParticleContainer.H:96
int m_grid
Definition AMReX_ParticleContainer.H:94
int m_tile
Definition AMReX_ParticleContainer.H:95
int m_lev
Definition AMReX_ParticleContainer.H:93
Box m_tilebox
Definition AMReX_ParticleContainer.H:98
Box m_gridbox
Definition AMReX_ParticleContainer.H:97
Definition AMReX_ParticleTile.H:721
ParticleTileData< StorageParticleType, NArrayReal, NArrayInt > ParticleTileDataType
Definition AMReX_ParticleTile.H:749
std::conditional_t< ParticleType::is_soa_particle, StructOfArrays< NArrayReal, NArrayInt, Allocator, true >, StructOfArrays< NArrayReal, NArrayInt, Allocator, false > > SoA
Definition AMReX_ParticleTile.H:744
ConstParticleTileData< StorageParticleType, NArrayReal, NArrayInt > ConstParticleTileDataType
Definition AMReX_ParticleTile.H:750
std::conditional_t< ParticleType::is_soa_particle, ThisParticleTileHasNoAoS, ArrayOfStructs< ParticleType, Allocator > > AoS
Definition AMReX_ParticleTile.H:738
The struct used to store particles.
Definition AMReX_Particle.H:404
ParticleReal RealType
The floating point type used for the particles.
Definition AMReX_Particle.H:416