IBAMR  IBAMR version 0.19.
Public Types | Public Member Functions | Protected Member Functions | Private Attributes | List of all members
SAMRAI::algs::HyperbolicLevelIntegrator< DIM > Class Template Reference

#include <ibamr/AdvectorPredictorCorrectorHyperbolicPatchOps.h>

Inheritance diagram for SAMRAI::algs::HyperbolicLevelIntegrator< DIM >:
Inheritance graph
[legend]

Public Types

enum  HYP_VAR_TYPE {
  TIME_DEP = 0, INPUT = 1, NO_FILL = 2, FLUX = 3,
  TEMPORARY = 4
}
 

Public Member Functions

 HyperbolicLevelIntegrator (const std::string &object_name, tbox::Pointer< tbox::Database > input_db, HyperbolicPatchStrategy< DIM > *patch_strategy, bool register_for_restart=true, bool use_time_refinement=true)
 
virtual ~HyperbolicLevelIntegrator ()
 
virtual void initializeLevelIntegrator (tbox::Pointer< mesh::BaseGriddingAlgorithm< DIM > > base_gridding_alg)
 
virtual double getLevelDt (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const double dt_time, const bool initial_time)
 
virtual double getMaxFinerLevelDt (const int finer_level_number, const double coarse_dt, const hier::IntVector< DIM > &ratio_to_coarser)
 
virtual double advanceLevel (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const double current_time, const double new_time, const bool first_step, const bool last_step, const bool regrid_advance=false)
 
virtual void standardLevelSynchronization (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int coarsest_level, const int finest_level, const double sync_time, const tbox::Array< double > &old_times)
 
virtual void standardLevelSynchronization (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int coarsest_level, const int finest_level, const double sync_time, const double old_time)
 
virtual void synchronizeNewLevels (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int coarsest_level, const int finest_level, const double sync_time, const bool initial_time)
 
virtual void resetTimeDependentData (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const double new_time, const bool can_be_refined)
 
virtual void resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< DIM > > level)
 
virtual void initializeLevelData (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double init_data_time, const bool can_be_refined, const bool initial_time, const tbox::Pointer< hier::BasePatchLevel< DIM > > old_level=tbox::Pointer< hier::BasePatchLevel< DIM > >(NULL), const bool allocate_data=true)
 
virtual void resetHierarchyConfiguration (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int coarsest_level, const int finest_level)
 
virtual void applyGradientDetector (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double error_data_time, const int tag_index, const bool initial_time, const bool uses_richardson_extrapolation_too)
 
virtual void applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< DIM > > level, const double error_data_time, const int tag_index, const double deltat, const int error_coarsen_ratio, const bool initial_time, const bool uses_gradient_detector_too)
 
virtual void coarsenDataForRichardsonExtrapolation (const tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< DIM > > coarse_level, const double coarsen_data_time, const bool before_advance)
 
virtual void registerVariable (const tbox::Pointer< hier::Variable< DIM > > var, const hier::IntVector< DIM > ghosts, const HYP_VAR_TYPE h_v_type, const tbox::Pointer< xfer::Geometry< DIM > > transfer_geom, const std::string &coarsen_name=std::string(), const std::string &refine_name=std::string())
 
virtual void printClassData (std::ostream &os) const
 
virtual void putToDatabase (tbox::Pointer< tbox::Database > db)
 
tbox::Pointer< hier::VariableContextgetCurrentContext () const
 
tbox::Pointer< hier::VariableContextgetNewContext () const
 
tbox::Pointer< hier::VariableContextgetOldContext () const
 
tbox::Pointer< hier::VariableContextgetScratchContext () const
 
tbox::Pointer< hier::VariableContextgetPlotContext () const
 
bool usingRefinedTimestepping () const
 
void printStatistics (std::ostream &s=tbox::plog) const
 

Protected Member Functions

virtual void getFromInput (tbox::Pointer< tbox::Database > db, bool is_from_restart)
 
virtual void getFromRestart ()
 
virtual void preprocessFluxData (const tbox::Pointer< hier::PatchLevel< DIM > > level, const double cur_time, const double new_time, const bool regrid_advance, const bool first_step, const bool last_step)
 
virtual void postprocessFluxData (const tbox::Pointer< hier::PatchLevel< DIM > > level, const bool regrid_advance, const bool first_step, const bool last_step)
 
virtual void copyTimeDependentData (const tbox::Pointer< hier::PatchLevel< DIM > > level, const tbox::Pointer< hier::VariableContext > src_context, const tbox::Pointer< hier::VariableContext > dst_context)
 
virtual void synchronizeLevelWithCoarser (const tbox::Pointer< hier::PatchLevel< DIM > > fine, const tbox::Pointer< hier::PatchLevel< DIM > > coarse, const double sync_time, const double coarse_sim_time)
 
void recordStatistics (const hier::PatchLevel< DIM > &patch_level, double current_time)
 

Private Attributes

HyperbolicPatchStrategy< DIM > * d_patch_strategy
 
std::string d_object_name
 
bool d_use_time_refinement
 
bool d_registered_for_restart
 
double d_cfl
 
double d_cfl_init
 
bool d_lag_dt_computation
 
bool d_use_ghosts_for_dt
 
bool d_flux_is_face
 
bool d_flux_face_registered
 
bool d_flux_side_registered
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_bdry_fill_advance
 
tbox::Array< tbox::Pointer< xfer::RefineSchedule< DIM > > > d_bdry_sched_advance
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_bdry_fill_advance_new
 
tbox::Array< tbox::Pointer< xfer::RefineSchedule< DIM > > > d_bdry_sched_advance_new
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_bdry_fill_advance_old
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_coarsen_fluxsum
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_coarsen_sync_data
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_sync_initial_data
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_coarsen_rich_extrap_init
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_coarsen_rich_extrap_final
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_fill_new_level
 
int d_number_time_data_levels
 
tbox::Pointer< hier::VariableContextd_scratch
 
tbox::Pointer< hier::VariableContextd_current
 
tbox::Pointer< hier::VariableContextd_new
 
tbox::Pointer< hier::VariableContextd_old
 
tbox::Pointer< hier::VariableContextd_plot_context
 
tbox::List< tbox::Pointer< hier::Variable< DIM > > > d_all_variables
 
tbox::List< tbox::Pointer< hier::Variable< DIM > > > d_time_dep_variables
 
tbox::List< tbox::Pointer< hier::Variable< DIM > > > d_flux_variables
 
tbox::List< tbox::Pointer< hier::Variable< DIM > > > d_fluxsum_variables
 
hier::ComponentSelector d_saved_var_scratch_data
 
hier::ComponentSelector d_temp_var_scratch_data
 
hier::ComponentSelector d_new_patch_init_data
 
hier::ComponentSelector d_new_time_dep_data
 
hier::ComponentSelector d_flux_var_data
 
hier::ComponentSelector d_fluxsum_data
 
bool d_have_flux_on_level_zero
 
hier::ComponentSelector d_old_time_dep_data
 
bool d_distinguish_mpi_reduction_costs
 
tbox::Pointer< tbox::Timert_advance_bdry_fill_comm
 
tbox::Pointer< tbox::Timert_error_bdry_fill_create
 
tbox::Pointer< tbox::Timert_error_bdry_fill_comm
 
tbox::Pointer< tbox::Timert_mpi_reductions
 
tbox::Pointer< tbox::Timert_initialize_level_data
 
tbox::Pointer< tbox::Timert_fill_new_level_create
 
tbox::Pointer< tbox::Timert_fill_new_level_comm
 
tbox::Pointer< tbox::Timert_advance_bdry_fill_create
 
tbox::Pointer< tbox::Timert_new_advance_bdry_fill_create
 
tbox::Pointer< tbox::Timert_apply_gradient_detector
 
tbox::Pointer< tbox::Timert_coarsen_rich_extrap
 
tbox::Pointer< tbox::Timert_get_level_dt
 
tbox::Pointer< tbox::Timert_get_level_dt_sync
 
tbox::Pointer< tbox::Timert_advance_level
 
tbox::Pointer< tbox::Timert_new_advance_bdry_fill_comm
 
tbox::Pointer< tbox::Timert_patch_num_kernel
 
tbox::Pointer< tbox::Timert_advance_level_sync
 
tbox::Pointer< tbox::Timert_std_level_sync
 
tbox::Pointer< tbox::Timert_sync_new_levels
 
tbox::Array< tbox::Pointer< tbox::Statistic > > d_boxes_stat
 
tbox::Array< tbox::Pointer< tbox::Statistic > > d_cells_stat
 
tbox::Array< tbox::Pointer< tbox::Statistic > > d_timestamp_stat
 

Detailed Description

template<int DIM>
class SAMRAI::algs::HyperbolicLevelIntegrator< DIM >

Class HyperbolicLevelIntegrator<DIM> provides routines needed to integrate a system of hyperbolic conservation laws on a structured AMR patch hierarchy using local time refinement. The routines include initializing a level, advance a level, and synchronize levels in a time-dependent AMR application. The AMR timestepping algorithm that cycles through the patch levels and calls these routines is provided by the TimeRefinementIntegrator<DIM> class. Together, that hierarchy integration class and this single level integration class produce the common AMR algorithm due to Berger, Colella and Oliger (see e.g., Berger and Colella, J. Comp. Phys. (82)1:64-84, 1989). The operations performed on single patches on each level are implemented in the user-defined, problem-specific class derived from the abstract base class HyperbolicPatchStrategy<DIM>.

It is important to note that the variable contexts used by the concrete patch strategy subclass must be consistent with those defined in this class which manages the data for the variables.

This class is derived from the abstract base class TimeRefinementLevelStrategy<DIM>, which defines routines needed by the time refinement integrator. There is an argument in the constructor that determines whether this class will be used by the time refinement integrator for refined timestepping or synchronized timestepping. The routines overloaded in TimeRefinementLevelStrategy<DIM> are: initializeLevelIntegrator(), getLevelDt(), getMaxFinerLevelDt(), advanceLevel(), standardLevelSynchronization(), synchronizeNewLevels(), resetTimeDependentData(), and resetDataToPreadvanceState(). This class is also derived from mesh::StandardTagAndInitStrategy<DIM>, which defines routines needed by the gridding algorithm classes. The routines overloaded in mesh::StandardTagAndInitStrategy<DIM> are: initializeLevelData(), resetHierarchyConfiguration(), applyGradientDetector(), applyRichardsonExtrapolation(), and coarsenDataForRichardsonExtrapolation().

An object of this class requires numerous parameters to be read from input. Also, data must be written to and read from files for restart. The input and restart data are summarized as follows.

Required input keys and data types: NONE

Optional input keys, data types, and defaults:

 - \b    cfl    
    double value for the CFL factor used for timestep selection 
    (dt used = CFL * max dt).  If no input value is given, a default 
    value of 0.9 is used.

 - \b    cfl_init   
    double value for CFL factor used for initial timestep.
    If no input value is given, a default value of 0.9 is used.

 - \b    lag_dt_computation   
    boolean value indicating whether dt is based on current
    solution or solution from previous step (possible optimization
    in communication for characteristic analysis).  If no input
    value is given, a default value of TRUE is used.


 - \b    use_ghosts_to_compute_dt   
    boolean value indicating whether ghost data must be filled before 
    timestep is computed on each patch (possible communication 
    optimization).  if no input value is given, a default value
    of TRUE is used.

 - \b    distinguish_mpi_reduction_costs   
    boolean specifying whether to separate reduction costs in tbox::MPI
    from costs of load imbalances.  By specifying it true, a
    barrier is put in place before the reduction call, so an extra
    operation is incurred.  For this reason, it is defaulted FALSE.

Note that when continuing from restart, the input values in the input file override all values read in from the restart database.

A sample input file entry might look like:

* 
*    cfl = 0.9
*    cfl_init = 0.9
*    lag_dt_computation = FALSE
*    use_ghosts_to_compute_dt = TRUE
*    distinguish_mpi_reduction_costs = TRUE
*
* 
See also
algs::TimeRefinementIntegrator
mesh::StandardTagAndInitStrategy
algs::HyperbolicPatchStrategy

Member Enumeration Documentation

◆ HYP_VAR_TYPE

Enumerated type for the different ways in which variable storage can be manipulated by the level integration algorithm. See registerVariable(...) function for more details.

  • TIME_DEP {Data that changes in time and needs more than one time level to be stored.}
  • INPUT {Data that is set once and do not change during the ghosts are never re-filled outside of user-defined routines.}
  • FLUX {Face-centered double values used in conservative difference and synchronization (i.e., refluxing) process. A corresponding variable to store flux integral information is created for each FLUX variable.}
  • TEMPORARY {Accessory values intended to live only for computation on a single patch (i.e., they cannot be assumed to exist between patch routine function calls.)}
Enumerator
TIME_DEP 
INPUT 
NO_FILL 
FLUX 
TEMPORARY 

Constructor & Destructor Documentation

◆ HyperbolicLevelIntegrator()

template<int DIM>
SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::HyperbolicLevelIntegrator ( const std::string &  object_name,
tbox::Pointer< tbox::Database input_db,
HyperbolicPatchStrategy< DIM > *  patch_strategy,
bool  register_for_restart = true,
bool  use_time_refinement = true 
)

Constructor for HyperbolicLevelIntegrator<DIM> initializes integration parameters to default values and constructs standard communication algorithms. Other data members are read in from the specified input database or the restart database corresponding to the specified object_name. The constructor also registers this object for restart using the specified object name when the boolean argument is true. Whether object will write its state to restart files during program execution is determined by this argument. Note that it has a default state of true. This class is used by the time refinement integrator for refined timestepping when the use_time_refinement argument is true, and for synchronized timestepping when the boolean is false.

When assertion checking is active, passing in any null pointer or an empty string will result in an unrecoverable exception.

◆ ~HyperbolicLevelIntegrator()

template<int DIM>
virtual SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::~HyperbolicLevelIntegrator ( )
virtual

The destructor for HyperbolicLevelIntegrator<DIM> unregisters the integrator object with the restart manager when so registered.

Member Function Documentation

◆ initializeLevelIntegrator()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::initializeLevelIntegrator ( tbox::Pointer< mesh::BaseGriddingAlgorithm< DIM > >  base_gridding_alg)
virtual

Initialize level integrator by by setting the number of time levels of data needed based on specifications of the gridding algorithm.

This routine also invokes variable registration in the patch strategy.

Assertion checking will throw unrecoverable exceptions if either pointer is null.

Implements SAMRAI::algs::TimeRefinementLevelStrategy< DIM >.

◆ getLevelDt()

template<int DIM>
virtual double SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getLevelDt ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level,
const double  dt_time,
const bool  initial_time 
)
virtual

Determine time increment to advance data on level and return that value. The double dt_time argument is the simulation time when the routine is called. The initial_time boolean is true if this routine is called during hierarchy initialization (i.e., at the initial simulation time). Otherwise, it is false. The recompute_dt option specifies whether to compute the timestep using the current level data or to return the value stored by the time integrator. The default true setting means the timestep will be computed if no value is supplied.

When assertion checking is active, an unrecoverable exception will result if the level pointer is null.

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ getMaxFinerLevelDt()

template<int DIM>
virtual double SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getMaxFinerLevelDt ( const int  finer_level_number,
const double  coarse_dt,
const hier::IntVector< DIM > &  ratio_to_coarser 
)
virtual

Return the maximum allowable time increment for the level with the specified level number based on the time increment for the next coarser level and the mesh refinement ratio between the two levels. For the common explicit integration methods for hyperbolic conservation laws (constrained by a CFL limit), the fine time increment is typically the coarse increment divided by the refinement ratio.

When assertion checking is active, an unrecoverable exception will result if the ratio vector is not acceptable (i.e., all values > 0).

Reimplemented from SAMRAI::algs::TimeRefinementLevelStrategy< DIM >.

◆ advanceLevel()

template<int DIM>
virtual double SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::advanceLevel ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level,
const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const double  current_time,
const double  new_time,
const bool  first_step,
const bool  last_step,
const bool  regrid_advance = false 
)
virtual

Integrate data on all patches on the given patch level from current time (current_time) to new time (new_time). This routine is used to advance the solution on each level in the hierarchy and during time-dependent regridding procedures, such as Richardson extrapolation. The boolean arguments are used to determine the state of the algorithm and the data when the advance routine is called. The first_step and last_step indicate whether the step is the first or last in the current timestep sequence on the level. Typically, the current timestep sequence means each step on the level between advance steps on a coarser level in the hierarchy, if one exists. The regrid_advance value is true when the advance is called as part of a time-dependent regridding procedure. Usually when this happens, the results of the colution advance will be discarded. So, for example, when this is true flux information is not maintained and flux integrals are not updated. The final boolean argument indicates whether or not the level resides in the hierarchy. For example, during time-dependent regridding, such as Richardson extrapolation, a temporary level that is not in the hierarchy is created and advanced. Then, a communication schedule must be generated for the level before the advance begins.

This routine is called at two different points during time integration: during the regular time advance sequence, and possibly at the initial simulation time. The second call advances the solution on a coarser level ahead in time to provide time-dependent boundary values for some finer level when time-dependent regridding is used. In the first case, the values of the boolean flags are:

  • first_step = true for first step in level time step sequence; else, false.
  • last_step = true for last step in level time step sequence; else, false.
  • regrid_advance = false.

In the second case, the values of the boolean flags are:

  • first_step = true.
  • last_step = false.
  • regrid_advance = true.

When time-dependent regridding (i.e., Richardson extrapolation) is used, the routine is called from two different points in addition to those described above: to advance a temporary level that is coarser than the hierarchy level on which error estimation is performed, and to advance the hierarchy level itself. In the first case, the values of the boolean flags are:

  • first_step = true.
  • last_step = true.
  • regrid_advance = true.

In the second case, the values of the boolean flags are:

  • first_step (when regridding during time integration sequence) = true when the level is not coarsest level to synchronize immediately before the regridding process; else, false. (when generating initial hierarchy construction) = true, even though there may be multiple advance steps.
  • last_step = true when the advance is the last in the Richardson extrapolation step sequence; else false.
  • regrid_advance = true.

When assertion checking is active, an unrecoverable exception will result if either the level or hierarchy pointer is null, or the new time is not greater than the given time.

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ standardLevelSynchronization() [1/2]

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::standardLevelSynchronization ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  coarsest_level,
const int  finest_level,
const double  sync_time,
const tbox::Array< double > &  old_times 
)
virtual

Synchronize data between given patch levels in patch hierarchy according to the standard hyperbolic AMR flux correction algorithm. This routine synchronizes data between two levels at a time from the level with index finest_level down to the level with index coarsest_level. The array of old time values are used in the
re-integration of the time-dependent data.

When assertion checking is active, an unrecoverable exception will result if the hierarchy pointer is null, the level numbers do not properly match existing levels in the hierarchy (either coarsest_level > finest_level or some level is null), or all of the old time values are less than the value of sync_time.

Reimplemented from SAMRAI::algs::TimeRefinementLevelStrategy< DIM >.

◆ standardLevelSynchronization() [2/2]

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::standardLevelSynchronization ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  coarsest_level,
const int  finest_level,
const double  sync_time,
const double  old_time 
)
virtual

This overloaded version of standardLevelSynchronization implements a routine used for synchronized timestepping. Only a single value for the old time is needed, since all levels would have the same old time.

Reimplemented from SAMRAI::algs::TimeRefinementLevelStrategy< DIM >.

◆ synchronizeNewLevels()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::synchronizeNewLevels ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  coarsest_level,
const int  finest_level,
const double  sync_time,
const bool  initial_time 
)
virtual

Coarsen current solution data from finest hierarchy level specified down through the coarsest hierarchy level specified, if initial_time is true. In this case, the hierarchy is being constructed at the initial simulation time, After data is coarsened, the application- specific initialization routine is called to set data before that solution is further coarsened to the next coarser level in the hierarchy. This operation makes the solution consistent between coarser levels and finer levels that did not exist when the coarse levels where created and initialized originally.

When initial_time is false, this routine does nothing since the standard hyperbolic AMR algorithm for conservation laws requires no data synchronization after regridding beyond interpolation of data from coarser levels in the hierarchy in some conservative fashion.

When assertion checking is active, an unrecoverable exception will result if the hierarchy pointer is null, the level numbers do not properly match existing levels in the hierarchy (either coarsest_level > finest_level or some level is null).

Implements SAMRAI::algs::TimeRefinementLevelStrategy< DIM >.

◆ resetTimeDependentData()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::resetTimeDependentData ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level,
const double  new_time,
const bool  can_be_refined 
)
virtual

Resets time-dependent data storage and update time for patch level.

When assertion checking is active, an unrecoverable exception will result if the level pointer is null.

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ resetDataToPreadvanceState()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::resetDataToPreadvanceState ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level)
virtual

Deallocate all new simulation data on the given level. This may be necessary during regridding, or setting up levels initially.

When assertion checking is active, an unrecoverable exception will result if the level pointer is null.

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ initializeLevelData()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::initializeLevelData ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  level_number,
const double  init_data_time,
const bool  can_be_refined,
const bool  initial_time,
const tbox::Pointer< hier::BasePatchLevel< DIM > >  old_level = tbox::Pointerhier::BasePatchLevel< DIM > >(NULL),
const bool  allocate_data = true 
)
virtual

Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm. The level number indicates that of the new level. The old_level pointer corresponds to the level that resided in the hierarchy before the level with the specified number was introduced. If the pointer is null, there was no level in the hierarchy prior to the call and the level data is set based on the user routines and the simulation time. Otherwise, the specified level replaces the old level and the new level receives data from the old level appropriately before it is destroyed.

Typically, when data is set, it is interpolated from coarser levels in the hierarchy. If the data is to be set, the level number must match that of the old level, if non-null. If the old level is non-null, then data is copied from the old level to the new level on regions of intersection between those levels before interpolation occurs. Then, user-supplied patch routines are called to further initialize the data if needed. The boolean argument initial_time is passed into the user's routines.

The boolean argument initial_time indicates whether the level is being introduced for the first time (i.e., at initialization time), or after some regrid process during the calculation beyond the initial hierarchy construction. This information is provided since the initialization of the data on a patch may be different in each of those circumstances. The can_be_refined boolean argument indicates whether the level is the finest level allowed in the hierarchy. This may or may not affect the data initialization process depending on the problem.

When assertion checking is active, an unrecoverable exception will result if the hierarchy pointer is null, the level number does not match any level in the hierarchy, or the old level number does not match the level number (if the old level pointer is non-null).

Implements SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ resetHierarchyConfiguration()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::resetHierarchyConfiguration ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  coarsest_level,
const int  finest_level 
)
virtual

Reset cached communication schedules after the hierarchy has changed (due to regidding, for example) and the data has been initialized on the new levels. The intent is that the cost of data movement on the hierarchy will be amortized across multiple communication cycles, if possible. The level numbers indicate the range of levels in the hierarchy that have changed. However, this routine updates communication schedules every level finer than and including that indexed by the coarsest level number given.

When assertion checking is active, an unrecoverable exception will result if the hierarchy pointer is null, any pointer to a level in the hierarchy that is coarser than the finest level is null, or the given level numbers not specified properly; e.g., coarsest_level > finest_level.

Implements SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ applyGradientDetector()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::applyGradientDetector ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  level_number,
const double  error_data_time,
const int  tag_index,
const bool  initial_time,
const bool  uses_richardson_extrapolation_too 
)
virtual

Set integer tags to "one" in cells where refinement of the given level should occur according to some user-supplied gradient criteria. The double time argument is the regrid time. The integer "tag_index" argument is the patch descriptor index of the cell-centered integer tag array on each patch in the hierarchy. The boolean argument initial_time indicates whether the level is being subject to refinement at the initial simulation time. If it is false, then the error estimation process is being invoked at some later time after the AMR hierarchy was initially constructed. The boolean argument uses_richardson_extrapolation_too is true when Richardson extrapolation error estimation is used in addition to the gradient detector, and false otherwise. This argument helps the user to manage multiple regridding criteria. This information is passed along to the user's patch tagging routines since the application of the gradient detector may be different in each case.

When assertion checking is active, an unrecoverable exception will result if the hierarchy pointer is null or the level number does not match any existing level in the hierarchy.

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ applyRichardsonExtrapolation()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::applyRichardsonExtrapolation ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const double  error_data_time,
const int  tag_index,
const double  deltat,
const int  error_coarsen_ratio,
const bool  initial_time,
const bool  uses_gradient_detector_too 
)
virtual

Set integer tags to "one" where refinement onf the given level should occur according to some user-supplied Richardson extrapolation criteria. The "error_data_time" argument is the regrid time. The "deltat" argument is the time increment to advance the solution on the level to be refined. Note that that level is finer than the level in the argument list, in general. The ratio between the argument level and the actual hierarchy level is given by the integer "coarsen ratio".

The integer "tag_index" argument is the patch descriptor index of the cell-centered integer tag array on each patch in the hierarchy.

The boolean argument initial_time indicates whether the level is being subject to refinement at the initial simulation time. If it is false, then the error estimation process is being invoked at some later time after the AMR hierarchy was initially constructed. Typically, this information is passed to the user's patch tagging routines since the application of the Richardson extrapolation process may be different in each case.

The boolean uses_gradient_detector_too is true when a gradient detector procedure is used in addition to Richardson extrapolation, and false otherwise. This argument helps the user to manage multiple regridding criteria.

When assertion checking is active, an unrecoverable exception will result if the level pointer is null.

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ coarsenDataForRichardsonExtrapolation()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::coarsenDataForRichardsonExtrapolation ( const tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy,
const int  level_number,
const tbox::Pointer< hier::PatchLevel< DIM > >  coarse_level,
const double  coarsen_data_time,
const bool  before_advance 
)
virtual

Coarsen solution data from level to coarse_level for Richardson extrapolation. Note that this routine will be called twice during the Richardson extrapolation error estimation process. The before_advance boolean argument indicates whether data is set on the coarse level by coarsening the "old" time level solution (i.e., before it has been advanced) or by coarsening the "new" solution on the fine level (i.e., after it has been advanced).

When assertion checking is active, an unrecoverable exception will result if either level pointer is null.

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< DIM >.

◆ registerVariable()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::registerVariable ( const tbox::Pointer< hier::Variable< DIM > >  var,
const hier::IntVector< DIM >  ghosts,
const HYP_VAR_TYPE  h_v_type,
const tbox::Pointer< xfer::Geometry< DIM > >  transfer_geom,
const std::string &  coarsen_name = std::string(),
const std::string &  refine_name = std::string() 
)
virtual

Register a variable with the hyperbolic integration algorithm. The variable type must be one of the options defined by the enumerated type defined above. Typically, this routine is called from the hyperbolic patch model when the variable registration process is invoked by calling the function initializeLevelIntegrator() above.
In fact, that function should be called before this routine is called.

When assertion checking is active, an unrecoverable exception will result if the variable pointer or geometry pointer is null.

◆ printClassData()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::printClassData ( std::ostream &  os) const
virtual

Print class data representation for hyperbolic level integrator object. This is done automatically, when an unrecoverable run-time exception is thrown within some member function of this class.

◆ putToDatabase()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::putToDatabase ( tbox::Pointer< tbox::Database db)
virtual

Write out object state to the given database.

When assertion checking is active, database point must be non-null.

Implements SAMRAI::tbox::Serializable.

◆ getCurrentContext()

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getCurrentContext ( ) const

Return pointer to "current" variable context used by integrator. Current data corresponds to state data at the beginning of a timestep, or when a new level is initialized.

◆ getNewContext()

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getNewContext ( ) const

Return pointer to "new" variable context used by integrator. New data corresponds to advanced state data at the end of a timestep. The data is one timestep later than the "current" data.

◆ getOldContext()

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getOldContext ( ) const

Return pointer to "old" variable context used by integrator. Old data corresponds to an extra time level of state data used for Richardson extrapolation error estimation. The data is one timestep earlier than the "current" data.

Note that only in certain cases when using time-dependent error estimation, such as Richardson extrapolation, is the returned pointer will non-null. See contructor for more information.

◆ getScratchContext()

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getScratchContext ( ) const

Return pointer to "scratch" variable context used by integrator. Scratch data typically corresponds to storage that user-routines in the concrete HyperbolicPatchStrategy<DIM> object manipulate; in particular, scratch data contains ghost cells.

◆ getPlotContext()

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getPlotContext ( ) const

Return pointer to variable context used for plotting. This context corresponds to the data storage that should be written to plot files. Typically, this is the same as the "current" context.

◆ usingRefinedTimestepping()

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::usingRefinedTimestepping ( ) const
virtual

Return true if this class has been constructed to use refined timestepping and false if it has been constructed to use synchronized timestepping.

Implements SAMRAI::algs::TimeRefinementLevelStrategy< DIM >.

◆ printStatistics()

template<int DIM>
void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::printStatistics ( std::ostream &  s = tbox::plog) const

◆ getFromInput()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getFromInput ( tbox::Pointer< tbox::Database db,
bool  is_from_restart 
)
protectedvirtual

Read input values, indicated above, from given database. The boolean argument is_from_restart should be set to true if the simulation is beginning from restart. Otherwise it should be set to false.

When assertion checking is active, the database pointer must be non-null.

◆ getFromRestart()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::getFromRestart ( )
protectedvirtual

Read object state from the restart file and initialize class data members. The database from which the restart data is read is determined by the object_name specified in the constructor.

Unrecoverable Errors:

  • The database corresponding to object_name is not found in the restart file.
  • The class version number and restart version number do not match.

◆ preprocessFluxData()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::preprocessFluxData ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const double  cur_time,
const double  new_time,
const bool  regrid_advance,
const bool  first_step,
const bool  last_step 
)
protectedvirtual

◆ postprocessFluxData()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::postprocessFluxData ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const bool  regrid_advance,
const bool  first_step,
const bool  last_step 
)
protectedvirtual

◆ copyTimeDependentData()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::copyTimeDependentData ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const tbox::Pointer< hier::VariableContext src_context,
const tbox::Pointer< hier::VariableContext dst_context 
)
protectedvirtual

◆ synchronizeLevelWithCoarser()

template<int DIM>
virtual void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::synchronizeLevelWithCoarser ( const tbox::Pointer< hier::PatchLevel< DIM > >  fine,
const tbox::Pointer< hier::PatchLevel< DIM > >  coarse,
const double  sync_time,
const double  coarse_sim_time 
)
protectedvirtual

Apply the standard AMR hyperbolic flux synchronization process preserve conservation properties in the solution between the fine level and the coarse level. The sync_time argument indicates the time at which the solution on the two levels is being synchronized. The variable coarse_sim_time indicates the previous simulation time on the coarse level (recall the conservative difference will be repeated on the coarse level during the synchronization process). After the synchronization, the flux and flux integral data storage is reset on the levels.

When assertion checking is turned on, an unrecoverable exception will result if either level pointer is null, the levels are not consecutive in the AMR hierarchy, or the coarse sim time is not less than the sync time.

◆ recordStatistics()

template<int DIM>
void SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::recordStatistics ( const hier::PatchLevel< DIM > &  patch_level,
double  current_time 
)
protected

Member Data Documentation

◆ d_patch_strategy

template<int DIM>
HyperbolicPatchStrategy<DIM>* SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_patch_strategy
private

◆ d_object_name

template<int DIM>
std::string SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_object_name
private

◆ d_use_time_refinement

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_use_time_refinement
private

◆ d_registered_for_restart

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_registered_for_restart
private

◆ d_cfl

template<int DIM>
double SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_cfl
private

◆ d_cfl_init

template<int DIM>
double SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_cfl_init
private

◆ d_lag_dt_computation

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_lag_dt_computation
private

◆ d_use_ghosts_for_dt

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_use_ghosts_for_dt
private

◆ d_flux_is_face

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_flux_is_face
private

◆ d_flux_face_registered

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_flux_face_registered
private

◆ d_flux_side_registered

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_flux_side_registered
private

◆ d_bdry_fill_advance

template<int DIM>
tbox::Pointer< xfer::RefineAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_bdry_fill_advance
private

◆ d_bdry_sched_advance

template<int DIM>
tbox::Array< tbox::Pointer< xfer::RefineSchedule<DIM> > > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_bdry_sched_advance
private

◆ d_bdry_fill_advance_new

template<int DIM>
tbox::Pointer< xfer::RefineAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_bdry_fill_advance_new
private

◆ d_bdry_sched_advance_new

template<int DIM>
tbox::Array< tbox::Pointer< xfer::RefineSchedule<DIM> > > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_bdry_sched_advance_new
private

◆ d_bdry_fill_advance_old

template<int DIM>
tbox::Pointer< xfer::RefineAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_bdry_fill_advance_old
private

◆ d_coarsen_fluxsum

template<int DIM>
tbox::Pointer< xfer::CoarsenAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_coarsen_fluxsum
private

◆ d_coarsen_sync_data

template<int DIM>
tbox::Pointer< xfer::CoarsenAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_coarsen_sync_data
private

◆ d_sync_initial_data

template<int DIM>
tbox::Pointer< xfer::CoarsenAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_sync_initial_data
private

◆ d_coarsen_rich_extrap_init

template<int DIM>
tbox::Pointer< xfer::CoarsenAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_coarsen_rich_extrap_init
private

◆ d_coarsen_rich_extrap_final

template<int DIM>
tbox::Pointer< xfer::CoarsenAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_coarsen_rich_extrap_final
private

◆ d_fill_new_level

template<int DIM>
tbox::Pointer< xfer::RefineAlgorithm<DIM> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_fill_new_level
private

◆ d_number_time_data_levels

template<int DIM>
int SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_number_time_data_levels
private

◆ d_scratch

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_scratch
private

◆ d_current

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_current
private

◆ d_new

◆ d_old

◆ d_plot_context

template<int DIM>
tbox::Pointer<hier::VariableContext> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_plot_context
private

◆ d_all_variables

template<int DIM>
tbox::List< tbox::Pointer< hier::Variable<DIM> > > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_all_variables
private

◆ d_time_dep_variables

template<int DIM>
tbox::List< tbox::Pointer< hier::Variable<DIM> > > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_time_dep_variables
private

◆ d_flux_variables

template<int DIM>
tbox::List< tbox::Pointer< hier::Variable<DIM> > > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_flux_variables
private

◆ d_fluxsum_variables

template<int DIM>
tbox::List< tbox::Pointer< hier::Variable<DIM> > > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_fluxsum_variables
private

◆ d_saved_var_scratch_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_saved_var_scratch_data
private

◆ d_temp_var_scratch_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_temp_var_scratch_data
private

◆ d_new_patch_init_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_new_patch_init_data
private

◆ d_new_time_dep_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_new_time_dep_data
private

◆ d_flux_var_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_flux_var_data
private

◆ d_fluxsum_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_fluxsum_data
private

◆ d_have_flux_on_level_zero

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_have_flux_on_level_zero
private

◆ d_old_time_dep_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_old_time_dep_data
private

◆ d_distinguish_mpi_reduction_costs

template<int DIM>
bool SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_distinguish_mpi_reduction_costs
private

◆ t_advance_bdry_fill_comm

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_advance_bdry_fill_comm
private

◆ t_error_bdry_fill_create

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_error_bdry_fill_create
private

◆ t_error_bdry_fill_comm

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_error_bdry_fill_comm
private

◆ t_mpi_reductions

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_mpi_reductions
private

◆ t_initialize_level_data

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_initialize_level_data
private

◆ t_fill_new_level_create

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_fill_new_level_create
private

◆ t_fill_new_level_comm

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_fill_new_level_comm
private

◆ t_advance_bdry_fill_create

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_advance_bdry_fill_create
private

◆ t_new_advance_bdry_fill_create

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_new_advance_bdry_fill_create
private

◆ t_apply_gradient_detector

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_apply_gradient_detector
private

◆ t_coarsen_rich_extrap

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_coarsen_rich_extrap
private

◆ t_get_level_dt

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_get_level_dt
private

◆ t_get_level_dt_sync

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_get_level_dt_sync
private

◆ t_advance_level

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_advance_level
private

◆ t_new_advance_bdry_fill_comm

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_new_advance_bdry_fill_comm
private

◆ t_patch_num_kernel

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_patch_num_kernel
private

◆ t_advance_level_sync

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_advance_level_sync
private

◆ t_std_level_sync

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_std_level_sync
private

◆ t_sync_new_levels

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::t_sync_new_levels
private

◆ d_boxes_stat

template<int DIM>
tbox::Array< tbox::Pointer<tbox::Statistic> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_boxes_stat
private

◆ d_cells_stat

template<int DIM>
tbox::Array< tbox::Pointer<tbox::Statistic> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_cells_stat
private

◆ d_timestamp_stat

template<int DIM>
tbox::Array< tbox::Pointer<tbox::Statistic> > SAMRAI::algs::HyperbolicLevelIntegrator< DIM >::d_timestamp_stat
private

The documentation for this class was generated from the following file: