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$f(\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 template<class Scalar>
41 class StepperTrapezoidal : virtual public Tempus::StepperImplicit<Scalar>
42 {
43 public:
44 
45  /** \brief Default constructor.
46  *
47  * Requires subsequent setModel(), setSolver() and initialize()
48  * calls before calling takeStep().
49  */
51 
52  /// Constructor
54  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
55  const Teuchos::RCP<StepperObserver<Scalar> >& obs,
56  const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
57  bool useFSAL,
58  std::string ICConsistency,
59  bool ICConsistencyCheck,
60  bool zeroInitialGuess);
61 
62  /// \name Basic stepper methods
63  //@{
64  virtual void setObserver(
65  Teuchos::RCP<StepperObserver<Scalar> > obs = Teuchos::null);
66 
67  virtual Teuchos::RCP<StepperObserver<Scalar> > getObserver() const
68  { return this->stepperTrapObserver_; }
69 
70  /// Initialize during construction and after changing input parameters.
71  virtual void initialize();
72 
73  /// Set the initial conditions and make them consistent.
74  virtual void setInitialConditions (
75  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
76 
77  /// Take the specified timestep, dt, and return true if successful.
78  virtual void takeStep(
79  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
80 
81  /// Get a default (initial) StepperState
82  virtual Teuchos::RCP<Tempus::StepperState<Scalar> > getDefaultStepperState();
83  virtual Scalar getOrder() const {return 2.0;}
84  virtual Scalar getOrderMin() const {return 2.0;}
85  virtual Scalar getOrderMax() const {return 2.0;}
86 
87  virtual bool isExplicit() const {return false;}
88  virtual bool isImplicit() const {return true;}
89  virtual bool isExplicitImplicit() const
90  {return isExplicit() and isImplicit();}
91  virtual bool isOneStepMethod() const {return true;}
92  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
93  virtual OrderODE getOrderODE() const {return FIRST_ORDER_ODE;}
94  //@}
95 
96  /// Return alpha = d(xDot)/dx.
97  virtual Scalar getAlpha(const Scalar dt) const { return Scalar(2.0)/dt; }
98  /// Return beta = d(x)/dx.
99  virtual Scalar getBeta (const Scalar ) const { return Scalar(1.0); }
100 
101  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
102 
103  /// \name Overridden from Teuchos::Describable
104  //@{
105  virtual void describe(Teuchos::FancyOStream & out,
106  const Teuchos::EVerbosityLevel verbLevel) const;
107  //@}
108 
109  virtual bool getUseFSALDefault() const { return true; }
110  virtual std::string getICConsistencyDefault() const { return "Consistent"; }
111 
112 private:
113 
114  Teuchos::RCP<StepperTrapezoidalObserver<Scalar> > stepperTrapObserver_;
115 
116 };
117 
118 /** \brief Time-derivative interface for Trapezoidal method.
119  *
120  * Given the state \f$x\f$, compute the Trapezoidal method time-derivative,
121  * \f[
122  * \dot{x}_{n} = \frac{(x_{n} - x_{n-1})}{(\Delta t_n/2)} - \dot{x}_{n-1}.
123  * \f]
124  * \f$\ddot{x}\f$ is not used and set to null.
125  */
126 template <typename Scalar>
128  : virtual public Tempus::TimeDerivative<Scalar>
129 {
130 public:
131 
132  /// Constructor
134  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xOld,
135  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xDotOld)
136  { initialize(s, xOld, xDotOld); }
137 
138  /// Destructor
140 
141  /// Compute the time derivative.
142  virtual void compute(
143  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x,
144  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDot,
145  Teuchos::RCP< Thyra::VectorBase<Scalar> > xDotDot = Teuchos::null)
146  {
147  xDotDot = Teuchos::null;
148  // Calculate the Trapezoidal method x dot vector
149  Thyra::V_StVpStV(xDot.ptr(),s_,*x,-s_,*xOld_);
150  Thyra::V_VpStV (xDot.ptr(),*xDot,Scalar(-1.0),*xDotOld_);
151  }
152 
153  virtual void initialize(Scalar s,
154  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xOld,
155  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xDotOld)
156  { s_ = s; xOld_ = xOld; xDotOld_ = xDotOld; }
157 
158 private:
159 
160  Scalar s_; // = 1.0/(dt/2)
161  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xOld_;
162  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xDotOld_;
163 };
164 
165 
166 } // namespace Tempus
167 
168 #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_
virtual void initialize()
Initialize during construction and after changing input parameters.
StepperObserver class for Stepper class.
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
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.
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