Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus Namespace Reference

Classes

class  AdjointAuxSensitivityModelEvaluator
 ModelEvaluator for forming adjoint sensitivity equations. More...
 
class  AdjointSensitivityModelEvaluator
 ModelEvaluator for forming adjoint sensitivity equations. More...
 
class  AuxiliaryIntegralModelEvaluator
 ModelEvaluator for integrating auxiliary equations. More...
 
class  CombinedForwardSensitivityModelEvaluator
 Transform a ModelEvaluator's sensitivity equations to its residual. More...
 
class  Stepper
 Thyra Base interface for time steppers. More...
 
class  SolutionHistory
 SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of SolutionStates for later retrival and reuse, such as checkpointing, restart, and undo operations. More...
 
class  TimeStepControl
 TimeStepControl manages the time step size. There several mechanicisms that effect the time step size and handled with this class: More...
 
class  Integrator
 Thyra Base interface for time integrators. Time integrators are designed to advance the solution from an initial time, $t_0$, to a final time, $t_f$. More...
 
class  IntegratorAdjointSensitivity
 Time integrator suitable for adjoint sensitivity analysis. More...
 
class  IntegratorBasic
 Basic time integrator. More...
 
class  IntegratorForwardSensitivity
 Time integrator implementing forward sensitivity analysis. More...
 
class  IntegratorObserver
 IntegratorObserver class for time integrators. More...
 
class  IntegratorObserverBasic
 IntegratorObserverBasic class for time integrators. This basic class has simple no-op functions, as all basic functionality should be handled through other methods. More...
 
class  IntegratorObserverComposite
 This observer is a composite observer,. More...
 
class  IntegratorObserverLogging
 This observer logs calls to observer functions. This observer simply logs and counts the calls to each of the observer functions. This is useful in monirtoring and debugging the time integration. More...
 
class  IntegratorObserverNoOp
 IntegratorObserverNoOp class for time integrators. This basic class has simple no-op functions, as all basic functionality should be handled through other methods. More...
 
class  IntegratorObserverSubcycling
 IntegratorObserverSubcycling class for time integrators. This basic class has simple no-op functions, as all basic functionality should be handled through other methods. More...
 
class  IntegratorPseudoTransientAdjointSensitivity
 Time integrator suitable for pseudotransient adjoint sensitivity analysis. More...
 
class  IntegratorPseudoTransientForwardSensitivity
 Time integrator suitable for pseudotransient forward sensitivity analysis. More...
 
class  Interpolator
 Base strategy class for interpolation functionality. More...
 
class  InterpolatorFactory
 Interpolator factory. More...
 
class  InterpolatorLagrange
 Concrete implemenation of Interpolator that does simple lagrange interpolation. More...
 
class  PhysicsState
 PhysicsState is a simple class to hold information about the physics. More...
 
class  RKButcherTableau
 Runge-Kutta methods. More...
 
class  SensitivityModelEvaluatorBase
 A ModelEvaluator decorator for sensitivity analysis. More...
 
class  SolutionState
 Solution state for integrators and steppers. SolutionState contains the metadata for solutions and the solutions themselves. More...
 
class  SolutionStateMetaData
 Solution state meta data. More...
 
class  StaggeredForwardSensitivityModelEvaluator
 Transform a ModelEvaluator's sensitivity equations to its residual. More...
 
class  StepperBackwardEuler
 Backward Euler time stepper. More...
 
class  StepperBackwardEulerTimeDerivative
 Time-derivative interface for Backward Euler. More...
 
class  StepperFactory
 Stepper factory. More...
 
class  StepperBackwardEulerAppAction
 Application Action for StepperBackwardEuler. More...
 
class  StepperBackwardEulerAppActionComposite
 This composite AppAction loops over added AppActions. More...
 
class  StepperBackwardEulerModifierBase
 Base modifier for StepperBackwardEuler. More...
 
class  StepperBackwardEulerModifierDefault
 Default modifier for StepperBackwardEuler. More...
 
class  StepperBackwardEulerModifierXBase
 Base ModifierX for StepperBackwardEuler. More...
 
class  StepperBackwardEulerModifierXDefault
 Default ModifierX for StepperBackwardEuler. More...
 
class  StepperBackwardEulerObserver
 StepperBackwardEulerObserver class for StepperBackwardEuler. More...
 
class  StepperBackwardEulerObserverBase
 Base observer for StepperBackwardEuler. More...
 
class  StepperBackwardEulerObserverDefault
 Default observer for StepperBackwardEuler. More...
 
class  StepperBDF2
 BDF2 (Backward-Difference-Formula-2) time stepper. More...
 
class  StepperBDF2TimeDerivative
 Time-derivative interface for BDF2. More...
 
class  StepperBDF2Observer
 StepperBDF2Observer class for StepperBDF2. More...
 
class  StepperDIRK
 Diagonally Implicit Runge-Kutta (DIRK) time stepper. More...
 
class  StepperDIRKTimeDerivative
 Time-derivative interface for DIRK. More...
 
class  StepperExplicit
 Thyra Base interface for implicit time steppers. More...
 
class  StepperExplicitRK
 Explicit Runge-Kutta time stepper. More...
 
class  StepperForwardEuler
 Forward Euler time stepper. More...
 
class  StepperForwardEulerAppAction
 Application Action for StepperForwardEuler. More...
 
class  StepperForwardEulerAppActionComposite
 This composite AppAction loops over added AppActions. More...
 
class  StepperForwardEulerModifierBase
 Base modifier for StepperBackwardEuler. More...
 
class  StepperForwardEulerModifierDefault
 Default modifier for StepperForwardEuler. More...
 
class  StepperForwardEulerModifierXBase
 Base ModifierX for StepperForwardEuler. More...
 
class  StepperForwardEulerModifierXDefault
 Default ModifierX for StepperForwardEuler. More...
 
class  StepperForwardEulerObserver
 StepperForwardEulerObserver class for StepperForwardEuler. More...
 
class  StepperForwardEulerObserverBase
 Base observer for StepperForwardEuler. More...
 
class  StepperForwardEulerObserverDefault
 Default observer for StepperForwardEuler. More...
 
class  StepperHHTAlpha
 HHT-Alpha time stepper. More...
 
class  StepperIMEX_RK
 Implicit-Explicit Runge-Kutta (IMEX-RK) time stepper. More...
 
class  StepperIMEX_RKTimeDerivative
 Time-derivative interface for IMEX RK. More...
 
class  StepperIMEX_RK_Partition
 Partitioned Implicit-Explicit Runge-Kutta (IMEX-RK) time stepper. More...
 
class  StepperIMEX_RKPartTimeDerivative
 Time-derivative interface for Partitioned IMEX RK. More...
 
class  ImplicitODEParameters
 
class  StepperImplicit
 Thyra Base interface for implicit time steppers. More...
 
class  StepperLeapfrog
 Leapfrog time stepper. More...
 
class  StepperLeapfrogObserver
 StepperLeapfrogObserver class for StepperLeapfrog. More...
 
class  StepperNewmarkExplicitAForm
 Newmark Explicit time stepper. More...
 
class  StepperNewmarkImplicitAForm
 Newmark time stepper in acceleration form (a-form). More...
 
class  StepperNewmarkImplicitDForm
 Newmark time stepper. More...
 
class  StepperObserver
 StepperObserver class for Stepper class. More...
 
class  StepperObserverBasic
 StepperObserverBasic class for Stepper class. More...
 
class  StepperObserverComposite
 This observer is a composite observer,. More...
 
class  StepperOperatorSplit
 OperatorSplit stepper loops through the Stepper list. More...
 
class  StepperOperatorSplitAppAction
 StepperOperatorSplitAppAction class for StepperOperatorSplit. More...
 
class  StepperOperatorSplitAppActionComposite
 This composite AppAction loops over added AppActions. More...
 
class  StepperOperatorSplitModifierBase
 Base modifier for OperatorSplit. More...
 
class  StepperOperatorSplitModifierDefault
 Default modifier for StepperOperatorSplit. More...
 
class  StepperOperatorSplitModifierXBase
 Base ModifierX for StepperOperatorSplit. More...
 
class  StepperOperatorSplitModifierXDefault
 Default ModifierX for StepperOperatorSplit. More...
 
class  StepperOperatorSplitObserver
 StepperOperatorSplitObserver class for StepperOperatorSplit. More...
 
class  StepperOperatorSplitObserverBase
 Base observer for StepperOperatorSplit. More...
 
class  StepperOperatorSplitObserverDefault
 Default observer for StepperOperatorSplit. More...
 
class  StepperOptimizationInterface
 Stepper interface to support full-space optimization. More...
 
class  StepperRKBase
 Base class for Runge-Kutta methods, ExplicitRK, DIRK and IMEX. More...
 
class  StepperRKAppAction
 Application Action for StepperRKBase. More...
 
class  StepperRKAppActionComposite
 This composite AppAction loops over added AppActions. More...
 
class  StepperERK_ForwardEuler
 Forward Euler Runge-Kutta Butcher Tableau. More...
 
class  StepperERK_4Stage4thOrder
 Runge-Kutta 4th order Butcher Tableau. More...
 
class  StepperERK_BogackiShampine32
 Explicit RK Bogacki-Shampine Butcher Tableau. More...
 
class  StepperERK_Merson45
 Explicit RK Merson Butcher Tableau. More...
 
class  StepperERK_3_8Rule
 Explicit RK 3/8th Rule Butcher Tableau. More...
 
class  StepperERK_4Stage3rdOrderRunge
 RK Explicit 4 Stage 3rd order by Runge. More...
 
class  StepperERK_5Stage3rdOrderKandG
 RK Explicit 5 Stage 3rd order by Kinnmark and Gray. More...
 
class  StepperERK_3Stage3rdOrder
 RK Explicit 3 Stage 3rd order. More...
 
class  StepperERK_3Stage3rdOrderTVD
 RK Explicit 3 Stage 3rd order TVD. More...
 
class  StepperERK_3Stage3rdOrderHeun
 RK Explicit 3 Stage 3rd order by Heun. More...
 
class  StepperERK_Midpoint
 RK Explicit Midpoint. More...
 
class  StepperERK_Ralston
 RK Explicit Ralston. More...
 
class  StepperERK_Trapezoidal
 RK Explicit Trapezoidal. More...
 
class  StepperERK_SSPERK54
 Strong Stability Preserving Explicit RK Butcher Tableau. More...
 
class  StepperERK_General
 General Explicit Runge-Kutta Butcher Tableau. More...
 
class  StepperDIRK_BackwardEuler
 Backward Euler Runge-Kutta Butcher Tableau. More...
 
class  StepperSDIRK_2Stage2ndOrder
 SDIRK 2 Stage 2nd order. More...
 
class  StepperSDIRK_3Stage2ndOrder
 SDIRK 3 Stage 2nd order. More...
 
class  StepperSDIRK_2Stage3rdOrder
 SDIRK 2 Stage 3rd order. More...
 
class  StepperEDIRK_2Stage3rdOrder
 EDIRK 2 Stage 3rd order. More...
 
class  StepperDIRK_1StageTheta
 SDIRK 1 Stage Theta. More...
 
class  StepperEDIRK_2StageTheta
 EDIRK 2 Stage Theta Method. More...
 
class  StepperEDIRK_TrapezoidalRule
 RK Trapezoidal Rule (A.K.A. RK Crank-Nicolson) More...
 
class  StepperSDIRK_ImplicitMidpoint
 SDIRK Implicit Midpoint. More...
 
class  StepperSDIRK_SSPDIRK22
 Strong Stability Preserving Diagonally-Implicit RK Butcher Tableau. More...
 
class  StepperSDIRK_SSPDIRK32
 Strong Stability Preserving Diagonally-Implicit RK Butcher Tableau. More...
 
class  StepperSDIRK_SSPDIRK23
 Strong Stability Preserving Diagonally-Implicit RK Butcher Tableau. More...
 
class  StepperSDIRK_SSPDIRK33
 Strong Stability Preserving Diagonally-Implicit RK Butcher Tableau. More...
 
class  StepperDIRK_1Stage1stOrderRadauIA
 RK Implicit 1 Stage 1st order Radau IA. More...
 
class  StepperDIRK_2Stage2ndOrderLobattoIIIB
 RK Implicit 2 Stage 2nd order Lobatto IIIB. More...
 
class  StepperSDIRK_5Stage4thOrder
 SDIRK 5 Stage 4th order. More...
 
class  StepperSDIRK_3Stage4thOrder
 SDIRK 3 Stage 4th order. More...
 
class  StepperSDIRK_5Stage5thOrder
 SDIRK 5 Stage 5th order. More...
 
class  StepperSDIRK_21Pair
 SDIRK 2(1) pair. More...
 
class  StepperDIRK_General
 General Implicit Runge-Kutta Butcher Tableau. More...
 
class  StepperRKModifierBase
 Base modifier for StepperRK. More...
 
class  StepperRKModifierDefault
 Default modifier for StepperRK. More...
 
class  StepperRKModifierXBase
 Base ModifierX for StepperRK. More...
 
class  StepperRKModifierXDefault
 Default ModifierX for StepperRK. More...
 
class  StepperRKObserver
 StepperRKObserver class for StepperRK. More...
 
class  StepperRKObserverBase
 Base observer for StepperRK. More...
 
class  StepperRKObserverComposite
 This observer is a composite observer,. More...
 
class  StepperRKObserverDefault
 Default observer for StepperRK. More...
 
class  StepperRKObserverLogging
 This observer logs calls to observer functions. This observer simply logs and counts the calls to each of the observer functions. This is useful in monirtoring and debugging the time integration. More...
 
class  StepperStaggeredForwardSensitivity
 A stepper implementing staggered forward sensitivity analysis. More...
 
class  StepperState
 StepperState is a simple class to hold state information about the stepper. More...
 
class  StepperSubcycling
 Subcycling time stepper. More...
 
class  StepperSubcyclingAppAction
 Application Action for StepperSubcycling. More...
 
class  StepperSubcyclingAppActionComposite
 This composite AppAction loops over added AppActions. More...
 
class  StepperSubcyclingModifierBase
 Base modifier for StepperSubcycling. More...
 
class  StepperSubcyclingModifierDefault
 Default modifier for StepperSubcycling. More...
 
class  StepperSubcyclingModifierXBase
 Base ModifierX for StepperSubcycling. More...
 
class  StepperSubcyclingModifierXDefault
 Default ModifierX for StepperSubcycling. More...
 
class  StepperSubcyclingObserver
 StepperSubcyclingObserver class for StepperSubcycling. More...
 
class  StepperSubcyclingObserverBase
 Base observer for StepperSubcycling. More...
 
class  StepperSubcyclingObserverDefault
 Default observer for StepperSubcycling. More...
 
class  StepperTrapezoidal
 Trapezoidal method time stepper. More...
 
class  StepperTrapezoidalTimeDerivative
 Time-derivative interface for Trapezoidal method. More...
 
class  StepperTrapezoidalObserver
 StepperTrapezoidalObserver class for StepperTrapezoidal. More...
 
class  TimeDerivative
 This interface defines the time derivative connection between an implicit Stepper and WrapperModelEvaluator. More...
 
class  TimeEventBase
 This class defines time events which can be used to "trigger" an action. Time events are points in time and/or timestep index where an an action should occur, such as, solution output (mesh and solution) diagnostic output, restart, screen dump, in-situ visualization, user-specified, or any other action. More...
 
class  TimeEventComposite
 This composite TimeEvent loops over added TimeEvents. More...
 
class  TimeEventList
 TimeEventList specifies a list of time events. More...
 
class  TimeEventListIndex
 TimeEventListIndex specifies a list of index events. More...
 
class  TimeEventRange
 TimeEventRange specifies a start, stop and stride time. More...
 
class  TimeEventRangeIndex
 TimeEventRangeIndex specifies a start, stop and stride index. More...
 
class  TimeStepControlStrategy
 StepControlStrategy class for TimeStepControl. More...
 
class  TimeStepControlStrategyBasicVS
 StepControlStrategy class for TimeStepControl. More...
 
class  TimeStepControlStrategyComposite
 StepControlStrategy class for TimeStepControl. More...
 
class  TimeStepControlStrategyConstant
 StepControlStrategy class for TimeStepControl. More...
 
class  TimeStepControlStrategyIntegralController
 StepControlStrategy class for TimeStepControl. More...
 
class  TimeStepControlStrategyPID
 StepControlStrategy class for TimeStepControl. More...
 
class  WrapperModelEvaluator
 A ModelEvaluator which wraps the application ModelEvaluator. More...
 
class  WrapperModelEvaluatorBasic
 A ModelEvaluator for residual evaluations given a state. This ModelEvaluator takes a state, x, and determines its residual, $ g(x) $, which is suitable for a nonlinear solve. This is accomplished by computing the time derivative of the state, x_dot, (through Lambda functions), supplying the current time, and calling the application application ModelEvaluator, $ f(x,\dot{x},t) $. More...
 
class  WrapperModelEvaluatorPairIMEX
 ModelEvaluator pair for implicit and explicit (IMEX) evaluations. More...
 
class  WrapperModelEvaluatorPairIMEX_Basic
 ModelEvaluator pair for implicit and explicit (IMEX) evaulations. More...
 
class  WrapperModelEvaluatorPairIMEX_CombinedFSA
 Specialization of IMEX ME for "combined" FSA method. More...
 
class  WrapperModelEvaluatorPairIMEX_StaggeredFSA
 Specialization of IMEX ME for "staggered" FSA method. More...
 
class  WrapperModelEvaluatorPairPartIMEX_Basic
 ModelEvaluator pair for implicit and explicit (IMEX) evaulations. More...
 
class  WrapperModelEvaluatorPairPartIMEX_CombinedFSA
 Specialization of IMEX-Part ME for "combined" FSA method. More...
 
class  WrapperModelEvaluatorPairPartIMEX_StaggeredFSA
 Specialization of IMEX-Part ME for "combined" FSA method. More...
 
class  WrapperModelEvaluatorSecondOrder
 A ModelEvaluator for residual evaluations given a state. This ModelEvaluator takes a state, x, and determines its residual, $ g(x) $, which is suitable for a nonlinear solve. This is accomplished by computing the time derivative of the state, x_dot, (through Lambda functions), supplying the current time, and calling the application application ModelEvaluator, $ f(\dot{x},x,t) $. More...
 

Enumerations

enum  StorageType {
  STORAGE_TYPE_INVALID = 0, STORAGE_TYPE_KEEP_NEWEST = 1, STORAGE_TYPE_UNDO = 2, STORAGE_TYPE_STATIC = 3,
  STORAGE_TYPE_UNLIMITED = 4
}
 
enum  OrderODE { FIRST_ORDER_ODE = 1, SECOND_ORDER_ODE = 2 }
 
enum  Status { PASSED, FAILED, WORKING }
 Status for the Integrator, the Stepper and the SolutionState. More...
 
enum  EVALUATION_TYPE { EVALUATE_RESIDUAL, SOLVE_FOR_X, SOLVE_FOR_XDOT_CONST_X }
 EVALUATION_TYPE indicates the evaluation to apply to the implicit ODE. More...
 

Functions

template<class Scalar >
Teuchos::RCP
< IntegratorAdjointSensitivity
< Scalar > > 
integratorAdjointSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< IntegratorAdjointSensitivity
< Scalar > > 
integratorAdjointSensitivity ()
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP< IntegratorBasic
< Scalar > > 
integratorBasic (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP< IntegratorBasic
< Scalar > > 
integratorBasic (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP< IntegratorBasic
< Scalar > > 
integratorBasic ()
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP< IntegratorBasic
< Scalar > > 
integratorBasic (Teuchos::RCP< Teuchos::ParameterList > pList, std::vector< Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > > models)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< IntegratorForwardSensitivity
< Scalar > > 
integratorForwardSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< IntegratorForwardSensitivity
< Scalar > > 
integratorForwardSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< IntegratorForwardSensitivity
< Scalar > > 
integratorForwardSensitivity ()
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< Tempus::IntegratorPseudoTransientAdjointSensitivity
< Scalar > > 
integratorPseudoTransientAdjointSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< Tempus::IntegratorPseudoTransientAdjointSensitivity
< Scalar > > 
integratorPseudoTransientAdjointSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< Tempus::IntegratorPseudoTransientAdjointSensitivity
< Scalar > > 
integratorPseudoTransientAdjointSensitivity ()
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< IntegratorPseudoTransientAdjointSensitivity
< Scalar > > 
integratorPseudoTransientAdjointSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< IntegratorPseudoTransientAdjointSensitivity
< Scalar > > 
integratorPseudoTransientAdjointSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< IntegratorPseudoTransientAdjointSensitivity
< Scalar > > 
integratorPseudoTransientAdjointSensitivity ()
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< Tempus::IntegratorPseudoTransientForwardSensitivity
< Scalar > > 
integratorPseudoTransientForwardSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< Tempus::IntegratorPseudoTransientForwardSensitivity
< Scalar > > 
integratorPseudoTransientForwardSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType)
 Non-member constructor. More...
 
template<class Scalar >
Teuchos::RCP
< Tempus::IntegratorPseudoTransientForwardSensitivity
< Scalar > > 
integratorPseudoTransientForwardSensitivity ()
 Non-member constructor. More...
 
template<class Scalar >
void interpolate (const Interpolator< Scalar > &interpolator, const Scalar &t, SolutionState< Scalar > *state_out)
 Non-member constructor. More...
 
template<class Scalar >
void interpolate (Interpolator< Scalar > &interpolator, const Teuchos::RCP< const std::vector< Teuchos::RCP< SolutionState< Scalar > > > > &nodes, const Scalar &t, SolutionState< Scalar > *state_out)
 Non-member constructor. More...
 
template<typename Scalar >
bool floating_compare_equals (const Scalar &a, const Scalar &b, const Scalar &scale)
 Helper function for comparing times. More...
 
template<class Scalar >
Teuchos::RCP
< InterpolatorLagrange< Scalar > > 
lagrangeInterpolator ()
 
template<class Scalar >
Teuchos::RCP< SolutionHistory
< Scalar > > 
createSolutionHistoryPL (Teuchos::RCP< Teuchos::ParameterList > pList)
 Nonmember constructor from a ParameterList. More...
 
template<class Scalar >
Teuchos::RCP< SolutionHistory
< Scalar > > 
createSolutionHistoryState (const Teuchos::RCP< SolutionState< Scalar > > &state)
 Nonmember contructor from a SolutionState. More...
 
template<class Scalar >
Teuchos::RCP< SolutionHistory
< Scalar > > 
createSolutionHistoryME (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
 Nonmember contructor from a Thyra ModelEvaluator. More...
 
template<class Scalar >
Teuchos::RCP< SolutionState
< Scalar > > 
createSolutionStateX (const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdot=Teuchos::null, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdotdot=Teuchos::null)
 Nonmember constructor from non-const solution vectors, x. More...
 
template<class Scalar >
Teuchos::RCP< SolutionState
< Scalar > > 
createSolutionStateX (const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &xdot=Teuchos::null, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &xdotdot=Teuchos::null)
 Nonmember constructor from const solution vectors, x. More...
 
template<class Scalar >
Teuchos::RCP< SolutionState
< Scalar > > 
createSolutionStateME (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< StepperState< Scalar > > &stepperState=Teuchos::null, const Teuchos::RCP< PhysicsState< Scalar > > &physicsState=Teuchos::null)
 Nonmember constructor from Thyra ModelEvaluator. More...
 
void trim (std::string &str)
 Removes whitespace at beginning and end of string. More...
 
void StringTokenizer (std::vector< std::string > &tokens, const std::string &str, const std::string delimiter=",", bool trim=false)
 Tokenize a string, put tokens in a vector. More...
 
void TokensToDoubles (std::vector< double > &values, const std::vector< std::string > &tokens)
 Turn a vector of tokens into a vector of doubles. More...
 
void TokensToInts (std::vector< int > &values, const std::vector< std::string > &tokens)
 Turn a vector of tokens into a vector of ints. More...
 
template<typename ScalarT >
ScalarT getScalarParameter (const std::string &field, const Teuchos::ParameterList &plist)
 
const std::string toString (const Status status)
 Convert Status to string. More...
 
std::string version ()
 
template<typename Scalar >
Teuchos::RCP
< SensitivityModelEvaluatorBase
< Scalar > > 
wrapCombinedFSAModelEvaluator (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
 
template<typename Scalar >
Teuchos::RCP
< SensitivityModelEvaluatorBase
< Scalar > > 
wrapCombinedFSAModelEvaluator (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
 
template<typename Scalar >
Teuchos::RCP
< SensitivityModelEvaluatorBase
< Scalar > > 
wrapStaggeredFSAModelEvaluator (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
 
template<typename Scalar >
Teuchos::RCP
< SensitivityModelEvaluatorBase
< Scalar > > 
wrapStaggeredFSAModelEvaluator (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
 
Helper functions
void getValidParametersBasic (Teuchos::RCP< Teuchos::ParameterList > pl, std::string stepperType)
 Provide basic parameters to Steppers. More...
 
template<class Scalar >
void validExplicitODE (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
 Validate that the model supports explicit ODE evaluation, f(x,t) [=xdot]. More...
 
template<class Scalar >
void validSecondOrderExplicitODE (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
 Validate that the model supports explicit second order ODE evaluation, f(x,xdot,t) [=xdotdot]. More...
 
template<class Scalar >
void validImplicitODE_DAE (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
 Validate ME supports implicit ODE/DAE evaluation, f(xdot,x,t) [= 0]. More...
 
template<class Scalar >
void validSecondOrderODE_DAE (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
 Validate ME supports 2nd order implicit ODE/DAE evaluation, f(xdotdot,xdot,x,t) [= 0]. More...
 
Teuchos::RCP
< Teuchos::ParameterList > 
defaultSolverParameters ()
 Returns the default solver ParameterList for implicit Steppers. More...
 

Enumeration Type Documentation

EVALUATION_TYPE indicates the evaluation to apply to the implicit ODE.

Enumerator
EVALUATE_RESIDUAL 

Evaluate residual for the implicit ODE.

SOLVE_FOR_X 

Solve for x and determine xDot from x.

SOLVE_FOR_XDOT_CONST_X 

Solve for xDot keeping x constant (for ICs).

Definition at line 18 of file Tempus_WrapperModelEvaluator.hpp.

Enumerator
FIRST_ORDER_ODE 

Stepper integrates first-order ODEs.

SECOND_ORDER_ODE 

Stepper integrates second-order ODEs.

Definition at line 27 of file Tempus_Stepper_decl.hpp.

Status for the Integrator, the Stepper and the SolutionState.

Enumerator
PASSED 
FAILED 
WORKING 

Definition at line 16 of file Tempus_Types.hpp.

Enumerator
STORAGE_TYPE_INVALID 

Invalid storage type.

STORAGE_TYPE_KEEP_NEWEST 

Keep the single newest state.

STORAGE_TYPE_UNDO 

Keep the 2 newest states for undo.

STORAGE_TYPE_STATIC 

Keep a fix number of states.

STORAGE_TYPE_UNLIMITED 

Grow the history as needed.

Definition at line 26 of file Tempus_SolutionHistory_decl.hpp.

Function Documentation

template<class Scalar >
Teuchos::RCP< SolutionHistory< Scalar > > Tempus::createSolutionHistoryME ( const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &  model)

Nonmember contructor from a Thyra ModelEvaluator.

Definition at line 615 of file Tempus_SolutionHistory_impl.hpp.

template<class Scalar >
Teuchos::RCP< SolutionHistory< Scalar > > Tempus::createSolutionHistoryPL ( Teuchos::RCP< Teuchos::ParameterList >  pList)

Nonmember constructor from a ParameterList.

Definition at line 594 of file Tempus_SolutionHistory_impl.hpp.

template<class Scalar >
Teuchos::RCP< SolutionHistory< Scalar > > Tempus::createSolutionHistoryState ( const Teuchos::RCP< SolutionState< Scalar > > &  state)

Nonmember contructor from a SolutionState.

Definition at line 604 of file Tempus_SolutionHistory_impl.hpp.

template<class Scalar >
Teuchos::RCP< SolutionState< Scalar > > Tempus::createSolutionStateME ( const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &  model,
const Teuchos::RCP< StepperState< Scalar > > &  stepperState = Teuchos::null,
const Teuchos::RCP< PhysicsState< Scalar > > &  physicsState = Teuchos::null 
)

Nonmember constructor from Thyra ModelEvaluator.

Definition at line 515 of file Tempus_SolutionState_impl.hpp.

template<class Scalar >
Teuchos::RCP< SolutionState< Scalar > > Tempus::createSolutionStateX ( const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  x,
const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  xdot = Teuchos::null,
const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  xdotdot = Teuchos::null 
)

Nonmember constructor from non-const solution vectors, x.

Definition at line 472 of file Tempus_SolutionState_impl.hpp.

template<class Scalar >
Teuchos::RCP< SolutionState< Scalar > > Tempus::createSolutionStateX ( const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  x,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  xdot = Teuchos::null,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  xdotdot = Teuchos::null 
)

Nonmember constructor from const solution vectors, x.

Definition at line 493 of file Tempus_SolutionState_impl.hpp.

Teuchos::RCP< Teuchos::ParameterList > Tempus::defaultSolverParameters ( )

Returns the default solver ParameterList for implicit Steppers.

Definition at line 78 of file Tempus_Stepper_impl.hpp.

template<typename Scalar >
bool Tempus::floating_compare_equals ( const Scalar &  a,
const Scalar &  b,
const Scalar &  scale 
)

Helper function for comparing times.

Definition at line 90 of file Tempus_Interpolator.hpp.

template<typename ScalarT >
ScalarT Tempus::getScalarParameter ( const std::string &  field,
const Teuchos::ParameterList &  plist 
)

Read in a parameter field and return the correct scalar field. This parses scalar type data

Definition at line 39 of file Tempus_String_Utilities.hpp.

void Tempus::getValidParametersBasic ( Teuchos::RCP< Teuchos::ParameterList >  pl,
std::string  stepperType 
)

Provide basic parameters to Steppers.

Definition at line 16 of file Tempus_Stepper_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorAdjointSensitivity< Scalar > > Tempus::integratorAdjointSensitivity ( Teuchos::RCP< Teuchos::ParameterList >  pList,
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model 
)

Non-member constructor.

Definition at line 552 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorAdjointSensitivity< Scalar > > Tempus::integratorAdjointSensitivity ( )

Non-member constructor.

Definition at line 564 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorBasic< Scalar > > Tempus::integratorBasic ( Teuchos::RCP< Teuchos::ParameterList >  pList,
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model 
)

Non-member constructor.

Definition at line 707 of file Tempus_IntegratorBasic_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorBasic< Scalar > > Tempus::integratorBasic ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model,
std::string  stepperType 
)

Non-member constructor.

Definition at line 718 of file Tempus_IntegratorBasic_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorBasic< Scalar > > Tempus::integratorBasic ( )

Non-member constructor.

Definition at line 729 of file Tempus_IntegratorBasic_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorBasic< Scalar > > Tempus::integratorBasic ( Teuchos::RCP< Teuchos::ParameterList >  pList,
std::vector< Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > >  models 
)

Non-member constructor.

Definition at line 738 of file Tempus_IntegratorBasic_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorForwardSensitivity< Scalar > > Tempus::integratorForwardSensitivity ( Teuchos::RCP< Teuchos::ParameterList >  pList,
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model 
)

Non-member constructor.

Definition at line 330 of file Tempus_IntegratorForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorForwardSensitivity< Scalar > > Tempus::integratorForwardSensitivity ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model,
std::string  stepperType 
)

Non-member constructor.

Definition at line 342 of file Tempus_IntegratorForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorForwardSensitivity< Scalar > > Tempus::integratorForwardSensitivity ( )

Non-member constructor.

Definition at line 354 of file Tempus_IntegratorForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> > Tempus::integratorPseudoTransientAdjointSensitivity ( Teuchos::RCP< Teuchos::ParameterList >  pList,
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model 
)

Non-member constructor.

Definition at line 461 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> > Tempus::integratorPseudoTransientAdjointSensitivity ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model,
std::string  stepperType 
)

Non-member constructor.

Definition at line 473 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> > Tempus::integratorPseudoTransientAdjointSensitivity ( )

Non-member constructor.

Definition at line 485 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> > Tempus::integratorPseudoTransientAdjointSensitivity ( Teuchos::RCP< Teuchos::ParameterList >  pList,
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model 
)

Non-member constructor.

Definition at line 461 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> > Tempus::integratorPseudoTransientAdjointSensitivity ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model,
std::string  stepperType 
)

Non-member constructor.

Definition at line 473 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> > Tempus::integratorPseudoTransientAdjointSensitivity ( )

Non-member constructor.

Definition at line 485 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar > > Tempus::integratorPseudoTransientForwardSensitivity ( Teuchos::RCP< Teuchos::ParameterList >  pList,
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model 
)

Non-member constructor.

Definition at line 550 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar > > Tempus::integratorPseudoTransientForwardSensitivity ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model,
std::string  stepperType 
)

Non-member constructor.

Definition at line 562 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar > > Tempus::integratorPseudoTransientForwardSensitivity ( )

Non-member constructor.

Definition at line 574 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
void Tempus::interpolate ( const Interpolator< Scalar > &  interpolator,
const Scalar &  t,
SolutionState< Scalar > *  state_out 
)

Non-member constructor.

Definition at line 70 of file Tempus_Interpolator.hpp.

template<class Scalar >
void Tempus::interpolate ( Interpolator< Scalar > &  interpolator,
const Teuchos::RCP< const std::vector< Teuchos::RCP< SolutionState< Scalar > > > > &  nodes,
const Scalar &  t,
SolutionState< Scalar > *  state_out 
)

Non-member constructor.

Definition at line 79 of file Tempus_Interpolator.hpp.

template<class Scalar >
Teuchos::RCP<InterpolatorLagrange<Scalar> > Tempus::lagrangeInterpolator ( )

Definition at line 74 of file Tempus_InterpolatorLagrange_decl.hpp.

void Tempus::StringTokenizer ( std::vector< std::string > &  tokens,
const std::string &  str,
const std::string  delimiters,
bool  trim 
)

Tokenize a string, put tokens in a vector.

Definition at line 31 of file Tempus_String_Utilities.cpp.

void Tempus::TokensToDoubles ( std::vector< double > &  values,
const std::vector< std::string > &  tokens 
)

Turn a vector of tokens into a vector of doubles.

Definition at line 63 of file Tempus_String_Utilities.cpp.

void Tempus::TokensToInts ( std::vector< int > &  values,
const std::vector< std::string > &  tokens 
)

Turn a vector of tokens into a vector of ints.

Definition at line 77 of file Tempus_String_Utilities.cpp.

const std::string Tempus::toString ( const Status  status)
inline

Convert Status to string.

Definition at line 25 of file Tempus_Types.hpp.

void Tempus::trim ( std::string &  str)

Removes whitespace at beginning and end of string.

Definition at line 14 of file Tempus_String_Utilities.cpp.

template<class Scalar >
void Tempus::validExplicitODE ( const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &  model)

Validate that the model supports explicit ODE evaluation, f(x,t) [=xdot].

Currently the convention to evaluate f(x,t) is to set xdot=null! There is no InArgs support to test if xdot is null, so we set xdot=null and hopefully the ModelEvaluator can handle it.

Definition at line 228 of file Tempus_Stepper_decl.hpp.

template<class Scalar >
void Tempus::validImplicitODE_DAE ( const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &  model)

Validate ME supports implicit ODE/DAE evaluation, f(xdot,x,t) [= 0].

Definition at line 292 of file Tempus_Stepper_decl.hpp.

template<class Scalar >
void Tempus::validSecondOrderExplicitODE ( const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &  model)

Validate that the model supports explicit second order ODE evaluation, f(x,xdot,t) [=xdotdot].

Currently the convention to evaluate f(x,xdot,t) is to set xdotdot=null! There is no InArgs support to test if xdotdot is null, so we set xdotdot=null and hopefully the ModelEvaluator can handle it.

Definition at line 260 of file Tempus_Stepper_decl.hpp.

template<class Scalar >
void Tempus::validSecondOrderODE_DAE ( const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &  model)

Validate ME supports 2nd order implicit ODE/DAE evaluation, f(xdotdot,xdot,x,t) [= 0].

Definition at line 338 of file Tempus_Stepper_decl.hpp.

std::string Tempus::version ( )

Definition at line 13 of file Tempus_Version.cpp.

template<typename Scalar >
Teuchos::RCP< SensitivityModelEvaluatorBase<Scalar> > Tempus::wrapCombinedFSAModelEvaluator ( const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &  model,
const Teuchos::RCP< const Teuchos::ParameterList > &  pList = Teuchos::null 
)

Helper function for creating a CombinedForwardSensitivityModelEvaluator from a given application model evaluator. It handles the complexity introducted by IMEX steppers where the sensitivity model evaluator needs to be put inside the IMEX pair model evaluators.

Definition at line 24 of file Tempus_WrapCombinedFSAModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP< SensitivityModelEvaluatorBase<Scalar> > Tempus::wrapCombinedFSAModelEvaluator ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model,
const Teuchos::RCP< const Teuchos::ParameterList > &  pList = Teuchos::null 
)

Definition at line 58 of file Tempus_WrapCombinedFSAModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP< SensitivityModelEvaluatorBase<Scalar> > Tempus::wrapStaggeredFSAModelEvaluator ( const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &  model,
const Teuchos::RCP< const Teuchos::ParameterList > &  pList = Teuchos::null 
)

Helper function for creating a StaggeredForwardSensitivityModelEvaluator from a given application model evaluator. It handles the complexity introducted by IMEX steppers where the sensitivity model evaluator needs to be put inside the IMEX pair model evaluators.

Definition at line 27 of file Tempus_WrapStaggeredFSAModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP< SensitivityModelEvaluatorBase<Scalar> > Tempus::wrapStaggeredFSAModelEvaluator ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &  model,
const Teuchos::RCP< const Teuchos::ParameterList > &  pList = Teuchos::null 
)

Definition at line 62 of file Tempus_WrapStaggeredFSAModelEvaluator.hpp.