9 #ifndef Tempus_StepperHHTAlpha_decl_hpp
10 #define Tempus_StepperHHTAlpha_decl_hpp
12 #include "Tempus_config.hpp"
13 #include "Tempus_StepperImplicit.hpp"
14 #include "Tempus_WrapperModelEvaluatorSecondOrder.hpp"
45 template<
class Scalar>
62 std::string ICConsistency,
63 bool ICConsistencyCheck,
64 bool zeroInitialGuess,
65 std::string schemeName,
94 if (
gamma_ == 0.5)
return 2.0;
111 {
return Scalar(1.0)/(
beta_*dt*dt); }
115 virtual Scalar
getBeta (
const Scalar )
const {
return Scalar(1.0); }
130 const Scalar dt)
const;
136 const Scalar dt)
const;
150 const Scalar dt)
const;
155 const Scalar dt)
const;
180 template<
class Scalar>
189 #endif // Tempus_StepperHHTAlpha_decl_hpp
void setBeta(Scalar beta)
void correctDisplacement(Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
virtual Scalar getW_xDotDot_coeff(const Scalar dt) const
Return W_xDotxDot_coeff = d(xDotDot)/d(x).
virtual Teuchos::RCP< StepperHHTAlphaAppAction< Scalar > > getAppAction() const
virtual bool isExplicitImplicit() const
virtual bool isValidSetup(Teuchos::FancyOStream &out) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
void correctVelocity(Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
virtual OrderODE getOrderODE() const
virtual Scalar getOrderMin() const
Application Action for HHT Alpha.
void correctAcceleration(Thyra::VectorBase< Scalar > &a_n_plus1, const Thyra::VectorBase< Scalar > &a_n) const
void setAlphaF(Scalar alpha_f)
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
StepperHHTAlpha()
Default constructor.
Stepper integrates second-order ODEs.
virtual Scalar getAlpha(const Scalar dt) const
Return alpha = d(xDot)/d(x).
Thyra Base interface for implicit time steppers.
void setAlphaM(Scalar alpha_m)
void predictVelocity_alpha_f(Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v) const
void setSchemeName(std::string schemeName)
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
void predictDisplacement_alpha_f(Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d) const
virtual Scalar getBeta(const Scalar) const
Return beta = d(x)/d(x).
virtual bool isOneStepMethod() const
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 Scalar getOrderMax() const
Teuchos::RCP< StepperHHTAlpha< Scalar > > createStepperHHTAlpha(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, Teuchos::RCP< Teuchos::ParameterList > pl)
Nonmember constructor - ModelEvaluator and ParameterList.
Teuchos::RCP< Teuchos::FancyOStream > out_
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 bool isImplicit() const
virtual void setAppAction(Teuchos::RCP< StepperHHTAlphaAppAction< Scalar > > appAction)
void setGamma(Scalar gamma)
virtual bool isMultiStepMethod() const
virtual Scalar getOrder() const
virtual bool isExplicit() const
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel)
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.
Teuchos::RCP< StepperHHTAlphaAppAction< Scalar > > stepperHHTAlphaAppAction_