Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_AmrParGDB.H
Go to the documentation of this file.
1#ifndef AMREX_AmrParGDB_H
2#define AMREX_AmrParGDB_H
3#include <AMReX_Config.H>
4
5#include <AMReX_ParGDB.H>
6#include <AMReX_AmrCore.H>
7
13namespace amrex {
14
16 : public ParGDBBase
17{
18 friend AmrCore;
19
20public:
21
27 explicit AmrParGDB (AmrCore* amr) noexcept
28 : m_amrcore(amr),
29 m_geom(amr->maxLevel()+1),
30 m_has_geom(amr->maxLevel()+1, 0),
31 m_dmap(amr->maxLevel()+1),
32 m_ba(amr->maxLevel()+1)
33 { }
34
36 [[nodiscard]] const Geometry& ParticleGeom (int level) const override;
38 [[nodiscard]] const Geometry& Geom (int level) const override;
39
41 [[nodiscard]] const Vector<Geometry>& ParticleGeom () const override;
43 [[nodiscard]] const Vector<Geometry>& Geom () const override;
44
46 [[nodiscard]] const DistributionMapping& ParticleDistributionMap (int level) const override;
48 [[nodiscard]] const DistributionMapping& DistributionMap (int level) const override;
49
51 [[nodiscard]] const Vector<DistributionMapping>& ParticleDistributionMap () const override;
53 [[nodiscard]] const Vector<DistributionMapping>& DistributionMap () const override;
54
56 [[nodiscard]] const BoxArray& ParticleBoxArray (int level) const override;
58 [[nodiscard]] const BoxArray& boxArray (int level) const override;
59
61 [[nodiscard]] const Vector<BoxArray>& ParticleBoxArray () const override;
63 [[nodiscard]] const Vector<BoxArray>& boxArray () const override;
64
66 void SetParticleBoxArray (int level, const BoxArray& new_ba) override;
68 void SetParticleDistributionMap (int level, const DistributionMapping& new_dm) override;
70 void SetParticleGeometry (int level, const Geometry& new_geom) override;
71
73 void ClearParticleBoxArray (int level) override;
75 void ClearParticleDistributionMap (int level) override;
77 void ClearParticleGeometry (int level) override;
78
80 [[nodiscard]] bool LevelDefined (int level) const override;
82 [[nodiscard]] int finestLevel () const override;
84 [[nodiscard]] int maxLevel () const override;
85
87 [[nodiscard]] IntVect refRatio (int level) const override;
89 [[nodiscard]] int MaxRefRatio (int level) const override;
90
92 [[nodiscard]] Vector<IntVect> refRatio () const override;
93
94protected:
95
101};
102
103inline
104const Geometry&
105AmrParGDB::ParticleGeom (int level) const
106{
107 if (! m_has_geom[level]) {
108 return m_amrcore->Geom(level);
109 } else {
110 return m_geom[level];
111 }
112}
113
114inline
115const Geometry&
116AmrParGDB::Geom (int level) const
117{
118 return m_amrcore->Geom(level);
119}
120
121inline
122const Vector<Geometry>&
124{
125 if (! m_has_geom[0]) {
126 return m_amrcore->Geom();
127 } else {
128 return m_geom;
129 }
130}
131
132inline
133const Vector<Geometry>&
135{
136 return m_amrcore->Geom();
137}
138
139inline
142{
143 if (m_dmap[level].empty()) {
144 return m_amrcore->DistributionMap(level);
145 } else {
146 return m_dmap[level];
147 }
148}
149
150inline
153{
154 return m_amrcore->DistributionMap(level);
155}
156
157inline
160{
161 if (m_dmap[0].empty()) {
162 return m_amrcore->DistributionMap();
163 } else {
164 return m_dmap;
165 }
166}
167
168inline
171{
172 return m_amrcore->DistributionMap();
173}
174
175inline
176const BoxArray&
178{
179 if (m_ba[level].empty()) {
180 return m_amrcore->boxArray(level);
181 } else {
182 return m_ba[level];
183 }
184}
185
186inline
187const BoxArray&
188AmrParGDB::boxArray (int level) const
189{
190 return m_amrcore->boxArray(level);
191}
192
193inline
194const Vector<BoxArray>&
196{
197 if (m_ba[0].empty()) {
198 return m_amrcore->boxArray();
199 } else {
200 return m_ba;
201 }
202}
203
204inline
205const Vector<BoxArray>&
207{
208 return m_amrcore->boxArray();
209}
210
211inline
212void AmrParGDB::SetParticleBoxArray (int level, const BoxArray& new_ba)
213{
214 m_ba[level] = new_ba;
215}
216
217inline
219{
220 m_dmap[level] = new_dmap;
221}
222
223inline
224void AmrParGDB::SetParticleGeometry (int level, const Geometry& new_geom)
225{
226 m_has_geom[level] = 1;
227 m_geom[level] = new_geom;
228}
229
230inline
232{
233 m_ba[level] = BoxArray();
234}
235
236inline
238{
239 m_dmap[level] = DistributionMapping();
240}
241
242inline
244{
245 m_geom[level] = Geometry();
246 m_has_geom[level] = 0;
247}
248
249inline
250bool
251AmrParGDB::LevelDefined (int level) const
252{
253 return m_amrcore->LevelDefined(level);
254}
255
256inline
257int
259{
260 return m_amrcore->finestLevel();
261}
262
263inline
264int
266{
267 return m_amrcore->maxLevel();
268}
269
270inline
272AmrParGDB::refRatio (int level) const
273{
274 return m_amrcore->refRatio(level);
275}
276
277inline
280{
281 return m_amrcore->refRatio();
282}
283
284inline
285int
286AmrParGDB::MaxRefRatio (int level) const
287{
288 return m_amrcore->MaxRefRatio(level);
289}
290
291}
292
293#endif
Base class for AMR applications that manages mesh hierarchy but not state data.
Provide basic functionalities to set up an AMR hierarchy.
Definition AMReX_AmrCore.H:31
const Vector< DistributionMapping > & DistributionMap() const noexcept
Definition AMReX_AmrMesh.H:150
int maxLevel() const noexcept
Return the max level.
Definition AMReX_AmrMesh.H:135
bool LevelDefined(int lev) const noexcept
Return true if valid grids/dmaps exist for level lev.
Definition AMReX_AmrMesh.cpp:443
int finestLevel() const noexcept
Return the finest level.
Definition AMReX_AmrMesh.H:138
IntVect refRatio(int lev) const noexcept
Return the refinement ratio for level lev.
Definition AMReX_AmrMesh.H:141
int MaxRefRatio(int lev) const noexcept
Return the maximum refinement ratio in any direction.
Definition AMReX_AmrMesh.cpp:391
const Vector< Geometry > & Geom() const noexcept
Definition AMReX_AmrMesh.H:149
const Vector< BoxArray > & boxArray() const noexcept
Definition AMReX_AmrMesh.H:151
Definition AMReX_AmrParGDB.H:17
AmrCore * m_amrcore
Definition AMReX_AmrParGDB.H:96
const Vector< DistributionMapping > & DistributionMap() const override
Return mesh DistributionMappings for all levels.
Definition AMReX_AmrParGDB.H:170
const Vector< BoxArray > & ParticleBoxArray() const override
Return particle BoxArrays for all levels.
Definition AMReX_AmrParGDB.H:195
Vector< IntVect > refRatio() const override
Return the refinement ratios for all levels.
Definition AMReX_AmrParGDB.H:279
Vector< int > m_has_geom
Definition AMReX_AmrParGDB.H:98
const Vector< Geometry > & Geom() const override
Return all mesh geometries from the underlying AmrCore.
Definition AMReX_AmrParGDB.H:134
Vector< BoxArray > m_ba
Definition AMReX_AmrParGDB.H:100
void ClearParticleGeometry(int level) override
Clear the cached particle Geometry on level level.
Definition AMReX_AmrParGDB.H:243
const Vector< BoxArray > & boxArray() const override
Return mesh BoxArrays for all levels.
Definition AMReX_AmrParGDB.H:206
int maxLevel() const override
Return the maximum level configured on the owning AmrCore.
Definition AMReX_AmrParGDB.H:265
void ClearParticleBoxArray(int level) override
Clear the cached particle BoxArray on level level.
Definition AMReX_AmrParGDB.H:231
void SetParticleDistributionMap(int level, const DistributionMapping &new_dm) override
Set a particle DistributionMapping for level level.
Definition AMReX_AmrParGDB.H:218
void SetParticleBoxArray(int level, const BoxArray &new_ba) override
Set a particle BoxArray for level level.
Definition AMReX_AmrParGDB.H:212
Vector< Geometry > m_geom
Definition AMReX_AmrParGDB.H:97
void ClearParticleDistributionMap(int level) override
Clear the cached particle DistributionMapping on level level.
Definition AMReX_AmrParGDB.H:237
const Vector< DistributionMapping > & ParticleDistributionMap() const override
Return particle DistributionMappings for all levels.
Definition AMReX_AmrParGDB.H:159
AmrParGDB(AmrCore *amr) noexcept
Construct a AmrParGDB connected to amr.
Definition AMReX_AmrParGDB.H:27
Vector< DistributionMapping > m_dmap
Definition AMReX_AmrParGDB.H:99
int finestLevel() const override
Return the current finest level known to the owning AmrCore.
Definition AMReX_AmrParGDB.H:258
bool LevelDefined(int level) const override
Return true if level level has been defined.
Definition AMReX_AmrParGDB.H:251
void SetParticleGeometry(int level, const Geometry &new_geom) override
Set a particle Geometry for level level.
Definition AMReX_AmrParGDB.H:224
int MaxRefRatio(int level) const override
Return the maximum scalar refinement ratio encountered up to level level.
Definition AMReX_AmrParGDB.H:286
const Vector< Geometry > & ParticleGeom() const override
Return all particle geometries.
Definition AMReX_AmrParGDB.H:123
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
Definition AMReX_ParGDB.H:13
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Definition AMReX_Amr.cpp:49