Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleContainerBase.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLECONTAINERBASE_H_
2#define AMREX_PARTICLECONTAINERBASE_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_Extension.H>
6#include <AMReX_INT.H>
7#include <AMReX_IntVect.H>
8#include <AMReX_ParGDB.H>
9#include <AMReX_Geometry.H>
11#include <AMReX_BoxArray.H>
12#include <AMReX_Vector.H>
13#include <AMReX_ParticleUtil.H>
14#include <AMReX_MultiFab.H>
15#include <AMReX_iMultiFab.H>
17#include <AMReX_DenseBins.H>
18
19#include <string>
20
21namespace amrex {
22
23#if defined(AMREX_USE_MPI)
40#endif
41
43{
44public:
45
47
49 :
50 m_gdb(gdb)
51 {}
52
54 const DistributionMapping & dmap,
55 const BoxArray & ba)
56 :
57 m_gdb_object(std::make_unique<ParGDB>(geom,dmap,ba)),
58 m_gdb(static_cast<ParGDBBase*>(m_gdb_object.get()))
59 {
60 }
61
63 const Vector<DistributionMapping> & dmap,
64 const Vector<BoxArray> & ba,
65 const Vector<int> & rr)
66 :
67 m_gdb_object(std::make_unique<ParGDB>(geom,dmap,ba,rr)),
68 m_gdb(static_cast<ParGDBBase*>(m_gdb_object.get()))
69 {
70 }
71
73 const Vector<DistributionMapping> & dmap,
74 const Vector<BoxArray> & ba,
75 const Vector<IntVect> & rr)
76 :
78 (std::make_unique<ParGDB>(geom,dmap,ba, [&]() -> Vector<int> {
79 Vector<int> ref_ratio;
80 for (auto const& iv : rr)
81 {
82#if AMREX_SPACEDIM > 1
83 AMREX_ASSERT(iv[0] == iv[1]);
84#endif
85#if AMREX_SPACEDIM > 2
86 AMREX_ASSERT(iv[0] == iv[2]);
87#endif
88 ref_ratio.push_back(iv[0]);
89 }
90 return ref_ratio; }() )),
91 m_gdb(static_cast<ParGDBBase*>(m_gdb_object.get()))
92 {
93 }
94
95 virtual ~ParticleContainerBase () = default;
96
99
101 ParticleContainerBase& operator= ( ParticleContainerBase && ) noexcept = default;
102
103 void Define (ParGDBBase* gdb) { m_gdb = gdb;}
104
105 void Define (const Geometry & geom,
106 const DistributionMapping & dmap,
107 const BoxArray & ba);
108
109 void Define (const Vector<Geometry> & geom,
110 const Vector<DistributionMapping> & dmap,
111 const Vector<BoxArray> & ba,
112 const Vector<int> & rr);
113
114 void Define (const Vector<Geometry> & geom,
115 const Vector<DistributionMapping> & dmap,
116 const Vector<BoxArray> & ba,
117 const Vector<IntVect> & rr);
118
119 bool isDefined () const { return m_gdb != nullptr; }
120
121 virtual void reserveData ();
122 virtual void resizeData ();
123 void RedefineDummyMF (int lev);
124
125 MFIter MakeMFIter (int lev, const MFItInfo& info) const {
126 AMREX_ASSERT(m_dummy_mf[lev] != nullptr);
127 return MFIter(*m_dummy_mf[lev], info);
128 }
129
130 MFIter MakeMFIter (int lev) const {
131 AMREX_ASSERT(m_dummy_mf[lev] != nullptr);
133 }
134
135 MFIter MakeMFIter (int lev, bool tile) const {
136 AMREX_ASSERT(m_dummy_mf[lev] != nullptr);
137 return MFIter(*m_dummy_mf[lev], tile ? tile_size : IntVect::TheZeroVector());
138 }
139
148 void SetParGDB (const Geometry & geom,
149 const DistributionMapping & dmap,
150 const BoxArray & ba);
151
161 void SetParGDB (const Vector<Geometry> & geom,
162 const Vector<DistributionMapping> & dmap,
163 const Vector<BoxArray> & ba,
164 const Vector<IntVect> & rr);
165
177 void SetParGDB (const Vector<Geometry> & geom,
178 const Vector<DistributionMapping> & dmap,
179 const Vector<BoxArray> & ba,
180 const Vector<int> & rr);
181
189 void SetParticleBoxArray (int lev, BoxArray new_ba);
190
198 void SetParticleDistributionMap (int lev, DistributionMapping new_dmap);
199
207 void SetParticleGeometry (int lev, Geometry new_geom);
208
213 const BoxArray& ParticleBoxArray (int lev) const
214 { return m_gdb->ParticleBoxArray(lev); }
215
221 { return m_gdb->ParticleDistributionMap(lev); }
222
227 const Geometry& Geom (int lev) const { return m_gdb->ParticleGeom(lev); }
228
233 const Geometry& ParticleGeom (int lev) const { return m_gdb->ParticleGeom(lev); }
234
236 int finestLevel () const { return m_gdb->finestLevel(); }
237
239 int maxLevel () const { return m_gdb->maxLevel(); }
240
242 int numLevels() const { return finestLevel() + 1; }
243
245 const ParGDBBase* GetParGDB () const { return m_gdb; }
246
248 ParGDBBase* GetParGDB () { return m_gdb; }
249
250 int Verbose () const { return m_verbose; }
251
252 void SetVerbose (int verbose) { m_verbose = verbose; }
253
254 [[nodiscard]] int stableRedistribute () const {return m_stable_redistribute; }
255
256 void setStableRedistribute (int stable) { m_stable_redistribute = stable; }
257
258 const ParticleBufferMap& BufferMap () const {return m_buffer_map;}
259
260#if defined(AMREX_USE_MPI)
261 void ensureParticleHandshakeWindow () const;
262 [[nodiscard]] Long* particleHandshakeBuffer () const
264 [[nodiscard]] MPI_Win particleHandshakeWindow () const
265 { return m_particle_handshake_window ? m_particle_handshake_window->win : MPI_WIN_NULL; }
266#endif
267
268 Vector<int> NeighborProcs(int ngrow) const
269 {
270 return computeNeighborProcs(this->GetParGDB(), ngrow);
271 }
272
273 template <class MF>
274 bool OnSameGrids (int level, const MF& mf) const { return m_gdb->OnSameGrids(level, mf); }
275
276 [[nodiscard]] Arena* arena () const {
277 return m_arena;
278 }
279
280 void SetArena (Arena* a) {
281 m_arena = a;
282 }
283
284 static const std::string& CheckpointVersion ();
285 static const std::string& PlotfileVersion ();
286 static const std::string& DataPrefix ();
287 static int MaxReaders ();
288 static Long MaxParticlesPerRead ();
289 static const std::string& AggregationType ();
290 static int AggregationBuffer ();
291
297
298protected:
299
300 void BuildRedistributeMask (int lev, int nghost=1) const;
301 void defineBufferMap () const;
302
303 int m_verbose{0};
305 std::unique_ptr<ParGDB> m_gdb_object = std::make_unique<ParGDB>();
306 ParGDBBase* m_gdb{nullptr};
308 Arena* m_arena = nullptr;
309
310 mutable std::unique_ptr<iMultiFab> redistribute_mask_ptr;
311 mutable int redistribute_mask_nghost = std::numeric_limits<int>::min();
314
315#if defined(AMREX_USE_MPI)
316 mutable std::unique_ptr<ParticleHandshakeWindow> m_particle_handshake_window;
317#endif
318
319};
320
321} // namespace amrex
322
323#endif // AMREX_PARTICLECONTAINERBASE_H_
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_EXPORT
Definition AMReX_Extension.H:191
Definition AMReX_ParticleLocator.H:282
A virtual base class for objects that manage their own dynamic memory allocation.
Definition AMReX_Arena.H:132
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
__host__ static __device__ constexpr IntVectND< dim > TheZeroVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:771
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
Definition AMReX_ParGDB.H:13
virtual int maxLevel() const =0
virtual const DistributionMapping & ParticleDistributionMap(int level) const =0
virtual const BoxArray & ParticleBoxArray(int level) const =0
bool OnSameGrids(int level, const MF &mf) const
Definition AMReX_ParGDB.H:133
virtual int finestLevel() const =0
virtual const Geometry & ParticleGeom(int level) const =0
we use this for non-Amr particle code
Definition AMReX_ParGDB.H:66
Definition AMReX_ParticleBufferMap.H:59
Definition AMReX_ParticleContainerBase.H:43
int stableRedistribute() const
Definition AMReX_ParticleContainerBase.H:254
ParticleContainerBase & operator=(const ParticleContainerBase &)=delete
amrex::Vector< int > neighbor_procs
Definition AMReX_ParticleContainerBase.H:312
bool isDefined() const
Definition AMReX_ParticleContainerBase.H:119
void BuildRedistributeMask(int lev, int nghost=1) const
Definition AMReX_ParticleContainerBase.cpp:332
std::unique_ptr< iMultiFab > redistribute_mask_ptr
Definition AMReX_ParticleContainerBase.H:310
ParticleContainerBase(const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
Definition AMReX_ParticleContainerBase.H:62
const DistributionMapping & ParticleDistributionMap(int lev) const
Get the DistributionMapping for a given level.
Definition AMReX_ParticleContainerBase.H:220
Vector< int > NeighborProcs(int ngrow) const
Definition AMReX_ParticleContainerBase.H:268
ParGDBBase * m_gdb
Definition AMReX_ParticleContainerBase.H:306
void SetVerbose(int verbose)
Definition AMReX_ParticleContainerBase.H:252
const BoxArray & ParticleBoxArray(int lev) const
Get the BoxArray for a given level.
Definition AMReX_ParticleContainerBase.H:213
virtual void reserveData()
Definition AMReX_ParticleContainerBase.cpp:41
ParticleBufferMap m_buffer_map
Definition AMReX_ParticleContainerBase.H:313
void SetArena(Arena *a)
Definition AMReX_ParticleContainerBase.H:280
MFIter MakeMFIter(int lev, bool tile) const
Definition AMReX_ParticleContainerBase.H:135
void setStableRedistribute(int stable)
Definition AMReX_ParticleContainerBase.H:256
std::unique_ptr< ParticleHandshakeWindow > m_particle_handshake_window
Definition AMReX_ParticleContainerBase.H:316
MFIter MakeMFIter(int lev) const
Definition AMReX_ParticleContainerBase.H:130
static const std::string & PlotfileVersion()
Definition AMReX_ParticleContainerBase.cpp:220
static const std::string & AggregationType()
Definition AMReX_ParticleContainerBase.cpp:292
const Geometry & Geom(int lev) const
Get the Geometry for a given level.
Definition AMReX_ParticleContainerBase.H:227
MPI_Win particleHandshakeWindow() const
Definition AMReX_ParticleContainerBase.H:264
ParticleContainerBase(ParticleContainerBase &&) noexcept=default
ParGDBBase * GetParGDB()
Get the ParGDB object used to define this container.
Definition AMReX_ParticleContainerBase.H:248
static const std::string & CheckpointVersion()
Definition AMReX_ParticleContainerBase.cpp:204
void Define(ParGDBBase *gdb)
Definition AMReX_ParticleContainerBase.H:103
static bool memEfficientSort
Definition AMReX_ParticleContainerBase.H:294
ParticleContainerBase(const ParticleContainerBase &)=delete
ParticleContainerBase(const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
Definition AMReX_ParticleContainerBase.H:53
int numLevels() const
the number of defined levels in the ParticleContainer
Definition AMReX_ParticleContainerBase.H:242
static const std::string & DataPrefix()
Definition AMReX_ParticleContainerBase.cpp:235
Arena * arena() const
Definition AMReX_ParticleContainerBase.H:276
ParticleContainerBase(ParGDBBase *gdb)
Definition AMReX_ParticleContainerBase.H:48
std::unique_ptr< ParGDB > m_gdb_object
Definition AMReX_ParticleContainerBase.H:305
static Long MaxParticlesPerRead()
Definition AMReX_ParticleContainerBase.cpp:267
int m_verbose
Definition AMReX_ParticleContainerBase.H:303
void SetParticleBoxArray(int lev, BoxArray new_ba)
Set the particle BoxArray. If the container was previously set to to track the AMR hierarchy of an Am...
Definition AMReX_ParticleContainerBase.cpp:166
const ParGDBBase * GetParGDB() const
Get the ParGDB object used to define this container (const version)
Definition AMReX_ParticleContainerBase.H:245
bool OnSameGrids(int level, const MF &mf) const
Definition AMReX_ParticleContainerBase.H:274
void ensureParticleHandshakeWindow() const
Definition AMReX_ParticleContainerBase.cpp:96
virtual ~ParticleContainerBase()=default
static bool do_tiling
Definition AMReX_ParticleContainerBase.H:292
AmrParticleLocator< DenseBins< Box > > m_particle_locator
Definition AMReX_ParticleContainerBase.H:296
int finestLevel() const
the finest level actually defined for the ParticleContainer
Definition AMReX_ParticleContainerBase.H:236
Long * particleHandshakeBuffer() const
Definition AMReX_ParticleContainerBase.H:262
int redistribute_mask_nghost
Definition AMReX_ParticleContainerBase.H:311
void defineBufferMap() const
Definition AMReX_ParticleContainerBase.cpp:73
static int AggregationBuffer()
Definition AMReX_ParticleContainerBase.cpp:312
void RedefineDummyMF(int lev)
Definition AMReX_ParticleContainerBase.cpp:55
int maxLevel() const
the finest allowed level in the ParticleContainer, whether it is defined or not.
Definition AMReX_ParticleContainerBase.H:239
static bool use_comms_arena
Definition AMReX_ParticleContainerBase.H:295
Vector< std::unique_ptr< MultiFab > > m_dummy_mf
Definition AMReX_ParticleContainerBase.H:307
MFIter MakeMFIter(int lev, const MFItInfo &info) const
Definition AMReX_ParticleContainerBase.H:125
const ParticleBufferMap & BufferMap() const
Definition AMReX_ParticleContainerBase.H:258
void SetParticleGeometry(int lev, Geometry new_geom)
Set the particle Geometry. If the container was previously set to to track the AMR hierarchy of an Am...
Definition AMReX_ParticleContainerBase.cpp:192
void SetParticleDistributionMap(int lev, DistributionMapping new_dmap)
Set the particle DistributionMapping. If the container was previously set to to track the AMR hierarc...
Definition AMReX_ParticleContainerBase.cpp:179
const Geometry & ParticleGeom(int lev) const
Get the particle Geometry for a given level.
Definition AMReX_ParticleContainerBase.H:233
static int MaxReaders()
Definition AMReX_ParticleContainerBase.cpp:245
void SetParGDB(const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
Set the particle Geometry, DistributionMapping, and BoxArray. If the container was previously set to ...
Definition AMReX_ParticleContainerBase.cpp:137
Arena * m_arena
Definition AMReX_ParticleContainerBase.H:308
int Verbose() const
Definition AMReX_ParticleContainerBase.H:250
int m_stable_redistribute
Definition AMReX_ParticleContainerBase.H:304
static IntVect tile_size
Definition AMReX_ParticleContainerBase.H:293
ParticleContainerBase(const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< IntVect > &rr)
Definition AMReX_ParticleContainerBase.H:72
virtual void resizeData()
Definition AMReX_ParticleContainerBase.cpp:46
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
amrex_long Long
Definition AMReX_INT.H:30
int MPI_Comm
Definition AMReX_ccse-mpi.H:51
static constexpr int MPI_COMM_NULL
Definition AMReX_ccse-mpi.H:59
Definition AMReX_Amr.cpp:49
Vector< int > computeNeighborProcs(const ParGDBBase *a_gdb, int ngrow)
Definition AMReX_ParticleUtil.cpp:22
const int[]
Definition AMReX_BLProfiler.cpp:1664
__host__ __device__ constexpr int get(IntVectND< dim > const &iv) noexcept
Get I'th element of IntVectND<dim>
Definition AMReX_IntVect.H:1335
Definition AMReX_MFIter.H:20
Definition AMReX_ParticleContainerBase.H:25
~ParticleHandshakeWindow()
Definition AMReX_ParticleContainerBase.cpp:84
ParticleHandshakeWindow(ParticleHandshakeWindow &&)=delete
Long * ptr
Definition AMReX_ParticleContainerBase.H:36
MPI_Win win
Definition AMReX_ParticleContainerBase.H:35
ParticleHandshakeWindow(const ParticleHandshakeWindow &)=delete
MPI_Comm comm
Definition AMReX_ParticleContainerBase.H:38
int nprocs
Definition AMReX_ParticleContainerBase.H:37
ParticleHandshakeWindow & operator=(const ParticleHandshakeWindow &)=delete