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: 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_config.hpp"
13 #include "Tempus_StepperImplicit.hpp"
14 #include "Tempus_WrapperModelEvaluatorSecondOrder.hpp"
16 
17 namespace Tempus {
18 
76 template <class Scalar>
78  : virtual public Tempus::StepperImplicit<Scalar> {
79  public:
86 
89  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel,
91  bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck,
92  bool zeroInitialGuess, std::string schemeName, Scalar beta, Scalar gamma,
94  stepperAppAction);
95 
97 
98  virtual void setModel(
99  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel);
100 
102  getAppAction() const
103  {
105  }
106 
107  virtual void setAppAction(
109 
111  virtual void setInitialConditions(
112  const Teuchos::RCP<SolutionHistory<Scalar>>& /* solutionHistory */)
113  {
114  }
115 
117  virtual void takeStep(
118  const Teuchos::RCP<SolutionHistory<Scalar>>& solutionHistory);
119 
122  virtual Scalar getOrder() const
123  {
124  if (gamma_ == 0.5)
125  return 2.0;
126  else
127  return 1.0;
128  }
129  virtual Scalar getOrderMin() const { return 1.0; }
130  virtual Scalar getOrderMax() const { return 2.0; }
131  virtual bool isExplicit() const { return false; }
132  virtual bool isImplicit() const { return true; }
133  virtual bool isExplicitImplicit() const
134  {
135  return isExplicit() && isImplicit();
136  }
137  virtual bool isOneStepMethod() const { return true; }
138  virtual bool isMultiStepMethod() const { return !isOneStepMethod(); }
139  virtual OrderODE getOrderODE() const { return SECOND_ORDER_ODE; }
141 
143  virtual Scalar getW_xDotDot_coeff(const Scalar dt) const
144  {
145  return Scalar(1.0) / (beta_ * dt * dt);
146  }
148  virtual Scalar getAlpha(const Scalar dt) const
149  {
150  return gamma_ / (beta_ * dt);
151  }
153  virtual Scalar getBeta(const Scalar) const { return Scalar(1.0); }
154 
156 
158 
159  virtual void describe(Teuchos::FancyOStream& out,
160  const Teuchos::EVerbosityLevel verbLevel) const;
162 
163  virtual bool isValidSetup(Teuchos::FancyOStream& out) const;
164 
166  const Thyra::VectorBase<Scalar>& v,
167  const Thyra::VectorBase<Scalar>& a,
168  const Scalar dt) const;
169 
171  const Thyra::VectorBase<Scalar>& d,
172  const Thyra::VectorBase<Scalar>& v,
173  const Thyra::VectorBase<Scalar>& a,
174  const Scalar dt) const;
175 
177  const Thyra::VectorBase<Scalar>& vPred,
178  const Thyra::VectorBase<Scalar>& a,
179  const Scalar dt) const;
180 
182  const Thyra::VectorBase<Scalar>& dPred,
183  const Thyra::VectorBase<Scalar>& a,
184  const Scalar dt) const;
185 
187  const Thyra::VectorBase<Scalar>& dPred,
188  const Thyra::VectorBase<Scalar>& d,
189  const Scalar dt) const;
190 
191  void setSchemeName(std::string schemeName);
192  void setBeta(Scalar beta);
193  void setGamma(Scalar gamma);
194 
195  protected:
196  std::string schemeName_;
197  Scalar beta_;
198  Scalar gamma_;
199 
203 };
204 
206 // ------------------------------------------------------------------------
207 template <class Scalar>
210  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& model,
212 
213 } // namespace Tempus
214 
215 #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