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