Rythmos - Transient Integration for Differential Equations
Version of the Day
|
Implementation of Rythmos::Stepper for explicit Taylor polynomial time integration of ODEs. More...
#include <Rythmos_ExplicitTaylorPolynomialStepper.hpp>
Public Types | |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | ScalarMag |
Typename of magnitude of scalars. More... | |
Public Types inherited from Rythmos::InterpolationBufferBase< Scalar > | |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | ScalarMag |
Public Member Functions | |
ExplicitTaylorPolynomialStepper () | |
Constructor. More... | |
~ExplicitTaylorPolynomialStepper () | |
Destructor. More... | |
RCP< const Thyra::VectorSpaceBase< Scalar > > | get_x_space () const |
Return the space for x and x_dot More... | |
void | setModel (const RCP< const Thyra::ModelEvaluator< Scalar > > &model) |
Set model. More... | |
void | setNonconstModel (const RCP< Thyra::ModelEvaluator< Scalar > > &model) |
Set model. More... | |
RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () const |
RCP< Thyra::ModelEvaluator < Scalar > > | getNonconstModel () |
void | setInitialCondition (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &initialCondition) |
Thyra::ModelEvaluatorBase::InArgs < Scalar > | getInitialCondition () const |
Scalar | takeStep (Scalar dt, StepSizeType flag) |
Take a time step of magnitude dt . More... | |
const StepStatus< Scalar > | getStepStatus () const |
void | setParameterList (RCP< Teuchos::ParameterList > const ¶mList) |
Redefined from Teuchos::ParameterListAcceptor. More... | |
RCP< Teuchos::ParameterList > | getNonconstParameterList () |
RCP< Teuchos::ParameterList > | unsetParameterList () |
RCP< const Teuchos::ParameterList > | getValidParameters () const |
std::string | description () const |
void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
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) |
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 |
Get values from buffer. More... | |
void | setRange (const TimeRange< Scalar > &range, const InterpolationBufferBase< Scalar > &IB) |
Fill data in from another interpolation buffer. More... | |
TimeRange< Scalar > | getTimeRange () const |
void | getNodes (Array< Scalar > *time_vec) const |
Get interpolation nodes. More... | |
void | removeNodes (Array< Scalar > &time_vec) |
Remove interpolation nodes. More... | |
int | getOrder () const |
Get order of interpolation. More... | |
Public Member Functions inherited from Rythmos::StepperBase< Scalar > | |
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 bool | modelIsConst () const |
Return of the model is only const or can be returned as a non-const object. More... | |
virtual void | setStepControlData (const StepperBase &stepper) |
Set step control data from another stepper. More... | |
Additional Inherited Members | |
Related Functions inherited from Rythmos::StepperBase< Scalar > | |
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... | |
Implementation of Rythmos::Stepper for explicit Taylor polynomial time integration of ODEs.
Let
be an ODE initial-value problem. This class implements a single time step of an explicit Taylor polynomial time integration method for computing numerical solutions to the IVP. The method consists of computing a local truncated Taylor series solution to the ODE (section Computing the Taylor Polynomial), estimating a step size within the radius of convergence of the Taylor series (section Computing a Step Size) and then summing the polynomial at that step to compute the next point in the numerical integration (section Summing the Polynomial). The algorithmic parameters to the method are controlled through the params
argument to the constructor which are described in section Parameters.
Let
be a power series solution to the IVP above. Then can be expaned in a power series along the solution curve :
where
By differentiating the power series for to compute a power series for and then comparing coefficients in the equation we find the following recurrence relationship for the Taylor coefficients :
where each coefficient is a function only of and can be efficiently computed using the Taylor polynomial mode of automatic differentation. This allows the Taylor coefficients to be iteratively computed starting with the initial point up to some fixed degree to yield a local truncated Taylor polynomial approximating the solution to the IVP.
With the truncated Taylor polynomial solution in hand, a step size is chosen by estimating the truncation error in the polynomial solution and forcing this error to be less than some prescribed tolerance. Let
so for . Assume for as well, then for any it can be shown that the truncation error is bounded by
A step size is then given by
for some error tolerance given an error of approximatly .
With a step size computed,
is used as the next integration point where a new Taylor series is calculated. Local error per step can also be controlled by computing . If this error is too large, the step size can be reduced to an appropriate size.
This method recognizes the following algorithmic parameters that can be set in the params
argument to the constructor:
Definition at line 163 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::ScalarMag |
Typename of magnitude of scalars.
Definition at line 168 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::ExplicitTaylorPolynomialStepper | ( | ) |
Constructor.
Definition at line 364 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::~ExplicitTaylorPolynomialStepper | ( | ) |
Destructor.
Definition at line 372 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Return the space for x
and x_dot
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 863 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Set model.
Implements Rythmos::StepperBase< Scalar >.
Definition at line 406 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Set model.
Implements Rythmos::StepperBase< Scalar >.
Definition at line 419 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 429 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 437 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 444 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 467 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Take a time step of magnitude dt
.
Implements Rythmos::StepperBase< Scalar >.
Definition at line 475 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 593 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
void Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::setParameterList | ( | RCP< Teuchos::ParameterList > const & | paramList | ) |
Redefined from Teuchos::ParameterListAcceptor.
Definition at line 626 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
RCP< Teuchos::ParameterList > Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::getNonconstParameterList | ( | ) |
Definition at line 661 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
RCP< Teuchos::ParameterList > Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::unsetParameterList | ( | ) |
Definition at line 669 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
RCP< const Teuchos::ParameterList > Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::getValidParameters | ( | ) | const |
Definition at line 679 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
std::string Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::description | ( | ) | const |
Definition at line 701 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
void Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::describe | ( | Teuchos::FancyOStream & | out, |
const Teuchos::EVerbosityLevel | verbLevel = Teuchos::Describable::verbLevel_default |
||
) | const |
Definition at line 709 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Redefined from InterpolationBufferBase Add points to buffer
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 741 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Get values from buffer.
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 752 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
void Rythmos::ExplicitTaylorPolynomialStepper< Scalar >::setRange | ( | const TimeRange< Scalar > & | range, |
const InterpolationBufferBase< Scalar > & | IB | ||
) |
Fill data in from another interpolation buffer.
|
virtual |
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 771 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Get interpolation nodes.
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 782 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Remove interpolation nodes.
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 798 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.
|
virtual |
Get order of interpolation.
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 805 of file Rythmos_ExplicitTaylorPolynomialStepper.hpp.