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

Manage hierarchy of levels for time-dependent AMR computations. More...

#include <AMReX_Amr.H>

Inheritance diagram for amrex::Amr:
amrex::AmrCore amrex::AmrMesh amrex::AmrInfo

Public Member Functions

 Amr (LevelBld *a_levelbld)
 Construct an AMR driver bound to a factory of AmrLevel implementations.
 
 Amr (const RealBox *rb, int max_level_in, const Vector< int > &n_cell_in, int coord, LevelBld *a_levelbld)
 Construct with explicit problem domain metadata and a level factory.
 
 Amr (const Amr &rhs)=delete
 Amr owns global driver state—disable copying to avoid duplicate hierarchies.
 
 Amr (Amr &&rhs)=delete
 Amr instances are not movable; they rely on pointers registered with static callbacks.
 
Amroperator= (const Amr &rhs)=delete
 Copy assignment is likewise disabled.
 
Amroperator= (Amr &&rhs)=delete
 Move assignment is disabled to keep ownership contracts simple.
 
void InitAmr ()
 Initialize Amr-wide defaults and parse the amr.* ParmParse options.
 
 ~Amr () override
 The destructor.
 
virtual void init (Real strt_time, Real stop_time)
 Initialize hierarchy state (grid creation, initial conditions, etc.).
 
void InitializeInit (Real strt_time, Real stop_time, const BoxArray *lev0_grids=nullptr, const Vector< int > *pmap=nullptr)
 First half of the two-phase initialization used by initialInit().
 
void FinalizeInit (Real strt_time, Real stop_time)
 Second half of initial initialization; pairs with InitializeInit() using times strt_time and stop_time.
 
void setDtLevel (const Vector< Real > &dt_lev) noexcept
 Set per-level time steps to the values in dt_lev.
 
void setDtLevel (Real dt, int lev) noexcept
 Set the time step stored at AMR level lev to dt.
 
void setDtMin (const Vector< Real > &dt_min_in) noexcept
 Override the maximum allowable dt (typically CFL-limited) for each level using dt_min_in.
 
void setNCycle (const Vector< int > &ns) noexcept
 Provide the number of time subcycles per level via ns.
 
int subCycle () const noexcept
 Subcycle in time?
 
const std::string & subcyclingMode () const noexcept
 How are we subcycling?
 
int level_being_advanced () const noexcept
 What is "level" in Amr::timeStep? This is only relevant if we are still in Amr::timeStep; it is set back to -1 on leaving Amr::timeStep.
 
Real cumTime () const noexcept
 Physical time (time advanced so far).
 
void setCumTime (Real t) noexcept
 Set the accumulated physical time to t.
 
Real startTime () const noexcept
 Physical time this simulation started.
 
void setStartTime (Real t) noexcept
 Set the recorded start time to t.
 
Real dtLevel (int level) const noexcept
 Time step currently stored for level level.
 
Real dtMin (int level) const noexcept
 Max time step (typically physics-limited) for level level.
 
const Vector< Real > & dtLevel () const noexcept
 Vector of time steps (dt) for all levels.
 
int nCycle (int level) const noexcept
 Number of subcycled time steps taken per coarse step on level level.
 
int levelSteps (int lev) const noexcept
 Number of time steps completed so far on level lev.
 
void setLevelSteps (int lev, int n) noexcept
 Override the recorded time-step count for level lev.
 
int levelCount (int lev) const noexcept
 Which cycle number are we on for level lev?
 
void setLevelCount (int lev, int n) noexcept
 Set the current cycle count for level lev.
 
int regridInt (int lev) const noexcept
 Interval (in steps) between regrids on level lev.
 
int checkInt () const noexcept
 Number of time steps between checkpoint files.
 
Real checkPer () const noexcept
 Time between checkpoint files.
 
int plotInt () const noexcept
 Number of time steps between plot files.
 
Real plotPer () const noexcept
 Time between plot files.
 
Real plotLogPer () const noexcept
 Spacing in log10(time) of logarithmically spaced plot files.
 
int plotMaxLevel () const noexcept
 Maximum level index to include when writing plot files.
 
int smallplotInt () const noexcept
 Number of time steps between small plot files.
 
Real smallplotPer () const noexcept
 Time between plot files.
 
Real smallplotLogPer () const noexcept
 Spacing in log10(time) of logarithmically spaced small plot files.
 
AmrLevelgetLevel (int lev) noexcept
 Access the AmrLevel object stored at index lev.
 
Vector< std::unique_ptr< AmrLevel > > & getAmrLevels () noexcept
 Return the container of level objects (index 0..finest).
 
Long cellCount () noexcept
 
Long cellCount (int lev) noexcept
 
int numGrids () noexcept
 
int numGrids (int lev) noexcept
 
int okToContinue () noexcept
 
void regrid (int lbase, Real time, bool initial=false) override
 Rebuild grid hierarchy finer than lbase.
 
void RegridOnly (Real time, bool do_io=true)
 Regrid the hierarchy at physical time time without advancing any levels.
 
bool okToRegrid (int level) noexcept
 Return true if level level is permitted to regrid on this step.
 
virtual void coarseTimeStep (Real stop_time)
 Advance the hierarchy by one coarse step.
 
Real coarseTimeStepDt (Real stop_time)
 Advance once and report the coarse level dt that was used.
 
std::unique_ptr< MultiFabderive (const std::string &name, Real time, int lev, int ngrow)
 Build a MultiFab containing the derived quantity named name.
 
void derive (const std::string &name, Real time, const Vector< MultiFab * > &mf, int dcomp)
 Compute derived data on all levels and store it in the supplied MultiFabs.
 
Vector< std::unique_ptr< MultiFab > > derive (const std::string &name, amrex::Real time, int ngrow)
 Compute derived data on all levels and return as owning MultiFabs.
 
const std::string & theRestartFile () const noexcept
 Name of the most recent restart checkpoint file.
 
const std::string & theRestartPlotFile () const noexcept
 Name of the restart plotfile.
 
std::ostream & DataLog (int i)
 Return the i th datalog stream (opened lazily) for writing diagnostics.
 
std::string DataLogName (int i) const noexcept
 The filename backing the i th datalog stream.
 
int NumDataLogs () noexcept
 
virtual void writePlotFile ()
 Write a full-resolution plot file to disk using the current hierarchy state.
 
int stepOfLastPlotFile () const noexcept
 Return the AMR step index associated with the most recent plot file.
 
virtual void writeSmallPlotFile ()
 Write the reduced “small plot” file (subset of variables/levels) to disk.
 
int stepOfLastSmallPlotFile () const noexcept
 Return the AMR step index associated with the most recent small plot file.
 
virtual void checkPoint ()
 Write current state into a chk* checkpoint directory.
 
int stepOfLastCheckPoint () const noexcept
 Return the AMR step index associated with the most recent checkpoint.
 
void setBoundaryGeometry (BoundaryPointList &IntersectLoX, BoundaryPointList &IntersectHiX, BoundaryPointList &IntersectLoY, BoundaryPointList &IntersectHiY) noexcept
 Specialized version: register intersection points for faces orthogonal to X/Y.
 
void setBoundaryGeometry (BoundaryPointList &IntersectLoX, BoundaryPointList &IntersectHiX, BoundaryPointList &IntersectLoY, BoundaryPointList &IntersectHiY, BoundaryPointList &IntersectLoZ, BoundaryPointList &IntersectHiZ) noexcept
 General version: register intersections for all Cartesian faces.
 
BoundaryPointList & getIntersectLoX () noexcept
 
BoundaryPointList & getIntersectHiX () noexcept
 
BoundaryPointList & getIntersectLoY () noexcept
 
BoundaryPointList & getIntersectHiY () noexcept
 
BoundaryPointList & getIntersectLoZ () noexcept
 
BoundaryPointList & getIntersectHiZ () noexcept
 
void RedistributeParticles ()
 Trigger a particle redistribution pass across all levels.
 
void InstallNewDistributionMap (int lev, const DistributionMapping &newdm)
 Replace the distribution map on level lev with newdm and update metadata.
 
- Public Member Functions inherited from amrex::AmrCore
 AmrCore ()
 Construct an empty AmrCore.
 
 AmrCore (const RealBox *rb, int max_level_in, const Vector< int > &n_cell_in, int coord=-1, Vector< IntVect > ref_ratios=Vector< IntVect >(), const int *is_per=nullptr)
 Construct using legacy inputs that mirror AmrMesh.
 
 AmrCore (const RealBox &rb, int max_level_in, const Vector< int > &n_cell_in, int coord, Vector< IntVect > const &ref_ratios, Array< int, 3 > const &is_per)
 Construct with an explicit domain and periodicity array.
 
 AmrCore (Geometry const &level_0_geom, AmrInfo const &amr_info)
 Construct from a ready-made Geometry and AmrInfo bundle.
 
 AmrCore (AmrCore &&rhs) noexcept
 
AmrCoreoperator= (AmrCore &&rhs) noexcept
 
 AmrCore (const AmrCore &rhs)=delete
 
AmrCoreoperator= (const AmrCore &rhs)=delete
 
 ~AmrCore () override
 
AmrParGDBGetParGDB () const noexcept
 Return AmrParGDB pointer for particle containers.
 
void InitFromScratch (Real time)
 Initialize BoxArray, DistributionMapping, and level data from scratch.
 
void printGridSummary (std::ostream &os, int min_lev, int max_lev) const noexcept
 Print a human-readable summary of grid sizes between min_lev and max_lev.
 
- Public Member Functions inherited from amrex::AmrMesh
 AmrMesh ()
 Construct an empty mesh.
 
 AmrMesh (const RealBox *rb, int max_level_in, const Vector< int > &n_cell_in, int coord=-1, Vector< IntVect > refrat=Vector< IntVect >(), const int *is_per=nullptr)
 Construct using legacy inputs (same signature as AmrCore).
 
 AmrMesh (const RealBox &rb, int max_level_in, const Vector< int > &n_cell_in, int coord, Vector< IntVect > const &a_refrat, Array< int, 3 > const &is_per)
 Construct with explicit domain and periodicity array.
 
 AmrMesh (Geometry const &level_0_geom, AmrInfo const &amr_info)
 Construct from a pre-built Geometry and AmrInfo bundle.
 
 AmrMesh (const AmrMesh &rhs)=delete
 
AmrMeshoperator= (const AmrMesh &rhs)=delete
 
 AmrMesh (AmrMesh &&rhs)=default
 
AmrMeshoperator= (AmrMesh &&rhs)=default
 
virtual ~AmrMesh ()=default
 
int Verbose () const noexcept
 
int maxLevel () const noexcept
 Return the max level.
 
int finestLevel () const noexcept
 Return the finest level.
 
IntVect refRatio (int lev) const noexcept
 Return the refinement ratio for level lev.
 
int MaxRefRatio (int lev) const noexcept
 Return the maximum refinement ratio in any direction.
 
const Vector< IntVect > & refRatio () const noexcept
 Return refinement ratios between all levels.
 
const Vector< Geometry > & Geom () const noexcept
 
const Vector< DistributionMapping > & DistributionMap () const noexcept
 
const Vector< BoxArray > & boxArray () const noexcept
 
const GeometryGeom (int lev) const noexcept
 
const DistributionMappingDistributionMap (int lev) const noexcept
 
const BoxArrayboxArray (int lev) const noexcept
 
Vector< GeometryGeom (int a_coarsest_lev, int a_finest_lev) const noexcept
 
Vector< BoxArrayboxArray (int a_coarsest_lev, int a_finest_lev) const noexcept
 
Vector< DistributionMappingDistributionMap (int a_coarsest_lev, int a_finest_lev) const noexcept
 
Vector< Geometry > & Geom () noexcept
 
GeometryGeom (int lev) noexcept
 
void SetMaxGridSize (int new_mgs) noexcept
 Override the max grid size constraint for all levels (scalar version).
 
void SetMaxGridSize (const IntVect &new_mgs) noexcept
 Override the max grid size constraint using explicit IntVect new_mgs.
 
void SetMaxGridSize (const Vector< int > &new_mgs) noexcept
 Override per-level max grid sizes using scalar-per-level vector new_mgs.
 
void SetMaxGridSize (const Vector< IntVect > &new_mgs) noexcept
 Override per-level max grid sizes using explicit IntVect-per-level vector new_mgs.
 
void SetBlockingFactor (int new_bf) noexcept
 Override the blocking factor constraint for all levels (scalar version).
 
void SetBlockingFactor (const IntVect &new_bf) noexcept
 Override blocking factor using explicit IntVect new_bf.
 
void SetBlockingFactor (const Vector< int > &new_bf) noexcept
 Override per-level blocking factors using scalar-per-level vector new_bf.
 
void SetBlockingFactor (const Vector< IntVect > &new_bf) noexcept
 Override per-level blocking factors using explicit IntVect vector new_bf.
 
void SetGridEff (Real eff) noexcept
 Set the minimum acceptable grid efficiency to eff.
 
void SetNProper (int n) noexcept
 Set the number of cell blocks required for proper nesting to n.
 
void SetFinestLevel (int new_finest_level) noexcept
 Update the recorded finest level index to new_finest_level (without reallocating).
 
void SetDistributionMap (int lev, const DistributionMapping &dmap_in) noexcept
 Install a DistributionMapping for level lev.
 
void SetBoxArray (int lev, const BoxArray &ba_in) noexcept
 Install a BoxArray for level lev.
 
void SetGeometry (int lev, const Geometry &geom_in) noexcept
 Install a Geometry object for level lev.
 
int GetLevel (Box const &domain) const noexcept
 Given domain box, return AMR level. Return -1 if there is no match.
 
void ClearDistributionMap (int lev) noexcept
 Clear the stored DistributionMapping for level lev.
 
void ClearBoxArray (int lev) noexcept
 Clear the stored BoxArray for level lev.
 
int nErrorBuf (int lev, int direction=0) const noexcept
 Return the number of buffer cells (as a single integer) in error estimator.
 
const IntVectnErrorBufVect (int lev) const noexcept
 Return the number of buffer cells (as an IntVect) in error estimator.
 
Real gridEff () const noexcept
 Return the minimum allowable grid efficiency.
 
int nProper () const noexcept
 Return the number of cell blocks to define proper nesting.
 
const IntVectblockingFactor (int lev) const noexcept
 Return the blocking factor at level lev.
 
const IntVectmaxGridSize (int lev) const noexcept
 Return the largest allowable grid.
 
bool LevelDefined (int lev) const noexcept
 Return true if valid grids/dmaps exist for level lev.
 
bool useFixedCoarseGrids () const noexcept
 Should we keep the coarser grids fixed (and not regrid those levels) at all?
 
int useFixedUpToLevel () const noexcept
 Up to what level should we keep the coarser grids fixed (and not regrid those levels)?
 
void ChopGrids (int lev, BoxArray &ba, int target_size) const
 Try to chop grids on level lev until BoxArray ba has at least target_size boxes.
 
BoxArray MakeBaseGrids () const
 Make a level-0 BoxArray covering the whole domain (without installing it).
 
void MakeNewGrids (int lbase, Real time, int &new_finest, Vector< BoxArray > &new_grids)
 Make new grids based on error estimates. This function expects that valid BoxArrays exist in this->grids from level lbase to level this->finest_level (the current finest level). new_grids. On return, the new finest level is stored in new_finest, and the new grids are stored in new_grids from Array element lbase+1 to new_finest_level (unless fixed grids are used). Note that this function adds at most one more level to the existing levels, and it may remove all levels above the base level. This function does not change the value of this->finest_level, nor does it modifies any BoxArrays stored in this->grids. It also does not modify new_grids's elements outside the range [lbase+1,new_finest_level].
 
void MakeNewGrids (Real time=0.0)
 Convenience overload that rebuilds every level (including level 0) at time time.
 
virtual void PostProcessBaseGrids (BoxArray &box_array) const
 Allow derived classes to edit the base-level grids before installing them.
 
Long CountCells (int lev) const noexcept
 Count the total number of cells at level lev.
 
virtual DistributionMapping MakeDistributionMap (int lev, BoxArray const &ba)
 Build a distribution map for level lev grid layout ba (override to customize load balancing).
 

Static Public Member Functions

static bool RegridOnRestart () noexcept
 Whether to regrid right after restart.
 
static const std::list< std::string > & statePlotVars () noexcept
 The names of state variables to output in the plotfile. They can be set using the amr.plot_vars variable in a ParmParse inputs file.
 
static const std::list< std::string > & stateSmallPlotVars () noexcept
 List of state variables written to the reduced “small” plot files.
 
static bool isStatePlotVar (const std::string &name)
 Return true if name is already registered for full plotfiles.
 
static bool isStateSmallPlotVar (const std::string &name)
 Return true if name is already registered for small plotfiles.
 
static void addStatePlotVar (const std::string &name)
 Register name as a state variable to include in full plotfiles (no-op if already present).
 
static void addStateSmallPlotVar (const std::string &name)
 Append name to the small-plot variable list if not already present.
 
static void deleteStatePlotVar (const std::string &name)
 Remove name from the full-plot list if present.
 
static void clearStatePlotVarList ()
 Clear the list of state_plot_vars.
 
static void clearStateSmallPlotVarList ()
 Clear the list of variables written to small plot files.
 
static void fillStatePlotVarList ()
 Fill the list of state_plot_vars with all of the state quantities.
 
static void fillStateSmallPlotVarList ()
 Populate state_small_plot_vars with every state quantity registered on AmrLevel 0.
 
static bool Plot_Files_Output ()
 Write out plotfiles (True/False)?
 
static const std::list< std::string > & derivePlotVars () noexcept
 The names of derived variables to output in the plotfile. They can be set using the amr.derive_plot_vars variable in a ParmParse inputs file.
 
static const std::list< std::string > & deriveSmallPlotVars () noexcept
 List of derived quantities included when writing small plot files.
 
static bool isDerivePlotVar (const std::string &name) noexcept
 Return true if name is already registered for full derived-output.
 
static bool isDeriveSmallPlotVar (const std::string &name) noexcept
 Return true if name is already registered for small derived-output.
 
static void addDerivePlotVar (const std::string &name)
 Register name as a derived quantity to include in full plotfiles (no-op if already present).
 
static void addDeriveSmallPlotVar (const std::string &name)
 Append name to the small-plot derived list if not already present.
 
static void deleteDerivePlotVar (const std::string &name)
 Remove name from the full derived list.
 
static void deleteDeriveSmallPlotVar (const std::string &name)
 Remove name from the small derived list.
 
static void clearDerivePlotVarList ()
 Clear the list of derive_plot_vars.
 
static void clearDeriveSmallPlotVarList ()
 Clear the list of derived quantities written to small plot files.
 
static void fillDerivePlotVarList ()
 Fill the list of derive_plot_vars with all derived quantities.
 
static void fillDeriveSmallPlotVarList ()
 Populate derive_small_plot_vars with every registered derived quantity.
 
static void setComputeNewDtOnRegrid (bool flag)
 Toggle whether new time steps are recomputed immediately after each regrid (flag=true recomputes).
 
static void Initialize ()
 Initialize global, Amr-wide static state (e.g., ParmParse defaults, logs).
 
static void Finalize ()
 Release global Amr static resources allocated by Initialize().
 
static const BoxArrayinitialBa (int level) noexcept
 Return the user-specified BoxArray that initializes level level (1-based).
 
static int initialBaLevels () noexcept
 
static Real computeOptimalSubcycling (int n, int *best, const Real *dt_max, const Real *est_work, const int *cycle_max)
 Compute the optimal subcycling pattern. This assumes that anything less than cycle_max[i] is a valid number of subcycles at level[i]. For example if ref_ratio[i] = cycle_max[i] = 4, then 1,2,3,4 are all valid values for n_cycles[i].
 
static const Vector< BoxArray > & getInitialBA () noexcept
 
static bool UsingPrecreateDirectories () noexcept
 Whether plot/checkpoint writers pre-create directories (helps avoid races on some filesystems).
 

Public Attributes

BoundaryPointList intersect_lox
 
BoundaryPointList intersect_loy
 
BoundaryPointList intersect_loz
 
BoundaryPointList intersect_hix
 
BoundaryPointList intersect_hiy
 
BoundaryPointList intersect_hiz
 

Static Public Attributes

static bool first_smallplotfile
 

Protected Member Functions

void initialInit (Real strt_time, Real stop_time, const BoxArray *lev0_grids=nullptr, const Vector< int > *pmap=nullptr)
 Initialize the grid hierarchy (called by Amr::init()).
 
void readProbinFile (int &init)
 Read the probin file, updating the flag init that tracks initialization iterations.
 
void checkInput ()
 Check for valid input.
 
void restart (const std::string &filename)
 Restart from a checkpoint file named filename.
 
void defBaseLevel (Real strt_time, const BoxArray *lev0_grids=nullptr, const Vector< int > *pmap=nullptr)
 Define and initialize the coarsest level using optional level-0 grids lev0_grids and map pmap.
 
void bldFineLevels (Real strt_time)
 Define and initialize refined levels at simulation time strt_time.
 
virtual void regrid_level_0_on_restart ()
 Regrid level 0 on restart.
 
void grid_places (int lbase, Real time, int &new_finest, Vector< BoxArray > &new_grids)
 Define new grid locations (called from regrid) and store them in new_grids.
 
DistributionMapping makeLoadBalanceDistributionMap (int lev, Real time, const BoxArray &ba) const
 Build a DistributionMapping for level lev at time time using grid layout ba.
 
void LoadBalanceLevel0 (Real time)
 Compute and install a new distribution map for level 0.
 
void ErrorEst (int lev, TagBoxArray &tags, Real time, int ngrow) override
 Default error-estimation driver; forwards to the levelbld-defined callbacks.
 
BoxArray GetAreaNotToTag (int lev) override
 Return cached regions that must not be tagged on level lev.
 
void ManualTagsPlacement (int lev, TagBoxArray &tags, const Vector< IntVect > &bf_lev) override
 Allow derived classes to override tagging manually on level lev.
 
virtual void timeStep (int level, Real time, int iteration, int niter, Real stop_time)
 Do a single time step on level level.
 
void MakeNewLevelFromScratch (int, Real, const BoxArray &, const DistributionMapping &) override
 Create and fill a new level from scratch.
 
void MakeNewLevelFromCoarse (int, Real, const BoxArray &, const DistributionMapping &) override
 Make a new level using provided metadata and populate it from the next coarser level.
 
void RemakeLevel (int, Real, const BoxArray &, const DistributionMapping &) override
 Remake an existing level using provided metadata and refill from old data.
 
void ClearLevel (int) override
 Delete level lev data structures owned by the derived class.
 
bool checkPointNow () noexcept
 Whether checkpoint criteria are satisfied on this step.
 
bool writePlotNow () noexcept
 Whether to write a plotfile now.
 
bool writeSmallPlotNow () noexcept
 Whether to emit a small plot file on this step.
 
void printGridInfo (std::ostream &os, int min_lev, int max_lev)
 Print per-level grid metadata to os for levels [min_lev, max_lev].
 
void setRecordGridInfo (const std::string &filename)
 Record the log filename filename used when dumping grid info metadata.
 
void setRecordRunInfo (const std::string &filename)
 Record the log filename filename used when writing detailed run information.
 
void setRecordRunInfoTerse (const std::string &filename)
 Record the log filename filename used for terse run information.
 
void setRecordDataInfo (int i, const std::string &filename)
 Configure the i th datalog stream to mirror file filename.
 
void initSubcycle ()
 Parse ParmParse state and initialize subcycling-related arrays/flags.
 
void initPltAndChk ()
 Initialize plot/checkpoint cadence, names, and counters from inputs.
 
int updateInSitu ()
 Drive any per-step in-situ actions; returns 0 on success.
 
- Protected Member Functions inherited from amrex::AmrMesh
void checkInput ()
 Validate AmrInfo inputs (throws if inconsistent).
 
void SetIterateToFalse () noexcept
 Disable the iterate-on-new-grids loop for the next MakeNewGrids call.
 
void SetUseNewChop () noexcept
 Force MakeNewGrids to use the newer chopping algorithm.
 

Static Protected Member Functions

static int initInSitu ()
 Set up in-situ analysis hooks; returns 0 on success.
 
static int finalizeInSitu ()
 Tear down in-situ infrastructure; returns 0 on success.
 

Protected Attributes

std::string regrid_grids_file
 Grids file that will bypass regridding.
 
std::string initial_grids_file
 Grids file that will bypass regridding only at initialization.
 
Vector< std::unique_ptr< AmrLevel > > amr_level
 Vector of levels.
 
Real cumtime = std::numeric_limits<Real>::lowest()
 Physical time variable.
 
Real start_time = std::numeric_limits<Real>::lowest()
 Physical time this simulation started.
 
Vector< Realdt_level
 Timestep at this level.
 
Vector< intlevel_steps
 Number of time steps at this level.
 
Vector< intlevel_count
 
Vector< intn_cycle
 
std::string subcycling_mode
 Type of subcycling to use.
 
Vector< Realdt_min
 
Vector< intregrid_int
 Interval between regridding.
 
int last_checkpoint
 Step number of previous checkpoint.
 
int check_int
 How often checkpoint (# time steps).
 
Real check_per
 How often checkpoint (units of time).
 
std::string check_file_root
 Root name of checkpoint file.
 
int last_plotfile
 Step number of previous plotfile.
 
int last_smallplotfile
 Step number of previous small plotfile.
 
int plot_int
 How often plotfile (# of time steps)
 
Real plot_per
 How often plotfile (in units of time)
 
Real plot_log_per
 How often plotfile (in units of log10(time))
 
int plot_max_level
 Maximum AMR level to write to a plotfile.
 
int small_plot_int
 How often small plotfile (# of time steps)
 
Real small_plot_per
 How often small plotfile (in units of time)
 
Real small_plot_log_per
 How often small plotfile (in units of log10(time))
 
bool write_plotfile_with_checkpoint
 Write out a plotfile whenever we checkpoint.
 
int file_name_digits
 How many digits to use in the plotfile and checkpoint names.
 
int message_int
 How often checking messages touched by user, such as "stop_run".
 
std::string plot_file_root
 Root name of plotfile.
 
std::string small_plot_file_root
 Root name of small plotfile.
 
int which_level_being_advanced = -1
 Only >=0 if we are in Amr::timeStep(level,...)
 
int record_grid_info
 
int record_run_info
 
int record_run_info_terse
 
std::ofstream gridlog
 
std::ofstream runlog
 
std::ofstream runlog_terse
 
Vector< std::unique_ptr< std::fstream > > datalog
 
Vector< std::string > datalogname
 
int sub_cycle
 
std::string restart_chkfile
 
std::string restart_pltfile
 
std::string probin_file
 
LevelBldlevelbld
 
bool abort_on_stream_retry_failure
 
int stream_max_tries
 
int loadbalance_with_workestimates
 
int loadbalance_level0_int
 
Real loadbalance_max_fac
 
bool bUserStopRequest
 
- Protected Attributes inherited from amrex::AmrCore
std::unique_ptr< AmrParGDBm_gdb
 Geometry, DistributionMapping and BoxArray for particle containers.
 
- Protected Attributes inherited from amrex::AmrMesh
int finest_level
 Current finest level.
 
Vector< Geometrygeom
 
Vector< DistributionMappingdmap
 
Vector< BoxArraygrids
 
unsigned int num_setdm = 0
 
unsigned int num_setba = 0
 
- Protected Attributes inherited from amrex::AmrInfo
int verbose = 0
 
int max_level = 0
 Maximum allowed level.
 
Vector< IntVectref_ratio {{IntVect(2)}}
 Refinement ratios.
 
Vector< IntVectblocking_factor {{IntVect(8)}}
 Blocking factor in grid generation (by level).
 
Vector< IntVectmax_grid_size {{IntVect( 64 )}}
 Maximum allowable grid size (by level).
 
Vector< IntVectn_error_buf {{IntVect(1)}}
 Buffer cells around each tagged cell.
 
Real grid_eff = static_cast<Real>(0.7)
 Grid efficiency.
 
int n_proper = 1
 Cells required for proper nesting.
 
int use_fixed_upto_level = 0
 
bool use_fixed_coarse_grids = false
 
bool refine_grid_layout = true
 
IntVect refine_grid_layout_dims = IntVect(1)
 
bool check_input = true
 
bool use_new_chop = false
 
bool iterate_on_new_grids = true
 
int max_grid_iterations = 4
 

Static Protected Attributes

static std::list< std::string > state_plot_vars
 State Vars to dump to plotfile.
 
static std::list< std::string > state_small_plot_vars
 State Vars to dump to small plotfile.
 
static std::list< std::string > derive_plot_vars
 Derived Vars to dump to plotfile.
 
static std::list< std::string > derive_small_plot_vars
 Derived Vars to dump to small plotfile.
 
static bool first_plotfile
 
static Vector< BoxArrayinitial_ba
 Array of BoxArrays read in to initially define grid hierarchy.
 
static Vector< BoxArrayregrid_ba
 Array of BoxArrays read in to externally define grid hierarchy at each regrid.
 
static bool compute_new_dt_on_regrid
 

Detailed Description

Manage hierarchy of levels for time-dependent AMR computations.

The Amr class is designed to manage parts of the computation which do not belong on a single level, like establishing and updating the hierarchy of levels, global timestepping, and managing the different AmrLevels

Constructor & Destructor Documentation

◆ Amr() [1/4]

amrex::Amr::Amr ( LevelBld a_levelbld)

Construct an AMR driver bound to a factory of AmrLevel implementations.

Parameters
a_levelbldFactory responsible for creating level objects at runtime (must be non-null).

◆ Amr() [2/4]

amrex::Amr::Amr ( const RealBox rb,
int  max_level_in,
const Vector< int > &  n_cell_in,
int  coord,
LevelBld a_levelbld 
)

Construct with explicit problem domain metadata and a level factory.

Parameters
rbOptional physical domain for level 0.
max_level_inMaximum AMR level allowed.
n_cell_inLevel-0 cell counts per direction.
coordCoordinate system index.
a_levelbldFactory responsible for creating AmrLevel instances.

◆ Amr() [3/4]

amrex::Amr::Amr ( const Amr rhs)
delete

Amr owns global driver state—disable copying to avoid duplicate hierarchies.

◆ Amr() [4/4]

amrex::Amr::Amr ( Amr &&  rhs)
delete

Amr instances are not movable; they rely on pointers registered with static callbacks.

◆ ~Amr()

amrex::Amr::~Amr ( )
override

The destructor.

Member Function Documentation

◆ addDerivePlotVar()

void amrex::Amr::addDerivePlotVar ( const std::string &  name)
static

Register name as a derived quantity to include in full plotfiles (no-op if already present).

Parameters
nameDerived quantity identifier.

◆ addDeriveSmallPlotVar()

void amrex::Amr::addDeriveSmallPlotVar ( const std::string &  name)
static

Append name to the small-plot derived list if not already present.

Parameters
nameDerived quantity identifier.

◆ addStatePlotVar()

void amrex::Amr::addStatePlotVar ( const std::string &  name)
static

Register name as a state variable to include in full plotfiles (no-op if already present).

Parameters
nameState variable identifier.

◆ addStateSmallPlotVar()

void amrex::Amr::addStateSmallPlotVar ( const std::string &  name)
static

Append name to the small-plot variable list if not already present.

Parameters
nameState variable identifier.

◆ bldFineLevels()

void amrex::Amr::bldFineLevels ( Real  strt_time)
protected

Define and initialize refined levels at simulation time strt_time.

◆ cellCount() [1/2]

Long amrex::Amr::cellCount ( )
noexcept
Returns
Total number of cells owned by all levels.

◆ cellCount() [2/2]

Long amrex::Amr::cellCount ( int  lev)
noexcept
Returns
Number of cells at level lev.

◆ checkInput()

void amrex::Amr::checkInput ( )
protected

Check for valid input.

◆ checkInt()

int amrex::Amr::checkInt ( ) const
inlinenoexcept

Number of time steps between checkpoint files.

◆ checkPer()

Real amrex::Amr::checkPer ( ) const
inlinenoexcept

Time between checkpoint files.

◆ checkPointNow()

bool amrex::Amr::checkPointNow ( )
protectednoexcept

Whether checkpoint criteria are satisfied on this step.

◆ clearDerivePlotVarList()

void amrex::Amr::clearDerivePlotVarList ( )
static

Clear the list of derive_plot_vars.

◆ clearDeriveSmallPlotVarList()

void amrex::Amr::clearDeriveSmallPlotVarList ( )
static

Clear the list of derived quantities written to small plot files.

◆ ClearLevel()

void amrex::Amr::ClearLevel ( int  lev)
inlineoverrideprotectedvirtual

Delete level lev data structures owned by the derived class.

Parameters
levLevel index whose data should be released.

Implements amrex::AmrCore.

◆ clearStatePlotVarList()

void amrex::Amr::clearStatePlotVarList ( )
static

Clear the list of state_plot_vars.

◆ clearStateSmallPlotVarList()

void amrex::Amr::clearStateSmallPlotVarList ( )
static

Clear the list of variables written to small plot files.

◆ coarseTimeStep()

void amrex::Amr::coarseTimeStep ( Real  stop_time)
virtual

Advance the hierarchy by one coarse step.

Parameters
stop_timeOptional time limit; time stepping will stop once this value is reached.

◆ coarseTimeStepDt()

Real amrex::Amr::coarseTimeStepDt ( Real  stop_time)

Advance once and report the coarse level dt that was used.

Parameters
stop_timeOptional time limit (same semantics as coarseTimeStep()).
Returns
The dt associated with level 0 for that step.

◆ computeOptimalSubcycling()

Real amrex::Amr::computeOptimalSubcycling ( int  n,
int best,
const Real dt_max,
const Real est_work,
const int cycle_max 
)
static

Compute the optimal subcycling pattern. This assumes that anything less than cycle_max[i] is a valid number of subcycles at level[i]. For example if ref_ratio[i] = cycle_max[i] = 4, then 1,2,3,4 are all valid values for n_cycles[i].

Parameters
nNumber of refinement levels under consideration.
bestOutput array receiving cycle counts per level (size n).
dt_maxMaximum dt allowed at each level.
est_workEstimated per-level work weights (same length as dt_max).
cycle_maxMaximum allowed subcycles per level (matches n).
Returns
The coarse-level dt associated with the best pattern.

◆ cumTime()

Real amrex::Amr::cumTime ( ) const
inlinenoexcept

Physical time (time advanced so far).

◆ DataLog()

std::ostream & amrex::Amr::DataLog ( int  i)

Return the i th datalog stream (opened lazily) for writing diagnostics.

Parameters
iDatalog index.

◆ DataLogName()

std::string amrex::Amr::DataLogName ( int  i) const
inlinenoexcept

The filename backing the i th datalog stream.

◆ defBaseLevel()

void amrex::Amr::defBaseLevel ( Real  strt_time,
const BoxArray lev0_grids = nullptr,
const Vector< int > *  pmap = nullptr 
)
protected

Define and initialize the coarsest level using optional level-0 grids lev0_grids and map pmap.

◆ deleteDerivePlotVar()

void amrex::Amr::deleteDerivePlotVar ( const std::string &  name)
static

Remove name from the full derived list.

Parameters
nameDerived quantity identifier.

◆ deleteDeriveSmallPlotVar()

void amrex::Amr::deleteDeriveSmallPlotVar ( const std::string &  name)
static

Remove name from the small derived list.

Parameters
nameDerived quantity identifier.

◆ deleteStatePlotVar()

void amrex::Amr::deleteStatePlotVar ( const std::string &  name)
static

Remove name from the full-plot list if present.

Parameters
nameState variable identifier.

◆ derive() [1/3]

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

Compute derived data on all levels and return as owning MultiFabs.

Parameters
nameIdentifier registered with AmrLevel::derive().
timePhysical time at which the derived results should be evaluated.
ngrowNumber of ghost cells to fill (applied uniformly to all levels).
Returns
Vector of MultiFabs, slot lev containing the data for that level.

◆ derive() [2/3]

void amrex::Amr::derive ( const std::string &  name,
Real  time,
const Vector< MultiFab * > &  mf,
int  dcomp 
)

Compute derived data on all levels and store it in the supplied MultiFabs.

Parameters
nameIdentifier registered with AmrLevel::derive().
timePhysical time at which to evaluate the derived quantity.
mfPer-level buffers that will receive the results.
dcompDestination component inside each MultiFab mf entry.

◆ derive() [3/3]

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

Build a MultiFab containing the derived quantity named name.

Parameters
nameIdentifier registered with AmrLevel::derive().
timePhysical time at which the quantity should be evaluated.
levAMR level on which to compute.
ngrowNumber of ghost cells to fill in the new MultiFab.
Returns
Owning pointer to a MultiFab populated with the requested data.

◆ derivePlotVars()

static const std::list< std::string > & amrex::Amr::derivePlotVars ( )
inlinestaticnoexcept

The names of derived variables to output in the plotfile. They can be set using the amr.derive_plot_vars variable in a ParmParse inputs file.

◆ deriveSmallPlotVars()

static const std::list< std::string > & amrex::Amr::deriveSmallPlotVars ( )
inlinestaticnoexcept

List of derived quantities included when writing small plot files.

◆ dtLevel() [1/2]

const Vector< Real > & amrex::Amr::dtLevel ( ) const
inlinenoexcept

Vector of time steps (dt) for all levels.

◆ dtLevel() [2/2]

Real amrex::Amr::dtLevel ( int  level) const
inlinenoexcept

Time step currently stored for level level.

◆ dtMin()

Real amrex::Amr::dtMin ( int  level) const
inlinenoexcept

Max time step (typically physics-limited) for level level.

◆ ErrorEst()

void amrex::Amr::ErrorEst ( int  lev,
TagBoxArray tags,
Real  time,
int  ngrow 
)
overrideprotectedvirtual

Default error-estimation driver; forwards to the levelbld-defined callbacks.

Parameters
levLevel being tagged.
tagsTag arrays populated by derived routines.
timePhysical time passed to tagging functors.
ngrowNumber of ghost cells available in tags.

Implements amrex::AmrCore.

◆ fillDerivePlotVarList()

void amrex::Amr::fillDerivePlotVarList ( )
static

Fill the list of derive_plot_vars with all derived quantities.

◆ fillDeriveSmallPlotVarList()

void amrex::Amr::fillDeriveSmallPlotVarList ( )
static

Populate derive_small_plot_vars with every registered derived quantity.

◆ fillStatePlotVarList()

void amrex::Amr::fillStatePlotVarList ( )
static

Fill the list of state_plot_vars with all of the state quantities.

◆ fillStateSmallPlotVarList()

void amrex::Amr::fillStateSmallPlotVarList ( )
static

Populate state_small_plot_vars with every state quantity registered on AmrLevel 0.

◆ Finalize()

void amrex::Amr::Finalize ( )
static

Release global Amr static resources allocated by Initialize().

◆ FinalizeInit()

void amrex::Amr::FinalizeInit ( Real  strt_time,
Real  stop_time 
)

Second half of initial initialization; pairs with InitializeInit() using times strt_time and stop_time.

◆ finalizeInSitu()

int amrex::Amr::finalizeInSitu ( )
staticprotected

Tear down in-situ infrastructure; returns 0 on success.

◆ getAmrLevels()

Vector< std::unique_ptr< AmrLevel > > & amrex::Amr::getAmrLevels ( )
noexcept

Return the container of level objects (index 0..finest).

◆ GetAreaNotToTag()

BoxArray amrex::Amr::GetAreaNotToTag ( int  lev)
overrideprotectedvirtual

Return cached regions that must not be tagged on level lev.

Reimplemented from amrex::AmrMesh.

◆ getInitialBA()

const Vector< BoxArray > & amrex::Amr::getInitialBA ( )
staticnoexcept
Returns
Collection of BoxArrays that define the initial grid hierarchy.

◆ getIntersectHiX()

BoundaryPointList & amrex::Amr::getIntersectHiX ( )
inlinenoexcept
Returns
Intersection points between the boundary and the high-X plane.

◆ getIntersectHiY()

BoundaryPointList & amrex::Amr::getIntersectHiY ( )
inlinenoexcept
Returns
Intersection points between the boundary and the high-Y plane.

◆ getIntersectHiZ()

BoundaryPointList & amrex::Amr::getIntersectHiZ ( )
inlinenoexcept
Returns
Intersection points between the boundary and the high-Z plane.

◆ getIntersectLoX()

BoundaryPointList & amrex::Amr::getIntersectLoX ( )
inlinenoexcept
Returns
Intersection points between the domain boundary and the low-X plane.

◆ getIntersectLoY()

BoundaryPointList & amrex::Amr::getIntersectLoY ( )
inlinenoexcept
Returns
Intersection points between the boundary and the low-Y plane.

◆ getIntersectLoZ()

BoundaryPointList & amrex::Amr::getIntersectLoZ ( )
inlinenoexcept
Returns
Intersection points between the boundary and the low-Z plane.

◆ getLevel()

AmrLevel & amrex::Amr::getLevel ( int  lev)
inlinenoexcept

Access the AmrLevel object stored at index lev.

Parameters
levLevel index (0 = coarsest).

◆ grid_places()

void amrex::Amr::grid_places ( int  lbase,
Real  time,
int new_finest,
Vector< BoxArray > &  new_grids 
)
protected

Define new grid locations (called from regrid) and store them in new_grids.

Parameters
lbaseCoarsest level that may change.
timePhysical time used when tagging/refining.
new_finestOutput: index of the resulting finest level.
new_gridsOutput: BoxArrays describing each new level.

◆ init()

void amrex::Amr::init ( Real  strt_time,
Real  stop_time 
)
virtual

Initialize hierarchy state (grid creation, initial conditions, etc.).

Parameters
strt_timeSimulation start time.
stop_timeOptional stop time for the initial solve (used by restart logic).

◆ InitAmr()

void amrex::Amr::InitAmr ( )

Initialize Amr-wide defaults and parse the amr.* ParmParse options.

This drives LevelBld::variableSetUp(), resets bookkeeping (plot/checkpoint cadence, logging, etc.), and seeds shared resources before grids are created.

◆ initialBa()

static const BoxArray & amrex::Amr::initialBa ( int  level)
inlinestaticnoexcept

Return the user-specified BoxArray that initializes level level (1-based).

Parameters
levelLogical level index (1-based).

◆ initialBaLevels()

static int amrex::Amr::initialBaLevels ( )
inlinestaticnoexcept
Returns
Number of levels for which grids were explicitly provided via inputs.

◆ initialInit()

void amrex::Amr::initialInit ( Real  strt_time,
Real  stop_time,
const BoxArray lev0_grids = nullptr,
const Vector< int > *  pmap = nullptr 
)
protected

Initialize the grid hierarchy (called by Amr::init()).

Parameters
strt_timeSimulation start time.
stop_timeOptional time limit for the initialization phase.
lev0_gridsOptional user-defined BoxArray for level 0.
pmapOptional processor map that pairs with lev0_grids.

◆ Initialize()

void amrex::Amr::Initialize ( )
static

Initialize global, Amr-wide static state (e.g., ParmParse defaults, logs).

◆ InitializeInit()

void amrex::Amr::InitializeInit ( Real  strt_time,
Real  stop_time,
const BoxArray lev0_grids = nullptr,
const Vector< int > *  pmap = nullptr 
)

First half of the two-phase initialization used by initialInit().

Parameters
strt_timeSimulation start time.
stop_timeOptional stop time for the initial solve.
lev0_gridsOptional user-supplied level-0 grids.
pmapOptional processor map matching lev0_grids.

◆ initInSitu()

int amrex::Amr::initInSitu ( )
staticprotected

Set up in-situ analysis hooks; returns 0 on success.

◆ initPltAndChk()

void amrex::Amr::initPltAndChk ( )
protected

Initialize plot/checkpoint cadence, names, and counters from inputs.

◆ initSubcycle()

void amrex::Amr::initSubcycle ( )
protected

Parse ParmParse state and initialize subcycling-related arrays/flags.

◆ InstallNewDistributionMap()

void amrex::Amr::InstallNewDistributionMap ( int  lev,
const DistributionMapping newdm 
)

Replace the distribution map on level lev with newdm and update metadata.

Parameters
levLevel index whose data will be remapped.
newdmDistributionMapping describing the new ownership of level lev grids.

◆ isDerivePlotVar()

bool amrex::Amr::isDerivePlotVar ( const std::string &  name)
staticnoexcept

Return true if name is already registered for full derived-output.

Parameters
nameDerived quantity identifier.

◆ isDeriveSmallPlotVar()

bool amrex::Amr::isDeriveSmallPlotVar ( const std::string &  name)
staticnoexcept

Return true if name is already registered for small derived-output.

Parameters
nameDerived quantity identifier.

◆ isStatePlotVar()

bool amrex::Amr::isStatePlotVar ( const std::string &  name)
static

Return true if name is already registered for full plotfiles.

Parameters
nameState variable identifier.

◆ isStateSmallPlotVar()

bool amrex::Amr::isStateSmallPlotVar ( const std::string &  name)
static

Return true if name is already registered for small plotfiles.

Parameters
nameState variable identifier.

◆ level_being_advanced()

int amrex::Amr::level_being_advanced ( ) const
inlinenoexcept

What is "level" in Amr::timeStep? This is only relevant if we are still in Amr::timeStep; it is set back to -1 on leaving Amr::timeStep.

◆ levelCount()

int amrex::Amr::levelCount ( int  lev) const
inlinenoexcept

Which cycle number are we on for level lev?

◆ levelSteps()

int amrex::Amr::levelSteps ( int  lev) const
inlinenoexcept

Number of time steps completed so far on level lev.

◆ LoadBalanceLevel0()

void amrex::Amr::LoadBalanceLevel0 ( Real  time)
protected

Compute and install a new distribution map for level 0.

Parameters
timeSimulation time at which load-balancing decisions are made.

◆ makeLoadBalanceDistributionMap()

DistributionMapping amrex::Amr::makeLoadBalanceDistributionMap ( int  lev,
Real  time,
const BoxArray ba 
) const
protected

Build a DistributionMapping for level lev at time time using grid layout ba.

Returns
Load-balanced DistributionMapping suitable for the supplied BoxArray.

◆ MakeNewLevelFromCoarse()

void amrex::Amr::MakeNewLevelFromCoarse ( int  lev,
Real  time,
const BoxArray ba,
const DistributionMapping dm 
)
inlineoverrideprotectedvirtual

Make a new level using provided metadata and populate it from the next coarser level.

Parameters
levLevel index being created.
timeSimulation time when the level is initialized.
baBoxArray describing the new grids.
dmDistribution mapping describing ownership of those grids.

Implements amrex::AmrCore.

◆ MakeNewLevelFromScratch()

void amrex::Amr::MakeNewLevelFromScratch ( int  lev,
Real  time,
const BoxArray ba,
const DistributionMapping dm 
)
inlineoverrideprotectedvirtual

Create and fill a new level from scratch.

Only invoked during initialization; derived classes must allocate and populate all data for level lev using the supplied metadata.

Implements amrex::AmrCore.

◆ ManualTagsPlacement()

void amrex::Amr::ManualTagsPlacement ( int  lev,
TagBoxArray tags,
const Vector< IntVect > &  bf_lev 
)
overrideprotectedvirtual

Allow derived classes to override tagging manually on level lev.

Parameters
levLevel being tagged.
tagsTag arrays that may be modified in-place.
bf_levBuffer sizes supplied by Amr::regrid().

Reimplemented from amrex::AmrMesh.

◆ nCycle()

int amrex::Amr::nCycle ( int  level) const
inlinenoexcept

Number of subcycled time steps taken per coarse step on level level.

◆ NumDataLogs()

int amrex::Amr::NumDataLogs ( )
noexcept
Returns
Number of datalog streams that have been opened so far.

◆ numGrids() [1/2]

int amrex::Amr::numGrids ( )
noexcept
Returns
Total number of grids currently allocated across all levels.

◆ numGrids() [2/2]

int amrex::Amr::numGrids ( int  lev)
noexcept
Returns
Number of grids at level lev.

◆ okToContinue()

int amrex::Amr::okToContinue ( )
noexcept
Returns
Nonzero if the driver should continue advancing.

◆ okToRegrid()

bool amrex::Amr::okToRegrid ( int  level)
noexcept

Return true if level level is permitted to regrid on this step.

Parameters
levelLevel index being queried.

◆ operator=() [1/2]

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

Move assignment is disabled to keep ownership contracts simple.

◆ operator=() [2/2]

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

Copy assignment is likewise disabled.

◆ Plot_Files_Output()

bool amrex::Amr::Plot_Files_Output ( )
static

Write out plotfiles (True/False)?

◆ plotInt()

int amrex::Amr::plotInt ( ) const
inlinenoexcept

Number of time steps between plot files.

◆ plotLogPer()

Real amrex::Amr::plotLogPer ( ) const
inlinenoexcept

Spacing in log10(time) of logarithmically spaced plot files.

◆ plotMaxLevel()

int amrex::Amr::plotMaxLevel ( ) const
inlinenoexcept

Maximum level index to include when writing plot files.

◆ plotPer()

Real amrex::Amr::plotPer ( ) const
inlinenoexcept

Time between plot files.

◆ printGridInfo()

void amrex::Amr::printGridInfo ( std::ostream &  os,
int  min_lev,
int  max_lev 
)
protected

Print per-level grid metadata to os for levels [min_lev, max_lev].

Parameters
osOutput stream receiving the report.
min_levSmallest AMR level to describe.
max_levLargest AMR level to describe.

◆ readProbinFile()

void amrex::Amr::readProbinFile ( int init)
protected

Read the probin file, updating the flag init that tracks initialization iterations.

◆ RedistributeParticles()

void amrex::Amr::RedistributeParticles ( )

Trigger a particle redistribution pass across all levels.

Mirrors AmrLevel::particle_redistribute() but drives the full hierarchy.

◆ regrid()

void amrex::Amr::regrid ( int  lbase,
Real  time,
bool  initial = false 
)
overridevirtual

Rebuild grid hierarchy finer than lbase.

Parameters
lbaseBase level (levels <= lbase are unchanged).
timePhysical time used when tagging/refining.
initialWhether this call happens during initialization (affects logging and I/O).

Reimplemented from amrex::AmrCore.

◆ regrid_level_0_on_restart()

void amrex::Amr::regrid_level_0_on_restart ( )
protectedvirtual

Regrid level 0 on restart.

◆ regridInt()

int amrex::Amr::regridInt ( int  lev) const
inlinenoexcept

Interval (in steps) between regrids on level lev.

◆ RegridOnly()

void amrex::Amr::RegridOnly ( Real  time,
bool  do_io = true 
)

Regrid the hierarchy at physical time time without advancing any levels.

Parameters
timePhysical time used during tagging/refinement.
do_ioWhen true, emit log messages describing the regrid.

◆ RegridOnRestart()

bool amrex::Amr::RegridOnRestart ( )
staticnoexcept

Whether to regrid right after restart.

◆ RemakeLevel()

void amrex::Amr::RemakeLevel ( int  lev,
Real  time,
const BoxArray ba,
const DistributionMapping dm 
)
inlineoverrideprotectedvirtual

Remake an existing level using provided metadata and refill from old data.

Parameters
levLevel index being regenerated.
timeSimulation time to stamp on the data.
baNew BoxArray layout for that level.
dmDistribution mapping describing ownership.

Implements amrex::AmrCore.

◆ restart()

void amrex::Amr::restart ( const std::string &  filename)
protected

Restart from a checkpoint file named filename.

◆ setBoundaryGeometry() [1/2]

void amrex::Amr::setBoundaryGeometry ( BoundaryPointList &  IntersectLoX,
BoundaryPointList &  IntersectHiX,
BoundaryPointList &  IntersectLoY,
BoundaryPointList &  IntersectHiY 
)
inlinenoexcept

Specialized version: register intersection points for faces orthogonal to X/Y.

Define BoundaryPointLists that give the intersections of the external geometry with constant (i,k) and (j,k) planes on the coarsest level.

Parameters
IntersectLoXIntersection list for the low-X face.
IntersectHiXIntersection list for the high-X face.
IntersectLoYIntersection list for the low-Y face.
IntersectHiYIntersection list for the high-Y face.

◆ setBoundaryGeometry() [2/2]

void amrex::Amr::setBoundaryGeometry ( BoundaryPointList &  IntersectLoX,
BoundaryPointList &  IntersectHiX,
BoundaryPointList &  IntersectLoY,
BoundaryPointList &  IntersectHiY,
BoundaryPointList &  IntersectLoZ,
BoundaryPointList &  IntersectHiZ 
)
inlinenoexcept

General version: register intersections for all Cartesian faces.

Define BoundaryPointLists that give the intersections of the external geometry with constant (i,k), (j,k), and (i,j) planes on the coarsest level.

Parameters
IntersectLoXIntersection list for the low-X face.
IntersectHiXIntersection list for the high-X face.
IntersectLoYIntersection list for the low-Y face.
IntersectHiYIntersection list for the high-Y face.
IntersectLoZIntersection list for the low-Z face.
IntersectHiZIntersection list for the high-Z face.

◆ setComputeNewDtOnRegrid()

static void amrex::Amr::setComputeNewDtOnRegrid ( bool  flag)
inlinestatic

Toggle whether new time steps are recomputed immediately after each regrid (flag=true recomputes).

◆ setCumTime()

void amrex::Amr::setCumTime ( Real  t)
inlinenoexcept

Set the accumulated physical time to t.

◆ setDtLevel() [1/2]

void amrex::Amr::setDtLevel ( const Vector< Real > &  dt_lev)
noexcept

Set per-level time steps to the values in dt_lev.

◆ setDtLevel() [2/2]

void amrex::Amr::setDtLevel ( Real  dt,
int  lev 
)
noexcept

Set the time step stored at AMR level lev to dt.

◆ setDtMin()

void amrex::Amr::setDtMin ( const Vector< Real > &  dt_min_in)
noexcept

Override the maximum allowable dt (typically CFL-limited) for each level using dt_min_in.

◆ setLevelCount()

void amrex::Amr::setLevelCount ( int  lev,
int  n 
)
inlinenoexcept

Set the current cycle count for level lev.

Parameters
levLevel index to update.
nNew cycle number.

◆ setLevelSteps()

void amrex::Amr::setLevelSteps ( int  lev,
int  n 
)
inlinenoexcept

Override the recorded time-step count for level lev.

Parameters
levLevel index to update.
nNew time-step count.

◆ setNCycle()

void amrex::Amr::setNCycle ( const Vector< int > &  ns)
noexcept

Provide the number of time subcycles per level via ns.

◆ setRecordDataInfo()

void amrex::Amr::setRecordDataInfo ( int  i,
const std::string &  filename 
)
protected

Configure the i th datalog stream to mirror file filename.

◆ setRecordGridInfo()

void amrex::Amr::setRecordGridInfo ( const std::string &  filename)
protected

Record the log filename filename used when dumping grid info metadata.

◆ setRecordRunInfo()

void amrex::Amr::setRecordRunInfo ( const std::string &  filename)
protected

Record the log filename filename used when writing detailed run information.

◆ setRecordRunInfoTerse()

void amrex::Amr::setRecordRunInfoTerse ( const std::string &  filename)
protected

Record the log filename filename used for terse run information.

◆ setStartTime()

void amrex::Amr::setStartTime ( Real  t)
inlinenoexcept

Set the recorded start time to t.

◆ smallplotInt()

int amrex::Amr::smallplotInt ( ) const
inlinenoexcept

Number of time steps between small plot files.

◆ smallplotLogPer()

Real amrex::Amr::smallplotLogPer ( ) const
inlinenoexcept

Spacing in log10(time) of logarithmically spaced small plot files.

◆ smallplotPer()

Real amrex::Amr::smallplotPer ( ) const
inlinenoexcept

Time between plot files.

◆ startTime()

Real amrex::Amr::startTime ( ) const
inlinenoexcept

Physical time this simulation started.

◆ statePlotVars()

static const std::list< std::string > & amrex::Amr::statePlotVars ( )
inlinestaticnoexcept

The names of state variables to output in the plotfile. They can be set using the amr.plot_vars variable in a ParmParse inputs file.

◆ stateSmallPlotVars()

static const std::list< std::string > & amrex::Amr::stateSmallPlotVars ( )
inlinestaticnoexcept

List of state variables written to the reduced “small” plot files.

◆ stepOfLastCheckPoint()

int amrex::Amr::stepOfLastCheckPoint ( ) const
inlinenoexcept

Return the AMR step index associated with the most recent checkpoint.

◆ stepOfLastPlotFile()

int amrex::Amr::stepOfLastPlotFile ( ) const
inlinenoexcept

Return the AMR step index associated with the most recent plot file.

◆ stepOfLastSmallPlotFile()

int amrex::Amr::stepOfLastSmallPlotFile ( ) const
inlinenoexcept

Return the AMR step index associated with the most recent small plot file.

◆ subCycle()

int amrex::Amr::subCycle ( ) const
inlinenoexcept

Subcycle in time?

◆ subcyclingMode()

const std::string & amrex::Amr::subcyclingMode ( ) const
inlinenoexcept

How are we subcycling?

◆ theRestartFile()

const std::string & amrex::Amr::theRestartFile ( ) const
inlinenoexcept

Name of the most recent restart checkpoint file.

◆ theRestartPlotFile()

const std::string & amrex::Amr::theRestartPlotFile ( ) const
inlinenoexcept

Name of the restart plotfile.

◆ timeStep()

void amrex::Amr::timeStep ( int  level,
Real  time,
int  iteration,
int  niter,
Real  stop_time 
)
protectedvirtual

Do a single time step on level level.

Parameters
levelLevel index to advance.
timePhysical time at the beginning of the step.
iterationWhich sub-iteration of this coarse step we are in.
niterTotal number of iterations that will run on this level.
stop_timeOptional physical time limit that halts advancement when reached.

◆ updateInSitu()

int amrex::Amr::updateInSitu ( )
protected

Drive any per-step in-situ actions; returns 0 on success.

◆ UsingPrecreateDirectories()

bool amrex::Amr::UsingPrecreateDirectories ( )
staticnoexcept

Whether plot/checkpoint writers pre-create directories (helps avoid races on some filesystems).

◆ writePlotNow()

bool amrex::Amr::writePlotNow ( )
protectednoexcept

Whether to write a plotfile now.

◆ writeSmallPlotFile()

void amrex::Amr::writeSmallPlotFile ( )
virtual

Write the reduced “small plot” file (subset of variables/levels) to disk.

◆ writeSmallPlotNow()

bool amrex::Amr::writeSmallPlotNow ( )
protectednoexcept

Whether to emit a small plot file on this step.

Member Data Documentation

◆ abort_on_stream_retry_failure

bool amrex::Amr::abort_on_stream_retry_failure
protected

◆ amr_level

Vector<std::unique_ptr<AmrLevel> > amrex::Amr::amr_level
protected

Vector of levels.

◆ bUserStopRequest

bool amrex::Amr::bUserStopRequest
protected

◆ check_file_root

std::string amrex::Amr::check_file_root
protected

Root name of checkpoint file.

◆ check_int

int amrex::Amr::check_int
protected

How often checkpoint (# time steps).

◆ check_per

Real amrex::Amr::check_per
protected

How often checkpoint (units of time).

◆ compute_new_dt_on_regrid

bool amrex::Amr::compute_new_dt_on_regrid
staticprotected

◆ cumtime

Real amrex::Amr::cumtime = std::numeric_limits<Real>::lowest()
protected

Physical time variable.

◆ datalog

Vector<std::unique_ptr<std::fstream> > amrex::Amr::datalog
protected

◆ datalogname

Vector<std::string> amrex::Amr::datalogname
protected

◆ derive_plot_vars

std::list< std::string > amrex::Amr::derive_plot_vars
staticprotected

Derived Vars to dump to plotfile.

◆ derive_small_plot_vars

std::list< std::string > amrex::Amr::derive_small_plot_vars
staticprotected

Derived Vars to dump to small plotfile.

◆ dt_level

Vector<Real> amrex::Amr::dt_level
protected

Timestep at this level.

◆ dt_min

Vector<Real> amrex::Amr::dt_min
protected

◆ file_name_digits

int amrex::Amr::file_name_digits
protected

How many digits to use in the plotfile and checkpoint names.

◆ first_plotfile

bool amrex::Amr::first_plotfile
staticprotected

◆ first_smallplotfile

bool amrex::Amr::first_smallplotfile
static

◆ gridlog

std::ofstream amrex::Amr::gridlog
protected

◆ initial_ba

Vector< BoxArray > amrex::Amr::initial_ba
staticprotected

Array of BoxArrays read in to initially define grid hierarchy.

◆ initial_grids_file

std::string amrex::Amr::initial_grids_file
protected

Grids file that will bypass regridding only at initialization.

◆ intersect_hix

BoundaryPointList amrex::Amr::intersect_hix

◆ intersect_hiy

BoundaryPointList amrex::Amr::intersect_hiy

◆ intersect_hiz

BoundaryPointList amrex::Amr::intersect_hiz

◆ intersect_lox

BoundaryPointList amrex::Amr::intersect_lox

◆ intersect_loy

BoundaryPointList amrex::Amr::intersect_loy

◆ intersect_loz

BoundaryPointList amrex::Amr::intersect_loz

◆ last_checkpoint

int amrex::Amr::last_checkpoint
protected

Step number of previous checkpoint.

◆ last_plotfile

int amrex::Amr::last_plotfile
protected

Step number of previous plotfile.

◆ last_smallplotfile

int amrex::Amr::last_smallplotfile
protected

Step number of previous small plotfile.

◆ level_count

Vector<int> amrex::Amr::level_count
protected

◆ level_steps

Vector<int> amrex::Amr::level_steps
protected

Number of time steps at this level.

◆ levelbld

LevelBld* amrex::Amr::levelbld
protected

◆ loadbalance_level0_int

int amrex::Amr::loadbalance_level0_int
protected

◆ loadbalance_max_fac

Real amrex::Amr::loadbalance_max_fac
protected

◆ loadbalance_with_workestimates

int amrex::Amr::loadbalance_with_workestimates
protected

◆ message_int

int amrex::Amr::message_int
protected

How often checking messages touched by user, such as "stop_run".

◆ n_cycle

Vector<int> amrex::Amr::n_cycle
protected

◆ plot_file_root

std::string amrex::Amr::plot_file_root
protected

Root name of plotfile.

◆ plot_int

int amrex::Amr::plot_int
protected

How often plotfile (# of time steps)

◆ plot_log_per

Real amrex::Amr::plot_log_per
protected

How often plotfile (in units of log10(time))

◆ plot_max_level

int amrex::Amr::plot_max_level
protected

Maximum AMR level to write to a plotfile.

◆ plot_per

Real amrex::Amr::plot_per
protected

How often plotfile (in units of time)

◆ probin_file

std::string amrex::Amr::probin_file
protected

◆ record_grid_info

int amrex::Amr::record_grid_info
protected

◆ record_run_info

int amrex::Amr::record_run_info
protected

◆ record_run_info_terse

int amrex::Amr::record_run_info_terse
protected

◆ regrid_ba

Vector< BoxArray > amrex::Amr::regrid_ba
staticprotected

Array of BoxArrays read in to externally define grid hierarchy at each regrid.

◆ regrid_grids_file

std::string amrex::Amr::regrid_grids_file
protected

Grids file that will bypass regridding.

◆ regrid_int

Vector<int> amrex::Amr::regrid_int
protected

Interval between regridding.

◆ restart_chkfile

std::string amrex::Amr::restart_chkfile
protected

◆ restart_pltfile

std::string amrex::Amr::restart_pltfile
protected

◆ runlog

std::ofstream amrex::Amr::runlog
protected

◆ runlog_terse

std::ofstream amrex::Amr::runlog_terse
protected

◆ small_plot_file_root

std::string amrex::Amr::small_plot_file_root
protected

Root name of small plotfile.

◆ small_plot_int

int amrex::Amr::small_plot_int
protected

How often small plotfile (# of time steps)

◆ small_plot_log_per

Real amrex::Amr::small_plot_log_per
protected

How often small plotfile (in units of log10(time))

◆ small_plot_per

Real amrex::Amr::small_plot_per
protected

How often small plotfile (in units of time)

◆ start_time

Real amrex::Amr::start_time = std::numeric_limits<Real>::lowest()
protected

Physical time this simulation started.

◆ state_plot_vars

std::list< std::string > amrex::Amr::state_plot_vars
staticprotected

State Vars to dump to plotfile.

◆ state_small_plot_vars

std::list< std::string > amrex::Amr::state_small_plot_vars
staticprotected

State Vars to dump to small plotfile.

◆ stream_max_tries

int amrex::Amr::stream_max_tries
protected

◆ sub_cycle

int amrex::Amr::sub_cycle
protected

◆ subcycling_mode

std::string amrex::Amr::subcycling_mode
protected

Type of subcycling to use.

◆ which_level_being_advanced

int amrex::Amr::which_level_being_advanced = -1
protected

Only >=0 if we are in Amr::timeStep(level,...)

◆ write_plotfile_with_checkpoint

bool amrex::Amr::write_plotfile_with_checkpoint
protected

Write out a plotfile whenever we checkpoint.


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