#include <Rythmos_IntegratorBase.hpp>
Inheritance diagram for Rythmos::IntegratorBase< Scalar >:

Public Types | |
| typedef Teuchos::ScalarTraits< Scalar >::magnitudeType | ScalarMag |
| | |
Public Member Functions | |
| virtual RCP< IntegratorBase< Scalar > > | cloneIntegrator () const |
| | |
| virtual void | setTrailingInterpolationBuffer (const RCP< InterpolationBufferBase< Scalar > > &trailingInterpBuffer) |
| | |
| virtual bool | acceptsTrailingInterpolationBuffer () 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. | |
| virtual Teuchos::RCP< const StepperBase< Scalar > > | getStepper () const =0 |
| Get the current stepper that is set. | |
| virtual RCP< StepperBase< Scalar > > | unSetStepper ()=0 |
Remove the stepper and set *this to an unitilaized state. | |
| 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. | |
| virtual TimeRange< Scalar > | getFwdTimeRange () const =0 |
| Return the valid range of points that the integrator can integrate over. | |
Related Functions | |
| (Note that these are not member functions.) | |
| 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. | |
| void | get_fwd_x_and_x_dot (IntegratorBase< Scalar > &integrator, const Scalar t, RCP< const Thyra::VectorBase< Scalar > > *x, RCP< const Thyra::VectorBase< Scalar > > *x_dot) |
| Nonmember helper function to get x and/or x_dot at s (forward) time t. | |
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 |
Reimplemented from Rythmos::InterpolationBufferBase< Scalar >.
Reimplemented in Rythmos::IntegratorDefault< Scalar >, and Rythmos::SimpleIntegrator< Scalar >.
Definition at line 67 of file Rythmos_IntegratorBase.hpp.
| virtual RCP<IntegratorBase<Scalar> > Rythmos::IntegratorBase< Scalar >::cloneIntegrator | ( | ) | const [inline, virtual] |
Reimplemented in Rythmos::SimpleIntegrator< Scalar >.
Definition at line 70 of file Rythmos_IntegratorBase.hpp.
| virtual void Rythmos::IntegratorBase< Scalar >::setTrailingInterpolationBuffer | ( | const RCP< InterpolationBufferBase< Scalar > > & | trailingInterpBuffer | ) | [inline, virtual] |
Reimplemented in Rythmos::IntegratorDefault< Scalar >.
Definition at line 74 of file Rythmos_IntegratorBase.hpp.
| virtual bool Rythmos::IntegratorBase< Scalar >::acceptsTrailingInterpolationBuffer | ( | ) | const [inline, virtual] |
Reimplemented in Rythmos::IntegratorDefault< Scalar >.
Definition at line 80 of file Rythmos_IntegratorBase.hpp.
| 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.
| stepper | [inout,persisting] Gives the stepper that will be used to advance the time solution. | |
| 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. |
!is_null(stepper) stepper->getTimeRange().size() >= 0.0 compareTimeValues(finalTime,stepper->getTimeRange().upper()) >= 0 compareTimeValues(trailingInterpBuffer_->getTimeRange().upper(),stepper->getTimeRange().lower())==0. Postconditions:
this->getStepper() == stepper compareTimeValues(this->getFwdTimeRange().lower(),stepper->getTimeRange().lower())==0 compareTimeValues(this->getFwdTimeRange().upper(),finalTime)==0 2007/08/24: rabartl: ToDo: We should add another enum argument that specifies if we should let the stepper step past finalTime or if it has to stop exactly (within some floating point error) on top of finalTime at the very end. In essense, we should pass in if finalTime is a soft or hard breakpoint.
Implemented in Rythmos::IntegratorDefault< Scalar >, and Rythmos::SimpleIntegrator< Scalar >.
| virtual Teuchos::RCP<const StepperBase<Scalar> > Rythmos::IntegratorBase< Scalar >::getStepper | ( | ) | const [pure virtual] |
Get the current stepper that is set.
returnVal==null which case *this is in an uninitialized state. Implemented in Rythmos::IntegratorDefault< Scalar >, and Rythmos::SimpleIntegrator< Scalar >.
| virtual RCP<StepperBase<Scalar> > Rythmos::IntegratorBase< Scalar >::unSetStepper | ( | ) | [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::IntegratorDefault< Scalar >, and Rythmos::SimpleIntegrator< 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.
| 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. |
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). |
getPoints() which allows the integrator class to step forward to get the points asked for.
Implemented in Rythmos::IntegratorDefault< Scalar >.
| virtual TimeRange<Scalar> Rythmos::IntegratorBase< Scalar >::getFwdTimeRange | ( | ) | const [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::IntegratorDefault< Scalar >, and Rythmos::SimpleIntegrator< 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 212 of file Rythmos_IntegratorBase.hpp.
| void get_fwd_x_and_x_dot | ( | IntegratorBase< Scalar > & | integrator, | |
| const Scalar | t, | |||
| RCP< const Thyra::VectorBase< Scalar > > * | x, | |||
| RCP< const Thyra::VectorBase< Scalar > > * | x_dot | |||
| ) | [related] |
Nonmember helper function to get x and/or x_dot at s (forward) time t.
Definition at line 228 of file Rythmos_IntegratorBase.hpp.
1.4.7