Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_StepperExplicitRK_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #ifndef Tempus_StepperExplicitRK_decl_hpp
10 #define Tempus_StepperExplicitRK_decl_hpp
11 
12 #include "Tempus_config.hpp"
13 #include "Tempus_StepperExplicit.hpp"
15 #include "Tempus_StepperRKObserverComposite.hpp"
16 
17 
18 namespace Tempus {
19 
20 /** \brief Explicit Runge-Kutta time stepper.
21  *
22  * For the explicit ODE system,
23  * \f[
24  * \dot{x} = \bar{f}(x,t),
25  * \f]
26  * the general explicit Runge-Kutta method for \f$s\f$-stages can be
27  * written as
28  * \f[
29  * X_{i} = x_{n-1}
30  * + \Delta t\,\sum_{j=1}^{i-1} a_{ij}\,\bar{f}(X_{j},t_{n-1}+c_{j}\Delta t)
31  * \f]
32  * \f[
33  * x_{n} = x_{n-1}
34  * + \Delta t\,\sum_{i=1}^{s}b_{i}\,\bar{f}(X_{i},t_{n-1}+c_{i}\Delta t)
35  * \f]
36  * where \f$X_{i}\f$ are intermediate approximations to the solution
37  * at times, \f$t_{n-1}+c_{i}\Delta t\f$, (stage solutions) which may
38  * be correct to a lower order of accuracy than the solution, \f$x_{n}\f$.
39  * We should note that these lower-order approximations are combined
40  * through \f$b_{i}\f$ so that error terms cancel out and produce a
41  * more accurate solution. Note for explicit RK that \f$a_{ij}=0\f$ for
42  * \f$j \leq i\f$ and does not require any solves.
43  * Note that the stage time derivatives are
44  * \f[
45  * \dot{X}_{i} = \bar{f}(X_{i},t_{n-1}+c_{i}\Delta t),
46  * \f]
47  * and the time derivative by definition is
48  * \f[
49  * \dot{x}_{n} = \bar{f}(x_{n},t_{n}),
50  * \f]
51  *
52  * <b> Algorithm </b>
53  * The single-timestep algorithm for Explicit RK is simply,
54  * - for \f$i = 1 \ldots s\f$ do
55  * - \f$X_i \leftarrow x_{n-1}
56  * + \Delta t\,\sum_{j=1}^{i-1} a_{ij}\,\dot{X}_j\f$
57  * - Evaluate \f$\bar{f}(X_{i},t_{n-1}+c_{i}\Delta t)\f$
58  * - \f$\dot{X}_i \leftarrow \bar{f}(X_i,t_{n-1}+c_i\Delta t)\f$
59  * - end for
60  * - \f$x_n \leftarrow x_{n-1} + \Delta t\,\sum_{i=1}^{s}b_i\,\dot{X}_i\f$
61  *
62  * When using the First-Step-As-Last (FSAL) priniciple, where one can
63  * reuse the last function evaulation as the first evaluation of the next
64  * time step, the algorithm is only slightly more complicated.
65  * - for \f$i = 1 \ldots s\f$ do
66  * - if ( i==1 && useFSAL && (previous step not failed) )
67  * - tmp = \f$\dot{X}_1\f$
68  * - \f$\dot{X}_1 = \dot{X}_s\f$
69  * - \f$\dot{X}_s\f$ = tmp
70  * - else
71  * - \f$X_i \leftarrow x_{n-1}
72  * + \Delta t\,\sum_{j=1}^{i-1} a_{ij}\,\dot{X}_j\f$
73  * - Evaluate \f$\bar{f}(X_{i},t_{n-1}+c_{i}\Delta t)\f$
74  * - \f$\dot{X}_i \leftarrow \bar{f}(X_i,t_{n-1}+c_i\Delta t)\f$
75  * - end for
76  * - \f$x_n \leftarrow x_{n-1} + \Delta t\,\sum_{i=1}^{s}b_i\,\dot{X}_i\f$
77  *
78  * For Explicit RK, FSAL requires \f$c_1 = 0\f$, \f$c_s = 1\f$, and
79  * be stiffly accurate (\f$a_{sj} = b_j\f$). An example of this is
80  * the Bogacki-Shampine 3(2) method.
81  * \f[
82  * \begin{array}{c|cccc} 0 & 0 & & & \\
83  * 1/3 & 1/2 & 0 & & \\
84  * 2/3 & 0 & 3/4 & 0 & \\
85  * 1 & 2/9 & 1/3 & 4/9 & 0 \\ \hline
86  * & 2/9 & 1/3 & 4/9 & 0 \\
87  * & 7/24 & 1/4 & 1/3 & 1/8 \end{array}
88  * \f]
89  */
90 template<class Scalar>
91 class StepperExplicitRK : virtual public Tempus::StepperExplicit<Scalar>
92 {
93 
94 public:
95 
96  /// \name Basic stepper methods
97  //@{
98  virtual void setObserver(
99  Teuchos::RCP<StepperObserver<Scalar> > obs = Teuchos::null);
100 
101  virtual Teuchos::RCP<const RKButcherTableau<Scalar> > getTableau()
102  { return tableau_; }
103 
104  virtual Teuchos::RCP<StepperObserver<Scalar> > getObserver() const
105  { return this->stepperObserver_; }
106 
107  /// Initialize during construction and after changing input parameters.
108  virtual void initialize();
109 
110  /// Set the initial conditions and make them consistent.
111  virtual void setInitialConditions (
112  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
113 
114  /// Take the specified timestep, dt, and return true if successful.
115  virtual void takeStep(
116  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
117 
118  /// Get a default (initial) StepperState
119  virtual Teuchos::RCP<Tempus::StepperState<Scalar> > getDefaultStepperState();
120  virtual Scalar getOrder() const {return tableau_->order();}
121  virtual Scalar getOrderMin() const {return tableau_->orderMin();}
122  virtual Scalar getOrderMax() const {return tableau_->orderMax();}
123  virtual Scalar getInitTimeStep(
124  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory) const;
125 
126  virtual Teuchos::RCP<Thyra::VectorBase<Scalar> > getStageX() {return stageX_;}
127 
128  virtual bool isExplicit() const {return true;}
129  virtual bool isImplicit() const {return false;}
130  virtual bool isExplicitImplicit() const
131  {return isExplicit() and isImplicit();}
132  virtual bool isOneStepMethod() const {return true;}
133  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
134 
135  virtual OrderODE getOrderODE() const {return FIRST_ORDER_ODE;}
136 
137  void getValidParametersBasicERK(Teuchos::RCP<Teuchos::ParameterList> pl) const;
138  virtual std::string getDescription() const = 0;
139  //@}
140 
141  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
142 
143  /// \name Overridden from Teuchos::Describable
144  //@{
145  virtual void describe(Teuchos::FancyOStream & out,
146  const Teuchos::EVerbosityLevel verbLevel) const;
147  //@}
148 
149  /// \name Accessors methods
150  //@{
151  /** \brief Use embedded if avialable. */
152  virtual void setUseEmbedded(bool a) { useEmbedded_ = a; }
153  virtual bool getUseEmbedded() const { return useEmbedded_; }
154  virtual bool getUseEmbeddedDefault() const { return false; }
155  //@}
156 
157 
158 protected:
159 
160  /// Default setup for constructor.
161  virtual void setupDefault();
162 
163  /// Setup for constructor.
164  virtual void setup(
165  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
166  const Teuchos::RCP<StepperRKObserverComposite<Scalar> >& obs,
167  bool useFSAL,
168  std::string ICConsistency,
169  bool ICConsistencyCheck,
170  bool useEmbedded);
171 
172  virtual void setupTableau() = 0;
173 
174 
175  Teuchos::RCP<RKButcherTableau<Scalar> > tableau_;
176 
177  std::vector<Teuchos::RCP<Thyra::VectorBase<Scalar> > > stageXDot_;
178  Teuchos::RCP<Thyra::VectorBase<Scalar> > stageX_;
179 
180  Teuchos::RCP<StepperRKObserverComposite<Scalar> > stepperObserver_;
181 
182  // For Embedded RK
184  Teuchos::RCP<Thyra::VectorBase<Scalar> > ee_;
185  Teuchos::RCP<Thyra::VectorBase<Scalar> > abs_u0;
186  Teuchos::RCP<Thyra::VectorBase<Scalar> > abs_u;
187  Teuchos::RCP<Thyra::VectorBase<Scalar> > sc;
188 
189 };
190 
191 } // namespace Tempus
192 
193 #endif // Tempus_StepperExplicitRK_decl_hpp
Teuchos::RCP< RKButcherTableau< Scalar > > tableau_
Explicit Runge-Kutta time stepper.
virtual std::string getDescription() const =0
virtual void setupDefault()
Default setup for constructor.
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > getStageX()
virtual Scalar getInitTimeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) const
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar > > > stageXDot_
Teuchos::RCP< Thyra::VectorBase< Scalar > > ee_
virtual void setup(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperRKObserverComposite< Scalar > > &obs, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded)
Setup for constructor.
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions and make them consistent.
virtual void setupTableau()=0
virtual Teuchos::RCP< const RKButcherTableau< Scalar > > getTableau()
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
virtual Teuchos::RCP< StepperObserver< Scalar > > getObserver() const
Get Observer.
virtual void setObserver(Teuchos::RCP< StepperObserver< Scalar > > obs=Teuchos::null)
Set Observer.
Teuchos::RCP< Thyra::VectorBase< Scalar > > sc
StepperObserver class for Stepper class.
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
Teuchos::RCP< Thyra::VectorBase< Scalar > > stageX_
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
void getValidParametersBasicERK(Teuchos::RCP< Teuchos::ParameterList > pl) const
Stepper integrates first-order ODEs.
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.
virtual void initialize()
Initialize during construction and after changing input parameters.
Teuchos::RCP< StepperRKObserverComposite< Scalar > > stepperObserver_
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
virtual void setUseEmbedded(bool a)
Use embedded if avialable.
Thyra Base interface for implicit time steppers.
Teuchos::RCP< Thyra::VectorBase< Scalar > > abs_u0
Teuchos::RCP< Thyra::VectorBase< Scalar > > abs_u