Block-Structured AMR Software Framework
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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>
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{
81 using RealType = ParticleReal;
82 int m_lev;
83 int m_grid;
85 RealType m_data[1 + AMREX_SPACEDIM];
86};
87
101
115template<int NStructReal, int NStructInt, int NArrayReal, int NArrayInt>
117{
118 std::array<double, NStructReal> real_struct_data;
119 std::array<int, NStructInt > int_struct_data;
120 std::array<double, NArrayReal > real_array_data;
121 std::array<int, NArrayInt > int_array_data;
122};
123
124template <bool is_const, typename T_ParticleType, int NArrayReal, int NArrayInt,
125 template<class> class Allocator, class CellAssignor>
126class ParIterBase_impl;
127
141template <typename T_ParticleType, int T_NArrayReal, int T_NArrayInt,
142 template<class> class Allocator=DefaultAllocator,
143 class T_CellAssignor=DefaultAssignor>
144
146{
147public:
148 using ParticleType = T_ParticleType;
149 using ConstParticleType = typename ParticleType::ConstType;
150 using CellAssignor = T_CellAssignor;
151
153 static constexpr int NStructReal = ParticleType::NReal;
155 static constexpr int NStructInt = ParticleType::NInt;
157 static constexpr int NArrayReal = T_NArrayReal;
159 static constexpr int NArrayInt = T_NArrayInt;
161
162private:
163 friend class ParIterBase_impl<true,ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>;
164 friend class ParIterBase_impl<false,ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>;
165
166public:
168 template <typename T>
169 using AllocatorType = Allocator<T>;
174
175#ifdef AMREX_SINGLE_PRECISION_PARTICLES
177#else
179#endif
180
184
187 using ParticleLevel = std::map<std::pair<int, int>, ParticleTileType>;
188 using AoS = typename ParticleTileType::AoS;
189 using SoA = typename ParticleTileType::SoA;
190
191 using RealVector = typename SoA::RealVector;
192 using IntVector = typename SoA::IntVector;
193 using ParticleVector = typename AoS::ParticleVector;
197
201 :
203 {
204 Initialize ();
205 }
206
222
231 const DistributionMapping & dmap,
232 const BoxArray & ba)
233 :
234 ParticleContainerBase(geom, dmap, ba)
235 {
236 Initialize ();
239 }
240
251 const Vector<DistributionMapping> & dmap,
252 const Vector<BoxArray> & ba,
253 const Vector<int> & rr)
254 :
255 ParticleContainerBase(geom, dmap, ba, rr)
256 {
257 Initialize ();
260 }
261
271 const Vector<DistributionMapping> & dmap,
272 const Vector<BoxArray> & ba,
273 const Vector<IntVect> & rr)
274 :
275 ParticleContainerBase(geom, dmap, ba, rr)
276 {
277 Initialize ();
280 }
281
282 ~ParticleContainer_impl () override = default;
283
286
288 ParticleContainer_impl& operator= ( ParticleContainer_impl && ) noexcept = default;
289
290
298 void Define (ParGDBBase* gdb)
299 {
301 reserveData();
302 resizeData();
303 }
304
312 void Define (const Geometry & geom,
313 const DistributionMapping & dmap,
314 const BoxArray & ba)
315 {
316 this->ParticleContainerBase::Define(geom, dmap, ba);
317 reserveData();
318 resizeData();
319 }
320
330 void Define (const Vector<Geometry> & geom,
331 const Vector<DistributionMapping> & dmap,
332 const Vector<BoxArray> & ba,
333 const Vector<int> & rr)
334 {
335 this->ParticleContainerBase::Define(geom, dmap, ba, rr);
336 reserveData();
337 resizeData();
338 }
339
349 void Define (const Vector<Geometry> & geom,
350 const Vector<DistributionMapping> & dmap,
351 const Vector<BoxArray> & ba,
352 const Vector<IntVect> & rr)
353 {
354 this->ParticleContainerBase::Define(geom, dmap, ba, rr);
355 reserveData();
356 resizeData();
357 }
358
360 int numLocalTilesAtLevel (int lev) const {
361 return (lev < m_particles.size()) ? m_particles[lev].size() : 0;
362 }
363
368 void reserveData () override;
369
376 void resizeData () override;
377
378 void InitFromAsciiFile (const std::string& file, int extradata,
379 const IntVect* Nrep = nullptr);
380
381 void InitFromBinaryFile (const std::string& file, int extradata);
382
383 void InitFromBinaryMetaFile (const std::string& file, int extradata);
384
401 void InitRandom (Long icount, ULong iseed,
402 const ParticleInitData& pdata,
403 bool serialize = false, RealBox bx = RealBox());
404
405
421 void InitRandomPerBox (Long icount, ULong iseed, const ParticleInitData& pdata);
422
423
438 void InitOnePerCell (Real x_off, Real y_off, Real z_off,
439 const ParticleInitData& pdata);
440
441
453 void InitNRandomPerCell (int n_per_cell, const ParticleInitData& pdata);
454
455 void Increment (MultiFab& mf, int level);
456
457 Long IncrementWithTotal (MultiFab& mf, int level, bool local = false);
458
495 void Redistribute (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
496 bool remove_negative=true);
497
498
510 template <class index_type>
511 void ReorderParticles (int lev, const MFIter& mfi, const index_type* permutations);
512
527 void SortParticlesForDeposition (IntVect idx_type);
528
532 void SortParticlesByCell ();
533
540 void SortParticlesByBin (IntVect bin_size);
541
552 bool OK (int lev_min = 0, int lev_max = -1, int nGrow = 0) const;
553
554 std::array<Long, 3> ByteSpread () const;
555
556 std::array<Long, 3> PrintCapacity () const;
557
558 void ShrinkToFit ();
559
570 Long NumberOfParticlesAtLevel (int level, bool only_valid = true, bool only_local = false) const;
571
572 Vector<Long> NumberOfParticlesInGrid (int level, bool only_valid = true, bool only_local = false) const;
573
582 Long TotalNumberOfParticles (bool only_valid=true, bool only_local=false) const;
583
584
592 void RemoveParticlesAtLevel (int level);
593
595
603 void CreateVirtualParticles (int level, AoS& virts) const;
604
613 void CreateGhostParticles (int level, int ngrow, AoS& ghosts) const;
614
622 void AddParticlesAtLevel (AoS& particles, int level, int nGrow=0);
623
631 void CreateVirtualParticles (int level, ParticleTileType& virts) const;
632
641 void CreateGhostParticles (int level, int ngrow, ParticleTileType& ghosts) const;
642
650 void AddParticlesAtLevel (ParticleTileType& particles, int level, int nGrow=0);
651
652
656 void clearParticles ();
657
658
667 template <class PCType,
668 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0>
669 void copyParticles (const PCType& other, bool local=false);
670
678 template <class PCType,
679 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0>
680 void addParticles (const PCType& other, bool local=false);
681
696 template <class F, class PCType,
697 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0,
698 std::enable_if_t<! std::is_integral_v<F>, int> bar = 0>
699 void copyParticles (const PCType& other, F&&f, bool local=false);
700
714 template <class F, class PCType,
715 std::enable_if_t<IsParticleContainer<PCType>::value, int> foo = 0,
716 std::enable_if_t<! std::is_integral_v<F>, int> bar = 0>
717 void addParticles (const PCType& other, F const& f, bool local=false);
718
726 void WriteParticleRealData (void* data, size_t size, std::ostream& os) const;
727
736 void ReadParticleRealData (void* data, size_t size, std::istream& is);
737
746 void Checkpoint (const std::string& dir, const std::string& name,
747 const Vector<std::string>& real_comp_names = Vector<std::string>(),
748 const Vector<std::string>& int_comp_names = Vector<std::string>()) const
749 {
750 Checkpoint(dir, name, true, real_comp_names, int_comp_names);
751 }
752
758 void Checkpoint (const std::string& dir, const std::string& name, bool is_checkpoint,
759 const Vector<std::string>& real_comp_names = Vector<std::string>(),
760 const Vector<std::string>& int_comp_names = Vector<std::string>()) const;
761
774 void Checkpoint (const std::string& dir, const std::string& name,
775 const Vector<int>& write_real_comp,
776 const Vector<int>& write_int_comp,
777 const Vector<std::string>& real_comp_names,
778 const Vector<std::string>& int_comp_names) const;
779
794 template <class F>
795 void WriteBinaryParticleData (const std::string& dir,
796 const std::string& name,
797 const Vector<int>& write_real_comp,
798 const Vector<int>& write_int_comp,
799 const Vector<std::string>& real_comp_names,
800 const Vector<std::string>& int_comp_names,
801 F&& f, bool is_checkpoint=false) const;
802
803 void CheckpointPre ();
804
805 void CheckpointPost ();
806
813 void Restart (const std::string& dir, const std::string& file);
814
822 void Restart (const std::string& dir, const std::string& file, bool is_checkpoint);
823
830 void WritePlotFile (const std::string& dir, const std::string& name) const;
831
843 template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>* = nullptr>
844 void WritePlotFile (const std::string& dir, const std::string& name, F&& f) const;
845
854 void WritePlotFile (const std::string& dir, const std::string& name,
855 const Vector<std::string>& real_comp_names,
856 const Vector<std::string>& int_comp_names) const;
857
871 template <class F>
872 void WritePlotFile (const std::string& dir, const std::string& name,
873 const Vector<std::string>& real_comp_names,
874 const Vector<std::string>& int_comp_names, F&& f) const;
875
884 void WritePlotFile (const std::string& dir, const std::string& name,
885 const Vector<std::string>& real_comp_names) const;
886
900 template <class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::string>>>* = nullptr>
901 void WritePlotFile (const std::string& dir, const std::string& name,
902 const Vector<std::string>& real_comp_names, F&& f) const;
903
913 void WritePlotFile (const std::string& dir,
914 const std::string& name,
915 const Vector<int>& write_real_comp,
916 const Vector<int>& write_int_comp) const;
917
932 template <class F>
933 void WritePlotFile (const std::string& dir,
934 const std::string& name,
935 const Vector<int>& write_real_comp,
936 const Vector<int>& write_int_comp, F&& f) const;
937
951 void WritePlotFile (const std::string& dir,
952 const std::string& name,
953 const Vector<int>& write_real_comp,
954 const Vector<int>& write_int_comp,
955 const Vector<std::string>& real_comp_names,
956 const Vector<std::string>& int_comp_names) const;
957
976 template <class F>
977 void WritePlotFile (const std::string& dir,
978 const std::string& name,
979 const Vector<int>& write_real_comp,
980 const Vector<int>& write_int_comp,
981 const Vector<std::string>& real_comp_names,
982 const Vector<std::string>& int_comp_names,
983 F&& f) const;
984
985 void WritePlotFilePre ();
986
987 void WritePlotFilePost ();
988
989 void WriteAsciiFile (const std::string& file);
990
996
1001 Vector <ParticleLevel>& GetParticles () { return m_particles; }
1002
1016 const ParticleLevel& GetParticles (int lev) const { return m_particles[lev]; }
1017
1031 ParticleLevel & GetParticles (int lev) { return m_particles[lev]; }
1032
1055 const ParticleTileType& ParticlesAt (int lev, int grid, int tile) const
1056 { return m_particles[lev].at(std::make_pair(grid, tile)); }
1057
1080 ParticleTileType& ParticlesAt (int lev, int grid, int tile)
1081 { return m_particles[lev].at(std::make_pair(grid, tile)); }
1082
1104 template <class Iterator>
1105 const ParticleTileType& ParticlesAt (int lev, const Iterator& iter) const
1106 { return ParticlesAt(lev, iter.index(), iter.LocalTileIndex()); }
1107
1129 template <class Iterator>
1130 ParticleTileType& ParticlesAt (int lev, const Iterator& iter)
1131 { return ParticlesAt(lev, iter.index(), iter.LocalTileIndex()); }
1132
1155 ParticleTileType& DefineAndReturnParticleTile (int lev, int grid, int tile)
1156 {
1157 m_particles[lev][std::make_pair(grid, tile)].define(NumRuntimeRealComps(), NumRuntimeIntComps(), &m_soa_rdata_names, &m_soa_idata_names);
1158
1159 return ParticlesAt(lev, grid, tile);
1160 }
1161
1183 template <class Iterator>
1184 ParticleTileType& DefineAndReturnParticleTile (int lev, const Iterator& iter)
1185 {
1186 auto index = std::make_pair(iter.index(), iter.LocalTileIndex());
1188 return ParticlesAt(lev, iter);
1189 }
1190
1201 void AssignDensity (int rho_index,
1202 Vector<std::unique_ptr<MultiFab> >& mf_to_be_filled,
1203 int lev_min, int ncomp, int finest_level, int ngrow=2) const;
1204
1205 void AssignCellDensitySingleLevel (int rho_index, MultiFab& mf, int level,
1206 int ncomp=1, int particle_lvl_offset = 0) const;
1207
1208 template <typename P, typename Assignor=CellAssignor>
1209 IntVect Index (const P& p, int lev) const;
1210
1220 ParticleLocData Reset (ParticleType& prt, bool update, bool verbose=true,
1221 ParticleLocData pld = ParticleLocData()) const;
1222
1228 template <typename P>
1229 bool PeriodicShift (P& p) const;
1230
1232
1234
1235 void SetUsePrePost (bool tf) const {
1236 usePrePost = tf;
1237 }
1238 bool GetUsePrePost () const {
1239 return usePrePost;
1240 }
1241
1242 int GetMaxNextIDPrePost () const { return maxnextidPrePost; }
1243 Long GetNParticlesPrePost () const { return nparticlesPrePost; }
1244
1245 void SetUseUnlink (bool tf) const {
1246 doUnlink = tf;
1247 }
1248
1249 bool GetUseUnlink () const {
1250 return doUnlink;
1251 }
1252
1253 void RedistributeCPU (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
1254 bool remove_negative=true);
1255
1256 void RedistributeGPU (int lev_min = 0, int lev_max = -1, int nGrow = 0, int local=0,
1257 bool remove_negative=true);
1258
1259 Long superParticleSize() const { return superparticle_size; }
1260
1261 void AddRealComp (std::string const & name, int communicate=1)
1262 {
1263 // names must be unique
1264 auto const it = std::find(m_soa_rdata_names.begin(), m_soa_rdata_names.end(), name);
1265 if (it != m_soa_rdata_names.end()) {
1266 throw std::runtime_error("AddRealComp: name '" + name + "' is already present in the SoA.");
1267 }
1268 m_soa_rdata_names.push_back(name);
1269
1272 h_redistribute_real_comp.push_back(communicate);
1274 this->resizeData();
1275
1276 // resize runtime SoA
1277 for (int lev = 0; lev < numLevels(); ++lev) {
1278 for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1279 auto& tile = DefineAndReturnParticleTile(lev, pti);
1280 auto np = tile.numParticles();
1281 if (np > 0) {
1282 auto& soa = tile.GetStructOfArrays();
1283 soa.resize(np);
1284 }
1285 }
1286 }
1287 }
1288
1289 void AddRealComp (int communicate=1)
1290 {
1291 AddRealComp(getDefaultCompNameReal<ParticleType>(NArrayReal+m_num_runtime_real), communicate);
1292 }
1293
1294 void AddIntComp (std::string const & name, int communicate=1)
1295 {
1296 // names must be unique
1297 auto const it = std::find(m_soa_idata_names.begin(), m_soa_idata_names.end(), name);
1298 if (it != m_soa_idata_names.end()) {
1299 throw std::runtime_error("AddIntComp: name '" + name + "' is already present in the SoA.");
1300 }
1301 m_soa_idata_names.push_back(name);
1302
1305 h_redistribute_int_comp.push_back(communicate);
1307 this->resizeData();
1308
1309 // resize runtime SoA
1310 for (int lev = 0; lev < numLevels(); ++lev) {
1311 for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1312 auto& tile = DefineAndReturnParticleTile(lev, pti);
1313 auto np = tile.numParticles();
1314 if (np > 0) {
1315 auto& soa = tile.GetStructOfArrays();
1316 soa.resize(np);
1317 }
1318 }
1319 }
1320 }
1321
1322 void AddIntComp (int communicate=1)
1323 {
1324 AddIntComp(getDefaultCompNameInt<ParticleType>(NArrayInt+m_num_runtime_int), communicate);
1325 }
1326
1328 int NumRuntimeIntComps () const { return m_num_runtime_int; }
1329
1330 int NumRealComps () const { return NArrayReal + NumRuntimeRealComps(); }
1331 int NumIntComps () const { return NArrayInt + NumRuntimeIntComps() ; }
1332
1338 void ResizeRuntimeRealComp (int new_size, bool communicate);
1339
1345 void ResizeRuntimeIntComp (int new_size, bool communicate);
1346
1348 template <template<class> class NewAllocator=amrex::DefaultAllocator>
1350
1361 template <template<class> class NewAllocator=amrex::DefaultAllocator>
1363 make_alike () const
1364 {
1366
1367 std::vector<std::string> const real_names = this->GetRealSoANames();
1368 std::vector<std::string> real_ct_names(NArrayReal);
1369 for (int ic = 0; ic < NArrayReal; ++ic) { real_ct_names.at(ic) = real_names[ic]; }
1370
1371 std::vector<std::string> const int_names = this->GetIntSoANames();
1372 std::vector<std::string> int_ct_names(NArrayInt);
1373 for (int ic = 0; ic < NArrayInt; ++ic) { int_ct_names.at(ic) = int_names[ic]; }
1374
1375 tmp.SetSoACompileTimeNames(real_ct_names, int_ct_names);
1376
1377 // add runtime real comps to tmp
1378 for (int ic = 0; ic < this->NumRuntimeRealComps(); ++ic) {
1379 tmp.AddRealComp(real_names.at(ic + NArrayReal));
1380 }
1381
1382 // add runtime int comps to tmp
1383 for (int ic = 0; ic < this->NumRuntimeIntComps(); ++ic) {
1384 tmp.AddIntComp(int_names.at(ic + NArrayInt));
1385 }
1386
1389
1390 return tmp;
1391 }
1392
1395
1398 mutable bool usePrePost;
1399 mutable bool doUnlink;
1401 mutable int nOutFilesPrePost;
1407 mutable std::string HdrFileNamePrePost;
1409
1410protected:
1411
1425 template <typename P>
1426 bool Where (const P& prt, ParticleLocData& pld,
1427 int lev_min = 0, int lev_max = -1, int nGrow=0, int local_grid=-1) const;
1428
1429
1440 template <typename P>
1441 bool EnforcePeriodicWhere (P& prt, ParticleLocData& pld,
1442 int lev_min = 0, int lev_max = -1, int local_grid=-1) const;
1443
1444public:
1445 void
1446 WriteParticles (int level, std::ofstream& ofs, int fnum,
1447 Vector<int>& which, Vector<int>& count, Vector<Long>& where,
1448 const Vector<int>& write_real_comp, const Vector<int>& write_int_comp,
1449 const Vector<std::map<std::pair<int, int>,IntVector>>& particle_io_flags, bool is_checkpoint) const;
1450#ifdef AMREX_USE_HDF5
1451#include "AMReX_ParticlesHDF5.H"
1452#endif
1453
1455 void SetSoACompileTimeNames (std::vector<std::string> const & rdata_name, std::vector<std::string> const & idata_name);
1456
1458 std::vector<std::string> GetRealSoANames () const {return m_soa_rdata_names;}
1459
1461 std::vector<std::string> GetIntSoANames () const {return m_soa_idata_names;}
1462
1468 bool HasRealComp (std::string const & name);
1469
1475 bool HasIntComp (std::string const & name);
1476
1484 int GetRealCompIndex (std::string const & name);
1485
1493 int GetIntCompIndex (std::string const & name);
1494
1495protected:
1496
1497 template <class RTYPE>
1498 void ReadParticles (int cnt, int grd, int lev, std::ifstream& ifs, int finest_level_in_file, bool convert_ids);
1499
1500 void SetParticleSize ();
1501
1503
1504private:
1505 virtual void particlePostLocate (ParticleType& /*p*/, const ParticleLocData& /*pld*/,
1506 const int /*lev*/) {}
1507
1508 virtual void correctCellVectors (int /*old_index*/, int /*new_index*/,
1509 int /*grid*/, const ParticleType& /*p*/) {}
1510
1511 void RedistributeMPI (std::map<int, Vector<char> >& not_ours,
1512 int lev_min = 0, int lev_max = 0, int nGrow = 0, int local=0);
1513
1514 template <typename P>
1515 void locateParticle (P& p, ParticleLocData& pld,
1516 int lev_min, int lev_max, int nGrow, int local_grid=-1) const;
1517
1518 void Initialize ();
1519
1523
1527
1528 // names of both compile-time and runtime Real and Int SoA data
1529 std::vector<std::string> m_soa_rdata_names;
1530 std::vector<std::string> m_soa_idata_names;
1531};
1532
1533template <int T_NStructReal, int T_NStructInt, int T_NArrayReal, int T_NArrayInt, template<class> class Allocator, class CellAssignor>
1534using ParticleContainer = ParticleContainer_impl<Particle<T_NStructReal, T_NStructInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
1535
1536template <int T_NArrayReal, int T_NArrayInt, template<class> class Allocator=DefaultAllocator, class CellAssignor=DefaultAssignor>
1537using ParticleContainerPureSoA = ParticleContainer_impl<SoAParticle<T_NArrayReal, T_NArrayInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
1538
1539}
1540
1541#include "AMReX_ParticleInit.H"
1543#include "AMReX_ParticleIO.H"
1544
1545#ifdef AMREX_USE_HDF5
1546#include "AMReX_ParticleHDF5.H"
1547#endif
1548
1549#endif /*_PARTICLES_H_*/
Definition AMReX_GpuAllocators.H:114
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:550
A container for storing items in a set of bins.
Definition AMReX_DenseBins.H:77
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
static const RealDescriptor & Native64RealDescriptor()
Definition AMReX_FPC.cpp:137
static const RealDescriptor & Native32RealDescriptor()
NativeRealDescriptor is equivalent to Native32RealDescriptor if BL_FLOAT is used. Otherwise,...
Definition AMReX_FPC.cpp:124
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:73
Definition AMReX_MFIter.H:57
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:141
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:38
Definition AMReX_PODVector.H:262
Definition AMReX_ParIter.H:142
Definition AMReX_ParGDB.H:13
Definition AMReX_ParIter.H:32
Definition AMReX_ParIter.H:113
Definition AMReX_ParticleContainerBase.H:23
ParGDBBase * m_gdb
Definition AMReX_ParticleContainerBase.H:270
void Define(ParGDBBase *gdb)
Definition AMReX_ParticleContainerBase.H:83
int numLevels() const
the number of defined levels in the ParticleContainer
Definition AMReX_ParticleContainerBase.H:222
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:146
ParticleTileType & DefineAndReturnParticleTile(int lev, const Iterator &iter)
Define and return the ParticleTile for level "lev", and Iterator "iter".
Definition AMReX_ParticleContainer.H:1184
void WritePlotFilePre()
Definition AMReX_ParticleIO.H:553
void AssignDensity(int rho_index, Vector< std::unique_ptr< MultiFab > > &mf_to_be_filled, int lev_min, int ncomp, int finest_level, int ngrow=2) const
Functions depending the layout of the data. Use with caution.
Definition AMReX_AmrParticles.H:17
void 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:1133
void SetUseUnlink(bool tf) const
Definition AMReX_ParticleContainer.H:1245
Vector< ParticleLevel > m_particles
Definition AMReX_ParticleContainer.H:1526
Vector< Long > nParticlesAtLevelPrePost
Definition AMReX_ParticleContainer.H:1403
bool GetLevelDirectoriesCreated() const
Definition AMReX_ParticleContainer.H:1233
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:1196
size_t superparticle_size
Definition AMReX_ParticleContainer.H:1524
void AddIntComp(std::string const &name, int communicate=1)
Definition AMReX_ParticleContainer.H:1294
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:1221
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:2679
bool doUnlink
Definition AMReX_ParticleContainer.H:1399
Vector< Vector< int > > whichPrePost
Definition AMReX_ParticleContainer.H:1404
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:270
std::string HdrFileNamePrePost
Definition AMReX_ParticleContainer.H:1407
Vector< Vector< Long > > wherePrePost
Definition AMReX_ParticleContainer.H:1406
RealDescriptor ParticleRealDescriptor
Definition AMReX_ParticleContainer.H:178
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:2452
ParticleTileType & ParticlesAt(int lev, const Iterator &iter)
Return the ParticleTile for level "lev" and Iterator "iter". Non-const version.
Definition AMReX_ParticleContainer.H:1130
int nOutFilesPrePost
Definition AMReX_ParticleContainer.H:1401
int maxnextidPrePost
Definition AMReX_ParticleContainer.H:1400
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:27
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:746
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:187
bool levelDirectoriesCreated
Variables for i/o optimization saved for pre and post checkpoint.
Definition AMReX_ParticleContainer.H:1397
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:1579
std::vector< std::string > m_soa_idata_names
Definition AMReX_ParticleContainer.H:1530
static constexpr int NArrayInt
Number of extra integer components stored in struct-of-array form.
Definition AMReX_ParticleContainer.H:159
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:157
int GetMaxNextIDPrePost() const
Definition AMReX_ParticleContainer.H:1242
const ParticleTileType & ParticlesAt(int lev, const Iterator &iter) const
Return the ParticleTile for level "lev" and Iterator "iter". Const version.
Definition AMReX_ParticleContainer.H:1105
void WritePlotFilePost()
Definition AMReX_ParticleIO.H:563
bool GetUsePrePost() const
Definition AMReX_ParticleContainer.H:1238
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:188
ParticleInitType< NStructReal, NStructInt, NArrayReal, NArrayInt > ParticleInitData
Definition AMReX_ParticleContainer.H:183
typename SoA::RealVector RealVector
Definition AMReX_ParticleContainer.H:191
typename ParticleType::ConstType ConstParticleType
Definition AMReX_ParticleContainer.H:149
ParticleContainer_impl(const ParticleContainer_impl &)=delete
void WriteAsciiFile(const std::string &file)
Definition AMReX_ParticleIO.H:1121
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:230
void locateParticle(P &p, ParticleLocData &pld, int lev_min, int lev_max, int nGrow, int local_grid=-1) const
Definition AMReX_ParticleContainerI.H:446
void SetParticleSize()
Definition AMReX_ParticleContainerI.H:16
void CheckpointPost()
Definition AMReX_ParticleIO.H:497
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:1231
void RemoveParticlesAtLevel(int level)
The Following methods are for managing Virtual and Ghost Particles.
Definition AMReX_ParticleContainerI.H:738
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:1055
typename AoS::ParticleVector ParticleVector
Definition AMReX_ParticleContainer.H:193
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:312
static constexpr int NStructInt
Number of extra integer components in the particle struct.
Definition AMReX_ParticleContainer.H:155
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:1353
ParticleContainer_impl()
Default constructor - construct an empty particle container that has no concept of a level hierarchy....
Definition AMReX_ParticleContainer.H:200
void SetUsePrePost(bool tf) const
Definition AMReX_ParticleContainer.H:1235
int NumRealComps() const
Definition AMReX_ParticleContainer.H:1330
void clearParticles()
Clear all the particles in this container. This does not free memory.
Definition AMReX_ParticleContainerI.H:1106
void ReadParticles(int cnt, int grd, int lev, std::ifstream &ifs, int finest_level_in_file, bool convert_ids)
Definition AMReX_ParticleIO.H:940
void Increment(MultiFab &mf, int level)
Definition AMReX_ParticleContainerI.H:701
Vector< int > h_redistribute_int_comp
Definition AMReX_ParticleContainer.H:1394
std::vector< std::string > m_soa_rdata_names
Definition AMReX_ParticleContainer.H:1529
typename ParticleTileType::SoA SoA
Definition AMReX_ParticleContainer.H:189
virtual void correctCellVectors(int, int, int, const ParticleType &)
Definition AMReX_ParticleContainer.H:1508
void AssignCellDensitySingleLevel(int rho_index, MultiFab &mf, int level, int ncomp=1, int particle_lvl_offset=0) const
Definition AMReX_ParticleContainerI.H:2513
int NumRuntimeIntComps() const
Definition AMReX_ParticleContainer.H:1328
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:1061
Long GetNParticlesPrePost() const
Definition AMReX_ParticleContainer.H:1243
void SortParticlesByBin(IntVect bin_size)
Sort the particles on each tile by groups of cells, given an IntVect bin_size.
Definition AMReX_ParticleContainerI.H:1320
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:150
const Vector< ParticleLevel > & GetParticles() const
Return the underlying Vector (over AMR levels) of ParticleLevels. Const version.
Definition AMReX_ParticleContainer.H:995
const ParticleLevel & GetParticles(int lev) const
Return the ParticleLevel for level "lev". Const version.
Definition AMReX_ParticleContainer.H:1016
std::vector< std::string > GetRealSoANames() const
Definition AMReX_ParticleContainer.H:1458
DenseBins< typename ParticleTileType::ParticleTileDataType > m_bins
Definition AMReX_ParticleContainer.H:1502
void CheckpointPre()
Definition AMReX_ParticleIO.H:440
void ShrinkToFit()
Definition AMReX_ParticleContainerI.H:681
int m_num_runtime_int
Definition AMReX_ParticleContainer.H:1522
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:214
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:349
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:1327
void ResizeRuntimeRealComp(int new_size, bool communicate)
Definition AMReX_ParticleContainerI.H:2653
void AddIntComp(int communicate=1)
Definition AMReX_ParticleContainer.H:1322
Long nparticlesPrePost
Definition AMReX_ParticleContainer.H:1402
int m_num_runtime_real
Definition AMReX_ParticleContainer.H:1521
Vector< std::string > filePrefixPrePost
Definition AMReX_ParticleContainer.H:1408
int numLocalTilesAtLevel(int lev) const
The total number of tiles on this rank on this level.
Definition AMReX_ParticleContainer.H:360
Long IncrementWithTotal(MultiFab &mf, int level, bool local=false)
Definition AMReX_ParticleContainerI.H:728
static constexpr int NStructReal
Number of extra Real components in the particle struct.
Definition AMReX_ParticleContainer.H:153
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:1311
Vector< int > h_redistribute_real_comp
Definition AMReX_ParticleContainer.H:1393
typename Particle< NStructReal, NStructInt >::RealType RealType
The type of the Real data.
Definition AMReX_ParticleContainer.H:173
ParticleTile< ParticleType, NArrayReal, NArrayInt, Allocator > ParticleTileType
Definition AMReX_ParticleContainer.H:182
bool GetUseUnlink() const
Definition AMReX_ParticleContainer.H:1249
int GetRealCompIndex(std::string const &name)
Definition AMReX_ParticleContainerI.H:161
void InitFromBinaryMetaFile(const std::string &file, int extradata)
Definition AMReX_ParticleInit.H:931
void Initialize()
Definition AMReX_ParticleContainerI.H:45
void RedistributeMPI(std::map< int, Vector< char > > &not_ours, int lev_min=0, int lev_max=0, int nGrow=0, int local=0)
Definition AMReX_ParticleContainerI.H:2099
int NumIntComps() const
Definition AMReX_ParticleContainer.H:1331
void RedistributeGPU(int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
Definition AMReX_ParticleContainerI.H:1383
ParticleContainer_impl & operator=(const ParticleContainer_impl &)=delete
Allocator< T > AllocatorType
The memory allocator in use.
Definition AMReX_ParticleContainer.H:169
ParticleLevel & GetParticles(int lev)
Return the ParticleLevel for level "lev". Non-const version.
Definition AMReX_ParticleContainer.H:1031
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:1155
std::vector< std::string > GetIntSoANames() const
Definition AMReX_ParticleContainer.H:1461
void RemoveParticlesNotAtFinestLevel()
Definition AMReX_ParticleContainerI.H:752
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:298
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:250
void Restart(const std::string &dir, const std::string &file)
Restart from checkpoint.
Definition AMReX_ParticleIO.H:637
Long superParticleSize() const
Definition AMReX_ParticleContainer.H:1259
typename SoA::IntVector IntVector
Definition AMReX_ParticleContainer.H:192
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:1122
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:573
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
int num_int_comm_comps
Definition AMReX_ParticleContainer.H:1525
size_t particle_size
Definition AMReX_ParticleContainer.H:1524
virtual void particlePostLocate(ParticleType &, const ParticleLocData &, const int)
Definition AMReX_ParticleContainer.H:1505
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:330
std::array< Long, 3 > PrintCapacity() const
Definition AMReX_ParticleContainerI.H:641
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:819
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:1080
void InitFromBinaryFile(const std::string &file, int extradata)
Definition AMReX_ParticleInit.H:485
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:110
void AddRealComp(std::string const &name, int communicate=1)
Definition AMReX_ParticleContainer.H:1261
int num_real_comm_comps
Definition AMReX_ParticleContainer.H:1525
std::array< Long, 3 > ByteSpread() const
Definition AMReX_ParticleContainerI.H:597
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:13
ContainerLike< NewAllocator > make_alike() const
Definition AMReX_ParticleContainer.H:1363
void AddRealComp(int communicate=1)
Definition AMReX_ParticleContainer.H:1289
void AddParticlesAtLevel(AoS &particles, int level, int nGrow=0)
Add particles from a pbox to the grid at this level.
Definition AMReX_ParticleContainerI.H:2467
~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:1001
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:416
T_ParticleType ParticleType
Definition AMReX_ParticleContainer.H:148
bool m_runtime_comps_defined
Definition AMReX_ParticleContainer.H:1520
Vector< Vector< int > > countPrePost
Definition AMReX_ParticleContainer.H:1405
bool usePrePost
Definition AMReX_ParticleContainer.H:1398
A Box with real dimensions. A RealBox is OK iff volume >= 0.
Definition AMReX_RealBox.H:21
A Descriptor of the Real Type.
Definition AMReX_FabConv.H:105
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:27
Long size() const noexcept
Definition AMReX_Vector.H:50
Definition AMReX_Amr.cpp:49
amrex::ArenaAllocator< T > DefaultAllocator
Definition AMReX_GpuAllocators.H:194
int verbose
Definition AMReX_DistributionMapping.cpp:36
Definition AMReX_ParticleUtil.H:432
A struct used for communicating particle data across processes during multi-level operations.
Definition AMReX_ParticleContainer.H:80
ParticleReal RealType
Definition AMReX_ParticleContainer.H:81
int m_grid
Definition AMReX_ParticleContainer.H:83
RealType m_data[1+AMREX_SPACEDIM]
Definition AMReX_ParticleContainer.H:85
IntVect m_cell
Definition AMReX_ParticleContainer.H:84
int m_lev
Definition AMReX_ParticleContainer.H:82
A struct used to pass initial data into the various Init methods. This struct is used to pass initial...
Definition AMReX_ParticleContainer.H:117
std::array< int, NStructInt > int_struct_data
Definition AMReX_ParticleContainer.H:119
std::array< int, NArrayInt > int_array_data
Definition AMReX_ParticleContainer.H:121
std::array< double, NArrayReal > real_array_data
Definition AMReX_ParticleContainer.H:120
std::array< double, NStructReal > real_struct_data
Definition AMReX_ParticleContainer.H:118
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:702
std::conditional_t< ParticleType::is_soa_particle, StructOfArrays< NArrayReal, NArrayInt, Allocator, true >, StructOfArrays< NArrayReal, NArrayInt, Allocator, false > > SoA
Definition AMReX_ParticleTile.H:725
std::conditional_t< ParticleType::is_soa_particle, ThisParticleTileHasNoAoS, ArrayOfStructs< ParticleType, Allocator > > AoS
Definition AMReX_ParticleTile.H:719
The struct used to store particles.
Definition AMReX_Particle.H:295
ParticleReal RealType
The floating point type used for the particles.
Definition AMReX_Particle.H:307