29 #ifndef Rythmos_INTERPOLATION_BUFFER_BASE_H
30 #define Rythmos_INTERPOLATION_BUFFER_BASE_H
32 #include "Rythmos_Types.hpp"
33 #include "Rythmos_TimeRange.hpp"
35 #include "Thyra_VectorBase.hpp"
37 #include "Teuchos_Describable.hpp"
38 #include "Teuchos_ParameterListAcceptor.hpp"
39 #include "Teuchos_VerboseObject.hpp"
40 #include "Teuchos_implicit_cast.hpp"
41 #include "Teuchos_Assert.hpp"
42 #include "Teuchos_as.hpp"
67 template<
class Scalar>
69 :
virtual public Teuchos::Describable
70 ,
virtual public Teuchos::ParameterListAcceptor
71 ,
virtual public Teuchos::VerboseObject<InterpolationBufferBase<Scalar> >
76 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
ScalarMag;
87 virtual RCP<const Thyra::VectorSpaceBase<Scalar> >
134 const Array<Scalar>& time_vec,
135 const Array<RCP<
const Thyra::VectorBase<Scalar> > >& x_vec,
136 const Array<RCP<
const Thyra::VectorBase<Scalar> > >& xdot_vec
193 const Array<Scalar>& time_vec,
194 Array<RCP<
const Thyra::VectorBase<Scalar> > >* x_vec,
195 Array<RCP<
const Thyra::VectorBase<Scalar> > >* xdot_vec,
196 Array<ScalarMag>* accuracy_vec
218 virtual void getNodes(Array<Scalar>* time_vec)
const = 0;
234 virtual void removeNodes(Array<Scalar>& time_vec) =0;
251 template<
class Scalar>
252 RCP<const Thyra::VectorBase<Scalar> >
255 using Teuchos::implicit_cast;
256 Array<Scalar> time_vec;
257 time_vec.push_back(t);
258 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
259 interpBuffer.
getPoints(time_vec,&x_vec,0,0);
260 TEUCHOS_ASSERT( 1 == implicit_cast<int>(x_vec.size()) );
275 template<
class Scalar>
276 RCP<const Thyra::VectorBase<Scalar> >
279 using Teuchos::implicit_cast;
280 Array<Scalar> time_vec;
281 time_vec.push_back(t);
282 Array<RCP<const Thyra::VectorBase<Scalar> > > xdot_vec;
283 interpBuffer.
getPoints(time_vec,0,&xdot_vec,0);
284 TEUCHOS_ASSERT( 1 == implicit_cast<int>(xdot_vec.size()) );
298 template<
class Scalar>
302 const Ptr<RCP<
const Thyra::VectorBase<Scalar> > > &x,
303 const Ptr<RCP<
const Thyra::VectorBase<Scalar> > > &x_dot
306 Array<Scalar> time_vec;
307 time_vec.push_back(t);
308 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
309 Array<RCP<const Thyra::VectorBase<Scalar> > > x_dot_vec;
312 nonnull(x) ? &x_vec : 0,
313 nonnull(x_dot) ? &x_dot_vec : 0,
316 if (nonnull(x)) *x = x_vec[0];
317 if (nonnull(x_dot)) *x_dot = x_dot_vec[0];
321 template<
class Scalar>
322 void get_x_and_x_dot(
325 RCP<
const Thyra::VectorBase<Scalar> > *x,
326 RCP<
const Thyra::VectorBase<Scalar> > *x_dot
329 get_x_and_x_dot(interpBuffer, t, Teuchos::ptr(x), Teuchos::ptr(x_dot));
344 #endif //Rythmos_INTERPOLATION_BUFFER_BASE_H
RCP< const Thyra::VectorBase< Scalar > > get_xdot(const InterpolationBufferBase< Scalar > &interpBuffer, const Scalar &t)
Get a single point xdot(t) from an interpolation buffer.
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.
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.
RCP< const Thyra::VectorBase< Scalar > > get_x(const InterpolationBufferBase< Scalar > &interpBuffer, const Scalar &t)
Get a single point x(t) from an interpolation buffer.
virtual void getNodes(Array< Scalar > *time_vec) const =0
Get interpolation nodes.
virtual void removeNodes(Array< Scalar > &time_vec)=0
Remove nodes from the interpolation buffer.
Base class for an interpolation buffer.
virtual RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const =0
Return the space for x and x_dot.
virtual TimeRange< Scalar > getTimeRange() const =0
Return the range of time values where interpolation calls can be performed.
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.
virtual int getOrder() const =0
Get order of interpolation.
Teuchos::ScalarTraits< Scalar >::magnitudeType ScalarMag