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 
73 template <class Scalar>
74 class StepperNewmarkImplicitDForm : virtual public Tempus::StepperImplicit<Scalar> {
75  public:
76 
83 
86  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
88  bool useFSAL,
89  std::string ICConsistency,
90  bool ICConsistencyCheck,
91  bool zeroInitialGuess,
92  std::string schemeName,
93  Scalar beta,
94  Scalar gamma,
96 
98 
99  virtual void
100  setModel(const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel);
101 
103  { return stepperNewmarkImpAppAction_; }
104 
105  virtual void setAppAction(
107 
109  virtual void setInitialConditions (
110  const Teuchos::RCP<SolutionHistory<Scalar> >& /* solutionHistory */){}
111 
113  virtual void
114  takeStep(const Teuchos::RCP<SolutionHistory<Scalar>>& solutionHistory);
115 
119  virtual Scalar
120  getOrder() const {
121  if (gamma_ == 0.5)
122  return 2.0;
123  else
124  return 1.0;
125  }
126  virtual Scalar
127  getOrderMin() const {
128  return 1.0;
129  }
130  virtual Scalar
131  getOrderMax() const {
132  return 2.0;
133  }
134  virtual bool isExplicit() const {return false;}
135  virtual bool isImplicit() const {return true;}
136  virtual bool isExplicitImplicit() const
137  {return isExplicit() && isImplicit();}
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  { return Scalar(1.0)/(beta_*dt*dt); }
147  virtual Scalar getAlpha(const Scalar dt) const { return gamma_/(beta_*dt); }
149  virtual Scalar getBeta (const Scalar ) const { return Scalar(1.0); }
150 
152 
154 
155  virtual void describe(Teuchos::FancyOStream& out,
156  const Teuchos::EVerbosityLevel verbLevel) const;
158 
159  virtual bool isValidSetup(Teuchos::FancyOStream & out) const;
160 
161  void
164  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const;
165 
166  void
170  const Scalar dt) const;
171 
172  void
175  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const;
176 
177  void
180  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const;
181 
182  void
185  const Thyra::VectorBase<Scalar>& d, const Scalar dt) const;
186 
187 
188  void setSchemeName(std::string schemeName);
189  void setBeta(Scalar beta);
190  void setGamma(Scalar gamma);
191 
192  protected:
193 
194  std::string schemeName_;
195  Scalar beta_;
196  Scalar gamma_;
197 
200 
201 };
202 
203 
205 // ------------------------------------------------------------------------
206 template<class Scalar>
209  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& model,
211 
212 
213 } // namespace Tempus
214 
215 #endif // Tempus_StepperNewmarkImplicitDForm_decl_hpp
Teuchos::RCP< StepperNewmarkImplicitDForm< Scalar > > createStepperNewmarkImplicitDForm(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, Teuchos::RCP< Teuchos::ParameterList > pl)
Nonmember constructor - ModelEvaluator and ParameterList.
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 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.
virtual void setAppAction(Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar > > appAction)
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.
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