LCOV - code coverage report
Current view: top level - src/Integrator - Integrator.H (source / functions) Coverage Total Hit
Test: coverage_merged.info Lines: 70.5 % 139 98
Test Date: 2025-02-27 04:17:48 Functions: 63.2 % 76 48

            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              :     /// This is the constructor for the intetgrator class, which reads timestep information,
      79              :     /// simulation output and AMR, initialized time substep, and creates a new directory.
      80              :     Integrator();
      81              : 
      82              :     /// Virtual destructure; make sure delete any pointers that you create here.
      83              :     virtual ~Integrator();
      84              : 
      85              :     /// Front-end method to initialize simulation on all levels
      86              :     void InitData();
      87              : 
      88              :     /// Read in output from previous simulation and start simulation at that point -
      89              :     /// Not currently tested
      90              :     void Restart(std::string restartfile, bool a_node = false);
      91              : 
      92              :     /// Front-end method to start simulation
      93              :     void Evolve();
      94              : 
      95              :     /// Simple setter to set filename
      96              :     void SetFilename(std::string _plot_file) { plot_file = _plot_file; };
      97              : 
      98              :     /// Simple getter to get filename
      99              :     std::string GetFilename() { return plot_file; };
     100              : 
     101              :     /// This overrides an AMReX method just to allow for explicit meshing when
     102              :     /// desired.
     103         7005 :     void regrid(int lbase, Set::Scalar time, bool initial = false) override
     104              :     {
     105         7005 :         if (!explicitmesh.on)
     106         6993 :             AmrCore::regrid(lbase, time, initial);
     107         7004 :     }
     108              : 
     109              :     /// This creates a new levels that have not previously been
     110              :     /// used.
     111           67 :     void InitFromScratch(Set::Scalar time)
     112              :     {
     113           67 :         if (!explicitmesh.on) AmrCore::InitFromScratch(time);
     114              :         else
     115              :         {
     116              :             // Generate the coarse level mesh.
     117              :             {
     118            9 :                 finest_level = 0;
     119            9 :                 const amrex::BoxArray& ba = MakeBaseGrids();
     120            9 :                 amrex::DistributionMapping dm(ba);
     121            9 :                 const auto old_num_setdm = num_setdm;
     122            9 :                 const auto old_num_setba = num_setba;
     123            9 :                 MakeNewLevelFromScratch(0, time, ba, dm);
     124            9 :                 if (old_num_setba == num_setba) SetBoxArray(0, ba);
     125            9 :                 if (old_num_setdm == num_setdm) SetDistributionMap(0, dm);
     126            9 :             }
     127              :             // Generate subsequent level meshes based on user input
     128           19 :             for (int ilev = 0; ilev < maxLevel(); ++ilev)
     129              :             {
     130           10 :                 finest_level = ilev + 1;
     131           10 :                 amrex::BoxArray grids(explicitmesh.box[ilev]);
     132           10 :                 ChopGrids(ilev + 1, grids, amrex::ParallelDescriptor::NProcs());
     133           10 :                 amrex::DistributionMapping dmap(grids);
     134           10 :                 SetBoxArray(ilev + 1, grids);
     135           10 :                 SetDistributionMap(ilev + 1, dmap);
     136           10 :                 MakeNewLevelFromScratch(ilev + 1, time, grids, dmap);
     137           10 :             }
     138              : 
     139              :         }
     140           66 :     }
     141              : 
     142              : protected:
     143              : 
     144              :     /// You must override this function to inherit this class;
     145              :     /// this function is called before the simulation begins, and is where
     146              :     /// initial conditions should be applied.
     147              :     virtual void Initialize(int lev ///<[in] AMR Level
     148              :     ) = 0;
     149              : 
     150              :     /// You must override this function to inherit this class;
     151              :     /// Advance is called every time(sub)step, and implements the evolution of 
     152              :     /// the system in time.
     153              :     virtual void Advance(int lev,          ///<[in] AMR Level
     154              :         amrex::Real time, ///< [in] System time
     155              :         amrex::Real dt    ///< [in] Timestep for this level
     156              :     ) = 0;
     157              : 
     158              :     /// You must override this function to inherit this class;
     159              :     /// Advance is called every time(sub)step, and implements the evolution of 
     160              :     /// the system in time.
     161              :     virtual void TagCellsForRefinement(int lev, amrex::TagBoxArray& tags, amrex::Real time,
     162              :         int ngrow) = 0;
     163              : 
     164              :     /// This optional function is  called at the beginning of every timestep, and can be used
     165              :     /// to complete additional global solves, e.g. a MLMG implicit solve.
     166         9740 :     virtual void TimeStepBegin(Set::Scalar /*time*/, int /*iter*/) {};
     167              : 
     168              :     /// This optional function is called at the end of every timestep, and can be used
     169              :     /// to complete additional global solves, e.g. a MLMG implicit solve.
     170        20162 :     virtual void TimeStepComplete(Set::Scalar /*time*/, int /*iter*/) {};
     171              : 
     172              :     /// This is a function that is called by `Integrator` to update the variables registered in
     173              :     /// RegisterIntegratedVariable; you can override this to do your integration.
     174            0 :     virtual void Integrate( int /*amrlev*/,
     175              :                             Set::Scalar /*time*/,
     176              :                             int /*iter*/,
     177              :                             const amrex::MFIter&/*mfi*/,
     178              :                             const amrex::Box&/*box*/
     179              :         )
     180              :     {
     181            0 :         if (thermo.number > 0)
     182            0 :             Util::Warning(INFO, "integrated variables registered, but no integration implemented!");
     183            0 :     }
     184              : 
     185              :     /// An optionally overridable method to trigger behavior whenver a regrid occurs.
     186         1070 :     virtual void Regrid(int /* amrlev */, Set::Scalar /* time */)
     187         1070 :     {}
     188              : 
     189              : 
     190              :     /// Add a new cell-based scalar field.
     191              :     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 = {});
     192              :     /// Add a new cell-based scalar field (with additional arguments).
     193              :     void RegisterNewFab(Set::Field<Set::Scalar>& new_fab, int ncomp, std::string name, bool writeout,std::vector<std::string> suffix = {});
     194              :     /// Add a new node-based scalar field
     195              :     void RegisterNodalFab(Set::Field<Set::Scalar>& new_fab, int ncomp, int nghost, std::string name, bool writeout,std::vector<std::string> suffix = {});
     196              :     /// Add a new node-based scalar field (wtih additional arguments)
     197              :     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 = {});
     198              :     template<class T>
     199              :     /// Add a templated nodal field
     200              :     void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool evolving = true);
     201              :     template<class T>
     202              :     /// Add a templated nodal field (additional arguments)
     203              :     void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, std::string a_name, bool evolving = true);
     204              :     template<class T>
     205              :     /// Add a templated nodal field (additional arguments)
     206              :     void RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool writeout, std::string a_name, bool evolving = true);
     207              : 
     208              :     /// Add a field with arbitrary type (templated with T) and grid location (templated with d).
     209              :     template<class T, int d>
     210              :     void AddField(  Set::Field<T>& new_field, BC::BC<T>* new_bc, int ncomp, int nghost,
     211              :                     std::string, bool writeout, bool evolving, std::vector<std::string> suffix = {});
     212              : 
     213              :     /// Utility to ensure that all fields know what the finest level is
     214       142817 :     void SetFinestLevel(const int a_finestlevel)
     215              :     {
     216       622095 :         for (unsigned int i = 0; i < cell.fab_array.size(); i++)
     217       479278 :             cell.fab_array[i]->finest_level = a_finestlevel;
     218       229248 :         for (unsigned int i = 0; i < node.fab_array.size(); i++)
     219        86431 :             node.fab_array[i]->finest_level = a_finestlevel;
     220       142817 :         for (unsigned int i = 0; i < m_basefields_cell.size(); i++)
     221            0 :             m_basefields_cell[i]->SetFinestLevel(finest_level);
     222       268398 :         for (unsigned int i = 0; i < m_basefields.size(); i++)
     223       125581 :             m_basefields[i]->SetFinestLevel(finest_level);
     224       142817 :     }
     225              : 
     226              :     /// Register a variable to be integrated over the spatial domain using
     227              :     /// the Integrate function.
     228              :     void RegisterIntegratedVariable(Set::Scalar* integrated_variable, std::string name, bool extensive=true);
     229              : 
     230              :     /// Utility to set the coarse-grid timestep
     231              :     void SetTimestep(Set::Scalar _timestep);
     232              :     
     233              :     /// Utility to set the frequency (in timesteps) of plotfile dumping
     234              :     void SetPlotInt(int plot_int);
     235              :     
     236              :     /// Utility to set the frequency (in timesteps) of thermo data calculation
     237            2 :     void SetThermoInt(int a_thermo_int) { thermo.interval = a_thermo_int; }
     238              :     /// Utility to set the frequency (in timesteps) of thermo data writing to file
     239              :     void SetThermoPlotInt(int a_thermo_plot_int) { thermo.plot_int = a_thermo_plot_int; }
     240              :     /// Utility to set the global stop time.
     241              :     void SetStopTime(Set::Scalar a_stop_time) { stop_time = a_stop_time; }
     242              : 
     243              : 
     244              :     // Dynamic timestep adjustment
     245              : 
     246              :     struct {
     247              :         // user params
     248              :         bool on = false;
     249              :         int verbose = -1;
     250              :         int nprevious = -1;
     251              :         Set::Scalar cfl = NAN;
     252              :         Set::Scalar min = NAN;
     253              :         Set::Scalar max = NAN;
     254              :         // internal variables
     255              :         std::vector<Set::Scalar> dt_limit_min;
     256              :         std::vector<Set::Scalar> previous_timesteps;
     257              :     } dynamictimestep; /// Params for the dynamic timestp
     258              :     void DynamicTimestep_SyncTimeStep(int lev, Set::Scalar dt_min)
     259              :     {
     260              :         if (!dynamictimestep.on) return;
     261              :         if (!dynamictimestep.dt_limit_min.size())
     262              :         {
     263              :             dynamictimestep.dt_limit_min.resize(max_level+1,
     264              :                                                 std::numeric_limits<Set::Scalar>::max());
     265              :         }
     266              :         dynamictimestep.dt_limit_min[lev] = std::min(dt_min,
     267              :                                                     dynamictimestep.dt_limit_min[lev]);
     268              :         
     269              :         amrex::ParallelDescriptor::ReduceRealMin(dynamictimestep.dt_limit_min[lev]);
     270              :     }
     271              :     void DynamicTimestep_Reset()
     272              :     {
     273              :         if (!dynamictimestep.on) return;
     274              :         dynamictimestep.previous_timesteps.clear();
     275              :     }
     276          657 :     void DynamicTimestep_Update()
     277              :     {
     278          657 :         if (!dynamictimestep.on) return;
     279            0 :         if (!dynamictimestep.dt_limit_min.size()) return;
     280              : 
     281            0 :         Set::Scalar final_timestep = NAN;
     282              : 
     283            0 :         if (amrex::ParallelDescriptor::IOProcessor())
     284              :         {
     285            0 :             Set::Scalar timestep_average = this->dt[0];
     286            0 :             if (dynamictimestep.previous_timesteps.size() > 0)
     287              :             {
     288            0 :                 timestep_average = 0.0;
     289            0 :                 for (unsigned int d = 0; d < dynamictimestep.previous_timesteps.size(); d++)
     290            0 :                     timestep_average += dynamictimestep.previous_timesteps[d];
     291            0 :                 timestep_average /= dynamictimestep.previous_timesteps.size();
     292              :             }
     293              : 
     294            0 :             Set::Scalar new_timestep = std::numeric_limits<Set::Scalar>::max();
     295            0 :             for (int lev = 0; lev <= this->max_level; lev++)
     296              :             {
     297              :                 //const Set::Scalar* DX = this->geom[lev].CellSize();
     298            0 :                 Set::Scalar dt_lev = dynamictimestep.dt_limit_min[lev];
     299              : 
     300            0 :                 Util::Message(INFO,"lev=",lev," ",dt_lev, " (",this->nsubsteps[lev],")");
     301              : 
     302            0 :                 for (int ilev = lev; ilev > 0; ilev--) dt_lev *= (Set::Scalar)(this->nsubsteps[ilev]);
     303              : 
     304            0 :                 Util::Message(INFO,"lev=",lev,"       -->    ",dt_lev);
     305              : 
     306            0 :                 new_timestep = std::min(new_timestep,dt_lev);
     307              :             }
     308              : 
     309            0 :             if (new_timestep < timestep_average)
     310              :             {
     311            0 :                 dynamictimestep.previous_timesteps.clear();
     312              : 
     313            0 :                 final_timestep = new_timestep;
     314            0 :                 final_timestep = std::max(final_timestep,dynamictimestep.min);
     315            0 :                 final_timestep = std::min(final_timestep,dynamictimestep.max);
     316              : 
     317            0 :                 dynamictimestep.previous_timesteps.push_back(new_timestep);
     318              :             }
     319              :             else
     320              :             {
     321            0 :                 final_timestep = timestep_average;
     322            0 :                 final_timestep = std::max(final_timestep,dynamictimestep.min);
     323            0 :                 final_timestep = std::min(final_timestep,dynamictimestep.max);
     324              : 
     325            0 :                 if ((int)(dynamictimestep.previous_timesteps.size()) > dynamictimestep.nprevious)
     326            0 :                     dynamictimestep.previous_timesteps.erase(dynamictimestep.previous_timesteps.begin()); // pop first
     327            0 :                 dynamictimestep.previous_timesteps.push_back(new_timestep); // push back new timestep
     328              :             }
     329              : 
     330              :         }
     331            0 :         amrex::ParallelDescriptor::Bcast(&final_timestep,1);
     332            0 :         this->SetTimestep(final_timestep);
     333            0 :         dynamictimestep.dt_limit_min.clear();
     334              :     }
     335              : 
     336              : 
     337              : 
     338              :     amrex::Vector<amrex::Real> t_new;         ///< Keep track of current old simulation time on each level
     339              :     amrex::Vector<int> istep;           ///< Keep track of where each level is
     340              :     // PLOT FILES
     341              :     std::string plot_file{ "plt" };   ///< Plotfile name
     342              : 
     343              : private:
     344              :     virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray& ba,
     345              :         const amrex::DistributionMapping& dm) override;
     346              :     virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray& ba,
     347              :         const amrex::DistributionMapping& dm) override;
     348              :     virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray& ba,
     349              :         const amrex::DistributionMapping& dm) override;
     350              :     virtual void ClearLevel(int lev) override;
     351              :     virtual void ErrorEst(int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
     352              : 
     353              : 
     354              :     /// This is the function that is responsible for updating patch data.
     355              :     void FillPatch(int lev, amrex::Real time,
     356              :         amrex::Vector<std::unique_ptr<amrex::MultiFab>>& source_mf,
     357              :         amrex::MultiFab& destination_multifab,
     358              :         BC::BC<Set::Scalar>& physbc,
     359              :         int icomp);
     360              :     /// Simple utility to count cells 
     361              :     long CountCells(int lev);
     362              :     /// Timestep marching
     363              :     void TimeStep(int lev, amrex::Real time, int iteration);
     364              :     void FillCoarsePatch(int lev, amrex::Real time, Set::Field<Set::Scalar>& mf, BC::BC<Set::Scalar>& physbc, int icomp, int ncomp);
     365              :     void GetData(const int lev, const amrex::Real time, amrex::Vector<amrex::MultiFab*>& data, amrex::Vector<amrex::Real>& datatime);
     366              : 
     367              :     std::vector<std::string> PlotFileName(int lev, std::string prefix = "") const;
     368              : protected:
     369              :     void IntegrateVariables(Set::Scalar cur_time, int step);
     370              :     void WritePlotFile(bool initial = false) const;
     371              :     void WritePlotFile(std::string prefix, Set::Scalar time, int step) const;
     372              :     void WritePlotFile(Set::Scalar time, amrex::Vector<int> iter, bool initial = false, std::string prefix = "") const;
     373              : 
     374              :     //
     375              :     // MEMBER VARIABLES
     376              :     //
     377              : 
     378              :     // TIME (STEP) KEEPINGamrex::Vector<std::unique_ptr<amrex::MultiFab> >
     379              : protected:
     380              :     amrex::Real timestep = NAN;   ///< Timestep for the base level of refinement
     381              :     amrex::Vector<amrex::Real> dt;  ///< Timesteps for each level of refinement
     382              :     amrex::Vector<int> nsubsteps;   ///< how many substeps on each level?
     383              : private:
     384              :     int max_plot_level = -1;
     385              : 
     386              :     amrex::Vector<amrex::Real> t_old;///< Keep track of current old simulation time on each level
     387              :     int max_step = std::numeric_limits<int>::max(); ///< Maximum allowable timestep
     388              :     amrex::Real tstart = 0; ///< Default start time (default: 0)
     389              :     amrex::Real stop_time = NAN; ///< Default stop time 
     390              : 
     391              : protected:
     392              :     bool integrate_variables_before_advance = true;
     393              :     bool integrate_variables_after_advance = false;
     394              : 
     395              : protected:
     396              :     struct {
     397              :         int number_of_fabs = 0;
     398              :         std::vector<Set::Field<Set::Scalar>*> fab_array;
     399              :         std::vector<int> ncomp_array;
     400              :         std::vector<int> nghost_array;
     401              :         std::vector<std::vector<std::string>> name_array;
     402              :         std::vector<BC::BC<Set::Scalar>*> physbc_array;
     403              :         std::vector<bool> writeout_array;
     404              :         bool any = true;
     405              :         bool all = false;
     406              :     } node;
     407              : 
     408              :     struct {
     409              :         int number_of_fabs = 0;
     410              :         std::vector<Set::Field<Set::Scalar>*> fab_array;
     411              :         std::vector<int> ncomp_array;
     412              :         std::vector<int> nghost_array;
     413              :         std::vector<std::vector<std::string>> name_array;
     414              :         std::vector<BC::BC<Set::Scalar>*> physbc_array;
     415              :         std::vector<bool> writeout_array;
     416              :         bool any = true;
     417              :         bool all = false;
     418              :     } cell;
     419              : 
     420              :     std::vector<BaseField*> m_basefields;
     421              :     std::vector<BaseField*> m_basefields_cell;
     422              : 
     423              :     BC::Nothing bcnothing;
     424              : 
     425              :     // KEEP TRACK OF ALL INTEGRATED VARIABLES
     426              :     struct {
     427              :         int interval = -1;
     428              :         Set::Scalar dt = NAN;
     429              :         int plot_int = -1;
     430              :         Set::Scalar plot_dt = NAN;
     431              :         int number = 0;
     432              :         std::vector<Set::Scalar*> vars;
     433              :         std::vector<std::string> names;
     434              :         std::vector<bool> extensives;
     435              :     } thermo;
     436              : 
     437              :     // REGRIDDING
     438              :     int regrid_int = -1;       ///< Determine how often to regrid (default: 2)
     439              :     int base_regrid_int = -1; ///< Determine how often to regrid based on coarse level only (default: 0)
     440              : 
     441              :     std::string restart_file_cell = "";
     442              :     std::string restart_file_node = "";
     443              : 
     444              :     struct {
     445              :         int on = 0;
     446              :         std::vector<amrex::Box> box;
     447              :     } explicitmesh;
     448              : 
     449              : protected:
     450              :     int plot_int = -1;               ///< How frequently to dump plot file (default: never)
     451              :     Set::Scalar plot_dt = -1.0;
     452              : 
     453              :     int abort_on_nan = true;
     454              : };
     455              : 
     456              : 
     457              : template<>
     458              : ALAMO_SINGLE_DEFINITION
     459           98 : void Integrator::AddField<Set::Scalar, Set::Hypercube::Cell>
     460              : (   Set::Field<Set::Scalar>& new_field,
     461              :     BC::BC<Set::Scalar>* new_bc,
     462              :     int ncomp,
     463              :     int nghost,
     464              :     std::string name,
     465              :     bool writeout,
     466              :     bool /*evolving*/,
     467              :     std::vector<std::string> suffix)
     468              : {
     469           98 :     int nlevs_max = maxLevel() + 1;
     470           98 :     new_field.resize(nlevs_max);
     471           98 :     cell.fab_array.push_back(&new_field);
     472           98 :     if (new_bc != nullptr) cell.physbc_array.push_back(new_bc);
     473            0 :     else                   cell.physbc_array.push_back(&bcnothing);
     474           98 :     cell.ncomp_array.push_back(ncomp);
     475           98 :     cell.nghost_array.push_back(nghost);
     476              :     //cell.name_array.push_back(name);
     477           98 :     cell.writeout_array.push_back(writeout);
     478           98 :     cell.number_of_fabs++;
     479              : 
     480          686 :     Util::Assert(INFO,TEST((int)suffix.size() == 0 || (int)suffix.size() == ncomp));
     481           98 :     std::vector<std::string> names;
     482           98 :     if (ncomp == 1)
     483           78 :         names.push_back(name);
     484              :     else
     485              :     {
     486           20 :         if (suffix.size() == 0)
     487          124 :             for (int j = 0; j < ncomp; j++)
     488          104 :                 names.push_back(amrex::Concatenate(name,j+1,3));
     489              :         else
     490            0 :             for (int j = 0; j < ncomp; j++)
     491            0 :                 names.push_back(name + suffix[j]);
     492              :     }
     493           98 :     cell.name_array.push_back(names);
     494           98 : }
     495              : 
     496              : template<>
     497              : ALAMO_SINGLE_DEFINITION
     498           47 : void Integrator::AddField<Set::Scalar, Set::Hypercube::Node>
     499              : (   Set::Field<Set::Scalar>& new_field,
     500              :     BC::BC<Set::Scalar>* new_bc,
     501              :     int ncomp,
     502              :     int nghost,
     503              :     std::string name,
     504              :     bool writeout,
     505              :     bool /*evolving*/,
     506              :     std::vector<std::string> suffix)
     507              : {
     508              :     BL_PROFILE("Integrator::RegisterNodalFab");
     509          329 :     Util::Assert(INFO, TEST(new_bc == nullptr));
     510           47 :     int nlevs_max = maxLevel() + 1;
     511           47 :     new_field.resize(nlevs_max);
     512           47 :     node.fab_array.push_back(&new_field);
     513           47 :     node.physbc_array.push_back(&bcnothing);
     514           47 :     node.ncomp_array.push_back(ncomp);
     515           47 :     node.nghost_array.push_back(nghost);
     516              :     //node.name_array.push_back(name);
     517           47 :     node.writeout_array.push_back(writeout);
     518           47 :     node.number_of_fabs++;
     519              : 
     520          329 :     Util::Assert(INFO,TEST((int)suffix.size() == 0 || (int)suffix.size() == ncomp));
     521           47 :     std::vector<std::string> names;
     522           47 :     if (ncomp == 1)
     523           35 :         names.push_back(name);
     524              :     else
     525              :     {
     526           12 :         if (suffix.size() == 0)
     527           44 :             for (int j = 0; j < ncomp; j++)
     528           32 :                 names.push_back(amrex::Concatenate(name,j+1,3));
     529              :         else
     530            0 :             for (int j = 0; j < ncomp; j++)
     531            0 :                 names.push_back(name + suffix[j]);
     532              :     }
     533           47 :     node.name_array.push_back(names);
     534           47 : }
     535              : 
     536              : template<class T, int d>
     537              : ALAMO_SINGLE_DEFINITION
     538          272 : void Integrator::AddField
     539              : (   Set::Field<T>& new_field, BC::BC<T>* new_bc, int ncomp,
     540              :     int nghost, std::string name, bool writeout, bool evolving,
     541              :     std::vector<std::string> /*suffix*/)
     542              : {
     543              :     if (d == Set::Hypercube::Node)
     544              :     {
     545         1904 :         Util::Assert(INFO, TEST(new_bc == nullptr));
     546          272 :         int nlevs_max = maxLevel() + 1;
     547          272 :         new_field.resize(nlevs_max);
     548          272 :         m_basefields.push_back(new Field<T>(new_field, geom, refRatio(), ncomp, nghost));
     549          272 :         m_basefields.back()->evolving = evolving;
     550          272 :         m_basefields.back()->writeout = writeout;
     551          272 :         m_basefields.back()->setName(name);
     552          272 :         m_basefields.back()->evolving = evolving;
     553          272 :         m_basefields.back()->m_gridtype = Set::Hypercube::Node;
     554              :     }
     555              :     else if (d == Set::Hypercube::Cell)
     556              :     {
     557              :         int nlevs_max = maxLevel() + 1;
     558              :         new_field.resize(nlevs_max);
     559              :         m_basefields_cell.push_back(new Field<T>(new_field, geom, refRatio(), ncomp, nghost));
     560              :         m_basefields_cell.back()->evolving = evolving;
     561              :         m_basefields_cell.back()->writeout = writeout;
     562              :         m_basefields_cell.back()->setName(name);
     563              :         m_basefields_cell.back()->evolving = evolving;
     564              :         if (new_bc) m_basefields_cell.back()->setBC(new_bc);
     565              :         else m_basefields_cell.back()->setBC(&bcnothing);
     566              :         m_basefields_cell.back()->m_gridtype = Set::Hypercube::Cell;
     567              :     }
     568              :     else
     569              :     {
     570              :         Util::Abort(INFO, "Only node and cell based fields can be added at this time");
     571              :     }
     572          272 : }
     573              : 
     574              : 
     575              : 
     576              : template<class T>
     577              : ALAMO_SINGLE_DEFINITION
     578           63 : void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool evolving) 
     579              : {
     580              :     //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
     581          189 :     AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, "", true, evolving);
     582           63 : }
     583              : template<class T>
     584              : ALAMO_SINGLE_DEFINITION
     585           17 : void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, std::string a_name, bool evolving)
     586              : {
     587              :     //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
     588           17 :     AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, a_name, true, evolving);
     589           17 : }
     590              : template<class T>
     591              : AMREX_ATTRIBUTE_WEAK
     592          192 : void Integrator::RegisterGeneralFab(Set::Field<T>& new_fab, int ncomp, int nghost, bool writeout, std::string a_name, bool evolving)
     593              : {
     594              :     //Util::Warning(INFO, "RegisterGeneralFab is depricated. Please replace with AddField");
     595          192 :     AddField<T, Set::Hypercube::Node>(new_fab, nullptr, ncomp, nghost, a_name, writeout, evolving);
     596          192 : }
     597              : 
     598              : }
     599              : #endif
        

Generated by: LCOV version 2.0-1