Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_StepperNewmarkImplicitAForm_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_StepperNewmarkImplicitAForm_decl_hpp
10 #define Tempus_StepperNewmarkImplicitAForm_decl_hpp
11 
12 #include "Tempus_StepperImplicit.hpp"
13 #include "Tempus_WrapperModelEvaluatorSecondOrder.hpp"
14 
15 namespace Tempus {
16 
17 
18 /** \brief Newmark time stepper in acceleration form (a-form).
19  *
20  * Here, we implement the Newmark scheme in predictor/corrector form;
21  * see equations (34)-(35) in: A. Mota, W. Klug, M. Ortiz, "Finite element
22  * simulation of firearm injury to the human cranium", Computational
23  * Mechanics 31(1) 115-121 (2003).
24  *
25  * Newmark has two parameters: \f$\beta\f$
26  * and \f$\gamma\f$, both of which need to be in the range \f$[0,1]\f$.
27  * Newmark can be an explicit or implicit method, depending on
28  * the value of the \f$\beta\f$ parameter. If \f$\beta = 0\f$, the method
29  * is explicit. Regardless of whether the method is implicit
30  * or explicit, a linear solve is required. This linear solve can be
31  * optimized, however, for the explicit case by lumping the mass matrix.
32  * This optimization can be invoked by running "Newmark Explicit d-Form"
33  * Stepper through the Piro::TempusSolver class.
34  *
35  * Newmark is second order accurate if \f$\gamma = 0.5\f$; otherwise it
36  * is first order accurate. Some additional properties about the Newmark
37  * scheme can be found
38  *<a href="http://opensees.berkeley.edu/wiki/index.php/Newmark_Method">here</a>.
39  *
40  * The governing equation solved by this stepper is
41  * \f[
42  * \mathbf{M}\, \ddot{\mathbf{x}} + \mathbf{C}\, \dot{\mathbf{x}}
43  * + \mathbf{K}\, \mathbf{x} + \mathbf{F}(t) = 0
44  * \f]
45  * For the A-form (i.e., solving for the acceleration,
46  * \f$\mathbf{a} = \ddot{\mathbf{x}}\f$), we have the following implicit ODE
47  * \f[
48  * \mathbf{M}\, \mathbf{a} + \mathbf{C}\, \mathbf{v}
49  * + \mathbf{K}\, \mathbf{d} + \mathbf{F}(t) =
50  * \mathbf{f}(\mathbf{d}, \mathbf{v}, \mathbf{a}, t) = 0
51  * \f]
52  * where \f$\mathbf{v} = \dot{\mathbf{x}}\f$ and \f$\mathbf{d} = \mathbf{x}\f$.
53  *
54  * <b> Algorithm </b>
55  * The algorithm for the Newmark implicit A-form with predictors and
56  * correctors is
57  * - \f$\mathbf{d}^{\ast} = \mathbf{d}^{n-1} + \Delta t \mathbf{v}^{n-1}
58  * + \Delta t^2 (1-2 \beta) \mathbf{a}^{n-1} / 2\f$
59  * - \f$\mathbf{v}^{\ast} =
60  * \mathbf{v}^{n-1} + \Delta t (1-\gamma) \mathbf{a}^{n-1}\f$
61  * - Solve
62  * \f$\mathbf{f}(\mathbf{d}^n, \mathbf{v}^n, \mathbf{a}^n, t^n) = 0\f$
63  * for \f$\mathbf{a}^n\f$ where
64  * - \f$\mathbf{d}^n = \mathbf{d}^{\ast} + \beta \Delta t^2 \mathbf{a}^n\f$
65  * - \f$\mathbf{v}^n = \mathbf{v}^{\ast} + \gamma \Delta t \mathbf{a}^n\f$
66  *
67  * The First-Step-As-Last (FSAL) principle is part of the Newmark
68  * implicit A-Form as the acceleration from the previous time step is
69  * used for the predictors. The default is to set useFSAL=true,
70  * and useFSAL=false will be ignored.
71  */
72 template<class Scalar>
74  : virtual public Tempus::StepperImplicit<Scalar>
75 {
76 public:
77 
78  /** \brief Default constructor.
79  *
80  * - Constructs with a default ParameterList.
81  * - Can reset ParameterList with setParameterList().
82  * - Requires subsequent setModel() and initialize() calls before calling
83  * takeStep().
84  */
86 
87  /// Constructor
89  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
90  Teuchos::RCP<Teuchos::ParameterList> pList = Teuchos::null);
91 
92  /// \name Basic stepper methods
93  //@{
94  virtual void setModel(
95  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel);
96 
97  virtual void setObserver(
98  Teuchos::RCP<StepperObserver<Scalar> > /* obs */ = Teuchos::null){}
99 
100  /// Initialize during construction and after changing input parameters.
101  virtual void initialize();
102 
103  /// Set the initial conditions and make them consistent.
104  virtual void setInitialConditions (
105  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
106 
107  /// Take the specified timestep, dt, and return true if successful.
108  virtual void takeStep(
109  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
110 
111  /// Get a default (initial) StepperState
112  virtual Teuchos::RCP<Tempus::StepperState<Scalar> >
114  virtual Scalar getOrder() const {
115  if (gamma_ == 0.5) return 2.0;
116  else return 1.0;
117  }
118  virtual Scalar getOrderMin() const {return 1.0;}
119  virtual Scalar getOrderMax() const {return 2.0;}
120 
121  virtual bool isExplicit() const {return false;}
122  virtual bool isImplicit() const {return true;}
123  virtual bool isExplicitImplicit() const
124  {return isExplicit() and isImplicit();}
125  virtual bool isOneStepMethod() const {return true;}
126  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
127 
128  virtual OrderODE getOrderODE() const {return SECOND_ORDER_ODE;}
129  //@}
130 
131  /// Return W_xDotxDot_coeff = d(xDotDot)/d(xDotDot).
132  virtual Scalar getW_xDotDot_coeff (const Scalar) const {return Scalar(1.0);}
133  /// Return alpha = d(xDot)/d(xDotDot).
134  virtual Scalar getAlpha(const Scalar dt) const { return gamma_*dt; }
135  /// Return beta = d(x)/d(xDotDot).
136  virtual Scalar getBeta (const Scalar dt) const { return beta_*dt*dt; }
137 
138  /// \name ParameterList methods
139  //@{
140  void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & pl);
141  Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList();
142  Teuchos::RCP<Teuchos::ParameterList> unsetParameterList();
143  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
144  Teuchos::RCP<Teuchos::ParameterList> getDefaultParameters() const;
145  //@}
146 
147  /// \name Overridden from Teuchos::Describable
148  //@{
149  virtual std::string description() const;
150  virtual void describe(Teuchos::FancyOStream & out,
151  const Teuchos::EVerbosityLevel verbLevel) const;
152  //@}
153 
154  void predictVelocity(Thyra::VectorBase<Scalar>& vPred,
155  const Thyra::VectorBase<Scalar>& v,
156  const Thyra::VectorBase<Scalar>& a,
157  const Scalar dt) const;
158 
159  void predictDisplacement(Thyra::VectorBase<Scalar>& dPred,
160  const Thyra::VectorBase<Scalar>& d,
161  const Thyra::VectorBase<Scalar>& v,
162  const Thyra::VectorBase<Scalar>& a,
163  const Scalar dt) const;
164 
165  void correctVelocity(Thyra::VectorBase<Scalar>& v,
166  const Thyra::VectorBase<Scalar>& vPred,
167  const Thyra::VectorBase<Scalar>& a,
168  const Scalar dt) const;
169 
170  void correctDisplacement(Thyra::VectorBase<Scalar>& d,
171  const Thyra::VectorBase<Scalar>& dPred,
172  const Thyra::VectorBase<Scalar>& a,
173  const Scalar dt) const;
174 private:
175 
176  Scalar beta_;
177  Scalar gamma_;
178 
179  Teuchos::RCP<Teuchos::FancyOStream> out_;
180 
181 };
182 } // namespace Tempus
183 
184 #endif // Tempus_StepperNewmarkImplicitAForm_decl_hpp
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl)
virtual Scalar getW_xDotDot_coeff(const Scalar) const
Return W_xDotxDot_coeff = d(xDotDot)/d(xDotDot).
void predictVelocity(Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
virtual Scalar getAlpha(const Scalar dt) const
Return alpha = d(xDot)/d(xDotDot).
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
virtual Scalar getBeta(const Scalar dt) const
Return beta = d(x)/d(xDotDot).
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Stepper integrates second-order ODEs.
Thyra Base interface for implicit time steppers.
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
StepperObserver class for Stepper class.
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions and make them consistent.
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
void predictDisplacement(Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel)
void correctDisplacement(Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
virtual void initialize()
Initialize during construction and after changing input parameters.
void correctVelocity(Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< Teuchos::ParameterList > getDefaultParameters() const
Newmark time stepper in acceleration form (a-form).
virtual void setObserver(Teuchos::RCP< StepperObserver< Scalar > >=Teuchos::null)
Set Observer.