Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_StepperNewmarkImplicitAForm_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_StepperNewmarkImplicitAForm_decl_hpp
11 #define Tempus_StepperNewmarkImplicitAForm_decl_hpp
12 
13 #include "Tempus_config.hpp"
14 #include "Tempus_StepperImplicit.hpp"
15 #include "Tempus_WrapperModelEvaluatorSecondOrder.hpp"
17 
18 namespace Tempus {
19 
91 template <class Scalar>
93  : virtual public Tempus::StepperImplicit<Scalar> {
94  public:
101 
104  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
106  bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck,
107  bool zeroInitialGuess, std::string schemeName, Scalar beta, Scalar gamma,
109  stepperAppAction);
110 
112 
113  virtual void setModel(
114  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel);
115 
116  virtual void setAppAction(
118 
120  getAppAction() const
121  {
123  }
124 
126  virtual void setInitialConditions(
127  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
128 
130  virtual void takeStep(
131  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
132 
135  virtual Scalar getOrder() const
136  {
137  if (gamma_ == 0.5)
138  return 2.0;
139  else
140  return 1.0;
141  }
142  virtual Scalar getOrderMin() const { return 1.0; }
143  virtual Scalar getOrderMax() const { return 2.0; }
144 
145  virtual bool isExplicit() const { return false; }
146  virtual bool isImplicit() const { return true; }
147  virtual bool isExplicitImplicit() const
148  {
149  return isExplicit() && isImplicit();
150  }
151  virtual bool isOneStepMethod() const { return true; }
152  virtual bool isMultiStepMethod() const { return !isOneStepMethod(); }
153  virtual void setUseFSAL(bool a) { this->setUseFSALTrueOnly(a); }
154  virtual OrderODE getOrderODE() const { return SECOND_ORDER_ODE; }
156 
158  virtual Scalar getW_xDotDot_coeff(const Scalar) const { return Scalar(1.0); }
160  virtual Scalar getAlpha(const Scalar dt) const { return gamma_ * dt; }
162  virtual Scalar getBeta(const Scalar dt) const { return beta_ * dt * dt; }
163 
165 
167 
168  virtual void describe(Teuchos::FancyOStream& out,
169  const Teuchos::EVerbosityLevel verbLevel) const;
171 
172  virtual bool isValidSetup(Teuchos::FancyOStream& out) const;
173 
175  const Thyra::VectorBase<Scalar>& v,
176  const Thyra::VectorBase<Scalar>& a,
177  const Scalar dt) const;
178 
180  const Thyra::VectorBase<Scalar>& d,
181  const Thyra::VectorBase<Scalar>& v,
182  const Thyra::VectorBase<Scalar>& a,
183  const Scalar dt) const;
184 
186  const Thyra::VectorBase<Scalar>& vPred,
187  const Thyra::VectorBase<Scalar>& a,
188  const Scalar dt) const;
189 
191  const Thyra::VectorBase<Scalar>& dPred,
192  const Thyra::VectorBase<Scalar>& a,
193  const Scalar dt) const;
194 
195  void setSchemeName(std::string schemeName);
196  void setBeta(Scalar beta);
197  void setGamma(Scalar gamma);
198 
199  private:
200  std::string schemeName_;
201  Scalar beta_;
202  Scalar gamma_;
203 
207 };
208 
210 // ------------------------------------------------------------------------
211 template <class Scalar>
214  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& model,
216 
217 } // namespace Tempus
218 
219 #endif // Tempus_StepperNewmarkImplicitAForm_decl_hpp
virtual Teuchos::RCP< StepperNewmarkImplicitAFormAppAction< Scalar > > getAppAction() const
Application Action for StepperNewmarkImplicitAForm.
Teuchos::RCP< StepperNewmarkImplicitAForm< Scalar > > createStepperNewmarkImplicitAForm(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, Teuchos::RCP< Teuchos::ParameterList > pl)
Nonmember constructor - ModelEvaluator and ParameterList.
virtual bool isValidSetup(Teuchos::FancyOStream &out) const
Teuchos::RCP< StepperNewmarkImplicitAFormAppAction< Scalar > > stepperNewmarkImpAppAction_
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).
virtual Scalar getBeta(const Scalar dt) const
Return beta = d(x)/d(xDotDot).
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual void setAppAction(Teuchos::RCP< StepperNewmarkImplicitAFormAppAction< Scalar > > appAction)
Stepper integrates second-order ODEs.
Thyra Base interface for implicit time steppers.
void setUseFSALTrueOnly(bool a)
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions and make them consistent.
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)
Set the model.
void correctDisplacement(Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
void correctVelocity(Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Newmark time stepper in acceleration form (a-form).