Rythmos - Transient Integration for Differential Equations
Version of the Day
|
Base class for defining stepper functionality. More...
#include <Rythmos_StepperBase_decl.hpp>
Public Member Functions | |
virtual bool | supportsCloning () const |
Return if this stepper supports cloning or not. More... | |
virtual RCP< StepperBase < Scalar > > | cloneStepperAlgorithm () const |
Clone the stepper object if supported. More... | |
virtual bool | isImplicit () const |
Return if this stepper is an implicit stepper. More... | |
virtual bool | acceptsModel () const |
Return if this stepper accepts a model. More... | |
virtual void | setModel (const RCP< const Thyra::ModelEvaluator< Scalar > > &model)=0 |
Specify the model problem to integrate. More... | |
virtual void | setNonconstModel (const RCP< Thyra::ModelEvaluator< Scalar > > &model)=0 |
Accept a nonconst model. More... | |
virtual bool | modelIsConst () const |
Return of the model is only const or can be returned as a non-const object. More... | |
virtual RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () const =0 |
Get the model. More... | |
virtual RCP < Thyra::ModelEvaluator < Scalar > > | getNonconstModel ()=0 |
Get the model nonconst. More... | |
virtual void | setInitialCondition (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &initialCondition)=0 |
Specify initial condition and re-initialize. More... | |
virtual Thyra::ModelEvaluatorBase::InArgs < Scalar > | getInitialCondition () const =0 |
Get the currently set initial condtion. More... | |
virtual Scalar | takeStep (Scalar dt, StepSizeType stepType)=0 |
Take a step. More... | |
virtual const StepStatus< Scalar > | getStepStatus () const =0 |
Get current stepper status after a step has been taken. More... | |
virtual void | setStepControlData (const StepperBase &stepper) |
Set step control data from another stepper. More... | |
Public Member Functions inherited from Rythmos::InterpolationBufferBase< Scalar > | |
virtual RCP< const Thyra::VectorSpaceBase< Scalar > > | get_x_space () const =0 |
Return the space for x and x_dot . More... | |
virtual void | addPoints (const Array< Scalar > &time_vec, const Array< RCP< const Thyra::VectorBase< Scalar > > > &x_vec, const Array< RCP< const Thyra::VectorBase< Scalar > > > &xdot_vec)=0 |
Add points to the buffer. More... | |
virtual TimeRange< Scalar > | getTimeRange () const =0 |
Return the range of time values where interpolation calls can be performed. More... | |
virtual void | getPoints (const Array< Scalar > &time_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *x_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *xdot_vec, Array< ScalarMag > *accuracy_vec) const =0 |
Get values from the buffer at different time points. More... | |
virtual void | getNodes (Array< Scalar > *time_vec) const =0 |
Get interpolation nodes. More... | |
virtual void | removeNodes (Array< Scalar > &time_vec)=0 |
Remove nodes from the interpolation buffer. More... | |
virtual int | getOrder () const =0 |
Get order of interpolation. More... | |
Related Functions | |
(Note that these are not member functions.) | |
template<class Scalar > | |
bool | isInitialized (const StepperBase< Scalar > &stepper) |
template<class Scalar > | |
bool | isInitialized (const StepperBase< Scalar > &stepper) |
Related Functions inherited from Rythmos::InterpolationBufferBase< Scalar > | |
template<class Scalar > | |
RCP< const Thyra::VectorBase < Scalar > > | get_x (const InterpolationBufferBase< Scalar > &interpBuffer, const Scalar &t) |
Get a single point x(t) from an interpolation buffer. More... | |
template<class Scalar > | |
RCP< const Thyra::VectorBase < Scalar > > | get_xdot (const InterpolationBufferBase< Scalar > &interpBuffer, const Scalar &t) |
Get a single point xdot(t) from an interpolation buffer. More... | |
template<class Scalar > | |
void | get_x_and_x_dot (const InterpolationBufferBase< Scalar > &interpBuffer, const Scalar t, const Ptr< RCP< const Thyra::VectorBase< Scalar > > > &x, const Ptr< RCP< const Thyra::VectorBase< Scalar > > > &x_dot) |
Nonmember helper function to get x and x_dot at t. More... | |
template<class Scalar > | |
void | assertTimePointsAreSorted (const Array< Scalar > &time_vec) |
Assert that a time point vector is sorted. More... | |
template<class Scalar > | |
void | assertNoTimePointsBeforeCurrentTimeRange (const InterpolationBufferBase< Scalar > &interpBuffer, const Array< Scalar > &time_vec, const int &startingTimePointIndex=0) |
Assert that none of the time points fall before the current time range for an interpolation buffer object. More... | |
template<class Scalar > | |
void | assertNoTimePointsInsideCurrentTimeRange (const InterpolationBufferBase< Scalar > &interpBuffer, const Array< Scalar > &time_vec) |
Assert that none of the time points fall inside the current time range for an interpolation buffer object. More... | |
template<class TimeType > | |
void | selectPointsInTimeRange (const Array< TimeType > &points_in, const TimeRange< TimeType > &range, const Ptr< Array< TimeType > > &points_out) |
Select points from an Array that sit in a TimeRange. More... | |
template<class TimeType > | |
void | removePointsInTimeRange (Array< TimeType > *points_in, const TimeRange< TimeType > &range) |
Remove points from an Array that sit in a TimeRange. More... | |
template<class Scalar > | |
bool | getCurrentPoints (const InterpolationBufferBase< Scalar > &interpBuffer, const Array< Scalar > &time_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *x_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *xdot_vec, int *nextTimePointIndex) |
Get time points in the current range of an interpolation buffer object. More... | |
Additional Inherited Members | |
Public Types inherited from Rythmos::InterpolationBufferBase< Scalar > | |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | ScalarMag |
Base class for defining stepper functionality.
A stepper object is only defined to step forward in time, never backward in time. Therefore a negative step length in the context of this interface is an invalid step length.
A stepper object can have one of three states of initialization:
Uninitialized: this->getTimeRange().isVaid()==false
Has initial condition: this->getTimeRange().lower() == this->getTimeRange().upper()
Has state buffer: this->getTimeRange().lower() < this->getTimeRange().upper()
ToDo: Finish documentation!
2007/05/17: rabartl: ToDo: Consider implementing a undoStep()
function that would erase the timestep that was just taken. This type of functionality would be needed for many different types of composed algorithms like operator split methods and for the staggered correct forward sensitivity stepper implementation with error control on the sensitivity variables.
Definition at line 16 of file Rythmos_IntegrationControlStrategyBase.hpp.
|
virtual |
Return if this stepper supports cloning or not.
If returnVal==true
, then cloneStepperAlgorithm()
will clone *this
object and return an non-null RCP.
The default implementation of this function simply returns false.
Reimplemented in Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ExplicitRKStepper< Scalar >.
Definition at line 66 of file Rythmos_StepperBase_def.hpp.
|
virtual |
Clone the stepper object if supported.
Postconditions:
supportsCloning()==true
] returnVal != Teuchos::null
supportsCloning()==false
] returnVal == Teuchos::null
Cloning a stepper in this case does not imply that the full state will be copied, shallow or deep. Instead, here cloning means to just clone the stepper algorithm and it will do a showllow copy of the model if a model is set. Since the model is stateless, this should be okay. Therefore, do not assume that the state of *returnVal
is exactly the same as the state of *this
. You have been warned!
The default implementation returns Teuchos::null
which is consistent with the default implementation of supportsCloning()==false
. If this function is overridden in a derived class to support cloning, then supportsCloning()
must be overridden to return true
.
Reimplemented in Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ExplicitRKStepper< Scalar >.
Definition at line 74 of file Rythmos_StepperBase_def.hpp.
|
virtual |
Return if this stepper is an implicit stepper.
The default implemntation returns false
and therefore, by default, a stepper is considered to be an excplicit stepper.
Reimplemented in Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ThetaStepper< Scalar >.
Definition at line 52 of file Rythmos_StepperBase_def.hpp.
|
virtual |
Return if this stepper accepts a model.
While it makes sense for most stepper objects to accept a compatible model, in some extreme cases it is not well defined what this model should be or how it relates to what is presented in the interpolation buffer interface from which this interface inherits from.
The default implemntation returns true and therefore, by default, a stepper must accept the model to be intergrated.
Reimplemented in Rythmos::ForwardSensitivityStepper< Scalar >.
Definition at line 59 of file Rythmos_StepperBase_def.hpp.
|
pure virtual |
Specify the model problem to integrate.
Preconditions:
acceptsModel()==true
!is_null(model)
model->createInArgs().supports(MEB::IN_ARG_t)==true
model->createInArgs().supports(MEB::IN_ARG_x)==true
model->createOutArgs().supports(MEB::OUT_ARG_f)==true
[isImplicit()
] model->createInArgs().supports(MEB::IN_ARG_x_dot)==true
[isImplicit()
] model->createInArgs().supports(MEB::IN_ARG_alpha)==true
[isImplicit()
] model->createInArgs().supports(MEB::IN_ARG_beta)==true
[isImplicit()
] model->createOutArgs().supports(MEB::OUT_ARG_W)==true
2007/06/10: rabartl : ToDo: Create helper macros that will assert these preconditions and call these macros in every stepper subclass that implements these function. We will have one for explicit steppers and one for implicit steppers.
Postconditions:
this->getModel() == model
this->modelIsConst() == true
Implemented in Rythmos::ForwardSensitivityStepper< Scalar >, Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ExplicitTaylorPolynomialStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, and Rythmos::ImplicitRKStepper< Scalar >.
|
pure virtual |
Accept a nonconst model.
See the full details on const version of this function above.
Postconditions:
this->getModel() == model
this->modelIsConst() == false
Implemented in Rythmos::ForwardSensitivityStepper< Scalar >, Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ExplicitTaylorPolynomialStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, and Rythmos::ImplicitRKStepper< Scalar >.
|
inlinevirtual |
Return of the model is only const or can be returned as a non-const object.
Definition at line 188 of file Rythmos_StepperBase_decl.hpp.
|
pure virtual |
Get the model.
Every stepper is expected to return the model that represents problem that it is integrating, even if acceptsModel()==false
. Exposing this model is necessary in order to get at the spaces and create the InArgs
object needed to set the initial condition.
Implemented in Rythmos::ForwardSensitivityStepper< Scalar >, Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ExplicitTaylorPolynomialStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ExplicitRKStepper< Scalar >.
|
pure virtual |
Get the model nonconst.
Implemented in Rythmos::ForwardSensitivityStepper< Scalar >, Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ExplicitTaylorPolynomialStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ExplicitRKStepper< Scalar >.
|
pure virtual |
Specify initial condition and re-initialize.
Preconditions:
!is_null(this->getModel())
The default implementation throws an exception.
Preconditions:
this->getModel() != null
ToDo: Remove this default implementation and make every concrete subclass implement this! Every stepper should except an initial condition that is separate from the model object.
Implemented in Rythmos::ForwardSensitivityStepper< Scalar >, Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ExplicitTaylorPolynomialStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ExplicitRKStepper< Scalar >.
|
pure virtual |
Get the currently set initial condtion.
Preconditions:
!is_null(this->getModel())
Implemented in Rythmos::ForwardSensitivityStepper< Scalar >, Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ExplicitTaylorPolynomialStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ExplicitRKStepper< Scalar >.
|
pure virtual |
Take a step.
dt | [in] The size of the step to take. |
stepType | [in] The type of step to take. |
Preconditions:
dt > 0.0
(i.e. only forward steps are allowed) !is_null(this->getModel())
isInitialized(*this) == true
Postconditions:
[returnVal > 0.0
] this->getTimeRange()
returns the time range [tk, tk + returnVal]
where tk
is the beginning of the timestep and tk + returnVal
is the end of the time step.
[returnVal > 0.0
] this->getPoints()
will return values of x(t)
and x_dot(t)
for all t
in this->getTimeRange()
.
returnVal > 0.0
, then a step of size returnVal
was taken. If returnVal == 0.0
, then the step could not be taken for some reason. If returnVal == 0.0
, then *this
is guaranteed to be in the same state after the function returns as it was before the function was called. This allows a client to try a different step size or make other adjustments.If stepType == STEP_TYPE_VARIABLE, and returnVal == 0.0 then no variable step will succeed in its current state.
Implemented in Rythmos::ForwardSensitivityStepper< Scalar >, Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ExplicitTaylorPolynomialStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ExplicitRKStepper< Scalar >.
|
pure virtual |
Get current stepper status after a step has been taken.
This function must have a low O(1)
complexity. In other words, it should be essentially free to call this function!
Preconditions:
isInitialized(*this) == true
Implemented in Rythmos::ForwardSensitivityStepper< Scalar >, Rythmos::BackwardEulerStepper< Scalar >, Rythmos::ForwardEulerStepper< Scalar >, Rythmos::ExplicitTaylorPolynomialStepper< Scalar >, Rythmos::ThetaStepper< Scalar >, Rythmos::ImplicitBDFStepper< Scalar >, Rythmos::ImplicitRKStepper< Scalar >, and Rythmos::ExplicitRKStepper< Scalar >.
|
virtual |
Set step control data from another stepper.
This is used to guarantee that you can re-use Jacobians from one stepper with another.
The default implementation simply does nothing so be warned!
Preconditions:
isInitialized(*this) == true
Reimplemented in Rythmos::ImplicitBDFStepper< Scalar >.
Definition at line 81 of file Rythmos_StepperBase_def.hpp.
|
related |
|
related |
Definition at line 41 of file Rythmos_StepperBase_def.hpp.