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
00126
00127
00128
00129
00130
00131
00132
00133 virtual void addPoints(
00134 const Array<Scalar>& time_vec,
00135 const Array<RCP<const Thyra::VectorBase<Scalar> > >& x_vec,
00136 const Array<RCP<const Thyra::VectorBase<Scalar> > >& xdot_vec
00137 ) = 0;
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00156 virtual TimeRange<Scalar> getTimeRange() const = 0;
00157
00184
00185
00186
00187
00188
00189
00190
00191
00192 virtual void getPoints(
00193 const Array<Scalar>& time_vec,
00194 Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00195 Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00196 Array<ScalarMag>* accuracy_vec
00197 ) const = 0;
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00217
00218 virtual void getNodes(Array<Scalar>* time_vec) const = 0;
00219
00233
00234 virtual void removeNodes(Array<Scalar>& time_vec) =0;
00235
00236
00242 virtual int getOrder() const = 0;
00243
00244 };
00245
00246
00251 template<class Scalar>
00252 RCP<const Thyra::VectorBase<Scalar> >
00253 get_x( const InterpolationBufferBase<Scalar> &interpBuffer, const Scalar &t )
00254 {
00255 using Teuchos::implicit_cast;
00256 Array<Scalar> time_vec;
00257 time_vec.push_back(t);
00258 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
00259 interpBuffer.getPoints(time_vec,&x_vec,0,0);
00260 TEUCHOS_ASSERT( 1 == implicit_cast<int>(x_vec.size()) );
00261 return x_vec[0];
00262 }
00263
00264
00265
00266
00267
00268
00269
00270
00275 template<class Scalar>
00276 RCP<const Thyra::VectorBase<Scalar> >
00277 get_xdot( const InterpolationBufferBase<Scalar> &interpBuffer, const Scalar &t )
00278 {
00279 using Teuchos::implicit_cast;
00280 Array<Scalar> time_vec;
00281 time_vec.push_back(t);
00282 Array<RCP<const Thyra::VectorBase<Scalar> > > xdot_vec;
00283 interpBuffer.getPoints(time_vec,0,&xdot_vec,0);
00284 TEUCHOS_ASSERT( 1 == implicit_cast<int>(xdot_vec.size()) );
00285 return xdot_vec[0];
00286 }
00287
00288
00289
00290
00291
00292
00293
00298 template<class Scalar>
00299 void get_x_and_x_dot(
00300 const InterpolationBufferBase<Scalar> &interpBuffer,
00301 const Scalar t,
00302 const Ptr<RCP<const Thyra::VectorBase<Scalar> > > &x,
00303 const Ptr<RCP<const Thyra::VectorBase<Scalar> > > &x_dot
00304 )
00305 {
00306 Array<Scalar> time_vec;
00307 time_vec.push_back(t);
00308 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
00309 Array<RCP<const Thyra::VectorBase<Scalar> > > x_dot_vec;
00310 interpBuffer.getPoints(
00311 time_vec,
00312 nonnull(x) ? &x_vec : 0,
00313 nonnull(x_dot) ? &x_dot_vec : 0,
00314 0
00315 );
00316 if (nonnull(x)) *x = x_vec[0];
00317 if (nonnull(x_dot)) *x_dot = x_dot_vec[0];
00318 }
00319
00321 template<class Scalar>
00322 void get_x_and_x_dot(
00323 const InterpolationBufferBase<Scalar> &interpBuffer,
00324 const Scalar t,
00325 RCP<const Thyra::VectorBase<Scalar> > *x,
00326 RCP<const Thyra::VectorBase<Scalar> > *x_dot
00327 )
00328 {
00329 get_x_and_x_dot(interpBuffer, t, Teuchos::ptr(x), Teuchos::ptr(x_dot));
00330 }
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341 }
00342
00343
00344 #endif //Rythmos_INTERPOLATION_BUFFER_BASE_H