Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_BndryRegister.H
Go to the documentation of this file.
1
2#ifndef AMREX_BNDRYREGISTER_H_
3#define AMREX_BNDRYREGISTER_H_
4#include <AMReX_Config.H>
5
6#include <AMReX_BoxArray.H>
7#include <AMReX_FabSet.H>
8#include <AMReX_LO_BCTYPES.H>
9#include <AMReX_Orientation.H>
10#include <AMReX_Utility.H>
11#include <utility>
12
13namespace amrex {
14
15class Orientation;
16
39template <typename MF>
41{
42public:
43
44 using value_type = typename MF::value_type;
45
47 BndryRegisterT () noexcept = default;
48
50 BndryRegisterT (const BoxArray& grids_, // NOLINT(modernize-pass-by-value)
51 const DistributionMapping& dmap,
52 int in_rad, int out_rad, int extent_rad, int ncomp);
53
55 ~BndryRegisterT () = default;
56
57 BndryRegisterT (BndryRegisterT<MF>&& rhs) noexcept = default;
58
59 BndryRegisterT (const BndryRegisterT<MF>& src) = delete;
60 BndryRegisterT& operator= (const BndryRegisterT<MF>& src) = delete;
61 BndryRegisterT& operator= (BndryRegisterT<MF>&& rhs) = delete;
62
63 void define (const BoxArray& grids_, const DistributionMapping& dmap,
64 int in_rad, int out_rad, int extent_rad, int ncomp);
65
67 void define (Orientation face, IndexType typ,
68 int in_rad, int out_rad, int extent_rad, int ncomp,
69 const DistributionMapping& dm);
70
71 void clear ();
72
74 const BoxArray& boxes () const noexcept { return grids; }
75
77 int size () const noexcept { return grids.size(); }
78
80 const FabSetT<MF>& operator[] (Orientation face) const noexcept { return bndry[face]; }
81
83 FabSetT<MF>& operator[] (Orientation face) noexcept { return bndry[face]; }
84
86 void setVal (value_type v);
87
91
93 BndryRegisterT<MF>& copyFrom (const MF& src, int nghost,
94 int src_comp, int dest_comp, int num_comp,
95 const Periodicity& period = Periodicity::NonPeriodic());
96
98 BndryRegisterT<MF>& plusFrom (const MF& src, int nghost,
99 int src_comp, int dest_comp, int num_comp,
100 const Periodicity& period = Periodicity::NonPeriodic());
101
103 BndryRegisterT<MF>& linComb (value_type a, const MF& mfa, int a_comp,
104 value_type b, const MF& mfb, int b_comp,
105 int dest_comp, int num_comp, int n_ghost = 0);
106
108 void setBoxes (const BoxArray& grids);
109
111 const DistributionMapping& DistributionMap () const noexcept { return bndry[0].DistributionMap(); }
112
114 void write (const std::string& name, std::ostream& os) const;
115
117 void read (const std::string& name, std::istream& is);
118
120 static void Copy (BndryRegisterT<MF>& dst, const BndryRegisterT<MF>& src);
121
122protected:
123
125 FabSetT<MF> bndry[2*AMREX_SPACEDIM];
127};
128
129template <typename MF>
130BndryRegisterT<MF>::BndryRegisterT (const BoxArray& grids_, // NOLINT(modernize-pass-by-value)
131 const DistributionMapping& dmap,
132 int in_rad, int out_rad,
133 int extent_rad, int ncomp)
134 : grids(grids_)
135{
136 BL_ASSERT(ncomp > 0);
137 BL_ASSERT(grids[0].cellCentered());
138
139 for (OrientationIter face; face; ++face)
140 {
141 define(face(),IndexType::TheCellType(),in_rad,out_rad,extent_rad,ncomp,dmap);
142 }
143}
144
145template <typename MF>
146void
148 const DistributionMapping& dmap,
149 int in_rad, int out_rad,
150 int extent_rad, int ncomp)
151{
152 grids = grids_;
153 for (OrientationIter face; face; ++face)
154 {
155 define(face(),IndexType::TheCellType(),in_rad,out_rad,extent_rad,ncomp,dmap);
156 }
157}
158
159template <typename MF>
160void
162{
163 for (auto & i : bndry) {
164 i.clear();
165 }
166 grids.clear();
167}
168
169template <typename MF>
170void
172 int _out_rad, int _extent_rad, int _ncomp,
173 const DistributionMapping& dmap)
174{
175 BoxArray fsBA(grids, BATransformer(_face,_typ,_in_rad,_out_rad,_extent_rad));
176
177 FabSetT<MF>& fabs = bndry[_face];
178
179 BL_ASSERT(fabs.size() == 0);
180
181 fabs.define(fsBA,dmap,_ncomp);
182 //
183 // Go ahead and assign values to the boundary register fabs
184 // since in some places APPLYBC (specifically in the tensor
185 // operator) the boundary registers are used for a few calculations
186 // before the masks are tested to see if you need them.
187 //
188 fabs.setVal(std::numeric_limits<value_type>::quiet_NaN());
189}
190
191template <typename MF>
192void
194{
195 BL_ASSERT(grids.empty());
196 BL_ASSERT(!_grids.empty());
197 BL_ASSERT(_grids[0].cellCentered());
198
199 grids = _grids;
200 //
201 // Check that bndry regions are not allocated.
202 //
203 for (auto const& k : bndry) {
205 BL_ASSERT(k.size() == 0);
206 }
207}
208
209template <typename MF>
211{
212 for (OrientationIter face; face; ++face)
213 {
214 bndry[face()].setVal(v);
215 }
216}
217
218template <typename MF>
221{
222 BL_ASSERT(grids == rhs.grids);
223 for (OrientationIter face; face; ++face) {
224 const auto f = face();
225 const int ncomp = bndry[f].nComp();
226#ifdef AMREX_USE_OMP
227#pragma omp parallel if (Gpu::notInLaunchRegion())
228#endif
229 for (FabSetIter bfsi(rhs[f]); bfsi.isValid(); ++bfsi) {
230 const Box& bx = bfsi.validbox();
231 auto const sfab = rhs[f].array(bfsi);
232 auto dfab = bndry[f].array(bfsi);
233 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
234 {
235 dfab(i,j,k,n) += sfab(i,j,k,n);
236 });
237 }
238 }
239 return *this;
240}
241
242template <typename MF>
245{
246 return operator+=(rhs);
247}
248
249template <typename MF>
251BndryRegisterT<MF>::linComb (value_type a, const MF& mfa, int a_comp,
252 value_type b, const MF& mfb, int b_comp,
253 int dest_comp, int num_comp, int n_ghost)
254{
255 for (OrientationIter face; face; ++face)
256 {
257 bndry[face()].linComb(a, mfa, a_comp,
258 b, mfb, b_comp,
259 dest_comp, num_comp, n_ghost);
260 }
261 return *this;
262}
263
264template <typename MF>
266BndryRegisterT<MF>::copyFrom (const MF& src, int nghost,
267 int src_comp, int dest_comp, int num_comp,
268 const Periodicity& period)
269{
270 for (OrientationIter face; face; ++face) {
271 bndry[face()].multiFab().ParallelCopy_nowait(src,src_comp,dest_comp,num_comp,
272 nghost,0,period);
273 }
274
275 for (OrientationIter face; face; ++face) {
276 bndry[face()].multiFab().ParallelCopy_finish();
277 }
278 return *this;
279}
280
281template <typename MF>
283BndryRegisterT<MF>::plusFrom (const MF& src, int nghost,
284 int src_comp, int dest_comp, int num_comp,
285 const Periodicity& period)
286{
287 for (OrientationIter face; face; ++face)
288 {
289 bndry[face()].plusFrom(src,nghost,src_comp,dest_comp,num_comp,period);
290 }
291 return *this;
292}
293
294template <typename MF>
295void
296BndryRegisterT<MF>::write (const std::string& name, std::ostream& os) const
297{
299 {
300 grids.writeOn(os);
301 os << '\n';
302 }
303
304 for (OrientationIter face; face; ++face)
305 {
306 //
307 // Take name here and make a "new" name unique to each face.
308 // Simplest thing would probably to append "_n" to the name,
309 // where n is the integer value of face().
310 //
311 const int i = face();
312 BL_ASSERT(i >= 0 && i <= 7);
313
314 std::string facename = amrex::Concatenate(name + '_', i, 1);
315
316 bndry[face()].write(facename);
317 }
318}
319
320template <typename MF>
321void
322BndryRegisterT<MF>::read (const std::string& name, std::istream& is)
323{
324 BoxArray grids_in;
325
326 grids_in.readFrom(is);
327
328 if (!amrex::match(grids,grids_in)) {
329 amrex::Abort("BndryRegisterT<MF>::read: grids do not match");
330 }
331
332 for (OrientationIter face; face; ++face)
333 {
334 //
335 // Take name here and make a "new" name unique to each face.
336 // Simplest thing would probably to append "_n" to the name,
337 // where n is the integer value of face().
338 //
339 const int i = face();
340 BL_ASSERT(i >= 0 && i <= 7);
341
342 std::string facename = amrex::Concatenate(name + '_', i, 1);
343
344 bndry[face()].read(facename);
345 }
346}
347
348// Local copy function
349template <typename MF>
350void
352{
353 for (OrientationIter face; face; ++face)
354 {
355 FabSetT<MF>::Copy(dst[face()], src[face()]);
356 }
357}
358
361
362}
363
364#endif /*_BNDRYREGISTER_H_*/
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
A BndryRegister organizes FabSets bounding each grid in a BoxArray. A FabSet is maintained for each b...
Definition AMReX_BndryRegister.H:41
void setVal(value_type v)
Set all boundary FABs to given value.
Definition AMReX_BndryRegister.H:210
int size() const noexcept
Return the number of grids in this domain.
Definition AMReX_BndryRegister.H:77
void write(const std::string &name, std::ostream &os) const
Write (used for writing to checkpoint)
Definition AMReX_BndryRegister.H:296
BndryRegisterT< MF > & copyFrom(const MF &src, int nghost, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic())
Fill the boundary FABs on intersection with given MF.
Definition AMReX_BndryRegister.H:266
BndryRegisterT< MF > & operator+=(const BndryRegisterT< MF > &rhs)
register += rhs
Definition AMReX_BndryRegister.H:220
void setBoxes(const BoxArray &grids)
Set box domain, if not set previously.
Definition AMReX_BndryRegister.H:193
BoxArray grids
Definition AMReX_BndryRegister.H:126
BndryRegisterT< MF > & plus(const BndryRegisterT< MF > &rhs)
Definition AMReX_BndryRegister.H:244
BndryRegisterT< MF > & plusFrom(const MF &src, int nghost, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic())
Increment the boundary FABs on intersect with given MF.
Definition AMReX_BndryRegister.H:283
const DistributionMapping & DistributionMap() const noexcept
Returns constant reference to associated DistributionMapping.
Definition AMReX_BndryRegister.H:111
void clear()
Definition AMReX_BndryRegister.H:161
const FabSetT< MF > & operator[](Orientation face) const noexcept
Return const set of FABs bounding the domain grid boxes on a given orientation.
Definition AMReX_BndryRegister.H:80
static void Copy(BndryRegisterT< MF > &dst, const BndryRegisterT< MF > &src)
Local copy function.
Definition AMReX_BndryRegister.H:351
void define(const BoxArray &grids_, const DistributionMapping &dmap, int in_rad, int out_rad, int extent_rad, int ncomp)
Definition AMReX_BndryRegister.H:147
typename MF::value_type value_type
Definition AMReX_BndryRegister.H:44
FabSetT< MF > bndry[2 *3]
The data.
Definition AMReX_BndryRegister.H:125
void read(const std::string &name, std::istream &is)
Read (used for reading from checkpoint)
Definition AMReX_BndryRegister.H:322
BndryRegisterT< MF > & linComb(value_type a, const MF &mfa, int a_comp, value_type b, const MF &mfb, int b_comp, int dest_comp, int num_comp, int n_ghost=0)
Linear combination: this := a*mfa + b*mfb on intersection of MFs with the boundary FABs.
Definition AMReX_BndryRegister.H:251
const BoxArray & boxes() const noexcept
Get box domain (as an array of boxes).
Definition AMReX_BndryRegister.H:74
BndryRegisterT() noexcept=default
The default constructor.
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:551
int readFrom(std::istream &is)
Initialize the BoxArray from the supplied istream. It is an error if the BoxArray has already been in...
Definition AMReX_BoxArray.cpp:461
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:598
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition AMReX_BoxArray.H:604
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
static void Copy(FabSetT< MF > &dst, const FabSetT< MF > &src)
Definition AMReX_FabSet.H:389
void setVal(value_type val)
Definition AMReX_FabSet.H:254
int size() const noexcept
Definition AMReX_FabSet.H:90
void define(const BoxArray &grids, const DistributionMapping &dmap, int ncomp)
Define a FabSetT<MF> constructed via default constructor.
Definition AMReX_FabSet.H:163
__host__ static __device__ constexpr IndexTypeND< dim > TheCellType() noexcept
This static member function returns an IndexTypeND object of value IndexTypeND::CELL....
Definition AMReX_IndexType.H:149
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:141
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
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
static const Periodicity & NonPeriodic() noexcept
Definition AMReX_Periodicity.cpp:52
bool IOProcessor() noexcept
Is this CPU the I/O Processor? To get the rank number, call IOProcessorNumber()
Definition AMReX_ParallelDescriptor.H:275
Definition AMReX_Amr.cpp:49
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:138
std::string Concatenate(const std::string &root, int num, int mindigits)
Returns rootNNNN where NNNN == num.
Definition AMReX_String.cpp:34
bool match(const BoxArray &x, const BoxArray &y)
Note that two BoxArrays that match are not necessarily equal.
Definition AMReX_BoxArray.cpp:1918
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
Definition AMReX_BoxArray.H:277