Block-Structured AMR Software Framework
amrex::FillPatcher< MF > Class Template Reference

FillPatcher is for filling a fine level MultiFab/FabArray. More...

#include <AMReX_FillPatcher.H>

Public Member Functions

 FillPatcher (BoxArray const &fba, DistributionMapping const &fdm, Geometry const &fgeom, BoxArray const &cba, DistributionMapping const &cdm, Geometry const &cgeom, IntVect const &nghost, int ncomp, InterpBase *interp, EB2::IndexSpace const *eb_index_space=nullptr)
 Constructor of FillPatcher. More...
 
template<typename BC , typename PreInterpHook = NullInterpHook<MF>, typename PostInterpHook = NullInterpHook<MF>>
void fill (MF &mf, IntVect const &nghost, Real time, Vector< MF * > const &cmf, Vector< Real > const &ct, Vector< MF * > const &fmf, Vector< Real > const &ft, int scomp, int dcomp, int ncomp, BC &cbc, int cbccomp, BC &fbc, int fbccomp, Vector< BCRec > const &bcs, int bcscomp, PreInterpHook const &pre_interp={}, PostInterpHook const &post_interp={})
 Function to fill data. More...
 
template<typename BC , typename PreInterpHook = NullInterpHook<MF>, typename PostInterpHook = NullInterpHook<MF>>
void fillCoarseFineBoundary (MF &mf, IntVect const &nghost, Real time, Vector< MF * > const &cmf, Vector< Real > const &ct, int scomp, int dcomp, int ncomp, BC &cbc, int cbccomp, Vector< BCRec > const &bcs, int bcscomp, PreInterpHook const &pre_interp={}, PostInterpHook const &post_interp={})
 Function to fill data at coarse/fine boundary only. More...
 
template<std::size_t order>
void storeRKCoarseData (Real time, Real dt, MF const &S_old, Array< MF, order > const &RK_k)
 Store coarse AMR level data for RK3 and RK4. More...
 
template<typename BC >
void fillRK (int stage, int iteration, int ncycle, MF &mf, Real time, BC &cbc, BC &fbc, Vector< BCRec > const &bcs)
 Fill ghost cells of fine AMR level for RK3 and RK4. More...
 

Private Member Functions

FabArrayBase::FPinfo const & getFPinfo ()
 

Private Attributes

BoxArray m_fba
 
BoxArray m_cba
 
DistributionMapping m_fdm
 
DistributionMapping m_cdm
 
Geometry m_fgeom
 
Geometry m_cgeom
 
IntVect m_nghost
 
int m_ncomp
 
InterpBasem_interp
 
EB2::IndexSpace const * m_eb_index_space = nullptr
 
MF m_sfine
 
IntVect m_ratio
 
Vector< std::pair< Real, std::unique_ptr< MF > > > m_cf_crse_data
 
std::unique_ptr< MF > m_cf_crse_data_tmp
 
std::unique_ptr< MF > m_cf_fine_data
 
Real m_dt_coarse = std::numeric_limits<Real>::lowest()
 

Detailed Description

template<class MF = MultiFab>
class amrex::FillPatcher< MF >

FillPatcher is for filling a fine level MultiFab/FabArray.

This class is not as general as the FillPatchTwoLevels functions. It fills the fine ghost cells not overlapping any fine level valid cells with interpolation of the coarse data. Then it fills the fine ghost cells overlapping fine level valid cells with the fine level data. If the valid cells of the destination need to be filled, it will be done as well. Finally, it will fill the physical boundary using the user provided functor. The fill member function can be used to do the operations just described. Alternatively, one can also use the fillCoarseFineBounary to fill the ghost cells at the coarse/fine boundary only. Then one can manually call FillBoundary to fill the other ghost cells, and use the physical BC functor to handle the physical boundary.

The communication of the coarse data needed for spatial interpolation is optimized at the cost of being error-prone. One must follow the following guidelines.

(1) This class is for filling data during time stepping, not during regrid. The fine level data passed as input must have the same BoxArray and DistributionMapping as the destination. It's OK they are the same MultiFab. For AmrLevel based codes, AmrLevel::FillPatcherFill will try to use FillPatcher if it can, and AmrLevel::FillPatch will use the fillpatch functions.

(2) When to build? It is recommended that one uses std::unique_ptr to store the FillPatcher object, and build it only when it is needed and it's a nullptr. For AmrLevel based codes, the AmrLevel class will build it for you as needed when you call the AmrLevel::FillPatcherFill function.

(3) When to destroy? Usually, we do time steppig on a coarse level first. Then we recursively do time stepping on fine levels. After the finer level finishes, we do reflux and average the fine data down to the coarse level. After that we should destroy the FillPatcher object associated with these two levels, because the coarse data stored in the object has become outdated. For AmrCore based codes, you could use Tests/Amr/Advection_AmrCore as an example. For AmrLevel based codes, you should do this in the post_timestep virtual function (see Tests/Amr/Advection_AmrLevel for an example).

(4) The source MultiFabs/FabArrays (i.e., the crse_data and fine_data arguments of the fill function) need to have exactly the same number of components as the ncomp argument of the constructor, even though it's allowed to fill only some of the components with the fill function.

(5) This only works for cell-centered and nodal data.

This class also provides support for RungeKutta::RK3 and RungeKutta::RK4. The storeRKCoarseData function can be used to store coarse AMR level data that are needed for filling fine level data's ghost cells in this class. The fillRK function can be used to fill ghost cells for fine AMR levels. This operation at the coarse/fine boundary is non-trivial for RK orders higher than 2. Note that it is expected that time stepping on the coarse level is perform before any fine level time stepping, and it's the user's responsibility to properly create and destroy this object. See AmrLevel::RK for an example of using the RungeKutta functions and FillPatcher together.

Constructor & Destructor Documentation

◆ FillPatcher()

template<class MF >
amrex::FillPatcher< MF >::FillPatcher ( BoxArray const &  fba,
DistributionMapping const &  fdm,
Geometry const &  fgeom,
BoxArray const &  cba,
DistributionMapping const &  cdm,
Geometry const &  cgeom,
IntVect const &  nghost,
int  ncomp,
InterpBase interp,
EB2::IndexSpace const *  eb_index_space = nullptr 
)

Constructor of FillPatcher.

Parameters
fbafine level BoxArray
fdmfine level DistributionMapping
fgeomfine level Geometry
cbacoarse level BoxArray
cdmcoarse level DistributionMapping
cgeomcoarse level Geometry
nghostmax number of ghost cells to be filled at coarse/fine boundary
ncompthe number of components
interpfor spatial interpolation
eb_index_spaceoptional argument for specifying EB IndexSpace

Member Function Documentation

◆ fill()

template<class MF >
template<typename BC , typename PreInterpHook , typename PostInterpHook >
void amrex::FillPatcher< MF >::fill ( MF &  mf,
IntVect const &  nghost,
Real  time,
Vector< MF * > const &  cmf,
Vector< Real > const &  ct,
Vector< MF * > const &  fmf,
Vector< Real > const &  ft,
int  scomp,
int  dcomp,
int  ncomp,
BC &  cbc,
int  cbccomp,
BC &  fbc,
int  fbccomp,
Vector< BCRec > const &  bcs,
int  bcscomp,
PreInterpHook const &  pre_interp = {},
PostInterpHook const &  post_interp = {} 
)

Function to fill data.

Parameters
mfdestination MultiFab/FabArray
nghostnumber of ghost cells to fill. This must be <= what's provided to the constructor
timetime associated with the destination
crse_datacoarse level data
crse_timetime associated with the coarse data
fine_datafine level data
fine_timetime associated with the fine data
scompstarting component of the source
dcompstarting component of the destination
ncompthe number of components to fill
cbcfor filling coarse level physical BC
cbccompstarting component of the coarse level BC functor
fbcfor filling fine level physical BC
fbccompstarting component of the fine level BC functor
bcsBCRec specifying physical boundary types
bcscompstarting component of the BCRec Vector.
pre_interpoptional pre-interpolation hook for modifying the coarse data
post_interpoptional post-interpolation hook for modifying the fine data

◆ fillCoarseFineBoundary()

template<class MF >
template<typename BC , typename PreInterpHook , typename PostInterpHook >
void amrex::FillPatcher< MF >::fillCoarseFineBoundary ( MF &  mf,
IntVect const &  nghost,
Real  time,
Vector< MF * > const &  cmf,
Vector< Real > const &  ct,
int  scomp,
int  dcomp,
int  ncomp,
BC &  cbc,
int  cbccomp,
Vector< BCRec > const &  bcs,
int  bcscomp,
PreInterpHook const &  pre_interp = {},
PostInterpHook const &  post_interp = {} 
)

Function to fill data at coarse/fine boundary only.

Parameters
mfdestination MultiFab/FabArray
nghostnumber of ghost cells to fill. This must be <= what's provided to the constructor
timetime associated with the destination
crse_datacoarse level data
crse_timetime associated with the coarse data
scompstarting component of the source
dcompstarting component of the destination
ncompthe number of components to fill
cbcfor filling coarse level physical BC
cbccompstarting component of the coarse level BC functor
bcsBCRec specifying physical boundary types
bcscompstarting component of the BCRec Vector.
pre_interpoptional pre-interpolation hook for modifying the coarse data
post_interpoptional post-interpolation hook for modifying the fine data

◆ fillRK()

template<typename MF >
template<typename BC >
void amrex::FillPatcher< MF >::fillRK ( int  stage,
int  iteration,
int  ncycle,
MF &  mf,
Real  time,
BC &  cbc,
BC &  fbc,
Vector< BCRec > const &  bcs 
)

Fill ghost cells of fine AMR level for RK3 and RK4.

Parameters
stageRK stage number starting from 1
iterationiteration number on fine level during a coarse time step. For an AMR simulation with subcycling and a refinement ratio of 2, the number is either 1 or 2, denoting the first and second substep, respectively.
ncyclenumber of subcyling steps. It's usually 2 or 4. Without subcycling, this will be 1.
cbcfilling physical boundary on coarse level
fbcfilling physical boundary on fine level
bcsphysical BC types

◆ getFPinfo()

template<class MF >
FabArrayBase::FPinfo const & amrex::FillPatcher< MF >::getFPinfo
private

◆ storeRKCoarseData()

template<typename MF >
template<std::size_t order>
void amrex::FillPatcher< MF >::storeRKCoarseData ( Real  time,
Real  dt,
MF const &  S_old,
Array< MF, order > const &  RK_k 
)

Store coarse AMR level data for RK3 and RK4.

Template Parameters
orderRK order. Must be 3 or 4.
Parameters
timetime at the beginning of the step
dttime step
S_olddata at time
RK_kright-hand side at RK stages

Member Data Documentation

◆ m_cba

template<class MF = MultiFab>
BoxArray amrex::FillPatcher< MF >::m_cba
private

◆ m_cdm

template<class MF = MultiFab>
DistributionMapping amrex::FillPatcher< MF >::m_cdm
private

◆ m_cf_crse_data

template<class MF = MultiFab>
Vector<std::pair<Real,std::unique_ptr<MF> > > amrex::FillPatcher< MF >::m_cf_crse_data
private

◆ m_cf_crse_data_tmp

template<class MF = MultiFab>
std::unique_ptr<MF> amrex::FillPatcher< MF >::m_cf_crse_data_tmp
private

◆ m_cf_fine_data

template<class MF = MultiFab>
std::unique_ptr<MF> amrex::FillPatcher< MF >::m_cf_fine_data
private

◆ m_cgeom

template<class MF = MultiFab>
Geometry amrex::FillPatcher< MF >::m_cgeom
private

◆ m_dt_coarse

template<class MF = MultiFab>
Real amrex::FillPatcher< MF >::m_dt_coarse = std::numeric_limits<Real>::lowest()
private

◆ m_eb_index_space

template<class MF = MultiFab>
EB2::IndexSpace const* amrex::FillPatcher< MF >::m_eb_index_space = nullptr
private

◆ m_fba

template<class MF = MultiFab>
BoxArray amrex::FillPatcher< MF >::m_fba
private

◆ m_fdm

template<class MF = MultiFab>
DistributionMapping amrex::FillPatcher< MF >::m_fdm
private

◆ m_fgeom

template<class MF = MultiFab>
Geometry amrex::FillPatcher< MF >::m_fgeom
private

◆ m_interp

template<class MF = MultiFab>
InterpBase* amrex::FillPatcher< MF >::m_interp
private

◆ m_ncomp

template<class MF = MultiFab>
int amrex::FillPatcher< MF >::m_ncomp
private

◆ m_nghost

template<class MF = MultiFab>
IntVect amrex::FillPatcher< MF >::m_nghost
private

◆ m_ratio

template<class MF = MultiFab>
IntVect amrex::FillPatcher< MF >::m_ratio
private

◆ m_sfine

template<class MF = MultiFab>
MF amrex::FillPatcher< MF >::m_sfine
private

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