|
Rythmos - Transient Integration for Differential Equations
Version of the Day
|
00001 //@HEADER 00002 // *********************************************************************** 00003 // 00004 // Rythmos Package 00005 // Copyright (2006) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // This library is free software; you can redistribute it and/or modify 00011 // it under the terms of the GNU Lesser General Public License as 00012 // published by the Free Software Foundation; either version 2.1 of the 00013 // License, or (at your option) any later version. 00014 // 00015 // This library is distributed in the hope that it will be useful, but 00016 // WITHOUT ANY WARRANTY; without even the implied warranty of 00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00018 // Lesser General Public License for more details. 00019 // 00020 // You should have received a copy of the GNU Lesser General Public 00021 // License along with this library; if not, write to the Free Software 00022 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 00023 // USA 00024 // Questions? Contact Todd S. Coffey (tscoffe@sandia.gov) 00025 // 00026 // *********************************************************************** 00027 //@HEADER 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 &my_what):ExceptionBase(my_what) {}}; 00048 00049 00050 } // namespace Exceptions 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 00103 virtual void setStepper( 00104 const RCP<StepperBase<Scalar> > &stepper, 00105 const Scalar &finalTime, 00106 const bool landOnFinalTime = true 00107 ) =0; 00108 00114 virtual Teuchos::RCP<const StepperBase<Scalar> > getStepper() const =0; 00115 00121 virtual Teuchos::RCP<StepperBase<Scalar> > getNonconstStepper() const =0; 00122 00131 virtual RCP<StepperBase<Scalar> > unSetStepper() =0; 00132 00175 virtual void getFwdPoints( 00176 const Array<Scalar>& time_vec, 00177 Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec, 00178 Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec, 00179 Array<ScalarMag>* accuracy_vec 00180 ) =0; 00181 00190 virtual TimeRange<Scalar> getFwdTimeRange() const =0; 00191 00192 }; 00193 00194 00195 // 2007/09/14: rabartl: ToDo: Below, Move these functions into a file 00196 // Rythmos_IntegratorBaseHelpers.hpp. 00197 00198 00203 template<class Scalar> 00204 RCP<const Thyra::VectorBase<Scalar> > 00205 get_fwd_x( IntegratorBase<Scalar>& integrator, const Scalar t ) 00206 { 00207 Array<Scalar> time_vec; 00208 time_vec.push_back(t); 00209 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec; 00210 integrator.getFwdPoints(time_vec,&x_vec,0,0); 00211 return x_vec[0]; 00212 } 00213 00214 00239 template<class Scalar> 00240 void get_fwd_x_and_x_dot( 00241 IntegratorBase<Scalar>& integrator, 00242 const Scalar t, 00243 const Ptr<RCP<const Thyra::VectorBase<Scalar> > > &x, 00244 const Ptr<RCP<const Thyra::VectorBase<Scalar> > > &x_dot 00245 ) 00246 { 00247 Array<Scalar> time_vec; 00248 time_vec.push_back(t); 00249 Array<RCP<const Thyra::VectorBase<Scalar> > > x_vec; 00250 Array<RCP<const Thyra::VectorBase<Scalar> > > x_dot_vec; 00251 integrator.getFwdPoints( 00252 time_vec, 00253 nonnull(x) ? &x_vec : 0, 00254 nonnull(x_dot) ? &x_dot_vec : 0, 00255 0 00256 ); 00257 if (nonnull(x)) 00258 *x = x_vec[0]; 00259 if (nonnull(x_dot)) 00260 *x_dot = x_dot_vec[0]; 00261 } 00262 00263 00265 template<class Scalar> 00266 void get_fwd_x_and_x_dot( 00267 IntegratorBase<Scalar>& integrator, 00268 const Scalar t, 00269 RCP<const Thyra::VectorBase<Scalar> > *x, 00270 RCP<const Thyra::VectorBase<Scalar> > *x_dot 00271 ) 00272 { 00273 get_fwd_x_and_x_dot(integrator, t, ptr(x), ptr(x_dot)); 00274 } 00275 00276 00277 } // namespace Rythmos 00278 00279 00280 #endif // Rythmos_INTEGRATOR_BASE_H 00281
1.7.6.1