Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_StepperOptimizationInterface.hpp
Go to the documentation of this file.
1//@HEADER
2// *****************************************************************************
3// Tempus: Time Integration and Sensitivity Analysis Package
4//
5// Copyright 2017 NTESS and the Tempus contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8//@HEADER
9
10#ifndef Tempus_Stepper_Optimization_Interface_hpp
11#define Tempus_Stepper_Optimization_Interface_hpp
12
13// Teuchos
14#include "Teuchos_Array.hpp"
15#include "Teuchos_RCP.hpp"
16
17// Thyra
18#include "Thyra_VectorBase.hpp"
19#include "Thyra_LinearOpBase.hpp"
20#include "Thyra_LinearOpWithSolveBase.hpp"
21#include "Tempus_config.hpp"
22
23namespace Tempus {
24
40template <class Scalar>
42 public:
44
46
48 virtual int stencilLength() const = 0;
49
51 virtual void computeStepResidual(
53 const Teuchos::Array<Teuchos::RCP<const Thyra::VectorBase<Scalar> > >& x,
54 const Teuchos::Array<Scalar>& t, const Thyra::VectorBase<Scalar>& p,
55 const int param_index) const = 0;
56
58
62 virtual void computeStepJacobian(
64 const Teuchos::Array<Teuchos::RCP<const Thyra::VectorBase<Scalar> > >& x,
65 const Teuchos::Array<Scalar>& t, const Thyra::VectorBase<Scalar>& p,
66 const int param_index, const int deriv_index) const = 0;
67
71 const Teuchos::Array<Teuchos::RCP<const Thyra::VectorBase<Scalar> > >& x,
72 const Teuchos::Array<Scalar>& t, const Thyra::VectorBase<Scalar>& p,
73 const int param_index) const = 0;
74
76
79 virtual void computeStepSolver(
81 const Teuchos::Array<Teuchos::RCP<const Thyra::VectorBase<Scalar> > >& x,
82 const Teuchos::Array<Scalar>& t, const Thyra::VectorBase<Scalar>& p,
83 const int param_index) const = 0;
84};
85
86} // namespace Tempus
87
88#endif // Tempus_Stepper_hpp
Stepper interface to support full-space optimization.
virtual int stencilLength() const =0
Return the number of solution vectors in the time step stencil.
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.
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.
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.
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.