Rythmos - Transient Integration for Differential Equations  Version of the Day
 All Classes Functions Variables Typedefs Pages
Rythmos_InterpolationBufferBase.hpp
1 //@HEADER
2 // ***********************************************************************
3 //
4 // Rythmos Package
5 // Copyright (2006) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // This library is free software; you can redistribute it and/or modify
11 // it under the terms of the GNU Lesser General Public License as
12 // published by the Free Software Foundation; either version 2.1 of the
13 // License, or (at your option) any later version.
14 //
15 // This library is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 // Lesser General Public License for more details.
19 //
20 // You should have received a copy of the GNU Lesser General Public
21 // License along with this library; if not, write to the Free Software
22 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
23 // USA
24 // Questions? Contact Todd S. Coffey (tscoffe@sandia.gov)
25 //
26 // ***********************************************************************
27 //@HEADER
28 
29 #ifndef Rythmos_INTERPOLATION_BUFFER_BASE_H
30 #define Rythmos_INTERPOLATION_BUFFER_BASE_H
31 
32 #include "Rythmos_Types.hpp"
33 #include "Rythmos_TimeRange.hpp"
34 
35 #include "Thyra_VectorBase.hpp"
36 
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"
43 
44 
45 namespace Rythmos {
46 
47 
67 template<class Scalar>
69  : virtual public Teuchos::Describable
70  , virtual public Teuchos::ParameterListAcceptor
71  , virtual public Teuchos::VerboseObject<InterpolationBufferBase<Scalar> >
72 {
73 public:
74 
76  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType ScalarMag;
77 
87  virtual RCP<const Thyra::VectorSpaceBase<Scalar> >
88  get_x_space() const =0;
89 
125  /* 11/24/08 tscoffe: Proposed new interface for addPoints
126  * virtual void addPoints(
127  * const ArrayView<const Scalar>& time_vec,
128  * const ArrayView<const Ptr<const Thyra::VectorBase<Scalar> > >& x_vec,
129  * const ArrayView<const Ptr<const Thyra::VectorBase<Scalar> > >& xdot_vec
130  ) = 0;
131  */
132 
133  virtual void addPoints(
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
137  ) = 0;
138 
139  /*
140  virtual void addPoints(
141  const ArrayView<const Scalar>& time_vec,
142  const ArrayView<const Ptr<const VectorBase<Scalar> > >& x_vec,
143  const ArrayView<const Ptr<const VectorBase<Scalar> > >& xdot_vec
144  ) =0;
145  */
146 
147 
156  virtual TimeRange<Scalar> getTimeRange() const = 0;
157 
184  /* 11/24/08 tscoffe: Proposed new interface for getPoints
185  * virtual void getPoints(
186  * const ArrayView<const Scalar>& time_vec,
187  * const ArrayView<const Ptr<Thyra::VectorBase<Scalar> > >& x_vec,
188  * const ArrayView<const Ptr<Thyra::VectorBase<Scalar> > >& xdot_vec,
189  * const ArrayView<ScalarMag>& accuracy_vec
190  * ) const = 0;
191  */
192  virtual void getPoints(
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
197  ) const = 0;
198 
199  /*
200  virtual void getPoints(
201  const ArrayView<const Scalar>& time_vec,
202  const ArrayView<const Ptr<VectorBase<Scalar> > >& x_vec,
203  const ArrayView<const Ptr<VectorBase<Scalar> > >& xdot_vec,
204  const ArrayView<ScalarMag>& accuracy_vec
205  ) const = 0;
206  */
207 
217  // 11/24/08 tscoffe: Proposal: get rid of "getNodes" in abstract base interface
218  virtual void getNodes(Array<Scalar>* time_vec) const = 0;
219 
233  // 11/24/08 tscoffe: Proposal: get rid of "removeNodes" in abstract base interface
234  virtual void removeNodes(Array<Scalar>& time_vec) =0;
235 
236 
242  virtual int getOrder() const = 0;
243 
244 };
245 
246 
251 template<class Scalar>
252 RCP<const Thyra::VectorBase<Scalar> >
253 get_x( const InterpolationBufferBase<Scalar> &interpBuffer, const Scalar &t )
254 {
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()) );
261  return x_vec[0];
262 }
263 
264 /* 11/24/08 tscoffe: Proposed new get_x function:
265  * template<class Scalar>
266  * void get_x( const InterpolationBufferBase<Scalar> &interpBuffer, const Scalar &t, const Ptr<Thyra::VectorBase<Scalar> >& x_out )
267  * This will copy data into your vector and it won't be responsible for allocating new memory
268  */
269 
270 
275 template<class Scalar>
276 RCP<const Thyra::VectorBase<Scalar> >
277 get_xdot( const InterpolationBufferBase<Scalar> &interpBuffer, const Scalar &t )
278 {
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()) );
285  return xdot_vec[0];
286 }
287 
288 /* 11/24/08 tscoffe: Proposed new get_xdot function
289  * template<class Scalar>
290  * void get_xdot( const InterpolationBufferBase<Scalar> &interpBuffer, const Scalar &t, const Ptr<Thyra::VectorBase<Scalar> >& xdot_out )
291  * This will copy data into your vector and it won't be responsible for allocating new memory
292  */
293 
298 template<class Scalar>
300  const InterpolationBufferBase<Scalar> &interpBuffer,
301  const Scalar t,
302  const Ptr<RCP<const Thyra::VectorBase<Scalar> > > &x,
303  const Ptr<RCP<const Thyra::VectorBase<Scalar> > > &x_dot
304  )
305 {
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;
310  interpBuffer.getPoints(
311  time_vec,
312  nonnull(x) ? &x_vec : 0,
313  nonnull(x_dot) ? &x_dot_vec : 0,
314  0
315  );
316  if (nonnull(x)) *x = x_vec[0];
317  if (nonnull(x_dot)) *x_dot = x_dot_vec[0];
318 }
319 
321 template<class Scalar>
322 void get_x_and_x_dot(
323  const InterpolationBufferBase<Scalar> &interpBuffer,
324  const Scalar t,
325  RCP<const Thyra::VectorBase<Scalar> > *x,
326  RCP<const Thyra::VectorBase<Scalar> > *x_dot
327  )
328 {
329  get_x_and_x_dot(interpBuffer, t, Teuchos::ptr(x), Teuchos::ptr(x_dot));
330 }
331 
332 /* 11/24/08 tscoffe: Proposed new get_x_and_xdot function:
333  * template<class Scalar>
334  * void get_x_and_xdot( const InterpolationBufferBase<Scalar> &interpBuffer,
335  * const Scalar &t,
336  * const Ptr<Thyra::VectorBase<Scalar> >& x_out,
337  * const Ptr<Thyra::VectorBase<Scalar> >& xdot_out )
338  * This will copy data into your vector and it won't be responsible for allocating new memory
339  */
340 
341 } // namespace Rythmos
342 
343 
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