Tempus  Version of the Day
Time Integration
Tempus_StepperHHTAlpha_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_StepperHHTAlpha_decl_hpp
10 #define Tempus_StepperHHTAlpha_decl_hpp
11 
12 #include "Tempus_StepperImplicit.hpp"
13 #include "Tempus_WrapperModelEvaluatorSecondOrder.hpp"
14 
15 namespace Tempus {
16 
17 
18 /** \brief HHT-Alpha time stepper.
19  *
20  * Here, we implement the HHT-Alpha scheme in predictor/corrector form;
21  * see equations (10) and (13)-(19) in: G.M. Hulbert, J. Chung,
22  * "Explicit time integration algorithms for structural dynamics with
23  * optimal numerical dissipation", Comput. Methods Appl. Mech. Engrg.
24  * 137 175-188 (1996).
25  *
26  * There are four parameters in the scheme: \f$\alpha_m\f$, \f$\alpha_f\f$,
27  * \f$\beta\f$ and \f$\gamma\f$, all of which must be in the range \f$[0,1]\f$.
28  * When \f$\alpha_m=\alpha_f = 0\f$, the scheme reduces to the Newmark Beta
29  * scheme (see Tempus::StepperNewmark for details). Like the Newmark Beta
30  * scheme, the HHT-Alpha scheme can be either first or second order accurate,
31  * and either explicit or implicit.
32  *
33  * Although the general form of the scheme has been implemented in Tempus,
34  * it has only been verified for the case when \f$\alpha_m=\alpha_f = 0\f$
35  * (corresponding to the Newmark Beta) scheme, so other values for these
36  * parameters are not allowed at the present time. Also, note that, like
37  * the Newmark Beta stepper, the linear solve for the explicit version of
38  * this scheme has not been optimized (the mass matrix is not lumped).
39  */
40 template<class Scalar>
41 class StepperHHTAlpha : virtual public Tempus::StepperImplicit<Scalar>
42 {
43 public:
44 
45  /// Constructor
47  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
48  Teuchos::RCP<Teuchos::ParameterList> pList = Teuchos::null);
49 
50  /// \name Basic stepper methods
51  //@{
52  virtual void setModel(
53  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel);
54 
55  virtual void setObserver(
56  Teuchos::RCP<StepperObserver<Scalar> > obs = Teuchos::null){}
57 
58  /// Initialize during construction and after changing input parameters.
59  virtual void initialize();
60 
61  /// Take the specified timestep, dt, and return true if successful.
62  virtual void takeStep(
63  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
64 
65  /// Pass initial guess to Newton solver (only relevant for implicit solvers)
66  virtual void setInitialGuess(Teuchos::RCP<const Thyra::VectorBase<Scalar> > initial_guess)
67  {initial_guess_ = initial_guess;}
68 
69  /// Get a default (initial) StepperState
70  virtual Teuchos::RCP<Tempus::StepperState<Scalar> > getDefaultStepperState();
71  virtual Scalar getOrder() const {
72  if (gamma_ == 0.5) return 2.0;
73  else return 1.0;
74  }
75  virtual Scalar getOrderMin() const {return 1.0;}
76  virtual Scalar getOrderMax() const {return 2.0;}
77 
78  virtual bool isExplicit() const {return false;}
79  virtual bool isImplicit() const {return true;}
80  virtual bool isExplicitImplicit() const
81  {return isExplicit() and isImplicit();}
82  virtual bool isOneStepMethod() const {return true;}
83  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
84  //@}
85 
86  /// \name ParameterList methods
87  //@{
88  void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & pl);
89  Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList();
90  Teuchos::RCP<Teuchos::ParameterList> unsetParameterList();
91  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
92  Teuchos::RCP<Teuchos::ParameterList> getDefaultParameters() const;
93  //@}
94 
95  /// \name Overridden from Teuchos::Describable
96  //@{
97  virtual std::string description() const;
98  virtual void describe(Teuchos::FancyOStream & out,
99  const Teuchos::EVerbosityLevel verbLevel) const;
100  //@}
101 
102  void predictVelocity(Thyra::VectorBase<Scalar>& vPred,
103  const Thyra::VectorBase<Scalar>& v,
104  const Thyra::VectorBase<Scalar>& a,
105  const Scalar dt) const;
106 
107  void predictDisplacement(Thyra::VectorBase<Scalar>& dPred,
108  const Thyra::VectorBase<Scalar>& d,
109  const Thyra::VectorBase<Scalar>& v,
110  const Thyra::VectorBase<Scalar>& a,
111  const Scalar dt) const;
112 
113  void predictVelocity_alpha_f(Thyra::VectorBase<Scalar>& vPred,
114  const Thyra::VectorBase<Scalar>& v) const;
115 
116  void predictDisplacement_alpha_f(Thyra::VectorBase<Scalar>& dPred,
117  const Thyra::VectorBase<Scalar>& d) const;
118 
119  void correctAcceleration(Thyra::VectorBase<Scalar>& a_n_plus1,
120  const Thyra::VectorBase<Scalar>& a_n) const;
121 
122  void correctVelocity(Thyra::VectorBase<Scalar>& v,
123  const Thyra::VectorBase<Scalar>& vPred,
124  const Thyra::VectorBase<Scalar>& a,
125  const Scalar dt) const;
126 
127  void correctDisplacement(Thyra::VectorBase<Scalar>& d,
128  const Thyra::VectorBase<Scalar>& dPred,
129  const Thyra::VectorBase<Scalar>& a,
130  const Scalar dt) const;
131 
132 private:
133 
134  /// Default Constructor -- not allowed
135  StepperHHTAlpha();
136 
137 private:
138 
139  Scalar beta_;
140  Scalar gamma_;
141  Scalar alpha_f_;
142  Scalar alpha_m_;
143 
144  Teuchos::RCP<Teuchos::FancyOStream> out_;
145 
146  Teuchos::RCP<const Thyra::VectorBase<Scalar> > initial_guess_;
147 
148 
149 };
150 } // namespace Tempus
151 
152 #endif // Tempus_StepperHHTAlpha_decl_hpp
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 predictDisplacement_alpha_f(Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
void correctVelocity(Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, 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 Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
StepperHHTAlpha()
Default Constructor – not allowed.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Thyra Base interface for implicit time steppers.
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl)
void correctDisplacement(Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
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...
virtual std::string description() const
Teuchos::RCP< Teuchos::FancyOStream > out_
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< Teuchos::ParameterList > getDefaultParameters() const
virtual void setInitialGuess(Teuchos::RCP< const Thyra::VectorBase< Scalar > > initial_guess)
Pass initial guess to Newton solver (only relevant for implicit solvers)
void correctAcceleration(Thyra::VectorBase< Scalar > &a_n_plus1, const Thyra::VectorBase< Scalar > &a_n) const
Teuchos::RCP< const Thyra::VectorBase< Scalar > > initial_guess_
virtual void setObserver(Teuchos::RCP< StepperObserver< Scalar > > obs=Teuchos::null)
Set Observer.
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
virtual void initialize()
Initialize during construction and after changing input parameters.
void predictVelocity_alpha_f(Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v) const
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel)
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.