Tempus
Version of the Day
Time Integration
|
ModelEvaluator for forming adjoint sensitivity equations. More...
#include <Tempus_AdjointSensitivityModelEvaluator_decl.hpp>
Public Types | |
typedef Thyra::VectorBase< Scalar > | Vector |
typedef Thyra::MultiVectorBase < Scalar > | MultiVector |
Public Types inherited from Thyra::ModelEvaluator< class > | |
enum | EInArgsMembers |
enum | EInArgs_p_mp |
enum | EEvalType |
enum | EDerivativeMultiVectorOrientation |
enum | EDerivativeLinearOp |
enum | EDerivativeLinearity |
enum | ERankStatus |
enum | EOutArgsMembers |
enum | EOutArgsDfDp |
enum | EOutArgsDgDx_dot |
enum | EOutArgsDgDx |
enum | EOutArgsDgDp |
enum | EOutArgsDfDp_mp |
enum | EOutArgs_g_mp |
enum | EOutArgsDgDx_dot_mp |
enum | EOutArgsDgDx_mp |
enum | EOutArgsDgDp_mp |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | ScalarMag |
Public Member Functions | |
AdjointSensitivityModelEvaluator (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &adjoint_residual_model, const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &adjoint_solve_model, const Scalar &t_init, const Scalar &t_final, const bool is_pseudotransient, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null) | |
Constructor. More... | |
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () const |
Get the underlying model 'f'. More... | |
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > | getAdjointResidualModel () const |
Get the underlying adjoint residual model. More... | |
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > | getAdjointSolveModel () const |
Get the underlying adjoint solve model. More... | |
void | setFinalTime (const Scalar t_final) |
Set the final time from the forward evaluation. More... | |
void | setForwardSolutionHistory (const Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > &sh) |
Set solution history from forward evaluation. More... | |
Public Member Functions inherited from Thyra::StateFuncModelEvaluatorBase< Scalar > | |
RCP< const VectorSpaceBase < Scalar > > | get_p_space (int l) const |
RCP< const Teuchos::Array < std::string > > | get_p_names (int l) const |
RCP< const VectorSpaceBase < Scalar > > | get_g_space (int j) const |
Teuchos::ArrayView< const std::string > | get_g_names (int j) const |
ModelEvaluatorBase::InArgs < Scalar > | getNominalValues () const |
ModelEvaluatorBase::InArgs < Scalar > | getLowerBounds () const |
ModelEvaluatorBase::InArgs < Scalar > | getUpperBounds () const |
RCP< LinearOpBase< Scalar > > | create_W_op () const |
RCP< PreconditionerBase< Scalar > > | create_W_prec () const |
RCP< const LinearOpWithSolveFactoryBase < Scalar > > | get_W_factory () const |
void | reportFinalPoint (const ModelEvaluatorBase::InArgs< Scalar > &finalPoint, const bool wasSolved) |
Public Member Functions inherited from Thyra::ModelEvaluator< class > | |
virtual int | Np () const =0 |
virtual int | Ng () const =0 |
virtual RCP< const VectorSpaceBase< Scalar > > | get_f_multiplier_space () const =0 |
virtual RCP< const VectorSpaceBase< Scalar > > | get_g_multiplier_space (int j) const =0 |
virtual RCP < LinearOpWithSolveBase < Scalar > > | create_W () const =0 |
virtual RCP< LinearOpBase < Scalar > > | create_DfDp_op (int l) const =0 |
virtual RCP< LinearOpBase < Scalar > > | create_DgDx_dot_op (int j) const =0 |
virtual RCP< LinearOpBase < Scalar > > | create_DgDx_op (int j) const =0 |
virtual RCP< LinearOpBase < Scalar > > | create_DgDp_op (int j, int l) const =0 |
virtual ModelEvaluatorBase::OutArgs < Scalar > | createOutArgs () const =0 |
virtual void | evalModel (const ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const =0 |
ModelEvaluatorBase () | |
std::string | toString (ModelEvaluatorBase::EInArgsMembers) |
std::string | toString (ModelEvaluatorBase::EOutArgsMembers) |
std::string | toString (ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
ModelEvaluatorBase::EDerivativeMultiVectorOrientation | getOtherDerivativeMultiVectorOrientation (ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
Public Member Functions inherited from Teuchos::Describable | |
virtual std::string | description () const |
virtual void | describe (FancyOStream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
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) |
Static Public Member Functions | |
static Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () |
Private Types | |
typedef Thyra::DefaultMultiVectorProductVectorSpace < Scalar > | DMVPVS |
typedef Thyra::DefaultMultiVectorProductVector < Scalar > | DMVPV |
Private Member Functions | |
Thyra::ModelEvaluatorBase::OutArgs < Scalar > | createOutArgsImpl () const |
void | evalModelImpl (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const |
Public functions overridden from ModelEvaulator. | |
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > | get_p_space (int p) const |
Teuchos::RCP< const Teuchos::Array< std::string > > | get_p_names (int p) const |
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > | get_x_space () const |
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > | get_f_space () const |
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > | get_g_space (int j) const |
Teuchos::RCP < Thyra::LinearOpBase< Scalar > > | create_W_op () const |
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase < Scalar > > | get_W_factory () const |
Thyra::ModelEvaluatorBase::InArgs < Scalar > | createInArgs () const |
Thyra::ModelEvaluatorBase::InArgs < Scalar > | getNominalValues () const |
Additional Inherited Members | |
Static Public Attributes inherited from Thyra::ModelEvaluator< class > | |
static const int | NUM_E_IN_ARGS_MEMBERS |
static const int | NUM_E_OUT_ARGS_MEMBERS |
Static Public Attributes inherited from Teuchos::Describable | |
static const EVerbosityLevel | verbLevel_default |
ModelEvaluator for forming adjoint sensitivity equations.
This class wraps a given ModelEvalutor encapsulating f(x_dot,x,p) and creates a new "residual" for the adjoint sensitivity equations: F(y) = d/dt( df/dx_dot^T*y ) - df/dx^T*y + dg/dx^T = 0 where y is the adjoint variable. Nominally y is a multi-vector belonging to the vector space of the residual f with number of columns equal to the dimension of g. To satisfy the model evaluator interface, y is converted to a product vector formed by its columns. This is the form of the adjoint equations suitable for computing pseudotransientsensitivities of a response function g(x(T),p) or adjoint sensitivities for a time-integrated response function int_0^T g(x(t),p).
To compute adjoint sensitivities, the equations f(x_dot,x,p) = 0 are integrated forward in time to some final time T, with the adjoint equations above integrated backward in time starting at T. Since the tempus integrator can only integrate forward in time, we define tau = T - t and transform the adjoint equations to: F(y) = d/dtau( df/dx_dot^T*y ) + df/dx^T*y - dg/dx^T = 0. The initial conditions for y are y(T) = 0. The sensitivity of g(T) is then int_0^T(dg/dp(T) - df/dp^T*y)dt + dx/dp(0)^T*df/dx_dot(0)^T*y(0) for transient adjoint sensitivites and is dg/dp(T) - df/dp^T*y for pseudotransient.
This model evaluator supports both implicit and explict forms of the governing equations. However it assumes df/dx_dot is constant with respect to x_dot, x, and t so that the adjoint equations become F(y) = df/dxdot^T*y_dot + df/dx^T*y - dg/dx^T = 0. Moving beyond this assumption will require modifications to the steppers in how they generate time-derivative terms.
Definition at line 55 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
typedef Thyra::VectorBase<Scalar> Tempus::AdjointSensitivityModelEvaluator< Scalar >::Vector |
Definition at line 58 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
typedef Thyra::MultiVectorBase<Scalar> Tempus::AdjointSensitivityModelEvaluator< Scalar >::MultiVector |
Definition at line 59 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 144 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 145 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
Tempus::AdjointSensitivityModelEvaluator< Scalar >::AdjointSensitivityModelEvaluator | ( | const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > & | model, |
const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > & | adjoint_residual_model, | ||
const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > & | adjoint_solve_model, | ||
const Scalar & | t_init, | ||
const Scalar & | t_final, | ||
const bool | is_pseudotransient, | ||
const Teuchos::RCP< const Teuchos::ParameterList > & | pList = Teuchos::null |
||
) |
Constructor.
t_final is the final integration time used to implement the time transformation described above. num_adjoint is the number of adjoint variables defining the number of columns in y, which is determined by the number of elements in the vector space for the response g. The optionally supplied parameter list supports the following options:
Definition at line 29 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
inline |
Get the underlying model 'f'.
Definition at line 91 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
inline |
Get the underlying adjoint residual model.
Definition at line 97 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
inline |
Get the underlying adjoint solve model.
Definition at line 104 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
void Tempus::AdjointSensitivityModelEvaluator< Scalar >::setFinalTime | ( | const Scalar | t_final | ) |
Set the final time from the forward evaluation.
Definition at line 123 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
void Tempus::AdjointSensitivityModelEvaluator< Scalar >::setForwardSolutionHistory | ( | const Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > & | sh | ) |
Set solution history from forward evaluation.
Definition at line 130 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 149 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 157 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 165 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 172 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 179 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 188 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 199 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 215 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< class >.
Definition at line 222 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
static |
Definition at line 541 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
private |
Definition at line 253 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
private |
Definition at line 259 of file Tempus_AdjointSensitivityModelEvaluator_impl.hpp.
|
private |
Definition at line 153 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 154 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 156 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 157 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 158 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 160 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 161 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 162 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 163 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 164 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 165 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 166 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 167 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 168 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 169 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 170 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
private |
Definition at line 171 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 173 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 174 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 175 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 176 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 177 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 178 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 179 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 180 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 181 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 182 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.
|
mutableprivate |
Definition at line 183 of file Tempus_AdjointSensitivityModelEvaluator_decl.hpp.