Tempus
Version of the Day
Time Integration
|
HHT-Alpha time stepper. More...
#include <Tempus_StepperHHTAlpha_decl.hpp>
Public Member Functions | |
StepperHHTAlpha () | |
Default constructor. More... | |
StepperHHTAlpha (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< StepperObserver< Scalar > > &obs, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool zeroInitialGuess, std::string schemeName, Scalar beta, Scalar gamma, Scalar alpha_f_, Scalar alpha_m_) | |
Constructor. More... | |
virtual Scalar | getW_xDotDot_coeff (const Scalar dt) const |
Return W_xDotxDot_coeff = d(xDotDot)/d(x). More... | |
virtual Scalar | getAlpha (const Scalar dt) const |
Return alpha = d(xDot)/d(x). More... | |
virtual Scalar | getBeta (const Scalar) const |
Return beta = d(x)/d(x). More... | |
Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const |
void | predictVelocity (Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) 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 |
void | predictVelocity_alpha_f (Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v) const |
void | predictDisplacement_alpha_f (Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d) const |
void | correctAcceleration (Thyra::VectorBase< Scalar > &a_n_plus1, const Thyra::VectorBase< Scalar > &a_n) const |
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 |
void | setSchemeName (std::string schemeName) |
void | setBeta (Scalar beta) |
void | setGamma (Scalar gamma) |
void | setAlphaF (Scalar alpha_f) |
void | setAlphaM (Scalar alpha_m) |
Basic stepper methods | |
virtual void | setModel (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel) |
virtual void | setObserver (Teuchos::RCP< StepperObserver< Scalar > >=Teuchos::null) |
Set Observer. More... | |
virtual Teuchos::RCP < StepperObserver< Scalar > > | getObserver () const |
Get Observer. More... | |
virtual void | initialize () |
Initialize during construction and after changing input parameters. More... | |
virtual void | setInitialConditions (const Teuchos::RCP< SolutionHistory< Scalar > > &) |
Set the initial conditions and make them consistent. More... | |
virtual void | takeStep (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) |
Take the specified timestep, dt, and return true if successful. More... | |
virtual Teuchos::RCP < Tempus::StepperState< Scalar > > | getDefaultStepperState () |
Get a default (initial) StepperState. More... | |
virtual Scalar | getOrder () const |
virtual Scalar | getOrderMin () const |
virtual Scalar | getOrderMax () const |
virtual bool | isExplicit () const |
virtual bool | isImplicit () const |
virtual bool | isExplicitImplicit () const |
virtual bool | isOneStepMethod () const |
virtual bool | isMultiStepMethod () const |
virtual OrderODE | getOrderODE () const |
Overridden from Teuchos::Describable | |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
Public Member Functions inherited from Tempus::StepperImplicit< Scalar > | |
virtual void | setNonConstModel (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &appModel) |
virtual Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () |
virtual Teuchos::RCP< const WrapperModelEvaluator< Scalar > > | getWrapperModel () |
virtual void | setSolver (Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > solver=Teuchos::null) |
Set solver. More... | |
virtual Teuchos::RCP < Thyra::NonlinearSolverBase < Scalar > > | getSolver () const |
Get solver. More... | |
const Thyra::SolveStatus< Scalar > | solveImplicitODE (const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x) |
Solve problem using x in-place. (Needs to be deprecated!) More... | |
const Thyra::SolveStatus< Scalar > | solveImplicitODE (const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xDot, const Scalar time, const Teuchos::RCP< ImplicitODEParameters< Scalar > > &p) |
Solve implicit ODE, f(x, xDot, t, p) = 0. More... | |
void | evaluateImplicitODE (Teuchos::RCP< Thyra::VectorBase< Scalar > > &f, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xDot, const Scalar time, const Teuchos::RCP< ImplicitODEParameters< Scalar > > &p) |
Evaluate implicit ODE residual, f(x, xDot, t, p). More... | |
virtual void | setInitialGuess (Teuchos::RCP< const Thyra::VectorBase< Scalar > > initial_guess) |
Pass initial guess to Newton solver (only relevant for implicit solvers) More... | |
virtual void | setZeroInitialGuess (bool zIG) |
Set parameter so that the initial guess is set to zero (=True) or use last timestep (=False). More... | |
virtual bool | getZeroInitialGuess () const |
virtual Scalar | getInitTimeStep (const Teuchos::RCP< SolutionHistory< Scalar > > &) const |
virtual void | setStepperXDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xDot) |
Set xDot for Stepper storage. More... | |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getStepperXDot (Teuchos::RCP< SolutionState< Scalar > > state) |
Get xDot from SolutionState or Stepper storage. More... | |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getStepperXDotDot (Teuchos::RCP< SolutionState< Scalar > > state) |
Get xDotDot from SolutionState or Stepper storage. More... | |
Public Member Functions inherited from Tempus::Stepper< Scalar > | |
void | setStepperType (std::string s) |
std::string | getStepperType () const |
void | setUseFSAL (bool a) |
bool | getUseFSAL () const |
virtual bool | getUseFSALDefault () const |
void | setICConsistency (std::string s) |
std::string | getICConsistency () const |
virtual std::string | getICConsistencyDefault () const |
void | setICConsistencyCheck (bool c) |
bool | getICConsistencyCheck () const |
virtual bool | getICConsistencyCheckDefault () const |
virtual std::string | description () const |
virtual void | createSubSteppers (std::vector< Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > >) |
Private Attributes | |
std::string | schemeName_ |
Scalar | beta_ |
Scalar | gamma_ |
Scalar | alpha_f_ |
Scalar | alpha_m_ |
Teuchos::RCP < Teuchos::FancyOStream > | out_ |
Additional Inherited Members | |
Protected Attributes inherited from Tempus::StepperImplicit< Scalar > | |
Teuchos::RCP < WrapperModelEvaluator < Scalar > > | wrapperModel_ |
Teuchos::RCP < Thyra::NonlinearSolverBase < Scalar > > | solver_ |
Teuchos::RCP< const Thyra::VectorBase< Scalar > > | initial_guess_ |
bool | zeroInitialGuess_ |
Teuchos::RCP< StepperObserver < Scalar > > | stepperObserver_ |
Teuchos::RCP < Thyra::VectorBase< Scalar > > | stepperXDot_ |
Teuchos::RCP < Thyra::VectorBase< Scalar > > | stepperXDotDot_ |
HHT-Alpha time stepper.
Here, we implement the HHT-Alpha scheme in predictor/corrector form; see equations (10) and (13)-(19) in: G.M. Hulbert, J. Chung, "Explicit time integration algorithms for structural dynamics with optimal numerical dissipation", Comput. Methods Appl. Mech. Engrg. 137 175-188 (1996).
There are four parameters in the scheme: , , and , all of which must be in the range . When , the scheme reduces to the Newmark Beta scheme (see Tempus::StepperNewmark for details). Like the Newmark Beta scheme, the HHT-Alpha scheme can be either first or second order accurate, and either explicit or implicit.
Although the general form of the scheme has been implemented in Tempus, it has only been verified for the case when (corresponding to the Newmark Beta) scheme, so other values for these parameters are not allowed at the present time. Also, note that, like the Newmark Beta stepper, the linear solve for the explicit version of this scheme has not been optimized (the mass matrix is not lumped).
The First-Step-As-Last (FSAL) principle is not used with the HHT-Alpha method.
Definition at line 44 of file Tempus_StepperHHTAlpha_decl.hpp.
Tempus::StepperHHTAlpha< Scalar >::StepperHHTAlpha | ( | ) |
Default constructor.
Requires subsequent setModel(), setSolver() and initialize() calls before calling takeStep().
Definition at line 233 of file Tempus_StepperHHTAlpha_impl.hpp.
Tempus::StepperHHTAlpha< Scalar >::StepperHHTAlpha | ( | const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > & | appModel, |
const Teuchos::RCP< StepperObserver< Scalar > > & | obs, | ||
const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > & | solver, | ||
bool | useFSAL, | ||
std::string | ICConsistency, | ||
bool | ICConsistencyCheck, | ||
bool | zeroInitialGuess, | ||
std::string | schemeName, | ||
Scalar | beta, | ||
Scalar | gamma, | ||
Scalar | alpha_f_, | ||
Scalar | alpha_m_ | ||
) |
Constructor.
Definition at line 254 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::correctAcceleration | ( | Thyra::VectorBase< Scalar > & | a_n_plus1, |
const Thyra::VectorBase< Scalar > & | a_n | ||
) | const |
Definition at line 88 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::correctDisplacement | ( | Thyra::VectorBase< Scalar > & | d, |
const Thyra::VectorBase< Scalar > & | dPred, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 117 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::correctVelocity | ( | Thyra::VectorBase< Scalar > & | v, |
const Thyra::VectorBase< Scalar > & | vPred, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 103 of file Tempus_StepperHHTAlpha_impl.hpp.
|
virtual |
Definition at line 456 of file Tempus_StepperHHTAlpha_impl.hpp.
|
inlinevirtual |
Return alpha = d(xDot)/d(x).
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 115 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Return beta = d(x)/d(x).
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 117 of file Tempus_StepperHHTAlpha_decl.hpp.
|
virtual |
Get a default (initial) StepperState.
Provide a StepperState to the SolutionState. This Stepper does not have any special state data, so just provide the base class StepperState with the Stepper description. This can be checked to ensure that the input StepperState can be used by this Stepper.
Implements Tempus::Stepper< Scalar >.
Definition at line 444 of file Tempus_StepperHHTAlpha_impl.hpp.
|
inlinevirtual |
Get Observer.
Implements Tempus::Stepper< Scalar >.
Definition at line 78 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 94 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 99 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 98 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 108 of file Tempus_StepperHHTAlpha_decl.hpp.
|
virtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 470 of file Tempus_StepperHHTAlpha_impl.hpp.
|
inlinevirtual |
Return W_xDotxDot_coeff = d(xDotDot)/d(x).
Definition at line 112 of file Tempus_StepperHHTAlpha_decl.hpp.
|
virtual |
Initialize during construction and after changing input parameters.
Implements Tempus::Stepper< Scalar >.
Definition at line 307 of file Tempus_StepperHHTAlpha_impl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 101 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 103 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 102 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 106 of file Tempus_StepperHHTAlpha_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 105 of file Tempus_StepperHHTAlpha_decl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::predictDisplacement | ( | Thyra::VectorBase< Scalar > & | dPred, |
const Thyra::VectorBase< Scalar > & | d, | ||
const Thyra::VectorBase< Scalar > & | v, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 42 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::predictDisplacement_alpha_f | ( | Thyra::VectorBase< Scalar > & | dPred, |
const Thyra::VectorBase< Scalar > & | d | ||
) | const |
Definition at line 76 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::predictVelocity | ( | Thyra::VectorBase< Scalar > & | vPred, |
const Thyra::VectorBase< Scalar > & | v, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 28 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::predictVelocity_alpha_f | ( | Thyra::VectorBase< Scalar > & | vPred, |
const Thyra::VectorBase< Scalar > & | v | ||
) | const |
Definition at line 63 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::setAlphaF | ( | Scalar | alpha_f | ) |
Definition at line 179 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::setAlphaM | ( | Scalar | alpha_m | ) |
Definition at line 191 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::setBeta | ( | Scalar | beta | ) |
Definition at line 132 of file Tempus_StepperHHTAlpha_impl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::setGamma | ( | Scalar | gamma | ) |
Definition at line 161 of file Tempus_StepperHHTAlpha_impl.hpp.
|
inlinevirtual |
Set the initial conditions and make them consistent.
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 85 of file Tempus_StepperHHTAlpha_decl.hpp.
|
virtual |
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 292 of file Tempus_StepperHHTAlpha_impl.hpp.
|
inlinevirtual |
Set Observer.
Implements Tempus::Stepper< Scalar >.
Definition at line 75 of file Tempus_StepperHHTAlpha_decl.hpp.
void Tempus::StepperHHTAlpha< Scalar >::setSchemeName | ( | std::string | schemeName | ) |
Definition at line 203 of file Tempus_StepperHHTAlpha_impl.hpp.
|
virtual |
Take the specified timestep, dt, and return true if successful.
Implements Tempus::Stepper< Scalar >.
Definition at line 321 of file Tempus_StepperHHTAlpha_impl.hpp.
|
private |
Definition at line 168 of file Tempus_StepperHHTAlpha_decl.hpp.
|
private |
Definition at line 169 of file Tempus_StepperHHTAlpha_decl.hpp.
|
private |
Definition at line 166 of file Tempus_StepperHHTAlpha_decl.hpp.
|
private |
Definition at line 167 of file Tempus_StepperHHTAlpha_decl.hpp.
|
private |
Definition at line 171 of file Tempus_StepperHHTAlpha_decl.hpp.
|
private |
Definition at line 165 of file Tempus_StepperHHTAlpha_decl.hpp.