![]() |
Block-Structured AMR Software Framework
|
Flux Register. More...
#include <AMReX_FluxRegister.H>
Public Types | |
| enum | FrOp { COPY = 0 , ADD = 1 } |
| An enum that says whether to add or copy src data to members. More... | |
Public Types inherited from amrex::BndryRegisterT< MF > | |
| using | value_type = typename MF::value_type |
Public Member Functions | |
| FluxRegister () | |
| The default constructor. | |
| FluxRegister (const BoxArray &fine_boxes, const DistributionMapping &dm, const IntVect &ref_ratio, int fine_lev, int nvar) | |
| Constructor. | |
| ~FluxRegister ()=default | |
| The destructor. | |
| FluxRegister (FluxRegister &&rhs) noexcept=default | |
| FluxRegister (const FluxRegister &rhs)=delete | |
| FluxRegister & | operator= (const FluxRegister &rhs)=delete |
| FluxRegister & | operator= (FluxRegister &&rhs)=delete |
| void | define (const BoxArray &fine_boxes, const DistributionMapping &dm, const IntVect &ref_ratio, int fine_lev, int nvar) |
| Initialize after default construction. | |
| void | clear () |
| Release all internal data structures and reset to the default-constructed state. | |
| const IntVect & | refRatio () const noexcept |
| Returns the refinement ratio. | |
| int | fineLevel () const noexcept |
| Returns the level number of the fine level. | |
| int | crseLevel () const noexcept |
| Returns the level number of the coarse level (fineLevel()-1). | |
| int | nComp () const noexcept |
| The number of components. | |
| const BoxArray & | coarsenedBoxes () const noexcept |
| The coarsened boxes. | |
| Real | SumReg (int comp) const |
| Returns the sum of the registers. | |
| void | CrseInit (const MultiFab &mflx, const MultiFab &area, int dir, int srccomp, int destcomp, int numcomp, Real mult=-1.0, FrOp op=FluxRegister::COPY) |
| Initialize flux correction with coarse data. | |
| void | CrseInit (const MultiFab &mflx, int dir, int srccomp, int destcomp, int numcomp, Real mult=-1.0, FrOp op=FluxRegister::COPY) |
| Initialize flux correction with coarse data (area already applied). | |
| void | CrseAdd (const MultiFab &mflx, const MultiFab &area, int dir, int srccomp, int destcomp, int numcomp, Real mult, const Geometry &geom) |
| Add coarse fluxes to the register without overwriting existing data. | |
| void | CrseAdd (const MultiFab &mflx, int dir, int srccomp, int destcomp, int numcomp, Real mult, const Geometry &geom) |
Like CrseAdd above but assumes areas are already multiplied into mflx. | |
| void | FineAdd (const MultiFab &mflx, int dir, int srccomp, int destcomp, int numcomp, Real mult) |
| Increment flux correction with fine data (areas already applied). | |
| void | FineAdd (const MultiFab &mflx, const MultiFab &area, int dir, int srccomp, int destcomp, int numcomp, Real mult) |
| Increment flux correction with fine data, supplying explicit areas. | |
| void | FineAdd (const FArrayBox &flux, int dir, int boxno, int srccomp, int destcomp, int numcomp, Real mult, RunOn runon) noexcept |
| Increment flux correction with fine data. | |
| void | FineAdd (const FArrayBox &flux, const FArrayBox &area, int dir, int boxno, int srccomp, int destcomp, int numcomp, Real mult, RunOn runon) noexcept |
| Increment flux correction with fine flux and area. | |
| void | FineSetVal (int dir, int boxno, int destcomp, int numcomp, Real val, RunOn runon) noexcept |
Set flux correction data for a fine box (given by boxno) to val. This routine used by FLASH does NOT run on GPU for safety. | |
| void | Reflux (MultiFab &mf, const MultiFab &volume, Real scale, int scomp, int dcomp, int nc, const Geometry &crse_geom) |
| Apply flux correction and divide by coarse cell volumes. | |
| void | Reflux (MultiFab &mf, const MultiFab &volume, int dir, Real scale, int scomp, int dcomp, int nc, const Geometry &crse_geom) |
Overload that only corrects faces in direction dir using explicit cell volumes. | |
| void | Reflux (MultiFab &mf, Real scale, int scomp, int dcomp, int nc, const Geometry &crse_geom) |
| Constant-volume version of Reflux() (no explicit volume MultiFab). | |
| void | Reflux (MultiFab &mf, int dir, Real scale, int scomp, int dcomp, int nc, const Geometry &crse_geom) |
Overload that only corrects faces in direction dir assuming unit volumes. | |
| void | OverwriteFlux (Array< MultiFab *, 3 > const &crse_fluxes, Real scale, int srccomp, int destcomp, int numcomp, const Geometry &crse_geom) |
| Overwrite the coarse flux at the coarse/fine interface (and the interface only) with the fine flux stored in the FluxRegister. | |
| void | ClearInternalBorders (const Geometry &crse_geom) |
| Set internal coarse/fine border fluxes to zero before reuse. | |
| void | write (const std::string &name, std::ostream &os) const |
| Write data to checkpoint files. | |
| void | read (const std::string &name, std::istream &is) |
| Read data from checkpoint files. | |
| void | Reflux (MultiFab &mf, const MultiFab &volume, Orientation face, Real scale, int scomp, int dcomp, int nc, const Geometry &geom) |
| Low-level helper used on GPUs: apply correction for a single face orientation. | |
Public Member Functions inherited from amrex::BndryRegisterT< MF > | |
| BndryRegisterT () noexcept=default | |
| The default constructor. | |
| BndryRegisterT (const BoxArray &grids_, const DistributionMapping &dmap, int in_rad, int out_rad, int extent_rad, int ncomp) | |
| The constructor, given number of cells in/out, extent and number of components (assumes cell-centered boxes, and allocates cell-centered FABs) | |
| ~BndryRegisterT ()=default | |
| The destructor. | |
| BndryRegisterT (BndryRegisterT< MF > &&rhs) noexcept=default | |
| BndryRegisterT (const BndryRegisterT< MF > &src)=delete | |
| BndryRegisterT & | operator= (const BndryRegisterT< MF > &src)=delete |
| BndryRegisterT & | operator= (BndryRegisterT< MF > &&rhs)=delete |
| 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. | |
| void | define (Orientation face, IndexType typ, int in_rad, int out_rad, int extent_rad, int ncomp, const DistributionMapping &dm) |
| Build FABs along a particular face orientation. | |
| void | clear () |
| Release all FABs and clear the associated BoxArray. | |
| const BoxArray & | boxes () const noexcept |
| Get box domain (as an array of boxes). | |
| int | size () const noexcept |
| Return the number of grids in this domain. | |
| const FabSetT< MF > & | operator[] (Orientation face) const noexcept |
Return const set of FABs bounding the domain grid boxes on orientation face. | |
| FabSetT< MF > & | operator[] (Orientation face) noexcept |
Return set of FABs bounding the domain grid boxes on orientation face. | |
| void | setVal (value_type v) |
Set all boundary FABs to the scalar value v across every component. | |
| BndryRegisterT< MF > & | operator+= (const BndryRegisterT< MF > &rhs) |
| Accumulate another register with identical layout into this one. | |
| BndryRegisterT< MF > & | plus (const BndryRegisterT< MF > &rhs) |
| Convenience alias for operator+=. | |
| 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, sampling nghost layers and honoring period. | |
| 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(). | |
| 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 writing to dest_comp for num_comp entries with n_ghost layers. | |
| void | setBoxes (const BoxArray &grids) |
Set the BoxArray domain to grids once prior to calling define(). | |
| const DistributionMapping & | DistributionMap () const noexcept |
| Returns the DistributionMapping describing which MPI ranks own each boundary FAB. | |
| void | write (const std::string &name, std::ostream &os) const |
| Write all FabSets to disk (used for checkpoint data). | |
| void | read (const std::string &name, std::istream &is) |
| Read all FabSets from disk (used for checkpoint data). | |
Additional Inherited Members | |
Static Public Member Functions inherited from amrex::BndryRegisterT< MF > | |
| static void | Copy (BndryRegisterT< MF > &dst, const BndryRegisterT< MF > &src) |
| Local copy function for duplicating boundary registers. | |
Protected Attributes inherited from amrex::BndryRegisterT< MF > | |
| FabSetT< MF > | bndry [2 *3] |
| The data. | |
| BoxArray | grids |
Flux Register.
Stores and manipulates fluxes at coarse-fine interfaces.
| amrex::FluxRegister::FluxRegister | ( | ) |
The default constructor.
| amrex::FluxRegister::FluxRegister | ( | const BoxArray & | fine_boxes, |
| const DistributionMapping & | dm, | ||
| const IntVect & | ref_ratio, | ||
| int | fine_lev, | ||
| int | nvar | ||
| ) |
Constructor.
| fine_boxes | Fine-level BoxArray. |
| dm | DistributionMapping describing ownership of fine boxes. |
| ref_ratio | Refinement ratio between coarse and fine. |
| fine_lev | Fine AMR level index. |
| nvar | Number of components carried in the registers. |
|
default |
The destructor.
|
defaultnoexcept |
|
delete |
| void amrex::FluxRegister::clear | ( | ) |
Release all internal data structures and reset to the default-constructed state.
| void amrex::FluxRegister::ClearInternalBorders | ( | const Geometry & | crse_geom | ) |
Set internal coarse/fine border fluxes to zero before reuse.
| crse_geom | Coarse Geometry describing the domain (needed for iteration). |
|
noexcept |
The coarsened boxes.
| void amrex::FluxRegister::CrseAdd | ( | const MultiFab & | mflx, |
| const MultiFab & | area, | ||
| int | dir, | ||
| int | srccomp, | ||
| int | destcomp, | ||
| int | numcomp, | ||
| Real | mult, | ||
| const Geometry & | geom | ||
| ) |
Add coarse fluxes to the register without overwriting existing data.
Unlike CrseInit(ADD), this path is used when fine-grid fluxes are unavailable on regions covered by fine levels (e.g., the FLASH code).
| void amrex::FluxRegister::CrseAdd | ( | const MultiFab & | mflx, |
| int | dir, | ||
| int | srccomp, | ||
| int | destcomp, | ||
| int | numcomp, | ||
| Real | mult, | ||
| const Geometry & | geom | ||
| ) |
Like CrseAdd above but assumes areas are already multiplied into mflx.
| mflx | Coarse flux with area weighting. |
| dir | Direction index. |
| srccomp | Starting component in mflx. |
| destcomp | Starting component in the register. |
| numcomp | Number of components. |
| mult | Scale factor applied to the flux. |
| geom | Coarse Geometry (needed for periodic wrap-around). |
| void amrex::FluxRegister::CrseInit | ( | const MultiFab & | mflx, |
| const MultiFab & | area, | ||
| int | dir, | ||
| int | srccomp, | ||
| int | destcomp, | ||
| int | numcomp, | ||
| Real | mult = -1.0, |
||
| FrOp | op = FluxRegister::COPY |
||
| ) |
Initialize flux correction with coarse data.
| mflx | Coarse flux MultiFab. |
| area | Area MultiFab on the coarse faces. |
| dir | Direction index of the flux. |
| srccomp | Starting component in mflx. |
| destcomp | Starting component stored in the register. |
| numcomp | Number of components to copy. |
| mult | Optional scale factor applied to the flux. |
| op | Whether to COPY or ADD the flux into the register. |
| void amrex::FluxRegister::CrseInit | ( | const MultiFab & | mflx, |
| int | dir, | ||
| int | srccomp, | ||
| int | destcomp, | ||
| int | numcomp, | ||
| Real | mult = -1.0, |
||
| FrOp | op = FluxRegister::COPY |
||
| ) |
Initialize flux correction with coarse data (area already applied).
| mflx | Coarse flux MultiFab with area weighting applied. |
| dir | Direction index. |
| srccomp | Starting component in mflx. |
| destcomp | Starting component in the register. |
| numcomp | Number of components. |
| mult | Optional scale factor. |
| op | Whether to COPY or ADD the flux into the register. |
|
noexcept |
Returns the level number of the coarse level (fineLevel()-1).
| void amrex::FluxRegister::define | ( | const BoxArray & | fine_boxes, |
| const DistributionMapping & | dm, | ||
| const IntVect & | ref_ratio, | ||
| int | fine_lev, | ||
| int | nvar | ||
| ) |
Initialize after default construction.
| fine_boxes | Fine-level BoxArray. |
| dm | DistributionMapping describing ownership. |
| ref_ratio | Refinement ratio between coarse and fine. |
| fine_lev | Fine AMR level index. |
| nvar | Number of components. |
|
noexcept |
Increment flux correction with fine flux and area.
| flux | Fine flux FAB. |
| area | Area FAB. |
| dir | Direction index. |
| boxno | Index of the owning FAB in the register. |
| srccomp | Starting component in flux. |
| destcomp | Starting component in the register. |
| numcomp | Number of components. |
| mult | Scale factor. |
| runon | Execution context (host/device) when copying into the register. |
|
noexcept |
Increment flux correction with fine data.
This overload assumes areas have already been multiplied into flux.
| flux | Fine flux FAB. |
| dir | Direction index. |
| boxno | Index of the owning FAB in the register. |
| srccomp | Starting component in flux. |
| destcomp | Starting component in the register. |
| numcomp | Number of components. |
| mult | Scale factor. |
| runon | Execution context (host/device) when copying into the register. |
| void amrex::FluxRegister::FineAdd | ( | const MultiFab & | mflx, |
| const MultiFab & | area, | ||
| int | dir, | ||
| int | srccomp, | ||
| int | destcomp, | ||
| int | numcomp, | ||
| Real | mult | ||
| ) |
Increment flux correction with fine data, supplying explicit areas.
| mflx | Fine flux MultiFab. |
| area | Fine face areas. |
| dir | Direction index. |
| srccomp | Starting component in mflx. |
| destcomp | Starting component in the register. |
| numcomp | Number of components. |
| mult | Scale factor. |
| void amrex::FluxRegister::FineAdd | ( | const MultiFab & | mflx, |
| int | dir, | ||
| int | srccomp, | ||
| int | destcomp, | ||
| int | numcomp, | ||
| Real | mult | ||
| ) |
Increment flux correction with fine data (areas already applied).
| mflx | Fine flux MultiFab. |
| dir | Direction index. |
| srccomp | Starting component in mflx. |
| destcomp | Starting component in the register. |
| numcomp | Number of components. |
| mult | Scale factor (typically +1). |
|
noexcept |
Returns the level number of the fine level.
|
noexcept |
Set flux correction data for a fine box (given by boxno) to val. This routine used by FLASH does NOT run on GPU for safety.
| dir | Direction index. |
| boxno | Index of the owning FAB in the register. |
| destcomp | Starting component to set. |
| numcomp | Number of components to set. |
| val | Scalar value to assign. |
| runon | Execution context (host/device). |
|
noexcept |
The number of components.
|
delete |
|
delete |
| void amrex::FluxRegister::OverwriteFlux | ( | Array< MultiFab *, 3 > const & | crse_fluxes, |
| Real | scale, | ||
| int | srccomp, | ||
| int | destcomp, | ||
| int | numcomp, | ||
| const Geometry & | crse_geom | ||
| ) |
Overwrite the coarse flux at the coarse/fine interface (and the interface only) with the fine flux stored in the FluxRegister.
| crse_fluxes | MultiFab pointers to coarse fluxes. |
| scale | scaling factor by which the fine flux is multiplied. |
| srccomp | starting component in fine flux in the FluxRegister |
| destcomp | starting component in coarse flux MultiFab |
| numcomp | number of components |
| crse_geom | coarse Geometry |
| void amrex::FluxRegister::read | ( | const std::string & | name, |
| std::istream & | is | ||
| ) |
Read data from checkpoint files.
| name | Base filename prefix. |
| is | Stream providing metadata. |
| void amrex::FluxRegister::Reflux | ( | MultiFab & | mf, |
| const MultiFab & | volume, | ||
| int | dir, | ||
| Real | scale, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | nc, | ||
| const Geometry & | crse_geom | ||
| ) |
Overload that only corrects faces in direction dir using explicit cell volumes.
| void amrex::FluxRegister::Reflux | ( | MultiFab & | mf, |
| const MultiFab & | volume, | ||
| Orientation | face, | ||
| Real | scale, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | nc, | ||
| const Geometry & | geom | ||
| ) |
Low-level helper used on GPUs: apply correction for a single face orientation.
| mf | Coarse MultiFab to correct. |
| volume | Coarse cell volumes. |
| face | Orientation being updated. |
| scale | Scaling factor applied to the accumulated flux. |
| scomp | Starting component within the register. |
| dcomp | Starting component within mf. |
| nc | Number of components to correct. |
| geom | Coarse Geometry (used for periodic wrap and face iteration). |
| void amrex::FluxRegister::Reflux | ( | MultiFab & | mf, |
| const MultiFab & | volume, | ||
| Real | scale, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | nc, | ||
| const Geometry & | crse_geom | ||
| ) |
Apply flux correction and divide by coarse cell volumes.
| mf | Coarse MultiFab to correct. |
| volume | Coarse cell volumes. |
| scale | Scaling factor applied to the correction. |
| scomp | Starting component within the register. |
| dcomp | Starting component within mf. |
| nc | Number of components to correct. |
| crse_geom | Coarse Geometry (used for periodic wrap and face orientation). |
| void amrex::FluxRegister::Reflux | ( | MultiFab & | mf, |
| int | dir, | ||
| Real | scale, | ||
| int | scomp, | ||
| int | dcomp, | ||
| int | nc, | ||
| const Geometry & | crse_geom | ||
| ) |
Overload that only corrects faces in direction dir assuming unit volumes.
|
noexcept |
Returns the refinement ratio.
Returns the sum of the registers.
| comp |
| void amrex::FluxRegister::write | ( | const std::string & | name, |
| std::ostream & | os | ||
| ) | const |
Write data to checkpoint files.
| name | Base filename prefix. |
| os | Stream storing metadata (e.g., BoxArray layout). |