Block-Structured AMR Software Framework
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 
6 #include <AMReX_BndryRegister.H>
7 #include <AMReX_BoundCond.H>
8 #include <AMReX_MultiMask.H>
9 
10 #include <memory>
11 #include <map>
12 
13 namespace amrex {
14 
38 template <typename MF>
40  : public BndryRegisterT<MF>
41 {
42 public:
44  enum MaskVal { covered = 0, not_covered = 1, outside_domain = 2, NumMaskVals = 3 };
45 
46  using value_type = typename MF::value_type;
47 
49  BndryDataT() noexcept = default;
50 
55  BndryDataT (const BoxArray& grids,
56  const DistributionMapping& dmap,
57  int ncomp,
58  const Geometry& geom);
59 
61  ~BndryDataT () = default;
62 
63  BndryDataT (const BndryDataT<MF>& src) = delete;
64  BndryDataT (BndryDataT<MF>&& rhs) = delete;
65  BndryDataT<MF>& operator= (const BndryDataT<MF>& src) = delete;
66  BndryDataT<MF>& operator= (BndryDataT<MF>&& rhs) = delete;
67 
69  void define (const BoxArray& grids,
70  const DistributionMapping& dmap,
71  int ncomp,
72  const Geometry& geom);
73  //
74  const MultiMask& bndryMasks (Orientation face) const noexcept { return masks[face]; }
75 
77  const FabSetT<MF>& bndryValues (Orientation face) const noexcept {
78  return this->bndry[face];
79  }
80 
83 
89  const RealTuple& bndryLocs (int igrid) const noexcept;
90  const RealTuple& bndryLocs (const MFIter& mfi) const noexcept;
91 
96  const Vector< Vector<BoundCond> >& bndryConds (int igrid) const noexcept;
97  const Vector< Vector<BoundCond> >& bndryConds (const MFIter& mfi) const noexcept;
98 
100  int nComp () const noexcept { return m_ncomp; }
101 
103  const Box& getDomain () const noexcept { return geom.Domain(); }
104 
106  const Geometry& getGeom () const noexcept { return geom; }
107 
109  void setValue (Orientation face, int n, Real val) noexcept;
110 
112  void setBoundCond (Orientation face,
113  int n,
114  int comp,
115  const BoundCond& bcn) noexcept;
116 
117  void setBoundCond (Orientation face,
118  const MFIter& mfi,
119  int comp,
120  const BoundCond& bcn) noexcept;
121 
123  void setBoundLoc (Orientation face,
124  int n,
125  value_type val) noexcept;
126 
127  void setBoundLoc (Orientation face,
128  const MFIter& mfi,
129  value_type val) noexcept;
130 
132  const FabSetT<MF>& operator[] (Orientation face) const noexcept { return BndryRegisterT<MF>::bndry[face]; }
133 
136 
137 protected:
143 
149  int m_ncomp{-1};
150  bool m_defined{false};
151 
152 private:
153  // Mask info required for this many cells past grid edge (here,
154  // e.g. ref=4, crse stencil width=3, and let stencil slide 2 past grid
155  // edge)
156  static constexpr int NTangHalfWidth = 5; // ref_ratio + 1, so won't work if ref_ratio > 4
157 };
158 
159 template <typename MF>
161  const DistributionMapping& _dmap,
162  int _ncomp,
163  const Geometry& _geom)
164  :
165  geom(_geom),
166  m_ncomp(_ncomp)
167 {
168  define(_grids,_dmap,_ncomp,_geom);
169 }
170 
171 template <typename MF>
172 void
174  int _n,
175  int _comp,
176  const BoundCond& _bcn) noexcept
177 {
178  bcond[_n][_face][_comp] = _bcn;
179 }
180 
181 template <typename MF>
182 void
184  const MFIter& mfi,
185  int _comp,
186  const BoundCond& _bcn) noexcept
187 {
188  bcond[mfi][_face][_comp] = _bcn;
189 }
190 
191 template <typename MF>
192 void
194  int _n,
195  value_type _val) noexcept
196 {
197  bcloc[_n][_face] = _val;
198 }
199 
200 template <typename MF>
201 void
203  const MFIter& mfi,
204  value_type _val) noexcept
205 {
206  bcloc[mfi][_face] = _val;
207 }
208 
209 template <typename MF>
211 BndryDataT<MF>::bndryConds (int igrid) const noexcept
212 {
213  return bcond[igrid];
214 }
215 
216 template <typename MF>
218 BndryDataT<MF>::bndryConds (const MFIter& mfi) const noexcept
219 {
220  return bcond[mfi];
221 }
222 
223 template <typename MF>
224 const typename BndryDataT<MF>::RealTuple&
225 BndryDataT<MF>::bndryLocs (int igrid) const noexcept
226 {
227  return bcloc[igrid];
228 }
229 
230 template <typename MF>
231 const typename BndryDataT<MF>::RealTuple&
232 BndryDataT<MF>::bndryLocs (const MFIter& mfi) const noexcept
233 {
234  return bcloc[mfi];
235 }
236 
237 template <typename MF>
238 void
240  const DistributionMapping& _dmap,
241  int _ncomp,
242  const Geometry& _geom)
243 {
244  BL_PROFILE("BndryData::define()");
245 
246  if (m_defined) {
247  if (_grids == this->boxes() && m_ncomp == _ncomp && _geom.Domain() == geom.Domain()) {
248  // We want to allow reuse of BndryDataT objects that were define()d exactly as a previous call.
249  return;
250  }
251  // Otherwise we'll just abort. We could make this work but it's just as easy to start with a fresh Bndrydata object.
252  amrex::Abort("BndryDataT<MF>::define(): object already built");
253  }
254  geom = _geom;
255  m_ncomp = _ncomp;
256 
258 
259  masks.clear();
260  masks.resize(2*AMREX_SPACEDIM);
261 
262  for (OrientationIter fi; fi; ++fi) {
263  Orientation face = fi();
264  BndryRegisterT<MF>::define(face,IndexType::TheCellType(),0,1,1,_ncomp,_dmap);
265  masks[face].define(_grids, _dmap, geom, face, 0, 2, NTangHalfWidth, 1, true);
266  }
267  //
268  // Define "bcond" and "bcloc".
269  //
270  // We note that all orientations of the FabSets have the same distribution.
271  // We'll use the low 0 side as the model.
272  //
273  //
274 
275  bcloc.define(_grids, _dmap);
276  bcond.define(_grids, _dmap);
277 
278  for (FabSetIter bfsi(this->bndry[Orientation(0,Orientation::low)]);
279  bfsi.isValid(); ++bfsi) {
280  Vector< Vector<BoundCond> >& abc = bcond[bfsi];
281  abc.resize(2*AMREX_SPACEDIM);
282  for (OrientationIter fi; fi; ++fi) {
283  abc[fi()].resize(_ncomp);
284  }
285  }
286 
287  m_defined = true;
288 }
289 
290 template <typename MF>
291 void
292 BndryDataT<MF>::setValue (Orientation face, int n, Real val) noexcept
293 {
294  auto& fab = this->bndry[face][n];
295  auto arr = fab.array();
296  const Box& bx = fab.box();
297  const int ncomp = m_ncomp;
298  AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, m,
299  {
300  arr(i,j,k,m) = val;
301  });
302 }
303 
306 
307 }
308 
309 #endif /*_BNDRYDATA_H_*/
#define BL_PROFILE(a)
Definition: AMReX_BLProfiler.H:551
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition: AMReX_GpuLaunch.nolint.H:55
A BndryData stores and manipulates boundary data information on each side of each box in a BoxArray.
Definition: AMReX_BndryData.H:41
LayoutData< Vector< Vector< BoundCond > > > bcond
Map of boundary condition type specifiers. The outer Array dimension is over Orientation.
Definition: AMReX_BndryData.H:142
int nComp() const noexcept
return number of components for which this object is intended
Definition: AMReX_BndryData.H:100
BndryDataT() noexcept=default
Default constructor.
bool m_defined
Definition: AMReX_BndryData.H:150
const MultiMask & bndryMasks(Orientation face) const noexcept
Definition: AMReX_BndryData.H:74
Array< value_type, 2 *AMREX_SPACEDIM > RealTuple
Some useful typedefs.
Definition: AMReX_BndryData.H:82
const Geometry & getGeom() const noexcept
return geometry used to define masks
Definition: AMReX_BndryData.H:106
void define(const BoxArray &grids, const DistributionMapping &dmap, int ncomp, const Geometry &geom)
allocate bndry fabs along given face
Definition: AMReX_BndryData.H:239
const FabSetT< MF > & operator[](Orientation face) const noexcept
implement public access to const BndryRegister::operator[]
Definition: AMReX_BndryData.H:132
void setBoundCond(Orientation face, int n, int comp, const BoundCond &bcn) noexcept
set boundary type specifier for given orientation on nth grid
Definition: AMReX_BndryData.H:173
const Box & getDomain() const noexcept
return domain used to define masks
Definition: AMReX_BndryData.H:103
Vector< MultiMask > masks
Boundary condition mask.
Definition: AMReX_BndryData.H:146
typename MF::value_type value_type
Definition: AMReX_BndryData.H:46
void setBoundLoc(Orientation face, int n, value_type val) noexcept
set boundary location for given orientation on nth grid
Definition: AMReX_BndryData.H:193
MaskVal
mask values enumeration
Definition: AMReX_BndryData.H:44
@ not_covered
Definition: AMReX_BndryData.H:44
@ NumMaskVals
Definition: AMReX_BndryData.H:44
@ covered
Definition: AMReX_BndryData.H:44
@ outside_domain
Definition: AMReX_BndryData.H:44
void setValue(Orientation face, int n, Real val) noexcept
set values of boundary Fab for given orientation on nth grid
Definition: AMReX_BndryData.H:292
int m_ncomp
Definition: AMReX_BndryData.H:149
const FabSetT< MF > & bndryValues(Orientation face) const noexcept
Return FabSet on given face.
Definition: AMReX_BndryData.H:77
const Vector< Vector< BoundCond > > & bndryConds(int igrid) const noexcept
Return boundary type specifier on given face for grids we own. It's an error if we don't own that gri...
Definition: AMReX_BndryData.H:211
Geometry geom
Domain used for mask definitions.
Definition: AMReX_BndryData.H:148
static constexpr int NTangHalfWidth
Definition: AMReX_BndryData.H:156
const RealTuple & bndryLocs(int igrid) const noexcept
Return boundary location on given face for grids we own. It's an error if we don't own that grid....
Definition: AMReX_BndryData.H:225
LayoutData< RealTuple > bcloc
Definition: AMReX_BndryData.H:144
A BndryRegister organizes FabSets bounding each grid in a BoxArray. A FabSet is maintained for each b...
Definition: AMReX_BndryRegister.H:41
void setBoxes(const BoxArray &grids)
Set box domain, if not set previously.
Definition: AMReX_BndryRegister.H:193
BoxArray grids
Definition: AMReX_BndryRegister.H:126
FabSetT< MF > bndry[2 *AMREX_SPACEDIM]
The data.
Definition: AMReX_BndryRegister.H:125
void define(const BoxArray &grids_, const DistributionMapping &dmap, int in_rad, int out_rad, int extent_rad, int ncomp)
Definition: AMReX_BndryRegister.H:147
Maintain an identifier for boundary condition types.
Definition: AMReX_BoundCond.H:20
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:530
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
Definition: AMReX_FabSet.H:149
A FabSet is a group of FArrayBox's. The grouping is designed specifically to represent regions along ...
Definition: AMReX_FabSet.H:46
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition: AMReX_Geometry.H:210
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IndexTypeND< dim > TheCellType() noexcept
This static member function returns an IndexTypeND object of value IndexTypeND::CELL....
Definition: AMReX_IndexType.H:149
a one-thingy-per-box distributed object
Definition: AMReX_LayoutData.H:13
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
Definition: AMReX_MultiMask.H:18
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:27
static int fi(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition: AMReX_SundialsIntegrator.H:49
Definition: AMReX_Amr.cpp:49
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition: AMReX.cpp:221
std::array< T, N > Array
Definition: AMReX_Array.H:23