| Tempus
    Version of the Day
    Time Integration | 
Stepper interface to support full-space optimization. More...
#include <Tempus_StepperOptimizationInterface.hpp>
 
  
 | Public Member Functions | |
| StepperOptimizationInterface () | |
| virtual | ~StepperOptimizationInterface () | 
| virtual int | stencilLength () const =0 | 
| Return the number of solution vectors in the time step stencil.  More... | |
| virtual void | computeStepResidual (Thyra::VectorBase< Scalar > &residual, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index) const =0 | 
| Compute time step residual.  More... | |
| virtual void | computeStepJacobian (Thyra::LinearOpBase< Scalar > &jacobian, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index, const int deriv_index) const =0 | 
| Compute time step Jacobian.  More... | |
| virtual void | computeStepParamDeriv (Thyra::LinearOpBase< Scalar > &deriv, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index) const =0 | 
| Compute time step derivative w.r.t. model parameters.  More... | |
| virtual void | computeStepSolver (Thyra::LinearOpWithSolveBase< Scalar > &jacobian_solver, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index) const =0 | 
| Compute time step Jacobian solver.  More... | |
Stepper interface to support full-space optimization.
This is a potential interface to support transient full-space optimizations methods such as those implemented by ROL through its DynamicConstraint interface. This interface is subject to major revision!
Design consideration: Take array of solution vectors as input, or solution history, or array of states?
The length of x is determined by the time step stencil, e.g., in an m-step BDF method x would contain x_n,...,x_{n-m} at time step n in x[0],...,x[m] with time values stored in t similarly. p is the vector of design parameters and param_index determines which model parameter this corresponds to.
Definition at line 41 of file Tempus_StepperOptimizationInterface.hpp.
| 
 | inline | 
Definition at line 43 of file Tempus_StepperOptimizationInterface.hpp.
| 
 | inlinevirtual | 
Definition at line 45 of file Tempus_StepperOptimizationInterface.hpp.
| 
 | pure virtual | 
Return the number of solution vectors in the time step stencil.
Implemented in Tempus::StepperBackwardEuler< Scalar >.
| 
 | pure virtual | 
Compute time step residual.
Implemented in Tempus::StepperBackwardEuler< Scalar >.
| 
 | pure virtual | 
Compute time step Jacobian.
deriv_index determines which component of x the derivative should be computed with respect to.
Implemented in Tempus::StepperBackwardEuler< Scalar >.
| 
 | pure virtual | 
Compute time step derivative w.r.t. model parameters.
Implemented in Tempus::StepperBackwardEuler< Scalar >.
| 
 | pure virtual | 
Compute time step Jacobian solver.
Derivative is always w.r.t. the most current solution vector
Implemented in Tempus::StepperBackwardEuler< Scalar >.
 1.8.5
 1.8.5