Tempus
Version of the Day
Time Integration
|
Solution state for integrators and steppers. More...
#include <Tempus_SolutionState_decl.hpp>
Public Member Functions | |
SolutionState () | |
Default Constructor – Not meant for immediate adding to SolutionHistory. This constructor does not set the solution vectors, x, xdot and xdotdot. which should be set via setX(), setXDot(), and/or setXDotDot() prior to being added to SolutionHistory. More... | |
SolutionState (const Teuchos::RCP< SolutionStateMetaData< Scalar > > ssmd, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdot, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdotdot, const Teuchos::RCP< StepperState< Scalar > > &stepperState, const Teuchos::RCP< PhysicsState< Scalar > > &physicsState) | |
SolutionState (const Teuchos::RCP< const SolutionStateMetaData< Scalar > > ssmd, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &xdot, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &xdotdot, const Teuchos::RCP< const StepperState< Scalar > > &stepperState, const Teuchos::RCP< const PhysicsState< Scalar > > &physicsState) | |
SolutionState (const SolutionState< Scalar > &ss) | |
This is a shallow copy constructor, use clone for a deep copy constructor. More... | |
virtual Teuchos::RCP < SolutionState< Scalar > > | clone () const |
This is a deep copy constructor. More... | |
virtual void | copy (const Teuchos::RCP< const SolutionState< Scalar > > &ss) |
This is a deep copy. More... | |
virtual void | copySolutionData (const Teuchos::RCP< const SolutionState< Scalar > > &s) |
Deep copy solution data, but keep metaData untouched. More... | |
virtual | ~SolutionState () |
Destructor. More... | |
virtual void | computeNorms (const Teuchos::RCP< const SolutionState< Scalar > > &ssIn=Teuchos::null) |
Compute the solution norms, and solution change from ssIn, if provided. More... | |
![]() | |
void | describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
virtual | ~Describable () |
LabeledObject () | |
virtual | ~LabeledObject () |
virtual void | setObjectLabel (const std::string &objectLabel) |
virtual std::string | getObjectLabel () const |
DescribableStreamManipulatorState | describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default) |
std::ostream & | operator<< (std::ostream &os, const DescribableStreamManipulatorState &d) |
![]() | |
VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual const VerboseObject & | setVerbLevel (const EVerbosityLevel verbLevel) const |
virtual const VerboseObject & | setOverridingVerbLevel (const EVerbosityLevel verbLevel) const |
virtual EVerbosityLevel | getVerbLevel () const |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT RCP< const ParameterList > | getValidVerboseObjectSublist () |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | setupVerboseObjectSublist (ParameterList *paramList) |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | readVerboseObjectSublist (ParameterList *paramList, RCP< FancyOStream > *oStream, EVerbosityLevel *verbLevel) |
void | readVerboseObjectSublist (ParameterList *paramList, VerboseObject< ObjectType > *verboseObject) |
![]() | |
virtual | ~VerboseObjectBase () |
VerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual const VerboseObjectBase & | setOStream (const RCP< FancyOStream > &oStream) const |
virtual const VerboseObjectBase & | setOverridingOStream (const RCP< FancyOStream > &oStream) const |
virtual VerboseObjectBase & | setLinePrefix (const std::string &linePrefix) |
virtual RCP< FancyOStream > | getOStream () const |
virtual RCP< FancyOStream > | getOverridingOStream () const |
virtual std::string | getLinePrefix () const |
virtual OSTab | getOSTab (const int tabs=1, const std::string &linePrefix="") const |
Private Attributes | |
Teuchos::RCP< const SolutionStateMetaData< Scalar > > | metaData_ |
Meta Data for the solution state. More... | |
Teuchos::RCP < SolutionStateMetaData < Scalar > > | metaData_nc_ |
Teuchos::RCP< const Thyra::VectorBase< Scalar > > | x_ |
Solution. More... | |
Teuchos::RCP < Thyra::VectorBase< Scalar > > | x_nc_ |
Teuchos::RCP< const Thyra::VectorBase< Scalar > > | xdot_ |
Time derivative of the solution. More... | |
Teuchos::RCP < Thyra::VectorBase< Scalar > > | xdot_nc_ |
Teuchos::RCP< const Thyra::VectorBase< Scalar > > | xdotdot_ |
Second time derivative of the solution. More... | |
Teuchos::RCP < Thyra::VectorBase< Scalar > > | xdotdot_nc_ |
Teuchos::RCP< const StepperState< Scalar > > | stepperState_ |
StepperState for this SolutionState. More... | |
Teuchos::RCP< StepperState < Scalar > > | stepperState_nc_ |
Teuchos::RCP< const PhysicsState< Scalar > > | physicsState_ |
PhysicsState for this SolutionState. More... | |
Teuchos::RCP< PhysicsState < Scalar > > | physicsState_nc_ |
Get MetaData values | |
virtual Teuchos::RCP< const SolutionStateMetaData< Scalar > > | getMetaData () const |
virtual Teuchos::RCP < SolutionStateMetaData < Scalar > > | getMetaData () |
virtual Scalar | getTime () const |
virtual int | getIndex () const |
virtual Scalar | getTimeStep () const |
virtual Scalar | getErrorAbs () const |
virtual Scalar | getErrorRel () const |
virtual Scalar | getErrorRelNm1 () const |
virtual Scalar | getErrorRelNm2 () const |
virtual int | getOrder () const |
virtual int | getNFailures () const |
virtual int | getNRunningFailures () const |
virtual int | getNConsecutiveFailures () const |
virtual Scalar | getTolAbs () const |
virtual Scalar | getTolRel () const |
virtual Scalar | getXNormL2 () const |
virtual Scalar | getDxNormL2Abs () const |
virtual Scalar | getDxNormL2Rel () const |
virtual bool | getComputeNorms () const |
virtual Status | getSolutionStatus () const |
virtual bool | getOutput () const |
virtual bool | getOutputScreen () const |
virtual bool | getIsSynced () const |
virtual bool | getIsInterpolated () const |
virtual bool | getAccuracy () const |
Set MetaData values | |
virtual void | setMetaData (Teuchos::RCP< const SolutionStateMetaData< Scalar > > md) |
virtual void | setMetaData (Teuchos::RCP< SolutionStateMetaData< Scalar > > md) |
virtual void | setTime (Scalar time) |
virtual void | setIndex (Scalar index) |
virtual void | setTimeStep (Scalar dt) |
virtual void | setErrorAbs (Scalar errorAbs) |
virtual void | setErrorRel (Scalar errorRel) |
virtual void | setOrder (int order) |
virtual void | setNFailures (int nFailures) |
virtual void | setNRunningFailures (int nFailures) |
virtual void | setNConsecutiveFailures (int nConsecutiveFailures) |
virtual void | setTolRel (Scalar tolRel) |
virtual void | setTolAbs (Scalar tolAbs) |
virtual void | setXNormL2 (Scalar xNormL2) |
virtual void | setDxNormL2Rel (Scalar dxNormL2Rel) |
virtual void | setDxNormL2Abs (Scalar dxNormL2Abs) |
virtual void | setComputeNorms (bool computeNorms) |
virtual void | setSolutionStatus (Status s) |
virtual void | setSolutionStatus (const Thyra::SolveStatus< Scalar > sStatus) |
virtual void | setOutput (bool output) |
virtual void | setOutputScreen (bool output) |
virtual void | setIsSynced (bool isSynced) |
virtual void | setIsInterpolated (bool isInterpolated) |
virtual void | setAccuracy (bool accuracy) |
Get State Data | |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getX () |
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > | getX () const |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getXDot () |
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > | getXDot () const |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getXDotDot () |
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > | getXDotDot () const |
virtual Teuchos::RCP < StepperState< Scalar > > | getStepperState () |
virtual Teuchos::RCP< const StepperState< Scalar > > | getStepperState () const |
virtual Teuchos::RCP < PhysicsState< Scalar > > | getPhysicsState () |
virtual Teuchos::RCP< const PhysicsState< Scalar > > | getPhysicsState () const |
Set State Data | |
virtual void | setX (Teuchos::RCP< Thyra::VectorBase< Scalar > > x) |
virtual void | setX (Teuchos::RCP< const Thyra::VectorBase< Scalar > > x) |
virtual void | setXDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xdot) |
virtual void | setXDot (Teuchos::RCP< const Thyra::VectorBase< Scalar > > xdot) |
virtual void | setXDotDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xdotdot) |
virtual void | setXDotDot (Teuchos::RCP< const Thyra::VectorBase< Scalar > > xdotdot) |
virtual void | setStepperState (Teuchos::RCP< StepperState< Scalar > > &ss) |
virtual void | setStepperState (const Teuchos::RCP< StepperState< Scalar > > &ss) |
virtual void | setPhysicsState (Teuchos::RCP< PhysicsState< Scalar > > &ps) |
virtual void | setPhysicsState (const Teuchos::RCP< PhysicsState< Scalar > > &ps) |
Comparison methods | |
bool | operator< (const SolutionState< Scalar > &ss) const |
Less than comparison for sorting based on time. More... | |
bool | operator<= (const SolutionState< Scalar > &ss) const |
Less than or equal to comparison for sorting based on time. More... | |
bool | operator< (const Scalar &t) const |
Less than comparison for sorting based on time. More... | |
bool | operator<= (const Scalar &t) const |
Less than or equal to comparison for sorting based on time. More... | |
bool | operator> (const SolutionState< Scalar > &ss) const |
Greater than comparison for sorting based on time. More... | |
bool | operator>= (const SolutionState< Scalar > &ss) const |
Greater than or equal to comparison for sorting based on time. More... | |
bool | operator> (const Scalar &t) const |
Greater than comparison for sorting based on time. More... | |
bool | operator>= (const Scalar &t) const |
Greater than or equal to comparison for sorting based on time. More... | |
bool | operator== (const SolutionState< Scalar > &ss) const |
Equality comparison for matching. More... | |
bool | operator== (const Scalar &t) const |
Equality comparison for matching. More... | |
Overridden from Teuchos::Describable | |
virtual std::string | description () const |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
Additional Inherited Members | |
![]() | |
static void | setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel) |
static EVerbosityLevel | getDefaultVerbLevel () |
![]() | |
static void | setDefaultOStream (const RCP< FancyOStream > &defaultOStream) |
static RCP< FancyOStream > | getDefaultOStream () |
![]() | |
static const EVerbosityLevel | verbLevel_default |
![]() | |
void | initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) |
![]() | |
void | initializeVerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) |
virtual void | informUpdatedVerbosityState () const |
Solution state for integrators and steppers.
The SolutionState object contains all the information and functions related to solution variables at a particular state or time, . The member data is simply the solution variables,
; optionally the solution time derivatives,
and
; and all associated metadata (more on this shortly). The primary driver for this definition is the requirement that the SolutionState has all the information necessary to restart the time integration, e.g., from a failed time step, a checkpointed solution, and/or for adjoint sensitivity analysis. Additionally, by having all the solution and metadata incapsulated, it allows easy passing via a single object, and many SolutionStates forming a history (see SolutionHistory) can be easily used for multi-step methods, solution interpolation, adjoint sensitivities, etc.
The solution variables, ,
, and which allows us to use all of its functionalities, e.g., norms. At a minimum, the solution,
, needs to have a valid RCP of a Thyra::Vector, and a Teuchos::null indicates the SolutionState is not fully constructed. The solution time derivatives,
and
, are optional. All Tempus steppers (Tempus::Stepper) will efficiently work with just
available (i.e., steppers will maintain an
and/or
as needed for efficiency). The primary decision on whether to have the time derivatives as part of the SolutionState is related to memory usage versus recomputing
and
. If the user needs readily available access to the time derivatives (e.g., in their ModelEvaluator, diagnostics or output), then they should likely be added to the SolutionState. Otherwise, the user should just store the solution,
, and recompute the time derivatives as needed.
There are two important concepts for the SolutionState related to the solution, and its time derivatives,
and
. The first concept is a solution and its time derivatives are said to be consistent with respect to its governing equation, if it satisfies its explicit ODE,
or
, or its implicit ODE
. Thus, there is a relationship between the solution and its time derivatives, such that they satisfy its governing equation. Obviously, there are times when the solution is not consistent with its time derivatives, e.g., initial guess for the nonlinear solve and interpolated solutions. Additionally, because of the discrete representation, this relationship is also dependent on the particular stepper, e.g., for Forward Euler, it is defined via
, so would not be consistent for another stepper such as the Trapezoidal Rule, where the time derivative is defined as
.
Consistency is a very important attribute of a solution to be able to use it for initial conditions and restart. At the cost of a solve, a SolutionState can be made consistent for the Stepper being used by resetting the time derivatives, and
, to match the solution,
(see Tempus::Stepper::setInitialConditions).
The second concept for the SolutionState is the idea of the solution, and its time derivatives,
and
of being synchronized. During the time-integration process (loop), the solution and its time derivatives may not be at the same time level, e.g.,
. This is common for Leapfrog time integration, where the solution,
; the first time derivative,
; and the second time derivative,
are all at different time levels, thus they are termed unsynchronized. Of course, there are methods to synchronize them to the same time level,
. Synchronization is important to ensure accuracy for output, visualization, and verification. All Tempus::Steppers can ensure the SolutionState is synchronized at the end of the time step, if needed. Otherwise, the unsynchronized SolutionState can be maintained for efficiency reasons. (see Tempus::SolutionState::getIsSynced()).
For more complex time integration where the physics has additional state information or the time integrator is not a one-step method (i.e., cannot accurately start from a single time step), this class can be inherited and the physics state (PhysicsState) or additional time-integration parameters can be managed.
SolutionStates can be interpolated to generate solutions at various times (see SolutionHistory). However not all metadata or state information can be interpolated. Thus interpolated solutions may not be suitable for checkpointing, restart and undo operations, but may be useful for diagnostics, output and/or adjoint sensitivities.
The solution vectors, ,
, and
, in SolutionState can be null pointers. This indicates that the application does not need them, so do not storage them. This can be a huge savings when saving many states in the solution history. Some Steppers will need temporary memory to store time derivative(s) (
, or
) for evaluation of the ODE/DAE (
), but each individual Stepper manages that.
Definition at line 137 of file Tempus_SolutionState_decl.hpp.
Tempus::SolutionState< Scalar >::SolutionState | ( | ) |
Default Constructor – Not meant for immediate adding to SolutionHistory. This constructor does not set the solution vectors, x, xdot and xdotdot. which should be set via setX(), setXDot(), and/or setXDotDot() prior to being added to SolutionHistory.
Definition at line 18 of file Tempus_SolutionState_impl.hpp.
Tempus::SolutionState< Scalar >::SolutionState | ( | const Teuchos::RCP< SolutionStateMetaData< Scalar > > | ssmd, |
const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | x, | ||
const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | xdot, | ||
const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | xdotdot, | ||
const Teuchos::RCP< StepperState< Scalar > > & | stepperState, | ||
const Teuchos::RCP< PhysicsState< Scalar > > & | physicsState | ||
) |
Definition at line 40 of file Tempus_SolutionState_impl.hpp.
Tempus::SolutionState< Scalar >::SolutionState | ( | const Teuchos::RCP< const SolutionStateMetaData< Scalar > > | ssmd, |
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > & | x, | ||
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > & | xdot, | ||
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > & | xdotdot, | ||
const Teuchos::RCP< const StepperState< Scalar > > & | stepperState, | ||
const Teuchos::RCP< const PhysicsState< Scalar > > & | physicsState | ||
) |
Definition at line 71 of file Tempus_SolutionState_impl.hpp.
Tempus::SolutionState< Scalar >::SolutionState | ( | const SolutionState< Scalar > & | ss | ) |
This is a shallow copy constructor, use clone for a deep copy constructor.
Definition at line 101 of file Tempus_SolutionState_impl.hpp.
|
inlinevirtual |
Destructor.
Definition at line 180 of file Tempus_SolutionState_decl.hpp.
|
virtual |
This is a deep copy constructor.
Definition at line 118 of file Tempus_SolutionState_impl.hpp.
|
virtual |
This is a deep copy.
Definition at line 149 of file Tempus_SolutionState_impl.hpp.
|
virtual |
Deep copy solution data, but keep metaData untouched.
Definition at line 158 of file Tempus_SolutionState_impl.hpp.
|
inlinevirtual |
Definition at line 185 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 186 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 190 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 191 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 192 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 193 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 194 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 195 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 196 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 197 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 198 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 199 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 200 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 201 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 202 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 203 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 204 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 205 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 206 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 207 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 208 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 209 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 210 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 211 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 212 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 217 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 220 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 223 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 225 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 227 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 229 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 231 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 233 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 236 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 238 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 240 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 243 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 245 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 248 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 250 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 252 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 254 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 257 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 259 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 267 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 269 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 271 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 273 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 275 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 281 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 284 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 286 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 288 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 290 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 292 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 295 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 298 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 301 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 303 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 309 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 311 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 313 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 315 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 317 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 319 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 322 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 324 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 327 of file Tempus_SolutionState_decl.hpp.
|
inlinevirtual |
Definition at line 329 of file Tempus_SolutionState_decl.hpp.
bool Tempus::SolutionState< Scalar >::operator< | ( | const SolutionState< Scalar > & | ss | ) | const |
Less than comparison for sorting based on time.
Definition at line 213 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator<= | ( | const SolutionState< Scalar > & | ss | ) | const |
Less than or equal to comparison for sorting based on time.
Definition at line 219 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator< | ( | const Scalar & | t | ) | const |
Less than comparison for sorting based on time.
Definition at line 225 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator<= | ( | const Scalar & | t | ) | const |
Less than or equal to comparison for sorting based on time.
Definition at line 231 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator> | ( | const SolutionState< Scalar > & | ss | ) | const |
Greater than comparison for sorting based on time.
Definition at line 237 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator>= | ( | const SolutionState< Scalar > & | ss | ) | const |
Greater than or equal to comparison for sorting based on time.
Definition at line 243 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator> | ( | const Scalar & | t | ) | const |
Greater than comparison for sorting based on time.
Definition at line 249 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator>= | ( | const Scalar & | t | ) | const |
Greater than or equal to comparison for sorting based on time.
Definition at line 255 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator== | ( | const SolutionState< Scalar > & | ss | ) | const |
Equality comparison for matching.
Definition at line 261 of file Tempus_SolutionState_impl.hpp.
bool Tempus::SolutionState< Scalar >::operator== | ( | const Scalar & | t | ) | const |
Equality comparison for matching.
Definition at line 267 of file Tempus_SolutionState_impl.hpp.
|
virtual |
Reimplemented from Teuchos::Describable.
Definition at line 273 of file Tempus_SolutionState_impl.hpp.
|
virtual |
Reimplemented from Teuchos::Describable.
Definition at line 285 of file Tempus_SolutionState_impl.hpp.
|
virtual |
Compute the solution norms, and solution change from ssIn, if provided.
Definition at line 321 of file Tempus_SolutionState_impl.hpp.
|
private |
Meta Data for the solution state.
Definition at line 382 of file Tempus_SolutionState_decl.hpp.
|
private |
Definition at line 383 of file Tempus_SolutionState_decl.hpp.
|
private |
Solution.
Definition at line 386 of file Tempus_SolutionState_decl.hpp.
|
private |
Definition at line 387 of file Tempus_SolutionState_decl.hpp.
|
private |
Time derivative of the solution.
Definition at line 390 of file Tempus_SolutionState_decl.hpp.
|
private |
Definition at line 391 of file Tempus_SolutionState_decl.hpp.
|
private |
Second time derivative of the solution.
Definition at line 394 of file Tempus_SolutionState_decl.hpp.
|
private |
Definition at line 395 of file Tempus_SolutionState_decl.hpp.
|
private |
StepperState for this SolutionState.
Definition at line 398 of file Tempus_SolutionState_decl.hpp.
|
private |
Definition at line 399 of file Tempus_SolutionState_decl.hpp.
|
private |
PhysicsState for this SolutionState.
Definition at line 402 of file Tempus_SolutionState_decl.hpp.
|
private |
Definition at line 403 of file Tempus_SolutionState_decl.hpp.