Rythmos - Transient Integration for Differential Equations Version of the Day
Rythmos_DefaultIntegrator_decl.hpp
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., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Todd S. Coffey (tscoffe@sandia.gov)
00025 //
00026 // ***********************************************************************
00027 //@HEADER
00028 
00029 #ifndef RYTHMOS_DEFAULT_INTEGRATOR_DECL_HPP
00030 #define RYTHMOS_DEFAULT_INTEGRATOR_DECL_HPP
00031 
00032 
00033 #include "Rythmos_IntegrationControlStrategyAcceptingIntegratorBase.hpp"
00034 #include "Rythmos_InterpolationBufferAppenderAcceptingIntegratorBase.hpp"
00035 #include "Rythmos_TrailingInterpolationBufferAcceptingIntegratorBase.hpp"
00036 #include "Rythmos_IntegrationObserverBase.hpp"
00037 #include "Rythmos_StepControlInfo.hpp"
00038 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
00039 
00040 
00041 
00042 
00043 namespace Rythmos {
00044 
00045 
00049 template<class Scalar> 
00050 class DefaultIntegrator
00051   : virtual public IntegrationControlStrategyAcceptingIntegratorBase<Scalar>,
00052     virtual public InterpolationBufferAppenderAcceptingIntegratorBase<Scalar>,
00053     virtual public TrailingInterpolationBufferAcceptingIntegratorBase<Scalar>,
00054     virtual public Teuchos::ParameterListAcceptorDefaultBase
00055 {
00056 public:
00057   
00059   typedef typename ScalarTraits<Scalar>::magnitudeType ScalarMag;
00060 
00063   
00065   DefaultIntegrator();
00066 
00068   void setIntegrationObserver(
00069     const RCP<IntegrationObserverBase<Scalar> > &integrationObserver
00070     );
00071 
00074 
00076   void setInterpolationBufferAppender(
00077     const RCP<InterpolationBufferAppenderBase<Scalar> > &interpBufferAppender
00078     );
00079 
00081   RCP<const InterpolationBufferAppenderBase<Scalar> >
00082     getInterpolationBufferAppender();
00083 
00085   RCP<InterpolationBufferAppenderBase<Scalar> >
00086     getNonconstInterpolationBufferAppender();
00087 
00089   RCP<InterpolationBufferAppenderBase<Scalar> >
00090     unSetInterpolationBufferAppender();
00091 
00093   
00096   
00098   void setIntegrationControlStrategy(
00099     const RCP<IntegrationControlStrategyBase<Scalar> > &integrationControlStrategy
00100     );
00101 
00103   RCP<IntegrationControlStrategyBase<Scalar> > 
00104     getNonconstIntegrationControlStrategy();
00105 
00107   RCP<const IntegrationControlStrategyBase<Scalar> > 
00108     getIntegrationControlStrategy() const;
00109 
00111 
00114 
00116   void setParameterList(RCP<ParameterList> const& paramList);
00117 
00119   RCP<const ParameterList> getValidParameters() const;
00120 
00122 
00125 
00127   RCP<IntegratorBase<Scalar> > cloneIntegrator() const;
00128   
00130   void setStepper(
00131     const RCP<StepperBase<Scalar> > &stepper,
00132     const Scalar &finalTime,
00133     const bool landOnFinalTime = true
00134     );
00135 
00137   RCP<StepperBase<Scalar> > unSetStepper();
00138 
00140   RCP<const StepperBase<Scalar> > getStepper() const;
00141 
00143   RCP<StepperBase<Scalar> > getNonconstStepper() const;
00144 
00147   
00149   void setTrailingInterpolationBuffer(
00150     const RCP<InterpolationBufferBase<Scalar> > &trailingInterpBuffer
00151     );
00152 
00154   RCP<InterpolationBufferBase<Scalar> >
00155     getNonconstTrailingInterpolationBuffer();
00156 
00158   RCP<const InterpolationBufferBase<Scalar> >
00159     getTrailingInterpolationBuffer() const;
00160 
00162   RCP<InterpolationBufferBase<Scalar> >
00163     unSetTrailingInterpolationBuffer();
00164 
00166 
00168   void getFwdPoints(
00169     const Array<Scalar>& time_vec,
00170     Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00171     Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00172     Array<ScalarMag>* accuracy_vec
00173     );
00174 
00176   TimeRange<Scalar> getFwdTimeRange() const;
00177 
00179 
00182 
00184   RCP<const Thyra::VectorSpaceBase<Scalar> > get_x_space() const;
00185     
00187   void addPoints(
00188     const Array<Scalar>& time_vec,
00189     const Array<RCP<const Thyra::VectorBase<Scalar> > >& x_vec,
00190     const Array<RCP<const Thyra::VectorBase<Scalar> > >& xdot_vec
00191     );
00192 
00194   void getPoints(
00195     const Array<Scalar>& time_vec,
00196     Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00197     Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00198     Array<ScalarMag>* accuracy_vec
00199     ) const;
00200 
00202   TimeRange<Scalar> getTimeRange() const;
00203 
00205   void getNodes(Array<Scalar>* time_vec) const;
00206 
00208   void removeNodes(Array<Scalar>& time_vec);
00209 
00211   int getOrder() const;
00212 
00214 
00215 private:
00216 
00217   // ////////////////////////
00218   // Private data members
00219 
00220   RCP<IntegrationControlStrategyBase<Scalar> > integrationControlStrategy_;
00221   RCP<IntegrationObserverBase<Scalar> > integrationObserver_;
00222 
00223   RCP<InterpolationBufferBase<Scalar> > trailingInterpBuffer_;
00224   RCP<InterpolationBufferAppenderBase<Scalar> > interpBufferAppender_;
00225   
00226   RCP<StepperBase<Scalar> > stepper_;
00227   TimeRange<Scalar> integrationTimeDomain_;
00228   bool landOnFinalTime_;
00229 
00230   int maxNumTimeSteps_;
00231   bool observeInitCond_;
00232 
00233   int currTimeStepIndex_;
00234   StepControlInfo<Scalar> stepCtrlInfoLast_;
00235 
00236   static const std::string maxNumTimeSteps_name_;
00237   static const std::string observeInitCond_name_;
00238   static const int maxNumTimeSteps_default_;
00239 
00240   // /////////////////////////
00241   // Private member functions
00242 
00243   void finalizeSetup();
00244 
00245   bool advanceStepperToTime( const Scalar& t );
00246 
00247 };
00248 
00249 
00254 template<class Scalar> 
00255 RCP<DefaultIntegrator<Scalar> >
00256 defaultIntegrator();
00257 
00258 
00263 template<class Scalar> 
00264 RCP<DefaultIntegrator<Scalar> >
00265 defaultIntegrator(
00266   const RCP<IntegrationControlStrategyBase<Scalar> > &integrationControlStrategy,
00267   const RCP<IntegrationObserverBase<Scalar> > &integrationObserver
00268   );
00269 
00270 
00275 template<class Scalar> 
00276 RCP<DefaultIntegrator<Scalar> >
00277 controlledDefaultIntegrator(
00278   const RCP<IntegrationControlStrategyBase<Scalar> > &integrationControlStrategy
00279   );
00280 
00281 
00286 template<class Scalar> 
00287 RCP<DefaultIntegrator<Scalar> >
00288 observedDefaultIntegrator(
00289   const RCP<IntegrationObserverBase<Scalar> > &integrationObserver
00290   );
00291 
00292 // 2007/08/30: rabartl: Above, note that I had to name the nonmember
00293 // constructors taking an single RCP argument different names from each other
00294 // in order to get around the classic ambiguity problem with implicit
00295 // conversions of smart pointers.
00296 
00297 
00298 } // namespace Rythmos
00299 
00300 
00301 #endif //RYTHMOS_DEFAULT_INTEGRATOR_DECL_HPP
 All Classes Functions Variables Typedefs Friends