Block-Structured AMR Software Framework
AMReX_Amr.H
Go to the documentation of this file.
1 
2 #ifndef AMREX_Amr_H_
3 #define AMREX_Amr_H_
4 #include <AMReX_Config.H>
5 
6 #include <AMReX_Box.H>
7 #include <AMReX_Geometry.H>
8 #include <AMReX_BoxArray.H>
9 #include <AMReX_Array.H>
10 #include <AMReX_Vector.H>
11 #include <AMReX_BCRec.H>
12 #include <AMReX_AmrCore.H>
13 
14 #include <iosfwd>
15 #include <list>
16 #include <memory>
17 
18 namespace amrex {
19 
20 class AmrLevel;
21 class LevelBld;
22 #if defined(AMREX_USE_SENSEI_INSITU) && !defined(AMREX_NO_SENSEI_AMR_INST)
23 class AmrInSituBridge;
24 #endif
25 
33 class Amr
34  : public AmrCore
35 {
36  using BoundaryPointList = std::multimap< std::pair<int, int>, double >;
37 
38 public:
40  Amr (LevelBld* a_levelbld /* One must pass LevelBld* as an argument now*/);
41 
42  Amr (const RealBox* rb, int max_level_in, const Vector<int>& n_cell_in, int coord,
43  LevelBld* a_levelbld /* One must pass LevelBld* as an argument now*/);
44 
45  Amr (const Amr& rhs) = delete;
46  Amr (Amr&& rhs) = delete;
47  Amr& operator= (const Amr& rhs) = delete;
48  Amr& operator= (Amr&& rhs) = delete;
49 
50  void InitAmr ();
51 
53  ~Amr () override;
54 
56  virtual void init (Real strt_time, Real stop_time);
57 
59  void InitializeInit (Real strt_time, Real stop_time,
60  const BoxArray* lev0_grids = nullptr, const Vector<int>* pmap = nullptr);
61 
63  void FinalizeInit (Real strt_time, Real stop_time);
64 
66  void setDtLevel (const Vector<Real>& dt_lev) noexcept;
67 
69  void setDtLevel (Real dt, int lev) noexcept;
70 
72  void setDtMin (const Vector<Real>& dt_min_in) noexcept;
73 
75  void setNCycle (const Vector<int>& ns) noexcept;
76 
78  int subCycle () const noexcept { return sub_cycle; }
79 
81  const std::string& subcyclingMode() const noexcept { return subcycling_mode; }
82 
87  int level_being_advanced () const noexcept { return which_level_being_advanced; }
89  Real cumTime () const noexcept { return cumtime; }
90  void setCumTime (Real t) noexcept {cumtime = t;}
92  Real startTime () const noexcept { return start_time; }
93  void setStartTime (Real t) noexcept {start_time = t;}
95  Real dtLevel (int level) const noexcept { return dt_level[level]; }
97  Real dtMin (int level) const noexcept { return dt_min[level]; }
99  const Vector<Real>& dtLevel () const noexcept { return dt_level; }
101  int nCycle (int level) const noexcept { return n_cycle[level]; }
103  int levelSteps (int lev) const noexcept { return level_steps[lev]; }
105  void setLevelSteps (int lev, int n) noexcept { level_steps[lev] = n; }
107  int levelCount (int lev) const noexcept { return level_count[lev]; }
109  void setLevelCount (int lev, int n) noexcept { level_count[lev] = n; }
111  static bool RegridOnRestart () noexcept;
113  int regridInt (int lev) const noexcept { return regrid_int[lev]; }
115  int checkInt () const noexcept { return check_int; }
117  Real checkPer() const noexcept { return check_per; }
119  int plotInt () const noexcept { return plot_int; }
121  Real plotPer () const noexcept { return plot_per; }
123  Real plotLogPer () const noexcept { return plot_log_per; }
125  int plotMaxLevel () const noexcept { return plot_max_level; }
127  int smallplotInt () const noexcept { return small_plot_int; }
129  Real smallplotPer () const noexcept { return small_plot_per; }
131  Real smallplotLogPer () const noexcept { return small_plot_log_per; }
137  static const std::list<std::string>& statePlotVars () noexcept { return state_plot_vars; }
138  static const std::list<std::string>& stateSmallPlotVars () noexcept { return state_small_plot_vars; }
140  static bool isStatePlotVar (const std::string& name);
141  static bool isStateSmallPlotVar (const std::string& name);
146  static void addStatePlotVar (const std::string& name);
147  static void addStateSmallPlotVar (const std::string& name);
149  static void deleteStatePlotVar (const std::string& name);
151  static void clearStatePlotVarList ();
152  static void clearStateSmallPlotVarList ();
154  static void fillStatePlotVarList ();
155  static void fillStateSmallPlotVarList ();
157  static bool Plot_Files_Output ();
163  static const std::list<std::string>& derivePlotVars () noexcept { return derive_plot_vars; }
164  static const std::list<std::string>& deriveSmallPlotVars () noexcept { return derive_small_plot_vars; }
166  static bool isDerivePlotVar (const std::string& name) noexcept;
167  static bool isDeriveSmallPlotVar (const std::string& name) noexcept;
172  static void addDerivePlotVar (const std::string& name);
173  static void addDeriveSmallPlotVar (const std::string& name);
175  static void deleteDerivePlotVar (const std::string& name);
176  static void deleteDeriveSmallPlotVar (const std::string& name);
178  static void clearDerivePlotVarList ();
179  static void clearDeriveSmallPlotVarList ();
181  static void fillDerivePlotVarList ();
182  static void fillDeriveSmallPlotVarList ();
183 
184  static void setComputeNewDtOnRegrid (bool flag) { compute_new_dt_on_regrid = flag; }
185 
186  static void Initialize ();
187  static void Finalize ();
189  AmrLevel& getLevel (int lev) noexcept { return *amr_level[lev]; }
193  Long cellCount () noexcept;
195  Long cellCount (int lev) noexcept;
197  int numGrids () noexcept;
199  int numGrids (int lev) noexcept;
201  int okToContinue () noexcept;
203  void regrid (int lbase,
204  Real time,
205  bool initial = false) override;
207  void RegridOnly (Real time, bool do_io = true);
209  bool okToRegrid (int level) noexcept;
211  static const BoxArray& initialBa (int level) noexcept
212  { BL_ASSERT(level-1 < initial_ba.size()); return initial_ba[level-1]; }
214  static int initialBaLevels () noexcept { return static_cast<int>(initial_ba.size()); }
216  virtual void coarseTimeStep (Real stop_time);
217 
219  Real coarseTimeStepDt (Real stop_time);
221  std::unique_ptr<MultiFab> derive (const std::string& name,
222  Real time,
223  int lev,
224  int ngrow);
226  const std::string& theRestartFile () const noexcept { return restart_chkfile; }
228  const std::string& theRestartPlotFile () const noexcept { return restart_pltfile; }
230  std::ostream& DataLog (int i);
232  std::string DataLogName (int i) const noexcept { return datalogname[i]; }
234  int NumDataLogs () noexcept;
242  static Real computeOptimalSubcycling (int n,
243  int* best,
244  const Real* dt_max,
245  const Real* est_work,
246  const int* cycle_max);
247 
249  virtual void writePlotFile ();
250  int stepOfLastPlotFile () const noexcept {return last_plotfile;}
252  virtual void writeSmallPlotFile ();
253  int stepOfLastSmallPlotFile () const noexcept {return last_smallplotfile;}
255  virtual void checkPoint ();
256  int stepOfLastCheckPoint () const noexcept {return last_checkpoint;}
257 
258  static const Vector<BoxArray>& getInitialBA() noexcept;
259 
267  BoundaryPointList& IntersectHiX,
268  BoundaryPointList& IntersectLoY,
269  BoundaryPointList& IntersectHiY) noexcept
270  {
271  intersect_lox = IntersectLoX;
272  intersect_hix = IntersectHiX;
273  intersect_loy = IntersectLoY;
274  intersect_hiy = IntersectHiY;
275  }
276 
285  BoundaryPointList& IntersectHiX,
286  BoundaryPointList& IntersectLoY,
287  BoundaryPointList& IntersectHiY,
288  BoundaryPointList& IntersectLoZ,
289  BoundaryPointList& IntersectHiZ) noexcept
290  {
291  intersect_lox = IntersectLoX;
292  intersect_hix = IntersectHiX;
293  intersect_loy = IntersectLoY;
294  intersect_hiy = IntersectHiY;
295  intersect_loz = IntersectLoZ;
296  intersect_hiz = IntersectHiZ;
297  }
298 
300  {
301  return intersect_lox;
302  }
304  {
305  return intersect_hix;
306  }
308  {
309  return intersect_loy;
310  }
312  {
313  return intersect_hiy;
314  }
316  {
317  return intersect_loz;
318  }
320  {
321  return intersect_hiz;
322  }
323 
324 #ifdef AMREX_PARTICLES
326  void RedistributeParticles ();
327 #endif
328 
329  void InstallNewDistributionMap (int lev, const DistributionMapping& newdm);
330 
331  static bool UsingPrecreateDirectories () noexcept;
332 
333 protected:
334 
336  void initialInit (Real strt_time, Real stop_time,
337  const BoxArray* lev0_grids = nullptr, const Vector<int>* pmap = nullptr);
338 #ifndef AMREX_NO_PROBINIT
340  void readProbinFile (int& init);
341 #endif
343  void checkInput ();
345  void restart (const std::string& filename);
347  void defBaseLevel (Real strt_time, const BoxArray* lev0_grids = nullptr, const Vector<int>* pmap = nullptr);
349  void bldFineLevels (Real strt_time);
351  virtual void regrid_level_0_on_restart ();
353  void grid_places (int lbase,
354  Real time,
355  int& new_finest,
356  Vector<BoxArray>& new_grids);
357 
358  DistributionMapping makeLoadBalanceDistributionMap (int lev, Real time, const BoxArray& ba) const;
359  void LoadBalanceLevel0 (Real time);
360 
361  void ErrorEst (int lev, TagBoxArray& tags, Real time, int ngrow) override;
362  BoxArray GetAreaNotToTag (int lev) override;
363  void ManualTagsPlacement (int lev, TagBoxArray& tags, const Vector<IntVect>& bf_lev) override;
364 
366  virtual void timeStep (int level,
367  Real time,
368  int iteration,
369  int niter,
370  Real stop_time);
371 
372  // pure virtual function in AmrCore
373  void MakeNewLevelFromScratch (int /*lev*/, Real /*time*/, const BoxArray& /*ba*/, const DistributionMapping& /*dm*/) override
374  { amrex::Abort("How did we get here!"); }
375  void MakeNewLevelFromCoarse (int /*lev*/, Real /*time*/, const BoxArray& /*ba*/, const DistributionMapping& /*dm*/) override
376  { amrex::Abort("How did we get here!"); }
377  void RemakeLevel (int /*lev*/, Real /*time*/, const BoxArray& /*ba*/, const DistributionMapping& /*dm*/) override
378  { amrex::Abort("How did we get here!"); }
379  void ClearLevel (int /*lev*/) override
380  { amrex::Abort("How did we get here!"); }
381 
383  bool writePlotNow () noexcept;
384  bool writeSmallPlotNow () noexcept;
385 
386  void printGridInfo (std::ostream& os,
387  int min_lev,
388  int max_lev);
389 
390  void setRecordGridInfo (const std::string&);
391 
392  void setRecordRunInfo (const std::string&);
393 
394  void setRecordRunInfoTerse (const std::string&);
395 
396  void setRecordDataInfo (int i, const std::string&);
397 
398  void initSubcycle();
399  void initPltAndChk();
400 
401  static int initInSitu();
402  int updateInSitu();
403  static int finalizeInSitu();
404 
405  //
406  // The data ...
407  //
408  std::string regrid_grids_file;
409  std::string initial_grids_file;
410  Vector<std::unique_ptr<AmrLevel> > amr_level;
411  Real cumtime = std::numeric_limits<Real>::lowest();
412  Real start_time = std::numeric_limits<Real>::lowest();
413  Vector<Real> dt_level;
417  std::string subcycling_mode;
418  Vector<Real> dt_min;
421  int check_int;
422  Real check_per;
423  std::string check_file_root;
426  int plot_int;
427  Real plot_per;
436  std::string plot_file_root;
437  std::string small_plot_file_root;
438 
440 
444  std::ofstream gridlog;
445  std::ofstream runlog;
446  std::ofstream runlog_terse;
447  Vector<std::unique_ptr<std::fstream> > datalog;
448  Vector<std::string> datalogname;
450  std::string restart_chkfile;
451  std::string restart_pltfile;
452 #ifndef AMREX_NO_PROBINIT
453  std::string probin_file;
454 #endif
461 
463 
464  //
465  // The static data ...
466  //
467  static std::list<std::string> state_plot_vars;
468  static std::list<std::string> state_small_plot_vars;
469  static std::list<std::string> derive_plot_vars;
470  static std::list<std::string> derive_small_plot_vars;
471  static bool first_plotfile;
477 
478 #if defined(AMREX_USE_SENSEI_INSITU) && !defined(AMREX_NO_SENSEI_AMR_INST)
479  static AmrInSituBridge *insitu_bridge;
480 #endif
481 
482 public:
489 
490  static bool first_smallplotfile;
491 
492 private:
493  void writePlotFileDoit (std::string const& pltfile, bool regular);
494 };
495 
496 }
497 
498 #endif /*_Amr_H_*/
#define BL_ASSERT(EX)
Definition: AMReX_BLassert.H:39
Provide basic functionalities to set up an AMR hierarchy.
Definition: AMReX_AmrCore.H:25
Contains the bridge code for simulations that use amrex::Amr.
Definition: AMReX_AmrInSituBridge.H:14
Virtual base class for managing individual levels. AmrLevel functions both as a container for state d...
Definition: AMReX_AmrLevel.H:38
Manage hierarchy of levels for time-dependent AMR computations.
Definition: AMReX_Amr.H:35
Real smallplotPer() const noexcept
Time between plot files.
Definition: AMReX_Amr.H:129
std::string check_file_root
Root name of checkpoint file.
Definition: AMReX_Amr.H:423
virtual void init(Real strt_time, Real stop_time)
Init data after construction. Must be called before timestepping.
Definition: AMReX_Amr.cpp:1120
std::string probin_file
Definition: AMReX_Amr.H:453
virtual void writeSmallPlotFile()
Write the small plot file to be used for visualization.
Definition: AMReX_Amr.cpp:880
BoundaryPointList & getIntersectHiY() noexcept
Definition: AMReX_Amr.H:311
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....
Definition: AMReX_Amr.H:137
bool bUserStopRequest
Definition: AMReX_Amr.H:462
void RemakeLevel(int, Real, const BoxArray &, const DistributionMapping &) override
Remake an existing level using provided BoxArray and DistributionMapping and fill with existing fine ...
Definition: AMReX_Amr.H:377
static void clearDeriveSmallPlotVarList()
Definition: AMReX_Amr.cpp:686
static void fillDerivePlotVarList()
Fill the list of derive_plot_vars with all derived quantities.
Definition: AMReX_Amr.cpp:650
Vector< int > n_cycle
Definition: AMReX_Amr.H:416
void checkInput()
Check for valid input.
Definition: AMReX_Amr.cpp:1042
void setRecordGridInfo(const std::string &)
Definition: AMReX_Amr.cpp:731
const std::string & theRestartPlotFile() const noexcept
Name of the restart plotfile.
Definition: AMReX_Amr.H:228
DistributionMapping makeLoadBalanceDistributionMap(int lev, Real time, const BoxArray &ba) const
Definition: AMReX_Amr.cpp:2755
void setRecordDataInfo(int i, const std::string &)
Definition: AMReX_Amr.cpp:773
void setNCycle(const Vector< int > &ns) noexcept
Set the cycle count on each level.
Definition: AMReX_Amr.cpp:801
virtual void writePlotFile()
Write the plot file to be used for visualization.
Definition: AMReX_Amr.cpp:842
Vector< int > level_count
Definition: AMReX_Amr.H:415
Vector< Real > dt_level
Timestep at this level.
Definition: AMReX_Amr.H:413
static const Vector< BoxArray > & getInitialBA() noexcept
Definition: AMReX_Amr.cpp:3382
void setLevelSteps(int lev, int n) noexcept
Number of time steps at specified level.
Definition: AMReX_Amr.H:105
int file_name_digits
How many digits to use in the plotfile and checkpoint names.
Definition: AMReX_Amr.H:434
static void deleteDerivePlotVar(const std::string &name)
Remove the string from derive_plot_vars.
Definition: AMReX_Amr.cpp:708
void writePlotFileDoit(std::string const &pltfile, bool regular)
Definition: AMReX_Amr.cpp:919
int sub_cycle
Definition: AMReX_Amr.H:449
void initPltAndChk()
Definition: AMReX_Amr.cpp:3221
int plot_max_level
Maximum AMR level to write to a plotfile.
Definition: AMReX_Amr.H:429
int stream_max_tries
Definition: AMReX_Amr.H:457
int NumDataLogs() noexcept
How many datalogs have been opened.
Definition: AMReX_Amr.cpp:162
int stepOfLastSmallPlotFile() const noexcept
Definition: AMReX_Amr.H:253
virtual void timeStep(int level, Real time, int iteration, int niter, Real stop_time)
Do a single timestep on level L.
Definition: AMReX_Amr.cpp:1921
void MakeNewLevelFromCoarse(int, Real, const BoxArray &, const DistributionMapping &) override
Make a new level using provided BoxArray and DistributionMapping and fill with interpolated coarse le...
Definition: AMReX_Amr.H:375
std::string restart_pltfile
Definition: AMReX_Amr.H:451
int smallplotInt() const noexcept
Number of time steps between small plot files.
Definition: AMReX_Amr.H:127
Vector< std::string > datalogname
Definition: AMReX_Amr.H:448
Real coarseTimeStepDt(Real stop_time)
Do a complete integration cycle and return the coarse dt.
Definition: AMReX_Amr.cpp:2088
Vector< int > regrid_int
Interval between regridding.
Definition: AMReX_Amr.H:419
void restart(const std::string &filename)
Restart from a checkpoint file.
Definition: AMReX_Amr.cpp:1377
int plotMaxLevel() const noexcept
Maximum level to plot.
Definition: AMReX_Amr.H:125
const std::string & subcyclingMode() const noexcept
How are we subcycling?
Definition: AMReX_Amr.H:81
Vector< std::unique_ptr< std::fstream > > datalog
Definition: AMReX_Amr.H:447
BoundaryPointList & getIntersectLoY() noexcept
Definition: AMReX_Amr.H:307
static void clearStatePlotVarList()
Clear the list of state_plot_vars.
Definition: AMReX_Amr.cpp:586
int message_int
How often checking messages touched by user, such as "stop_run".
Definition: AMReX_Amr.H:435
Vector< int > level_steps
Number of time steps at this level.
Definition: AMReX_Amr.H:414
void setDtMin(const Vector< Real > &dt_min_in) noexcept
Set the dtmin on each level.
Definition: AMReX_Amr.cpp:174
Real plotLogPer() const noexcept
Spacing in log10(time) of logarithmically spaced plot files.
Definition: AMReX_Amr.H:123
std::unique_ptr< MultiFab > derive(const std::string &name, Real time, int lev, int ngrow)
Retrieve derived data. User is responsible for deleting pointer.
Definition: AMReX_Amr.cpp:200
static bool compute_new_dt_on_regrid
Definition: AMReX_Amr.H:476
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 ...
Definition: AMReX_Amr.cpp:3329
Vector< Real > dt_min
Definition: AMReX_Amr.H:418
static bool isStateSmallPlotVar(const std::string &name)
Definition: AMReX_Amr.cpp:565
Vector< std::unique_ptr< AmrLevel > > amr_level
Vector of levels.
Definition: AMReX_Amr.H:410
static Vector< BoxArray > initial_ba
Array of BoxArrays read in to initially define grid hierarchy.
Definition: AMReX_Amr.H:473
void printGridInfo(std::ostream &os, int min_lev, int max_lev)
Definition: AMReX_Amr.cpp:2883
BoundaryPointList intersect_hiz
Definition: AMReX_Amr.H:488
int loadbalance_level0_int
Definition: AMReX_Amr.H:459
void InitializeInit(Real strt_time, Real stop_time, const BoxArray *lev0_grids=nullptr, const Vector< int > *pmap=nullptr)
First part of initialInit.
Definition: AMReX_Amr.cpp:1269
static void clearDerivePlotVarList()
Clear the list of derive_plot_vars.
Definition: AMReX_Amr.cpp:680
static std::list< std::string > derive_small_plot_vars
Derived Vars to dump to small plotfile.
Definition: AMReX_Amr.H:470
static bool isDerivePlotVar(const std::string &name) noexcept
Is the string the name of a variable in derive_plot_vars?
Definition: AMReX_Amr.cpp:636
std::string small_plot_file_root
Root name of small plotfile.
Definition: AMReX_Amr.H:437
static int finalizeInSitu()
Definition: AMReX_Amr.cpp:542
int nCycle(int level) const noexcept
Number of subcycled time steps.
Definition: AMReX_Amr.H:101
static void addDerivePlotVar(const std::string &name)
If the string is not the name of a variable in derive_plot_vars, add it to derive_plot_vars.
Definition: AMReX_Amr.cpp:692
void defBaseLevel(Real strt_time, const BoxArray *lev0_grids=nullptr, const Vector< int > *pmap=nullptr)
Define and initialize coarsest level.
Definition: AMReX_Amr.cpp:2509
Real plotPer() const noexcept
Time between plot files.
Definition: AMReX_Amr.H:121
BoundaryPointList & getIntersectLoZ() noexcept
Definition: AMReX_Amr.H:315
virtual void regrid_level_0_on_restart()
Regrid level 0 on restart.
Definition: AMReX_Amr.cpp:2821
AmrLevel & getLevel(int lev) noexcept
AmrLevel lev.
Definition: AMReX_Amr.H:189
static void clearStateSmallPlotVarList()
Definition: AMReX_Amr.cpp:606
int loadbalance_with_workestimates
Definition: AMReX_Amr.H:458
int stepOfLastPlotFile() const noexcept
Definition: AMReX_Amr.H:250
int small_plot_int
How often small plotfile (# of time steps)
Definition: AMReX_Amr.H:430
void RegridOnly(Real time, bool do_io=true)
Regrid only!
Definition: AMReX_Amr.cpp:1886
static std::list< std::string > state_small_plot_vars
State Vars to dump to small plotfile.
Definition: AMReX_Amr.H:468
BoundaryPointList intersect_loy
Definition: AMReX_Amr.H:484
Amr(Amr &&rhs)=delete
static std::list< std::string > state_plot_vars
State Vars to dump to plotfile.
Definition: AMReX_Amr.H:467
std::ofstream runlog_terse
Definition: AMReX_Amr.H:446
static bool isDeriveSmallPlotVar(const std::string &name) noexcept
Definition: AMReX_Amr.cpp:643
int plotInt() const noexcept
Number of time steps between plot files.
Definition: AMReX_Amr.H:119
int regridInt(int lev) const noexcept
Interval between regridding.
Definition: AMReX_Amr.H:113
BoundaryPointList & getIntersectHiZ() noexcept
Definition: AMReX_Amr.H:319
std::string plot_file_root
Root name of plotfile.
Definition: AMReX_Amr.H:436
void setStartTime(Real t) noexcept
Definition: AMReX_Amr.H:93
Vector< std::unique_ptr< AmrLevel > > & getAmrLevels() noexcept
Array of AmrLevels.
Definition: AMReX_Amr.cpp:182
Real startTime() const noexcept
Physical time this simulation started.
Definition: AMReX_Amr.H:92
int which_level_being_advanced
Only >=0 if we are in Amr::timeStep(level,...)
Definition: AMReX_Amr.H:439
std::string subcycling_mode
Type of subcycling to use.
Definition: AMReX_Amr.H:417
BoundaryPointList intersect_lox
Definition: AMReX_Amr.H:483
static void setComputeNewDtOnRegrid(bool flag)
Definition: AMReX_Amr.H:184
std::ostream & DataLog(int i)
The ith datalog file. Do with it what you want.
Definition: AMReX_Amr.cpp:156
bool writeSmallPlotNow() noexcept
Definition: AMReX_Amr.cpp:2436
static std::list< std::string > derive_plot_vars
Derived Vars to dump to plotfile.
Definition: AMReX_Amr.H:469
void ErrorEst(int lev, TagBoxArray &tags, Real time, int ngrow) override
Tag cells for refinement. TagBoxArray tags is built on level lev grids.
Definition: AMReX_Amr.cpp:3029
int levelSteps(int lev) const noexcept
Number of time steps at specified level.
Definition: AMReX_Amr.H:103
void setCumTime(Real t) noexcept
Definition: AMReX_Amr.H:90
int stepOfLastCheckPoint() const noexcept
Definition: AMReX_Amr.H:256
Amr(const Amr &rhs)=delete
static const std::list< std::string > & stateSmallPlotVars() noexcept
Definition: AMReX_Amr.H:138
Amr(LevelBld *a_levelbld)
The constructor.
Definition: AMReX_Amr.cpp:208
int updateInSitu()
Definition: AMReX_Amr.cpp:529
virtual void checkPoint()
Write current state into a chk* file.
Definition: AMReX_Amr.cpp:1701
int last_plotfile
Step number of previous plotfile.
Definition: AMReX_Amr.H:424
~Amr() override
The destructor.
Definition: AMReX_Amr.cpp:723
void ClearLevel(int) override
Delete level data.
Definition: AMReX_Amr.H:379
static bool isStatePlotVar(const std::string &name)
Is the string the name of a variable in state_plot_vars?
Definition: AMReX_Amr.cpp:558
int plot_int
How often plotfile (# of time steps)
Definition: AMReX_Amr.H:426
std::ofstream gridlog
Definition: AMReX_Amr.H:444
void grid_places(int lbase, Real time, int &new_finest, Vector< BoxArray > &new_grids)
Define new grid locations (called from regrid) and put into new_grids.
Definition: AMReX_Amr.cpp:2928
const Vector< Real > & dtLevel() const noexcept
Array of time steps at all levels.
Definition: AMReX_Amr.H:99
void setBoundaryGeometry(BoundaryPointList &IntersectLoX, BoundaryPointList &IntersectHiX, BoundaryPointList &IntersectLoY, BoundaryPointList &IntersectHiY) noexcept
Specialized version: Define BoundaryPointLists that give the intersections of the external geometry w...
Definition: AMReX_Amr.H:266
void InitAmr()
Definition: AMReX_Amr.cpp:227
int okToContinue() noexcept
More work to be done?
Definition: AMReX_Amr.cpp:829
void readProbinFile(int &init)
Read the probin file.
Definition: AMReX_Amr.cpp:1159
static int initInSitu()
Definition: AMReX_Amr.cpp:515
bool abort_on_stream_retry_failure
Definition: AMReX_Amr.H:456
LevelBld * levelbld
Definition: AMReX_Amr.H:455
static void deleteDeriveSmallPlotVar(const std::string &name)
Definition: AMReX_Amr.cpp:716
Real cumtime
Physical time variable.
Definition: AMReX_Amr.H:411
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....
Definition: AMReX_Amr.H:163
void setRecordRunInfo(const std::string &)
Definition: AMReX_Amr.cpp:745
void ManualTagsPlacement(int lev, TagBoxArray &tags, const Vector< IntVect > &bf_lev) override
Manually tag. Note that tags is built on level lev grids coarsened by bf_lev[lev].
Definition: AMReX_Amr.cpp:3041
void regrid(int lbase, Real time, bool initial=false) override
Rebuild grid hierarchy finer than lbase.
Definition: AMReX_Amr.cpp:2570
Real check_per
How often checkpoint (units of time).
Definition: AMReX_Amr.H:422
Amr & operator=(const Amr &rhs)=delete
static int initialBaLevels() noexcept
Number of levels at which the grids are initially specified.
Definition: AMReX_Amr.H:214
Real checkPer() const noexcept
Time between checkpoint files.
Definition: AMReX_Amr.H:117
static void Initialize()
Definition: AMReX_Amr.cpp:105
static void fillStatePlotVarList()
Fill the list of state_plot_vars with all of the state quantities.
Definition: AMReX_Amr.cpp:572
BoundaryPointList intersect_hix
Definition: AMReX_Amr.H:486
static bool Plot_Files_Output()
Write out plotfiles (True/False)?
Definition: AMReX_Amr.cpp:153
Real small_plot_per
How often small plotfile (in units of time)
Definition: AMReX_Amr.H:431
std::string regrid_grids_file
Grids file that will bypass regridding.
Definition: AMReX_Amr.H:408
std::multimap< std::pair< int, int >, double > BoundaryPointList
Definition: AMReX_Amr.H:36
int check_int
How often checkpoint (# time steps).
Definition: AMReX_Amr.H:421
Real dtLevel(int level) const noexcept
Time step at specified level.
Definition: AMReX_Amr.H:95
static void fillDeriveSmallPlotVarList()
Definition: AMReX_Amr.cpp:665
static bool first_smallplotfile
Definition: AMReX_Amr.H:490
std::string initial_grids_file
Grids file that will bypass regridding only at initialization.
Definition: AMReX_Amr.H:409
Long cellCount() noexcept
Total number of cells.
Definition: AMReX_Amr.cpp:809
void LoadBalanceLevel0(Real time)
Definition: AMReX_Amr.cpp:2799
int record_run_info_terse
Definition: AMReX_Amr.H:443
virtual void coarseTimeStep(Real stop_time)
Do a complete integration cycle.
Definition: AMReX_Amr.cpp:2095
static void deleteStatePlotVar(const std::string &name)
Remove the string from state_plot_vars.
Definition: AMReX_Amr.cpp:628
Real small_plot_log_per
How often small plotfile (in units of log10(time))
Definition: AMReX_Amr.H:432
void bldFineLevels(Real strt_time)
Define and initialize refined levels.
Definition: AMReX_Amr.cpp:3047
static const std::list< std::string > & deriveSmallPlotVars() noexcept
Definition: AMReX_Amr.H:164
std::ofstream runlog
Definition: AMReX_Amr.H:445
int record_grid_info
Definition: AMReX_Amr.H:441
std::string restart_chkfile
Definition: AMReX_Amr.H:450
Real plot_log_per
How often plotfile (in units of log10(time))
Definition: AMReX_Amr.H:428
int checkInt() const noexcept
Number of time steps between checkpoint files.
Definition: AMReX_Amr.H:115
static void addStatePlotVar(const std::string &name)
If the string is not the name of a variable in state_plot_vars, add it to state_plot_vars.
Definition: AMReX_Amr.cpp:612
int last_smallplotfile
Step number of previous small plotfile.
Definition: AMReX_Amr.H:425
void InstallNewDistributionMap(int lev, const DistributionMapping &newdm)
Definition: AMReX_Amr.cpp:2808
const std::string & theRestartFile() const noexcept
Name of the restart chkpoint file.
Definition: AMReX_Amr.H:226
bool writePlotNow() noexcept
Whether to write a plotfile now.
Definition: AMReX_Amr.cpp:2363
void setLevelCount(int lev, int n) noexcept
Which step are we at for the specified level?
Definition: AMReX_Amr.H:109
void initSubcycle()
Definition: AMReX_Amr.cpp:3119
BoundaryPointList intersect_loz
Definition: AMReX_Amr.H:485
std::string DataLogName(int i) const noexcept
The filename of the ith datalog file.
Definition: AMReX_Amr.H:232
static const BoxArray & initialBa(int level) noexcept
Array of BoxArrays read in to initially define grid hierarchy.
Definition: AMReX_Amr.H:211
bool write_plotfile_with_checkpoint
Write out a plotfile whenever we checkpoint.
Definition: AMReX_Amr.H:433
void FinalizeInit(Real strt_time, Real stop_time)
Second part of initialInit.
Definition: AMReX_Amr.cpp:1301
static void fillStateSmallPlotVarList()
Definition: AMReX_Amr.cpp:592
Real dtMin(int level) const noexcept
Max time step (typically based on physics) at specified level.
Definition: AMReX_Amr.H:97
static Vector< BoxArray > regrid_ba
Array of BoxArrays read in to externally define grid hierarchy at each regrid.
Definition: AMReX_Amr.H:475
int levelCount(int lev) const noexcept
Which step are we at for the specified level?
Definition: AMReX_Amr.H:107
static bool RegridOnRestart() noexcept
Whether to regrid right after restart.
Definition: AMReX_Amr.cpp:168
BoundaryPointList intersect_hiy
Definition: AMReX_Amr.H:487
int numGrids() noexcept
Total number of grids.
Definition: AMReX_Amr.cpp:819
int record_run_info
Definition: AMReX_Amr.H:442
Real start_time
Physical time this simulation started.
Definition: AMReX_Amr.H:412
Real loadbalance_max_fac
Definition: AMReX_Amr.H:460
static void Finalize()
Definition: AMReX_Amr.cpp:141
static void addStateSmallPlotVar(const std::string &name)
Definition: AMReX_Amr.cpp:620
BoundaryPointList & getIntersectHiX() noexcept
Definition: AMReX_Amr.H:303
void setBoundaryGeometry(BoundaryPointList &IntersectLoX, BoundaryPointList &IntersectHiX, BoundaryPointList &IntersectLoY, BoundaryPointList &IntersectHiY, BoundaryPointList &IntersectLoZ, BoundaryPointList &IntersectHiZ) noexcept
More general version: Define BoundaryPointLists that give the intersections of the external geometry ...
Definition: AMReX_Amr.H:284
Real cumTime() const noexcept
Physical time.
Definition: AMReX_Amr.H:89
int subCycle() const noexcept
Subcycle in time?
Definition: AMReX_Amr.H:78
void setDtLevel(const Vector< Real > &dt_lev) noexcept
Set the timestep on each level.
Definition: AMReX_Amr.cpp:787
Real smallplotLogPer() const noexcept
Spacing in log10(time) of logarithmically spaced small plot files.
Definition: AMReX_Amr.H:131
void MakeNewLevelFromScratch(int, Real, const BoxArray &, const DistributionMapping &) override
Definition: AMReX_Amr.H:373
BoundaryPointList & getIntersectLoX() noexcept
Definition: AMReX_Amr.H:299
int last_checkpoint
Step number of previous checkpoint.
Definition: AMReX_Amr.H:420
bool okToRegrid(int level) noexcept
Should we regrid this level?
Definition: AMReX_Amr.cpp:3319
Real plot_per
How often plotfile (in units of time)
Definition: AMReX_Amr.H:427
static bool UsingPrecreateDirectories() noexcept
Definition: AMReX_Amr.cpp:99
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 b...
Definition: AMReX_Amr.H:87
static bool first_plotfile
Definition: AMReX_Amr.H:471
void setRecordRunInfoTerse(const std::string &)
Definition: AMReX_Amr.cpp:759
void initialInit(Real strt_time, Real stop_time, const BoxArray *lev0_grids=nullptr, const Vector< int > *pmap=nullptr)
Initialize grid hierarchy – called by Amr::init.
Definition: AMReX_Amr.cpp:1251
static void addDeriveSmallPlotVar(const std::string &name)
Definition: AMReX_Amr.cpp:700
BoxArray GetAreaNotToTag(int lev) override
Definition: AMReX_Amr.cpp:3035
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:549
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
Builds problem-specific AmrLevels.
Definition: AMReX_LevelBld.H:22
A Box with real dimensions. A RealBox is OK iff volume >= 0.
Definition: AMReX_RealBox.H:21
An array of TagBoxes.
Definition: AMReX_TagBox.H:150
Definition: AMReX_Amr.cpp:49
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition: AMReX.cpp:225