Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_BndryData.H
Go to the documentation of this file.
1
2#ifndef AMREX_BNDRYDATA_H_
3#define AMREX_BNDRYDATA_H_
4#include <AMReX_Config.H>
5
7#include <AMReX_BoundCond.H>
8#include <AMReX_MultiMask.H>
9
15#include <memory>
16#include <map>
17
18namespace amrex {
19
43template <typename MF>
45 : public BndryRegisterT<MF>
46{
47public:
49 enum MaskVal { covered = 0,
53
54 using value_type = typename MF::value_type;
55
57 BndryDataT() noexcept = default;
58
64 const DistributionMapping& dmap,
65 int ncomp,
66 const Geometry& geom);
67
69 ~BndryDataT () = default;
70
71 BndryDataT (const BndryDataT<MF>& src) = delete;
72 BndryDataT (BndryDataT<MF>&& rhs) = delete;
73 BndryDataT<MF>& operator= (const BndryDataT<MF>& src) = delete;
74 BndryDataT<MF>& operator= (BndryDataT<MF>&& rhs) = delete;
75
84 void define (const BoxArray& grids,
85 const DistributionMapping& dmap,
86 int ncomp,
87 const Geometry& geom);
88
90 const MultiMask& bndryMasks (Orientation face) const noexcept { return masks[face]; }
91
93 const FabSetT<MF>& bndryValues (Orientation face) const noexcept {
94 return this->bndry[face];
95 }
96
99
106 const RealTuple& bndryLocs (int igrid) const noexcept;
107
109 const RealTuple& bndryLocs (const MFIter& mfi) const noexcept;
110
117 const Vector< Vector<BoundCond> >& bndryConds (int igrid) const noexcept;
118
120 const Vector< Vector<BoundCond> >& bndryConds (const MFIter& mfi) const noexcept;
121
123 int nComp () const noexcept { return m_ncomp; }
124
126 const Box& getDomain () const noexcept { return geom.Domain(); }
127
129 const Geometry& getGeom () const noexcept { return geom; }
130
132 void setValue (Orientation face, int n, Real val) noexcept;
133
135 void setBoundCond (Orientation face,
136 int n,
137 int comp,
138 const BoundCond& bcn) noexcept;
139
141 void setBoundCond (Orientation face,
142 const MFIter& mfi,
143 int comp,
144 const BoundCond& bcn) noexcept;
145
147 void setBoundLoc (Orientation face,
148 int n,
149 value_type val) noexcept;
150
152 void setBoundLoc (Orientation face,
153 const MFIter& mfi,
154 value_type val) noexcept;
155
157 const FabSetT<MF>& operator[] (Orientation face) const noexcept { return BndryRegisterT<MF>::bndry[face]; }
158
161
162protected:
168
174 int m_ncomp{-1};
175 bool m_defined{false};
176
177private:
178 // Mask info required for this many cells past grid edge (here,
179 // e.g. ref=4, crse stencil width=3, and let stencil slide 2 past grid
180 // edge)
181 static constexpr int NTangHalfWidth = 5; // ref_ratio + 1, so won't work if ref_ratio > 4
182};
183
184template <typename MF>
186 const DistributionMapping& _dmap,
187 int _ncomp,
188 const Geometry& _geom)
189 :
190 geom(_geom),
191 m_ncomp(_ncomp)
192{
193 define(_grids,_dmap,_ncomp,_geom);
194}
195
196template <typename MF>
197void
199 int _n,
200 int _comp,
201 const BoundCond& _bcn) noexcept
202{
203 bcond[_n][_face][_comp] = _bcn;
204}
205
206template <typename MF>
207void
209 const MFIter& mfi,
210 int _comp,
211 const BoundCond& _bcn) noexcept
212{
213 bcond[mfi][_face][_comp] = _bcn;
214}
215
216template <typename MF>
217void
219 int _n,
220 value_type _val) noexcept
221{
222 bcloc[_n][_face] = _val;
223}
224
225template <typename MF>
226void
228 const MFIter& mfi,
229 value_type _val) noexcept
230{
231 bcloc[mfi][_face] = _val;
232}
233
234template <typename MF>
236BndryDataT<MF>::bndryConds (int igrid) const noexcept
237{
238 return bcond[igrid];
239}
240
241template <typename MF>
243BndryDataT<MF>::bndryConds (const MFIter& mfi) const noexcept
244{
245 return bcond[mfi];
246}
247
248template <typename MF>
249const typename BndryDataT<MF>::RealTuple&
250BndryDataT<MF>::bndryLocs (int igrid) const noexcept
251{
252 return bcloc[igrid];
253}
254
255template <typename MF>
256const typename BndryDataT<MF>::RealTuple&
257BndryDataT<MF>::bndryLocs (const MFIter& mfi) const noexcept
258{
259 return bcloc[mfi];
260}
261
262template <typename MF>
263void
265 const DistributionMapping& _dmap,
266 int _ncomp,
267 const Geometry& _geom)
268{
269 BL_PROFILE("BndryData::define()");
270
271 if (m_defined) {
272 if (m_ncomp == _ncomp && _geom.Domain() == geom.Domain() &&
273 _geom.periodicity() == geom.periodicity() &&
274 _grids == this->boxes() && _dmap == this->DistributionMap()) {
275 // We want to allow reuse of BndryDataT objects that were define()d exactly as a previous call.
276 return;
277 }
278 // Otherwise we'll just abort. We could make this work but it's just as easy to start with a fresh Bndrydata object.
279 amrex::Abort("BndryDataT<MF>::define(): object already built");
280 }
281 geom = _geom;
282 m_ncomp = _ncomp;
283
285
286 masks.clear();
287 masks.resize(2*AMREX_SPACEDIM);
288
289 for (OrientationIter fi; fi; ++fi) {
290 Orientation face = fi();
291 BndryRegisterT<MF>::define(face,IndexType::TheCellType(),0,1,1,_ncomp,_dmap);
292 masks[face].define(_grids, _dmap, geom, face, 0, 2, NTangHalfWidth, 1, true);
293 }
294 //
295 // Define "bcond" and "bcloc".
296 //
297 // We note that all orientations of the FabSets have the same distribution.
298 // We'll use the low 0 side as the model.
299 //
300 //
301
302 bcloc.define(_grids, _dmap);
303 bcond.define(_grids, _dmap);
304
305 for (FabSetIter bfsi(this->bndry[Orientation(0,Orientation::low)]);
306 bfsi.isValid(); ++bfsi) {
307 Vector< Vector<BoundCond> >& abc = bcond[bfsi];
308 abc.resize(2*AMREX_SPACEDIM);
309 for (OrientationIter fi; fi; ++fi) {
310 abc[fi()].resize(_ncomp);
311 }
312 }
313
314 m_defined = true;
315}
316
317template <typename MF>
318void
319BndryDataT<MF>::setValue (Orientation face, int n, Real val) noexcept
320{
321 auto& fab = this->bndry[face][n];
322 auto arr = fab.array();
323 const Box& bx = fab.box();
324 const int ncomp = m_ncomp;
325 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, m,
326 {
327 arr(i,j,k,m) = val;
328 });
329}
330
333
334}
335
336#endif /*_BNDRYDATA_H_*/
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
Infrastructure for storing per-face boundary data in FabSets.
Lightweight wrapper for integer boundary-condition identifiers.
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
MultiFab-style container that stores Mask objects per grid.
A BndryData stores and manipulates boundary data information on each side of each box in a BoxArray.
Definition AMReX_BndryData.H:46
LayoutData< Vector< Vector< BoundCond > > > bcond
Map of boundary condition type specifiers. The outer Array dimension is over Orientation.
Definition AMReX_BndryData.H:167
int nComp() const noexcept
Return number of components this boundary data manages.
Definition AMReX_BndryData.H:123
BndryDataT() noexcept=default
Default constructor.
bool m_defined
Definition AMReX_BndryData.H:175
const FabSetT< MF > & bndryValues(Orientation face) const noexcept
Access the boundary values stored on orientation face.
Definition AMReX_BndryData.H:93
const Geometry & getGeom() const noexcept
Return geometry used to define masks.
Definition AMReX_BndryData.H:129
const MultiMask & bndryMasks(Orientation face) const noexcept
Return the per-face mask built during define() for orientation face.
Definition AMReX_BndryData.H:90
void define(const BoxArray &grids, const DistributionMapping &dmap, int ncomp, const Geometry &geom)
Allocate masks and boundary registers for the supplied grids.
Definition AMReX_BndryData.H:264
const Box & getDomain() const noexcept
Return domain used to define masks.
Definition AMReX_BndryData.H:126
void setBoundCond(Orientation face, int n, int comp, const BoundCond &bcn) noexcept
Set boundary type at grid index n, component comp, and face face to bcn.
Definition AMReX_BndryData.H:198
Vector< MultiMask > masks
Boundary condition mask.
Definition AMReX_BndryData.H:171
typename MF::value_type value_type
Definition AMReX_BndryData.H:54
void setBoundLoc(Orientation face, int n, value_type val) noexcept
Set the physical boundary location val for grid index n and face face.
Definition AMReX_BndryData.H:218
MaskVal
Mask classification used for each boundary cell.
Definition AMReX_BndryData.H:49
@ not_covered
Location is inside the domain but does not coincide with any valid grid points.
Definition AMReX_BndryData.H:50
@ NumMaskVals
Definition AMReX_BndryData.H:52
@ covered
Location coincides with a valid grid point at the same AMR level.
Definition AMReX_BndryData.H:49
@ outside_domain
Location lies outside the physical domain.
Definition AMReX_BndryData.H:51
void setValue(Orientation face, int n, Real val) noexcept
Set all components on boundary face face of grid n to scalar value val.
Definition AMReX_BndryData.H:319
Array< value_type, 2 *3 > RealTuple
Some useful typedefs.
Definition AMReX_BndryData.H:98
int m_ncomp
Definition AMReX_BndryData.H:174
const Vector< Vector< BoundCond > > & bndryConds(int igrid) const noexcept
Return boundary type specifiers stored for a grid.
Definition AMReX_BndryData.H:236
Geometry geom
Domain used for mask definitions.
Definition AMReX_BndryData.H:173
const RealTuple & bndryLocs(int igrid) const noexcept
Return physical boundary locations for all faces owned by a grid.
Definition AMReX_BndryData.H:250
LayoutData< RealTuple > bcloc
Definition AMReX_BndryData.H:169
const FabSetT< MF > & operator[](Orientation face) const noexcept
Provide access to boundary registers for orientation face.
Definition AMReX_BndryData.H:157
A BndryRegister organizes FabSets bounding each grid in a BoxArray. A FabSet is maintained for each b...
Definition AMReX_BndryRegister.H:46
void setBoxes(const BoxArray &grids)
Set the BoxArray domain to grids once prior to calling define().
Definition AMReX_BndryRegister.H:239
BoxArray grids
Definition AMReX_BndryRegister.H:173
void define(const BoxArray &grids_, const DistributionMapping &dmap, int in_rad, int out_rad, int extent_rad, int ncomp)
Allocate boundary FabSets for all faces.
Definition AMReX_BndryRegister.H:193
FabSetT< MF > bndry[2 *3]
The data.
Definition AMReX_BndryRegister.H:172
Maintain an identifier for boundary condition types.
Definition AMReX_BoundCond.H:25
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
Convenience iterator that reuses MFIter semantics for FabSet traversal.
Definition AMReX_FabSet.H:281
A FabSet is a group of FArrayBox's. The grouping is designed specifically to represent regions along ...
Definition AMReX_FabSet.H:51
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:211
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:356
__host__ static __device__ constexpr IndexTypeND< dim > TheCellType() noexcept
This static member function returns an IndexTypeND object of value IndexTypeND::CELL....
Definition AMReX_IndexType.H:152
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
Definition AMReX_MultiMask.H:23
An Iterator over the Orientation of Faces of a Box.
Definition AMReX_Orientation.H:135
Encapsulation of the Orientation of the Faces of a Box.
Definition AMReX_Orientation.H:29
@ low
Definition AMReX_Orientation.H:34
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
std::array< T, N > Array
Definition AMReX_Array.H:26
Definition AMReX_Amr.cpp:49
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240