Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::AmrLevel Class Referenceabstract

Virtual base class for managing individual levels. AmrLevel functions both as a container for state data on a level and also manages the advancement of data in time. More...

#include <AMReX_AmrLevel.H>

Public Types

enum  TimeLevel {
  AmrOldTime , AmrHalfTime , AmrNewTime , Amr1QtrTime ,
  Amr3QtrTime , AmrOtherTime
}
 What time are we at? More...
 

Public Member Functions

virtual ~AmrLevel ()
 The destructor.
 
 AmrLevel (const AmrLevel &)=delete
 Disable copying because AmrLevel owns state data and geometry caches.
 
 AmrLevel (AmrLevel &&)=delete
 Disable moving; derived classes manage non-movable FillPatch caches.
 
AmrLeveloperator= (const AmrLevel &)=delete
 Copy assignment is undefined for the same ownership reasons.
 
AmrLeveloperator= (AmrLevel &&)=delete
 Move assignment is also disabled—use factory methods instead.
 
void LevelDirectoryNames (const std::string &dir, std::string &LevelDir, std::string &FullPath) const
 Derive the directory name for this level inside a plot/checkpoint tree.
 
virtual void CreateLevelDirectory (const std::string &dir)
 Create the Level_* directory inside dir (if needed) before writing plots/checkpoints.
 
void SetLevelDirectoryCreated (bool ldc) noexcept
 Record whether CreateLevelDirectory() has already created the folder.
 
virtual std::string thePlotFileType () const
 A string written as the first item in writePlotFile() at level zero. It is so we can distinguish between different types of plot files. This default "HyperCLaw-V1.1" is for VisIt software and some of our internal postprocessing routines.
 
virtual void writePlotFile (const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles)
 Write plot file data for this level into directory dir.
 
virtual void writePlotFilePre (const std::string &dir, std::ostream &os)
 Perform pre-plotfile operations (e.g., gather metadata) for directory dir using os.
 
virtual void writePlotFilePost (const std::string &dir, std::ostream &os)
 Perform post-plotfile cleanup for directory dir once writes complete.
 
virtual void writeSmallPlotFile (const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles)
 Write a reduced plot file (same parameters as writePlotFile()) to dir.
 
virtual void checkPoint (const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles, bool dump_old=true)
 Write current state to a checkpoint directory.
 
virtual void checkPointPre (const std::string &dir, std::ostream &os)
 Do pre-checkpoint work (e.g., directory creation) for dir using os.
 
virtual void checkPointPost (const std::string &dir, std::ostream &os)
 Do post-checkpoint cleanup for dir once writes complete.
 
virtual void restart (Amr &papa, std::istream &is, bool bReadSpecial=false)
 Populate this level from checkpoint input stream is.
 
virtual void set_state_in_checkpoint (Vector< int > &state_in_checkpoint)
 Hook to adjust which state descriptors are expected in checkpoint streams.
 
virtual void computeInitialDt (int finest_level, int sub_cycle, Vector< int > &n_cycle, const Vector< IntVect > &ref_ratio, Vector< Real > &dt_level, Real stop_time)=0
 Compute the initial time step size for each level.
 
virtual void computeNewDt (int finest_level, int sub_cycle, Vector< int > &n_cycle, const Vector< IntVect > &ref_ratio, Vector< Real > &dt_min, Vector< Real > &dt_level, Real stop_time, int post_regrid_flag)=0
 Compute the next time step sizes after an advance.
 
virtual Real advance (Real time, Real dt, int iteration, int ncycle)=0
 Advance solution on this level by dt and return the next stable dt.
 
virtual void post_timestep (int iteration)
 Perform level-specific updates after each subcycle.
 
virtual void postCoarseTimeStep (Real time)
 Called once per coarse step after fine levels finish refluxing.
 
virtual void post_restart ()
 Override to perform any fixups immediately after restart (default: no-op).
 
virtual void post_regrid (int lbase, int new_finest)=0
 Operations to be performed immediately after regridding completes.
 
virtual void post_init (Real stop_time)=0
 Operations to run after initialization finishes.
 
virtual int okToContinue ()
 Return nonzero to continue advancing; zero to halt the run.
 
virtual int okToRegrid ()
 Return nonzero if refinement starting from this level is allowed this cycle.
 
virtual void initData ()=0
 Init grid data at problem start-up. This is a pure virtual function and hence MUST be implemented by derived classes.
 
virtual void setTimeLevel (Real time, Real dt_old, Real dt_new)
 Set the time metadata stored in the state data arrays (old/new).
 
virtual void allocOldData ()
 Allocate storage for "old" time slices before an advance.
 
virtual void removeOldData ()
 Delete the old-time slots once they are no longer needed.
 
virtual void init (AmrLevel &old)=0
 Init data on this level from another AmrLevel (during regrid). This is a pure virtual function and hence MUST be implemented by derived classes.
 
virtual void init ()=0
 
void reset ()
 Reset data to initial time by swapping new and old time slots.
 
int Level () const noexcept
 Returns this AmrLevel index (0 = coarsest).
 
const BoxArrayboxArray () const noexcept
 List of grids at this level.
 
const BoxArraygetEdgeBoxArray (int dir) const noexcept
 Return the edge-centered grid layout for direction dir (cached).
 
const BoxArraygetNodalBoxArray () const noexcept
 Return the nodal grid layout for this level (cached).
 
const DistributionMappingDistributionMap () const noexcept
 
const FabFactory< FArrayBox > & Factory () const noexcept
 
const EBFArrayBoxFactoryEBFactory () const noexcept
 
int numGrids () const noexcept
 Number of grids at this level.
 
int numStates () const noexcept
 Number of states at this level.
 
const BoxDomain () const noexcept
 Returns the indices defining physical domain.
 
int nStep () const noexcept
 Timestep n at this level.
 
const GeometryGeom () const noexcept
 Returns the geometry object.
 
const IntVectfineRatio () const noexcept
 
Long countCells () const noexcept
 Returns number of valid cells on this level.
 
const BoxArraygetAreaNotToTag () noexcept
 Get the cached BoxArray describing regions that should not be tagged.
 
const BoxgetAreaToTag () noexcept
 Get the Box (single region) that is permitted to be tagged.
 
void constructAreaNotToTag ()
 Build internal data structures defining the area not to tag.
 
void setAreaNotToTag (BoxArray &ba) noexcept
 Override the area-not-to-tag cache with user-supplied BoxArray ba.
 
void resetFillPatcher ()
 Clear cached FillPatcher objects (call after reflux/averaging).
 
virtual void errorEst (TagBoxArray &tb, int clearval, int tagval, Real time, int n_error_buf=0, int ngrow=0)=0
 Error estimation for regridding. This is a pure virtual function and hence MUST be implemented by derived classes.
 
void FillCoarsePatch (MultiFab &mf, int dcomp, Real time, int state_idx, int scomp, int ncomp, int nghost=0)
 Interpolate from coarse level to the valid area in mf.
 
virtual void setPhysBoundaryValues (FArrayBox &dest, int state_indx, Real time, int dest_comp, int src_comp, int num_comp)
 Apply physical boundary conditions to cell-centered data.
 
virtual std::unique_ptr< MultiFabderive (const std::string &name, Real time, int ngrow)
 Returns a MultiFab containing the derived data for this level. If ngrow>0 the MultiFab is built on the appropriately grown BoxArray.
 
virtual void derive (const std::string &name, Real time, MultiFab &mf, int dcomp)
 This version of derive() fills the dcomp'th component of mf with the derived quantity.
 
StateDataget_state_data (int state_indx) noexcept
 State data object accessor for descriptor state_indx.
 
MultiFabget_old_data (int state_indx) noexcept
 Mutable access to the MultiFab at the old time slot for state state_indx.
 
const MultiFabget_old_data (int state_indx) const noexcept
 Const access to the MultiFab at the old time slot for state state_indx.
 
MultiFabget_new_data (int state_indx) noexcept
 Mutable access to the MultiFab at the new time slot for state state_indx.
 
const MultiFabget_new_data (int state_indx) const noexcept
 Const access to the MultiFab at the new time slot for state state_indx.
 
int postStepRegrid () const noexcept
 Returns whether a post-timestep regrid is currently requested (nonzero for true).
 
void setPostStepRegrid (int new_val) noexcept
 Set the post-timestep regrid trigger (nonzero requests a regrid).
 
void UpdateDistributionMaps (DistributionMapping &dmap)
 Update every StateData's DistributionMapping to match dmap (e.g., after rebalancing).
 
Vector< intgetBCArray (int State_Type, int gridno, int strt_comp, int ncomp)
 Return an array describing boundary types for components [strt_comp,strt_comp+ncomp).
 
MultiFabget_data (int state_indx, Real time) noexcept
 Fetch the MultiFab (old or new) that best matches time for state state_indx.
 
virtual void set_preferred_boundary_values (MultiFab &S, int state_index, int scomp, int dcomp, int ncomp, Real time) const
 Hook to override fillpatched boundary values before applying physical BCs (default no-op).
 
virtual void manual_tags_placement (TagBoxArray &tags, const Vector< IntVect > &bf_lev)
 Adjust tags after core tagging routines run (default: no-op).
 
virtual void setPlotVariables ()
 Override to customize which variables appear in full plotfiles on this level (default uses global list).
 
virtual void setSmallPlotVariables ()
 Override to customize which variables appear in small plotfiles on this level.
 
virtual Real estimateWork ()
 Estimate the amount of work required to advance just this level, typically based on cell counts.
 
virtual int WorkEstType ()
 Select which state descriptor drives work estimates.
 
TimeLevel which_time (int state_indx, Real time) const noexcept
 Returns one the TimeLevel enums. Asserts that time is between AmrOldTime and AmrNewTime.
 
virtual bool checkPointNow ()
 Return true if this level requests an immediate checkpoint (default false).
 
virtual bool writePlotNow ()
 Return true if this level requests an immediate plotfile (default false).
 
virtual bool writeSmallPlotNow ()
 Return true if this level requests an immediate small-plot file (default false).
 
virtual void particle_redistribute (int lbase=0, bool a_init=false)
 Optional particle redistribution hook called by the parent Amr object.
 
void FillPatcherFill (amrex::MultiFab &mf, int dcomp, int ncomp, int nghost, amrex::Real time, int state_index, int scomp)
 Fill ghost cells using cached FillPatcher (level>0) or FillPatch (level 0).
 
template<typename F , typename P = RungeKutta::PostStageNoOp>
void RK (int order, int state_type, Real time, Real dt, int iteration, int ncycle, F &&f, P &&p=RungeKutta::PostStageNoOp())
 Evolve one step with Runge-Kutta (2, 3, or 4)
 

Static Public Member Functions

static bool isStateVariable (const std::string &name, int &state_indx, int &ncomp)
 Return true if name matches a registered state variable, also providing descriptor info.
 
static void FlushFPICache ()
 Resize/clear cached FillPatchIterator state in preparation for new fill operations.
 
static const DescriptorListget_desc_lst () noexcept
 Returns list of Descriptors.
 
static DeriveListget_derive_lst () noexcept
 Returns list of derived variables.
 
static void FillPatch (AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
 Convenience wrapper around FillPatchSingleLevel/TwoLevels.
 
static void FillPatchAdd (AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
 Same as FillPatch(), but adds into leveldata instead of overwriting.
 
static void SetEBMaxGrowCells (int nbasic, int nvolume, int nfull) noexcept
 Configure the ghost-cell requirements (basic/volume/full) for EB FillPatch operations.
 
static void SetEBSupportLevel (EBSupport ebs)
 Select which level of EB data to build (basic, volume, full) via ebs.
 
static IntVect ProperBlockingFactor (AmrLevel const &amr_level, int boxGrow, IndexType const &boxType, StateDescriptor const &desc, int SComp)
 Recommendation of a proper blocking factor for filling.
 

Static Public Attributes

static int m_eb_basic_grow_cells = 5
 
static int m_eb_volume_grow_cells = 4
 
static int m_eb_full_grow_cells = 2
 
static EBSupport m_eb_support_level = EBSupport::volume
 

Protected Member Functions

 AmrLevel () noexcept
 The constructors – for derived classes.
 
 AmrLevel (Amr &papa, int lev, const Geometry &level_geom, const BoxArray &ba, const DistributionMapping &dm, Real time)
 
void finishConstructor ()
 Common code used by all constructors.
 

Protected Attributes

int level {-1}
 
Geometry geom
 
BoxArray grids
 
DistributionMapping dmap
 
Amrparent {nullptr}
 
IntVect crse_ratio
 
IntVect fine_ratio
 
Vector< StateDatastate
 
BoxArray m_AreaNotToTag
 
Box m_AreaToTag
 
int post_step_regrid {0}
 
bool levelDirectoryCreated {false}
 
std::unique_ptr< FabFactory< FArrayBox > > m_factory
 
Vector< std::unique_ptr< FillPatcher< MultiFab > > > m_fillpatcher
 

Static Protected Attributes

static DeriveList derive_lst
 
static DescriptorList desc_lst
 

Friends

class Amr
 
class FillPatchIterator
 
class FillPatchIteratorHelper
 

Detailed Description

Virtual base class for managing individual levels. AmrLevel functions both as a container for state data on a level and also manages the advancement of data in time.

Member Enumeration Documentation

◆ TimeLevel

What time are we at?

Enumerator
AmrOldTime 
AmrHalfTime 
AmrNewTime 
Amr1QtrTime 
Amr3QtrTime 
AmrOtherTime 

Constructor & Destructor Documentation

◆ ~AmrLevel()

amrex::AmrLevel::~AmrLevel ( )
virtual

The destructor.

◆ AmrLevel() [1/4]

amrex::AmrLevel::AmrLevel ( const AmrLevel )
delete

Disable copying because AmrLevel owns state data and geometry caches.

◆ AmrLevel() [2/4]

amrex::AmrLevel::AmrLevel ( AmrLevel &&  )
delete

Disable moving; derived classes manage non-movable FillPatch caches.

◆ AmrLevel() [3/4]

amrex::AmrLevel::AmrLevel ( )
inlineprotectednoexcept

The constructors – for derived classes.

◆ AmrLevel() [4/4]

amrex::AmrLevel::AmrLevel ( Amr papa,
int  lev,
const Geometry level_geom,
const BoxArray ba,
const DistributionMapping dm,
Real  time 
)
protected

Member Function Documentation

◆ advance()

virtual Real amrex::AmrLevel::advance ( Real  time,
Real  dt,
int  iteration,
int  ncycle 
)
pure virtual

Advance solution on this level by dt and return the next stable dt.

Parameters
timeStart time of the advance.
dtStep size to take.
iterationSub-step counter within the coarse step.
ncycleNumber of fine subcycles per coarse step.
Returns
Recommended dt for the next call (based on CFL, etc.).

◆ allocOldData()

void amrex::AmrLevel::allocOldData ( )
virtual

Allocate storage for "old" time slices before an advance.

◆ boxArray()

const BoxArray & amrex::AmrLevel::boxArray ( ) const
inlinenoexcept

List of grids at this level.

◆ checkPoint()

void amrex::AmrLevel::checkPoint ( const std::string &  dir,
std::ostream &  os,
VisMF::How  how = VisMF::NFiles,
bool  dump_old = true 
)
virtual

Write current state to a checkpoint directory.

Parameters
dirDestination directory (e.g., chk00010).
osMetadata stream shared across levels.
howVisMF write mode.
dump_oldWhether to include the "old" state data.

◆ checkPointNow()

bool amrex::AmrLevel::checkPointNow ( )
virtual

Return true if this level requests an immediate checkpoint (default false).

Overrides can trigger checkpointing after conditional events.

◆ checkPointPost()

void amrex::AmrLevel::checkPointPost ( const std::string &  dir,
std::ostream &  os 
)
virtual

Do post-checkpoint cleanup for dir once writes complete.

Parameters
dirCheckpoint directory that was just written.
osMetadata stream used during writing.

◆ checkPointPre()

void amrex::AmrLevel::checkPointPre ( const std::string &  dir,
std::ostream &  os 
)
virtual

Do pre-checkpoint work (e.g., directory creation) for dir using os.

Parameters
dirCheckpoint directory being created.
osMetadata stream for header output.

◆ computeInitialDt()

virtual void amrex::AmrLevel::computeInitialDt ( int  finest_level,
int  sub_cycle,
Vector< int > &  n_cycle,
const Vector< IntVect > &  ref_ratio,
Vector< Real > &  dt_level,
Real  stop_time 
)
pure virtual

Compute the initial time step size for each level.

Parameters
finest_levelFinest AMR level existing so far.
sub_cycleWhether subcycling is enabled.
n_cycleIn/out: number of subcycles per level.
ref_ratioLevel-to-level refinement ratios.
dt_levelOutput per-level dt.
stop_timeSimulation stop time (used to clamp dt if needed).

◆ computeNewDt()

virtual void amrex::AmrLevel::computeNewDt ( int  finest_level,
int  sub_cycle,
Vector< int > &  n_cycle,
const Vector< IntVect > &  ref_ratio,
Vector< Real > &  dt_min,
Vector< Real > &  dt_level,
Real  stop_time,
int  post_regrid_flag 
)
pure virtual

Compute the next time step sizes after an advance.

Parameters
finest_levelFinest AMR level existing.
sub_cycleWhether subcycling is enabled.
n_cycleIn/out: number of subcycles per level.
ref_ratioLevel-to-level refinement ratios.
dt_minOutput per-level dt upper bounds (physics-limited).
dt_levelOutput per-level dt to use.
stop_timeSimulation stop time (used to clamp dt).
post_regrid_flagFlag describing regrid state (0 normal, 1 immediate post-regrid).

◆ constructAreaNotToTag()

void amrex::AmrLevel::constructAreaNotToTag ( )

Build internal data structures defining the area not to tag.

◆ countCells()

Long amrex::AmrLevel::countCells ( ) const
noexcept

Returns number of valid cells on this level.

◆ CreateLevelDirectory()

void amrex::AmrLevel::CreateLevelDirectory ( const std::string &  dir)
virtual

Create the Level_* directory inside dir (if needed) before writing plots/checkpoints.

Parameters
dirBase output directory (plot or checkpoint root).

◆ derive() [1/2]

std::unique_ptr< MultiFab > amrex::AmrLevel::derive ( const std::string &  name,
Real  time,
int  ngrow 
)
virtual

Returns a MultiFab containing the derived data for this level. If ngrow>0 the MultiFab is built on the appropriately grown BoxArray.

Parameters
nameDerived quantity identifier.
timePhysical time at which to evaluate the quantity.
ngrowHow many ghost cells to include.
Returns
Owning pointer to a freshly allocated MultiFab.

◆ derive() [2/2]

void amrex::AmrLevel::derive ( const std::string &  name,
Real  time,
MultiFab mf,
int  dcomp 
)
virtual

This version of derive() fills the dcomp'th component of mf with the derived quantity.

Parameters
nameDerived quantity identifier.
timePhysical time at which to evaluate the quantity.
mfDestination MultiFab.
dcompDestination component index.

◆ DistributionMap()

const DistributionMapping & amrex::AmrLevel::DistributionMap ( ) const
inlinenoexcept

◆ Domain()

const Box & amrex::AmrLevel::Domain ( ) const
inlinenoexcept

Returns the indices defining physical domain.

◆ EBFactory()

const EBFArrayBoxFactory & amrex::AmrLevel::EBFactory ( ) const
inlinenoexcept

◆ errorEst()

virtual void amrex::AmrLevel::errorEst ( TagBoxArray tb,
int  clearval,
int  tagval,
Real  time,
int  n_error_buf = 0,
int  ngrow = 0 
)
pure virtual

Error estimation for regridding. This is a pure virtual function and hence MUST be implemented by derived classes.

Parameters
tbTag array to mark for refinement (output).
clearvalValue that clears a tag entry.
tagvalValue that marks a cell for refinement.
timePhysical time at which tagging occurs.
n_error_bufNumber of grow cells added around tagged features (default 0).
ngrowNumber of ghost layers available in tb (default 0).

◆ estimateWork()

Real amrex::AmrLevel::estimateWork ( )
virtual

Estimate the amount of work required to advance just this level, typically based on cell counts.

Returns
Work units proportional to the expected compute effort (higher means more expensive).

◆ Factory()

const FabFactory< FArrayBox > & amrex::AmrLevel::Factory ( ) const
inlinenoexcept

◆ FillCoarsePatch()

void amrex::AmrLevel::FillCoarsePatch ( MultiFab mf,
int  dcomp,
Real  time,
int  state_idx,
int  scomp,
int  ncomp,
int  nghost = 0 
)

Interpolate from coarse level to the valid area in mf.

Fill mf on this level using data interpolated from the coarser level.

Parameters
mfDestination MultiFab.
dcompDestination component offset in mf.
timeTime at which to sample coarse data.
state_idxState descriptor index selecting which state to fill.
scompSource component index within the coarse state.
ncompNumber of components to fill.
nghostNumber of ghost cells to fill (default 0).

◆ FillPatch()

void amrex::AmrLevel::FillPatch ( AmrLevel amrlevel,
MultiFab leveldata,
int  boxGrow,
Real  time,
int  index,
int  scomp,
int  ncomp,
int  dcomp = 0 
)
static

Convenience wrapper around FillPatchSingleLevel/TwoLevels.

Parameters
amrlevelLevel providing state data.
leveldataDestination MultiFab to fill.
boxGrowNumber of grow cells to request.
timeTime at which to sample.
indexState descriptor index.
scompSource component.
ncompNumber of components.
dcompDestination component (default 0).

◆ FillPatchAdd()

void amrex::AmrLevel::FillPatchAdd ( AmrLevel amrlevel,
MultiFab leveldata,
int  boxGrow,
Real  time,
int  index,
int  scomp,
int  ncomp,
int  dcomp = 0 
)
static

Same as FillPatch(), but adds into leveldata instead of overwriting.

Parameters
amrlevelLevel providing state data.
leveldataDestination MultiFab to accumulate into.
boxGrowNumber of grow cells to request.
timeTime at which to sample.
indexState descriptor index.
scompSource component.
ncompNumber of components.
dcompDestination component (default 0).

◆ FillPatcherFill()

void amrex::AmrLevel::FillPatcherFill ( amrex::MultiFab mf,
int  dcomp,
int  ncomp,
int  nghost,
amrex::Real  time,
int  state_index,
int  scomp 
)

Fill ghost cells using cached FillPatcher (level>0) or FillPatch (level 0).

Parameters
mfDestination MultiFab.
dcompDestination component offset.
ncompNumber of components to fill.
nghostNumber of ghost cells to fill.
timeTime at which to sample data.
state_indexState descriptor index to pull from.
scompSource component offset within the state.

◆ fineRatio()

const IntVect & amrex::AmrLevel::fineRatio ( ) const
inlinenoexcept

◆ finishConstructor()

void amrex::AmrLevel::finishConstructor ( )
protected

Common code used by all constructors.

◆ FlushFPICache()

static void amrex::AmrLevel::FlushFPICache ( )
static

Resize/clear cached FillPatchIterator state in preparation for new fill operations.

◆ Geom()

const Geometry & amrex::AmrLevel::Geom ( ) const
inlinenoexcept

Returns the geometry object.

◆ get_data()

MultiFab & amrex::AmrLevel::get_data ( int  state_indx,
Real  time 
)
noexcept

Fetch the MultiFab (old or new) that best matches time for state state_indx.

Parameters
state_indxState descriptor index.
timeTarget physical time.
Returns
Reference to either the old or new MultiFab, whichever brackets time.

◆ get_derive_lst()

DeriveList & amrex::AmrLevel::get_derive_lst ( )
staticnoexcept

Returns list of derived variables.

◆ get_desc_lst()

static const DescriptorList & amrex::AmrLevel::get_desc_lst ( )
inlinestaticnoexcept

Returns list of Descriptors.

◆ get_new_data() [1/2]

const MultiFab & amrex::AmrLevel::get_new_data ( int  state_indx) const
inlinenoexcept

Const access to the MultiFab at the new time slot for state state_indx.

◆ get_new_data() [2/2]

MultiFab & amrex::AmrLevel::get_new_data ( int  state_indx)
inlinenoexcept

Mutable access to the MultiFab at the new time slot for state state_indx.

◆ get_old_data() [1/2]

const MultiFab & amrex::AmrLevel::get_old_data ( int  state_indx) const
inlinenoexcept

Const access to the MultiFab at the old time slot for state state_indx.

◆ get_old_data() [2/2]

MultiFab & amrex::AmrLevel::get_old_data ( int  state_indx)
inlinenoexcept

Mutable access to the MultiFab at the old time slot for state state_indx.

◆ get_state_data()

StateData & amrex::AmrLevel::get_state_data ( int  state_indx)
inlinenoexcept

State data object accessor for descriptor state_indx.

◆ getAreaNotToTag()

const BoxArray & amrex::AmrLevel::getAreaNotToTag ( )
noexcept

Get the cached BoxArray describing regions that should not be tagged.

◆ getAreaToTag()

const Box & amrex::AmrLevel::getAreaToTag ( )
noexcept

Get the Box (single region) that is permitted to be tagged.

◆ getBCArray()

Vector< int > amrex::AmrLevel::getBCArray ( int  State_Type,
int  gridno,
int  strt_comp,
int  ncomp 
)

Return an array describing boundary types for components [strt_comp,strt_comp+ncomp).

Parameters
State_TypeState descriptor index.
gridnoGrid index on this level.
strt_compStarting component.
ncompNumber of components.

◆ getEdgeBoxArray()

const BoxArray & amrex::AmrLevel::getEdgeBoxArray ( int  dir) const
noexcept

Return the edge-centered grid layout for direction dir (cached).

◆ getNodalBoxArray()

const BoxArray & amrex::AmrLevel::getNodalBoxArray ( ) const
noexcept

Return the nodal grid layout for this level (cached).

◆ init() [1/2]

virtual void amrex::AmrLevel::init ( )
pure virtual

Init data on this level after regridding if old AmrLevel did not previously exist. This is a pure virtual function and hence MUST be implemented by derived classes.

◆ init() [2/2]

virtual void amrex::AmrLevel::init ( AmrLevel old)
pure virtual

Init data on this level from another AmrLevel (during regrid). This is a pure virtual function and hence MUST be implemented by derived classes.

Parameters
oldSource level to interpolate from.

◆ initData()

virtual void amrex::AmrLevel::initData ( )
pure virtual

Init grid data at problem start-up. This is a pure virtual function and hence MUST be implemented by derived classes.

◆ isStateVariable()

bool amrex::AmrLevel::isStateVariable ( const std::string &  name,
int state_indx,
int ncomp 
)
static

Return true if name matches a registered state variable, also providing descriptor info.

Parameters
nameState variable identifier.
state_indxOutput: descriptor index when found.
ncompOutput: number of components associated with that state.

◆ Level()

int amrex::AmrLevel::Level ( ) const
inlinenoexcept

Returns this AmrLevel index (0 = coarsest).

◆ LevelDirectoryNames()

void amrex::AmrLevel::LevelDirectoryNames ( const std::string &  dir,
std::string &  LevelDir,
std::string &  FullPath 
) const

Derive the directory name for this level inside a plot/checkpoint tree.

Parameters
dirRoot directory (typically plot or checkpoint folder).
LevelDirOutput: relative level directory (e.g., Level_3).
FullPathOutput: absolute path combining dir and LevelDir.

◆ manual_tags_placement()

void amrex::AmrLevel::manual_tags_placement ( TagBoxArray tags,
const Vector< IntVect > &  bf_lev 
)
virtual

Adjust tags after core tagging routines run (default: no-op).

Parameters
tagsTag arrays built on this level.
bf_levBuffer sizes per level from Amr.

◆ nStep()

int amrex::AmrLevel::nStep ( ) const
inlinenoexcept

Timestep n at this level.

◆ numGrids()

int amrex::AmrLevel::numGrids ( ) const
inlinenoexcept

Number of grids at this level.

◆ numStates()

int amrex::AmrLevel::numStates ( ) const
inlinenoexcept

Number of states at this level.

◆ okToContinue()

virtual int amrex::AmrLevel::okToContinue ( )
inlinevirtual

Return nonzero to continue advancing; zero to halt the run.

◆ okToRegrid()

int amrex::AmrLevel::okToRegrid ( )
virtual

Return nonzero if refinement starting from this level is allowed this cycle.

Only evaluated when regrid_int > 0 and level_count >= regrid_int. Default returns true.

◆ operator=() [1/2]

AmrLevel & amrex::AmrLevel::operator= ( AmrLevel &&  )
delete

Move assignment is also disabled—use factory methods instead.

◆ operator=() [2/2]

AmrLevel & amrex::AmrLevel::operator= ( const AmrLevel )
delete

Copy assignment is undefined for the same ownership reasons.

◆ particle_redistribute()

virtual void amrex::AmrLevel::particle_redistribute ( int  lbase = 0,
bool  a_init = false 
)
inlinevirtual

Optional particle redistribution hook called by the parent Amr object.

Parameters
lbaseCoarsest level whose particles should be considered.
a_initWhether the call happens during initialization.

◆ post_init()

virtual void amrex::AmrLevel::post_init ( Real  stop_time)
pure virtual

Operations to run after initialization finishes.

Parameters
stop_timeSimulation stop time (useful for multi-stage inits).

◆ post_regrid()

virtual void amrex::AmrLevel::post_regrid ( int  lbase,
int  new_finest 
)
pure virtual

Operations to be performed immediately after regridding completes.

Parameters
lbaseCoarsest level that triggered the regrid.
new_finestNew finest level index after regridding.

◆ post_restart()

virtual void amrex::AmrLevel::post_restart ( )
inlinevirtual

Override to perform any fixups immediately after restart (default: no-op).

◆ post_timestep()

void amrex::AmrLevel::post_timestep ( int  iteration)
virtual

Perform level-specific updates after each subcycle.

Parameters
iterationSub-step index within the coarse step.
Note
Override implementations should reset any FillPatcher caches they own.

◆ postCoarseTimeStep()

void amrex::AmrLevel::postCoarseTimeStep ( Real  time)
virtual

Called once per coarse step after fine levels finish refluxing.

Parameters
timeCurrent coarse-level time.

◆ postStepRegrid()

int amrex::AmrLevel::postStepRegrid ( ) const
inlinenoexcept

Returns whether a post-timestep regrid is currently requested (nonzero for true).

◆ ProperBlockingFactor()

IntVect amrex::AmrLevel::ProperBlockingFactor ( AmrLevel const &  amr_level,
int  boxGrow,
IndexType const &  boxType,
StateDescriptor const &  desc,
int  SComp 
)
static

Recommendation of a proper blocking factor for filling.

Parameters
amr_levelOwning level supplying geometry.
boxGrowGrow cells requested by the caller.
boxTypeIndex type (cell-centered, nodal, etc.).
descState descriptor tied to the data.
SCompStarting component index.

◆ removeOldData()

void amrex::AmrLevel::removeOldData ( )
virtual

Delete the old-time slots once they are no longer needed.

◆ reset()

void amrex::AmrLevel::reset ( )

Reset data to initial time by swapping new and old time slots.

◆ resetFillPatcher()

void amrex::AmrLevel::resetFillPatcher ( )

Clear cached FillPatcher objects (call after reflux/averaging).

◆ restart()

void amrex::AmrLevel::restart ( Amr papa,
std::istream &  is,
bool  bReadSpecial = false 
)
virtual

Populate this level from checkpoint input stream is.

Parameters
papaOwning Amr driver (provides geometry and config).
isInput stream pointing at the level's checkpoint data.
bReadSpecialWhether to parse additional application-specific fields.

◆ RK()

template<typename F , typename P >
void amrex::AmrLevel::RK ( int  order,
int  state_type,
Real  time,
Real  dt,
int  iteration,
int  ncycle,
F &&  f,
P &&  p = RungeKutta::PostStageNoOp() 
)

Evolve one step with Runge-Kutta (2, 3, or 4)

To use RK, the StateData must have all the ghost cells needed. See namespace RungeKutta for expected function signatures of the callable parameters.

Parameters
orderorder of RK
state_typeindex of StateData
timetime at the beginning of the step.
dttime step
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.
fcomputing right-hand side for evolving the StateData. One can also register data for flux registers in this.
poptionally post-processing RK stage results

◆ set_preferred_boundary_values()

void amrex::AmrLevel::set_preferred_boundary_values ( MultiFab S,
int  state_index,
int  scomp,
int  dcomp,
int  ncomp,
Real  time 
) const
virtual

Hook to override fillpatched boundary values before applying physical BCs (default no-op).

Parameters
SMultiFab to edit.
state_indexState descriptor index.
scompSource component offset.
dcompDestination component offset.
ncompNumber of components to adjust.
timePhysical time at which BCs are evaluated.

◆ set_state_in_checkpoint()

void amrex::AmrLevel::set_state_in_checkpoint ( Vector< int > &  state_in_checkpoint)
virtual

Hook to adjust which state descriptors are expected in checkpoint streams.

Parameters
state_in_checkpointArray of flags (one per descriptor) that this level may modify.

◆ setAreaNotToTag()

void amrex::AmrLevel::setAreaNotToTag ( BoxArray ba)
noexcept

Override the area-not-to-tag cache with user-supplied BoxArray ba.

◆ SetEBMaxGrowCells()

static void amrex::AmrLevel::SetEBMaxGrowCells ( int  nbasic,
int  nvolume,
int  nfull 
)
inlinestaticnoexcept

Configure the ghost-cell requirements (basic/volume/full) for EB FillPatch operations.

◆ SetEBSupportLevel()

static void amrex::AmrLevel::SetEBSupportLevel ( EBSupport  ebs)
inlinestatic

Select which level of EB data to build (basic, volume, full) via ebs.

◆ SetLevelDirectoryCreated()

void amrex::AmrLevel::SetLevelDirectoryCreated ( bool  ldc)
inlinenoexcept

Record whether CreateLevelDirectory() has already created the folder.

Parameters
ldcTrue if the level directory exists; false to clear the flag.

◆ setPhysBoundaryValues()

void amrex::AmrLevel::setPhysBoundaryValues ( FArrayBox dest,
int  state_indx,
Real  time,
int  dest_comp,
int  src_comp,
int  num_comp 
)
virtual

Apply physical boundary conditions to cell-centered data.

Parameters
destDestination FArrayBox.
state_indxState descriptor index associated with dest.
timePhysical time at which BCs should be evaluated.
dest_compDestination component offset in dest.
src_compComponent offset in the state data.
num_compNumber of components to update.

◆ setPlotVariables()

void amrex::AmrLevel::setPlotVariables ( )
virtual

Override to customize which variables appear in full plotfiles on this level (default uses global list).

◆ setPostStepRegrid()

void amrex::AmrLevel::setPostStepRegrid ( int  new_val)
inlinenoexcept

Set the post-timestep regrid trigger (nonzero requests a regrid).

Parameters
new_valTrigger value applied after each step.

◆ setSmallPlotVariables()

void amrex::AmrLevel::setSmallPlotVariables ( )
virtual

Override to customize which variables appear in small plotfiles on this level.

◆ setTimeLevel()

void amrex::AmrLevel::setTimeLevel ( Real  time,
Real  dt_old,
Real  dt_new 
)
virtual

Set the time metadata stored in the state data arrays (old/new).

Parameters
timeNew-time value to record.
dt_oldDt associated with the old slot.
dt_newDt associated with the new slot.

◆ thePlotFileType()

virtual std::string amrex::AmrLevel::thePlotFileType ( ) const
inlinevirtual

A string written as the first item in writePlotFile() at level zero. It is so we can distinguish between different types of plot files. This default "HyperCLaw-V1.1" is for VisIt software and some of our internal postprocessing routines.

< Default identifier consumed by VisIt and AMReX tools.

◆ UpdateDistributionMaps()

void amrex::AmrLevel::UpdateDistributionMaps ( DistributionMapping dmap)

Update every StateData's DistributionMapping to match dmap (e.g., after rebalancing).

Update the distribution maps in StateData based on the size of the map.

Parameters
dmapDistribution mapping describing the new ownership.

◆ which_time()

AmrLevel::TimeLevel amrex::AmrLevel::which_time ( int  state_indx,
Real  time 
) const
noexcept

Returns one the TimeLevel enums. Asserts that time is between AmrOldTime and AmrNewTime.

Parameters
state_indxState descriptor index providing time metadata.
timeTarget physical time.

◆ WorkEstType()

virtual int amrex::AmrLevel::WorkEstType ( )
inlinevirtual

Select which state descriptor drives work estimates.

Returns
State descriptor index used for estimateWork() sampling, or -1 to disable.

◆ writePlotFile()

void amrex::AmrLevel::writePlotFile ( const std::string &  dir,
std::ostream &  os,
VisMF::How  how = VisMF::NFiles 
)
virtual

Write plot file data for this level into directory dir.

Parameters
dirOutput root (e.g., plt00010).
osMetadata stream shared across levels.
howVisMF write mode (NFiles or OneFilePerCPU).

◆ writePlotFilePost()

void amrex::AmrLevel::writePlotFilePost ( const std::string &  dir,
std::ostream &  os 
)
virtual

Perform post-plotfile cleanup for directory dir once writes complete.

Parameters
dirPlotfile directory that was just written.
osMetadata stream (same one passed to writePlotFilePre()).

◆ writePlotFilePre()

void amrex::AmrLevel::writePlotFilePre ( const std::string &  dir,
std::ostream &  os 
)
virtual

Perform pre-plotfile operations (e.g., gather metadata) for directory dir using os.

Parameters
dirPlotfile directory about to be written.
osStream receiving metadata (often the Header).

◆ writePlotNow()

bool amrex::AmrLevel::writePlotNow ( )
virtual

Return true if this level requests an immediate plotfile (default false).

◆ writeSmallPlotFile()

virtual void amrex::AmrLevel::writeSmallPlotFile ( const std::string &  dir,
std::ostream &  os,
VisMF::How  how = VisMF::NFiles 
)
inlinevirtual

Write a reduced plot file (same parameters as writePlotFile()) to dir.

Parameters
dirTarget directory.
osMetadata stream.
howVisMF write mode.

◆ writeSmallPlotNow()

bool amrex::AmrLevel::writeSmallPlotNow ( )
virtual

Return true if this level requests an immediate small-plot file (default false).

Friends And Related Symbol Documentation

◆ Amr

friend class Amr
friend

◆ FillPatchIterator

friend class FillPatchIterator
friend

◆ FillPatchIteratorHelper

friend class FillPatchIteratorHelper
friend

Member Data Documentation

◆ crse_ratio

IntVect amrex::AmrLevel::crse_ratio
protected

◆ derive_lst

DeriveList amrex::AmrLevel::derive_lst
staticprotected

◆ desc_lst

DescriptorList amrex::AmrLevel::desc_lst
staticprotected

◆ dmap

DistributionMapping amrex::AmrLevel::dmap
protected

◆ fine_ratio

IntVect amrex::AmrLevel::fine_ratio
protected

◆ geom

Geometry amrex::AmrLevel::geom
protected

◆ grids

BoxArray amrex::AmrLevel::grids
protected

◆ level

int amrex::AmrLevel::level {-1}
protected

◆ levelDirectoryCreated

bool amrex::AmrLevel::levelDirectoryCreated {false}
protected

◆ m_AreaNotToTag

BoxArray amrex::AmrLevel::m_AreaNotToTag
protected

◆ m_AreaToTag

Box amrex::AmrLevel::m_AreaToTag
protected

◆ m_eb_basic_grow_cells

int amrex::AmrLevel::m_eb_basic_grow_cells = 5
static

◆ m_eb_full_grow_cells

int amrex::AmrLevel::m_eb_full_grow_cells = 2
static

◆ m_eb_support_level

EBSupport amrex::AmrLevel::m_eb_support_level = EBSupport::volume
static

◆ m_eb_volume_grow_cells

int amrex::AmrLevel::m_eb_volume_grow_cells = 4
static

◆ m_factory

std::unique_ptr<FabFactory<FArrayBox> > amrex::AmrLevel::m_factory
protected

◆ m_fillpatcher

Vector<std::unique_ptr<FillPatcher<MultiFab> > > amrex::AmrLevel::m_fillpatcher
protected

◆ parent

Amr* amrex::AmrLevel::parent {nullptr}
protected

◆ post_step_regrid

int amrex::AmrLevel::post_step_regrid {0}
protected

◆ state

Vector<StateData> amrex::AmrLevel::state
protected

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