40#ifndef INTEGRATOR_INTEGRATOR_H
41#define INTEGRATOR_INTEGRATOR_H
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>
75 :
public amrex::AmrCore
93 void Restart(std::string restartfile,
bool a_node =
false);
99 void SetFilename(std::string _plot_file) { plot_file = _plot_file; };
108 if (!explicitmesh.on)
109 AmrCore::regrid(lbase, time, initial);
116 if (!explicitmesh.on) AmrCore::InitFromScratch(time);
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);
131 for (
int ilev = 0; ilev < maxLevel(); ++ilev)
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);
180 const amrex::MFIter&,
184 if (thermo.number > 0)
185 Util::Warning(
INFO,
"integrated variables registered, but no integration implemented!");
195 clock_running =
false;
199 clock_running =
true;
206 void RegisterNewFab(
Set::Field<Set::Scalar>& new_fab,
int ncomp, std::string name,
bool writeout,
bool evolving=
true, std::vector<std::string> suffix = {});
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 = {});
222 template<
class T,
int d>
224 std::string,
bool writeout,
bool evolving, std::vector<std::string> suffix = {});
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);
241 void RegisterIntegratedVariable(
Set::Scalar* integrated_variable, std::string name,
bool extensive=
true);
247 void SetPlotInt(
int plot_int);
250 void SetThermoInt(
int a_thermo_int) { thermo.interval = a_thermo_int; }
273 if (!dynamictimestep.on)
return;
274 if (!dynamictimestep.dt_limit_min.size())
276 dynamictimestep.dt_limit_min.resize(max_level+1,
277 std::numeric_limits<Set::Scalar>::max());
279 dynamictimestep.dt_limit_min[lev] = std::min(dt_min,
280 dynamictimestep.dt_limit_min[lev]);
282 amrex::ParallelDescriptor::ReduceRealMin(dynamictimestep.dt_limit_min[lev]);
286 if (!dynamictimestep.on)
return;
287 dynamictimestep.previous_timesteps.clear();
291 if (!dynamictimestep.on)
return;
292 if (!dynamictimestep.dt_limit_min.size())
return;
296 if (amrex::ParallelDescriptor::IOProcessor())
299 if (dynamictimestep.previous_timesteps.size() > 0)
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();
307 Set::Scalar new_timestep = std::numeric_limits<Set::Scalar>::max();
308 for (
int lev = 0; lev <= this->max_level; lev++)
311 Set::Scalar dt_lev = dynamictimestep.dt_limit_min[lev];
315 for (
int ilev = lev; ilev > 0; ilev--) dt_lev *= (
Set::Scalar)(this->nsubsteps[ilev]);
319 new_timestep = std::min(new_timestep,dt_lev);
322 if (new_timestep < timestep_average)
324 dynamictimestep.previous_timesteps.clear();
326 final_timestep = new_timestep;
327 final_timestep = std::max(final_timestep,dynamictimestep.min);
328 final_timestep = std::min(final_timestep,dynamictimestep.max);
330 dynamictimestep.previous_timesteps.push_back(new_timestep);
334 final_timestep = timestep_average;
335 final_timestep = std::max(final_timestep,dynamictimestep.min);
336 final_timestep = std::min(final_timestep,dynamictimestep.max);
338 if ((
int)(dynamictimestep.previous_timesteps.size()) > dynamictimestep.nprevious)
339 dynamictimestep.previous_timesteps.erase(dynamictimestep.previous_timesteps.begin());
340 dynamictimestep.previous_timesteps.push_back(new_timestep);
344 amrex::ParallelDescriptor::Bcast(&final_timestep,1);
345 this->SetTimestep(final_timestep);
346 dynamictimestep.dt_limit_min.clear();
354 std::string plot_file{
"plt" };
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;
368 void FillPatch(
int lev, amrex::Real time,
369 amrex::Vector<std::unique_ptr<amrex::MultiFab>>& source_mf,
370 amrex::MultiFab& destination_multifab,
374 long CountCells(
int lev);
376 void TimeStep(
int lev, amrex::Real time,
int iteration);
378 void GetData(
const int lev,
const amrex::Real time, amrex::Vector<amrex::MultiFab*>& data, amrex::Vector<amrex::Real>& datatime);
380 std::vector<std::string> PlotFileName(
int lev, std::string prefix =
"")
const;
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;
393 amrex::Real timestep = NAN;
394 amrex::Vector<amrex::Real>
dt;
397 int max_plot_level = -1;
398 int print_ghost_nodes = 0;
399 int print_ghost_cells = 0;
402 int max_step = std::numeric_limits<int>::max();
403 amrex::Real tstart = 0;
404 amrex::Real stop_time = NAN;
407 bool integrate_variables_before_advance =
true;
408 bool integrate_variables_after_advance =
false;
410 bool clock_running =
true;
414 int number_of_fabs = 0;
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;
451 std::vector<Set::Scalar*>
vars;
458 int base_regrid_int = -1;
460 std::string restart_file_cell =
"";
461 std::string restart_file_node =
"";
465 std::vector<amrex::Box>
box;
472 int abort_on_nan =
true;
486 std::vector<std::string> suffix)
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);
496 cell.writeout_array.push_back(writeout);
497 cell.evolving_array.push_back(evolving);
498 cell.number_of_fabs++;
501 std::vector<std::string> names;
503 names.push_back(name);
506 if (suffix.size() == 0)
507 for (
int j = 0; j < ncomp; j++)
508 names.push_back(amrex::Concatenate(name,j+1,3));
510 for (
int j = 0; j < ncomp; j++)
511 names.push_back(name + suffix[j]);
513 cell.name_array.push_back(names);
526 std::vector<std::string> suffix)
528 BL_PROFILE(
"Integrator::RegisterNodalFab");
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++;
541 std::vector<std::string> names;
543 names.push_back(name);
546 if (suffix.size() == 0)
547 for (
int j = 0; j < ncomp; j++)
548 names.push_back(amrex::Concatenate(name,j+1,3));
550 for (
int j = 0; j < ncomp; j++)
551 names.push_back(name + suffix[j]);
553 node.name_array.push_back(names);
556template<
class T,
int d>
560 int nghost, std::string name,
bool writeout,
bool evolving,
561 std::vector<std::string> )
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;
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);
590 Util::Abort(
INFO,
"Only node and cell based fields can be added at this time");
601 AddField<T, Set::Hypercube::Node>(new_fab,
nullptr, ncomp, nghost,
"",
true, evolving);
608 AddField<T, Set::Hypercube::Node>(new_fab,
nullptr, ncomp, nghost, a_name,
true, evolving);
615 AddField<T, Set::Hypercube::Node>(new_fab,
nullptr, ncomp, nghost, a_name, writeout, evolving);
#define ALAMO_SINGLE_DEFINITION
std::vector< bool > evolving_array
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::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)
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)
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)