Rythmos - Transient Integration for Differential Equations
Version of the Day
|
Abstract interface for time integrators. More...
#include <Rythmos_IntegratorBase.hpp>
Public Types | |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | ScalarMag |
Public Types inherited from Rythmos::InterpolationBufferBase< Scalar > | |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | ScalarMag |
Public Member Functions | |
virtual RCP< IntegratorBase < Scalar > > | cloneIntegrator () const |
virtual void | setStepper (const RCP< StepperBase< Scalar > > &stepper, const Scalar &finalTime, const bool landOnFinalTime=true)=0 |
Specify the stepper to use for integration which effectively reinitializes the intergrator. More... | |
virtual Teuchos::RCP< const StepperBase< Scalar > > | getStepper () const =0 |
Get the current stepper that is set. More... | |
virtual Teuchos::RCP < StepperBase< Scalar > > | getNonconstStepper () const =0 |
Get the current stepper that is set. More... | |
virtual RCP< StepperBase < Scalar > > | unSetStepper ()=0 |
Remove the stepper and set *this to an unitilaized state. More... | |
virtual void | getFwdPoints (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)=0 |
Get values at time points both inside and outside (forward) of current TimeRange. More... | |
virtual TimeRange< Scalar > | getFwdTimeRange () const =0 |
Return the valid range of points that the integrator can integrate over. 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 > | |
RCP< const Thyra::VectorBase < Scalar > > | get_fwd_x (IntegratorBase< Scalar > &integrator, const Scalar t) |
Nonmember helper function to get x at a (forward) time t. More... | |
template<class Scalar > | |
void | get_fwd_x_and_x_dot (IntegratorBase< Scalar > &integrator, 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/or x_dot at s (forward) time t. More... | |
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... | |
Abstract interface for time integrators.
A time integrator accepts a fully initialized stepper object (and a final time) and then carries out the time integration in some fasion. The client drives the integrator by requesting value of the state at different points in time. If possible, the client should request time points only forward in time if possible since.
Definition at line 62 of file Rythmos_IntegratorBase.hpp.
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType Rythmos::IntegratorBase< Scalar >::ScalarMag |
Definition at line 67 of file Rythmos_IntegratorBase.hpp.
|
inlinevirtual |
Reimplemented in Rythmos::DefaultIntegrator< Scalar >.
Definition at line 70 of file Rythmos_IntegratorBase.hpp.
|
pure virtual |
Specify the stepper to use for integration which effectively reinitializes the intergrator.
stepper | [inout,persisting] Gives the stepper that will be used to advance the time solution. Note that it is expected that the stepper will be fully initialized and ready to start taking time steps. This includes having an initial condition which defines the initial time. |
finalTime | [in] Gives the final time that the integrator will allow itself to integrate too. |
landOnFinalTime | [in] If true , then the integrator should stop exactly (within roundoff) on finalTime . If false , then the integrator can step over the final time if the stepper desires it. |
Preconditions:
!is_null(stepper)
stepper->getTimeRange().size() >= 0.0
compareTimeValues(finalTime,stepper->getTimeRange().upper()) >= 0
compareTimeValues(getTrailingInterpBuffer()->getTimeRange().upper(),stepper->getTimeRange().lower())==0
. Postconditions:
this->getStepper() == stepper
compareTimeValues(this->getFwdTimeRange().lower(),stepper->getTimeRange().lower())==0
compareTimeValues(this->getFwdTimeRange().upper(),finalTime)==0
Implemented in Rythmos::DefaultIntegrator< Scalar >.
|
pure virtual |
Get the current stepper that is set.
returnVal==null
which case *this
is in an uninitialized state. Implemented in Rythmos::DefaultIntegrator< Scalar >.
|
pure virtual |
Get the current stepper that is set.
returnVal==null
which case *this
is in an uninitialized state. Implemented in Rythmos::DefaultIntegrator< Scalar >.
|
pure virtual |
Remove the stepper and set *this
to an unitilaized state.
Postconditions:
is_null(this->getStepper()) == true
this->getTimeRange().isValid() == false
Implemented in Rythmos::DefaultIntegrator< Scalar >.
|
pure virtual |
Get values at time points both inside and outside (forward) of current TimeRange.
time_vec | [in] Array (length n ) of time points to get. |
x_vec | [out] On output, if x_vec != 0 , *x_vec will be resized to n = time_vec.size() and (*x_vec)[i] will be the state vector at time time_vec[i] , for i=0...n-1 . This argument can be left NULL in which case it will not be filled. |
xdot_vec | [out] On output, if xdot_vec != 0 , *xdot_vec will be resized to n = time_vec.size() and (*xdot_vec)[i] will be the state derivative vector at time time_vec[i] , for i=0...n-1 . This argument can be left NULL in which case it will not be filled. |
accuracy_vec | [out] This contains an estimate of the accuracy of the interpolation. This argument can be left NULL in which case it will not be filled. If you asked for a node, this should be zero. |
Preconditions:
range.lower() <= time_vec[i] <= range.upper()
, for i=0...n-1
, where range = this->getFwdTimeRange()
. time_vec
must have unique and sorted values in ascending order Postconditions:
this->getTimeRange().lower()
may be greater after this function returns than before this function was called! That is why this is a non-const function! Throwns | Exceptions::GetFwdPointsFailed if all of the time points could not be reached for some reason (e.g. the max number of time-step iterations was exceeded). |
This is a non-const version of the const function getPoints()
which allows the integrator class to step forward to get the points asked for.
|
pure virtual |
Return the valid range of points that the integrator can integrate over.
Postconditions:
this->getFwdTimeRange().lower() == this->getTimeRange().lower()
this->getFwdTimeRange().upper() >= this->getTimeRange().upper()
Implemented in Rythmos::DefaultIntegrator< Scalar >.
|
related |
Nonmember helper function to get x at a (forward) time t.
Definition at line 205 of file Rythmos_IntegratorBase.hpp.
|
related |
Nonmember helper function to get x and/or x_dot at s (forward) time t.
integrator | [modified] This is the integrator that the x and x_dot will be gotten from. This is declared non-const since the state of the integrator may change when getting forward points in time. |
t | [in] The time point to get x and x_dot at. |
x | [out] RCP to state vector to get if x!=0 . |
x_dot | [out] RCP to state derivative vector to get if x_dot!=0 . |
Note that *x
(if x!=0
) and *x_dot
(if x_dot!=0
) should be null on input since the RCP will get set internally. Note that on output *x
and *x_dot
is a const vector so the client can not change them. The client should destroy (or set to null) the RCP objects as soon as it is finished reading the vectors. The input values of x==0
and/or x_dot==0
are valid which means that these vectors will not be returned.
Definition at line 240 of file Rythmos_IntegratorBase.hpp.