Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
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
18namespace amrex {
19
20class AmrLevel;
21class LevelBld;
22#if defined(AMREX_USE_SENSEI_INSITU) && !defined(AMREX_NO_SENSEI_AMR_INST)
23class AmrInSituBridge;
24#endif
25
33class Amr
34 : public AmrCore
35{
36 using BoundaryPointList = std::multimap< std::pair<int, int>, double >;
37
38public:
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);
227 void derive (const std::string& name,
228 Real time,
229 const Vector<MultiFab*>& mf,
230 int dcomp);
232 Vector<std::unique_ptr<MultiFab>> derive (const std::string& name,
233 amrex::Real time,
234 int ngrow);
236 const std::string& theRestartFile () const noexcept { return restart_chkfile; }
238 const std::string& theRestartPlotFile () const noexcept { return restart_pltfile; }
240 std::ostream& DataLog (int i);
242 std::string DataLogName (int i) const noexcept { return datalogname[i]; }
244 int NumDataLogs () noexcept;
252 static Real computeOptimalSubcycling (int n,
253 int* best,
254 const Real* dt_max,
255 const Real* est_work,
256 const int* cycle_max);
257
259 virtual void writePlotFile ();
260 int stepOfLastPlotFile () const noexcept {return last_plotfile;}
262 virtual void writeSmallPlotFile ();
263 int stepOfLastSmallPlotFile () const noexcept {return last_smallplotfile;}
265 virtual void checkPoint ();
266 int stepOfLastCheckPoint () const noexcept {return last_checkpoint;}
267
268 static const Vector<BoxArray>& getInitialBA() noexcept;
269
277 BoundaryPointList& IntersectHiX,
278 BoundaryPointList& IntersectLoY,
279 BoundaryPointList& IntersectHiY) noexcept
280 {
281 intersect_lox = IntersectLoX;
282 intersect_hix = IntersectHiX;
283 intersect_loy = IntersectLoY;
284 intersect_hiy = IntersectHiY;
285 }
286
295 BoundaryPointList& IntersectHiX,
296 BoundaryPointList& IntersectLoY,
297 BoundaryPointList& IntersectHiY,
298 BoundaryPointList& IntersectLoZ,
299 BoundaryPointList& IntersectHiZ) noexcept
300 {
301 intersect_lox = IntersectLoX;
302 intersect_hix = IntersectHiX;
303 intersect_loy = IntersectLoY;
304 intersect_hiy = IntersectHiY;
305 intersect_loz = IntersectLoZ;
306 intersect_hiz = IntersectHiZ;
307 }
308
310 {
311 return intersect_lox;
312 }
314 {
315 return intersect_hix;
316 }
318 {
319 return intersect_loy;
320 }
322 {
323 return intersect_hiy;
324 }
326 {
327 return intersect_loz;
328 }
330 {
331 return intersect_hiz;
332 }
333
334#ifdef AMREX_PARTICLES
336 void RedistributeParticles ();
337#endif
338
339 void InstallNewDistributionMap (int lev, const DistributionMapping& newdm);
340
341 static bool UsingPrecreateDirectories () noexcept;
342
343protected:
344
346 void initialInit (Real strt_time, Real stop_time,
347 const BoxArray* lev0_grids = nullptr, const Vector<int>* pmap = nullptr);
348#ifndef AMREX_NO_PROBINIT
350 void readProbinFile (int& init);
351#endif
353 void checkInput ();
355 void restart (const std::string& filename);
357 void defBaseLevel (Real strt_time, const BoxArray* lev0_grids = nullptr, const Vector<int>* pmap = nullptr);
359 void bldFineLevels (Real strt_time);
361 virtual void regrid_level_0_on_restart ();
363 void grid_places (int lbase,
364 Real time,
365 int& new_finest,
366 Vector<BoxArray>& new_grids);
367
368 DistributionMapping makeLoadBalanceDistributionMap (int lev, Real time, const BoxArray& ba) const;
369 void LoadBalanceLevel0 (Real time);
370
371 void ErrorEst (int lev, TagBoxArray& tags, Real time, int ngrow) override;
372 BoxArray GetAreaNotToTag (int lev) override;
373 void ManualTagsPlacement (int lev, TagBoxArray& tags, const Vector<IntVect>& bf_lev) override;
374
376 virtual void timeStep (int level,
377 Real time,
378 int iteration,
379 int niter,
380 Real stop_time);
381
382 // pure virtual function in AmrCore
383 void MakeNewLevelFromScratch (int /*lev*/, Real /*time*/, const BoxArray& /*ba*/, const DistributionMapping& /*dm*/) override
384 { amrex::Abort("How did we get here!"); }
385 void MakeNewLevelFromCoarse (int /*lev*/, Real /*time*/, const BoxArray& /*ba*/, const DistributionMapping& /*dm*/) override
386 { amrex::Abort("How did we get here!"); }
387 void RemakeLevel (int /*lev*/, Real /*time*/, const BoxArray& /*ba*/, const DistributionMapping& /*dm*/) override
388 { amrex::Abort("How did we get here!"); }
389 void ClearLevel (int /*lev*/) override
390 { amrex::Abort("How did we get here!"); }
391
392 bool checkPointNow () noexcept;
394 bool writePlotNow () noexcept;
395 bool writeSmallPlotNow () noexcept;
396
397 void printGridInfo (std::ostream& os,
398 int min_lev,
399 int max_lev);
400
401 void setRecordGridInfo (const std::string&);
402
403 void setRecordRunInfo (const std::string&);
404
405 void setRecordRunInfoTerse (const std::string&);
406
407 void setRecordDataInfo (int i, const std::string&);
408
409 void initSubcycle();
410 void initPltAndChk();
411
412 static int initInSitu();
413 int updateInSitu();
414 static int finalizeInSitu();
415
416 //
417 // The data ...
418 //
419 std::string regrid_grids_file;
420 std::string initial_grids_file;
421 Vector<std::unique_ptr<AmrLevel> > amr_level;
422 Real cumtime = std::numeric_limits<Real>::lowest();
423 Real start_time = std::numeric_limits<Real>::lowest();
428 std::string subcycling_mode;
434 std::string check_file_root;
438 Real plot_per;
447 std::string plot_file_root;
449
451
455 std::ofstream gridlog;
456 std::ofstream runlog;
457 std::ofstream runlog_terse;
458 Vector<std::unique_ptr<std::fstream> > datalog;
459 Vector<std::string> datalogname;
461 std::string restart_chkfile;
462 std::string restart_pltfile;
463#ifndef AMREX_NO_PROBINIT
464 std::string probin_file;
465#endif
472
474
475 //
476 // The static data ...
477 //
478 static std::list<std::string> state_plot_vars;
479 static std::list<std::string> state_small_plot_vars;
480 static std::list<std::string> derive_plot_vars;
481 static std::list<std::string> derive_small_plot_vars;
482 static bool first_plotfile;
488
489#if defined(AMREX_USE_SENSEI_INSITU) && !defined(AMREX_NO_SENSEI_AMR_INST)
490 static AmrInSituBridge *insitu_bridge;
491#endif
492
493public:
500
502
503private:
504 void writePlotFileDoit (std::string const& pltfile, bool regular);
505};
506
507}
508
509#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
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:434
virtual void init(Real strt_time, Real stop_time)
Init data after construction. Must be called before timestepping.
Definition AMReX_Amr.cpp:1160
std::string probin_file
Definition AMReX_Amr.H:464
virtual void writeSmallPlotFile()
Write the small plot file to be used for visualization.
Definition AMReX_Amr.cpp:920
bool bUserStopRequest
Definition AMReX_Amr.H:473
BoundaryPointList & getIntersectHiX() noexcept
Definition AMReX_Amr.H:313
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:387
static void clearDeriveSmallPlotVarList()
Definition AMReX_Amr.cpp:726
static const std::list< std::string > & deriveSmallPlotVars() noexcept
Definition AMReX_Amr.H:164
static void fillDerivePlotVarList()
Fill the list of derive_plot_vars with all derived quantities.
Definition AMReX_Amr.cpp:690
Vector< int > n_cycle
Definition AMReX_Amr.H:427
void checkInput()
Check for valid input.
Definition AMReX_Amr.cpp:1082
void setRecordGridInfo(const std::string &)
Definition AMReX_Amr.cpp:771
DistributionMapping makeLoadBalanceDistributionMap(int lev, Real time, const BoxArray &ba) const
Definition AMReX_Amr.cpp:2801
void setRecordDataInfo(int i, const std::string &)
Definition AMReX_Amr.cpp:813
void setNCycle(const Vector< int > &ns) noexcept
Set the cycle count on each level.
Definition AMReX_Amr.cpp:841
virtual void writePlotFile()
Write the plot file to be used for visualization.
Definition AMReX_Amr.cpp:882
Vector< int > level_count
Definition AMReX_Amr.H:426
Vector< Real > dt_level
Timestep at this level.
Definition AMReX_Amr.H:424
static const Vector< BoxArray > & getInitialBA() noexcept
Definition AMReX_Amr.cpp:3429
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:445
static const BoxArray & initialBa(int level) noexcept
Array of BoxArrays read in to initially define grid hierarchy.
Definition AMReX_Amr.H:211
static void deleteDerivePlotVar(const std::string &name)
Remove the string from derive_plot_vars.
Definition AMReX_Amr.cpp:748
void writePlotFileDoit(std::string const &pltfile, bool regular)
Definition AMReX_Amr.cpp:959
int sub_cycle
Definition AMReX_Amr.H:460
const Vector< Real > & dtLevel() const noexcept
Array of time steps at all levels.
Definition AMReX_Amr.H:99
void initPltAndChk()
Definition AMReX_Amr.cpp:3268
int plot_max_level
Maximum AMR level to write to a plotfile.
Definition AMReX_Amr.H:440
int stream_max_tries
Definition AMReX_Amr.H:468
int NumDataLogs() noexcept
How many datalogs have been opened.
Definition AMReX_Amr.cpp:162
int stepOfLastSmallPlotFile() const noexcept
Definition AMReX_Amr.H:263
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:1961
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:385
std::string restart_pltfile
Definition AMReX_Amr.H:462
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:459
Real coarseTimeStepDt(Real stop_time)
Do a complete integration cycle and return the coarse dt.
Definition AMReX_Amr.cpp:2128
Vector< int > regrid_int
Interval between regridding.
Definition AMReX_Amr.H:430
void restart(const std::string &filename)
Restart from a checkpoint file.
Definition AMReX_Amr.cpp:1417
int plotMaxLevel() const noexcept
Maximum level to plot.
Definition AMReX_Amr.H:125
Vector< std::unique_ptr< std::fstream > > datalog
Definition AMReX_Amr.H:458
static void clearStatePlotVarList()
Clear the list of state_plot_vars.
Definition AMReX_Amr.cpp:626
int message_int
How often checking messages touched by user, such as "stop_run".
Definition AMReX_Amr.H:446
Vector< int > level_steps
Number of time steps at this level.
Definition AMReX_Amr.H:425
AmrLevel & getLevel(int lev) noexcept
AmrLevel lev.
Definition AMReX_Amr.H:189
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:487
BoundaryPointList & getIntersectLoY() noexcept
Definition AMReX_Amr.H:317
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:3376
Vector< Real > dt_min
Definition AMReX_Amr.H:429
static bool isStateSmallPlotVar(const std::string &name)
Definition AMReX_Amr.cpp:605
const std::string & theRestartPlotFile() const noexcept
Name of the restart plotfile.
Definition AMReX_Amr.H:238
Vector< std::unique_ptr< AmrLevel > > amr_level
Vector of levels.
Definition AMReX_Amr.H:421
static Vector< BoxArray > initial_ba
Array of BoxArrays read in to initially define grid hierarchy.
Definition AMReX_Amr.H:484
void printGridInfo(std::ostream &os, int min_lev, int max_lev)
Definition AMReX_Amr.cpp:2929
BoundaryPointList intersect_hiz
Definition AMReX_Amr.H:499
int loadbalance_level0_int
Definition AMReX_Amr.H:470
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:1309
static void clearDerivePlotVarList()
Clear the list of derive_plot_vars.
Definition AMReX_Amr.cpp:720
static std::list< std::string > derive_small_plot_vars
Derived Vars to dump to small plotfile.
Definition AMReX_Amr.H:481
BoundaryPointList & getIntersectHiY() noexcept
Definition AMReX_Amr.H:321
static bool isDerivePlotVar(const std::string &name) noexcept
Is the string the name of a variable in derive_plot_vars?
Definition AMReX_Amr.cpp:676
std::string small_plot_file_root
Root name of small plotfile.
Definition AMReX_Amr.H:448
static int finalizeInSitu()
Definition AMReX_Amr.cpp:582
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:732
void defBaseLevel(Real strt_time, const BoxArray *lev0_grids=nullptr, const Vector< int > *pmap=nullptr)
Define and initialize coarsest level.
Definition AMReX_Amr.cpp:2555
Real plotPer() const noexcept
Time between plot files.
Definition AMReX_Amr.H:121
virtual void regrid_level_0_on_restart()
Regrid level 0 on restart.
Definition AMReX_Amr.cpp:2867
static void clearStateSmallPlotVarList()
Definition AMReX_Amr.cpp:646
int loadbalance_with_workestimates
Definition AMReX_Amr.H:469
int stepOfLastPlotFile() const noexcept
Definition AMReX_Amr.H:260
int small_plot_int
How often small plotfile (# of time steps)
Definition AMReX_Amr.H:441
void RegridOnly(Real time, bool do_io=true)
Regrid only!
Definition AMReX_Amr.cpp:1926
BoundaryPointList & getIntersectLoX() noexcept
Definition AMReX_Amr.H:309
static std::list< std::string > state_small_plot_vars
State Vars to dump to small plotfile.
Definition AMReX_Amr.H:479
BoundaryPointList intersect_loy
Definition AMReX_Amr.H:495
Amr(Amr &&rhs)=delete
const std::string & theRestartFile() const noexcept
Name of the restart chkpoint file.
Definition AMReX_Amr.H:236
static std::list< std::string > state_plot_vars
State Vars to dump to plotfile.
Definition AMReX_Amr.H:478
std::ofstream runlog_terse
Definition AMReX_Amr.H:457
static bool isDeriveSmallPlotVar(const std::string &name) noexcept
Definition AMReX_Amr.cpp:683
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
std::string plot_file_root
Root name of plotfile.
Definition AMReX_Amr.H:447
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:450
std::string subcycling_mode
Type of subcycling to use.
Definition AMReX_Amr.H:428
BoundaryPointList intersect_lox
Definition AMReX_Amr.H:494
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:2482
void RedistributeParticles()
Redistribute particles.
Definition AMReX_Amr.cpp:3436
static std::list< std::string > derive_plot_vars
Derived Vars to dump to plotfile.
Definition AMReX_Amr.H:480
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:3076
int levelSteps(int lev) const noexcept
Number of time steps at specified level.
Definition AMReX_Amr.H:103
BoundaryPointList & getIntersectLoZ() noexcept
Definition AMReX_Amr.H:325
void setCumTime(Real t) noexcept
Definition AMReX_Amr.H:90
int stepOfLastCheckPoint() const noexcept
Definition AMReX_Amr.H:266
Amr(const Amr &rhs)=delete
BoundaryPointList & getIntersectHiZ() noexcept
Definition AMReX_Amr.H:329
int updateInSitu()
Definition AMReX_Amr.cpp:569
virtual void checkPoint()
Write current state into a chk* file.
Definition AMReX_Amr.cpp:1741
int last_plotfile
Step number of previous plotfile.
Definition AMReX_Amr.H:435
~Amr() override
The destructor.
Definition AMReX_Amr.cpp:763
void ClearLevel(int) override
Delete level data.
Definition AMReX_Amr.H:389
static bool isStatePlotVar(const std::string &name)
Is the string the name of a variable in state_plot_vars?
Definition AMReX_Amr.cpp:598
int plot_int
How often plotfile (# of time steps)
Definition AMReX_Amr.H:437
std::ofstream gridlog
Definition AMReX_Amr.H:455
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:2975
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:276
void InitAmr()
Definition AMReX_Amr.cpp:267
int okToContinue() noexcept
More work to be done?
Definition AMReX_Amr.cpp:869
void readProbinFile(int &init)
Read the probin file.
Definition AMReX_Amr.cpp:1199
static int initInSitu()
Definition AMReX_Amr.cpp:555
bool abort_on_stream_retry_failure
Definition AMReX_Amr.H:467
LevelBld * levelbld
Definition AMReX_Amr.H:466
static void deleteDeriveSmallPlotVar(const std::string &name)
Definition AMReX_Amr.cpp:756
Real cumtime
Physical time variable.
Definition AMReX_Amr.H:422
void setRecordRunInfo(const std::string &)
Definition AMReX_Amr.cpp:785
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:3088
void regrid(int lbase, Real time, bool initial=false) override
Rebuild grid hierarchy finer than lbase.
Definition AMReX_Amr.cpp:2616
Real check_per
How often checkpoint (units of time).
Definition AMReX_Amr.H:433
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:612
BoundaryPointList intersect_hix
Definition AMReX_Amr.H:497
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:442
std::string regrid_grids_file
Grids file that will bypass regridding.
Definition AMReX_Amr.H:419
std::multimap< std::pair< int, int >, double > BoundaryPointList
Definition AMReX_Amr.H:36
Amr & operator=(const Amr &rhs)=delete
int check_int
How often checkpoint (# time steps).
Definition AMReX_Amr.H:432
Real dtLevel(int level) const noexcept
Time step at specified level.
Definition AMReX_Amr.H:95
static const std::list< std::string > & stateSmallPlotVars() noexcept
Definition AMReX_Amr.H:138
static void fillDeriveSmallPlotVarList()
Definition AMReX_Amr.cpp:705
static bool first_smallplotfile
Definition AMReX_Amr.H:501
bool checkPointNow() noexcept
Definition AMReX_Amr.cpp:2361
std::string initial_grids_file
Grids file that will bypass regridding only at initialization.
Definition AMReX_Amr.H:420
const std::string & subcyclingMode() const noexcept
How are we subcycling?
Definition AMReX_Amr.H:81
Long cellCount() noexcept
Total number of cells.
Definition AMReX_Amr.cpp:849
void LoadBalanceLevel0(Real time)
Definition AMReX_Amr.cpp:2845
int record_run_info_terse
Definition AMReX_Amr.H:454
virtual void coarseTimeStep(Real stop_time)
Do a complete integration cycle.
Definition AMReX_Amr.cpp:2135
static void deleteStatePlotVar(const std::string &name)
Remove the string from state_plot_vars.
Definition AMReX_Amr.cpp:668
Real small_plot_log_per
How often small plotfile (in units of log10(time))
Definition AMReX_Amr.H:443
void bldFineLevels(Real strt_time)
Define and initialize refined levels.
Definition AMReX_Amr.cpp:3094
std::ofstream runlog
Definition AMReX_Amr.H:456
int record_grid_info
Definition AMReX_Amr.H:452
std::string restart_chkfile
Definition AMReX_Amr.H:461
Real plot_log_per
How often plotfile (in units of log10(time))
Definition AMReX_Amr.H:439
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:652
int last_smallplotfile
Step number of previous small plotfile.
Definition AMReX_Amr.H:436
void InstallNewDistributionMap(int lev, const DistributionMapping &newdm)
Definition AMReX_Amr.cpp:2854
bool writePlotNow() noexcept
Whether to write a plotfile now.
Definition AMReX_Amr.cpp:2409
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:3166
BoundaryPointList intersect_loz
Definition AMReX_Amr.H:496
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
std::string DataLogName(int i) const noexcept
The filename of the ith datalog file.
Definition AMReX_Amr.H:242
bool write_plotfile_with_checkpoint
Write out a plotfile whenever we checkpoint.
Definition AMReX_Amr.H:444
void FinalizeInit(Real strt_time, Real stop_time)
Second part of initialInit.
Definition AMReX_Amr.cpp:1341
static void fillStateSmallPlotVarList()
Definition AMReX_Amr.cpp:632
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:486
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:498
int numGrids() noexcept
Total number of grids.
Definition AMReX_Amr.cpp:859
int record_run_info
Definition AMReX_Amr.H:453
Real start_time
Physical time this simulation started.
Definition AMReX_Amr.H:423
Real loadbalance_max_fac
Definition AMReX_Amr.H:471
static void Finalize()
Definition AMReX_Amr.cpp:141
static void addStateSmallPlotVar(const std::string &name)
Definition AMReX_Amr.cpp:660
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:294
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:827
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:383
int last_checkpoint
Step number of previous checkpoint.
Definition AMReX_Amr.H:431
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
bool okToRegrid(int level) noexcept
Should we regrid this level?
Definition AMReX_Amr.cpp:3366
Real plot_per
How often plotfile (in units of time)
Definition AMReX_Amr.H:438
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:482
void setRecordRunInfoTerse(const std::string &)
Definition AMReX_Amr.cpp:799
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:1291
static void addDeriveSmallPlotVar(const std::string &name)
Definition AMReX_Amr.cpp:740
BoxArray GetAreaNotToTag(int lev) override
Definition AMReX_Amr.cpp:3082
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:551
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
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Definition AMReX_Amr.cpp:49
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230