Tempus  Version of the Day
Time Integration
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"
13 #include "Tempus_StepperRKBase.hpp"
14 #include "Tempus_StepperImplicit.hpp"
16 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
17  #include "Tempus_StepperRKObserverComposite.hpp"
18 #endif
19 
20 
21 namespace Tempus {
22 
23 /** \brief Diagonally Implicit Runge-Kutta (DIRK) time stepper.
24  *
25  * For the implicit ODE system, \f$\mathcal{F}(\dot{x},x,t) = 0\f$,
26  * the general DIRK method for \f$s\f$-stages, can be written as
27  * \f[
28  * X_{i} = x_{n-1}
29  * + \Delta t\,\sum_{j=1}^{i-1} a_{ij}\,\bar{f}(X_{j},t_{n-1}+c_{j}\Delta t)
30  * + \Delta t\, a_{ii}\,\bar{f}(X_{i},t_{n-1}+c_{i}\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$\dot{x}=\bar{f}(x,t)\f$ is the explicit form of the
37  * ODE, \f$X_{i}\f$ are intermediate approximations to the solution
38  * at times, \f$t_{n-1}+c_{i}\Delta t\f$, (stage solutions) which may
39  * be correct to a lower order of accuracy than the solution, \f$x_{n}\f$.
40  * We should note that these lower-order approximations are combined
41  * through \f$b_{i}\f$ so that error terms cancel out and produce a
42  * more accurate solution. Note for DIRK methods that \f$a_{ii}=a\f$
43  * for all the diagonal components is referred to as
44  * Singly Diagonal Implicit Runge-Kutta (SDIRK) methods.
45  *
46  * Note that the stage time derivatives,
47  * \f[
48  * \dot{X}_{i} = \bar{f}(X_{i},t_{n-1}+c_{i}\Delta t),
49  * \f]
50  * can be found via
51  * \f{eqnarray*}{
52  * \dot{X}_{i} & = & \frac{1}{a_{ii} \Delta t} [ X_{i} - x_{n-1}
53  * - \Delta t\,\sum_{j=1}^{i-1} a_{ij}\,\dot{X}_{j} ] \\
54  * \dot{X}_{i} & = & \frac{X_{i} - \tilde{X}}{a_{ii} \Delta t}
55  * \f}
56  * where
57  * \f[
58  * \tilde{X} = x_{n-1} + \Delta t \sum_{j=1}^{i-1} a_{ij}\, \dot{X}_{j}
59  * \f]
60  * Recalling that the definition for a DIRK is that for \f$j>i\f$,
61  * \f$a_{ij} = 0\f$ and \f$a_{ii} \neq 0\f$ for at least one \f$i\f$.
62  * Thus for stages where \f$a_{ii} = 0\f$, we can use the explicit RK
63  * methods (see StepperExplicitRK for additional details).
64  *
65  * <b> Algorithm </b>
66  * The single-timestep algorithm for DIRK is
67  *
68  * \f{algorithm}{
69  * \renewcommand{\thealgorithm}{}
70  * \caption{DIRK with the application-action locations indicated.}
71  * \begin{algorithmic}[1]
72  * \State {\it appAction.execute(solutionHistory, stepper, BEGIN\_STEP)}
73  * \If {``Reset initial guess.''}
74  * \State $X \leftarrow x_{n-1}$
75  * \Comment{Reset initial guess to last timestep.}
76  * \EndIf
77  * \For {$i = 0 \ldots s-1$}
78  * \If { $a_{k,i} = 0 \;\forall k = (i+1,\ldots, s-1)$, $b(i) = 0$, $b^\ast(i) = 0$}
79  * \State $\dot{X}_i \leftarrow 0$
80  * \Comment{Not needed for later calculations.}
81  * \State {\bf continue}
82  * \EndIf
83  * \State $\tilde{X} \leftarrow
84  * x_{n-1} +\Delta t \sum_{j=1}^{i-1} a_{ij}\,\dot{X}_{j}$
85  * \State {\it appAction.execute(solutionHistory, stepper, BEGIN\_STAGE)}
86  * \If {$a_{ii} = 0$} \Comment{Explicit stage.}
87  * \If {$i=0$ and ``Use FSAL''} \Comment{Save an evaluation?}
88  * \State $\dot{X}_0 \leftarrow \dot{X}_{s-1}$
89  * \Comment{Use $\dot{X}_{s-1}$ from $n-1$ time step.}
90  * \Else
91  * \State $\dot{X}_i \leftarrow \bar{f}(\tilde{X},t_{n-1}+c_i\Delta t)$
92  * \EndIf
93  * \Else \Comment{Implicit stage.}
94  * \State {\it appAction.execute(solutionHistory, stepper, BEFORE\_SOLVE)}
95  * \If {``Zero initial guess.''}
96  * \State $X \leftarrow 0$
97  * \Comment{Else use previous stage value as initial guess.}
98  * \EndIf
99  * \State Solve $\mathcal{F}_i(
100  * \dot{X}_i = \frac{X - \tilde{X}}{a_{ii} \Delta t},
101  * X, t_{n-1}+c_{i}\Delta t) = 0$ for $X$
102  * \State {\it appAction.execute(solutionHistory, stepper, AFTER\_SOLVE)}
103  * \State $\dot{X}_i \leftarrow \frac{X - \tilde{X}}{a_{ii} \Delta t}$
104  * \EndIf
105  * \State {\it appAction.execute(solutionHistory, stepper, BEFORE\_EXPLICIT\_EVAL)}
106  * \State {\it appAction.execute(solutionHistory, stepper, END\_STAGE)}
107  * \EndFor
108  * \State $x_n \leftarrow x_{n-1} + \Delta t\,\sum_{i=0}^{s-1}b_i\,\dot{X}_i$
109  * \If {``Embedded''} \Comment{Compute the local truncation error estimate.}
110  * \State $\mathbf{e} \leftarrow
111  * \sum_{i=0}^{s-1} (b_i-b^\ast_i)\Delta t\,\dot{X}_i$
112  * \State $\tilde{\mathbf{e}} \leftarrow
113  * \mathbf{e}/(a_{tol} + \max(\|x_n\|, \|x_{n-1}\|)r_{tol})$
114  * \State $e_n \leftarrow \|\tilde{\mathbf{e}}\|_\infty$
115  * \EndIf
116  * \State {\it appAction.execute(solutionHistory, stepper, END\_STEP)}
117  * \end{algorithmic}
118  * \f}
119  *
120  * The First-Step-As-Last (FSAL) principle is not needed with DIRK, but
121  * maybe useful if the first stage is explicit (EDIRK) (e.g., Trapezoidal
122  * Method). The default is to set useFSAL=false.
123  *
124  * <b> Iteration Matrix, \f$W\f$.</b>
125  * Recalling that the definition of the iteration matrix, \f$W\f$, is
126  * \f[
127  * W = \alpha \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n}
128  * + \beta \frac{\partial \mathcal{F}_n}{\partial x_n},
129  * \f]
130  * where \f$ \alpha \equiv \frac{\partial \dot{x}_n(x_n) }{\partial x_n}, \f$
131  * and \f$ \beta \equiv \frac{\partial x_n}{\partial x_n} = 1\f$. For the stage
132  * solutions, we have
133  * \f[
134  * \mathcal{F}_i = \dot{X}_{i} - \bar{f}(X_{i},t_{n-1}+c_{i}\Delta t) =0.
135  * \f]
136  * where \f$\mathcal{F}_n \rightarrow \mathcal{F}_i\f$,
137  * \f$x_n \rightarrow X_{i}\f$, and
138  * \f$\dot{x}_n(x_n) \rightarrow \dot{X}_{i}(X_{i})\f$.
139  * The time derivative for the DIRK stages is
140  * \f[
141  * \dot{X}_{i}(X_{i}) = \frac{X_{i} - \tilde{X}}{a_{ii} \Delta t},
142  * \f]
143  * and we can determine that
144  * \f$ \alpha = \frac{1}{a_{ii} \Delta t} \f$
145  * and \f$ \beta = 1 \f$, and therefore write
146  * \f[
147  * W = \frac{1}{a_{ii} \Delta t}
148  * \frac{\partial \mathcal{F}_i}{\partial \dot{X}_i}
149  * + \frac{\partial \mathcal{F}_i}{\partial X_i}.
150  * \f]
151  */
152 template<class Scalar>
153 class StepperDIRK : virtual public Tempus::StepperImplicit<Scalar>,
154  virtual public Tempus::StepperRKBase<Scalar>
155 {
156 public:
157 
158  /// \name Basic stepper methods
159  //@{
160 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
161  virtual void setObserver(
162  Teuchos::RCP<StepperObserver<Scalar> > obs = Teuchos::null);
163 
164  virtual Teuchos::RCP<StepperObserver<Scalar> > getObserver() const
165  { return this->stepperObserver_; }
166 #endif
167  /// Initialize after construction and changing input parameters.
168  virtual void initialize();
169 
170  /// Set the initial conditions and make them consistent.
171  virtual void setInitialConditions (
172  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
173 
174  /// Set parameter so that the initial guess is reset at the beginning of each timestep.
175  virtual void setResetInitialGuess(bool reset_guess)
176  { resetGuess_ = reset_guess; }
177  virtual bool getResetInitialGuess() const
178  { return resetGuess_; }
179 
180  /// Take the specified timestep, dt, and return true if successful.
181  virtual void takeStep(
182  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
183 
184  /// Get a default (initial) StepperState
185  virtual Teuchos::RCP<Tempus::StepperState<Scalar> >getDefaultStepperState();
186 
187  virtual bool isExplicit() const
188  {
189  const int numStages = this->tableau_->numStages();
190  Teuchos::SerialDenseMatrix<int,Scalar> A = this->tableau_->A();
191  bool isExplicit = false;
192  for (int i=0; i<numStages; ++i) if (A(i,i) == 0.0) isExplicit = true;
193  return isExplicit;
194  }
195  virtual bool isImplicit() const {return true;}
196  virtual bool isExplicitImplicit() const
197  {return isExplicit() and isImplicit();}
198  virtual bool isOneStepMethod() const {return true;}
199  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
200 
201  virtual OrderODE getOrderODE() const {return FIRST_ORDER_ODE;}
202 
204  Teuchos::RCP<Teuchos::ParameterList> pl) const;
205 
206  virtual std::string getDescription() const = 0;
207  //@}
208 
209  std::vector<Teuchos::RCP<Thyra::VectorBase<Scalar> > >& getStageXDot() {return stageXDot_;};
210  Teuchos::RCP<Thyra::VectorBase<Scalar> >& getXTilde() {return xTilde_;};
211 
212  /// Return alpha = d(xDot)/dx.
213  virtual Scalar getAlpha(const Scalar dt) const
214  {
215  const Teuchos::SerialDenseMatrix<int,Scalar> & A=this->tableau_->A();
216  return Scalar(1.0)/(dt*A(0,0)); // Getting the first diagonal coeff!
217  }
218  /// Return beta = d(x)/dx.
219  virtual Scalar getBeta (const Scalar ) const { return Scalar(1.0); }
220 
221  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
222 
223  /// \name Overridden from Teuchos::Describable
224  //@{
225  virtual void describe(Teuchos::FancyOStream & out,
226  const Teuchos::EVerbosityLevel verbLevel) const;
227  //@}
228 
229  virtual bool isValidSetup(Teuchos::FancyOStream & out) const;
230 
231  /// \name Accessors methods
232  //@{
233  /** \brief Use embedded if avialable. */
234  virtual void setUseEmbedded(bool a) { useEmbedded_ = a; }
235  virtual bool getUseEmbedded() const { return useEmbedded_; }
236  virtual bool getUseEmbeddedDefault() const { return false; }
237  //@}
238 
239 
240 protected:
241 
242  /// Default setup for constructor.
243  virtual void setupDefault();
244 
245  /// Setup for constructor.
246 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
247  virtual void setup(
248  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& wrapperModel,
249  const Teuchos::RCP<StepperRKObserver<Scalar> >& obs,
250  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
251  bool useFSAL,
252  std::string ICConsistency,
253  bool ICConsistencyCheck,
254  bool useEmbedded,
255  bool zeroInitialGuess);
256 #endif
257  virtual void setup(
258  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& wrapperModel,
259  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
260  bool useFSAL,
261  std::string ICConsistency,
262  bool ICConsistencyCheck,
263  bool useEmbedded,
264  bool zeroInitialGuess,
265  const Teuchos::RCP<StepperRKAppAction<Scalar> >& stepperRKAppAction);
266 
267  virtual void setupTableau() = 0;
268 
269  std::vector<Teuchos::RCP<Thyra::VectorBase<Scalar> > > stageXDot_;
270  Teuchos::RCP<Thyra::VectorBase<Scalar> > xTilde_;
271 
272 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
273  Teuchos::RCP<StepperRKObserverComposite<Scalar> > stepperObserver_;
274 #endif
275 
276  // For Embedded RK
278  Teuchos::RCP<Thyra::VectorBase<Scalar> > ee_;
279  Teuchos::RCP<Thyra::VectorBase<Scalar> > abs_u0;
280  Teuchos::RCP<Thyra::VectorBase<Scalar> > abs_u;
281  Teuchos::RCP<Thyra::VectorBase<Scalar> > sc;
282 
283  bool resetGuess_ = true;
284 };
285 
286 
287 /** \brief Time-derivative interface for DIRK.
288  *
289  * Given the stage state \f$X_i\f$ and
290  * \f[
291  * \tilde{X} = x_{n-1} +\Delta t \sum_{j=1}^{i-1} a_{ij}\,\dot{X}_{j},
292  * \f]
293  * compute the DIRK stage time-derivative,
294  * \f[
295  * \dot{X}_i = \frac{X_{i} - \tilde{X}}{a_{ii} \Delta t}
296  * \f]
297  * \f$\ddot{x}\f$ is not used and set to null.
298  */
299 template <typename Scalar>
301  : virtual public Tempus::TimeDerivative<Scalar>
302 {
303 public:
304 
305  /// Constructor
307  Scalar s, Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde)
308  { initialize(s, xTilde); }
309 
310  /// Destructor
312 
313  /// Compute the time derivative.
314  virtual void compute(
315  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x,
316  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDot,
317  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDotDot = Teuchos::null)
318  {
319  xDotDot = Teuchos::null;
320  Thyra::V_StVpStV(xDot.ptr(),s_,*x,-s_,*xTilde_);
321  }
322 
323  virtual void initialize(Scalar s,
324  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde)
325  { s_ = s; xTilde_ = xTilde; }
326 
327 private:
328 
329  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde_;
330  Scalar s_; // = 1/(dt*a_ii)
331 };
332 
333 
334 } // namespace Tempus
335 
336 #endif // Tempus_StepperDIRK_decl_hpp
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Time-derivative interface for DIRK.
StepperDIRKTimeDerivative(Scalar s, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xTilde)
Constructor.
virtual void initialize(Scalar s, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xTilde)
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< const Thyra::VectorBase< Scalar > > xTilde_
Diagonally Implicit Runge-Kutta (DIRK) time stepper.
Teuchos::RCP< Thyra::VectorBase< Scalar > > xTilde_
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual void setupDefault()
Default setup for constructor.
virtual bool isImplicit() const
virtual bool getUseEmbeddedDefault() const
virtual OrderODE getOrderODE() const
virtual bool isExplicitImplicit() const
virtual bool isOneStepMethod() const
virtual Scalar getAlpha(const Scalar dt) const
Return alpha = d(xDot)/dx.
virtual void setObserver(Teuchos::RCP< StepperObserver< Scalar > > obs=Teuchos::null)
Set Observer.
Teuchos::RCP< Thyra::VectorBase< Scalar > > abs_u0
virtual bool getResetInitialGuess() const
Teuchos::RCP< Thyra::VectorBase< Scalar > > sc
virtual std::string getDescription() const =0
std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar > > > & getStageXDot()
virtual bool getUseEmbedded() const
virtual Teuchos::RCP< StepperObserver< Scalar > > getObserver() const
Get Observer.
virtual Scalar getBeta(const Scalar) const
Return beta = d(x)/dx.
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions and make them consistent.
virtual void initialize()
Initialize after construction and changing input parameters.
virtual bool isValidSetup(Teuchos::FancyOStream &out) const
virtual void setResetInitialGuess(bool reset_guess)
Set parameter so that the initial guess is reset at the beginning of each timestep.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
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.
void getValidParametersBasicDIRK(Teuchos::RCP< Teuchos::ParameterList > pl) const
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
virtual void setupTableau()=0
Teuchos::RCP< Thyra::VectorBase< Scalar > > & getXTilde()
Teuchos::RCP< Thyra::VectorBase< Scalar > > ee_
Teuchos::RCP< Thyra::VectorBase< Scalar > > abs_u
virtual bool isExplicit() const
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.
virtual bool isMultiStepMethod() const
Teuchos::RCP< StepperRKObserverComposite< Scalar > > stepperObserver_
virtual void setUseEmbedded(bool a)
Use embedded if avialable.
std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar > > > stageXDot_
Thyra Base interface for implicit time steppers.
StepperObserver class for Stepper class.
Application Action for StepperRKBase.
Base class for Runge-Kutta methods, ExplicitRK, DIRK and IMEX.
Teuchos::RCP< RKButcherTableau< Scalar > > tableau_
StepperRKObserver class for StepperRK.
This interface defines the time derivative connection between an implicit Stepper and WrapperModelEva...
@ FIRST_ORDER_ODE
Stepper integrates first-order ODEs.