Rythmos - Transient Integration for Differential Equations  Version of the Day
 All Classes Functions Variables Typedefs Pages
Related Functions | List of all members
Rythmos::AdjointModelEvaluator< Scalar > Class Template Reference

Standard concrete adjoint ModelEvaluator for time-constant mass matrix models. More...

#include <Rythmos_AdjointModelEvaluator.hpp>

Inherits StateFuncModelEvaluatorBase< Scalar >.

Related Functions

(Note that these are not member functions.)

template<class Scalar >
RCP< AdjointModelEvaluator
< Scalar > > 
adjointModelEvaluator (const RCP< const Thyra::ModelEvaluator< Scalar > > &fwdStateModel, const TimeRange< Scalar > &fwdTimeRange)
 Nonmember constructor. More...
 

Constructors/Intializers/Accessors

 AdjointModelEvaluator ()
 
void setFwdStateModel (const RCP< const Thyra::ModelEvaluator< Scalar > > &fwdStateModel, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &basePoint)
 Set the underlying forward model and base point. More...
 
void setFwdTimeRange (const TimeRange< Scalar > &fwdTimeRange)
 Set the forward time range that this adjoint model will be defined over. More...
 
void setFwdStateSolutionBuffer (const RCP< const InterpolationBufferBase< Scalar > > &fwdStateSolutionBuffer)
 Set the interpolation buffer that will return values of the state solution x and x_dot at various points t as the adjoint is integrated backwards in time. More...
 

Public functions overridden from ModelEvaulator.

RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_x_space () const
 
RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_f_space () const
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
getNominalValues () const
 
RCP
< Thyra::LinearOpWithSolveBase
< Scalar > > 
create_W () const
 
RCP< Thyra::LinearOpBase
< Scalar > > 
create_W_op () const
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
createInArgs () const
 

Detailed Description

template<class Scalar>
class Rythmos::AdjointModelEvaluator< Scalar >

Standard concrete adjoint ModelEvaluator for time-constant mass matrix models.

Overview

This concrete ModelEvalautor subclass takes any suitable ModelEvalautor object and creates the adjoint model for use by any appropriate time integration method..

When the mass matrix d(f)/d(x_dot) is not a function of t (which this class assumes), then the adjoint model can be represented as:

    d(f)/d(x_dot)^T * lambda_dot - d(f)/d(x)^T * lambda + d(g)/d(x)^T

This model is stated in reverse time t_bar <: [0,t_final-t_initial] (with d/d(t) = -d/d(t_bar)) which results in the new adjoint formuation

  f_bar(x_bar_dot, x_bar, t_bar)
    = d(f)/d(x_dot)^T * lambda_rev_dot + d(f)/d(x)^T * lambda - d(g)/d(x)^T

Where:

WARNING! When interacting with this interface you must take note that reverse time is being used as defined above! This is especially important if you are going to use lambda_dot for anything. You have been warned!

Implementation Notes

Here, we describe how the residual of the adjoint system f_bar(...) is actually computed from the capabilities of the underlying forward model.

First, note that

  W_bar = alpha_bar * d(f_bar)/d(x_bar_dot) + beta_bar * d(f_bar)/d(x_bar)

        = alpha_bar * d(f)/d(x_dot)^T + beta_bar * d(f)/d(x)^T

This means that W_bar can be computed directly as W_bar_adj on the underlying forward state ModelEvaluator object as:

  W_bar_adj = alpha_bar * d(f)/d(x_dot) + beta_bar * d(f)/d(x)

by passing in alpha = alpha_bar and beta = beta_bar. We then use the subclass Thyra::DefaultAdjointLinearOpWithSolve to create W_bar = adjoint(W_bar_adj) and that is it.

Now, given that the client will request the form of W_bar = adjoint(W_bar_adj) described above, we want to use this W_bar_adj object in computing the adjoint equation residual f_bar. To see how to do this, note that from the above definition of W_bar that we have:

  d(f)/d(x)^T = 1/beta_bar * W_bar_adj^T
    - alpha_bar/beta_bar * d(f)/d(x_dot)^T

By using the above equation for d(f)/d(x)^T, we can eliminate d(f)/d(x) from f_bar to yield:

  f_bar = d(f)/d(x_dot)^T * lambda_hat + 1/beta_bar * W_bar_adj^T * lambda
          - d(g)/d(x)^T

     where:

        lambda_hat = lambda_rev_dot - alpha_bar/beta_bar * lambda

Above, we see we need to compute d(f)/d(x_dot) sperately from W_bar_adj from the underlying forward ModelEvaluator. Note that for many forward models, that d(f)/d(x_dot) will actually be constant and can be computed up front and reused throughout.

Todo:
Add support to response function derivative source d(g)/d(x)^T.

Todo:
Add support for more than one adjoint through the DefaultMultiVectorProductVector[Space] sublasses.

Todo:
Add functionality to the Thyra::ModelEvaluator::OutArgs class to communicate the dependence of a function on its input arguments. We need to know the exact dependance of f(...) on x_dot, x, and t to know if this class can be used and what shortcuts can be used with it.

Definition at line 173 of file Rythmos_AdjointModelEvaluator.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Rythmos::AdjointModelEvaluator< Scalar >::AdjointModelEvaluator ( )

Definition at line 297 of file Rythmos_AdjointModelEvaluator.hpp.

Member Function Documentation

template<class Scalar >
void Rythmos::AdjointModelEvaluator< Scalar >::setFwdStateModel ( const RCP< const Thyra::ModelEvaluator< Scalar > > &  fwdStateModel,
const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  basePoint 
)

Set the underlying forward model and base point.

Definition at line 303 of file Rythmos_AdjointModelEvaluator.hpp.

template<class Scalar >
void Rythmos::AdjointModelEvaluator< Scalar >::setFwdTimeRange ( const TimeRange< Scalar > &  fwdTimeRange)

Set the forward time range that this adjoint model will be defined over.

Definition at line 316 of file Rythmos_AdjointModelEvaluator.hpp.

template<class Scalar >
void Rythmos::AdjointModelEvaluator< Scalar >::setFwdStateSolutionBuffer ( const RCP< const InterpolationBufferBase< Scalar > > &  fwdStateSolutionBuffer)

Set the interpolation buffer that will return values of the state solution x and x_dot at various points t as the adjoint is integrated backwards in time.

NOTE: If the model is linear in x and x_dot, then the client can avoid setting this interpolation buffer since it will never be called.

NOTE: This object need be fully initialized at this point. It only needs to be fully initialized before this->evalModel(...) is called. This just sets up the basic link to this object.

Definition at line 324 of file Rythmos_AdjointModelEvaluator.hpp.

template<class Scalar >
RCP< const Thyra::VectorSpaceBase< Scalar > > Rythmos::AdjointModelEvaluator< Scalar >::get_x_space ( ) const

Definition at line 337 of file Rythmos_AdjointModelEvaluator.hpp.

template<class Scalar >
RCP< const Thyra::VectorSpaceBase< Scalar > > Rythmos::AdjointModelEvaluator< Scalar >::get_f_space ( ) const

Definition at line 346 of file Rythmos_AdjointModelEvaluator.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > Rythmos::AdjointModelEvaluator< Scalar >::getNominalValues ( ) const

Definition at line 355 of file Rythmos_AdjointModelEvaluator.hpp.

template<class Scalar >
RCP< Thyra::LinearOpWithSolveBase< Scalar > > Rythmos::AdjointModelEvaluator< Scalar >::create_W ( ) const

Definition at line 364 of file Rythmos_AdjointModelEvaluator.hpp.

template<class Scalar >
RCP< Thyra::LinearOpBase< Scalar > > Rythmos::AdjointModelEvaluator< Scalar >::create_W_op ( ) const

Definition at line 373 of file Rythmos_AdjointModelEvaluator.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > Rythmos::AdjointModelEvaluator< Scalar >::createInArgs ( ) const

Definition at line 382 of file Rythmos_AdjointModelEvaluator.hpp.

Friends And Related Function Documentation

template<class Scalar >
RCP< AdjointModelEvaluator< Scalar > > adjointModelEvaluator ( const RCP< const Thyra::ModelEvaluator< Scalar > > &  fwdStateModel,
const TimeRange< Scalar > &  fwdTimeRange 
)
related

Nonmember constructor.

Definition at line 276 of file Rythmos_AdjointModelEvaluator.hpp.


The documentation for this class was generated from the following file: