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

Class MethodOfLinesIntegrator<DIM> implements a spatially adaptive version of the Strong Stability Preserving (SSP) Runge-Kutta time integration algorithm. More...

#include <MethodOfLinesPatchStrategy.h>

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

Public Types

enum  MOL_VAR_TYPE { SOLN = 0, RHS = 1 }
 

Public Member Functions

 MethodOfLinesIntegrator (const std::string &object_name, tbox::Pointer< tbox::Database > input_db, MethodOfLinesPatchStrategy< DIM > *patch_strategy, bool register_for_restart=true)
 
virtual ~MethodOfLinesIntegrator ()
 
void initializeIntegrator (tbox::Pointer< mesh::GriddingAlgorithm< DIM > > gridding_alg)
 
double getTimestep (const tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy, const double time) const
 
void advanceHierarchy (const tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy, const double time, const double dt)
 
void registerVariable (const tbox::Pointer< hier::Variable< DIM > > variable, const hier::IntVector< DIM > &ghosts, const MOL_VAR_TYPE m_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
 
void initializeLevelData (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double init_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)
 
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 time, const int tag_index, const bool initial_time, const bool uses_richardson_extrapolation_too)
 
void putToDatabase (tbox::Pointer< tbox::Database > db)
 
virtual double getLevelDt (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const double dt_time, const bool initial_time)
 
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 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 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 > > coarser_level, const double coarsen_data_time, const bool before_advance)
 

Private Member Functions

void copyCurrentToScratch (const tbox::Pointer< hier::PatchLevel< DIM > > level) const
 
void copyScratchToCurrent (const tbox::Pointer< hier::PatchLevel< DIM > > level) const
 
void getFromInput (tbox::Pointer< tbox::Database > db, bool is_from_restart)
 
void getFromRestart ()
 

Private Attributes

std::string d_object_name
 
bool d_registered_for_restart
 
int d_order
 
tbox::Array< doubled_alpha_1
 
tbox::Array< doubled_alpha_2
 
tbox::Array< doubled_beta
 
MethodOfLinesPatchStrategy< DIM > * d_patch_strategy
 
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_fill_after_regrid
 
tbox::Pointer< xfer::RefineAlgorithm< DIM > > d_fill_before_tagging
 
tbox::Pointer< xfer::CoarsenAlgorithm< DIM > > d_coarsen_algorithm
 
tbox::Array< tbox::Pointer< xfer::CoarsenSchedule< DIM > > > d_coarsen_schedule
 
tbox::Pointer< hier::VariableContextd_current
 
tbox::Pointer< hier::VariableContextd_scratch
 
tbox::List< tbox::Pointer< hier::Variable< DIM > > > d_soln_variables
 
tbox::List< tbox::Pointer< hier::Variable< DIM > > > d_rhs_variables
 
hier::ComponentSelector d_current_data
 
hier::ComponentSelector d_scratch_data
 
hier::ComponentSelector d_rhs_data
 

Detailed Description

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

The original non-adaptive version of the algorithm is described in S. Gottlieb, C.W. Shu, E. Tadmor, SIAM Review, Vol. 43, No. 1, pp. 89-112. The advanceHierarchy() method integrates all levels of an AMR hierarchy through a specified timestep. See this method for details of the time-stepping process. Application-specific numerical routines that are necessary for these operations are provided by the MethodOfLinesPatchStrategy<DIM> data member. The collaboration between this class and the patch strategy follows the the Strategy design pattern. A concrete patch strategy object is derived from the base class to provide those routines for a specific problem.

This class is derived from the mesh::StandardTagAndInitStrategy<DIM> abstract base class which defines an interface for routines required by the dynamic adaptive mesh refinement routines in the mesh::GriddingAlgorithm<DIM> class. This collaboration also follows the Strategy design pattern.

Initialization of an MethodOfLinesIntegrator<DIM> object is performed by first setting default values, then reading from input. All input values may override values read from restart. Data read from input is summarized as follows:

Required input keys and data types: NONE

Optional input keys, data types, and defaults:

The following represents a sample input entry:

*  MethodOfLinesIntegrator{
*     order                 = 3         
*     alpha_1               = 1., 0.75, 0.33333
*     alpha_2               = 0., 0.25, 0.66666
*     beta                  = 1., 0.25, 0.66666
*  }
*  
See also
mesh::StandardTagAndInitStrategy

Member Enumeration Documentation

◆ MOL_VAR_TYPE

Enumerated type for the different categories of variable quantities allowed by the method of lines integration algorithm. See registerVariable(...) function for more details.

  • SOLN {Solution quantity for time-dependent ODE problem solved by RK time-stepping algorithm.}
  • RHS {Right-hand-side of ODE problem solved; i.e., du/dt = RHS.}
Enumerator
SOLN 
RHS 

Constructor & Destructor Documentation

◆ MethodOfLinesIntegrator()

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

The constructor for MethodOfLinesIntegrator<DIM> configures the method of lines integration algorithm with the concrete patch strategy object (containing problem-specific numerical routines) and initializes integration algorithm parameters provided in the specified input database and in 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.

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

◆ ~MethodOfLinesIntegrator()

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

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

Member Function Documentation

◆ initializeIntegrator()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< DIM >::initializeIntegrator ( tbox::Pointer< mesh::GriddingAlgorithm< DIM > >  gridding_alg)

Initialize integrator 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.

◆ getTimestep()

template<int DIM>
double SAMRAI::algs::MethodOfLinesIntegrator< DIM >::getTimestep ( const tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy,
const double  time 
) const

Return a suitable time increment over which to integrate the ODE problem. A minimum is taken over the increment computed on each patch in the hierarchy.

◆ advanceHierarchy()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< DIM >::advanceHierarchy ( const tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy,
const double  time,
const double  dt 
)

Advance the solution through the specified dt, which is assumed for the problem and state of the solution. Advances all patches in the hierarchy passed in.

◆ registerVariable()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< DIM >::registerVariable ( const tbox::Pointer< hier::Variable< DIM > >  variable,
const hier::IntVector< DIM > &  ghosts,
const MOL_VAR_TYPE  m_v_type,
const tbox::Pointer< xfer::Geometry< DIM > > &  transfer_geom,
const std::string &  coarsen_name = std::string(),
const std::string &  refine_name = std::string() 
)

Register variable quantity defined in the patch strategy with the method of lines integrator which manipulates its storage.

◆ printClassData()

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

Print all data members of MethodOfLinesIntegrator<DIM> object.

◆ initializeLevelData()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< DIM >::initializeLevelData ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  level_number,
const double  init_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 after_regrid is passed into the user's routines.

The boolean argument initial_time indicates whether the integration time corresponds to the initial simulation time. If true, the level should be initialized with initial simulation values. Otherwise, it should be assumed that the simulation time is at some point after the start of the simulation. 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 allowable level in the hierarchy.

Note: This function is overloaded from the base class mesh::StandardTagAndInitStrategy<DIM>.

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

◆ resetHierarchyConfiguration()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< 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 regridding, 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. Note, that whenever this routine is called, communication schedules are updated for every level finer than and including that indexed by coarsest_level.

Note: This function is overloaded from the base class mesh::StandardTagAndInitStrategy<DIM>.

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

◆ applyGradientDetector()

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

Set integer tags to "one" on the given level where refinement of that level should occur using the user-supplied gradient detector. The boolean argument initial_time is true when the level is being subject to error estimation at initialization time. If it is false, 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 data tagging routines since the application of the error estimator may be different in each of those circumstances.

Note: This function is overloaded from the base class mesh::StandardTagAndInitStrategy<DIM>.

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

◆ putToDatabase()

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

Writes object state out to the given database.

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

Implements SAMRAI::tbox::Serializable.

◆ copyCurrentToScratch()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< DIM >::copyCurrentToScratch ( const tbox::Pointer< hier::PatchLevel< DIM > >  level) const
private

◆ copyScratchToCurrent()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< DIM >::copyScratchToCurrent ( const tbox::Pointer< hier::PatchLevel< DIM > >  level) const
private

◆ getFromInput()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< DIM >::getFromInput ( tbox::Pointer< tbox::Database db,
bool  is_from_restart 
)
private

◆ getFromRestart()

template<int DIM>
void SAMRAI::algs::MethodOfLinesIntegrator< DIM >::getFromRestart ( )
private

◆ getLevelDt()

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

Determine time increment to advance data on level. 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.

This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

Reimplemented in SAMRAI::algs::HyperbolicLevelIntegrator< DIM >, and SAMRAI::algs::HyperbolicLevelIntegrator< NDIM >.

◆ advanceLevel()

template<int DIM>
virtual double SAMRAI::mesh::StandardTagAndInitStrategy< 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 
)
virtualinherited

Advance data on all patches on specified patch level from current time (current_time) to new time (new_time). This routine is called only during time-dependent regridding procedures, such as Richardson extrapolation. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed. The boolean arguments are used to determine the state of the algorithm and the data when the advance routine is called. Note that this advance function is also used during normal time integration steps.

When this function is called, the level data required to begin the advance must be allocated and be defined appropriately. Typically, this is equivalent to what is needed to initialize a new level after regridding. Upon exiting this routine, both current and new data may exist on the level. This data is needed until level synchronization occurs, in general. Current and new data may be reset by calling the member function resetTimeDependentData().

This routine is called from two different points within the Richardson exptrapolation process: 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.

Reimplemented in SAMRAI::algs::HyperbolicLevelIntegrator< DIM >, and SAMRAI::algs::HyperbolicLevelIntegrator< NDIM >.

◆ resetTimeDependentData()

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

Reset time-dependent data storage for the specified patch level.

This routine only applies when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

Reimplemented in SAMRAI::algs::HyperbolicLevelIntegrator< DIM >, and SAMRAI::algs::HyperbolicLevelIntegrator< NDIM >.

◆ resetDataToPreadvanceState()

template<int DIM>
virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::resetDataToPreadvanceState ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level)
virtualinherited

Reset data on the patch level by destroying all patch data other than that which is needed to initialize the solution on that level. In other words, this is the data needed to begin a time integration step on the level.

This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

Reimplemented in SAMRAI::algs::HyperbolicLevelIntegrator< DIM >, and SAMRAI::algs::HyperbolicLevelIntegrator< NDIM >.

◆ applyRichardsonExtrapolation()

template<int DIM>
virtual void SAMRAI::mesh::StandardTagAndInitStrategy< 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 
)
virtualinherited

Set integer tags to "one" in cells where refinement of 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.

This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

Reimplemented in SAMRAI::algs::HyperbolicLevelIntegrator< DIM >, and SAMRAI::algs::HyperbolicLevelIntegrator< NDIM >.

◆ coarsenDataForRichardsonExtrapolation()

template<int DIM>
virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::coarsenDataForRichardsonExtrapolation ( const tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy,
const int  level_number,
const tbox::Pointer< hier::PatchLevel< DIM > >  coarser_level,
const double  coarsen_data_time,
const bool  before_advance 
)
virtualinherited

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, once to set data on the coarser level and once to coarsen data from after advancing the fine level. The init_coarse_level boolean argument indicates whether data is set on the coarse level by coarsening the "old" time level solution or by coarsening the "new" solution on the fine level (i.e., after it has been advanced).

This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

Reimplemented in SAMRAI::algs::HyperbolicLevelIntegrator< DIM >, and SAMRAI::algs::HyperbolicLevelIntegrator< NDIM >.

Member Data Documentation

◆ d_object_name

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

◆ d_registered_for_restart

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

◆ d_order

template<int DIM>
int SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_order
private

◆ d_alpha_1

template<int DIM>
tbox::Array<double> SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_alpha_1
private

◆ d_alpha_2

template<int DIM>
tbox::Array<double> SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_alpha_2
private

◆ d_beta

template<int DIM>
tbox::Array<double> SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_beta
private

◆ d_patch_strategy

template<int DIM>
MethodOfLinesPatchStrategy<DIM>* SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_patch_strategy
private

◆ d_bdry_fill_advance

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

◆ d_bdry_sched_advance

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

◆ d_fill_after_regrid

template<int DIM>
tbox::Pointer< xfer::RefineAlgorithm<DIM> > SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_fill_after_regrid
private

◆ d_fill_before_tagging

template<int DIM>
tbox::Pointer< xfer::RefineAlgorithm<DIM> > SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_fill_before_tagging
private

◆ d_coarsen_algorithm

template<int DIM>
tbox::Pointer< xfer::CoarsenAlgorithm<DIM> > SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_coarsen_algorithm
private

◆ d_coarsen_schedule

template<int DIM>
tbox::Array< tbox::Pointer< xfer::CoarsenSchedule<DIM> > > SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_coarsen_schedule
private

◆ d_current

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

◆ d_scratch

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

◆ d_soln_variables

template<int DIM>
tbox::List< tbox::Pointer< hier::Variable<DIM> > > SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_soln_variables
private

◆ d_rhs_variables

template<int DIM>
tbox::List< tbox::Pointer< hier::Variable<DIM> > > SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_rhs_variables
private

◆ d_current_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_current_data
private

◆ d_scratch_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_scratch_data
private

◆ d_rhs_data

template<int DIM>
hier::ComponentSelector SAMRAI::algs::MethodOfLinesIntegrator< DIM >::d_rhs_data
private

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