Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_StepperDIRK_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_StepperDIRK_decl_hpp
10 #define Tempus_StepperDIRK_decl_hpp
11 
12 #include "Tempus_config.hpp"
14 #include "Tempus_StepperImplicit.hpp"
16 #include "Tempus_StepperRKObserverComposite.hpp"
17 
18 
19 namespace Tempus {
20 
21 /** \brief Diagonally Implicit Runge-Kutta (DIRK) time stepper.
22  *
23  * For the implicit ODE system, \f$\mathcal{F}(\dot{x},x,t) = 0\f$,
24  * the general DIRK method for \f$s\f$-stages, can be written as
25  * \f[
26  * X_{i} = x_{n-1}
27  * + \Delta t\,\sum_{j=1}^{i-1} a_{ij}\,\bar{f}(X_{j},t_{n-1}+c_{j}\Delta t)
28  * + \Delta t\, a_{ii}\,\bar{f}(X_{i},t_{n-1}+c_{i}\Delta t)
29  * \f]
30  * \f[
31  * x_{n} = x_{n-1}
32  * + \Delta t\,\sum_{i=1}^{s}b_{i}\,\bar{f}(X_{i},t_{n-1}+c_{i}\Delta t)
33  * \f]
34  * where \f$\dot{x}=\bar{f}(x,t)\f$ is the explicit form of the
35  * ODE, \f$X_{i}\f$ are intermediate approximations to the solution
36  * at times, \f$t_{n-1}+c_{i}\Delta t\f$, (stage solutions) which may
37  * be correct to a lower order of accuracy than the solution, \f$x_{n}\f$.
38  * We should note that these lower-order approximations are combined
39  * through \f$b_{i}\f$ so that error terms cancel out and produce a
40  * more accurate solution. Note for DIRK methods that \f$a_{ii}=a\f$
41  * for all the diagonal components is referred to as
42  * Singly Diagonal Implicit Runge-Kutta (SDIRK) methods.
43  *
44  * Note that the stage time derivatives,
45  * \f[
46  * \dot{X}_{i} = \bar{f}(X_{i},t_{n-1}+c_{i}\Delta t),
47  * \f]
48  * can be found via
49  * \f{eqnarray*}{
50  * \dot{X}_{i} & = & \frac{1}{a_{ii} \Delta t} [ X_{i} - x_{n-1}
51  * - \Delta t\,\sum_{j=1}^{i-1} a_{ij}\,\dot{X}_{j} ] \\
52  * \dot{X}_{i} & = & \frac{X_{i} - \tilde{X}}{a_{ii} \Delta t}
53  * \f}
54  * where
55  * \f[
56  * \tilde{X} = x_{n-1} + \Delta t \sum_{j=1}^{i-1} a_{ij}\, \dot{X}_{j}
57  * \f]
58  * Recalling that the definition for a DIRK is that for \f$j>i\f$,
59  * \f$a_{ij} = 0\f$ and \f$a_{ii} \neq 0\f$ for at least one \f$i\f$.
60  * Thus for stages where \f$a_{ii} = 0\f$, we can use the explicit RK
61  * methods (see StepperExplicitRK for additional details).
62  *
63  * <b> Algorithm </b>
64  * The single-timestep algorithm for DIRK is,
65  * - for \f$i = 1 \ldots s\f$ do
66  * - if \f$a_{ii} = 0\f$
67  * - \f$X_i \leftarrow x_{n-1}
68  * + \Delta t\,\sum_{j=1}^{i-1} a_{ij}\,\dot{X}_j\f$
69  * - Evaluate \f$\bar{f}(X_{i},t_{n-1}+c_{i}\Delta t)\f$
70  * - \f$\dot{X}_i \leftarrow \bar{f}(X_i,t_{n-1}+c_i\Delta t)\f$
71  * - else
72  * - \f$\tilde{X} =
73  * x_{n-1} +\Delta t \sum_{j=1}^{i-1} a_{ij}\,\dot{X}_{j}\f$
74  * - Define \f$\dot{X}_i \leftarrow
75  * \frac{X_{i} - \tilde{X}}{a_{ii} \Delta t}\f$
76  * - Solve \f$f(\dot{x} =
77  * \dot{X}_i,X_i,t_{n-1}+c_{i}\Delta t)=0\f$ for \f$X_i\f$
78  * - \f$\dot{X}_i \leftarrow \frac{X_{i} - \tilde{X}}{a_{ii} \Delta t}\f$
79  * - end for
80  * - \f$x_n \leftarrow x_{n-1} + \Delta t\,\sum_{i=1}^{s}b_i\,\dot{X}_i\f$
81  *
82  * The First-Step-As-Last (FSAL) principle is not needed with DIRK, but
83  * maybe useful if the first stage is explicit (i.e, EDIRK).
84  * The default is to set useFSAL=false.
85  */
86 template<class Scalar>
87 class StepperDIRK : virtual public Tempus::StepperImplicit<Scalar>
88 {
89 public:
90 
91  /// \name Basic stepper methods
92  //@{
93  virtual void setObserver(
94  Teuchos::RCP<StepperObserver<Scalar> > obs = Teuchos::null);
95 
96  virtual Teuchos::RCP<const RKButcherTableau<Scalar> > getTableau()
97  { return tableau_; }
98 
99  virtual Teuchos::RCP<StepperObserver<Scalar> > getObserver() const
100  { return this->stepperObserver_; }
101 
102  /// Initialize during construction and after changing input parameters.
103  virtual void initialize();
104 
105  /// Set the initial conditions and make them consistent.
106  virtual void setInitialConditions (
107  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
108 
109  /// Set parameter so that the initial guess is reset at the beginning of each timestep.
110  virtual void setResetInitialGuess(bool reset_guess)
111  { resetGuess_ = reset_guess; }
112  virtual bool getResetInitialGuess() const
113  { return resetGuess_; }
114 
115  /// Take the specified timestep, dt, and return true if successful.
116  virtual void takeStep(
117  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
118 
119  /// Get a default (initial) StepperState
120  virtual Teuchos::RCP<Tempus::StepperState<Scalar> >getDefaultStepperState();
121  virtual Scalar getOrder() const{return tableau_->order();}
122  virtual Scalar getOrderMin() const{return tableau_->orderMin();}
123  virtual Scalar getOrderMax() const{return tableau_->orderMax();}
124 
125  virtual bool isExplicit() const
126  {
127  const int numStages = tableau_->numStages();
128  Teuchos::SerialDenseMatrix<int,Scalar> A = tableau_->A();
129  bool isExplicit = false;
130  for (int i=0; i<numStages; ++i) if (A(i,i) == 0.0) isExplicit = true;
131  return isExplicit;
132  }
133  virtual bool isImplicit() const {return true;}
134  virtual bool isExplicitImplicit() const
135  {return isExplicit() and isImplicit();}
136  virtual bool isOneStepMethod() const {return true;}
137  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
138 
139  virtual OrderODE getOrderODE() const {return FIRST_ORDER_ODE;}
140 
142  Teuchos::RCP<Teuchos::ParameterList> pl) const;
143 
144  virtual std::string getDescription() const = 0;
145  //@}
146 
147  /// Return alpha = d(xDot)/dx.
148  virtual Scalar getAlpha(const Scalar dt) const
149  {
150  const Teuchos::SerialDenseMatrix<int,Scalar> & A=tableau_->A();
151  return Scalar(1.0)/(dt*A(0,0)); // Getting the first diagonal coeff!
152  }
153  /// Return beta = d(x)/dx.
154  virtual Scalar getBeta (const Scalar ) const { return Scalar(1.0); }
155 
156  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
157 
158  /// \name Overridden from Teuchos::Describable
159  //@{
160  virtual void describe(Teuchos::FancyOStream & out,
161  const Teuchos::EVerbosityLevel verbLevel) const;
162  //@}
163 
164  /// \name Accessors methods
165  //@{
166  /** \brief Use embedded if avialable. */
167  virtual void setUseEmbedded(bool a) { useEmbedded_ = a; }
168  virtual bool getUseEmbedded() const { return useEmbedded_; }
169  virtual bool getUseEmbeddedDefault() const { return false; }
170  //@}
171 
172 
173 protected:
174 
175  /// Default setup for constructor.
176  virtual void setupDefault();
177 
178  /// Setup for constructor.
179  virtual void setup(
180  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& wrapperModel,
181  const Teuchos::RCP<StepperRKObserver<Scalar> >& obs,
182  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
183  bool useFSAL,
184  std::string ICConsistency,
185  bool ICConsistencyCheck,
186  bool useEmbedded,
187  bool zeroInitialGuess);
188 
189  virtual void setupTableau() = 0;
190 
191  Teuchos::RCP<RKButcherTableau<Scalar> > tableau_;
192 
193  std::vector<Teuchos::RCP<Thyra::VectorBase<Scalar> > > stageXDot_;
194  Teuchos::RCP<Thyra::VectorBase<Scalar> > stageX_;
195  Teuchos::RCP<Thyra::VectorBase<Scalar> > xTilde_;
196 
197  Teuchos::RCP<StepperRKObserverComposite<Scalar> > stepperObserver_;
198 
199  // For Embedded RK
201  Teuchos::RCP<Thyra::VectorBase<Scalar> > ee_;
202  Teuchos::RCP<Thyra::VectorBase<Scalar> > abs_u0;
203  Teuchos::RCP<Thyra::VectorBase<Scalar> > abs_u;
204  Teuchos::RCP<Thyra::VectorBase<Scalar> > sc;
205 
206  bool resetGuess_ = true;
207 };
208 
209 
210 /** \brief Time-derivative interface for DIRK.
211  *
212  * Given the stage state \f$X_i\f$ and
213  * \f[
214  * \tilde{X} = x_{n-1} +\Delta t \sum_{j=1}^{i-1} a_{ij}\,\dot{X}_{j},
215  * \f]
216  * compute the DIRK stage time-derivative,
217  * \f[
218  * \dot{X}_i = \frac{X_{i} - \tilde{X}}{a_{ii} \Delta t}
219  * \f]
220  * \f$\ddot{x}\f$ is not used and set to null.
221  */
222 template <typename Scalar>
224  : virtual public Tempus::TimeDerivative<Scalar>
225 {
226 public:
227 
228  /// Constructor
230  Scalar s, Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde)
231  { initialize(s, xTilde); }
232 
233  /// Destructor
235 
236  /// Compute the time derivative.
237  virtual void compute(
238  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x,
239  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDot,
240  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDotDot = Teuchos::null)
241  {
242  xDotDot = Teuchos::null;
243  Thyra::V_StVpStV(xDot.ptr(),s_,*x,-s_,*xTilde_);
244  }
245 
246  virtual void initialize(Scalar s,
247  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde)
248  { s_ = s; xTilde_ = xTilde; }
249 
250 private:
251 
252  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde_;
253  Scalar s_; // = 1/(dt*a_ii)
254 };
255 
256 
257 } // namespace Tempus
258 
259 #endif // Tempus_StepperDIRK_decl_hpp
virtual Scalar getBeta(const Scalar) const
Return beta = d(x)/dx.
Teuchos::RCP< Thyra::VectorBase< Scalar > > xTilde_
virtual bool isMultiStepMethod() const
Teuchos::RCP< const Thyra::VectorBase< Scalar > > xTilde_
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
Teuchos::RCP< Thyra::VectorBase< Scalar > > ee_
virtual void setup(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &wrapperModel, const Teuchos::RCP< StepperRKObserver< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, bool zeroInitialGuess)
Setup for constructor.
virtual void setupTableau()=0
void getValidParametersBasicDIRK(Teuchos::RCP< Teuchos::ParameterList > pl) const
virtual Teuchos::RCP< const RKButcherTableau< Scalar > > getTableau()
virtual bool getUseEmbedded() const
virtual void compute(Teuchos::RCP< const Thyra::VectorBase< Scalar > > x, Teuchos::RCP< Thyra::VectorBase< Scalar > > xDot, Teuchos::RCP< Thyra::VectorBase< Scalar > > xDotDot=Teuchos::null)
Compute the time derivative.
Teuchos::RCP< Thyra::VectorBase< Scalar > > stageX_
virtual void setResetInitialGuess(bool reset_guess)
Set parameter so that the initial guess is reset at the beginning of each timestep.
virtual Teuchos::RCP< StepperObserver< Scalar > > getObserver() const
Get Observer.
virtual bool getUseEmbeddedDefault() const
Diagonally Implicit Runge-Kutta (DIRK) time stepper.
virtual bool isExplicitImplicit() const
virtual Scalar getOrderMin() const
Thyra Base interface for implicit time steppers.
virtual bool isExplicit() const
Teuchos::RCP< Thyra::VectorBase< Scalar > > abs_u
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions and make them consistent.
Teuchos::RCP< Thyra::VectorBase< Scalar > > abs_u0
virtual void initialize()
Initialize during construction and after changing input parameters.
virtual bool isImplicit() const
Time-derivative interface for DIRK.
virtual void setupDefault()
Default setup for constructor.
StepperObserver class for Stepper class.
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Teuchos::RCP< StepperRKObserverComposite< Scalar > > stepperObserver_
virtual Scalar getOrderMax() const
Stepper integrates first-order ODEs.
virtual bool isOneStepMethod() const
virtual std::string getDescription() const =0
Teuchos::RCP< RKButcherTableau< Scalar > > tableau_
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
virtual void setObserver(Teuchos::RCP< StepperObserver< Scalar > > obs=Teuchos::null)
Set Observer.
virtual Scalar getAlpha(const Scalar dt) const
Return alpha = d(xDot)/dx.
virtual OrderODE getOrderODE() const
This interface defines the time derivative connection between an implicit Stepper and WrapperModelEva...
virtual Scalar getOrder() const
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
StepperRKObserver class for StepperRK.
Teuchos::RCP< Thyra::VectorBase< Scalar > > sc
std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar > > > stageXDot_
virtual void initialize(Scalar s, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xTilde)
virtual void setUseEmbedded(bool a)
Use embedded if avialable.
virtual bool getResetInitialGuess() const
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.
StepperDIRKTimeDerivative(Scalar s, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xTilde)
Constructor.