Alamo
Integrator.H
Go to the documentation of this file.
1//
2// Pure abstract class for managing data structures, time integration (with substepping),
3// mesh refinement, and I/O.
4//
5// Native input file parameters:
6//
7// .. code-block:: make
8//
9// max_step = [maximum number of timesteps]
10// stop_time = [maximum simulation time]
11// timestep = [time step for coarsest level]
12//
13// amr.regrid_int = [number of timesteps between regridding]
14// amr.plot_int = [number of timesteps between dumping output]
15// amr.plot_file = [base name of output directory]
16//
17// amr.nsubsteps = [number of temporal substeps at each level. This can be
18// either a single int (which is then applied to every refinement
19// level) or an array of ints (equal to amr.max_level)
20// corresponding to the refinement for each level.]
21//
22// Inherited input file parameters (from amrex AmrMesh class):
23//
24// .. code-block:: make
25//
26// amr.v = [verbosity level]
27// amr.max_level = [maximum level of refinement]
28// amr.n_proper =
29// amr.grid_eff =
30// amr.n_error_buff =
31// amr.ref_ratio_vect = [refinement ratios in each direction]
32// amr.ref_ratio = [refinement ratio in all directions (cannot be used with ref_ratio_vect)]
33// amr.max_grid_x =
34// amr.blocking_factor =
35// amr.n_cell = [number of cells on coarsest level]
36// amr.refine_grid_layout =
37// amr.check_input =
38//
39
40#ifndef INTEGRATOR_INTEGRATOR_H
41#define INTEGRATOR_INTEGRATOR_H
42
43#include "IO/ParmParse.H"
44#include <chrono>
45#include <ctime>
46#include <string>
47#include <limits>
48#include <memory>
49
50#ifdef _OPENMP
51#include <omp.h>
52#endif
53
54#include <AMReX_ParallelDescriptor.H>
55#include <AMReX_ParmParse.H>
56#include <AMReX_MultiFabUtil.H>
57#include <AMReX_FillPatchUtil.H>
58#include <AMReX_BC_TYPES.H>
59#include <AMReX_AmrCore.H>
60#include <AMReX_FluxRegister.H>
61#include <AMReX_Utility.H>
62#include <AMReX_PlotFileUtil.H>
63
64#include "Set/Set.H"
65#include "BC/BC.H"
66#include "BC/Nothing.H"
67#include "IO/WriteMetaData.H"
68#include "BaseField.H"
69
70/// \brief Collection of numerical integrator objects
71namespace Integrator
72{
73
75 : public amrex::AmrCore
76{
77public:
78
79 /// This is the constructor for the intetgrator class, which reads timestep information,
80 /// simulation output and AMR, initialized time substep, and creates a new directory.
81 Integrator();
82
83 static void Parse(Integrator &, IO::ParmParse &);
84
85 /// Virtual destructure; make sure delete any pointers that you create here.
86 virtual ~Integrator();
87
88 /// Front-end method to initialize simulation on all levels
89 void InitData();
90
91 /// Read in output from previous simulation and start simulation at that point -
92 /// Not currently tested
93 void Restart(std::string restartfile, bool a_node = false);
94
95 /// Front-end method to start simulation
96 void Evolve();
97
98 /// Simple setter to set filename
99 void SetFilename(std::string _plot_file) { plot_file = _plot_file; };
100
101 /// Simple getter to get filename
102 std::string GetFilename() { return plot_file; };
103
104 /// This overrides an AMReX method just to allow for explicit meshing when
105 /// desired.
106 void regrid(int lbase, Set::Scalar time, bool initial = false) override
107 {
108 if (!explicitmesh.on)
109 AmrCore::regrid(lbase, time, initial);
110 }
111
112 /// This creates a new levels that have not previously been
113 /// used.
115 {
116 if (!explicitmesh.on) AmrCore::InitFromScratch(time);
117 else
118 {
119 // Generate the coarse level mesh.
120 {
121 finest_level = 0;
122 const amrex::BoxArray& ba = MakeBaseGrids();
123 amrex::DistributionMapping dm(ba);
124 const auto old_num_setdm = num_setdm;
125 const auto old_num_setba = num_setba;
126 MakeNewLevelFromScratch(0, time, ba, dm);
127 if (old_num_setba == num_setba) SetBoxArray(0, ba);
128 if (old_num_setdm == num_setdm) SetDistributionMap(0, dm);
129 }
130 // Generate subsequent level meshes based on user input
131 for (int ilev = 0; ilev < maxLevel(); ++ilev)
132 {
133 finest_level = ilev + 1;
134 amrex::BoxArray grids(explicitmesh.box[ilev]);
135 ChopGrids(ilev + 1, grids, amrex::ParallelDescriptor::NProcs());
136 amrex::DistributionMapping dmap(grids);
137 SetBoxArray(ilev + 1, grids);
138 SetDistributionMap(ilev + 1, dmap);
139 MakeNewLevelFromScratch(ilev + 1, time, grids, dmap);
140 }
141
142 }
143 }
144
145protected:
146
147 /// You must override this function to inherit this class;
148 /// this function is called before the simulation begins, and is where
149 /// initial conditions should be applied.
150 virtual void Initialize(int lev ///<[in] AMR Level
151 ) = 0;
152
153 /// You must override this function to inherit this class;
154 /// Advance is called every time(sub)step, and implements the evolution of
155 /// the system in time.
156 virtual void Advance(int lev, ///<[in] AMR Level
157 amrex::Real time, ///< [in] System time
158 amrex::Real dt ///< [in] Timestep for this level
159 ) = 0;
160
161 /// You must override this function to inherit this class;
162 /// Advance is called every time(sub)step, and implements the evolution of
163 /// the system in time.
164 virtual void TagCellsForRefinement(int lev, amrex::TagBoxArray& tags, amrex::Real time,
165 int ngrow) = 0;
166
167 /// This optional function is called at the beginning of every timestep, and can be used
168 /// to complete additional global solves, e.g. a MLMG implicit solve.
169 virtual void TimeStepBegin(Set::Scalar /*time*/, int /*iter*/) {};
170
171 /// This optional function is called at the end of every timestep, and can be used
172 /// to complete additional global solves, e.g. a MLMG implicit solve.
173 virtual void TimeStepComplete(Set::Scalar /*time*/, int /*iter*/) {};
174
175 /// This is a function that is called by `Integrator` to update the variables registered in
176 /// RegisterIntegratedVariable; you can override this to do your integration.
177 virtual void Integrate( int /*amrlev*/,
178 Set::Scalar /*time*/,
179 int /*iter*/,
180 const amrex::MFIter&/*mfi*/,
181 const amrex::Box&/*box*/
182 )
183 {
184 if (thermo.number > 0)
185 Util::Warning(INFO, "integrated variables registered, but no integration implemented!");
186 }
187
188 /// An optionally overridable method to trigger behavior whenver a regrid occurs.
189 virtual void Regrid(int /* amrlev */, Set::Scalar /* time */)
190 {}
191
192
193 inline void StopClock()
194 {
195 clock_running = false;
196 }
197 inline void RestartClock()
198 {
199 clock_running = true;
200 }
201
202
203 /// Add a new cell-based scalar field.
204 void RegisterNewFab(Set::Field<Set::Scalar>& new_fab, BC::BC<Set::Scalar>* new_bc, int ncomp, int nghost, std::string name, bool writeout,bool evolving = true,std::vector<std::string> suffix = {});
205 /// Add a new cell-based scalar field (with additional arguments).
206 void RegisterNewFab(Set::Field<Set::Scalar>& new_fab, int ncomp, std::string name, bool writeout,bool evolving=true, std::vector<std::string> suffix = {});
207 /// Add a new node-based scalar field
208 void RegisterNodalFab(Set::Field<Set::Scalar>& new_fab, int ncomp, int nghost, std::string name, bool writeout,bool evolving=true,std::vector<std::string> suffix = {});
209 /// Add a new node-based scalar field (wtih additional arguments)
210 void RegisterNodalFab(Set::Field<Set::Scalar>& new_fab, BC::BC<Set::Scalar>* new_bc, int ncomp, int nghost, std::string name, bool writeout,bool evolving=true, std::vector<std::string> suffix = {});
211 template<class T>
212 /// Add a templated nodal field
213 void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool evolving = true);
214 template<class T>
215 /// Add a templated nodal field (additional arguments)
216 void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, std::string a_name, bool evolving = true);
217 template<class T>
218 /// Add a templated nodal field (additional arguments)
219 void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool writeout, std::string a_name, bool evolving = true);
220
221 /// Add a field with arbitrary type (templated with T) and grid location (templated with d).
222 template<class T, int d>
223 void AddField( Set::Field<T>& new_field, BC::BC<T>* new_bc, int ncomp, int nghost,
224 std::string, bool writeout, bool evolving, std::vector<std::string> suffix = {});
225
226 /// Utility to ensure that all fields know what the finest level is
227 void SetFinestLevel(const int a_finestlevel)
228 {
229 for (unsigned int i = 0; i < cell.fab_array.size(); i++)
230 cell.fab_array[i]->finest_level = a_finestlevel;
231 for (unsigned int i = 0; i < node.fab_array.size(); i++)
232 node.fab_array[i]->finest_level = a_finestlevel;
233 for (unsigned int i = 0; i < m_basefields_cell.size(); i++)
234 m_basefields_cell[i]->SetFinestLevel(finest_level);
235 for (unsigned int i = 0; i < m_basefields.size(); i++)
236 m_basefields[i]->SetFinestLevel(finest_level);
237 }
238
239 /// Register a variable to be integrated over the spatial domain using
240 /// the Integrate function.
241 void RegisterIntegratedVariable(Set::Scalar* integrated_variable, std::string name, bool extensive=true);
242
243 /// Utility to set the coarse-grid timestep
244 void SetTimestep(Set::Scalar _timestep);
245
246 /// Utility to set the frequency (in timesteps) of plotfile dumping
247 void SetPlotInt(int plot_int);
248
249 /// Utility to set the frequency (in timesteps) of thermo data calculation
250 void SetThermoInt(int a_thermo_int) { thermo.interval = a_thermo_int; }
251 /// Utility to set the frequency (in timesteps) of thermo data writing to file
252 void SetThermoPlotInt(int a_thermo_plot_int) { thermo.plot_int = a_thermo_plot_int; }
253 /// Utility to set the global stop time.
254 void SetStopTime(Set::Scalar a_stop_time) { stop_time = a_stop_time; }
255
256
257 // Dynamic timestep adjustment
258
259 struct {
260 // user params
261 bool on = false;
262 int verbose = -1;
263 int nprevious = -1;
264 Set::Scalar cfl = NAN;
265 Set::Scalar min = NAN;
266 Set::Scalar max = NAN;
267 // internal variables
268 std::vector<Set::Scalar> dt_limit_min;
269 std::vector<Set::Scalar> previous_timesteps;
270 } dynamictimestep; /// Params for the dynamic timestp
272 {
273 if (!dynamictimestep.on) return;
274 if (!dynamictimestep.dt_limit_min.size())
275 {
276 dynamictimestep.dt_limit_min.resize(max_level+1,
277 std::numeric_limits<Set::Scalar>::max());
278 }
279 dynamictimestep.dt_limit_min[lev] = std::min(dt_min,
280 dynamictimestep.dt_limit_min[lev]);
281
282 amrex::ParallelDescriptor::ReduceRealMin(dynamictimestep.dt_limit_min[lev]);
283 }
285 {
286 if (!dynamictimestep.on) return;
287 dynamictimestep.previous_timesteps.clear();
288 }
290 {
291 if (!dynamictimestep.on) return;
292 if (!dynamictimestep.dt_limit_min.size()) return;
293
294 Set::Scalar final_timestep = NAN;
295
296 if (amrex::ParallelDescriptor::IOProcessor())
297 {
298 Set::Scalar timestep_average = this->dt[0];
299 if (dynamictimestep.previous_timesteps.size() > 0)
300 {
301 timestep_average = 0.0;
302 for (unsigned int d = 0; d < dynamictimestep.previous_timesteps.size(); d++)
303 timestep_average += dynamictimestep.previous_timesteps[d];
304 timestep_average /= dynamictimestep.previous_timesteps.size();
305 }
306
307 Set::Scalar new_timestep = std::numeric_limits<Set::Scalar>::max();
308 for (int lev = 0; lev <= this->max_level; lev++)
309 {
310 //const Set::Scalar* DX = this->geom[lev].CellSize();
311 Set::Scalar dt_lev = dynamictimestep.dt_limit_min[lev];
312
313 Util::Message(INFO,"lev=",lev," ",dt_lev, " (",this->nsubsteps[lev],")");
314
315 for (int ilev = lev; ilev > 0; ilev--) dt_lev *= (Set::Scalar)(this->nsubsteps[ilev]);
316
317 Util::Message(INFO,"lev=",lev," --> ",dt_lev);
318
319 new_timestep = std::min(new_timestep,dt_lev);
320 }
321
322 if (new_timestep < timestep_average)
323 {
324 dynamictimestep.previous_timesteps.clear();
325
326 final_timestep = new_timestep;
327 final_timestep = std::max(final_timestep,dynamictimestep.min);
328 final_timestep = std::min(final_timestep,dynamictimestep.max);
329
330 dynamictimestep.previous_timesteps.push_back(new_timestep);
331 }
332 else
333 {
334 final_timestep = timestep_average;
335 final_timestep = std::max(final_timestep,dynamictimestep.min);
336 final_timestep = std::min(final_timestep,dynamictimestep.max);
337
338 if ((int)(dynamictimestep.previous_timesteps.size()) > dynamictimestep.nprevious)
339 dynamictimestep.previous_timesteps.erase(dynamictimestep.previous_timesteps.begin()); // pop first
340 dynamictimestep.previous_timesteps.push_back(new_timestep); // push back new timestep
341 }
342
343 }
344 amrex::ParallelDescriptor::Bcast(&final_timestep,1);
345 this->SetTimestep(final_timestep);
346 dynamictimestep.dt_limit_min.clear();
347 }
348
349
350
351 amrex::Vector<amrex::Real> t_new; ///< Keep track of current old simulation time on each level
352 amrex::Vector<int> istep; ///< Keep track of where each level is
353 // PLOT FILES
354 std::string plot_file{ "plt" }; ///< Plotfile name
355
356private:
357 virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray& ba,
358 const amrex::DistributionMapping& dm) override;
359 virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray& ba,
360 const amrex::DistributionMapping& dm) override;
361 virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray& ba,
362 const amrex::DistributionMapping& dm) override;
363 virtual void ClearLevel(int lev) override;
364 virtual void ErrorEst(int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
365
366
367 /// This is the function that is responsible for updating patch data.
368 void FillPatch(int lev, amrex::Real time,
369 amrex::Vector<std::unique_ptr<amrex::MultiFab>>& source_mf,
370 amrex::MultiFab& destination_multifab,
371 BC::BC<Set::Scalar>& physbc,
372 int icomp);
373 /// Simple utility to count cells
374 long CountCells(int lev);
375 /// Timestep marching
376 void TimeStep(int lev, amrex::Real time, int iteration);
377 void FillCoarsePatch(int lev, amrex::Real time, Set::Field<Set::Scalar>& mf, BC::BC<Set::Scalar>& physbc, int icomp, int ncomp);
378 void GetData(const int lev, const amrex::Real time, amrex::Vector<amrex::MultiFab*>& data, amrex::Vector<amrex::Real>& datatime);
379
380 std::vector<std::string> PlotFileName(int lev, std::string prefix = "") const;
381protected:
382 void IntegrateVariables(Set::Scalar cur_time, int step);
383 void WritePlotFile(bool initial = false) const;
384 void WritePlotFile(std::string prefix, Set::Scalar time, int step) const;
385 void WritePlotFile(Set::Scalar time, amrex::Vector<int> iter, bool initial = false, std::string prefix = "") const;
386
387 //
388 // MEMBER VARIABLES
389 //
390
391 // TIME (STEP) KEEPINGamrex::Vector<std::unique_ptr<amrex::MultiFab> >
392protected:
393 amrex::Real timestep = NAN; ///< Timestep for the base level of refinement
394 amrex::Vector<amrex::Real> dt; ///< Timesteps for each level of refinement
395 amrex::Vector<int> nsubsteps; ///< how many substeps on each level?
396private:
397 int max_plot_level = -1;
398 int print_ghost_nodes = 0;
399 int print_ghost_cells = 0;
400
401 amrex::Vector<amrex::Real> t_old;///< Keep track of current old simulation time on each level
402 int max_step = std::numeric_limits<int>::max(); ///< Maximum allowable timestep
403 amrex::Real tstart = 0; ///< Default start time (default: 0)
404 amrex::Real stop_time = NAN; ///< Default stop time
405
406protected:
407 bool integrate_variables_before_advance = true;
408 bool integrate_variables_after_advance = false;
409
410 bool clock_running = true;
411
412protected:
413 struct {
414 int number_of_fabs = 0;
415 std::vector<Set::Field<Set::Scalar>*> fab_array;
416 std::vector<int> ncomp_array;
417 std::vector<int> nghost_array;
418 std::vector<std::vector<std::string>> name_array;
419 std::vector<BC::BC<Set::Scalar>*> physbc_array;
420 std::vector<bool> writeout_array;
421 std::vector<bool> evolving_array;
422 bool any = true;
423 bool all = false;
424 } node;
425
426 struct {
427 int number_of_fabs = 0;
428 std::vector<Set::Field<Set::Scalar>*> fab_array;
429 std::vector<int> ncomp_array;
430 std::vector<int> nghost_array;
431 std::vector<std::vector<std::string>> name_array;
432 std::vector<BC::BC<Set::Scalar>*> physbc_array;
433 std::vector<bool> writeout_array;
434 std::vector<bool> evolving_array;
435 bool any = true;
436 bool all = false;
437 } cell;
438
439 std::vector<BaseField*> m_basefields;
440 std::vector<BaseField*> m_basefields_cell;
441
443
444 // KEEP TRACK OF ALL INTEGRATED VARIABLES
445 struct {
446 int interval = -1;
447 Set::Scalar dt = NAN;
448 int plot_int = -1;
449 Set::Scalar plot_dt = NAN;
450 int number = 0;
451 std::vector<Set::Scalar*> vars;
452 std::vector<std::string> names;
453 std::vector<bool> extensives;
454 } thermo;
455
456 // REGRIDDING
457 int regrid_int = -1; ///< Determine how often to regrid (default: 2)
458 int base_regrid_int = -1; ///< Determine how often to regrid based on coarse level only (default: 0)
459
460 std::string restart_file_cell = "";
461 std::string restart_file_node = "";
462
463 struct {
464 int on = 0;
465 std::vector<amrex::Box> box;
466 } explicitmesh;
467
468protected:
469 int plot_int = -1; ///< How frequently to dump plot file (default: never)
470 Set::Scalar plot_dt = -1.0;
471
472 int abort_on_nan = true;
473};
474
475
476template<>
479( Set::Field<Set::Scalar>& new_field,
480 BC::BC<Set::Scalar>* new_bc,
481 int ncomp,
482 int nghost,
483 std::string name,
484 bool writeout,
485 bool evolving,
486 std::vector<std::string> suffix)
487{
488 int nlevs_max = maxLevel() + 1;
489 new_field.resize(nlevs_max);
490 cell.fab_array.push_back(&new_field);
491 if (new_bc != nullptr) cell.physbc_array.push_back(new_bc);
492 else cell.physbc_array.push_back(&bcnothing);
493 cell.ncomp_array.push_back(ncomp);
494 cell.nghost_array.push_back(nghost);
495 //cell.name_array.push_back(name);
496 cell.writeout_array.push_back(writeout);
497 cell.evolving_array.push_back(evolving);
498 cell.number_of_fabs++;
499
500 Util::Assert(INFO,TEST((int)suffix.size() == 0 || (int)suffix.size() == ncomp));
501 std::vector<std::string> names;
502 if (ncomp == 1)
503 names.push_back(name);
504 else
505 {
506 if (suffix.size() == 0)
507 for (int j = 0; j < ncomp; j++)
508 names.push_back(amrex::Concatenate(name,j+1,3));
509 else
510 for (int j = 0; j < ncomp; j++)
511 names.push_back(name + suffix[j]);
512 }
513 cell.name_array.push_back(names);
514}
515
516template<>
519( Set::Field<Set::Scalar>& new_field,
520 BC::BC<Set::Scalar>* new_bc,
521 int ncomp,
522 int nghost,
523 std::string name,
524 bool writeout,
525 bool evolving,
526 std::vector<std::string> suffix)
527{
528 BL_PROFILE("Integrator::RegisterNodalFab");
529 Util::Assert(INFO, TEST(new_bc == nullptr));
530 int nlevs_max = maxLevel() + 1;
531 new_field.resize(nlevs_max);
532 node.fab_array.push_back(&new_field);
533 node.physbc_array.push_back(&bcnothing);
534 node.ncomp_array.push_back(ncomp);
535 node.nghost_array.push_back(nghost);
536 node.evolving_array.push_back(evolving);
537 node.writeout_array.push_back(writeout);
538 node.number_of_fabs++;
539
540 Util::Assert(INFO,TEST((int)suffix.size() == 0 || (int)suffix.size() == ncomp));
541 std::vector<std::string> names;
542 if (ncomp == 1)
543 names.push_back(name);
544 else
545 {
546 if (suffix.size() == 0)
547 for (int j = 0; j < ncomp; j++)
548 names.push_back(amrex::Concatenate(name,j+1,3));
549 else
550 for (int j = 0; j < ncomp; j++)
551 names.push_back(name + suffix[j]);
552 }
553 node.name_array.push_back(names);
554}
555
556template<class T, int d>
559( Set::Field<T>& new_field, BC::BC<T>* new_bc, int ncomp,
560 int nghost, std::string name, bool writeout, bool evolving,
561 std::vector<std::string> /*suffix*/)
562{
563 if (d == Set::Hypercube::Node)
564 {
565 Util::Assert(INFO, TEST(new_bc == nullptr));
566 int nlevs_max = maxLevel() + 1;
567 new_field.resize(nlevs_max);
568 m_basefields.push_back(new Field<T>(new_field, geom, refRatio(), ncomp, nghost));
569 m_basefields.back()->evolving = evolving;
570 m_basefields.back()->writeout = writeout;
571 m_basefields.back()->setName(name);
572 m_basefields.back()->evolving = evolving;
573 m_basefields.back()->m_gridtype = Set::Hypercube::Node;
574 }
575 else if (d == Set::Hypercube::Cell)
576 {
577 int nlevs_max = maxLevel() + 1;
578 new_field.resize(nlevs_max);
579 m_basefields_cell.push_back(new Field<T>(new_field, geom, refRatio(), ncomp, nghost));
580 m_basefields_cell.back()->evolving = evolving;
581 m_basefields_cell.back()->writeout = writeout;
582 m_basefields_cell.back()->setName(name);
583 m_basefields_cell.back()->evolving = evolving;
584 if (new_bc) m_basefields_cell.back()->setBC(new_bc);
585 else m_basefields_cell.back()->setBC(&bcnothing);
586 m_basefields_cell.back()->m_gridtype = Set::Hypercube::Cell;
587 }
588 else
589 {
590 Util::Abort(INFO, "Only node and cell based fields can be added at this time");
591 }
592}
593
594
595
596template<class T>
598void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool evolving)
599{
600 //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
601 AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, "", true, evolving);
602}
603template<class T>
605void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, std::string a_name, bool evolving)
606{
607 //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
608 AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, a_name, true, evolving);
609}
610template<class T>
611AMREX_ATTRIBUTE_WEAK
612void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool writeout, std::string a_name, bool evolving)
613{
614 //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
615 AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, a_name, writeout, evolving);
616}
617
618}
619#endif
#define TEST(x)
Definition Util.H:24
#define ALAMO_SINGLE_DEFINITION
Definition Util.H:28
#define INFO
Definition Util.H:23
Definition BC.H:43
std::vector< bool > evolving_array
Definition Integrator.H:421
std::vector< BaseField * > m_basefields
Definition Integrator.H:439
virtual void Initialize(int lev)=0
You must override this function to inherit this class; this function is called before the simulation ...
void SetFinestLevel(const int a_finestlevel)
Utility to ensure that all fields know what the finest level is.
Definition Integrator.H:227
void SetStopTime(Set::Scalar a_stop_time)
Utility to set the global stop time.
Definition Integrator.H:254
std::string GetFilename()
Simple getter to get filename.
Definition Integrator.H:102
void DynamicTimestep_SyncTimeStep(int lev, Set::Scalar dt_min)
Params for the dynamic timestp.
Definition Integrator.H:271
void RegisterGeneralFab(Set::Field< T > &new_fab, int ncomp, int nghost, bool writeout, std::string a_name, bool evolving=true)
Add a templated nodal field (additional arguments)
std::vector< int > nghost_array
Definition Integrator.H:417
void regrid(int lbase, Set::Scalar time, bool initial=false) override
This overrides an AMReX method just to allow for explicit meshing when desired.
Definition Integrator.H:106
void SetThermoInt(int a_thermo_int)
Utility to set the frequency (in timesteps) of thermo data calculation.
Definition Integrator.H:250
std::vector< bool > extensives
Definition Integrator.H:453
virtual void TimeStepBegin(Set::Scalar, int)
This optional function is called at the beginning of every timestep, and can be used to complete addi...
Definition Integrator.H:169
std::vector< Set::Scalar * > vars
Definition Integrator.H:451
void RegisterGeneralFab(Set::Field< T > &new_fab, int ncomp, int nghost, std::string a_name, bool evolving=true)
Add a templated nodal field (additional arguments)
std::vector< amrex::Box > box
Definition Integrator.H:465
std::vector< BaseField * > m_basefields_cell
Definition Integrator.H:440
std::vector< BC::BC< Set::Scalar > * > physbc_array
Definition Integrator.H:419
virtual void Advance(int lev, amrex::Real time, amrex::Real dt)=0
You must override this function to inherit this class; Advance is called every time(sub)step,...
std::vector< int > ncomp_array
Definition Integrator.H:416
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition Integrator.H:394
virtual void TagCellsForRefinement(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow)=0
You must override this function to inherit this class; Advance is called every time(sub)step,...
amrex::Vector< amrex::Real > t_old
Keep track of current old simulation time on each level.
Definition Integrator.H:401
virtual void Integrate(int, Set::Scalar, int, const amrex::MFIter &, const amrex::Box &)
This is a function that is called by Integrator to update the variables registered in RegisterIntegra...
Definition Integrator.H:177
std::vector< Set::Scalar > previous_timesteps
Definition Integrator.H:269
void AddField(Set::Field< T > &new_field, BC::BC< T > *new_bc, int ncomp, int nghost, std::string, bool writeout, bool evolving, std::vector< std::string > suffix={})
Add a field with arbitrary type (templated with T) and grid location (templated with d).
amrex::Vector< int > nsubsteps
how many substeps on each level?
Definition Integrator.H:395
std::vector< std::vector< std::string > > name_array
Definition Integrator.H:418
std::vector< std::string > names
Definition Integrator.H:452
void SetThermoPlotInt(int a_thermo_plot_int)
Utility to set the frequency (in timesteps) of thermo data writing to file.
Definition Integrator.H:252
amrex::Vector< amrex::Real > t_new
Keep track of current old simulation time on each level.
Definition Integrator.H:351
amrex::Vector< int > istep
Keep track of where each level is.
Definition Integrator.H:352
std::vector< Set::Scalar > dt_limit_min
Definition Integrator.H:268
void GetData(const int lev, const amrex::Real time, amrex::Vector< amrex::MultiFab * > &data, amrex::Vector< amrex::Real > &datatime)
void InitFromScratch(Set::Scalar time)
This creates a new levels that have not previously been used.
Definition Integrator.H:114
void RegisterGeneralFab(Set::Field< T > &new_fab, int ncomp, int nghost, bool evolving=true)
Add a templated nodal field.
virtual void Regrid(int, Set::Scalar)
An optionally overridable method to trigger behavior whenver a regrid occurs.
Definition Integrator.H:189
std::vector< bool > writeout_array
Definition Integrator.H:420
std::vector< Set::Field< Set::Scalar > * > fab_array
Definition Integrator.H:415
void SetFilename(std::string _plot_file)
Simple setter to set filename.
Definition Integrator.H:99
virtual void TimeStepComplete(Set::Scalar, int)
This optional function is called at the end of every timestep, and can be used to complete additional...
Definition Integrator.H:173
Collection of numerical integrator objects.
Definition AllenCahn.H:42
ALAMO_SINGLE_DEFINITION void Integrator::AddField< Set::Scalar, Set::Hypercube::Node >(Set::Field< Set::Scalar > &new_field, BC::BC< Set::Scalar > *new_bc, int ncomp, int nghost, std::string name, bool writeout, bool evolving, std::vector< std::string > suffix)
Definition Integrator.H:519
ALAMO_SINGLE_DEFINITION void Integrator::AddField< Set::Scalar, Set::Hypercube::Cell >(Set::Field< Set::Scalar > &new_field, BC::BC< Set::Scalar > *new_bc, int ncomp, int nghost, std::string name, bool writeout, bool evolving, std::vector< std::string > suffix)
Definition Integrator.H:479
amrex::Real Scalar
Definition Base.H:18
@ Cell
Definition Set.H:32
@ Node
Definition Set.H:32
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
void Abort(const char *msg)
Definition Util.cpp:263
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:57
void Warning(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:201
void Message(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:128