Rythmos - Transient Integration for Differential Equations  Version of the Day
 All Classes Functions Variables Typedefs Pages
Public Types | Public Member Functions | Related Functions | List of all members
Rythmos::IntegratorBase< Scalar > Class Template Referenceabstract

Abstract interface for time integrators. More...

#include <Rythmos_IntegratorBase.hpp>

Inheritance diagram for Rythmos::IntegratorBase< Scalar >:
Inheritance graph
[legend]

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

Detailed Description

template<class Scalar>
class Rythmos::IntegratorBase< Scalar >

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.

Member Typedef Documentation

template<class Scalar>
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType Rythmos::IntegratorBase< Scalar >::ScalarMag

Definition at line 67 of file Rythmos_IntegratorBase.hpp.

Member Function Documentation

template<class Scalar>
virtual RCP<IntegratorBase<Scalar> > Rythmos::IntegratorBase< Scalar >::cloneIntegrator ( ) const
inlinevirtual

Reimplemented in Rythmos::DefaultIntegrator< Scalar >.

Definition at line 70 of file Rythmos_IntegratorBase.hpp.

template<class Scalar>
virtual void Rythmos::IntegratorBase< Scalar >::setStepper ( const RCP< StepperBase< Scalar > > &  stepper,
const Scalar &  finalTime,
const bool  landOnFinalTime = true 
)
pure virtual

Specify the stepper to use for integration which effectively reinitializes the intergrator.

Parameters
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:

Postconditions:

Implemented in Rythmos::DefaultIntegrator< Scalar >.

template<class Scalar>
virtual Teuchos::RCP<const StepperBase<Scalar> > Rythmos::IntegratorBase< Scalar >::getStepper ( ) const
pure virtual

Get the current stepper that is set.

Returns
This function can return returnVal==null which case *this is in an uninitialized state.

Implemented in Rythmos::DefaultIntegrator< Scalar >.

template<class Scalar>
virtual Teuchos::RCP<StepperBase<Scalar> > Rythmos::IntegratorBase< Scalar >::getNonconstStepper ( ) const
pure virtual

Get the current stepper that is set.

Returns
This function can return returnVal==null which case *this is in an uninitialized state.

Implemented in Rythmos::DefaultIntegrator< Scalar >.

template<class Scalar>
virtual RCP<StepperBase<Scalar> > Rythmos::IntegratorBase< Scalar >::unSetStepper ( )
pure virtual

Remove the stepper and set *this to an unitilaized state.

Postconditions:

Implemented in Rythmos::DefaultIntegrator< Scalar >.

template<class Scalar>
virtual void Rythmos::IntegratorBase< Scalar >::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 
)
pure virtual

Get values at time points both inside and outside (forward) of current TimeRange.

Parameters
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:

  • Returns all of requested time points if no exception is thrown.
  • Note that 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!
Exceptions
ThrownsExceptions::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.

template<class Scalar>
virtual TimeRange<Scalar> Rythmos::IntegratorBase< Scalar >::getFwdTimeRange ( ) const
pure virtual

Return the valid range of points that the integrator can integrate over.

Postconditions:

Implemented in Rythmos::DefaultIntegrator< Scalar >.

Friends And Related Function Documentation

template<class Scalar >
RCP< const Thyra::VectorBase< Scalar > > get_fwd_x ( IntegratorBase< Scalar > &  integrator,
const Scalar  t 
)
related

Nonmember helper function to get x at a (forward) time t.

Definition at line 205 of file Rythmos_IntegratorBase.hpp.

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 
)
related

Nonmember helper function to get x and/or x_dot at s (forward) time t.

Parameters
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.


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