Tempus Version of the Day
Time Integration
|
Newmark time stepper. More...
#include <Tempus_StepperNewmarkImplicitDForm_decl.hpp>
Public Member Functions | |
StepperNewmarkImplicitDForm () | |
Default constructor. | |
StepperNewmarkImplicitDForm (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool zeroInitialGuess, std::string schemeName, Scalar beta, Scalar gamma, const Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar > > &stepperAppAction) | |
Constructor. | |
![]() | |
Teuchos::RCP< Teuchos::ParameterList > | getValidParametersBasicImplicit () const |
void | setStepperImplicitValues (Teuchos::RCP< Teuchos::ParameterList > pl) |
Set StepperImplicit member data from the ParameterList. | |
void | setStepperSolverValues (Teuchos::RCP< Teuchos::ParameterList > pl) |
Set solver from ParameterList. | |
void | setSolverName (std::string i) |
Set the Solver Name. | |
std::string | getSolverName () const |
Get the Solver Name. | |
virtual Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () const override |
virtual Teuchos::RCP< const WrapperModelEvaluator< Scalar > > | getWrapperModel () |
virtual void | setDefaultSolver () |
virtual void | setSolver (Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > solver) override |
Set solver. | |
virtual Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > | getSolver () const override |
Get solver. | |
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, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &y=Teuchos::null, const int index=0) |
Solve implicit ODE, f(x, xDot, t, p) = 0. | |
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). | |
virtual void | setInitialGuess (Teuchos::RCP< const Thyra::VectorBase< Scalar > > initialGuess) override |
Pass initial guess to Newton solver (only relevant for implicit solvers) | |
virtual void | setZeroInitialGuess (bool zIG) |
Set parameter so that the initial guess is set to zero (=True) or use last timestep (=False). | |
virtual bool | getZeroInitialGuess () const |
virtual Scalar | getInitTimeStep (const Teuchos::RCP< SolutionHistory< Scalar > > &) const override |
![]() | |
virtual std::string | description () const |
void | setStepperValues (const Teuchos::RCP< Teuchos::ParameterList > pl) |
Set Stepper member data from ParameterList. | |
Teuchos::RCP< Teuchos::ParameterList > | getValidParametersBasic () const |
Add basic parameters to Steppers ParameterList. | |
virtual void | initialize () |
Initialize after construction and changing input parameters. | |
virtual bool | isInitialized () |
True if stepper's member data is initialized. | |
virtual void | checkInitialized () |
Check initialization, and error out on failure. | |
void | setStepperName (std::string s) |
Set the stepper name. | |
std::string | getStepperName () const |
Get the stepper name. | |
std::string | getStepperType () const |
Get the stepper type. The stepper type is used as an identifier for the stepper, and can only be set by the derived Stepper class. | |
virtual void | setUseFSAL (bool a) |
void | setUseFSALTrueOnly (bool a) |
void | setUseFSALFalseOnly (bool a) |
bool | getUseFSAL () const |
void | setICConsistency (std::string s) |
std::string | getICConsistency () const |
void | setICConsistencyCheck (bool c) |
bool | getICConsistencyCheck () const |
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > | getStepperX () |
Get Stepper x. | |
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > | getStepperXDot () |
Get Stepper xDot. | |
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > | getStepperXDotDot () |
Get Stepper xDotDot. | |
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > | getStepperXDotDot (Teuchos::RCP< SolutionState< Scalar > > state) |
Get xDotDot from SolutionState or Stepper storage. | |
Overridden from Teuchos::Describable | |
std::string | schemeName_ |
Scalar | beta_ |
Scalar | gamma_ |
Teuchos::RCP< Teuchos::FancyOStream > | out_ |
Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar > > | stepperNewmarkImpAppAction_ |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
virtual bool | isValidSetup (Teuchos::FancyOStream &out) 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 | 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 | correctAcceleration (Thyra::VectorBase< Scalar > &a, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d, const Scalar dt) const |
void | setSchemeName (std::string schemeName) |
void | setBeta (Scalar beta) |
void | setGamma (Scalar gamma) |
Basic stepper methods | |
virtual void | setModel (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel) |
Set the model. | |
virtual Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar > > | getAppAction () const |
virtual void | setAppAction (Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar > > appAction) |
virtual void | setInitialConditions (const Teuchos::RCP< SolutionHistory< Scalar > > &) |
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. | |
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > | getDefaultStepperState () |
Get a default (initial) StepperState. | |
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 |
virtual Scalar | getW_xDotDot_coeff (const Scalar dt) const |
Return W_xDotxDot_coeff = d(xDotDot)/d(x). | |
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). | |
Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const |
Additional Inherited Members | |
![]() | |
virtual void | setStepperX (Teuchos::RCP< Thyra::VectorBase< Scalar > > x) |
Set x for Stepper storage. | |
virtual void | setStepperXDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xDot) |
Set xDot for Stepper storage. | |
virtual void | setStepperXDotDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xDotDot) |
Set x for Stepper storage. | |
void | setStepperType (std::string s) |
Set the stepper type. | |
![]() | |
Teuchos::RCP< WrapperModelEvaluator< Scalar > > | wrapperModel_ |
Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > | solver_ |
Teuchos::RCP< const Thyra::VectorBase< Scalar > > | initialGuess_ |
bool | zeroInitialGuess_ |
std::string | solverName_ |
![]() | |
bool | useFSAL_ = false |
Use First-Same-As-Last (FSAL) principle. | |
bool | isInitialized_ = false |
True if stepper's member data is initialized. | |
Newmark time stepper.
Here, we implement the Newmark scheme in displacement predictor/corrector form; see equations (34)-(35) in: A. Mota, W. Klug, M. Ortiz, "Finite element simulation of firearm injury to the human cranium", Computational Mechanics 31(1) 115-121 (2003).
Newmark has two parameters:
Newmark is second order accurate if
Algorithm The algorithm for the Newmark implicit D-form with predictors and correctors is
The First-Same-As-Last (FSAL) principle is not used with the Newmark implicit D-Form method. The default is to set useFSAL=false, however useFSAL=true will also work but have no affect (i.e., no-op).
Definition at line 74 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
Tempus::StepperNewmarkImplicitDForm< Scalar >::StepperNewmarkImplicitDForm |
Default constructor.
Requires subsequent setModel(), setSolver() and initialize() calls before calling takeStep().
Definition at line 172 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
Tempus::StepperNewmarkImplicitDForm< Scalar >::StepperNewmarkImplicitDForm | ( | const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > & | appModel, |
const Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > & | solver, | ||
bool | useFSAL, | ||
std::string | ICConsistency, | ||
bool | ICConsistencyCheck, | ||
bool | zeroInitialGuess, | ||
std::string | schemeName, | ||
Scalar | beta, | ||
Scalar | gamma, | ||
const Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar > > & | stepperAppAction | ||
) |
Constructor.
Definition at line 191 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
|
virtual |
Set the model.
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 224 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
|
inlinevirtual |
Definition at line 102 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
virtual |
Definition at line 244 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
|
inlinevirtual |
Set the initial conditions and make them consistent.
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 109 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
virtual |
Take the specified timestep, dt, and return true if successful.
Implements Tempus::Stepper< Scalar >.
Definition at line 262 of file Tempus_StepperNewmarkImplicitDForm_impl.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 439 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 120 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 127 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 131 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 134 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 135 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 136 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 138 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 139 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 140 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Return W_xDotxDot_coeff = d(xDotDot)/d(x).
Definition at line 144 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Return alpha = d(xDot)/d(x).
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 147 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
inlinevirtual |
Return beta = d(x)/d(x).
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 149 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
virtual |
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 500 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
|
virtual |
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 450 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
|
virtual |
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 471 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
void Tempus::StepperNewmarkImplicitDForm< Scalar >::predictVelocity | ( | Thyra::VectorBase< Scalar > & | vPred, |
const Thyra::VectorBase< Scalar > & | v, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 23 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
void Tempus::StepperNewmarkImplicitDForm< 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 35 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
void Tempus::StepperNewmarkImplicitDForm< Scalar >::correctVelocity | ( | Thyra::VectorBase< Scalar > & | v, |
const Thyra::VectorBase< Scalar > & | vPred, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 53 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
void Tempus::StepperNewmarkImplicitDForm< Scalar >::correctDisplacement | ( | Thyra::VectorBase< Scalar > & | d, |
const Thyra::VectorBase< Scalar > & | dPred, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 65 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
void Tempus::StepperNewmarkImplicitDForm< Scalar >::correctAcceleration | ( | Thyra::VectorBase< Scalar > & | a, |
const Thyra::VectorBase< Scalar > & | dPred, | ||
const Thyra::VectorBase< Scalar > & | d, | ||
const Scalar | dt | ||
) | const |
Definition at line 77 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
void Tempus::StepperNewmarkImplicitDForm< Scalar >::setSchemeName | ( | std::string | schemeName | ) |
Definition at line 140 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
void Tempus::StepperNewmarkImplicitDForm< Scalar >::setBeta | ( | Scalar | beta | ) |
Definition at line 89 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
void Tempus::StepperNewmarkImplicitDForm< Scalar >::setGamma | ( | Scalar | gamma | ) |
Definition at line 120 of file Tempus_StepperNewmarkImplicitDForm_impl.hpp.
|
protected |
Definition at line 194 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
protected |
Definition at line 195 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
protected |
Definition at line 196 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
protected |
Definition at line 198 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.
|
protected |
Definition at line 199 of file Tempus_StepperNewmarkImplicitDForm_decl.hpp.