Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_StepperNewmarkImplicitDForm_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_StepperNewmarkImplicitDForm_decl_hpp
10 #define Tempus_StepperNewmarkImplicitDForm_decl_hpp
11 
12 #include "Tempus_WrapperModelEvaluatorSecondOrder.hpp"
13 #include "Tempus_StepperImplicit.hpp"
14 
15 namespace Tempus {
16 
17 /** \brief Newmark time stepper.
18  *
19  * Here, we implement the Newmark scheme in displacement predictor/corrector
20  * form; see equations (34)-(35) in: A. Mota, W. Klug, M. Ortiz,
21  * "Finite element simulation of firearm injury to the human cranium",
22  * Computational Mechanics 31(1) 115-121 (2003).
23  *
24  * Newmark has two parameters: \f$\beta\f$
25  * and \f$\gamma\f$, both of which need to be in the range \f$[0,1]\f$.
26  * Newmark can be an explicit or implicit method, depending on
27  * the value of the \f$\beta\f$ parameter. If \f$\beta = 0\f$, the method
28  * is explicit; but note that the d-form of the Newmark scheme is not defined
29  * for the explicit case.
30  *
31  * Newmark is second order accurate if \f$\gamma = 0.5\f$; otherwise it
32  * is first order accurate. Some additional properties about the Newmark
33  * Beta scheme can be found
34  * <a href="http://opensees.berkeley.edu/wiki/index.php/Newmark_Method">here</a>.
35  *
36  * The First-Step-As-Last (FSAL) principle is not used with the
37  * Newmark implicit D-Form method.
38  */
39 template <class Scalar>
40 class StepperNewmarkImplicitDForm : virtual public Tempus::StepperImplicit<Scalar> {
41  public:
42 
43  /** \brief Default constructor.
44  *
45  * - Constructs with a default ParameterList.
46  * - Can reset ParameterList with setParameterList().
47  * - Requires subsequent setModel() and initialize() calls before calling
48  * takeStep().
49  */
51 
52  /// Constructor
54  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel,
55  Teuchos::RCP<Teuchos::ParameterList> pList = Teuchos::null);
56 
57  /// \name Basic stepper methods
58  //@{
59  virtual void
60  setModel(const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel);
61 
62  virtual void setObserver(
63  Teuchos::RCP<StepperObserver<Scalar> > /* obs */ = Teuchos::null){}
64 
65  /// Initialize during construction and after changing input parameters.
66  virtual void
67  initialize();
68 
69  /// Set the initial conditions and make them consistent.
70  virtual void setInitialConditions (
71  const Teuchos::RCP<SolutionHistory<Scalar> >& /* solutionHistory */){}
72 
73  /// Take the specified timestep, dt, and return true if successful.
74  virtual void
75  takeStep(const Teuchos::RCP<SolutionHistory<Scalar>>& solutionHistory);
76 
77  /// Get a default (initial) StepperState
78  virtual Teuchos::RCP<Tempus::StepperState<Scalar>>
80  virtual Scalar
81  getOrder() const {
82  if (gamma_ == 0.5)
83  return 2.0;
84  else
85  return 1.0;
86  }
87  virtual Scalar
88  getOrderMin() const {
89  return 1.0;
90  }
91  virtual Scalar
92  getOrderMax() const {
93  return 2.0;
94  }
95  virtual bool isExplicit() const {return false;}
96  virtual bool isImplicit() const {return true;}
97  virtual bool isExplicitImplicit() const
98  {return isExplicit() and isImplicit();}
99  virtual bool isOneStepMethod() const {return true;}
100  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
101 
102  virtual OrderODE getOrderODE() const {return SECOND_ORDER_ODE;}
103  //@}
104 
105  /// Return W_xDotxDot_coeff = d(xDotDot)/d(x).
106  virtual Scalar getW_xDotDot_coeff (const Scalar dt) const
107  { return Scalar(1.0)/(beta_*dt*dt); }
108  /// Return alpha = d(xDot)/d(x).
109  virtual Scalar getAlpha(const Scalar dt) const { return gamma_/(beta_*dt); }
110  /// Return beta = d(x)/d(x).
111  virtual Scalar getBeta (const Scalar ) const { return Scalar(1.0); }
112 
113  /// \name ParameterList methods
114  //@{
115  void
116  setParameterList(const Teuchos::RCP<Teuchos::ParameterList>& pl);
117  Teuchos::RCP<Teuchos::ParameterList>
119  Teuchos::RCP<Teuchos::ParameterList>
121  Teuchos::RCP<const Teuchos::ParameterList>
122  getValidParameters() const;
123  Teuchos::RCP<Teuchos::ParameterList>
124  getDefaultParameters() const;
125  //@}
126 
127  /// \name Overridden from Teuchos::Describable
128  //@{
129  virtual std::string
130  description() const;
131  virtual void
132  describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel)
133  const;
134  //@}
135 
136  void
138  Thyra::VectorBase<Scalar>& vPred, const Thyra::VectorBase<Scalar>& v,
139  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const;
140 
141  void
143  Thyra::VectorBase<Scalar>& dPred, const Thyra::VectorBase<Scalar>& d,
144  const Thyra::VectorBase<Scalar>& v, const Thyra::VectorBase<Scalar>& a,
145  const Scalar dt) const;
146 
147  void
149  Thyra::VectorBase<Scalar>& v, const Thyra::VectorBase<Scalar>& vPred,
150  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const;
151 
152  void
154  Thyra::VectorBase<Scalar>& d, const Thyra::VectorBase<Scalar>& dPred,
155  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const;
156 
157  void
159  Thyra::VectorBase<Scalar>& a, const Thyra::VectorBase<Scalar>& dPred,
160  const Thyra::VectorBase<Scalar>& d, const Scalar dt) const;
161 
162  protected:
163 
164  Thyra::ModelEvaluatorBase::InArgs<Scalar> inArgs_;
165  Thyra::ModelEvaluatorBase::OutArgs<Scalar> outArgs_;
166 
167  Scalar beta_;
168  Scalar gamma_;
169 
170  Teuchos::RCP<Teuchos::FancyOStream> out_;
171 
172 };
173 } // namespace Tempus
174 
175 #endif // Tempus_StepperNewmarkImplicitDForm_decl_hpp
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Thyra::ModelEvaluatorBase::OutArgs< Scalar > outArgs_
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl)
void correctVelocity(Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
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.
Thyra::ModelEvaluatorBase::InArgs< Scalar > inArgs_
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
Stepper integrates second-order ODEs.
Thyra Base interface for implicit time steppers.
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar >> &solutionHistory)
Take the specified timestep, dt, and return true if successful.
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
void predictVelocity(Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
virtual void setObserver(Teuchos::RCP< StepperObserver< Scalar > >=Teuchos::null)
Set Observer.
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &)
Set the initial conditions and make them consistent.
StepperObserver class for Stepper class.
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
virtual Scalar getW_xDotDot_coeff(const Scalar dt) const
Return W_xDotxDot_coeff = d(xDotDot)/d(x).
Teuchos::RCP< Teuchos::ParameterList > getDefaultParameters() const
virtual Scalar getAlpha(const Scalar dt) const
Return alpha = d(xDot)/d(x).
virtual Scalar getBeta(const Scalar) const
Return beta = d(x)/d(x).
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar >> &appModel)
void correctAcceleration(Thyra::VectorBase< Scalar > &a, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d, const Scalar dt) const