![]() |
Block-Structured AMR Software Framework
|
A FabSet is a group of FArrayBox's. The grouping is designed specifically to represent regions along the boundary of Box's, and are used to implement boundary conditions to discretized partial differential equations. More...
#include <AMReX_FabSet.H>
Public Types | |
| using | value_type = typename FabDataType< MF >::value_type |
| using | FAB = typename FabDataType< MF >::fab_type |
Public Member Functions | |
| FabSetT () noexcept=default | |
| The default constructor – you must later call define(). | |
| FabSetT (const BoxArray &grids, const DistributionMapping &dmap, int ncomp) | |
Construct a FabSetT<MF> covering grids with mapping dmap and ncomp components. | |
| ~FabSetT ()=default | |
| FabSetT (FabSetT< MF > &&rhs) noexcept=default | |
| FabSetT (const FabSetT< MF > &rhs)=delete | |
| FabSetT< MF > & | operator= (const FabSetT< MF > &rhs)=delete |
| FabSetT< MF > & | operator= (FabSetT< MF > &&rhs)=delete |
| void | define (const BoxArray &grids, const DistributionMapping &dmap, int ncomp) |
Define a FabSetT<MF> constructed via default constructor on grids/ ncomp components. | |
| FAB const & | operator[] (const MFIter &mfi) const noexcept |
Access the FAB referenced by iterator mfi (const). | |
| FAB & | operator[] (const MFIter &mfi) noexcept |
Access the FAB referenced by iterator mfi (mutable). | |
| FAB const & | operator[] (int i) const noexcept |
Access the FAB with index i (const). | |
| FAB & | operator[] (int i) noexcept |
Access the FAB with index i (mutable). | |
| Array4< value_type const > | array (const MFIter &mfi) const noexcept |
Return an Array4 view for iterator mfi (const). | |
| Array4< value_type > | array (const MFIter &mfi) noexcept |
Return an Array4 view for iterator mfi (mutable). | |
| Array4< value_type const > | array (int i) const noexcept |
Return an Array4 view for FAB index i (const). | |
| Array4< value_type > | array (int i) noexcept |
Return an Array4 view for FAB index i (mutable). | |
| Array4< value_type const > | const_array (const MFIter &mfi) const noexcept |
Return a const Array4 view for iterator mfi. | |
| Array4< value_type const > | const_array (int i) const noexcept |
Return a const Array4 view for FAB index i. | |
| MultiArray4< value_type const > | arrays () const noexcept |
| Return multi-fab const Array4 views to every FAB. | |
| MultiArray4< value_type > | arrays () noexcept |
| Return multi-fab mutable Array4 views to every FAB. | |
| MultiArray4< value_type const > | const_arrays () const noexcept |
| Return multi-fab const Array4 views (alias of arrays()). | |
| Box | fabbox (int K) const noexcept |
Get the bounding Box of FAB index K. | |
| int | size () const noexcept |
| Return number of FABs stored in this set. | |
| const BoxArray & | boxArray () const noexcept |
| Return the BoxArray defining FAB locations. | |
| const DistributionMapping & | DistributionMap () const noexcept |
| Return the DistributionMapping describing FAB ownership. | |
| MF & | multiFab () noexcept |
| Mutable access to the underlying MultiFab. | |
| MF const & | multiFab () const noexcept |
| Const access to the underlying MultiFab. | |
| int | nComp () const noexcept |
| Number of components carried by each FAB. | |
| void | clear () |
| Release storage and reset the internal MultiFab. | |
| FabSetT< MF > & | copyFrom (const FabSetT< MF > &src, int scomp, int dcomp, int ncomp) |
| Copy components from another FabSet with identical layout. | |
| FabSetT< MF > & | copyFrom (const MF &src, int ngrow, int scomp, int dcomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic()) |
| Copy boundary data from a MultiFab intersecting the same boxes. | |
| FabSetT< MF > & | plusFrom (const FabSetT< MF > &src, int scomp, int dcomp, int ncomp) |
| Accumulate data from another FabSet with identical layout. | |
| FabSetT< MF > & | plusFrom (const MF &src, int ngrow, int scomp, int dcomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic(), bool deterministic=false) |
| Accumulate data from a MultiFab into this boundary set. | |
| void | copyTo (MF &dest, int ngrow, int scomp, int dcomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic()) const |
| Copy data from this boundary set back into a MultiFab. | |
| void | plusTo (MF &dest, int ngrow, int scomp, int dcomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic()) const |
| Add this boundary data back into a MultiFab. | |
| void | setVal (value_type val) |
Fill every component of every FAB with scalar value val. | |
| void | setVal (value_type val, int comp, int num_comp) |
Fill num_comp components starting at comp with scalar value val. | |
| FabSetT< MF > & | linComb (value_type a, value_type b, const FabSetT< MF > &src, int scomp, int dcomp, int ncomp) |
| Linear combination: this := a*this + b*src. | |
| FabSetT< MF > & | linComb (value_type a, const MF &mfa, int a_comp, value_type b, const MF &mfb, int b_comp, int dcomp, int ncomp, int ngrow) |
| Boundary linear combination from two MultiFabs: this := a*mfa + b*mfb. | |
| void | write (const std::string &name) const |
| Write (used for writing to checkpoint). | |
| void | read (const std::string &name) |
| Read (used for reading from checkpoint). | |
Static Public Member Functions | |
| static void | Copy (FabSetT< MF > &dst, const FabSetT< MF > &src) |
| Local copy function. | |
Friends | |
| class | FabSetIter |
| class | FluxRegister |
A FabSet is a group of FArrayBox's. The grouping is designed specifically to represent regions along the boundary of Box's, and are used to implement boundary conditions to discretized partial differential equations.
A FabSet is an array of pointers to FABs. The standard FAB operators, however, have been modified to be more useful for maintaining boundary conditions for partial differential equations discretized on boxes. Under normal circumstances, a FAB will be created for each face of a box. For a group of boxes, a FabSet will be the group of FABs at a particular orientation (ie. the lo-i side of each grid in a list).
Since a FabSet FAB will likely be used to bound a grid box, FArrayBox::resize() operations are disallowed. Also, to preserve flexibility in applicable boundary scenarios, intersecting FABs in the FabSet are not guaranteed to contain identical data–thus copy operations from a FabSet to any FAB-like structure may be order-dependent.
FabSets are used primarily as a data storage mechanism, and are manipulated by more sophisticated control classes.
| using amrex::FabSetT< MF >::FAB = typename FabDataType<MF>::fab_type |
| using amrex::FabSetT< MF >::value_type = typename FabDataType<MF>::value_type |
|
defaultnoexcept |
The default constructor – you must later call define().
| amrex::FabSetT< MF >::FabSetT | ( | const BoxArray & | grids, |
| const DistributionMapping & | dmap, | ||
| int | ncomp | ||
| ) |
Construct a FabSetT<MF> covering grids with mapping dmap and ncomp components.
|
default |
|
defaultnoexcept |
|
delete |
|
inlinenoexcept |
Return an Array4 view for iterator mfi (const).
|
inlinenoexcept |
Return an Array4 view for iterator mfi (mutable).
|
inlinenoexcept |
Return an Array4 view for FAB index i (const).
|
inlinenoexcept |
Return an Array4 view for FAB index i (mutable).
|
inlinenoexcept |
Return multi-fab const Array4 views to every FAB.
|
inlinenoexcept |
Return multi-fab mutable Array4 views to every FAB.
|
inlinenoexcept |
Return the BoxArray defining FAB locations.
|
inline |
Release storage and reset the internal MultiFab.
|
inlinenoexcept |
Return a const Array4 view for iterator mfi.
|
inlinenoexcept |
Return a const Array4 view for FAB index i.
|
inlinenoexcept |
Return multi-fab const Array4 views (alias of arrays()).
|
static |
Local copy function.
| dst | Destination FabSet. |
| src | Source FabSet. |
| FabSetT< MF > & amrex::FabSetT< MF >::copyFrom | ( | const FabSetT< MF > & | src, |
| int | scomp, | ||
| int | dcomp, | ||
| int | ncomp | ||
| ) |
Copy components from another FabSet with identical layout.
| src | Source FabSet. |
| scomp | Starting component in src. |
| dcomp | Destination component index. |
| ncomp | Number of components to copy. |
| FabSetT< MF > & amrex::FabSetT< MF >::copyFrom | ( | const MF & | src, |
| int | ngrow, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | ncomp, | ||
| const Periodicity & | period = Periodicity::NonPeriodic() |
||
| ) |
Copy boundary data from a MultiFab intersecting the same boxes.
| src | Source MultiFab. |
| ngrow | Number of ghost cells pulled from src. |
| scomp | Source component index. |
| dcomp | Destination component index. |
| ncomp | Number of components. |
| period | Periodicity description. |
| void amrex::FabSetT< MF >::copyTo | ( | MF & | dest, |
| int | ngrow, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | ncomp, | ||
| const Periodicity & | period = Periodicity::NonPeriodic() |
||
| ) | const |
Copy data from this boundary set back into a MultiFab.
| dest | Destination MultiFab. |
| ngrow | Number of grow cells on dest to fill. |
| scomp | Source component index inside this set. |
| dcomp | Destination component index. |
| ncomp | Number of components. |
| period | Periodicity description. |
| void amrex::FabSetT< MF >::define | ( | const BoxArray & | grids, |
| const DistributionMapping & | dmap, | ||
| int | ncomp | ||
| ) |
Define a FabSetT<MF> constructed via default constructor on grids/ dmap withncomp components.
|
inlinenoexcept |
Return the DistributionMapping describing FAB ownership.
|
inlinenoexcept |
Get the bounding Box of FAB index K.
| FabSetT< MF > & amrex::FabSetT< MF >::linComb | ( | value_type | a, |
| const MF & | mfa, | ||
| int | a_comp, | ||
| value_type | b, | ||
| const MF & | mfb, | ||
| int | b_comp, | ||
| int | dcomp, | ||
| int | ncomp, | ||
| int | ngrow | ||
| ) |
Boundary linear combination from two MultiFabs: this := a*mfa + b*mfb.
| FabSetT< MF > & amrex::FabSetT< MF >::linComb | ( | value_type | a, |
| value_type | b, | ||
| const FabSetT< MF > & | src, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | ncomp | ||
| ) |
Linear combination: this := a*this + b*src.
| a | Scaling applied to the current FabSet. |
| b | Scaling applied to src. |
| src | Source FabSet (must match layout). |
| scomp | Source component offset. |
| dcomp | Destination component offset. |
| ncomp | Number of components. |
|
inlinenoexcept |
Const access to the underlying MultiFab.
|
inlinenoexcept |
Mutable access to the underlying MultiFab.
|
inlinenoexcept |
Number of components carried by each FAB.
|
delete |
|
delete |
|
inlinenoexcept |
Access the FAB referenced by iterator mfi (const).
|
inlinenoexcept |
Access the FAB referenced by iterator mfi (mutable).
|
inlinenoexcept |
Access the FAB with index i (const).
|
inlinenoexcept |
Access the FAB with index i (mutable).
| FabSetT< MF > & amrex::FabSetT< MF >::plusFrom | ( | const FabSetT< MF > & | src, |
| int | scomp, | ||
| int | dcomp, | ||
| int | ncomp | ||
| ) |
Accumulate data from another FabSet with identical layout.
| src | Source FabSet. |
| scomp | Source component index. |
| dcomp | Destination component index. |
| ncomp | Number of components to accumulate. |
| FabSetT< MF > & amrex::FabSetT< MF >::plusFrom | ( | const MF & | src, |
| int | ngrow, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | ncomp, | ||
| const Periodicity & | period = Periodicity::NonPeriodic(), |
||
| bool | deterministic = false |
||
| ) |
Accumulate data from a MultiFab into this boundary set.
| src | Source MultiFab. |
| ngrow | Number of ghost cells read from src. |
| scomp | Source component index. |
| dcomp | Destination component index. |
| ncomp | Number of components. |
| period | Periodicity description. |
| deterministic | Whether GPU operations must be deterministic. |
| void amrex::FabSetT< MF >::plusTo | ( | MF & | dest, |
| int | ngrow, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | ncomp, | ||
| const Periodicity & | period = Periodicity::NonPeriodic() |
||
| ) | const |
Add this boundary data back into a MultiFab.
| dest | Destination MultiFab. |
| ngrow | Number of grow cells to cover on dest. |
| scomp | Source component index inside this set. |
| dcomp | Destination component index. |
| ncomp | Number of components. |
| period | Periodicity description. |
| void amrex::FabSetT< MF >::read | ( | const std::string & | name | ) |
Read (used for reading from checkpoint).
| name | Base filename previously written by write(). |
| void amrex::FabSetT< MF >::setVal | ( | value_type | val | ) |
Fill every component of every FAB with scalar value val.
| void amrex::FabSetT< MF >::setVal | ( | value_type | val, |
| int | comp, | ||
| int | num_comp | ||
| ) |
Fill num_comp components starting at comp with scalar value val.
|
inlinenoexcept |
Return number of FABs stored in this set.
| void amrex::FabSetT< MF >::write | ( | const std::string & | name | ) | const |
Write (used for writing to checkpoint).
| name | Base filename to write. |
|
friend |
|
friend |