40#ifndef INTEGRATOR_INTEGRATOR_H
41#define INTEGRATOR_INTEGRATOR_H
53#include <AMReX_ParallelDescriptor.H>
54#include <AMReX_ParmParse.H>
55#include <AMReX_MultiFabUtil.H>
56#include <AMReX_FillPatchUtil.H>
57#include <AMReX_BC_TYPES.H>
58#include <AMReX_AmrCore.H>
59#include <AMReX_FluxRegister.H>
60#include <AMReX_Utility.H>
61#include <AMReX_PlotFileUtil.H>
74 :
public amrex::AmrCore
90 void Restart(std::string restartfile,
bool a_node =
false);
96 void SetFilename(std::string _plot_file) { plot_file = _plot_file; };
105 if (!explicitmesh.on)
106 AmrCore::regrid(lbase, time, initial);
113 if (!explicitmesh.on) AmrCore::InitFromScratch(time);
119 const amrex::BoxArray& ba = MakeBaseGrids();
120 amrex::DistributionMapping dm(ba);
121 const auto old_num_setdm = num_setdm;
122 const auto old_num_setba = num_setba;
123 MakeNewLevelFromScratch(0, time, ba, dm);
124 if (old_num_setba == num_setba) SetBoxArray(0, ba);
125 if (old_num_setdm == num_setdm) SetDistributionMap(0, dm);
128 for (
int ilev = 0; ilev < maxLevel(); ++ilev)
130 finest_level = ilev + 1;
131 amrex::BoxArray grids(explicitmesh.box[ilev]);
132 ChopGrids(ilev + 1, grids, amrex::ParallelDescriptor::NProcs());
133 amrex::DistributionMapping dmap(grids);
134 SetBoxArray(ilev + 1, grids);
135 SetDistributionMap(ilev + 1, dmap);
136 MakeNewLevelFromScratch(ilev + 1, time, grids, dmap);
177 const amrex::MFIter&,
181 if (thermo.number > 0)
182 Util::Warning(
INFO,
"integrated variables registered, but no integration implemented!");
193 void RegisterNewFab(
Set::Field<Set::Scalar>& new_fab,
int ncomp, std::string name,
bool writeout,std::vector<std::string> suffix = {});
195 void RegisterNodalFab(
Set::Field<Set::Scalar>& new_fab,
int ncomp,
int nghost, std::string name,
bool writeout,std::vector<std::string> suffix = {});
209 template<
class T,
int d>
211 std::string,
bool writeout,
bool evolving, std::vector<std::string> suffix = {});
216 for (
unsigned int i = 0; i < cell.fab_array.size(); i++)
217 cell.fab_array[i]->finest_level = a_finestlevel;
218 for (
unsigned int i = 0; i < node.fab_array.size(); i++)
219 node.fab_array[i]->finest_level = a_finestlevel;
220 for (
unsigned int i = 0; i < m_basefields_cell.size(); i++)
221 m_basefields_cell[i]->SetFinestLevel(finest_level);
222 for (
unsigned int i = 0; i < m_basefields.size(); i++)
223 m_basefields[i]->SetFinestLevel(finest_level);
228 void RegisterIntegratedVariable(
Set::Scalar* integrated_variable, std::string name,
bool extensive=
true);
234 void SetPlotInt(
int plot_int);
237 void SetThermoInt(
int a_thermo_int) { thermo.interval = a_thermo_int; }
260 if (!dynamictimestep.on)
return;
261 if (!dynamictimestep.dt_limit_min.size())
263 dynamictimestep.dt_limit_min.resize(max_level+1,
264 std::numeric_limits<Set::Scalar>::max());
266 dynamictimestep.dt_limit_min[lev] = std::min(dt_min,
267 dynamictimestep.dt_limit_min[lev]);
269 amrex::ParallelDescriptor::ReduceRealMin(dynamictimestep.dt_limit_min[lev]);
273 if (!dynamictimestep.on)
return;
274 dynamictimestep.previous_timesteps.clear();
278 if (!dynamictimestep.on)
return;
279 if (!dynamictimestep.dt_limit_min.size())
return;
283 if (amrex::ParallelDescriptor::IOProcessor())
286 if (dynamictimestep.previous_timesteps.size() > 0)
288 timestep_average = 0.0;
289 for (
unsigned int d = 0; d < dynamictimestep.previous_timesteps.size(); d++)
290 timestep_average += dynamictimestep.previous_timesteps[d];
291 timestep_average /= dynamictimestep.previous_timesteps.size();
294 Set::Scalar new_timestep = std::numeric_limits<Set::Scalar>::max();
295 for (
int lev = 0; lev <= this->max_level; lev++)
298 Set::Scalar dt_lev = dynamictimestep.dt_limit_min[lev];
302 for (
int ilev = lev; ilev > 0; ilev--) dt_lev *= (
Set::Scalar)(this->nsubsteps[ilev]);
306 new_timestep = std::min(new_timestep,dt_lev);
309 if (new_timestep < timestep_average)
311 dynamictimestep.previous_timesteps.clear();
313 final_timestep = new_timestep;
314 final_timestep = std::max(final_timestep,dynamictimestep.min);
315 final_timestep = std::min(final_timestep,dynamictimestep.max);
317 dynamictimestep.previous_timesteps.push_back(new_timestep);
321 final_timestep = timestep_average;
322 final_timestep = std::max(final_timestep,dynamictimestep.min);
323 final_timestep = std::min(final_timestep,dynamictimestep.max);
325 if ((
int)(dynamictimestep.previous_timesteps.size()) > dynamictimestep.nprevious)
326 dynamictimestep.previous_timesteps.erase(dynamictimestep.previous_timesteps.begin());
327 dynamictimestep.previous_timesteps.push_back(new_timestep);
331 amrex::ParallelDescriptor::Bcast(&final_timestep,1);
332 this->SetTimestep(final_timestep);
333 dynamictimestep.dt_limit_min.clear();
341 std::string plot_file{
"plt" };
344 virtual void MakeNewLevelFromScratch(
int lev, amrex::Real time,
const amrex::BoxArray& ba,
345 const amrex::DistributionMapping& dm)
override;
346 virtual void MakeNewLevelFromCoarse(
int lev, amrex::Real time,
const amrex::BoxArray& ba,
347 const amrex::DistributionMapping& dm)
override;
348 virtual void RemakeLevel(
int lev, amrex::Real time,
const amrex::BoxArray& ba,
349 const amrex::DistributionMapping& dm)
override;
350 virtual void ClearLevel(
int lev)
override;
351 virtual void ErrorEst(
int lev, amrex::TagBoxArray& tags, amrex::Real time,
int ngrow)
override;
355 void FillPatch(
int lev, amrex::Real time,
356 amrex::Vector<std::unique_ptr<amrex::MultiFab>>& source_mf,
357 amrex::MultiFab& destination_multifab,
361 long CountCells(
int lev);
363 void TimeStep(
int lev, amrex::Real time,
int iteration);
365 void GetData(
const int lev,
const amrex::Real time, amrex::Vector<amrex::MultiFab*>& data, amrex::Vector<amrex::Real>& datatime);
367 std::vector<std::string> PlotFileName(
int lev, std::string prefix =
"")
const;
369 void IntegrateVariables(
Set::Scalar cur_time,
int step);
370 void WritePlotFile(
bool initial =
false)
const;
371 void WritePlotFile(std::string prefix,
Set::Scalar time,
int step)
const;
372 void WritePlotFile(
Set::Scalar time, amrex::Vector<int> iter,
bool initial =
false, std::string prefix =
"")
const;
380 amrex::Real timestep = NAN;
381 amrex::Vector<amrex::Real>
dt;
384 int max_plot_level = -1;
387 int max_step = std::numeric_limits<int>::max();
388 amrex::Real tstart = 0;
389 amrex::Real stop_time = NAN;
392 bool integrate_variables_before_advance =
true;
393 bool integrate_variables_after_advance =
false;
397 int number_of_fabs = 0;
409 int number_of_fabs = 0;
410 std::vector<Set::Field<Set::Scalar>*> fab_array;
411 std::vector<int> ncomp_array;
412 std::vector<int> nghost_array;
413 std::vector<std::vector<std::string>> name_array;
414 std::vector<BC::BC<Set::Scalar>*> physbc_array;
415 std::vector<bool> writeout_array;
432 std::vector<Set::Scalar*>
vars;
439 int base_regrid_int = -1;
441 std::string restart_file_cell =
"";
442 std::string restart_file_node =
"";
446 std::vector<amrex::Box>
box;
453 int abort_on_nan =
true;
467 std::vector<std::string> suffix)
469 int nlevs_max = maxLevel() + 1;
470 new_field.resize(nlevs_max);
471 cell.fab_array.push_back(&new_field);
472 if (new_bc !=
nullptr) cell.physbc_array.push_back(new_bc);
473 else cell.physbc_array.push_back(&bcnothing);
474 cell.ncomp_array.push_back(ncomp);
475 cell.nghost_array.push_back(nghost);
477 cell.writeout_array.push_back(writeout);
478 cell.number_of_fabs++;
481 std::vector<std::string> names;
483 names.push_back(name);
486 if (suffix.size() == 0)
487 for (
int j = 0; j < ncomp; j++)
488 names.push_back(amrex::Concatenate(name,j+1,3));
490 for (
int j = 0; j < ncomp; j++)
491 names.push_back(name + suffix[j]);
493 cell.name_array.push_back(names);
506 std::vector<std::string> suffix)
508 BL_PROFILE(
"Integrator::RegisterNodalFab");
510 int nlevs_max = maxLevel() + 1;
511 new_field.resize(nlevs_max);
512 node.fab_array.push_back(&new_field);
513 node.physbc_array.push_back(&bcnothing);
514 node.ncomp_array.push_back(ncomp);
515 node.nghost_array.push_back(nghost);
517 node.writeout_array.push_back(writeout);
518 node.number_of_fabs++;
521 std::vector<std::string> names;
523 names.push_back(name);
526 if (suffix.size() == 0)
527 for (
int j = 0; j < ncomp; j++)
528 names.push_back(amrex::Concatenate(name,j+1,3));
530 for (
int j = 0; j < ncomp; j++)
531 names.push_back(name + suffix[j]);
533 node.name_array.push_back(names);
536template<
class T,
int d>
540 int nghost, std::string name,
bool writeout,
bool evolving,
541 std::vector<std::string> )
546 int nlevs_max = maxLevel() + 1;
547 new_field.resize(nlevs_max);
548 m_basefields.push_back(
new Field<T>(new_field, geom, refRatio(), ncomp, nghost));
549 m_basefields.back()->evolving = evolving;
550 m_basefields.back()->writeout = writeout;
551 m_basefields.back()->setName(name);
552 m_basefields.back()->evolving = evolving;
557 int nlevs_max = maxLevel() + 1;
558 new_field.resize(nlevs_max);
559 m_basefields_cell.push_back(
new Field<T>(new_field, geom, refRatio(), ncomp, nghost));
560 m_basefields_cell.back()->evolving = evolving;
561 m_basefields_cell.back()->writeout = writeout;
562 m_basefields_cell.back()->setName(name);
563 m_basefields_cell.back()->evolving = evolving;
564 if (new_bc) m_basefields_cell.back()->setBC(new_bc);
565 else m_basefields_cell.back()->setBC(&bcnothing);
570 Util::Abort(
INFO,
"Only node and cell based fields can be added at this time");
581 AddField<T, Set::Hypercube::Node>(new_fab,
nullptr, ncomp, nghost,
"",
true, evolving);
588 AddField<T, Set::Hypercube::Node>(new_fab,
nullptr, ncomp, nghost, a_name,
true, evolving);
595 AddField<T, Set::Hypercube::Node>(new_fab,
nullptr, ncomp, nghost, a_name, writeout, evolving);
#define ALAMO_SINGLE_DEFINITION
std::vector< BaseField * > m_basefields
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.
void SetStopTime(Set::Scalar a_stop_time)
Utility to set the global stop time.
std::string GetFilename()
Simple getter to get filename.
void DynamicTimestep_SyncTimeStep(int lev, Set::Scalar dt_min)
Params for the dynamic timestp.
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
void regrid(int lbase, Set::Scalar time, bool initial=false) override
This overrides an AMReX method just to allow for explicit meshing when desired.
void SetThermoInt(int a_thermo_int)
Utility to set the frequency (in timesteps) of thermo data calculation.
std::vector< bool > extensives
virtual void TimeStepBegin(Set::Scalar, int)
This optional function is called at the beginning of every timestep, and can be used to complete addi...
std::vector< Set::Scalar * > vars
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
std::vector< BaseField * > m_basefields_cell
std::vector< BC::BC< Set::Scalar > * > physbc_array
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
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
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.
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...
std::vector< Set::Scalar > previous_timesteps
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?
std::vector< std::vector< std::string > > name_array
std::vector< std::string > names
void SetThermoPlotInt(int a_thermo_plot_int)
Utility to set the frequency (in timesteps) of thermo data writing to file.
amrex::Vector< amrex::Real > t_new
Keep track of current old simulation time on each level.
amrex::Vector< int > istep
Keep track of where each level is.
void DynamicTimestep_Reset()
std::vector< Set::Scalar > dt_limit_min
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.
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.
void DynamicTimestep_Update()
std::vector< bool > writeout_array
std::vector< Set::Field< Set::Scalar > * > fab_array
void SetFilename(std::string _plot_file)
Simple setter to set filename.
virtual void TimeStepComplete(Set::Scalar, int)
This optional function is called at the end of every timestep, and can be used to complete additional...
Collection of numerical integrator objects.
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, std::vector< std::string > suffix)
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, std::vector< std::string > suffix)
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
void Abort(const char *msg)
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
void Warning(std::string file, std::string func, int line, Args const &... args)
void Message(std::string file, std::string func, int line, Args const &... args)