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