Rythmos_InterpolationBufferBase.hpp
00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef Rythmos_INTERPOLATION_BUFFER_BASE_H
00030 #define Rythmos_INTERPOLATION_BUFFER_BASE_H
00031
00032 #include "Rythmos_Types.hpp"
00033 #include "Rythmos_TimeRange.hpp"
00034
00035 #include "Thyra_VectorBase.hpp"
00036
00037 #include "Teuchos_Describable.hpp"
00038 #include "Teuchos_ParameterListAcceptor.hpp"
00039 #include "Teuchos_VerboseObject.hpp"
00040 #include "Teuchos_implicit_cast.hpp"
00041 #include "Teuchos_Assert.hpp"
00042 #include "Teuchos_as.hpp"
00043
00044
00045 namespace Rythmos {
00046
00047
00067 template<class Scalar>
00068 class InterpolationBufferBase
00069 : virtual public Teuchos::Describable
00070 , virtual public Teuchos::ParameterListAcceptor
00071 , virtual public Teuchos::VerboseObject<InterpolationBufferBase<Scalar> >
00072 {
00073 public:
00074
00076 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType ScalarMag;
00077
00087 virtual RCP<const Thyra::VectorSpaceBase<Scalar> >
00088 get_x_space() const =0;
00089
00125 virtual void addPoints(
00126 const Array<Scalar>& time_vec,
00127 const Array<RCP<const Thyra::VectorBase<Scalar> > >& x_vec,
00128 const Array<RCP<const Thyra::VectorBase<Scalar> > >& xdot_vec
00129 ) = 0;
00130
00131
00140 virtual TimeRange<Scalar> getTimeRange() const = 0;
00141
00168 virtual void getPoints(
00169 const Array<Scalar>& time_vec,
00170 Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00171 Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00172 Array<ScalarMag>* accuracy_vec
00173 ) const = 0;
00174
00184 virtual void getNodes(Array<Scalar>* time_vec) const = 0;
00185
00196 virtual void removeNodes(Array<Scalar>& time_vec) =0;
00197
00198
00204 virtual int getOrder() const = 0;
00205
00206 };
00207
00208
00213 template<class Scalar>
00214 RCP<const Thyra::VectorBase<Scalar> >
00215 get_x( const InterpolationBufferBase<Scalar> &interpBuffer, const Scalar &t )
00216 {
00217 using Teuchos::implicit_cast;
00218 Array<Scalar> time_vec;
00219 time_vec.push_back(t);
00220 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
00221 interpBuffer.getPoints(time_vec,&x_vec,0,0);
00222 TEUCHOS_ASSERT( 1 == implicit_cast<int>(x_vec.size()) );
00223 return x_vec[0];
00224 }
00225
00226
00231 template<class Scalar>
00232 RCP<const Thyra::VectorBase<Scalar> >
00233 get_xdot( const InterpolationBufferBase<Scalar> &interpBuffer, const Scalar &t )
00234 {
00235 using Teuchos::implicit_cast;
00236 Array<Scalar> time_vec;
00237 time_vec.push_back(t);
00238 Array<RCP<const Thyra::VectorBase<Scalar> > > xdot_vec;
00239 interpBuffer.getPoints(time_vec,0,&xdot_vec,0);
00240 TEUCHOS_ASSERT( 1 == implicit_cast<int>(xdot_vec.size()) );
00241 return xdot_vec[0];
00242 }
00243
00244
00249 template<class Scalar>
00250 void get_x_and_x_dot(
00251 const InterpolationBufferBase<Scalar> &interpBuffer,
00252 const Scalar t,
00253 RCP<const Thyra::VectorBase<Scalar> > *x,
00254 RCP<const Thyra::VectorBase<Scalar> > *x_dot
00255 )
00256 {
00257 Array<Scalar> time_vec;
00258 time_vec.push_back(t);
00259 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
00260 Array<RCP<const Thyra::VectorBase<Scalar> > > x_dot_vec;
00261 interpBuffer.getPoints(
00262 time_vec,
00263 x ? &x_vec : 0,
00264 x_dot ? &x_dot_vec : 0,
00265 0
00266 );
00267 if (x) *x = x_vec[0];
00268 if (x_dot) *x_dot = x_dot_vec[0];
00269 }
00270
00271
00272 }
00273
00274
00275 #endif //Rythmos_INTERPOLATION_BUFFER_BASE_H