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: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #ifndef Tempus_StepperNewmarkImplicitAForm_decl_hpp
10 #define Tempus_StepperNewmarkImplicitAForm_decl_hpp
11 
12 #include "Tempus_config.hpp"
13 #include "Tempus_StepperImplicit.hpp"
14 #include "Tempus_WrapperModelEvaluatorSecondOrder.hpp"
16 
17 namespace Tempus {
18 
19 
91 template<class Scalar>
93  : virtual public Tempus::StepperImplicit<Scalar>
94 {
95 public:
96 
103 
106  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
108  bool useFSAL,
109  std::string ICConsistency,
110  bool ICConsistencyCheck,
111  bool zeroInitialGuess,
112  std::string schemeName,
113  Scalar beta,
114  Scalar gamma,
116 
118 
119  virtual void setModel(
120  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel);
121 
122  virtual void setAppAction(
124 
126  { return stepperNewmarkImpAppAction_; }
127 
129  virtual void setInitialConditions (
130  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
131 
133  virtual void takeStep(
134  const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory);
135 
139  virtual Scalar getOrder() const {
140  if (gamma_ == 0.5) return 2.0;
141  else return 1.0;
142  }
143  virtual Scalar getOrderMin() const {return 1.0;}
144  virtual Scalar getOrderMax() const {return 2.0;}
145 
146  virtual bool isExplicit() const {return false;}
147  virtual bool isImplicit() const {return true;}
148  virtual bool isExplicitImplicit() const
149  {return isExplicit() && isImplicit();}
150  virtual bool isOneStepMethod() const {return true;}
151  virtual bool isMultiStepMethod() const {return !isOneStepMethod();}
152  virtual void setUseFSAL(bool a) { this->setUseFSALTrueOnly(a); }
153  virtual OrderODE getOrderODE() const {return SECOND_ORDER_ODE;}
155 
157  virtual Scalar getW_xDotDot_coeff (const Scalar) const {return Scalar(1.0);}
159  virtual Scalar getAlpha(const Scalar dt) const { return gamma_*dt; }
161  virtual Scalar getBeta (const Scalar dt) const { return beta_*dt*dt; }
162 
164 
166 
167  virtual void describe(Teuchos::FancyOStream & out,
168  const Teuchos::EVerbosityLevel verbLevel) const;
170 
171  virtual bool isValidSetup(Teuchos::FancyOStream & out) const;
172 
174  const Thyra::VectorBase<Scalar>& v,
175  const Thyra::VectorBase<Scalar>& a,
176  const Scalar dt) const;
177 
179  const Thyra::VectorBase<Scalar>& d,
180  const Thyra::VectorBase<Scalar>& v,
181  const Thyra::VectorBase<Scalar>& a,
182  const Scalar dt) const;
183 
185  const Thyra::VectorBase<Scalar>& vPred,
186  const Thyra::VectorBase<Scalar>& a,
187  const Scalar dt) const;
188 
190  const Thyra::VectorBase<Scalar>& dPred,
191  const Thyra::VectorBase<Scalar>& a,
192  const Scalar dt) const;
193 
194  void setSchemeName(std::string schemeName);
195  void setBeta(Scalar beta);
196  void setGamma(Scalar gamma);
197 
198 private:
199 
200  std::string schemeName_;
201  Scalar beta_;
202  Scalar gamma_;
203 
206 
207 };
208 
209 
211 // ------------------------------------------------------------------------
212 template<class Scalar>
215  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& model,
217 
218 
219 } // namespace Tempus
220 
221 #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)
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).