Line data Source code
1 : //
2 : // Pure abstract class for managing data structures, time integration (with substepping),
3 : // mesh refinement, and I/O.
4 : //
5 : // Native input file parameters:
6 : //
7 : // .. code-block:: make
8 : //
9 : // max_step = [maximum number of timesteps]
10 : // stop_time = [maximum simulation time]
11 : // timestep = [time step for coarsest level]
12 : //
13 : // amr.regrid_int = [number of timesteps between regridding]
14 : // amr.plot_int = [number of timesteps between dumping output]
15 : // amr.plot_file = [base name of output directory]
16 : //
17 : // amr.nsubsteps = [number of temporal substeps at each level. This can be
18 : // either a single int (which is then applied to every refinement
19 : // level) or an array of ints (equal to amr.max_level)
20 : // corresponding to the refinement for each level.]
21 : //
22 : // Inherited input file parameters (from amrex AmrMesh class):
23 : //
24 : // .. code-block:: make
25 : //
26 : // amr.v = [verbosity level]
27 : // amr.max_level = [maximum level of refinement]
28 : // amr.n_proper =
29 : // amr.grid_eff =
30 : // amr.n_error_buff =
31 : // amr.ref_ratio_vect = [refinement ratios in each direction]
32 : // amr.ref_ratio = [refinement ratio in all directions (cannot be used with ref_ratio_vect)]
33 : // amr.max_grid_x =
34 : // amr.blocking_factor =
35 : // amr.n_cell = [number of cells on coarsest level]
36 : // amr.refine_grid_layout =
37 : // amr.check_input =
38 : //
39 :
40 : #ifndef INTEGRATOR_INTEGRATOR_H
41 : #define INTEGRATOR_INTEGRATOR_H
42 :
43 : #include <chrono>
44 : #include <ctime>
45 : #include <string>
46 : #include <limits>
47 : #include <memory>
48 :
49 : #ifdef _OPENMP
50 : #include <omp.h>
51 : #endif
52 :
53 : #include <AMReX_ParallelDescriptor.H>
54 : #include <AMReX_ParmParse.H>
55 : #include <AMReX_MultiFabUtil.H>
56 : #include <AMReX_FillPatchUtil.H>
57 : #include <AMReX_BC_TYPES.H>
58 : #include <AMReX_AmrCore.H>
59 : #include <AMReX_FluxRegister.H>
60 : #include <AMReX_Utility.H>
61 : #include <AMReX_PlotFileUtil.H>
62 :
63 : #include "Set/Set.H"
64 : #include "BC/BC.H"
65 : #include "BC/Nothing.H"
66 : #include "IO/WriteMetaData.H"
67 : #include "BaseField.H"
68 :
69 : /// \brief Collection of numerical integrator objects
70 : namespace Integrator
71 : {
72 :
73 : class Integrator
74 : : public amrex::AmrCore
75 : {
76 : public:
77 :
78 : /// \function Integrator
79 : /// \brief Constructor
80 : ///
81 : /// Does the following things:
82 : /// - Read in simulation TIME(STEP) information
83 : /// - Read in simulation output and AMR information
84 : /// - Initalize timestep substep information
85 : /// - Create a clean directory
86 : /// For derived classes this **must** be called for the derived constructor. For instance: `code`
87 : /// ```cpp
88 : /// class MyDerivedClass : Integrator
89 : /// {
90 : /// MyDerivedClass() : Integrator() { ... }
91 : /// ...
92 : /// }
93 : /// ```
94 : Integrator();
95 :
96 : virtual ~Integrator();
97 :
98 : /// \fn FrontData
99 : /// \brief Front-end method to initialize simulation
100 : void InitData();
101 :
102 : void Restart(std::string restartfile, bool a_node = false);
103 :
104 : /// \fn Evolve
105 : /// \brief Front-end method to start simulation
106 : void Evolve();
107 :
108 : void SetFilename(std::string _plot_file) { plot_file = _plot_file; };
109 : std::string GetFilename() { return plot_file; };
110 :
111 2196 : void regrid(int lbase, Set::Scalar time, bool initial = false) override
112 : {
113 2196 : if (!explicitmesh.on)
114 2188 : AmrCore::regrid(lbase, time, initial);
115 2196 : }
116 :
117 31 : void InitFromScratch(Set::Scalar time)
118 : {
119 31 : if (!explicitmesh.on) AmrCore::InitFromScratch(time);
120 : else
121 : {
122 : // Generate the coarse level mesh.
123 : {
124 7 : finest_level = 0;
125 14 : const amrex::BoxArray& ba = MakeBaseGrids();
126 14 : amrex::DistributionMapping dm(ba);
127 7 : const auto old_num_setdm = num_setdm;
128 7 : const auto old_num_setba = num_setba;
129 7 : MakeNewLevelFromScratch(0, time, ba, dm);
130 7 : if (old_num_setba == num_setba) SetBoxArray(0, ba);
131 7 : if (old_num_setdm == num_setdm) SetDistributionMap(0, dm);
132 : }
133 : // Generate subsequent level meshes based on user input
134 13 : for (int ilev = 0; ilev < maxLevel(); ++ilev)
135 : {
136 6 : finest_level = ilev + 1;
137 12 : amrex::BoxArray grids(explicitmesh.box[ilev]);
138 6 : ChopGrids(ilev + 1, grids, amrex::ParallelDescriptor::NProcs());
139 12 : amrex::DistributionMapping dmap(grids);
140 6 : SetBoxArray(ilev + 1, grids);
141 6 : SetDistributionMap(ilev + 1, dmap);
142 6 : MakeNewLevelFromScratch(ilev + 1, time, grids, dmap);
143 : }
144 :
145 : }
146 31 : }
147 :
148 : protected:
149 :
150 : /// \fn Initialize
151 : /// \brief Apply initial conditions
152 : ///
153 : /// You **must** override this function to inherit this class.
154 : /// This function is called before the simulation begins, and is where
155 : /// initial conditions should be applied.
156 : virtual void Initialize(int lev ///<[in] AMR Level
157 : ) = 0;
158 :
159 : /// \fn Advance
160 : /// \brief Perform computation
161 : ///
162 : /// You **must** override this function to inherit this class.
163 : /// Advance is called every time(sub)step, and implements the evolution of
164 : /// the system in time.
165 : ///
166 : virtual void Advance(int lev, ///<[in] AMR Level
167 : amrex::Real time, ///< [in] System time
168 : amrex::Real dt ///< [in] Timestep for this level
169 : ) = 0;
170 :
171 : /// \fn TagCellsForRefinement
172 : /// \brief Tag cells where mesh refinement is needed
173 : ///
174 : /// You **must** override this function to inherit this class.
175 : /// Advance is called every time(sub)step, and implements the evolution of
176 : /// the system in time.
177 : ///
178 : virtual void TagCellsForRefinement(int lev, amrex::TagBoxArray& tags, amrex::Real time,
179 : int ngrow) = 0;
180 :
181 : /// \fn TimeStepBegin
182 : /// \brief Run another system calculation (e.g. implicit solve) before integration step
183 : ///
184 : /// This function is called at the beginning of every timestep. This function can be used
185 : /// to complete additional global solves, e.g. a MLMG implicit solve.
186 : ///
187 : /// Overriding is optional; the default is to do nothing.
188 : ///
189 1606 : virtual void TimeStepBegin(Set::Scalar /*time*/, int /*iter*/) {};
190 :
191 : /// \fn TimeStepComplete
192 : /// \brief Run another system calculation (e.g. implicit solve) after integration step
193 : ///
194 : /// This function is called at the end of every timestep. This function can be used
195 : /// to complete additional global solves, e.g. a MLMG implicit solve.
196 : ///
197 : /// Overriding is optional; the default is to do nothing.
198 : ///
199 2032 : virtual void TimeStepComplete(Set::Scalar /*time*/, int /*iter*/) {};
200 :
201 : /// \fn Integrate
202 : /// \brief Perform an integration to compute integrated quantities
203 : ///
204 : /// This is a function that is called by `Integrator` to update the variables registered in
205 : /// RegisterIntegratedVariable.
206 : /// The following variables are used:
207 : /// - amrlev: current amr level
208 : /// - time: current simulation time
209 : /// - iter: current simulation iteration
210 : /// - mfi: current MFIter object (used to get FArrayBox from MultiFab)
211 : /// - box: Use this box (not mfi.tilebox). This box covers only cells on this level that are
212 : /// not also on a finer level.
213 0 : virtual void Integrate(int /*amrlev*/, Set::Scalar /*time*/, int /*iter*/,
214 : const amrex::MFIter&/*mfi*/, const amrex::Box&/*box*/)
215 : {
216 0 : if (thermo.number > 0)
217 0 : Util::Warning(INFO, "integrated variables registered, but no integration implemented!");
218 0 : }
219 :
220 409 : virtual void Regrid(int /* amrlev */, Set::Scalar /* time */)
221 409 : {}
222 :
223 :
224 :
225 :
226 : /// \fn RegisterNewFab
227 : /// \brief Register a field variable for AMR with this class
228 : void RegisterNewFab(Set::Field<Set::Scalar>& new_fab, BC::BC<Set::Scalar>* new_bc, int ncomp, int nghost, std::string name, bool writeout,std::vector<std::string> suffix = {});
229 : void RegisterNewFab(Set::Field<Set::Scalar>& new_fab, int ncomp, std::string name, bool writeout,std::vector<std::string> suffix = {});
230 : void RegisterNodalFab(Set::Field<Set::Scalar>& new_fab, int ncomp, int nghost, std::string name, bool writeout,std::vector<std::string> suffix = {});
231 : void RegisterNodalFab(Set::Field<Set::Scalar>& new_fab, BC::BC<Set::Scalar>* new_bc, int ncomp, int nghost, std::string name, bool writeout,std::vector<std::string> suffix = {});
232 : template<class T>
233 : void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool evolving = true);
234 : template<class T>
235 : void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, std::string a_name, bool evolving = true);
236 : template<class T>
237 : void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool writeout, std::string a_name, bool evolving = true);
238 :
239 : template<class T, int d>
240 : void AddField( Set::Field<T>& new_field, BC::BC<T>* new_bc, int ncomp, int nghost,
241 : std::string, bool writeout, bool evolving, std::vector<std::string> suffix = {});
242 :
243 70407 : void SetFinestLevel(const int a_finestlevel)
244 : {
245 232949 : for (unsigned int i = 0; i < cell.fab_array.size(); i++)
246 162542 : cell.fab_array[i]->finest_level = a_finestlevel;
247 125200 : for (unsigned int i = 0; i < node.fab_array.size(); i++)
248 54793 : node.fab_array[i]->finest_level = a_finestlevel;
249 70407 : for (unsigned int i = 0; i < m_basefields_cell.size(); i++)
250 0 : m_basefields_cell[i]->SetFinestLevel(finest_level);
251 74227 : for (unsigned int i = 0; i < m_basefields.size(); i++)
252 3820 : m_basefields[i]->SetFinestLevel(finest_level);
253 70407 : }
254 :
255 : void RegisterIntegratedVariable(Set::Scalar* integrated_variable, std::string name, bool extensive=true);
256 :
257 : void SetTimestep(Set::Scalar _timestep);
258 : void SetPlotInt(int plot_int);
259 1 : void SetThermoInt(int a_thermo_int) { thermo.interval = a_thermo_int; }
260 : void SetThermoPlotInt(int a_thermo_plot_int) { thermo.plot_int = a_thermo_plot_int; }
261 : void SetStopTime(Set::Scalar a_stop_time) { stop_time = a_stop_time; }
262 :
263 :
264 : // Dynamic timestep adjustment
265 :
266 : struct {
267 : // user params
268 : bool on = false;
269 : int verbose = -1;
270 : int nprevious = -1;
271 : Set::Scalar cfl = NAN;
272 : Set::Scalar min = NAN;
273 : Set::Scalar max = NAN;
274 : // internal variables
275 : std::vector<Set::Scalar> dt_limit_min;
276 : std::vector<Set::Scalar> previous_timesteps;
277 : } dynamictimestep;
278 : void DynamicTimestep_SyncTimeStep(int lev, Set::Scalar dt_min)
279 : {
280 : if (!dynamictimestep.on) return;
281 : if (!dynamictimestep.dt_limit_min.size())
282 : {
283 : dynamictimestep.dt_limit_min.resize(max_level+1,
284 : std::numeric_limits<Set::Scalar>::max());
285 : }
286 : dynamictimestep.dt_limit_min[lev] = std::min(dt_min,
287 : dynamictimestep.dt_limit_min[lev]);
288 :
289 : amrex::ParallelDescriptor::ReduceRealMin(dynamictimestep.dt_limit_min[lev]);
290 : }
291 : void DynamicTimestep_Reset()
292 : {
293 : if (!dynamictimestep.on) return;
294 : dynamictimestep.previous_timesteps.clear();
295 : }
296 : void DynamicTimestep_Update()
297 : {
298 : if (!dynamictimestep.on) return;
299 : if (!dynamictimestep.dt_limit_min.size()) return;
300 :
301 : Set::Scalar final_timestep = NAN;
302 :
303 : if (amrex::ParallelDescriptor::IOProcessor())
304 : {
305 : Set::Scalar timestep_average = this->dt[0];
306 : if (dynamictimestep.previous_timesteps.size() > 0)
307 : {
308 : timestep_average = 0.0;
309 : for (unsigned int d = 0; d < dynamictimestep.previous_timesteps.size(); d++)
310 : timestep_average += dynamictimestep.previous_timesteps[d];
311 : timestep_average /= dynamictimestep.previous_timesteps.size();
312 : }
313 :
314 : Set::Scalar new_timestep = std::numeric_limits<Set::Scalar>::max();
315 : for (int lev = 0; lev <= this->max_level; lev++)
316 : {
317 : //const Set::Scalar* DX = this->geom[lev].CellSize();
318 : Set::Scalar dt_lev = dynamictimestep.dt_limit_min[lev];
319 :
320 : Util::Message(INFO,"lev=",lev," ",dt_lev, " (",this->nsubsteps[lev],")");
321 :
322 : for (int ilev = lev; ilev > 0; ilev--) dt_lev *= (Set::Scalar)(this->nsubsteps[ilev]);
323 :
324 : Util::Message(INFO,"lev=",lev," --> ",dt_lev);
325 :
326 : new_timestep = std::min(new_timestep,dt_lev);
327 : }
328 :
329 : if (new_timestep < timestep_average)
330 : {
331 : dynamictimestep.previous_timesteps.clear();
332 :
333 : final_timestep = new_timestep;
334 : final_timestep = std::max(final_timestep,dynamictimestep.min);
335 : final_timestep = std::min(final_timestep,dynamictimestep.max);
336 :
337 : dynamictimestep.previous_timesteps.push_back(new_timestep);
338 : }
339 : else
340 : {
341 : final_timestep = timestep_average;
342 : final_timestep = std::max(final_timestep,dynamictimestep.min);
343 : final_timestep = std::min(final_timestep,dynamictimestep.max);
344 :
345 : if ((int)(dynamictimestep.previous_timesteps.size()) > dynamictimestep.nprevious)
346 : dynamictimestep.previous_timesteps.erase(dynamictimestep.previous_timesteps.begin()); // pop first
347 : dynamictimestep.previous_timesteps.push_back(new_timestep); // push back new timestep
348 : }
349 :
350 : }
351 : amrex::ParallelDescriptor::Bcast(&final_timestep,1);
352 : this->SetTimestep(final_timestep);
353 : dynamictimestep.dt_limit_min.clear();
354 : }
355 :
356 :
357 :
358 : amrex::Vector<amrex::Real> t_new; ///< Keep track of current old simulation time on each level
359 : amrex::Vector<int> istep; ///< Keep track of where each level is
360 : // PLOT FILES
361 : std::string plot_file{ "plt" }; ///< Plotfile name
362 :
363 : private:
364 : virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray& ba,
365 : const amrex::DistributionMapping& dm) override;
366 : virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray& ba,
367 : const amrex::DistributionMapping& dm) override;
368 : virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray& ba,
369 : const amrex::DistributionMapping& dm) override;
370 : virtual void ClearLevel(int lev) override;
371 : virtual void ErrorEst(int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
372 :
373 :
374 : void FillPatch(int lev, amrex::Real time,
375 : amrex::Vector<std::unique_ptr<amrex::MultiFab>>& source_mf,
376 : amrex::MultiFab& destination_multifab,
377 : BC::BC<Set::Scalar>& physbc,
378 : int icomp);
379 : long CountCells(int lev);
380 : void TimeStep(int lev, amrex::Real time, int iteration);
381 : void FillCoarsePatch(int lev, amrex::Real time, Set::Field<Set::Scalar>& mf, BC::BC<Set::Scalar>& physbc, int icomp, int ncomp);
382 : void GetData(const int lev, const amrex::Real time, amrex::Vector<amrex::MultiFab*>& data, amrex::Vector<amrex::Real>& datatime);
383 :
384 : std::vector<std::string> PlotFileName(int lev, std::string prefix = "") const;
385 : protected:
386 : void IntegrateVariables(Set::Scalar cur_time, int step);
387 : void WritePlotFile(bool initial = false) const;
388 : void WritePlotFile(std::string prefix, Set::Scalar time, int step) const;
389 : void WritePlotFile(Set::Scalar time, amrex::Vector<int> iter, bool initial = false, std::string prefix = "") const;
390 :
391 : //
392 : // MEMBER VARIABLES
393 : //
394 :
395 : // TIME (STEP) KEEPINGamrex::Vector<std::unique_ptr<amrex::MultiFab> >
396 : protected:
397 : amrex::Real timestep = NAN; ///< Timestep for the base level of refinement
398 : amrex::Vector<amrex::Real> dt; ///< Timesteps for each level of refinement
399 : amrex::Vector<int> nsubsteps; ///< how many substeps on each level?
400 : private:
401 : int max_plot_level = -1;
402 :
403 : amrex::Vector<amrex::Real> t_old;///< Keep track of current old simulation time on each level
404 : int max_step = std::numeric_limits<int>::max(); ///< Maximum allowable timestep
405 : amrex::Real tstart = 0; ///< Default start time (default: 0)
406 : amrex::Real stop_time = NAN; ///< Default stop time
407 :
408 : protected:
409 : bool integrate_variables_before_advance = true;
410 : bool integrate_variables_after_advance = false;
411 :
412 : protected:
413 : struct {
414 : int number_of_fabs = 0;
415 : std::vector<Set::Field<Set::Scalar>*> fab_array;
416 : std::vector<int> ncomp_array;
417 : std::vector<int> nghost_array;
418 : std::vector<std::vector<std::string>> name_array;
419 : std::vector<BC::BC<Set::Scalar>*> physbc_array;
420 : std::vector<bool> writeout_array;
421 : bool any = true;
422 : bool all = false;
423 : } node;
424 :
425 : struct {
426 : int number_of_fabs = 0;
427 : std::vector<Set::Field<Set::Scalar>*> fab_array;
428 : std::vector<int> ncomp_array;
429 : std::vector<int> nghost_array;
430 : std::vector<std::vector<std::string>> name_array;
431 : std::vector<BC::BC<Set::Scalar>*> physbc_array;
432 : std::vector<bool> writeout_array;
433 : bool any = true;
434 : bool all = false;
435 : } cell;
436 :
437 : std::vector<BaseField*> m_basefields;
438 : std::vector<BaseField*> m_basefields_cell;
439 :
440 : BC::Nothing bcnothing;
441 :
442 : // KEEP TRACK OF ALL INTEGRATED VARIABLES
443 : struct {
444 : int interval = -1;
445 : Set::Scalar dt = NAN;
446 : int plot_int = -1;
447 : Set::Scalar plot_dt = NAN;
448 : int number = 0;
449 : std::vector<Set::Scalar*> vars;
450 : std::vector<std::string> names;
451 : std::vector<bool> extensives;
452 : } thermo;
453 :
454 : // REGRIDDING
455 : int regrid_int = -1; ///< Determine how often to regrid (default: 2)
456 : int base_regrid_int = -1; ///< Determine how often to regrid based on coarse level only (default: 0)
457 :
458 : std::string restart_file_cell = "";
459 : std::string restart_file_node = "";
460 :
461 : struct {
462 : int on = 0;
463 : std::vector<amrex::Box> box;
464 : } explicitmesh;
465 :
466 : protected:
467 : int plot_int = -1; ///< How frequently to dump plot file (default: never)
468 : Set::Scalar plot_dt = -1.0;
469 :
470 : int abort_on_nan = true;
471 : };
472 :
473 :
474 : template<>
475 : ALAMO_SINGLE_DEFINITION
476 40 : void Integrator::AddField<Set::Scalar, Set::Hypercube::Cell>
477 : ( Set::Field<Set::Scalar>& new_field,
478 : BC::BC<Set::Scalar>* new_bc,
479 : int ncomp,
480 : int nghost,
481 : std::string name,
482 : bool writeout,
483 : bool /*evolving*/,
484 : std::vector<std::string> suffix)
485 : {
486 40 : int nlevs_max = maxLevel() + 1;
487 40 : new_field.resize(nlevs_max);
488 40 : cell.fab_array.push_back(&new_field);
489 40 : if (new_bc != nullptr) cell.physbc_array.push_back(new_bc);
490 0 : else cell.physbc_array.push_back(&bcnothing);
491 40 : cell.ncomp_array.push_back(ncomp);
492 40 : cell.nghost_array.push_back(nghost);
493 : //cell.name_array.push_back(name);
494 40 : cell.writeout_array.push_back(writeout);
495 40 : cell.number_of_fabs++;
496 :
497 80 : Util::Assert(INFO,TEST((int)suffix.size() == 0 || (int)suffix.size() == ncomp));
498 80 : std::vector<std::string> names;
499 40 : if (ncomp == 1)
500 34 : names.push_back(name);
501 : else
502 : {
503 6 : if (suffix.size() == 0)
504 34 : for (int j = 0; j < ncomp; j++)
505 28 : names.push_back(amrex::Concatenate(name,j+1,3));
506 : else
507 0 : for (int j = 0; j < ncomp; j++)
508 0 : names.push_back(name + suffix[j]);
509 : }
510 40 : cell.name_array.push_back(names);
511 40 : }
512 :
513 : template<>
514 : ALAMO_SINGLE_DEFINITION
515 21 : void Integrator::AddField<Set::Scalar, Set::Hypercube::Node>
516 : ( Set::Field<Set::Scalar>& new_field,
517 : BC::BC<Set::Scalar>* new_bc,
518 : int ncomp,
519 : int nghost,
520 : std::string name,
521 : bool writeout,
522 : bool /*evolving*/,
523 : std::vector<std::string> suffix)
524 : {
525 : BL_PROFILE("Integrator::RegisterNodalFab");
526 42 : Util::Assert(INFO, TEST(new_bc == nullptr));
527 21 : int nlevs_max = maxLevel() + 1;
528 21 : new_field.resize(nlevs_max);
529 21 : node.fab_array.push_back(&new_field);
530 21 : node.physbc_array.push_back(&bcnothing);
531 21 : node.ncomp_array.push_back(ncomp);
532 21 : node.nghost_array.push_back(nghost);
533 : //node.name_array.push_back(name);
534 21 : node.writeout_array.push_back(writeout);
535 21 : node.number_of_fabs++;
536 :
537 42 : Util::Assert(INFO,TEST((int)suffix.size() == 0 || (int)suffix.size() == ncomp));
538 42 : std::vector<std::string> names;
539 21 : if (ncomp == 1)
540 18 : names.push_back(name);
541 : else
542 : {
543 3 : if (suffix.size() == 0)
544 9 : for (int j = 0; j < ncomp; j++)
545 6 : names.push_back(amrex::Concatenate(name,j+1,3));
546 : else
547 0 : for (int j = 0; j < ncomp; j++)
548 0 : names.push_back(name + suffix[j]);
549 : }
550 21 : node.name_array.push_back(names);
551 21 : }
552 :
553 : template<class T, int d>
554 : ALAMO_SINGLE_DEFINITION
555 98 : void Integrator::AddField
556 : ( Set::Field<T>& new_field, BC::BC<T>* new_bc, int ncomp,
557 : int nghost, std::string name, bool writeout, bool evolving,
558 : std::vector<std::string> /*suffix*/)
559 : {
560 : if (d == Set::Hypercube::Node)
561 : {
562 196 : Util::Assert(INFO, TEST(new_bc == nullptr));
563 98 : int nlevs_max = maxLevel() + 1;
564 98 : new_field.resize(nlevs_max);
565 98 : m_basefields.push_back(new Field<T>(new_field, geom, refRatio(), ncomp, nghost));
566 98 : m_basefields.back()->evolving = evolving;
567 98 : m_basefields.back()->writeout = writeout;
568 98 : m_basefields.back()->setName(name);
569 98 : m_basefields.back()->evolving = evolving;
570 98 : m_basefields.back()->m_gridtype = Set::Hypercube::Node;
571 : }
572 : else if (d == Set::Hypercube::Cell)
573 : {
574 : int nlevs_max = maxLevel() + 1;
575 : new_field.resize(nlevs_max);
576 : m_basefields_cell.push_back(new Field<T>(new_field, geom, refRatio(), ncomp, nghost));
577 : m_basefields_cell.back()->evolving = evolving;
578 : m_basefields_cell.back()->writeout = writeout;
579 : m_basefields_cell.back()->setName(name);
580 : m_basefields_cell.back()->evolving = evolving;
581 : if (new_bc) m_basefields_cell.back()->setBC(new_bc);
582 : else m_basefields_cell.back()->setBC(&bcnothing);
583 : m_basefields_cell.back()->m_gridtype = Set::Hypercube::Cell;
584 : }
585 : else
586 : {
587 : Util::Abort(INFO, "Only node and cell based fields can be added at this time");
588 : }
589 98 : }
590 :
591 :
592 :
593 : template<class T>
594 : ALAMO_SINGLE_DEFINITION
595 19 : void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool evolving)
596 : {
597 : //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
598 19 : AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, "", true, evolving);
599 19 : }
600 : template<class T>
601 : ALAMO_SINGLE_DEFINITION
602 3 : void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, std::string a_name, bool evolving)
603 : {
604 : //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
605 3 : AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, a_name, true, evolving);
606 3 : }
607 : template<class T>
608 : AMREX_ATTRIBUTE_WEAK
609 76 : void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool writeout, std::string a_name, bool evolving)
610 : {
611 : //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
612 76 : AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, a_name, writeout, evolving);
613 76 : }
614 :
615 : }
616 : #endif
|