Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::FluxRegister Class Reference

Flux Register. More...

#include <AMReX_FluxRegister.H>

Inheritance diagram for amrex::FluxRegister:
amrex::BndryRegisterT< MF >

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
 
FluxRegisteroperator= (const FluxRegister &rhs)=delete
 
FluxRegisteroperator= (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 IntVectrefRatio () 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 BoxArraycoarsenedBoxes () 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
 
BndryRegisterToperator= (const BndryRegisterT< MF > &src)=delete
 
BndryRegisterToperator= (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 BoxArrayboxes () 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 DistributionMappingDistributionMap () 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
 

Detailed Description

Flux Register.

Stores and manipulates fluxes at coarse-fine interfaces.

Member Enumeration Documentation

◆ FrOp

An enum that says whether to add or copy src data to members.

Enumerator
COPY 
ADD 

Constructor & Destructor Documentation

◆ FluxRegister() [1/4]

amrex::FluxRegister::FluxRegister ( )

The default constructor.

◆ FluxRegister() [2/4]

amrex::FluxRegister::FluxRegister ( const BoxArray fine_boxes,
const DistributionMapping dm,
const IntVect ref_ratio,
int  fine_lev,
int  nvar 
)

Constructor.

Parameters
fine_boxesFine-level BoxArray.
dmDistributionMapping describing ownership of fine boxes.
ref_ratioRefinement ratio between coarse and fine.
fine_levFine AMR level index.
nvarNumber of components carried in the registers.

◆ ~FluxRegister()

amrex::FluxRegister::~FluxRegister ( )
default

The destructor.

◆ FluxRegister() [3/4]

amrex::FluxRegister::FluxRegister ( FluxRegister &&  rhs)
defaultnoexcept

◆ FluxRegister() [4/4]

amrex::FluxRegister::FluxRegister ( const FluxRegister rhs)
delete

Member Function Documentation

◆ clear()

void amrex::FluxRegister::clear ( )

Release all internal data structures and reset to the default-constructed state.

◆ ClearInternalBorders()

void amrex::FluxRegister::ClearInternalBorders ( const Geometry crse_geom)

Set internal coarse/fine border fluxes to zero before reuse.

Parameters
crse_geomCoarse Geometry describing the domain (needed for iteration).

◆ coarsenedBoxes()

const BoxArray & amrex::FluxRegister::coarsenedBoxes ( ) const
noexcept

The coarsened boxes.

◆ CrseAdd() [1/2]

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).

Parameters
mflxCoarse flux MultiFab.
areaCoarse face areas.
dirDirection index.
srccompStarting component in mflx.
destcompStarting component in the register.
numcompNumber of components.
multScale factor applied to the flux.
geomCoarse Geometry (needed for periodic wrap-around).

◆ CrseAdd() [2/2]

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.

Parameters
mflxCoarse flux with area weighting.
dirDirection index.
srccompStarting component in mflx.
destcompStarting component in the register.
numcompNumber of components.
multScale factor applied to the flux.
geomCoarse Geometry (needed for periodic wrap-around).

◆ CrseInit() [1/2]

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.

Parameters
mflxCoarse flux MultiFab.
areaArea MultiFab on the coarse faces.
dirDirection index of the flux.
srccompStarting component in mflx.
destcompStarting component stored in the register.
numcompNumber of components to copy.
multOptional scale factor applied to the flux.
opWhether to COPY or ADD the flux into the register.

◆ CrseInit() [2/2]

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).

Parameters
mflxCoarse flux MultiFab with area weighting applied.
dirDirection index.
srccompStarting component in mflx.
destcompStarting component in the register.
numcompNumber of components.
multOptional scale factor.
opWhether to COPY or ADD the flux into the register.

◆ crseLevel()

int amrex::FluxRegister::crseLevel ( ) const
noexcept

Returns the level number of the coarse level (fineLevel()-1).

◆ define()

void amrex::FluxRegister::define ( const BoxArray fine_boxes,
const DistributionMapping dm,
const IntVect ref_ratio,
int  fine_lev,
int  nvar 
)

Initialize after default construction.

Parameters
fine_boxesFine-level BoxArray.
dmDistributionMapping describing ownership.
ref_ratioRefinement ratio between coarse and fine.
fine_levFine AMR level index.
nvarNumber of components.

◆ FineAdd() [1/4]

void amrex::FluxRegister::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.

Parameters
fluxFine flux FAB.
areaArea FAB.
dirDirection index.
boxnoIndex of the owning FAB in the register.
srccompStarting component in flux.
destcompStarting component in the register.
numcompNumber of components.
multScale factor.
runonExecution context (host/device) when copying into the register.

◆ FineAdd() [2/4]

void amrex::FluxRegister::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.

This overload assumes areas have already been multiplied into flux.

Parameters
fluxFine flux FAB.
dirDirection index.
boxnoIndex of the owning FAB in the register.
srccompStarting component in flux.
destcompStarting component in the register.
numcompNumber of components.
multScale factor.
runonExecution context (host/device) when copying into the register.

◆ FineAdd() [3/4]

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.

Parameters
mflxFine flux MultiFab.
areaFine face areas.
dirDirection index.
srccompStarting component in mflx.
destcompStarting component in the register.
numcompNumber of components.
multScale factor.

◆ FineAdd() [4/4]

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).

Parameters
mflxFine flux MultiFab.
dirDirection index.
srccompStarting component in mflx.
destcompStarting component in the register.
numcompNumber of components.
multScale factor (typically +1).

◆ fineLevel()

int amrex::FluxRegister::fineLevel ( ) const
noexcept

Returns the level number of the fine level.

◆ FineSetVal()

void amrex::FluxRegister::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.

Parameters
dirDirection index.
boxnoIndex of the owning FAB in the register.
destcompStarting component to set.
numcompNumber of components to set.
valScalar value to assign.
runonExecution context (host/device).

◆ nComp()

int amrex::FluxRegister::nComp ( ) const
noexcept

The number of components.

◆ operator=() [1/2]

FluxRegister & amrex::FluxRegister::operator= ( const FluxRegister rhs)
delete

◆ operator=() [2/2]

FluxRegister & amrex::FluxRegister::operator= ( FluxRegister &&  rhs)
delete

◆ OverwriteFlux()

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.

Parameters
crse_fluxesMultiFab pointers to coarse fluxes.
scalescaling factor by which the fine flux is multiplied.
srccompstarting component in fine flux in the FluxRegister
destcompstarting component in coarse flux MultiFab
numcompnumber of components
crse_geomcoarse Geometry

◆ read()

void amrex::FluxRegister::read ( const std::string &  name,
std::istream &  is 
)

Read data from checkpoint files.

Parameters
nameBase filename prefix.
isStream providing metadata.

◆ Reflux() [1/5]

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.

◆ Reflux() [2/5]

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.

Parameters
mfCoarse MultiFab to correct.
volumeCoarse cell volumes.
faceOrientation being updated.
scaleScaling factor applied to the accumulated flux.
scompStarting component within the register.
dcompStarting component within mf.
ncNumber of components to correct.
geomCoarse Geometry (used for periodic wrap and face iteration).

◆ Reflux() [3/5]

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.

Parameters
mfCoarse MultiFab to correct.
volumeCoarse cell volumes.
scaleScaling factor applied to the correction.
scompStarting component within the register.
dcompStarting component within mf.
ncNumber of components to correct.
crse_geomCoarse Geometry (used for periodic wrap and face orientation).

◆ Reflux() [4/5]

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.

◆ Reflux() [5/5]

void amrex::FluxRegister::Reflux ( MultiFab mf,
Real  scale,
int  scomp,
int  dcomp,
int  nc,
const Geometry crse_geom 
)

Constant-volume version of Reflux() (no explicit volume MultiFab).

Parameters
mfCoarse MultiFab to correct.
scaleScaling factor applied to the correction.
scompStarting component within the register.
dcompStarting component within mf.
ncNumber of components to correct.
crse_geomCoarse Geometry.

◆ refRatio()

const IntVect & amrex::FluxRegister::refRatio ( ) const
noexcept

Returns the refinement ratio.

◆ SumReg()

Real amrex::FluxRegister::SumReg ( int  comp) const

Returns the sum of the registers.

Parameters
comp

◆ write()

void amrex::FluxRegister::write ( const std::string &  name,
std::ostream &  os 
) const

Write data to checkpoint files.

Parameters
nameBase filename prefix.
osStream storing metadata (e.g., BoxArray layout).

The documentation for this class was generated from the following files: