2#ifndef AMREX_BNDRYREGISTER_H_
3#define AMREX_BNDRYREGISTER_H_
4#include <AMReX_Config.H>
57 int in_rad,
int out_rad,
int extent_rad,
int ncomp);
79 int in_rad,
int out_rad,
int extent_rad,
int ncomp);
93 int in_rad,
int out_rad,
int extent_rad,
int ncomp,
122 int src_comp,
int dest_comp,
int num_comp,
127 int src_comp,
int dest_comp,
int num_comp,
133 int dest_comp,
int num_comp,
int n_ghost = 0);
151 void write (
const std::string& name, std::ostream& os)
const;
159 void read (
const std::string& name, std::istream& is);
176template <
typename MF>
179 int in_rad,
int out_rad,
180 int extent_rad,
int ncomp)
191template <
typename MF>
195 int in_rad,
int out_rad,
196 int extent_rad,
int ncomp)
205template <
typename MF>
209 for (
auto & i : bndry) {
215template <
typename MF>
218 int _out_rad,
int _extent_rad,
int _ncomp,
221 BoxArray fsBA(grids, BATransformer(_face,_typ,_in_rad,_out_rad,_extent_rad));
227 fabs.
define(fsBA,dmap,_ncomp);
234 fabs.
setVal(std::numeric_limits<value_type>::quiet_NaN());
237template <
typename MF>
247 for (
auto const& k : bndry) {
253template <
typename MF>
258 bndry[face()].setVal(v);
262template <
typename MF>
268 const auto f = face();
269 const int ncomp = bndry[f].nComp();
271#pragma omp parallel if (Gpu::notInLaunchRegion())
274 const Box& bx = bfsi.validbox();
275 auto const sfab = rhs[f].array(bfsi);
276 auto dfab = bndry[f].array(bfsi);
279 dfab(i,j,k,n) += sfab(i,j,k,n);
286template <
typename MF>
290 return operator+=(rhs);
293template <
typename MF>
297 int dest_comp,
int num_comp,
int n_ghost)
301 bndry[face()].
linComb(a, mfa, a_comp,
303 dest_comp, num_comp, n_ghost);
308template <
typename MF>
311 int src_comp,
int dest_comp,
int num_comp,
315 bndry[face()].multiFab().ParallelCopy_nowait(src,src_comp,dest_comp,num_comp,
320 bndry[face()].multiFab().ParallelCopy_finish();
325template <
typename MF>
328 int src_comp,
int dest_comp,
int num_comp,
333 bndry[face()].
plusFrom(src,nghost,src_comp,dest_comp,num_comp,period);
338template <
typename MF>
355 const int i = face();
360 bndry[face()].write(facename);
364template <
typename MF>
373 amrex::Abort(
"BndryRegisterT<MF>::read: grids do not match");
383 const int i = face();
388 bndry[face()].read(facename);
393template <
typename MF>
#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