Tempus
Version of the Day
Time Integration
Toggle main menu visibility
Loading...
Searching...
No Matches
src
Tempus_StepperOptimizationInterface.hpp
Go to the documentation of this file.
1
// @HEADER
2
// ****************************************************************************
3
// Tempus: Copyright (2017) Sandia Corporation
4
//
5
// Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6
// ****************************************************************************
7
// @HEADER
8
9
#ifndef Tempus_Stepper_Optimization_Interface_hpp
10
#define Tempus_Stepper_Optimization_Interface_hpp
11
12
//Teuchos
13
#include "Teuchos_Array.hpp"
14
#include "Teuchos_RCP.hpp"
15
16
// Thyra
17
#include "Thyra_VectorBase.hpp"
18
#include "Thyra_LinearOpBase.hpp"
19
#include "Thyra_LinearOpWithSolveBase.hpp"
20
#include "Tempus_config.hpp"
21
22
namespace
Tempus
{
23
39
template
<
class
Scalar>
40
class
StepperOptimizationInterface
41
{
42
public
:
43
44
StepperOptimizationInterface
() {}
45
46
virtual
~StepperOptimizationInterface
() {}
47
49
virtual
int
stencilLength
()
const
= 0;
50
52
virtual
void
computeStepResidual
(
53
Thyra::VectorBase<Scalar>
& residual,
54
const
Teuchos::Array< Teuchos::RCP<
const
Thyra::VectorBase<Scalar>
> >& x,
55
const
Teuchos::Array<Scalar>& t,
56
const
Thyra::VectorBase<Scalar>
& p,
57
const
int
param_index)
const
= 0;
58
60
64
virtual
void
computeStepJacobian
(
65
Thyra::LinearOpBase<Scalar>
& jacobian,
66
const
Teuchos::Array< Teuchos::RCP<
const
Thyra::VectorBase<Scalar>
> >& x,
67
const
Teuchos::Array<Scalar>& t,
68
const
Thyra::VectorBase<Scalar>
& p,
69
const
int
param_index,
70
const
int
deriv_index)
const
= 0;
71
73
virtual
void
computeStepParamDeriv
(
74
Thyra::LinearOpBase<Scalar>
& deriv,
75
const
Teuchos::Array< Teuchos::RCP<
const
Thyra::VectorBase<Scalar>
> >& x,
76
const
Teuchos::Array<Scalar>& t,
77
const
Thyra::VectorBase<Scalar>
& p,
78
const
int
param_index)
const
= 0;
79
81
84
virtual
void
computeStepSolver
(
85
Thyra::LinearOpWithSolveBase<Scalar>
& jacobian_solver,
86
const
Teuchos::Array< Teuchos::RCP<
const
Thyra::VectorBase<Scalar>
> >& x,
87
const
Teuchos::Array<Scalar>& t,
88
const
Thyra::VectorBase<Scalar>
& p,
89
const
int
param_index)
const
= 0;
90
};
91
92
}
// namespace Tempus
93
94
#endif
// Tempus_Stepper_hpp
Tempus::StepperOptimizationInterface::stencilLength
virtual int stencilLength() const =0
Return the number of solution vectors in the time step stencil.
Tempus::StepperOptimizationInterface::computeStepResidual
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.
Tempus::StepperOptimizationInterface::computeStepParamDeriv
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.
Tempus::StepperOptimizationInterface::StepperOptimizationInterface
StepperOptimizationInterface()
Definition
Tempus_StepperOptimizationInterface.hpp:44
Tempus::StepperOptimizationInterface::~StepperOptimizationInterface
virtual ~StepperOptimizationInterface()
Definition
Tempus_StepperOptimizationInterface.hpp:46
Tempus::StepperOptimizationInterface::computeStepJacobian
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.
Tempus::StepperOptimizationInterface::computeStepSolver
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.
Thyra::LinearOpBase
Thyra::LinearOpWithSolveBase
Thyra::VectorBase
Tempus
Definition
Tempus_AdjointAuxSensitivityModelEvaluator_decl.hpp:21
Generated by
1.17.0