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"
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  /** \brief Default constructor.
92  *
93  * - Constructs with a default ParameterList.
94  * - Can reset ParameterList with setParameterList().
95  * - Requires subsequent setModel() and initialize() calls before calling
96  * takeStep().
97  */
98  StepperDIRK();
99 
100  /// Constructor to specialize Stepper parameters.
101  StepperDIRK(
102  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
103  Teuchos::RCP<Teuchos::ParameterList> pList);
104 
105  /// Constructor to use default Stepper parameters.
106  StepperDIRK(
107  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
108  std::string stepperType = "SDIRK 2 Stage 2nd order");
109 
110  /// Constructor for StepperFactory.
111  StepperDIRK(
112  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
113  std::string stepperType, Teuchos::RCP<Teuchos::ParameterList> pList);
114 
115  /// \name Basic stepper methods
116  //@{
117  virtual void setObserver(
118  Teuchos::RCP<StepperObserver<Scalar> > obs = Teuchos::null);
119 
120  virtual void setTableau(std::string stepperType);
121 
122  virtual void setTableau(
123  Teuchos::RCP<Teuchos::ParameterList> pList = Teuchos::null);
124 
125  virtual void setTableau(
126  Teuchos::RCP<const RKButcherTableau<Scalar> > DIRK_ButcherTableau);
127 
128  /// Initialize during construction and after changing input parameters.
129  virtual void initialize();
130 
131  /// Set the initial conditions and make them consistent.
132  virtual void setInitialConditions (
133  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
134 
135  /// Take the specified timestep, dt, and return true if successful.
136  virtual void takeStep(
137  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
138 
139  /// Get a default (initial) StepperState
140  virtual Teuchos::RCP<Tempus::StepperState<Scalar> >getDefaultStepperState();
141  virtual Scalar getOrder() const{return DIRK_ButcherTableau_->order();}
142  virtual Scalar getOrderMin() const{return DIRK_ButcherTableau_->orderMin();}
143  virtual Scalar getOrderMax() const{return DIRK_ButcherTableau_->orderMax();}
144 
145  virtual bool isExplicit() const
146  {
147  const int numStages = DIRK_ButcherTableau_->numStages();
148  Teuchos::SerialDenseMatrix<int,Scalar> A = DIRK_ButcherTableau_->A();
149  bool isExplicit = false;
150  for (int i=0; i<numStages; ++i) if (A(i,i) == 0.0) isExplicit = true;
151  return isExplicit;
152  }
153  virtual bool isImplicit() const {return true;}
154  virtual bool isExplicitImplicit() const
155  {return isExplicit() and isImplicit();}
156  virtual bool isOneStepMethod() const {return true;}
157  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
158 
159  virtual OrderODE getOrderODE() const {return FIRST_ORDER_ODE;}
160  //@}
161 
162  /// Return alpha = d(xDot)/dx.
163  virtual Scalar getAlpha(const Scalar dt) const
164  {
165  const Teuchos::SerialDenseMatrix<int,Scalar> & A=DIRK_ButcherTableau_->A();
166  return Scalar(1.0)/(dt*A(0,0)); // Getting the first diagonal coeff!
167  }
168  /// Return beta = d(x)/dx.
169  virtual Scalar getBeta (const Scalar ) const { return Scalar(1.0); }
170 
171  /// \name ParameterList methods
172  //@{
173  void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & pl);
174  Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList();
175  Teuchos::RCP<Teuchos::ParameterList> unsetParameterList();
176  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
177  Teuchos::RCP<Teuchos::ParameterList> getDefaultParameters() const;
178  //@}
179 
180  /// \name Overridden from Teuchos::Describable
181  //@{
182  virtual std::string description() const;
183  virtual void describe(Teuchos::FancyOStream & out,
184  const Teuchos::EVerbosityLevel verbLevel) const;
185  //@}
186 
187 protected:
188 
189  Teuchos::RCP<const RKButcherTableau<Scalar> > DIRK_ButcherTableau_;
190 
191  std::vector<Teuchos::RCP<Thyra::VectorBase<Scalar> > > stageXDot_;
192  Teuchos::RCP<Thyra::VectorBase<Scalar> > stageX_;
193  Teuchos::RCP<Thyra::VectorBase<Scalar> > xTilde_;
194 
195  Teuchos::RCP<StepperDIRKObserver<Scalar> > stepperDIRKObserver_;
196 
197  Teuchos::RCP<Thyra::VectorBase<Scalar> > ee_;
198 
199  // For Embedded RK
200  Teuchos::RCP<Thyra::VectorBase<Scalar> > abs_u0;
201  Teuchos::RCP<Thyra::VectorBase<Scalar> > abs_u;
202  Teuchos::RCP<Thyra::VectorBase<Scalar> > sc;
203 
204 };
205 
206 
207 /** \brief Time-derivative interface for DIRK.
208  *
209  * Given the stage state \f$X_i\f$ and
210  * \f[
211  * \tilde{X} = x_{n-1} +\Delta t \sum_{j=1}^{i-1} a_{ij}\,\dot{X}_{j},
212  * \f]
213  * compute the DIRK stage time-derivative,
214  * \f[
215  * \dot{X}_i = \frac{X_{i} - \tilde{X}}{a_{ii} \Delta t}
216  * \f]
217  * \f$\ddot{x}\f$ is not used and set to null.
218  */
219 template <typename Scalar>
221  : virtual public Tempus::TimeDerivative<Scalar>
222 {
223 public:
224 
225  /// Constructor
227  Scalar s, Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde)
228  { initialize(s, xTilde); }
229 
230  /// Destructor
232 
233  /// Compute the time derivative.
234  virtual void compute(
235  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x,
236  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDot,
237  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDotDot = Teuchos::null)
238  {
239  xDotDot = Teuchos::null;
240  Thyra::V_StVpStV(xDot.ptr(),s_,*x,-s_,*xTilde_);
241  }
242 
243  virtual void initialize(Scalar s,
244  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde)
245  { s_ = s; xTilde_ = xTilde; }
246 
247 private:
248 
249  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xTilde_;
250  Scalar s_; // = 1/(dt*a_ii)
251 };
252 
253 
254 } // namespace Tempus
255 
256 #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 void setTableau(std::string stepperType)
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
Teuchos::RCP< Thyra::VectorBase< Scalar > > ee_
StepperDIRK()
Default constructor.
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_
Teuchos::RCP< const RKButcherTableau< Scalar > > DIRK_ButcherTableau_
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 std::string description() const
StepperObserver class for Stepper class.
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl)
virtual Scalar getOrderMax() const
Teuchos::RCP< StepperDIRKObserver< Scalar > > stepperDIRKObserver_
Stepper integrates first-order ODEs.
virtual bool isOneStepMethod() const
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
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)
Teuchos::RCP< Teuchos::ParameterList > getDefaultParameters() 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.