Rythmos - Transient Integration for Differential Equations  Version of the Day
 All Classes Functions Variables Typedefs Pages
Rythmos_IntegratorBase.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_INTEGRATOR_BASE_H
30 #define Rythmos_INTEGRATOR_BASE_H
31 
32 #include "Rythmos_InterpolationBufferBase.hpp"
33 #include "Rythmos_StepperBase.hpp"
34 #include "Teuchos_as.hpp"
35 
36 
37 namespace Rythmos {
38 
39 
40 namespace Exceptions {
41 
42 
47 {public: GetFwdPointsFailed(const std::string &my_what):ExceptionBase(my_what) {}};
48 
49 
50 } // namespace Exceptions
51 
52 
61 template<class Scalar>
62 class IntegratorBase : virtual public InterpolationBufferBase<Scalar>
63 {
64 public:
65 
67  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType ScalarMag;
68 
70  virtual RCP<IntegratorBase<Scalar> > cloneIntegrator() const
71  { return Teuchos::null; }
72 
103  virtual void setStepper(
104  const RCP<StepperBase<Scalar> > &stepper,
105  const Scalar &finalTime,
106  const bool landOnFinalTime = true
107  ) =0;
108 
114  virtual Teuchos::RCP<const StepperBase<Scalar> > getStepper() const =0;
115 
121  virtual Teuchos::RCP<StepperBase<Scalar> > getNonconstStepper() const =0;
122 
131  virtual RCP<StepperBase<Scalar> > unSetStepper() =0;
132 
175  virtual void getFwdPoints(
176  const Array<Scalar>& time_vec,
177  Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
178  Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
179  Array<ScalarMag>* accuracy_vec
180  ) =0;
181 
190  virtual TimeRange<Scalar> getFwdTimeRange() const =0;
191 
192 };
193 
194 
195 // 2007/09/14: rabartl: ToDo: Below, Move these functions into a file
196 // Rythmos_IntegratorBaseHelpers.hpp.
197 
198 
203 template<class Scalar>
204 RCP<const Thyra::VectorBase<Scalar> >
205 get_fwd_x( IntegratorBase<Scalar>& integrator, const Scalar t )
206 {
207  Array<Scalar> time_vec;
208  time_vec.push_back(t);
209  Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
210  integrator.getFwdPoints(time_vec,&x_vec,0,0);
211  return x_vec[0];
212 }
213 
214 
239 template<class Scalar>
241  IntegratorBase<Scalar>& integrator,
242  const Scalar t,
243  const Ptr<RCP<const Thyra::VectorBase<Scalar> > > &x,
244  const Ptr<RCP<const Thyra::VectorBase<Scalar> > > &x_dot
245  )
246 {
247  Array<Scalar> time_vec;
248  time_vec.push_back(t);
249  Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
250  Array<RCP<const Thyra::VectorBase<Scalar> > > x_dot_vec;
251  integrator.getFwdPoints(
252  time_vec,
253  nonnull(x) ? &x_vec : 0,
254  nonnull(x_dot) ? &x_dot_vec : 0,
255  0
256  );
257  if (nonnull(x))
258  *x = x_vec[0];
259  if (nonnull(x_dot))
260  *x_dot = x_dot_vec[0];
261 }
262 
263 
265 template<class Scalar>
266 void get_fwd_x_and_x_dot(
267  IntegratorBase<Scalar>& integrator,
268  const Scalar t,
269  RCP<const Thyra::VectorBase<Scalar> > *x,
270  RCP<const Thyra::VectorBase<Scalar> > *x_dot
271  )
272 {
273  get_fwd_x_and_x_dot(integrator, t, ptr(x), ptr(x_dot));
274 }
275 
276 
277 } // namespace Rythmos
278 
279 
280 #endif // Rythmos_INTEGRATOR_BASE_H
281 
RCP< const Thyra::VectorBase< Scalar > > get_fwd_x(IntegratorBase< Scalar > &integrator, const Scalar t)
Nonmember helper function to get x at a (forward) time t.
virtual Teuchos::RCP< const StepperBase< Scalar > > getStepper() const =0
Get the current stepper that is set.
Base class for defining stepper functionality.
Abstract interface for time integrators.
Teuchos::ScalarTraits< Scalar >::magnitudeType ScalarMag
virtual RCP< StepperBase< Scalar > > unSetStepper()=0
Remove the stepper and set *this to an unitilaized state.
virtual void setStepper(const RCP< StepperBase< Scalar > > &stepper, const Scalar &finalTime, const bool landOnFinalTime=true)=0
Specify the stepper to use for integration which effectively reinitializes the intergrator.
virtual RCP< IntegratorBase< Scalar > > cloneIntegrator() const
virtual TimeRange< Scalar > getFwdTimeRange() const =0
Return the valid range of points that the integrator can integrate over.
Base class for an interpolation buffer.
void get_fwd_x_and_x_dot(IntegratorBase< Scalar > &integrator, 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/or x_dot at s (forward) time t.
virtual void getFwdPoints(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)=0
Get values at time points both inside and outside (forward) of current TimeRange. ...
virtual Teuchos::RCP< StepperBase< Scalar > > getNonconstStepper() const =0
Get the current stepper that is set.