Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_StepperNewmarkImplicitDForm_decl.hpp
Go to the documentation of this file.
1 //@HEADER
2 // *****************************************************************************
3 // Tempus: Time Integration and Sensitivity Analysis Package
4 //
5 // Copyright 2017 NTESS and the Tempus contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 //@HEADER
9 
10 #ifndef Tempus_StepperNewmarkImplicitDForm_decl_hpp
11 #define Tempus_StepperNewmarkImplicitDForm_decl_hpp
12 
13 #include "Tempus_config.hpp"
14 #include "Tempus_StepperImplicit.hpp"
15 #include "Tempus_WrapperModelEvaluatorSecondOrder.hpp"
17 
18 namespace Tempus {
19 
77 template <class Scalar>
79  : virtual public Tempus::StepperImplicit<Scalar> {
80  public:
87 
90  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel,
92  bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck,
93  bool zeroInitialGuess, std::string schemeName, Scalar beta, Scalar gamma,
95  stepperAppAction);
96 
98 
99  virtual void setModel(
100  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel);
101 
103  getAppAction() const
104  {
106  }
107 
108  virtual void setAppAction(
110 
112  virtual void setInitialConditions(
113  const Teuchos::RCP<SolutionHistory<Scalar>>& /* solutionHistory */)
114  {
115  }
116 
118  virtual void takeStep(
119  const Teuchos::RCP<SolutionHistory<Scalar>>& solutionHistory);
120 
123  virtual Scalar getOrder() const
124  {
125  if (gamma_ == 0.5)
126  return 2.0;
127  else
128  return 1.0;
129  }
130  virtual Scalar getOrderMin() const { return 1.0; }
131  virtual Scalar getOrderMax() const { return 2.0; }
132  virtual bool isExplicit() const { return false; }
133  virtual bool isImplicit() const { return true; }
134  virtual bool isExplicitImplicit() const
135  {
136  return isExplicit() && isImplicit();
137  }
138  virtual bool isOneStepMethod() const { return true; }
139  virtual bool isMultiStepMethod() const { return !isOneStepMethod(); }
140  virtual OrderODE getOrderODE() const { return SECOND_ORDER_ODE; }
142 
144  virtual Scalar getW_xDotDot_coeff(const Scalar dt) const
145  {
146  return Scalar(1.0) / (beta_ * dt * dt);
147  }
149  virtual Scalar getAlpha(const Scalar dt) const
150  {
151  return gamma_ / (beta_ * dt);
152  }
154  virtual Scalar getBeta(const Scalar) const { return Scalar(1.0); }
155 
157 
159 
160  virtual void describe(Teuchos::FancyOStream& out,
161  const Teuchos::EVerbosityLevel verbLevel) const;
163 
164  virtual bool isValidSetup(Teuchos::FancyOStream& out) const;
165 
167  const Thyra::VectorBase<Scalar>& v,
168  const Thyra::VectorBase<Scalar>& a,
169  const Scalar dt) const;
170 
172  const Thyra::VectorBase<Scalar>& d,
173  const Thyra::VectorBase<Scalar>& v,
174  const Thyra::VectorBase<Scalar>& a,
175  const Scalar dt) const;
176 
178  const Thyra::VectorBase<Scalar>& vPred,
179  const Thyra::VectorBase<Scalar>& a,
180  const Scalar dt) const;
181 
183  const Thyra::VectorBase<Scalar>& dPred,
184  const Thyra::VectorBase<Scalar>& a,
185  const Scalar dt) const;
186 
188  const Thyra::VectorBase<Scalar>& dPred,
189  const Thyra::VectorBase<Scalar>& d,
190  const Scalar dt) const;
191 
192  void setSchemeName(std::string schemeName);
193  void setBeta(Scalar beta);
194  void setGamma(Scalar gamma);
195 
196  protected:
197  std::string schemeName_;
198  Scalar beta_;
199  Scalar gamma_;
200 
204 };
205 
207 // ------------------------------------------------------------------------
208 template <class Scalar>
211  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& model,
213 
214 } // namespace Tempus
215 
216 #endif // Tempus_StepperNewmarkImplicitDForm_decl_hpp
void correctVelocity(Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
Teuchos::RCP< StepperNewmarkImplicitDForm< Scalar > > createStepperNewmarkImplicitDForm(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar >> &model, Teuchos::RCP< Teuchos::ParameterList > pl)
Nonmember constructor - ModelEvaluator and ParameterList.
void correctDisplacement(Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
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 setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar >> &)
Set the initial conditions and make them consistent.
virtual void setAppAction(Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar >> appAction)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Application Action for StepperNewmarkImplicitDForm.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
virtual Scalar getW_xDotDot_coeff(const Scalar dt) const
Return W_xDotxDot_coeff = d(xDotDot)/d(x).
virtual bool isValidSetup(Teuchos::FancyOStream &out) const
Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar > > stepperNewmarkImpAppAction_
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) 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
virtual Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar > > getAppAction() const