Block-Structured AMR Software Framework
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. More...
 
 AmrLevel (const AmrLevel &)=delete
 
 AmrLevel (AmrLevel &&)=delete
 
AmrLeveloperator= (const AmrLevel &)=delete
 
AmrLeveloperator= (AmrLevel &&)=delete
 
void LevelDirectoryNames (const std::string &dir, std::string &LevelDir, std::string &FullPath) const
 Get the level directory names. More...
 
virtual void CreateLevelDirectory (const std::string &dir)
 Create the Level_ directory for checkpoint and plot files. More...
 
void SetLevelDirectoryCreated (bool ldc) noexcept
 Set if the Level_ directory was created or to clear the value. CreateLevelDirectory sets levelDirectoryCreated = true. More...
 
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. More...
 
virtual void writePlotFile (const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles)
 Write plot file stuff to specified directory. More...
 
virtual void writePlotFilePre (const std::string &dir, std::ostream &os)
 Do pre-plotfile work to avoid synchronizations while writing the amr hierarchy. More...
 
virtual void writePlotFilePost (const std::string &dir, std::ostream &os)
 Do post-plotfile work to avoid synchronizations while writing the amr hierarchy. More...
 
virtual void writeSmallPlotFile (const std::string &, std::ostream &, VisMF::How=VisMF::NFiles)
 Write small plot file stuff to specified directory. More...
 
virtual void checkPoint (const std::string &dir, std::ostream &os, VisMF::How how=VisMF::NFiles, bool dump_old=true)
 Write current state to checkpoint file. More...
 
virtual void checkPointPre (const std::string &dir, std::ostream &os)
 Do pre-checkpoint work to avoid synchronizations while writing the amr hierarchy. More...
 
virtual void checkPointPost (const std::string &dir, std::ostream &os)
 Do post-checkpoint work to avoid synchronizations while writing the amr hierarchy. More...
 
virtual void restart (Amr &papa, std::istream &is, bool bReadSpecial=false)
 Restart from a checkpoint file. More...
 
virtual void set_state_in_checkpoint (Vector< int > &state_in_checkpoint)
 Old checkpoint may have different number of states than the new source code. More...
 
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. This is a pure virtual function and hence MUST be implemented by derived classes. More...
 
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. This is a pure virtual function and hence MUST be implemented by derived classes. More...
 
virtual Real advance (Real time, Real dt, int iteration, int ncycle)=0
 Do an integration step on this level. Returns maximum safe time step. This is a pure virtual function and hence MUST be implemented by derived classes. More...
 
virtual void post_timestep (int iteration)
 Contains operations to be done after a timestep. If this function is overridden, don't forget to reset FillPatcher. More...
 
virtual void postCoarseTimeStep (Real time)
 Contains operations to be done only after a full coarse timestep. The default implementation does nothing. More...
 
virtual void post_restart ()
 Operations to be done after restart. More...
 
virtual void post_regrid (int lbase, int new_finest)=0
 Operations to be done after regridding This is a pure virtual function and hence MUST be implemented by derived classes. More...
 
virtual void post_init (Real stop_time)=0
 Operations to be done after initialization. This is a pure virtual function and hence MUST be implemented by derived classes. More...
 
virtual int okToContinue ()
 Is it ok to continue the calculation? More...
 
virtual int okToRegrid ()
 Should I regrid with this level as base level? This test is only evaluated if regrid_int > 0 and level_count >= regrid_int as well. Defaults to true. More...
 
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. More...
 
virtual void setTimeLevel (Real time, Real dt_old, Real dt_new)
 Set the time levels of state data. More...
 
virtual void allocOldData ()
 Alloc space for old time data. More...
 
virtual void removeOldData ()
 Delete old-time data. More...
 
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. More...
 
virtual void init ()=0
 
void reset ()
 Reset data to initial time by swapping new and old time data. More...
 
int Level () const noexcept
 Returns this AmrLevel. More...
 
const BoxArrayboxArray () const noexcept
 List of grids at this level. More...
 
const BoxArraygetEdgeBoxArray (int dir) const noexcept
 
const BoxArraygetNodalBoxArray () const noexcept
 
const DistributionMappingDistributionMap () const noexcept
 
const FabFactory< FArrayBox > & Factory () const noexcept
 
int numGrids () const noexcept
 Number of grids at this level. More...
 
int numStates () const noexcept
 Number of states at this level. More...
 
const BoxDomain () const noexcept
 Returns the indices defining physical domain. More...
 
int nStep () const noexcept
 Timestep n at this level. More...
 
const GeometryGeom () const noexcept
 Returns the geometry object. More...
 
const IntVectfineRatio () const noexcept
 
Long countCells () const noexcept
 Returns number of cells on level. More...
 
const BoxArraygetAreaNotToTag () noexcept
 Get the area not to tag. More...
 
const BoxgetAreaToTag () noexcept
 
void constructAreaNotToTag ()
 Construct the area not to tag. More...
 
void setAreaNotToTag (BoxArray &ba) noexcept
 Set the area not to tag. More...
 
void resetFillPatcher ()
 
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. More...
 
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. More...
 
virtual void setPhysBoundaryValues (FArrayBox &dest, int state_indx, Real time, int dest_comp, int src_comp, int num_comp)
 Function to set physical boundary conditions. More...
 
virtual std::unique_ptr< MultiFabderive (const std::string &name, Real time, int ngrow)
 Returns a MultiFab containing the derived data for this level. The user is responsible for deleting this pointer when done with it. If ngrow>0 the MultiFab is built on the appropriately grown BoxArray. More...
 
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. More...
 
StateDataget_state_data (int state_indx) noexcept
 State data object. More...
 
MultiFabget_old_data (int state_indx) noexcept
 State data at old time. More...
 
const MultiFabget_old_data (int state_indx) const noexcept
 State data at old time. More...
 
MultiFabget_new_data (int state_indx) noexcept
 State data at new time. More...
 
const MultiFabget_new_data (int state_indx) const noexcept
 State data at new time. More...
 
int postStepRegrid () const noexcept
 Returns whether or not we want a post-timestep regrid. More...
 
void setPostStepRegrid (int new_val) noexcept
 Sets a new value for the post-timestep regrid trigger. More...
 
void UpdateDistributionMaps (DistributionMapping &dmap)
 Update the distribution maps in StateData based on the size of the map. More...
 
Vector< intgetBCArray (int State_Type, int gridno, int strt_comp, int ncomp)
 Boundary condition access function. More...
 
MultiFabget_data (int state_indx, Real time) noexcept
 Get state data at specified index and time. More...
 
virtual void set_preferred_boundary_values (MultiFab &S, int state_index, int scomp, int dcomp, int ncomp, Real time) const
 Hack to allow override of (non-fine-fine) fillpatched boundary data. More...
 
virtual void manual_tags_placement (TagBoxArray &tags, const Vector< IntVect > &bf_lev)
 Called in grid_places after other tagging routines to modify the list of tagged points. Default implementation does nothing. More...
 
virtual void setPlotVariables ()
 Modify list of variables to be plotted. More...
 
virtual void setSmallPlotVariables ()
 Modify list of variables to be plotted. More...
 
virtual Real estimateWork ()
 Estimate the amount of work required to advance Just this level based on the number of cells. This estimate can be overwritten with different methods. More...
 
virtual int WorkEstType ()
 Which state data type is for work estimates? -1 means none. More...
 
TimeLevel which_time (int state_indx, Real time) const noexcept
 Returns one the TimeLevel enums. Asserts that time is between AmrOldTime and AmrNewTime. More...
 
virtual bool writePlotNow ()
 Does the AmrLevel want Amr to write a plotfile now? More...
 
virtual bool writeSmallPlotNow ()
 Does the AmrLevel want Amr to write a small plotfile now? More...
 
void FillPatcherFill (amrex::MultiFab &mf, int dcomp, int ncomp, int nghost, amrex::Real time, int state_index, int scomp)
 Fill with FillPatcher on level > 0 and AmrLevel::FillPatch on level 0. More...
 
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) More...
 

Static Public Member Functions

static bool isStateVariable (const std::string &name, int &state_indx, int &ncomp)
 Is name a state variable? More...
 
static void FlushFPICache ()
 
static const DescriptorListget_desc_lst () noexcept
 Returns list of Descriptors. More...
 
static DeriveListget_derive_lst () noexcept
 Returns list of derived variables. More...
 
static void FillPatch (AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
 
static void FillPatchAdd (AmrLevel &amrlevel, MultiFab &leveldata, int boxGrow, Real time, int index, int scomp, int ncomp, int dcomp=0)
 
static IntVect ProperBlockingFactor (AmrLevel const &amr_level, int boxGrow, IndexType const &boxType, StateDescriptor const &desc, int SComp)
 Recommendation of a proper blocking factor. More...
 

Protected Member Functions

 AmrLevel () noexcept
 The constructors – for derived classes. More...
 
 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. More...
 

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
 

Private Member Functions

template<std::size_t order>
void storeRKCoarseData (int state_type, Real time, Real dt, MultiFab const &S_old, Array< MultiFab, order > const &rkk)
 
void FillRKPatch (int state_index, MultiFab &S, Real time, int stage, int iteration, int ncycle)
 

Private Attributes

BoxArray edge_grids [AMREX_SPACEDIM]
 
BoxArray nodal_grids
 

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

◆ AmrLevel() [2/4]

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

◆ 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

Do an integration step on this level. Returns maximum safe time step. This is a pure virtual function and hence MUST be implemented by derived classes.

◆ allocOldData()

void amrex::AmrLevel::allocOldData ( )
virtual

Alloc space for old time data.

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

◆ checkPointPost()

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

Do post-checkpoint work to avoid synchronizations while writing the amr hierarchy.

◆ checkPointPre()

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

Do pre-checkpoint work to avoid synchronizations while writing the amr hierarchy.

◆ 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. This is a pure virtual function and hence MUST be implemented by derived classes.

◆ 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. This is a pure virtual function and hence MUST be implemented by derived classes.

◆ constructAreaNotToTag()

void amrex::AmrLevel::constructAreaNotToTag ( )

Construct the area not to tag.

◆ countCells()

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

Returns number of cells on level.

◆ CreateLevelDirectory()

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

Create the Level_ directory for checkpoint and plot files.

◆ 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. The user is responsible for deleting this pointer when done with it. If ngrow>0 the MultiFab is built on the appropriately grown BoxArray.

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

◆ DistributionMap()

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

◆ Domain()

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

Returns the indices defining physical domain.

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

◆ estimateWork()

Real amrex::AmrLevel::estimateWork ( )
virtual

Estimate the amount of work required to advance Just this level based on the number of cells. This estimate can be overwritten with different methods.

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

◆ FillPatch()

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

◆ FillPatchAdd()

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

◆ FillPatcherFill()

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

Fill with FillPatcher on level > 0 and AmrLevel::FillPatch on level 0.

Parameters
mfdestination MultiFab
dcompstarting component for the destination
ncompnumber of component to fill
nghostnumber of ghost cells to fill
timetime
state_indexStateData index
scompstarting component in the StateData

◆ FillRKPatch()

void amrex::AmrLevel::FillRKPatch ( int  state_index,
MultiFab S,
Real  time,
int  stage,
int  iteration,
int  ncycle 
)
private

◆ 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

◆ 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

Get state data at specified index and 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

State data at new time.

◆ get_new_data() [2/2]

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

State data at new time.

◆ get_old_data() [1/2]

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

State data at old time.

◆ get_old_data() [2/2]

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

State data at old time.

◆ get_state_data()

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

State data object.

◆ getAreaNotToTag()

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

Get the area not to tag.

◆ getAreaToTag()

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

◆ getBCArray()

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

Boundary condition access function.

◆ getEdgeBoxArray()

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

◆ getNodalBoxArray()

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

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

◆ 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

Is name a state variable?

◆ Level()

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

Returns this AmrLevel.

◆ LevelDirectoryNames()

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

Get the level directory names.

◆ manual_tags_placement()

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

Called in grid_places after other tagging routines to modify the list of tagged points. Default implementation does nothing.

◆ 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

Is it ok to continue the calculation?

◆ okToRegrid()

int amrex::AmrLevel::okToRegrid ( )
virtual

Should I regrid with this level as base level? This test is only evaluated if regrid_int > 0 and level_count >= regrid_int as well. Defaults to true.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ post_init()

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

Operations to be done after initialization. This is a pure virtual function and hence MUST be implemented by derived classes.

◆ post_regrid()

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

Operations to be done after regridding This is a pure virtual function and hence MUST be implemented by derived classes.

◆ post_restart()

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

Operations to be done after restart.

◆ post_timestep()

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

Contains operations to be done after a timestep. If this function is overridden, don't forget to reset FillPatcher.

◆ postCoarseTimeStep()

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

Contains operations to be done only after a full coarse timestep. The default implementation does nothing.

◆ postStepRegrid()

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

Returns whether or not we want a post-timestep regrid.

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

◆ removeOldData()

void amrex::AmrLevel::removeOldData ( )
virtual

Delete old-time data.

◆ reset()

void amrex::AmrLevel::reset ( )

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

◆ resetFillPatcher()

void amrex::AmrLevel::resetFillPatcher ( )

◆ restart()

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

Restart from a checkpoint file.

◆ 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

Hack to allow override of (non-fine-fine) fillpatched boundary data.

◆ set_state_in_checkpoint()

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

Old checkpoint may have different number of states than the new source code.

◆ setAreaNotToTag()

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

Set the area not to tag.

◆ SetLevelDirectoryCreated()

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

Set if the Level_ directory was created or to clear the value. CreateLevelDirectory sets levelDirectoryCreated = true.

◆ setPhysBoundaryValues()

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

Function to set physical boundary conditions.

◆ setPlotVariables()

void amrex::AmrLevel::setPlotVariables ( )
virtual

Modify list of variables to be plotted.

◆ setPostStepRegrid()

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

Sets a new value for the post-timestep regrid trigger.

◆ setSmallPlotVariables()

void amrex::AmrLevel::setSmallPlotVariables ( )
virtual

Modify list of variables to be plotted.

◆ setTimeLevel()

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

Set the time levels of state data.

◆ storeRKCoarseData()

template<std::size_t order>
void amrex::AmrLevel::storeRKCoarseData ( int  state_type,
Real  time,
Real  dt,
MultiFab const &  S_old,
Array< MultiFab, order > const &  rkk 
)
private

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

◆ UpdateDistributionMaps()

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

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

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

◆ WorkEstType()

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

Which state data type is for work estimates? -1 means none.

◆ writePlotFile()

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

Write plot file stuff to specified directory.

◆ writePlotFilePost()

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

Do post-plotfile work to avoid synchronizations while writing the amr hierarchy.

◆ writePlotFilePre()

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

Do pre-plotfile work to avoid synchronizations while writing the amr hierarchy.

◆ writePlotNow()

bool amrex::AmrLevel::writePlotNow ( )
virtual

Does the AmrLevel want Amr to write a plotfile now?

◆ writeSmallPlotFile()

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

Write small plot file stuff to specified directory.

◆ writeSmallPlotNow()

bool amrex::AmrLevel::writeSmallPlotNow ( )
virtual

Does the AmrLevel want Amr to write a small plotfile now?

Friends And Related Function 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

◆ edge_grids

BoxArray amrex::AmrLevel::edge_grids[AMREX_SPACEDIM]
mutableprivate

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

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

◆ m_fillpatcher

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

◆ nodal_grids

BoxArray amrex::AmrLevel::nodal_grids
mutableprivate

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