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
18namespace amrex {
19
20class Orientation;
21
44template <typename MF>
46{
47public:
48
49 using value_type = typename MF::value_type;
50
52 BndryRegisterT () noexcept = default;
53
55 BndryRegisterT (const BoxArray& grids_, // NOLINT(modernize-pass-by-value)
56 const DistributionMapping& dmap,
57 int in_rad, int out_rad, int extent_rad, int ncomp);
58
60 ~BndryRegisterT () = default;
61
62 BndryRegisterT (BndryRegisterT<MF>&& rhs) noexcept = default;
63
64 BndryRegisterT (const BndryRegisterT<MF>& src) = delete;
65 BndryRegisterT& operator= (const BndryRegisterT<MF>& src) = delete;
66 BndryRegisterT& operator= (BndryRegisterT<MF>&& rhs) = delete;
67
78 void define (const BoxArray& grids_, const DistributionMapping& dmap,
79 int in_rad, int out_rad, int extent_rad, int ncomp);
80
92 void define (Orientation face, IndexType typ,
93 int in_rad, int out_rad, int extent_rad, int ncomp,
94 const DistributionMapping& dm);
95
97 void clear ();
98
100 const BoxArray& boxes () const noexcept { return grids; }
101
103 int size () const noexcept { return grids.size(); }
104
106 const FabSetT<MF>& operator[] (Orientation face) const noexcept { return bndry[face]; }
107
109 FabSetT<MF>& operator[] (Orientation face) noexcept { return bndry[face]; }
110
112 void setVal (value_type v);
113
116
119
121 BndryRegisterT<MF>& copyFrom (const MF& src, int nghost,
122 int src_comp, int dest_comp, int num_comp,
123 const Periodicity& period = Periodicity::NonPeriodic());
124
126 BndryRegisterT<MF>& plusFrom (const MF& src, int nghost,
127 int src_comp, int dest_comp, int num_comp,
128 const Periodicity& period = Periodicity::NonPeriodic());
129
131 BndryRegisterT<MF>& linComb (value_type a, const MF& mfa, int a_comp,
132 value_type b, const MF& mfb, int b_comp,
133 int dest_comp, int num_comp, int n_ghost = 0);
134
140 void setBoxes (const BoxArray& grids);
141
143 const DistributionMapping& DistributionMap () const noexcept { return bndry[0].DistributionMap(); }
144
151 void write (const std::string& name, std::ostream& os) const;
152
159 void read (const std::string& name, std::istream& is);
160
167 static void Copy (BndryRegisterT<MF>& dst, const BndryRegisterT<MF>& src);
168
169protected:
170
172 FabSetT<MF> bndry[2*AMREX_SPACEDIM];
174};
175
176template <typename MF>
177BndryRegisterT<MF>::BndryRegisterT (const BoxArray& grids_, // NOLINT(modernize-pass-by-value)
178 const DistributionMapping& dmap,
179 int in_rad, int out_rad,
180 int extent_rad, int ncomp)
181 : grids(grids_)
182{
183 AMREX_ASSERT(ncomp > 0 && !grids.empty() && grids[0].cellCentered());
184
185 for (OrientationIter face; face; ++face)
186 {
187 define(face(),IndexType::TheCellType(),in_rad,out_rad,extent_rad,ncomp,dmap);
188 }
189}
190
191template <typename MF>
192void
194 const DistributionMapping& dmap,
195 int in_rad, int out_rad,
196 int extent_rad, int ncomp)
197{
198 grids = grids_;
199 for (OrientationIter face; face; ++face)
200 {
201 define(face(),IndexType::TheCellType(),in_rad,out_rad,extent_rad,ncomp,dmap);
202 }
203}
204
205template <typename MF>
206void
208{
209 for (auto & i : bndry) {
210 i.clear();
211 }
212 grids.clear();
213}
214
215template <typename MF>
216void
218 int _out_rad, int _extent_rad, int _ncomp,
219 const DistributionMapping& dmap)
220{
221 BoxArray fsBA(grids, BATransformer(_face,_typ,_in_rad,_out_rad,_extent_rad));
222
223 FabSetT<MF>& fabs = bndry[_face];
224
225 AMREX_ASSERT(fabs.size() == 0);
226
227 fabs.define(fsBA,dmap,_ncomp);
228 //
229 // Go ahead and assign values to the boundary register fabs
230 // since in some places APPLYBC (specifically in the tensor
231 // operator) the boundary registers are used for a few calculations
232 // before the masks are tested to see if you need them.
233 //
234 fabs.setVal(std::numeric_limits<value_type>::quiet_NaN());
235}
236
237template <typename MF>
238void
240{
241 AMREX_ASSERT(grids.empty() && !_grids.empty() && _grids[0].cellCentered());
242
243 grids = _grids;
244 //
245 // Check that bndry regions are not allocated.
246 //
247 for (auto const& k : bndry) {
249 AMREX_ASSERT(k.size() == 0);
250 }
251}
252
253template <typename MF>
255{
256 for (OrientationIter face; face; ++face)
257 {
258 bndry[face()].setVal(v);
259 }
260}
261
262template <typename MF>
265{
266 AMREX_ASSERT(grids == rhs.grids);
267 for (OrientationIter face; face; ++face) {
268 const auto f = face();
269 const int ncomp = bndry[f].nComp();
270#ifdef AMREX_USE_OMP
271#pragma omp parallel if (Gpu::notInLaunchRegion())
272#endif
273 for (FabSetIter bfsi(rhs[f]); bfsi.isValid(); ++bfsi) {
274 const Box& bx = bfsi.validbox();
275 auto const sfab = rhs[f].array(bfsi);
276 auto dfab = bndry[f].array(bfsi);
277 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
278 {
279 dfab(i,j,k,n) += sfab(i,j,k,n);
280 });
281 }
282 }
283 return *this;
284}
285
286template <typename MF>
289{
290 return operator+=(rhs);
291}
292
293template <typename MF>
295BndryRegisterT<MF>::linComb (value_type a, const MF& mfa, int a_comp,
296 value_type b, const MF& mfb, int b_comp,
297 int dest_comp, int num_comp, int n_ghost)
298{
299 for (OrientationIter face; face; ++face)
300 {
301 bndry[face()].linComb(a, mfa, a_comp,
302 b, mfb, b_comp,
303 dest_comp, num_comp, n_ghost);
304 }
305 return *this;
306}
307
308template <typename MF>
310BndryRegisterT<MF>::copyFrom (const MF& src, int nghost,
311 int src_comp, int dest_comp, int num_comp,
312 const Periodicity& period)
313{
314 for (OrientationIter face; face; ++face) {
315 bndry[face()].multiFab().ParallelCopy_nowait(src,src_comp,dest_comp,num_comp,
316 nghost,0,period);
317 }
318
319 for (OrientationIter face; face; ++face) {
320 bndry[face()].multiFab().ParallelCopy_finish();
321 }
322 return *this;
323}
324
325template <typename MF>
327BndryRegisterT<MF>::plusFrom (const MF& src, int nghost,
328 int src_comp, int dest_comp, int num_comp,
329 const Periodicity& period)
330{
331 for (OrientationIter face; face; ++face)
332 {
333 bndry[face()].plusFrom(src,nghost,src_comp,dest_comp,num_comp,period);
334 }
335 return *this;
336}
337
338template <typename MF>
339void
340BndryRegisterT<MF>::write (const std::string& name, std::ostream& os) const
341{
343 {
344 grids.writeOn(os);
345 os << '\n';
346 }
347
348 for (OrientationIter face; face; ++face)
349 {
350 //
351 // Take name here and make a "new" name unique to each face.
352 // Simplest thing would probably to append "_n" to the name,
353 // where n is the integer value of face().
354 //
355 const int i = face();
356 AMREX_ASSERT(i >= 0 && i <= 7);
357
358 std::string facename = amrex::Concatenate(name + '_', i, 1);
359
360 bndry[face()].write(facename);
361 }
362}
363
364template <typename MF>
365void
366BndryRegisterT<MF>::read (const std::string& name, std::istream& is)
367{
368 BoxArray grids_in;
369
370 grids_in.readFrom(is);
371
372 if (!amrex::match(grids,grids_in)) {
373 amrex::Abort("BndryRegisterT<MF>::read: grids do not match");
374 }
375
376 for (OrientationIter face; face; ++face)
377 {
378 //
379 // Take name here and make a "new" name unique to each face.
380 // Simplest thing would probably to append "_n" to the name,
381 // where n is the integer value of face().
382 //
383 const int i = face();
384 AMREX_ASSERT(i >= 0 && i <= 7);
385
386 std::string facename = amrex::Concatenate(name + '_', i, 1);
387
388 bndry[face()].read(facename);
389 }
390}
391
392// Local copy function
393template <typename MF>
394void
396{
397 for (OrientationIter face; face; ++face)
398 {
399 FabSetT<MF>::Copy(dst[face()], src[face()]);
400 }
401}
402
405
406}
407
408#endif /*_BNDRYREGISTER_H_*/
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
Convenience wrapper for groups of FArrayBoxes located on grid boundaries.
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
Constants describing linear-operator boundary condition types.
A BndryRegister organizes FabSets bounding each grid in a BoxArray. A FabSet is maintained for each b...
Definition AMReX_BndryRegister.H:46
void setVal(value_type v)
Set all boundary FABs to the scalar value v across every component.
Definition AMReX_BndryRegister.H:254
int size() const noexcept
Return the number of grids in this domain.
Definition AMReX_BndryRegister.H:103
void write(const std::string &name, std::ostream &os) const
Write all FabSets to disk (used for checkpoint data).
Definition AMReX_BndryRegister.H:340
BndryRegisterT< MF > & copyFrom(const MF &src, int nghost, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic())
Fill boundary FABs with num_comp components copied from src (component src_comp) into dest_comp,...
Definition AMReX_BndryRegister.H:310
BndryRegisterT< MF > & operator+=(const BndryRegisterT< MF > &rhs)
Accumulate another register with identical layout into this one.
Definition AMReX_BndryRegister.H:264
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
BndryRegisterT< MF > & plus(const BndryRegisterT< MF > &rhs)
Convenience alias for operator+=.
Definition AMReX_BndryRegister.H:288
BndryRegisterT< MF > & plusFrom(const MF &src, int nghost, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic())
Increment boundary FABs with num_comp components from src, using the same arguments as copyFrom().
Definition AMReX_BndryRegister.H:327
const DistributionMapping & DistributionMap() const noexcept
Returns the DistributionMapping describing which MPI ranks own each boundary FAB.
Definition AMReX_BndryRegister.H:143
void clear()
Release all FABs and clear the associated BoxArray.
Definition AMReX_BndryRegister.H:207
const FabSetT< MF > & operator[](Orientation face) const noexcept
Return const set of FABs bounding the domain grid boxes on orientation face.
Definition AMReX_BndryRegister.H:106
static void Copy(BndryRegisterT< MF > &dst, const BndryRegisterT< MF > &src)
Local copy function for duplicating boundary registers.
Definition AMReX_BndryRegister.H:395
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
typename MF::value_type value_type
Definition AMReX_BndryRegister.H:49
FabSetT< MF > bndry[2 *3]
The data.
Definition AMReX_BndryRegister.H:172
void read(const std::string &name, std::istream &is)
Read all FabSets from disk (used for checkpoint data).
Definition AMReX_BndryRegister.H:366
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 the boundary, using components a_comp/b_comp and wri...
Definition AMReX_BndryRegister.H:295
const BoxArray & boxes() const noexcept
Get box domain (as an array of boxes).
Definition AMReX_BndryRegister.H:100
BndryRegisterT() noexcept=default
The default constructor.
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
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:468
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:615
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition AMReX_BoxArray.H:621
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
static void Copy(FabSetT< MF > &dst, const FabSetT< MF > &src)
Local copy function.
Definition AMReX_FabSet.H:522
void setVal(value_type val)
Fill every component of every FAB with scalar value val.
Definition AMReX_FabSet.H:387
int size() const noexcept
Return number of FABs stored in this set.
Definition AMReX_FabSet.H:120
void define(const BoxArray &grids, const DistributionMapping &dmap, int ncomp)
Define a FabSetT<MF> constructed via default constructor on grids/dmap with ncomp components.
Definition AMReX_FabSet.H:295
__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
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
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:289
Definition AMReX_Amr.cpp:49
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
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:1934
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240