Block-Structured AMR Software Framework
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. More...
 
 FluxRegister (const BoxArray &fine_boxes, const DistributionMapping &dm, const IntVect &ref_ratio, int fine_lev, int nvar)
 The constructor. This version allows setting the DistributionMapping. More...
 
 ~FluxRegister ()=default
 The destructor. More...
 
 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 using default constructor. This version allows setting the DistributionMapping. More...
 
void clear ()
 
const IntVectrefRatio () const noexcept
 Returns the refinement ratio. More...
 
int fineLevel () const noexcept
 Returns the level number of the fine level. More...
 
int crseLevel () const noexcept
 Returns the level number of the coarse level (fineLevel()-1). More...
 
int nComp () const noexcept
 The number of components. More...
 
const BoxArraycoarsenedBoxes () const noexcept
 The coarsened boxes. More...
 
Real SumReg (int comp) const
 Returns the sum of the registers. More...
 
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. More...
 
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. More...
 
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 flux register. This is different from CrseInit with FluxRegister::ADD. This is used for cases in which the grids covered by fine do not have fluxes computed (e.g., FLASH). More...
 
void CrseAdd (const MultiFab &mflx, int dir, int srccomp, int destcomp, int numcomp, Real mult, const Geometry &geom)
 /in this version the area is assumed to multiplied into the flux (if not, use scale to fix) More...
 
void FineAdd (const MultiFab &mflx, int dir, int srccomp, int destcomp, int numcomp, Real mult)
 Increment flux correction with fine data. More...
 
void FineAdd (const MultiFab &mflx, const MultiFab &area, int dir, int srccomp, int destcomp, int numcomp, Real mult)
 Increment flux correction with fine data. More...
 
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. More...
 
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 data. More...
 
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 a given value. This routine used by FLASH does NOT run on gpu for safety. More...
 
void Reflux (MultiFab &mf, const MultiFab &volume, Real scale, int scomp, int dcomp, int nc, const Geometry &crse_geom)
 Apply flux correction. Note that this takes the coarse Geometry. More...
 
void Reflux (MultiFab &mf, const MultiFab &volume, int dir, Real scale, int scomp, int dcomp, int nc, const Geometry &crse_geom)
 
void Reflux (MultiFab &mf, Real scale, int scomp, int dcomp, int nc, const Geometry &crse_geom)
 Constant volume version of Reflux(). Note that this takes the coarse Geometry. More...
 
void Reflux (MultiFab &mf, int dir, Real scale, int scomp, int dcomp, int nc, const Geometry &crse_geom)
 
void OverwriteFlux (Array< MultiFab *, AMREX_SPACEDIM > 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. More...
 
void ClearInternalBorders (const Geometry &crse_geom)
 Set internal borders to zero. More...
 
void write (const std::string &name, std::ostream &os) const
 Write (used for writing to checkpoint) More...
 
void read (const std::string &name, std::istream &is)
 Read (used for reading from checkpoint) More...
 
void Reflux (MultiFab &mf, const MultiFab &volume, Orientation face, Real scale, int scomp, int dcomp, int nc, const Geometry &geom)
 
- Public Member Functions inherited from amrex::BndryRegisterT< MF >
 BndryRegisterT () noexcept=default
 The default constructor. More...
 
 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) More...
 
 ~BndryRegisterT ()=default
 The destructor. More...
 
 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)
 
void define (Orientation face, IndexType typ, int in_rad, int out_rad, int extent_rad, int ncomp, const DistributionMapping &dm)
 Build FABs along given face, specifying the DistributionMapping. More...
 
void clear ()
 
const BoxArrayboxes () const noexcept
 Get box domain (as an array of boxes). More...
 
int size () const noexcept
 Return the number of grids in this domain. More...
 
const FabSetT< MF > & operator[] (Orientation face) const noexcept
 Return const set of FABs bounding the domain grid boxes on a given orientation. More...
 
FabSetT< MF > & operator[] (Orientation face) noexcept
 Return set of FABs bounding the domain grid boxes on a given orientation. More...
 
void setVal (value_type v)
 Set all boundary FABs to given value. More...
 
BndryRegisterT< MF > & operator+= (const BndryRegisterT< MF > &rhs)
 register += rhs More...
 
BndryRegisterT< MF > & plus (const BndryRegisterT< MF > &rhs)
 
BndryRegisterT< MF > & copyFrom (const MF &src, int nghost, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic())
 Fill the boundary FABs on intersection with given MF. More...
 
BndryRegisterT< MF > & plusFrom (const MF &src, int nghost, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic())
 Increment the boundary FABs on intersect with given MF. More...
 
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 intersection of MFs with the boundary FABs. More...
 
void setBoxes (const BoxArray &grids)
 Set box domain, if not set previously. More...
 
const DistributionMappingDistributionMap () const noexcept
 Returns constant reference to associated DistributionMapping. More...
 
void write (const std::string &name, std::ostream &os) const
 Write (used for writing to checkpoint) More...
 
void read (const std::string &name, std::istream &is)
 Read (used for reading from checkpoint) More...
 

Private Attributes

IntVect ratio
 Refinement ratio. More...
 
int fine_level
 Current level + 1. More...
 
int ncomp
 Number of state components. More...
 

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. More...
 
- Protected Attributes inherited from amrex::BndryRegisterT< MF >
FabSetT< MF > bndry [2 *AMREX_SPACEDIM]
 The data. More...
 
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 
)

The constructor. This version allows setting the DistributionMapping.

Parameters
fine_boxes
dm
ref_ratio
fine_lev
nvar

◆ ~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 ( )

◆ ClearInternalBorders()

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

Set internal borders to zero.

Parameters
crse_geom

◆ 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 flux register. This is different from CrseInit with FluxRegister::ADD. This is used for cases in which the grids covered by fine do not have fluxes computed (e.g., FLASH).

Parameters
mflx
area
dir
srccomp
destcomp
numcomp
mult
geom

◆ CrseAdd() [2/2]

void amrex::FluxRegister::CrseAdd ( const MultiFab mflx,
int  dir,
int  srccomp,
int  destcomp,
int  numcomp,
Real  mult,
const Geometry geom 
)

/in this version the area is assumed to multiplied into the flux (if not, use scale to fix)

Parameters
mflx
dir
srccomp
destcomp
numcomp
mult
geom

◆ 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
mflx
area
dir
srccomp
destcomp
numcomp
mult
op

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

Parameters
mflx
dir
srccomp
destcomp
numcomp
mult
op

◆ 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 using default constructor. This version allows setting the DistributionMapping.

Parameters
fine_boxes
dm
ref_ratio
fine_lev
nvar

◆ 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 data.

Parameters
flux
area
dir
boxno
srccomp
destcomp
numcomp
mult

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

in this version the area is assumed to multiplied into the flux (if not, use scale to fix)

Parameters
flux
dir
boxno
srccomp
destcomp
numcomp
mult

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

Parameters
mflx
area
dir
srccomp
destcomp
numcomp
mult

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

/in this version the area is assumed to multiplied into the flux (if not, use scale to fix)

Parameters
mflx
dir
srccomp
destcomp
numcomp
mult

◆ 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 a given value. This routine used by FLASH does NOT run on gpu for safety.

Parameters
dir
boxno
destcomp
numcomp
val

◆ 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 *, AMREX_SPACEDIM > 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.
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 (used for reading from checkpoint)

Parameters
name
is

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

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

◆ 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. Note that this takes the coarse Geometry.

Parameters
mf
volume
scale
srccomp
destcomp
numcomp
crse_geom

◆ Reflux() [4/5]

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

◆ 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(). Note that this takes the coarse Geometry.

Parameters
mf
scale
srccomp
destcomp
numcomp
crse_geom

◆ 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 (used for writing to checkpoint)

Parameters
name
os

Member Data Documentation

◆ fine_level

int amrex::FluxRegister::fine_level
private

Current level + 1.

◆ ncomp

int amrex::FluxRegister::ncomp
private

Number of state components.

◆ ratio

IntVect amrex::FluxRegister::ratio
private

Refinement ratio.


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