Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_StepperTrapezoidal_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_StepperTrapezoidal_decl_hpp
10 #define Tempus_StepperTrapezoidal_decl_hpp
11 
12 #include "Tempus_StepperImplicit.hpp"
15 
16 
17 namespace Tempus {
18 
19 /** \brief Trapezoidal method time stepper.
20  *
21  * For the implicit ODE system, \f$\mathcal{F}(\dot{x},x,t) = 0\f$,
22  * the solution, \f$\dot{x}\f$ and \f$x\f$, is determined using a
23  * solver (e.g., a non-linear solver, like NOX).
24  *
25  * <b> Algorithm </b>
26  * The single-timestep algorithm for Trapezoidal method is simply,
27  * - Solve \f$\mathcal{F}_n(\dot{x}=(x_n-x_{n-1})/(\Delta t_n/2) - \dot{x}_{n-1}, x_n, t_n)=0\f$ for \f$x_n\f$
28  * - \f$\dot{x}_n \leftarrow (x_n-x_{n-1})/(\Delta t_n/2) - \dot{x}_{n-1}\f$
29  *
30  * The First-Step-As-Last (FSAL) principle is required for the Trapezoidal
31  * Stepper (i.e., useFSAL=true)! There are at least two ways around this,
32  * but are not implemented.
33  * - Do a solve for xDotOld, xDot_{n-1}, at each time step as for the
34  * initial conditions. This is expensive since you would be doing
35  * two solves every time step.
36  * - Use evaluateExplicitODE to get xDot_{n-1} if the application
37  * provides it. Explicit evaluations are cheaper but requires the
38  * application to implement xDot = f(x,t).
39  *
40  * <b> Iteration Matrix, \f$W\f$.</b>
41  * Recalling that the definition of the iteration matrix, \f$W\f$, is
42  * \f[
43  * W = \alpha \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n}
44  * + \beta \frac{\partial \mathcal{F}_n}{\partial x_n},
45  * \f]
46  * where \f$ \alpha \equiv \frac{\partial \dot{x}_n(x_n) }{\partial x_n}, \f$
47  * and \f$ \beta \equiv \frac{\partial x_n}{\partial x_n} = 1\f$, and
48  * the time derivative for Trapezoidal method is
49  * \f[
50  * \dot{x}_n = (x_n-x_{n-1})/(\Delta t/2) - \dot{x}_{n-1},
51  * \f]
52  * we can determine that
53  * \f$ \alpha = \frac{2}{\Delta t} \f$
54  * and \f$ \beta = 1 \f$, and therefore write
55  * \f[
56  * W = \frac{2}{\Delta t}
57  * \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n}
58  * + \frac{\partial \mathcal{F}_n}{\partial x_n}.
59  * \f]
60  */
61 template<class Scalar>
62 class StepperTrapezoidal : virtual public Tempus::StepperImplicit<Scalar>
63 {
64 public:
65 
66  /** \brief Default constructor.
67  *
68  * Requires subsequent setModel(), setSolver() and initialize()
69  * calls before calling takeStep().
70  */
72 
73  /// Constructor
75  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
76  const Teuchos::RCP<StepperObserver<Scalar> >& obs,
77  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
78  bool useFSAL,
79  std::string ICConsistency,
80  bool ICConsistencyCheck,
81  bool zeroInitialGuess);
82 
83  /// \name Basic stepper methods
84  //@{
85  virtual void setObserver(
86  Teuchos::RCP<StepperObserver<Scalar> > obs = Teuchos::null);
87 
88  virtual Teuchos::RCP<StepperObserver<Scalar> > getObserver() const
89  { return this->stepperTrapObserver_; }
90 
91  /// Set the initial conditions and make them consistent.
92  virtual void setInitialConditions (
93  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
94 
95  /// Take the specified timestep, dt, and return true if successful.
96  virtual void takeStep(
97  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
98 
99  /// Get a default (initial) StepperState
100  virtual Teuchos::RCP<Tempus::StepperState<Scalar> > getDefaultStepperState();
101  virtual Scalar getOrder() const {return 2.0;}
102  virtual Scalar getOrderMin() const {return 2.0;}
103  virtual Scalar getOrderMax() const {return 2.0;}
104 
105  virtual bool isExplicit() const {return false;}
106  virtual bool isImplicit() const {return true;}
107  virtual bool isExplicitImplicit() const
108  {return isExplicit() and isImplicit();}
109  virtual bool isOneStepMethod() const {return true;}
110  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
111  virtual OrderODE getOrderODE() const {return FIRST_ORDER_ODE;}
112  //@}
113 
114  /// Return alpha = d(xDot)/dx.
115  virtual Scalar getAlpha(const Scalar dt) const { return Scalar(2.0)/dt; }
116  /// Return beta = d(x)/dx.
117  virtual Scalar getBeta (const Scalar ) const { return Scalar(1.0); }
118 
119  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
120 
121  /// \name Overridden from Teuchos::Describable
122  //@{
123  virtual void describe(Teuchos::FancyOStream & out,
124  const Teuchos::EVerbosityLevel verbLevel) const;
125  //@}
126 
127  virtual bool isValidSetup(Teuchos::FancyOStream & out) const;
128 
129  virtual bool getUseFSALDefault() const { return true; }
130  virtual std::string getICConsistencyDefault() const { return "Consistent"; }
131 
132 private:
133 
134  Teuchos::RCP<StepperTrapezoidalObserver<Scalar> > stepperTrapObserver_;
135 
136 };
137 
138 /** \brief Time-derivative interface for Trapezoidal method.
139  *
140  * Given the state \f$x\f$, compute the Trapezoidal method time-derivative,
141  * \f[
142  * \dot{x}_{n} = \frac{(x_{n} - x_{n-1})}{(\Delta t_n/2)} - \dot{x}_{n-1}.
143  * \f]
144  * \f$\ddot{x}\f$ is not used and set to null.
145  */
146 template <typename Scalar>
148  : virtual public Tempus::TimeDerivative<Scalar>
149 {
150 public:
151 
152  /// Constructor
154  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xOld,
155  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xDotOld)
156  { initialize(s, xOld, xDotOld); }
157 
158  /// Destructor
160 
161  /// Compute the time derivative.
162  virtual void compute(
163  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x,
164  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDot,
165  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDotDot = Teuchos::null)
166  {
167  xDotDot = Teuchos::null;
168  // Calculate the Trapezoidal method x dot vector
169  Thyra::V_StVpStV(xDot.ptr(),s_,*x,-s_,*xOld_);
170  Thyra::V_VpStV (xDot.ptr(),*xDot,Scalar(-1.0),*xDotOld_);
171  }
172 
173  virtual void initialize(Scalar s,
174  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xOld,
175  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xDotOld)
176  { s_ = s; xOld_ = xOld; xDotOld_ = xDotOld; }
177 
178 private:
179 
180  Scalar s_; // = 1.0/(dt/2)
181  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xOld_;
182  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xDotOld_;
183 };
184 
185 
186 } // namespace Tempus
187 
188 #endif // Tempus_StepperTrapezoidal_decl_hpp
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
Trapezoidal method time stepper.
StepperTrapezoidalTimeDerivative(Scalar s, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xOld, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xDotOld)
Constructor.
virtual void initialize(Scalar s, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xOld, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xDotOld)
Teuchos::RCP< StepperTrapezoidalObserver< Scalar > > stepperTrapObserver_
virtual Teuchos::RCP< StepperObserver< Scalar > > getObserver() const
Get Observer.
Teuchos::RCP< const Thyra::VectorBase< Scalar > > xDotOld_
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions and make them consistent.
Thyra Base interface for implicit time steppers.
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 > > xOld_
StepperObserver class for Stepper class.
virtual std::string getICConsistencyDefault() const
Time-derivative interface for Trapezoidal method.
virtual Scalar getAlpha(const Scalar dt) const
Return alpha = d(xDot)/dx.
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
virtual Scalar getBeta(const Scalar) const
Return beta = d(x)/dx.
Stepper integrates first-order ODEs.
virtual bool isValidSetup(Teuchos::FancyOStream &out) const
This interface defines the time derivative connection between an implicit Stepper and WrapperModelEva...
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual void setObserver(Teuchos::RCP< StepperObserver< Scalar > > obs=Teuchos::null)
Set Observer.
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const