Alamo
Integrator.cpp
Go to the documentation of this file.
1///
2/// \file Integrator.cpp
3/// \brief Compute the volume integral of two multiplied Fourier series
4///
5
6#include "Integrator.H"
7#include "IO/FileNameParse.H"
8#include "IO/ParmParse.H"
9#include "Util/Util.H"
10#include "Unit/Unit.H"
11#include <numeric>
12
13
14
15namespace Integrator
16{
17
18Integrator::Integrator() : amrex::AmrCore()
19{
20 BL_PROFILE("Integrator::Integrator()");
21 {
22 // These are basic parameters that are, in
23 // general, common to all Alamo simulations.
25 // Number of iterations before ending (default is maximum possible int)
26 pp_query_default("max_step", max_step, 2147483647);
27 // Simulation time before ending
28 pp_query_required("stop_time", stop_time, Unit::Time());
29 // Nominal timestep on amrlev = 0
30 pp_query_required("timestep", timestep, Unit::Time());
31 pp_query("restart", restart_file_cell); // Name of restart file to READ from
32 pp_query("restart_cell", restart_file_cell); // Name of cell-fab restart file to read from
33 pp_query("restart_node", restart_file_node); // Name of node-fab restart file to read from
34 }
35
36 {
37 // This allows the user to ignore certain arguments that
38 // would otherwise cause problems.
39 // Most generally this is used in the event of a "above inputs
40 // specified but not used" error.
41 // The primary purpose of this was to fix those errors that arise
42 // in regression tests.
44 std::vector<std::string> ignore;
45 if (pp.contains("ignore")) Util::Message(INFO, "Ignore directive detected");
46 pp_queryarr("ignore", ignore); // Space-separated list of entries to ignore
47 for (unsigned int i = 0; i < ignore.size(); i++)
48 {
49 Util::Message(INFO, "ignoring ", ignore[i]);
50 pp.remove(ignore[i].c_str());
51 }
52 }
53 {
54 // These are parameters that are specific to
55 // the AMR/regridding part of the code.
56 IO::ParmParse pp("amr");
57 pp_query_default("regrid_int", regrid_int, 2); // Regridding interval in step numbers
58 pp_query_default("base_regrid_int", base_regrid_int, 0); // Regridding interval based on coarse level only
59 pp_query_default("plot_int", plot_int, -1); // Interval (in timesteps) between plotfiles (Default negative value will cause the plot interval to be ignored.)
60 pp.query_default("plot_dt", plot_dt, "-1.0", Unit::Time()); // Interval (in simulation time) between plotfiles (Default negative value will cause the plot dt to be ignored.)
61
62
63 // Output file: see IO::FileNameParse for wildcards and variable substitution
64 pp_query_default("plot_file", plot_file, "output");
65
66 pp_query_default("cell.all", cell.all, false); // Turn on to write all output in cell fabs (default: off)
67 pp_query_default("cell.any", cell.any, true); // Turn off to prevent any cell based output (default: on)
68 pp_query_default("node.all", node.all, false); // Turn on to write all output in node fabs (default: off)
69 pp_query_default("node.any", node.any, true); // Turn off to prevent any node based output (default: on)
70
71 pp_query_default("abort_on_nan",abort_on_nan, true); // Abort if a plotfile contains nan or inf.
72
73 Util::Assert(INFO, TEST(!(!cell.any && cell.all)));
74 Util::Assert(INFO, TEST(!(!node.any && node.all)));
75
76 pp_query_default("max_plot_level", max_plot_level, -1); // Specify a maximum level of refinement for output files (NO REFINEMENT)
77
78 IO::FileNameParse(plot_file);
79
80 nsubsteps.resize(maxLevel() + 1, 1);
81 int cnt = pp.countval("nsubsteps");
82 if (cnt != 0)
83 if (cnt == maxLevel()) {
84 pp_queryarr("nsubsteps", nsubsteps); // Number of substeps to take on each level (default: 2)
85 nsubsteps.insert(nsubsteps.begin(), 1);
86 nsubsteps.pop_back();
87 }
88 else if (cnt == 1)
89 {
90 int nsubsteps_all;
91 pp_query_required("nsubsteps", nsubsteps_all);// Number of substeps to take on each level (set all levels to this value)
92 for (int lev = 1; lev <= maxLevel(); ++lev) nsubsteps[lev] = nsubsteps_all;
93 }
94 else
95 Util::Abort(INFO, "number of nsubsteps input must equal either 1 or amr.max_level");
96 else
97 for (int lev = 1; lev <= maxLevel(); ++lev)
98 nsubsteps[lev] = MaxRefRatio(lev - 1);
99 }
100 {
101 IO::ParmParse pp("dynamictimestep");
102 // activate dynamic CFL-based timestep
103 pp_query("on",dynamictimestep.on);
104 if (dynamictimestep.on)
105 {
106 // how much information to print
107 pp_query_validate("verbose",dynamictimestep.verbose,{0,1});
108 // number of previous timesteps for rolling average
109 pp_query_default("nprevious",dynamictimestep.nprevious,5);
110 // dynamic teimstep CFL condition
111 pp_query_default("cfl",dynamictimestep.cfl,1.0);
112 // minimum timestep size allowed shen stepping dynamically
113 pp_query_default("min",dynamictimestep.min,timestep);
114 // maximum timestep size allowed shen stepping dynamically
115 pp_query_default("max",dynamictimestep.max,timestep);
116
117 Util::AssertException(INFO,TEST(dynamictimestep.max >= dynamictimestep.min));
118 }
119 }
120 {
121 // Information on how to generate thermodynamic
122 // data (to show up in thermo.dat)
123 IO::ParmParse pp("amr.thermo");
124 thermo.interval = 1; // Default: integrate every time.
125 pp_query_default("int", thermo.interval, 1); // Integration interval (1)
126 pp_query_default("plot_int", thermo.plot_int, -1); // Interval (in timesteps) between writing (Default negative value will cause the plot interval to be ignored.)
127 pp_query_default("plot_dt", thermo.plot_dt, "-1.0", Unit::Time()); // Interval (in simulation time) between writing (Default negative value will cause the plot dt to be ignored.)
128 }
129
130 {
131 // Instead of using AMR, prescribe an explicit, user-defined
132 // set of grids to work on. This is pretty much always used
133 // for testing purposes only.
134 IO::ParmParse pp("explicitmesh");
135 pp_query_default("on", explicitmesh.on, 0); // Use explicit mesh instead of AMR
136 if (explicitmesh.on)
137 {
138 for (int ilev = 0; ilev < maxLevel(); ++ilev)
139 {
140 std::string strlo = "lo" + std::to_string(ilev + 1);
141 std::string strhi = "hi" + std::to_string(ilev + 1);
142
143 Util::Assert(INFO, TEST(pp.contains(strlo.c_str())));
144 Util::Assert(INFO, TEST(pp.contains(strhi.c_str())));
145
146 amrex::Vector<int> lodata, hidata;
147 pp_queryarr(strlo.c_str(), lodata);
148 pp_queryarr(strhi.c_str(), hidata);
149 amrex::IntVect lo(AMREX_D_DECL(lodata[0], lodata[1], lodata[2]));
150 amrex::IntVect hi(AMREX_D_DECL(hidata[0], hidata[1], hidata[2]));
151
152 explicitmesh.box.push_back(amrex::Box(lo, hi));
153 }
154 }
155 }
156
157
158 {
159 //
160 // These parameters are NOT USED!
161 // They are shadow paramters for AMReX::TimeIntegration.
162 // The purpose here is:
163 // (1) to set default values
164 // (2) to prevent "unused input" warnings when integration is parsed, since it is parsed in the
165 // AMReX convention, not the Alamo convention.
166 //
167
168 IO::ParmParse pp("integration");
169 std::string str;
170 // Type of time integration to use (see amrex::TimeIntegrator for more details)
171 pp.query_validate("type", str, {"ForwardEuler","RungeKutta"});
172 if (str == "RungeKutta")
173 {
174 int type;
175 // If RungeKutta specified, which order to use (3=SSPRK3, 4=RK4)
176 pp.query_validate("rk.type", type, {1,2,3,4});
177 }
178 }
179
180 int nlevs_max = maxLevel() + 1;
181
182 istep.resize(nlevs_max, 0);
183
184 t_new.resize(nlevs_max, 0.0);
185 t_old.resize(nlevs_max, -1.e100);
186 SetTimestep(timestep);
187
188 plot_file = Util::GetFileName();
190}
191
192// Destructor
194{
195 if (Util::finalized)
196 {
197 std::cout << "!! ERROR !! Integrator destructor called after alamo has been finalized." << std::endl;
198 std::cout << " Behavior occurring after this is undefined." << std::endl;
199 std::abort();
200 }
201
202 // Close out the metadata file and mark completed.
204
205 // De-initialize all of the base fields and clear the arrays.
206 for (unsigned int i = 0; i < m_basefields.size(); i++) delete m_basefields[i];
207 for (unsigned int i = 0; i < m_basefields_cell.size(); i++) delete m_basefields_cell[i];
208 m_basefields.clear();
209 m_basefields_cell.clear();
210}
211
213{
214 BL_PROFILE("Integrator::SetTimestep");
215 int nlevs_max = maxLevel() + 1;
216 timestep = _timestep;
217 dt.resize(nlevs_max, 1.e100);
218 dt[0] = timestep;
219 for (int i = 1; i < nlevs_max; i++)
220 dt[i] = dt[i - 1] / (amrex::Real)nsubsteps[i];
221}
222void Integrator::SetPlotInt(int a_plot_int)
223{
224 BL_PROFILE("Integrator::SetPlotInt");
225 plot_int = a_plot_int;
226}
227
228/// \fn Integrator::MakeNewLevelFromCoarse
229/// \brief Wrapper to call FillCoarsePatch
230/// \note **THIS OVERRIDES A PURE VIRTUAL METHOD - DO NOT CHANGE**
231///
232void
233Integrator::MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray& cgrids, const amrex::DistributionMapping& dm)
234{
235 BL_PROFILE("Integrator::MakeNewLevelFromCoarse");
236
237 for (int n = 0; n < cell.number_of_fabs; n++)
238 {
239 const int ncomp = (*cell.fab_array[n])[lev - 1]->nComp();
240 const int nghost = (*cell.fab_array[n])[lev - 1]->nGrow();
241
242 (*cell.fab_array[n])[lev].reset(new amrex::MultiFab(cgrids, dm, ncomp, nghost));
243
244 (*cell.fab_array[n])[lev]->setVal(0.0);
245
246 FillCoarsePatch(lev, time, *cell.fab_array[n], *cell.physbc_array[n], 0, ncomp);
247 }
248
249 amrex::BoxArray ngrids = cgrids;
250 ngrids.convert(amrex::IntVect::TheNodeVector());
251
252 for (int n = 0; n < node.number_of_fabs; n++)
253 {
254 const int ncomp = (*node.fab_array[n])[lev - 1]->nComp();
255 const int nghost = (*node.fab_array[n])[lev - 1]->nGrow();
256
257 (*node.fab_array[n])[lev].reset(new amrex::MultiFab(ngrids, dm, ncomp, nghost));
258 (*node.fab_array[n])[lev]->setVal(0.0);
259
260 FillCoarsePatch(lev, time, *node.fab_array[n], *node.physbc_array[n], 0, ncomp);
261 }
262
263 for (unsigned int n = 0; n < m_basefields.size(); n++)
264 {
265 m_basefields[n]->MakeNewLevelFromCoarse(lev, time, cgrids, dm);
266 }
267 for (unsigned int n = 0; n < m_basefields_cell.size(); n++)
268 {
269 m_basefields_cell[n]->MakeNewLevelFromCoarse(lev, time, cgrids, dm);
270 }
271
272 Regrid(lev, time);
273}
274
275
276///
277/// RESETS ALL MULTIFABS AT A GIVEN LEVEL
278///
279/// (OVERRIDES PURE VIRTUAL METHOD - DO NOT CHANGE)
280///
281void
282Integrator::RemakeLevel(int lev, ///<[in] AMR Level
283 amrex::Real time, ///<[in] Simulation time
284 const amrex::BoxArray &cgrids, ///<[in] Coarse grids
285 const amrex::DistributionMapping &dm ///[in] Distribution mapping
286 )
287{
288 BL_PROFILE("Integrator::RemakeLevel");
289 for (int n = 0; n < cell.number_of_fabs; n++)
290 {
291 const int ncomp = (*cell.fab_array[n])[lev]->nComp();
292 const int nghost = (*cell.fab_array[n])[lev]->nGrow();
293
294 amrex::MultiFab new_state(cgrids, dm, ncomp, nghost);
295
296 new_state.setVal(0.0);
297 FillPatch(lev, time, *cell.fab_array[n], new_state, *cell.physbc_array[n], 0);
298 std::swap(new_state, *(*cell.fab_array[n])[lev]);
299 }
300
301 amrex::BoxArray ngrids = cgrids;
302 ngrids.convert(amrex::IntVect::TheNodeVector());
303
304 for (int n = 0; n < node.number_of_fabs; n++)
305 {
306 const int ncomp = (*node.fab_array[n])[lev]->nComp();
307 const int nghost = (*node.fab_array[n])[lev]->nGrow();
308
309 amrex::MultiFab new_state(ngrids, dm, ncomp, nghost);
310
311 new_state.setVal(0.0);
312 FillPatch(lev, time, *node.fab_array[n], new_state, *node.physbc_array[n], 0);
313 std::swap(new_state, *(*node.fab_array[n])[lev]);
314 }
315
316 for (unsigned int n = 0; n < m_basefields_cell.size(); n++)
317 {
318 m_basefields_cell[n]->RemakeLevel(lev, time, cgrids, dm);
319 }
320 for (unsigned int n = 0; n < m_basefields.size(); n++)
321 {
322 m_basefields[n]->RemakeLevel(lev, time, cgrids, dm);
323 }
324 Regrid(lev, time);
325}
326
327//
328// DELETE EVERYTHING
329//
330// (OVERRIDES PURE VIRTUAL METHOD - DO NOT CHANGE)
331//
332void
334{
335 BL_PROFILE("Integrator::ClearLevel");
336 for (int n = 0; n < cell.number_of_fabs; n++)
337 {
338 (*cell.fab_array[n])[lev].reset(nullptr);
339 }
340 for (int n = 0; n < node.number_of_fabs; n++)
341 {
342 (*node.fab_array[n])[lev].reset(nullptr);
343 }
344}
345
346//
347//
348//
349
350
351
352
353
354void
355Integrator::RegisterNewFab(Set::Field<Set::Scalar>& new_fab, BC::BC<Set::Scalar>* new_bc, int ncomp, int nghost, std::string name, bool writeout, bool evolving, std::vector<std::string> suffix)
356{
357 //Util::Warning(INFO, "RegisterNewFab is depricated. Please replace with AddField");
358 AddField<Set::Scalar, Set::Hypercube::Cell>(new_fab, new_bc, ncomp, nghost, name, writeout, evolving, suffix);
359}
360void
361Integrator::RegisterNewFab(Set::Field<Set::Scalar>& new_fab, int ncomp, std::string name, bool writeout, bool evolving, std::vector<std::string> suffix)
362{
363 //Util::Warning(INFO, "RegisterNewFab is depricated. Please replace with AddField");
364 AddField<Set::Scalar, Set::Hypercube::Cell>(new_fab, nullptr, ncomp, 0, name, writeout, evolving, suffix);
365}
366void
367Integrator::RegisterNodalFab(Set::Field<Set::Scalar>& new_fab, BC::BC<Set::Scalar>* new_bc, int ncomp, int nghost, std::string name, bool writeout, bool evolving, std::vector<std::string> suffix)
368{
369 //Util::Warning(INFO, "RegisterNodalFab is depricated. Please replace with AddField");
370 AddField<Set::Scalar, Set::Hypercube::Node>(new_fab, new_bc, ncomp, nghost, name, writeout, evolving,suffix);
371}
372void
373Integrator::RegisterNodalFab(Set::Field<Set::Scalar>& new_fab, int ncomp, int nghost, std::string name, bool writeout, bool evolving, std::vector<std::string> suffix)
374{
375 //Util::Warning(INFO, "RegisterNodalFab is depricated. Please replace with AddField");
376 AddField<Set::Scalar, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, name, writeout, evolving,suffix);
377}
378
379
380
381
382void // CUSTOM METHOD - CHANGEABLE
383Integrator::RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive)
384{
385 BL_PROFILE("Integrator::RegisterIntegratedVariable");
386 thermo.vars.push_back(integrated_variable);
387 thermo.names.push_back(name);
388 thermo.extensives.push_back(extensive);
389 thermo.number++;
390}
391
392long // CUSTOM METHOD - CHANGEABLE
394{
395 BL_PROFILE("Integrator::CountCells");
396 const int N = grids[lev].size();
397
398 long cnt = 0;
399
400 for (int i = 0; i < N; ++i)
401 {
402 cnt += grids[lev][i].numPts();
403 }
404
405 return cnt;
406}
407
408void // CUSTOM METHOD - CHANGEABLE
409Integrator::FillPatch(int lev, amrex::Real time,
410 amrex::Vector<std::unique_ptr<amrex::MultiFab>>& source_mf,
411 amrex::MultiFab& destination_mf,
412 BC::BC<Set::Scalar>& physbc, int icomp)
413{
414 BL_PROFILE("Integrator::FillPatch");
415 if (lev == 0)
416 {
417
418 amrex::Vector<amrex::MultiFab*> smf;
419 smf.push_back(source_mf[lev].get());
420 amrex::Vector<amrex::Real> stime;
421 stime.push_back(time);
422
423 physbc.define(geom[lev]);
424 amrex::FillPatchSingleLevel(destination_mf, // Multifab
425 time, // time
426 smf, // Vector<MultiFab*> &smf (CONST)
427 stime, // Vector<Real> &stime (CONST)
428 0, // scomp - Source component
429 icomp, // dcomp - Destination component
430 destination_mf.nComp(), // ncomp - Number of components
431 geom[lev], // Geometry (CONST)
432 physbc,
433 0); // BC
434 }
435 else
436 {
437 amrex::Vector<amrex::MultiFab*> cmf, fmf;
438 cmf.push_back(source_mf[lev - 1].get());
439 fmf.push_back(source_mf[lev].get());
440 amrex::Vector<amrex::Real> ctime, ftime;
441 ctime.push_back(time);
442 ftime.push_back(time);
443
444 physbc.define(geom[lev]);
445
446 amrex::Interpolater* mapper;
447
448 if (destination_mf.boxArray().ixType() == amrex::IndexType::TheNodeType())
449 mapper = &amrex::node_bilinear_interp;
450 else
451 mapper = &amrex::cell_cons_interp;
452
453 amrex::Vector<amrex::BCRec> bcs(destination_mf.nComp(), physbc.GetBCRec()); // todo
454 amrex::FillPatchTwoLevels(destination_mf, time, cmf, ctime, fmf, ftime,
455 0, icomp, destination_mf.nComp(), geom[lev - 1], geom[lev],
456 physbc, 0,
457 physbc, 0,
458 refRatio(lev - 1),
459 mapper, bcs, 0);
460 }
461}
462
463/// \fn Integrator::FillCoarsePatch
464/// \brief Fill a fab at current level with the data from one level up
465///
466/// \note This is a custom method and is changeable
467void
468Integrator::FillCoarsePatch(int lev, ///<[in] AMR level
469 amrex::Real time, ///<[in] Simulatinon time
470 Set::Field<Set::Scalar>& mf, ///<[in] Fab to fill
471 BC::BC<Set::Scalar>& physbc, ///<[in] BC object applying to Fab
472 int icomp, ///<[in] start component
473 int ncomp) ///<[in] end component (i.e. applies to components `icomp`...`ncomp`)
474{
475 BL_PROFILE("Integrator::FillCoarsePatch");
476 AMREX_ASSERT(lev > 0);
477 amrex::Vector<amrex::MultiFab*> cmf;
478 cmf.push_back(mf[lev - 1].get());
479 amrex::Vector<amrex::Real> ctime;
480 ctime.push_back(time);
481
482 physbc.define(geom[lev]);
483
484 amrex::Interpolater* mapper;
485 if (mf[lev]->boxArray().ixType() == amrex::IndexType::TheNodeType())
486 mapper = &amrex::node_bilinear_interp;
487 else
488 mapper = &amrex::cell_cons_interp;
489
490 amrex::Vector<amrex::BCRec> bcs(ncomp, physbc.GetBCRec());
491 amrex::InterpFromCoarseLevel(*mf[lev], time, *cmf[0], 0, icomp, ncomp, geom[lev - 1], geom[lev],
492 physbc, 0,
493 physbc, 0,
494 refRatio(lev - 1),
495 mapper, bcs, 0);
496}
497
498void
499Integrator::ErrorEst(int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow)
500{
501 BL_PROFILE("Integrator::ErrorEst");
502 TagCellsForRefinement(lev, tags, time, ngrow);
503}
504
505
506void
508{
509 BL_PROFILE("Integrator::InitData");
510
511 if (restart_file_cell == "" && restart_file_node == "")
512 {
513 const amrex::Real time = 0.0;
514 InitFromScratch(time);
515
516 for (int lev = finest_level - 1; lev >= 0; --lev)
517 {
518 if (lev < max_level) regrid(lev, 0.0);
519 for (int n = 0; n < cell.number_of_fabs; n++)
520 amrex::average_down(*(*cell.fab_array[n])[lev + 1], *(*cell.fab_array[n])[lev],
521 geom[lev + 1], geom[lev],
522 0, (*cell.fab_array[n])[lev]->nComp(), refRatio(lev));
523 }
524 SetFinestLevel(finest_level);
525 }
526 if (restart_file_cell != "")
527 {
528 Restart(restart_file_cell, false);
529 }
530 if (restart_file_node != "")
531 {
532 Restart(restart_file_node, true);
533 }
534
535 if (plot_int > 0 || plot_dt > 0.0) {
536 WritePlotFile();
537 }
538}
539
540void
541Integrator::Restart(const std::string dirname, bool a_nodal)
542{
543 BL_PROFILE("Integrator::Restart");
544
545 if (a_nodal && node.fab_array.size() == 0)
546 {
547 Util::Message(INFO, "Nothing here for nodal fabs");
548 return;
549 }
550 if (!a_nodal && cell.fab_array.size() == 0)
551 {
552 Util::Message(INFO, "Nothing here for cell-based fabs");
553 return;
554 }
555
556 std::string filename = dirname + "/Header";
557 std::string chkptfilename = dirname + "/Checkpoint";
558 amrex::VisMF::IO_Buffer io_buffer(amrex::VisMF::GetIOBufferSize());
559 amrex::Vector<char> fileCharPtr, chkptfileCharPtr;
560 amrex::ParallelDescriptor::ReadAndBcastFile(filename, fileCharPtr);
561 amrex::ParallelDescriptor::ReadAndBcastFile(chkptfilename, chkptfileCharPtr);
562 std::string fileCharPtrString(fileCharPtr.dataPtr());
563 std::string chkptfileCharPtrString(chkptfileCharPtr.dataPtr());
564 std::istringstream is(fileCharPtrString, std::istringstream::in);
565 std::istringstream chkpt_is(chkptfileCharPtrString, std::istringstream::in);
566
567 std::string line, word;
568
569 // Get version
570 std::getline(is, line);
571 Util::Message(INFO, "Version: ", line);
572
573 // Get number of fabs
574 int tmp_numfabs;
575 std::getline(is, line); tmp_numfabs = std::stoi(line);
576 Util::Message(INFO, "number of fabs:", tmp_numfabs);
577 std::vector<std::string> tmp_name_array;
578
579 for (int i = 0; i < tmp_numfabs; i++)
580 {
581 std::getline(is, line);
582 tmp_name_array.push_back(line);
583 }
584
585 // Dimension?
586 std::getline(is, line);
587 Util::Warning(INFO, "Dimension: " + line);
588
589 // Current time
590 Set::Scalar tmp_time = 0.0;
591 std::getline(is, line); tmp_time = std::stof(line); Util::Message(INFO, "Current time: ", tmp_time);
592 for (int i = 0; i < max_level + 1; i++)
593 {
594 t_new[i] = tmp_time; t_old[i] = tmp_time;
595 }
596
597 // AMR level
598 int tmp_max_level;
599 std::getline(is, line); tmp_max_level = std::stoi(line); Util::Message(INFO, "Max AMR level: ", line);
600 if (tmp_max_level > max_level)
601 Util::Abort(INFO, "The max level specified (", max_level, ") is smaller than the finest level in the restart file (", tmp_max_level, ")");
602 finest_level = tmp_max_level;
603 // Geometry ?
604 std::getline(is, line); Util::Message(INFO, "Input geometry: ", line);
605 std::getline(is, line); Util::Message(INFO, " ", line);
606
607 // Mesh refinement ratio?
608 std::getline(is, line); Util::Message(INFO, "Mesh refinement ratio: ", line);
609
610 // Domain
611 std::getline(is, line); Util::Warning(INFO, "Domain: ", line);
612
613 // Domain
614 std::getline(is, line);
615 std::vector<std::string> tmp_iters = Util::String::Split(line);
616 if (tmp_iters.size() != (unsigned int)(finest_level + 1)) Util::Abort(INFO, "Error reading in interation counts: line = ", line);
617 for (int lev = 0; lev <= finest_level; lev++) { istep[lev] = std::stoi(tmp_iters[lev]); Util::Message(INFO, "Iter on level ", lev, " = ", istep[lev]); }
618
619 amrex::Vector<amrex::MultiFab> tmpdata(tmp_max_level + 1);
620 int total_ncomp = 0;
621
622 if (a_nodal)
623 {
624 for (unsigned int i = 0; i < node.fab_array.size(); i++)
625 if (node.writeout_array[i])
626 total_ncomp += node.ncomp_array[i];
627 }
628 else
629 {
630 for (unsigned int i = 0; i < cell.fab_array.size(); i++)
631 if (cell.writeout_array[i])
632 total_ncomp += cell.ncomp_array[i];
633 }
634
635 int total_nghost = a_nodal ? 0 : cell.nghost_array[0];
636
637 for (int lev = 0; lev <= finest_level; lev++)
638 {
639 amrex::BoxArray tmp_ba;
640 tmp_ba.readFrom(chkpt_is);
641 SetBoxArray(lev, tmp_ba);
642 amrex::DistributionMapping tmp_dm(tmp_ba, amrex::ParallelDescriptor::NProcs());
643 SetDistributionMap(lev, tmp_dm);
644
645 if (a_nodal)
646 {
647 amrex::BoxArray ngrids = grids[lev];
648 ngrids.convert(amrex::IntVect::TheNodeVector());
649 //tmpdata[lev].define(ngrids, dmap[lev], total_ncomp, total_nghost);
650 }
651 else
652 {
653 //tmpdata[lev].define(grids[lev], dmap[lev], total_ncomp, total_nghost);
654 }
655 Util::Message(INFO,max_level);
656 Util::Message(INFO,finest_level);
657 Util::Message(INFO,lev,dirname);
658 Util::Message(INFO,lev,grids[lev]);
659 Util::Message(INFO,amrex::MultiFabFileFullPrefix(lev, dirname, "Level_", "Cell"));
660 Util::Message(INFO,total_ncomp);
661 Util::Message(INFO,total_nghost);
662 amrex::VisMF::Read(tmpdata[lev],
663 amrex::MultiFabFileFullPrefix(lev, dirname, "Level_", "Cell"));
664
665 if (a_nodal)
666 for (int i = 0; i < node.number_of_fabs; i++)
667 {
668 amrex::BoxArray ngrids = grids[lev];
669 ngrids.convert(amrex::IntVect::TheNodeVector());
670 (*node.fab_array[i])[lev].reset(new amrex::MultiFab(ngrids, dmap[lev], node.ncomp_array[i], node.nghost_array[i]));
671 (*node.fab_array[i])[lev]->setVal(0.);
672 }
673 else
674 for (int i = 0; i < cell.number_of_fabs; i++)
675 (*cell.fab_array[i])[lev].reset(new amrex::MultiFab(grids[lev], dmap[lev], cell.ncomp_array[i], cell.nghost_array[i]));
676 for (int i = 0; i < tmp_numfabs; i++)
677 {
678 bool match = false;
679 if (a_nodal)
680 {
681 for (int j = 0; j < node.number_of_fabs; j++)
682 {
683 if (tmp_name_array[i] == node.name_array[i][j])
684 {
685 match = true;
686 Util::Message(INFO, "Initializing ", node.name_array[i][j], "; nghost=", node.nghost_array[j], " with ", tmp_name_array[i]);
687 amrex::MultiFab::Copy(*((*node.fab_array[j])[lev]).get(), tmpdata[lev], i, 0, 1, total_nghost);
688 }
689 for (int k = 0; k < node.ncomp_array[j]; k++)
690 {
691 if (tmp_name_array[i] == node.name_array[j][k])
692 {
693 match = true;
694 Util::Message(INFO, "Initializing ", node.name_array[j][k], "; ncomp=", node.ncomp_array[j], "; nghost=", node.nghost_array[j], " with ", tmp_name_array[i]);
695 amrex::MultiFab::Copy(*((*node.fab_array[j])[lev]).get(), tmpdata[lev], i, k, 1, total_nghost);
696 }
697 }
698 Util::RealFillBoundary(*((*node.fab_array[j])[lev]).get(), geom[lev]);
699 }
700 }
701 else
702 {
703 for (int j = 0; j < cell.number_of_fabs; j++)
704 {
705 for (int k = 0; k < cell.ncomp_array[j]; k++)
706 {
707 if (tmp_name_array[i] == cell.name_array[j][k])
708 {
709 match = true;
710 Util::Message(INFO, "Initializing ", cell.name_array[j][k], "; ncomp=", cell.ncomp_array[j], "; nghost=", cell.nghost_array[j], " with ", tmp_name_array[i]);
711 amrex::MultiFab::Copy(*((*cell.fab_array[j])[lev]).get(), tmpdata[lev], i, k, 1, 0 /*cell.nghost_array[j]*/);
712 }
713 }
714 Util::RealFillBoundary(*(*cell.fab_array[j])[lev].get(), geom[lev]);
715 }
716 }
717 if (!match) Util::Warning(INFO, "Fab ", tmp_name_array[i], " is in the restart file, but there is no fab with that name here.");
718 }
719
720 for (unsigned int n = 0; n < m_basefields_cell.size(); n++)
721 {
722 m_basefields_cell[n]->MakeNewLevelFromScratch(lev, t_new[lev], grids[lev], dmap[lev]);
723 }
724 for (unsigned int n = 0; n < m_basefields.size(); n++)
725 {
726 m_basefields[n]->MakeNewLevelFromScratch(lev, t_new[lev], grids[lev], dmap[lev]);
727 }
728
729
730 for (int n = 0; n < cell.number_of_fabs; n++)
731 {
732 if (cell.writeout_array[n])
733 FillPatch(lev, t_new[lev], *cell.fab_array[n], *(*cell.fab_array[n])[lev], *cell.physbc_array[n], 0);
734 }
735
736 }
737
738 SetFinestLevel(max_level);
739}
740
741void
742Integrator::MakeNewLevelFromScratch(int lev, amrex::Real t, const amrex::BoxArray& cgrids,
743 const amrex::DistributionMapping& dm)
744{
745 BL_PROFILE("Integrator::MakeNewLevelFromScratch");
746 for (int n = 0; n < cell.number_of_fabs; n++)
747 {
748 (*cell.fab_array[n])[lev].reset(new amrex::MultiFab(cgrids, dm, cell.ncomp_array[n], cell.nghost_array[n]));
749 (*cell.fab_array[n])[lev]->setVal(0.0);
750 }
751
752 amrex::BoxArray ngrids = cgrids;
753 ngrids.convert(amrex::IntVect::TheNodeVector());
754 for (int n = 0; n < node.number_of_fabs; n++)
755 {
756 (*node.fab_array[n])[lev].reset(new amrex::MultiFab(ngrids, dm, node.ncomp_array[n], node.nghost_array[n]));
757 (*node.fab_array[n])[lev]->setVal(0.0);
758 }
759 for (unsigned int n = 0; n < m_basefields_cell.size(); n++)
760 {
761 m_basefields_cell[n]->MakeNewLevelFromScratch(lev, t, cgrids, dm);
762 }
763 for (unsigned int n = 0; n < m_basefields.size(); n++)
764 {
765 m_basefields[n]->MakeNewLevelFromScratch(lev, t, cgrids, dm);
766 }
767
768 t_new[lev] = t;
769 t_old[lev] = t - dt[lev];
770
771 Initialize(lev);
772
773 for (int n = 0; n < cell.number_of_fabs; n++)
774 {
775 cell.physbc_array[n]->define(geom[lev]);
776 cell.physbc_array[n]->FillBoundary(*(*cell.fab_array[n])[lev], 0, 0, t, 0);
777 }
778
779 //for (int n = 0 ; n < node.number_of_fabs; n++)
780 //{
781 // bcnothing->define(geom[lev]);
782 // for (amrex::MFIter mfi(*(*node.fab_array[n])[lev],true); mfi.isValid(); ++mfi)
783 // {
784 // amrex::BaseFab<Set::Scalar> &patch = (*(*node.fab_array[n])[lev])[mfi];
785 // const amrex::Box& box = mfi.tilebox();
786 // bcnothing->FillBoundary(patch,box,0,0,0,t);
787 // }
788 //}
789
790 for (unsigned int n = 0; n < m_basefields_cell.size(); n++)
791 {
792 m_basefields_cell[n]->FillBoundary(lev, t);
793 }
794 for (unsigned int n = 0; n < m_basefields.size(); n++)
795 {
796 m_basefields[n]->FillBoundary(lev, t);
797 }
798}
799
800std::vector<std::string>
801Integrator::PlotFileName(int lev, std::string prefix) const
802{
803 BL_PROFILE("Integrator::PlotFileName");
804 std::vector<std::string> name;
805 name.push_back(plot_file + "/" + prefix + "/");
806 name.push_back(amrex::Concatenate("", lev, 5));
807 return name;
808}
809
810void
811Integrator::WritePlotFile(bool initial) const
812{
813 WritePlotFile(t_new[0], istep, initial, "");
814}
815void
816Integrator::WritePlotFile(std::string prefix, Set::Scalar time, int step) const
817{
818 amrex::Vector<int> istep(max_level + 1, step);
819 WritePlotFile(time, istep, false, prefix);
820}
821
822void
823Integrator::WritePlotFile(Set::Scalar time, amrex::Vector<int> iter, bool initial, std::string prefix) const
824{
825 BL_PROFILE("Integrator::WritePlotFile");
826 int nlevels = finest_level + 1;
827 if (max_plot_level >= 0) nlevels = std::min(nlevels, max_plot_level);
828
829 int ccomponents = 0, ncomponents = 0, bfcomponents_cell = 0, bfcomponents = 0;
830 amrex::Vector<std::string> cnames, nnames, bfnames_cell, bfnames;
831 for (int i = 0; i < cell.number_of_fabs; i++)
832 {
833 if (!cell.writeout_array[i]) continue;
834 ccomponents += cell.ncomp_array[i];
835 if (cell.ncomp_array[i] > 1)
836 for (int j = 0; j < cell.ncomp_array[i]; j++)
837 cnames.push_back(cell.name_array[i][j]);
838 else
839 cnames.push_back(cell.name_array[i][0]);
840 }
841 for (int i = 0; i < node.number_of_fabs; i++)
842 {
843 if (!node.writeout_array[i]) continue;
844 ncomponents += node.ncomp_array[i];
845 if (node.ncomp_array[i] > 1)
846 for (int j = 0; j < node.ncomp_array[i]; j++)
847 nnames.push_back(node.name_array[i][j]);
848 else
849 nnames.push_back(node.name_array[i][0]);
850 }
851 for (unsigned int i = 0; i < m_basefields_cell.size(); i++)
852 {
853 if (m_basefields_cell[i]->writeout)
854 {
855 bfcomponents_cell += m_basefields_cell[i]->NComp();
856 for (int j = 0; j < m_basefields_cell[i]->NComp(); j++)
857 bfnames_cell.push_back(m_basefields_cell[i]->Name(j));
858 }
859 }
860 for (unsigned int i = 0; i < m_basefields.size(); i++)
861 {
862 if (m_basefields[i]->writeout)
863 {
864 bfcomponents += m_basefields[i]->NComp();
865 for (int j = 0; j < m_basefields[i]->NComp(); j++)
866 bfnames.push_back(m_basefields[i]->Name(j));
867 }
868 }
869
870 amrex::Vector<amrex::MultiFab> cplotmf(nlevels), nplotmf(nlevels);
871
872 bool do_cell_plotfile = (ccomponents + bfcomponents_cell > 0 || (ncomponents + bfcomponents > 0 && cell.all)) && cell.any;
873 bool do_node_plotfile = (ncomponents + bfcomponents > 0 || (ccomponents + bfcomponents_cell > 0 && node.all)) && node.any;
874
875 for (int ilev = 0; ilev < nlevels; ++ilev)
876 {
877 if (do_cell_plotfile)
878 {
879 int ncomp = ccomponents + bfcomponents_cell;
880 if (cell.all) ncomp += ncomponents + bfcomponents;
881 cplotmf[ilev].define(grids[ilev], dmap[ilev], ncomp, 0);
882
883 int n = 0;
884 int cnames_cnt = 0;
885 for (int i = 0; i < cell.number_of_fabs; i++)
886 {
887 if (!cell.writeout_array[i]) continue;
888 if ((*cell.fab_array[i])[ilev]->contains_nan())
889 {
890 if (abort_on_nan) Util::Abort(INFO, cnames[cnames_cnt], " contains nan (i=", i, ")");
891 else Util::Warning(INFO, cnames[cnames_cnt], " contains nan (i=", i, ")");
892 }
893 if ((*cell.fab_array[i])[ilev]->contains_inf())
894 {
895 if (abort_on_nan) Util::Abort(INFO, cnames[cnames_cnt], " contains inf (i=", i, ")");
896 else Util::Warning(INFO, cnames[cnames_cnt], " contains inf (i=", i, ")");
897 }
898 cnames_cnt++;
899 amrex::MultiFab::Copy(cplotmf[ilev], *(*cell.fab_array[i])[ilev], 0, n, cell.ncomp_array[i], 0);
900 n += cell.ncomp_array[i];
901 }
902 for (unsigned int i = 0; i < m_basefields_cell.size(); i++)
903 {
904 if (m_basefields_cell[i]->writeout)
905 {
906 m_basefields_cell[i]->Copy(ilev, cplotmf[ilev], n, 0);
907 n += m_basefields_cell[i]->NComp();
908 }
909 }
910
911 if (cell.all)
912 {
913 int nnames_cnt = 0;
914 for (int i = 0; i < node.number_of_fabs; i++)
915 {
916 if (!node.writeout_array[i]) continue;
917 if ((*node.fab_array[i])[ilev]->contains_nan())
918 {
919 if (abort_on_nan) Util::Abort(INFO, nnames[nnames_cnt], " contains nan (i=", i, ")");
920 else Util::Warning(INFO, nnames[nnames_cnt], " contains nan (i=", i, ")");
921 }
922 if ((*node.fab_array[i])[ilev]->contains_inf())
923 {
924 if (abort_on_nan) Util::Abort(INFO, nnames[nnames_cnt], " contains inf (i=", i, ")");
925 else Util::Warning(INFO, nnames[nnames_cnt], " contains inf (i=", i, ")");
926 }
927 nnames_cnt++;
928 amrex::average_node_to_cellcenter(cplotmf[ilev], n, *(*node.fab_array[i])[ilev], 0, node.ncomp_array[i], 0);
929 n += node.ncomp_array[i];
930 }
931 if (bfcomponents > 0)
932 {
933 amrex::BoxArray ngrids = grids[ilev];
934 ngrids.convert(amrex::IntVect::TheNodeVector());
935 amrex::MultiFab bfplotmf(ngrids, dmap[ilev], bfcomponents, 0);
936 int ctr = 0;
937 for (unsigned int i = 0; i < m_basefields.size(); i++)
938 {
939 if (m_basefields[i]->writeout)
940 {
941 m_basefields[i]->Copy(ilev, bfplotmf, ctr, 0);
942 ctr += m_basefields[i]->NComp();
943 }
944 }
945 amrex::average_node_to_cellcenter(cplotmf[ilev], n, bfplotmf, 0, bfcomponents);
946 n += bfcomponents;
947 }
948 }
949 }
950
951 if (do_node_plotfile)
952 {
953 amrex::BoxArray ngrids = grids[ilev];
954 ngrids.convert(amrex::IntVect::TheNodeVector());
955 int ncomp = ncomponents + bfcomponents;
956 if (node.all) ncomp += ccomponents + bfcomponents_cell;
957 nplotmf[ilev].define(ngrids, dmap[ilev], ncomp, 0);
958
959 int n = 0;
960 for (int i = 0; i < node.number_of_fabs; i++)
961 {
962 if (!node.writeout_array[i]) continue;
963 if ((*node.fab_array[i])[ilev]->contains_nan()) Util::Warning(INFO, nnames[i], " contains nan (i=", i, "). Resetting to zero.");
964 if ((*node.fab_array[i])[ilev]->contains_inf()) Util::Abort(INFO, nnames[i], " contains inf (i=", i, ")");
965 amrex::MultiFab::Copy(nplotmf[ilev], *(*node.fab_array[i])[ilev], 0, n, node.ncomp_array[i], 0);
966 n += node.ncomp_array[i];
967 }
968 for (unsigned int i = 0; i < m_basefields.size(); i++)
969 {
970 if (m_basefields[i]->writeout)
971 {
972 m_basefields[i]->Copy(ilev, nplotmf[ilev], n, 0);
973 n += m_basefields[i]->NComp();
974 }
975 }
976
977 if (node.all)
978 {
979 for (int i = 0; i < cell.number_of_fabs; i++)
980 {
981 if (!cell.writeout_array[i]) continue;
982 if ((*cell.fab_array[i])[ilev]->contains_nan()) Util::Warning(INFO, cnames[i], " contains nan (i=", i, "). Resetting to zero.");
983 if ((*cell.fab_array[i])[ilev]->contains_inf()) Util::Abort(INFO, cnames[i], " contains inf (i=", i, ")");
984 if ((*cell.fab_array[i])[ilev]->nGrow() == 0)
985 {
986 if (initial) Util::Warning(INFO, cnames[i], " has no ghost cells and will not be included in nodal output");
987 continue;
988 }
989 Util::AverageCellcenterToNode(nplotmf[ilev], n, *(*cell.fab_array[i])[ilev], 0, cell.ncomp_array[i]);
990 n += cell.ncomp_array[i];
991 }
992
993 if (bfcomponents_cell > 0)
994 {
995 amrex::BoxArray cgrids = grids[ilev];
996 amrex::MultiFab bfplotmf(cgrids, dmap[ilev], bfcomponents_cell, 0);
997 int ctr = 0;
998 for (unsigned int i = 0; i < m_basefields_cell.size(); i++)
999 {
1000 if (m_basefields_cell[i]->writeout)
1001 {
1002 m_basefields_cell[i]->Copy(ilev, bfplotmf, ctr, 0);
1003 ctr += m_basefields_cell[i]->NComp();
1004 }
1005 }
1006 Util::AverageCellcenterToNode(nplotmf[ilev], n, bfplotmf, 0, bfcomponents_cell);
1007 n += bfcomponents_cell;
1008 }
1009 }
1010 }
1011 }
1012
1013 std::vector<std::string> plotfilename = PlotFileName(istep[0], prefix);
1014 if (initial) plotfilename[1] = plotfilename[1] + "init";
1015
1016 if (do_cell_plotfile)
1017 {
1018 amrex::Vector<std::string> allnames = cnames;
1019 allnames.insert(allnames.end(), bfnames_cell.begin(), bfnames_cell.end());
1020 if (cell.all) {
1021 allnames.insert(allnames.end(), nnames.begin(), nnames.end());
1022 allnames.insert(allnames.end(), bfnames.begin(), bfnames.end());
1023 }
1024 WriteMultiLevelPlotfile(plotfilename[0] + plotfilename[1] + "cell", nlevels, amrex::GetVecOfConstPtrs(cplotmf), allnames,
1025 Geom(), time, iter, refRatio());
1026
1027 std::ofstream chkptfile;
1028 chkptfile.open(plotfilename[0] + plotfilename[1] + "cell/Checkpoint");
1029 for (int i = 0; i <= max_level; i++) boxArray(i).writeOn(chkptfile);
1030 chkptfile.close();
1031 }
1032
1033 if (do_node_plotfile)
1034 {
1035 amrex::Vector<std::string> allnames = nnames;
1036 allnames.insert(allnames.end(), bfnames.begin(), bfnames.end());
1037 if (node.all) allnames.insert(allnames.end(), cnames.begin(), cnames.end());
1038 WriteMultiLevelPlotfile(plotfilename[0] + plotfilename[1] + "node", nlevels, amrex::GetVecOfConstPtrs(nplotmf), allnames,
1039 Geom(), time, iter, refRatio());
1040
1041 std::ofstream chkptfile;
1042 chkptfile.open(plotfilename[0] + plotfilename[1] + "node/Checkpoint");
1043 for (int i = 0; i <= max_level; i++) boxArray(i).writeOn(chkptfile);
1044 chkptfile.close();
1045 }
1046
1047 if (amrex::ParallelDescriptor::IOProcessor())
1048 {
1049 std::ofstream coutfile, noutfile;
1050 if (istep[0] == 0)
1051 {
1052 if (do_cell_plotfile) coutfile.open(plot_file + "/celloutput.visit", std::ios_base::out);
1053 if (do_node_plotfile) noutfile.open(plot_file + "/nodeoutput.visit", std::ios_base::out);
1054 }
1055 else
1056 {
1057 if (do_cell_plotfile) coutfile.open(plot_file + "/celloutput.visit", std::ios_base::app);
1058 if (do_node_plotfile) noutfile.open(plot_file + "/nodeoutput.visit", std::ios_base::app);
1059 }
1060 if (do_cell_plotfile) coutfile << plotfilename[1] + "cell" + "/Header" << std::endl;
1061 if (do_node_plotfile) noutfile << plotfilename[1] + "node" + "/Header" << std::endl;
1062 }
1063}
1064
1065void
1067{
1068 BL_PROFILE("Integrator::Evolve");
1069 amrex::Real cur_time = t_new[0];
1070 int last_plot_file_step = 0;
1071
1072 for (int step = istep[0]; step < max_step && cur_time < stop_time; ++step)
1073 {
1074 if (amrex::ParallelDescriptor::IOProcessor()) {
1075 std::cout << "\nSTEP " << step + 1 << " starts ..." << std::endl;
1076 }
1077 int lev = 0;
1078 int iteration = 1;
1079 TimeStepBegin(cur_time, step);
1080 if (integrate_variables_before_advance) IntegrateVariables(cur_time, step);
1081 TimeStep(lev, cur_time, iteration);
1082 if (integrate_variables_after_advance) IntegrateVariables(cur_time, step);
1083 TimeStepComplete(cur_time, step);
1084 cur_time += dt[0];
1085
1086 if (amrex::ParallelDescriptor::IOProcessor()) {
1087 std::cout << "STEP " << step + 1 << " ends."
1088 << " TIME = " << cur_time << " DT = " << dt[0]
1089 << std::endl;
1090 }
1091
1092 // sync up time
1093 for (int lev = 0; lev <= finest_level; ++lev) {
1094 t_new[lev] = cur_time;
1095 }
1096
1097 if (plot_int > 0 && (step + 1) % plot_int == 0) {
1098 last_plot_file_step = step + 1;
1099 WritePlotFile();
1100 IO::WriteMetaData(plot_file, IO::Status::Running, (int)(100.0 * cur_time / stop_time));
1101 }
1102 else if (std::fabs(std::remainder(cur_time, plot_dt)) < 0.5 * dt[0])
1103 {
1104 last_plot_file_step = step + 1;
1105 WritePlotFile();
1106 IO::WriteMetaData(plot_file, IO::Status::Running, (int)(100.0 * cur_time / stop_time));
1107 }
1108
1109 if (cur_time >= stop_time - 1.e-6 * dt[0]) break;
1110 }
1111 if (plot_int > 0 && istep[0] > last_plot_file_step) {
1112 WritePlotFile();
1113 }
1114}
1115
1116void
1117Integrator::IntegrateVariables(amrex::Real time, int step)
1118{
1119 BL_PROFILE("Integrator::IntegrateVariables");
1120 if (!thermo.number) return;
1121
1122 if ((thermo.interval > 0 && (step) % thermo.interval == 0) ||
1123 ((thermo.dt > 0.0) && (std::fabs(std::remainder(time, plot_dt)) < 0.5 * dt[0])))
1124 {
1125 // Zero out all variables
1126 for (int i = 0; i < thermo.number; i++)
1127 {
1128 if (thermo.extensives[i]) *thermo.vars[i] = 0;
1129 }
1130
1131 // All levels except the finest
1132 for (int ilev = 0; ilev < max_level; ilev++)
1133 {
1134 const amrex::BoxArray& cfba = amrex::coarsen(grids[ilev + 1], refRatio(ilev));
1135
1136#ifdef OMP
1137#pragma omp parallel
1138#endif
1139 for (amrex::MFIter mfi(grids[ilev], dmap[ilev], true); mfi.isValid(); ++mfi)
1140 {
1141 const amrex::Box& box = mfi.tilebox();
1142 const amrex::BoxArray& comp = amrex::complementIn(box, cfba);
1143
1144 for (int i = 0; i < comp.size(); i++)
1145 {
1146 Integrate(ilev, time, step,
1147 mfi, comp[i]);
1148 }
1149 }
1150 }
1151 // Now do the finest level
1152 {
1153#ifdef OMP
1154#pragma omp parallel
1155#endif
1156 for (amrex::MFIter mfi(grids[max_level], dmap[max_level], true); mfi.isValid(); ++mfi)
1157 {
1158 const amrex::Box& box = mfi.tilebox();
1159 Integrate(max_level, time, step, mfi, box);
1160 }
1161 }
1162
1163 // Sum up across all processors
1164 for (int i = 0; i < thermo.number; i++)
1165 {
1166 if (thermo.extensives[i])
1167 amrex::ParallelDescriptor::ReduceRealSum(*thermo.vars[i]);
1168 }
1169 }
1170 if (amrex::ParallelDescriptor::IOProcessor() &&
1171 (
1172 (thermo.plot_int > 0 && step % thermo.plot_int == 0) ||
1173 (thermo.plot_dt > 0.0 && std::fabs(std::remainder(time, thermo.plot_dt)) < 0.5 * dt[0])
1174 ))
1175 {
1176 std::ofstream outfile;
1177 if (step == 0)
1178 {
1179 outfile.open(plot_file + "/thermo.dat", std::ios_base::out);
1180 outfile << "time";
1181 for (int i = 0; i < thermo.number; i++)
1182 outfile << "\t" << thermo.names[i];
1183 outfile << std::endl;
1184 }
1185 else outfile.open(plot_file + "/thermo.dat", std::ios_base::app);
1186 outfile << time;
1187 for (int i = 0; i < thermo.number; i++)
1188 outfile << "\t" << *thermo.vars[i];
1189 outfile << std::endl;
1190 outfile.close();
1191 }
1192
1193}
1194
1195
1196void
1197Integrator::TimeStep(int lev, amrex::Real time, int /*iteration*/)
1198{
1199 BL_PROFILE("Integrator::TimeStep");
1200 if (base_regrid_int <= 0 || istep[0] % base_regrid_int == 0)
1201 {
1202 if (regrid_int > 0 || base_regrid_int > 0) // We may need to regrid
1203 {
1204 static amrex::Vector<int> last_regrid_step(max_level + 1, 0);
1205
1206 // regrid doesn't change the base level, so we don't regrid on max_level
1207 if (lev < max_level && istep[lev] > last_regrid_step[lev])
1208 {
1209 if (istep[lev] % regrid_int == 0)
1210 {
1211 regrid(lev, time, false);
1212 }
1213 }
1214 }
1215 }
1216 SetFinestLevel(finest_level);
1217
1218 if (Verbose() && amrex::ParallelDescriptor::IOProcessor()) {
1219 std::cout << "[Level " << lev
1220 << " step " << istep[lev] + 1 << "] ";
1221 std::cout << "ADVANCE with dt = "
1222 << dt[lev]
1223 << std::endl;
1224 }
1225
1226 for (int n = 0; n < cell.number_of_fabs; n++)
1227 if (cell.evolving_array[n]) FillPatch(lev, time, *cell.fab_array[n], *(*cell.fab_array[n])[lev], *cell.physbc_array[n], 0);
1228 for (int n = 0; n < node.number_of_fabs; n++)
1229 if (node.evolving_array[n]) FillPatch(lev, time, *node.fab_array[n], *(*node.fab_array[n])[lev], *node.physbc_array[n], 0);
1230 for (unsigned int n = 0; n < m_basefields_cell.size(); n++)
1231 if (m_basefields_cell[n]->evolving) m_basefields_cell[n]->FillPatch(lev, time);
1232 for (unsigned int n = 0; n < m_basefields.size(); n++)
1233 if (m_basefields[n]->evolving) m_basefields[n]->FillPatch(lev, time);
1234
1235 Advance(lev, time, dt[lev]);
1236 ++istep[lev];
1237
1238 if (Verbose() && amrex::ParallelDescriptor::IOProcessor())
1239 {
1240 std::cout << "[Level " << lev
1241 << " step " << istep[lev] << "] ";
1242 std::cout << "Advanced "
1243 << CountCells(lev)
1244 << " cells"
1245 << std::endl;
1246 }
1247
1248 if (lev < finest_level)
1249 {
1250 for (int i = 1; i <= nsubsteps[lev + 1]; ++i)
1251 TimeStep(lev + 1, time + (i - 1) * dt[lev + 1], i);
1252
1253 for (int n = 0; n < cell.number_of_fabs; n++)
1254 {
1255 amrex::average_down(*(*cell.fab_array[n])[lev + 1], *(*cell.fab_array[n])[lev],
1256 geom[lev + 1], geom[lev],
1257 0, (*cell.fab_array[n])[lev]->nComp(), refRatio(lev));
1258 }
1259 for (int n = 0; n < node.number_of_fabs; n++)
1260 {
1261 amrex::average_down(*(*node.fab_array[n])[lev + 1], *(*node.fab_array[n])[lev],
1262 0, (*node.fab_array[n])[lev]->nComp(), refRatio(lev));
1263 }
1264 for (unsigned int n = 0; n < m_basefields_cell.size(); n++)
1265 {
1266 if (m_basefields_cell[n]->evolving)
1267 m_basefields_cell[n]->AverageDown(lev, refRatio(lev));
1268 }
1269 for (unsigned int n = 0; n < m_basefields.size(); n++)
1270 {
1271 if (m_basefields[n]->evolving)
1272 m_basefields[n]->AverageDown(lev, refRatio(lev));
1273 }
1274
1275 }
1276}
1277}
std::time_t t
#define pp_query_validate(...)
Definition ParmParse.H:103
#define pp_queryarr(...)
Definition ParmParse.H:105
#define pp_query_required(...)
Definition ParmParse.H:101
#define pp_query(...)
Definition ParmParse.H:108
#define pp_query_default(...)
Definition ParmParse.H:102
#define TEST(x)
Definition Util.H:22
#define INFO
Definition Util.H:21
Definition BC.H:43
void define(const amrex::Geometry &a_geom)
Definition BC.H:48
virtual amrex::BCRec GetBCRec()=0
bool contains(std::string name)
Definition ParmParse.H:173
int query_default(std::string name, T &value, T defaultvalue)
Definition ParmParse.H:293
int query_validate(std::string name, int &value, std::vector< int > possibleintvals)
Definition ParmParse.H:336
void WritePlotFile(bool initial=false) const
Integrator()
This is the constructor for the intetgrator class, which reads timestep information,...
std::vector< std::string > PlotFileName(int lev, std::string prefix="") const
void FillCoarsePatch(int lev, amrex::Real time, Set::Field< Set::Scalar > &mf, BC::BC< Set::Scalar > &physbc, int icomp, int ncomp)
Fill a fab at current level with the data from one level up.
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={})
Add a new node-based scalar field.
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={})
Add a new cell-based scalar field.
void Restart(std::string restartfile, bool a_node=false)
Read in output from previous simulation and start simulation at that point - Not currently tested.
void TimeStep(int lev, amrex::Real time, int iteration)
Timestep marching.
virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
virtual void ClearLevel(int lev) override
long CountCells(int lev)
Simple utility to count cells.
void FillPatch(int lev, amrex::Real time, amrex::Vector< std::unique_ptr< amrex::MultiFab > > &source_mf, amrex::MultiFab &destination_multifab, BC::BC< Set::Scalar > &physbc, int icomp)
This is the function that is responsible for updating patch data.
virtual ~Integrator()
Virtual destructure; make sure delete any pointers that you create here.
void InitData()
Front-end method to initialize simulation on all levels.
void Evolve()
Front-end method to start simulation.
void SetPlotInt(int plot_int)
Utility to set the frequency (in timesteps) of plotfile dumping.
virtual void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
void IntegrateVariables(Set::Scalar cur_time, int step)
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Register a variable to be integrated over the spatial domain using the Integrate function.
void SetTimestep(Set::Scalar _timestep)
Utility to set the coarse-grid timestep.
virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Wrapper to call FillCoarsePatch.
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
RESETS ALL MULTIFABS AT A GIVEN LEVEL.
void WriteMetaData(std::string plot_file, Status status, int per)
void FileNameParse(std::string &filename)
Internal function to do processing of the file name.
@ Complete
@ Running
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:514
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:474
amrex::Real Scalar
Definition Base.H:18
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
AMREX_FORCE_INLINE std::vector< std::string > Split(std::string &str, const char delim=' ')
Definition String.H:138
std::string GetFileName()
Definition Util.cpp:27
AMREX_FORCE_INLINE void AssertException(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:251
bool finalized
Definition Util.cpp:25
void Abort(const char *msg)
Definition Util.cpp:243
void AverageCellcenterToNode(amrex::MultiFab &node_mf, const int &dcomp, const amrex::MultiFab &cell_mf, const int &scomp, const int &ncomp)
Definition Util.cpp:365
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:55
void Warning(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:199
void Message(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:126
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T > > &a_mf, const amrex::Geometry &, const int nghost=2)
Definition Util.H:335
static Unit Time()
Definition Unit.H:199