Rythmos_IntegratorBase.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_INTEGRATOR_BASE_H
00030 #define Rythmos_INTEGRATOR_BASE_H
00031
00032 #include "Rythmos_InterpolationBufferBase.hpp"
00033 #include "Rythmos_StepperBase.hpp"
00034 #include "Teuchos_as.hpp"
00035
00036
00037 namespace Rythmos {
00038
00039
00040 namespace Exceptions {
00041
00042
00046 class GetFwdPointsFailed : public ::Rythmos::Exceptions::ExceptionBase
00047 {public: GetFwdPointsFailed(const std::string &what):ExceptionBase(what) {}};
00048
00049
00050 }
00051
00052
00061 template<class Scalar>
00062 class IntegratorBase : virtual public InterpolationBufferBase<Scalar>
00063 {
00064 public:
00065
00067 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType ScalarMag;
00068
00070 virtual RCP<IntegratorBase<Scalar> > cloneIntegrator() const
00071 { return Teuchos::null; }
00072
00074 virtual bool acceptsTrailingInterpolationBuffer() const
00075 { return false; }
00076
00078 virtual void setTrailingInterpolationBuffer(
00079 const RCP<InterpolationBufferBase<Scalar> > &trailingInterpBuffer
00080 )
00081 { TEST_FOR_EXCEPT(true); }
00082
00084 virtual RCP<InterpolationBufferBase<Scalar> >
00085 getTrailingInterpolationBuffer()
00086 { return Teuchos::null; }
00087
00089 virtual RCP<const InterpolationBufferBase<Scalar> >
00090 getTrailingInterpolationBuffer() const
00091 { return Teuchos::null; }
00092
00122 virtual void setStepper(
00123 const RCP<StepperBase<Scalar> > &stepper,
00124 const Scalar &finalTime,
00125 const bool landOnFinalTime = true
00126 ) =0;
00127
00133 virtual Teuchos::RCP<const StepperBase<Scalar> > getStepper() const =0;
00134
00143 virtual RCP<StepperBase<Scalar> > unSetStepper() =0;
00144
00187 virtual void getFwdPoints(
00188 const Array<Scalar>& time_vec,
00189 Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00190 Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00191 Array<ScalarMag>* accuracy_vec
00192 ) =0;
00193
00202 virtual TimeRange<Scalar> getFwdTimeRange() const =0;
00203
00204 };
00205
00206
00207
00208
00209
00210
00215 template<class Scalar>
00216 RCP<const Thyra::VectorBase<Scalar> >
00217 get_fwd_x( IntegratorBase<Scalar>& integrator, const Scalar t )
00218 {
00219 Array<Scalar> time_vec;
00220 time_vec.push_back(t);
00221 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
00222 integrator.getFwdPoints(time_vec,&x_vec,0,0);
00223 return x_vec[0];
00224 }
00225
00226
00251 template<class Scalar>
00252 void get_fwd_x_and_x_dot(
00253 IntegratorBase<Scalar>& integrator,
00254 const Scalar t,
00255 RCP<const Thyra::VectorBase<Scalar> > *x,
00256 RCP<const Thyra::VectorBase<Scalar> > *x_dot
00257 )
00258 {
00259 Array<Scalar> time_vec;
00260 time_vec.push_back(t);
00261 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec;
00262 Array<RCP<const Thyra::VectorBase<Scalar> > > x_dot_vec;
00263 integrator.getFwdPoints(
00264 time_vec,
00265 x ? &x_vec : 0,
00266 x_dot ? &x_dot_vec : 0,
00267 0
00268 );
00269 if (x) *x = x_vec[0];
00270 if (x_dot) *x_dot = x_dot_vec[0];
00271 }
00272
00273
00274 }
00275
00276
00277 #endif // Rythmos_INTEGRATOR_BASE_H
00278