Rythmos - Transient Integration for Differential Equations
Version of the Day
|
#include <Rythmos_HermiteInterpolator_decl.hpp>
Public Member Functions | |
~HermiteInterpolator () | |
Destructor. More... | |
HermiteInterpolator () | |
Constructor. More... | |
void | setNodes (const RCP< const typename DataStore< Scalar >::DataStoreVector_t > &nodes) |
Store pointer to interpolation nodes. More... | |
void | interpolate (const Array< Scalar > &t_values, typename DataStore< Scalar >::DataStoreVector_t *data_out) const |
Interpolation: More... | |
int | order () const |
Order of interpolation: More... | |
std::string | description () const |
Inherited from Teuchos::Describable. More... | |
void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
void | setParameterList (RCP< ParameterList > const ¶mList) |
Redefined from ParameterListAcceptor. More... | |
RCP< ParameterList > | getNonconstParameterList () |
RCP< ParameterList > | unsetParameterList () |
RCP< const Teuchos::ParameterList > | getValidParameters () const |
Public Member Functions inherited from Rythmos::InterpolatorBase< Scalar > | |
virtual bool | supportsCloning () const |
Return if this interpolator supports cloning or not. More... | |
virtual Teuchos::RCP < InterpolatorBase< Scalar > > | cloneInterpolator () const |
Clone the interpolator object if supported. More... | |
Additional Inherited Members | |
Related Functions inherited from Rythmos::InterpolatorBase< Scalar > | |
template<class Scalar > | |
void | assertBaseInterpolatePreconditions (const typename DataStore< Scalar >::DataStoreVector_t &data_in, const Array< Scalar > &t_values, typename DataStore< Scalar >::DataStoreVector_t *data_out) |
This class implements piecewise Hermite interpolation on each interval where the data is: (t0,x(t0)), (t1,x(t1)), (t0,x'(t0)), (t1,x'(t1)) The Hermite Interpolation polynomial is: H_3(t) = x[z0] + x[z0,z1](t-t0) + x[z0,z1,z2](t-t0)^2 + x[z0,z1,z2,z3](t-t0)^2(t-t1) where z0 = z1 = t0 and z2 = z3 = t1 and x[z0,z1] = x'(t0) and x[z2,z3] = x'(t1) This reduces to: H_3(t) = x(t0) + x'(t0)(t-t0) + ((x(t1)-x(t0))/(t1-t0) - x'(t0))(t-t0)^2/(t1-t0) +(x'(t1) - 2(x(t1)-x(t0))/(t1-t0) + x'(t0))(t-t0)^2(t-t1)/(t1-t0)^2 With derivative: H_3'(t) = x'(t0) + 2*((x(t1)-x(t0))/(t1-t0) - x'(t0))(t-t0)/(t1-t0) +(x'(t1) - 2(x(t1)-x(t0))/(t1-t0) + x'(t0))[2*(t-t0)(t-t1) + (t-t0)^2]/(t1-t0)^2 With the error expression: x(t) - H_3(t) = (x^{(3)}((t))/(4!))(t-t0)^2(t-t1)^2 The Hermite Interpolant will match 3rd degree polynomials exactly with both function values and derivatives.
Definition at line 54 of file Rythmos_HermiteInterpolator_decl.hpp.
|
inline |
Destructor.
Definition at line 59 of file Rythmos_HermiteInterpolator_decl.hpp.
Rythmos::HermiteInterpolator< Scalar >::HermiteInterpolator | ( | ) |
Constructor.
Definition at line 40 of file Rythmos_HermiteInterpolator_def.hpp.
|
virtual |
Store pointer to interpolation nodes.
This function represent a persisting relationship between the interpolation nodes and the interpolator. For simple interpolators like linear and Hermite, this is not needed, but for interpolators like cubic splines where there is some computational work in assembling the interpolant, this is important.
Preconditions:
nodes
must have unique time values and be sorted in ascending time order Postconditions:
nodes
is a different pointer than was previously called, then it is possible that the interpolant will be recomputed when interpolate
is next called. Implements Rythmos::InterpolatorBase< Scalar >.
Definition at line 47 of file Rythmos_HermiteInterpolator_def.hpp.
|
virtual |
Interpolation:
Hermite interpolation function.
Preconditions:
(*nodes_)[i].xdot != Teuchos::null
for all i=0..(*nodes_).size()-1
Implements Rythmos::InterpolatorBase< Scalar >.
Definition at line 55 of file Rythmos_HermiteInterpolator_def.hpp.
|
virtual |
Order of interpolation:
Implements Rythmos::InterpolatorBase< Scalar >.
Definition at line 181 of file Rythmos_HermiteInterpolator_def.hpp.
std::string Rythmos::HermiteInterpolator< Scalar >::description | ( | ) | const |
Inherited from Teuchos::Describable.
Definition at line 187 of file Rythmos_HermiteInterpolator_def.hpp.
void Rythmos::HermiteInterpolator< Scalar >::describe | ( | Teuchos::FancyOStream & | out, |
const Teuchos::EVerbosityLevel | verbLevel | ||
) | const |
Definition at line 194 of file Rythmos_HermiteInterpolator_def.hpp.
void Rythmos::HermiteInterpolator< Scalar >::setParameterList | ( | RCP< ParameterList > const & | paramList | ) |
Redefined from ParameterListAcceptor.
Definition at line 214 of file Rythmos_HermiteInterpolator_def.hpp.
RCP< ParameterList > Rythmos::HermiteInterpolator< Scalar >::getNonconstParameterList | ( | ) |
Definition at line 223 of file Rythmos_HermiteInterpolator_def.hpp.
RCP< ParameterList > Rythmos::HermiteInterpolator< Scalar >::unsetParameterList | ( | ) |
Definition at line 229 of file Rythmos_HermiteInterpolator_def.hpp.
RCP< const Teuchos::ParameterList > Rythmos::HermiteInterpolator< Scalar >::getValidParameters | ( | ) | const |
Definition at line 237 of file Rythmos_HermiteInterpolator_def.hpp.