ROL
Public Member Functions | List of all members
ROL::LinearOperator< Real > Class Template Referenceabstract

Provides the interface to apply a linear operator. More...

#include <ROL_LinearOperator.hpp>

+ Inheritance diagram for ROL::LinearOperator< Real >:

Public Member Functions

virtual ~LinearOperator ()
 
virtual void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update linear operator. More...
 
virtual void apply (Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const =0
 Apply linear operator. More...
 
virtual void applyInverse (Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
 Apply inverse of linear operator. More...
 
virtual void applyAdjoint (Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
 Apply adjoint of linear operator. More...
 
virtual void applyAdjointInverse (Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
 Apply adjoint of the inverse linear operator. More...
 

Detailed Description

template<class Real>
class ROL::LinearOperator< Real >

Provides the interface to apply a linear operator.

ROL's linear operator interface is designed to interface with ROL's Krylov methods. These linear operators often represent projected Hessians or preconditioners. The basic operator interace, to be implemented by the user, requires:

The user may also implement:

Definition at line 71 of file ROL_LinearOperator.hpp.

Constructor & Destructor Documentation

template<class Real>
virtual ROL::LinearOperator< Real >::~LinearOperator ( )
inlinevirtual

Definition at line 74 of file ROL_LinearOperator.hpp.

Member Function Documentation

template<class Real>
virtual void ROL::LinearOperator< Real >::update ( const Vector< Real > &  x,
bool  flag = true,
int  iter = -1 
)
inlinevirtual

Update linear operator.

This function updates the linear operator at new iterations.

Parameters
[in]xis the new iterate.
[in]flagis true if the iterate has changed.
[in]iteris the outer algorithm iterations count.

Reimplemented in ROL::InteriorPoint::PrimalDualSymmetrizer< Real >, ROL::InteriorPoint::PrimalDualSymmetrizer< Real >, ROL::PrimalDualInteriorPointBlock22< Real >, ROL::PrimalDualInteriorPointBlock21< Real >, ROL::StdLinearOperator< Real >, ROL::LinearOperatorSum< Real >, ROL::LinearOperatorProduct< Real >, ROL::PrimalDualInteriorPointBlock11< Real >, ROL::BlockOperator2Determinant< Real >, and ROL::DiagonalOperator< Real >.

Definition at line 83 of file ROL_LinearOperator.hpp.

template<class Real>
virtual void ROL::LinearOperator< Real >::apply ( Vector< Real > &  Hv,
const Vector< Real > &  v,
Real &  tol 
) const
pure virtual

Apply linear operator.

This function applies the linear operator to a vector.

Parameters
[out]Hvis the output vector.
[in]vis the input vector.
[in]tolis a tolerance for inexact linear operator application.

Implemented in ROL::InteriorPoint::PrimalDualSymmetrizer< Real >, ROL::InteriorPoint::PrimalDualSymmetrizer< Real >, ROL::PrimalDualInteriorPointBlock22< Real >, ROL::BoundFletcher< Real >::AugSystemPrecond, ROL::BoundFletcher< Real >::AugSystemNonSym, ROL::PrimalDualActiveSetStep< Real >::PrecondPD, ROL::BoundFletcher< Real >::AugSystemSym, ROL::PrimalDualInteriorPointBlock21< Real >, ROL::PrimalDualActiveSetStep< Real >::HessianPD, ROL::Secant< Real >, ROL::NullSpaceOperator< Real >, ROL::Fletcher< Real >::AugSystemPrecond, ROL::PrimalDualInteriorPointBlock12< Real >, ROL::RangeSpaceOperator< Real >, ROL::ProjectedNewtonKrylovStep< Real >::PrecondPNK, ROL::Fletcher< Real >::AugSystem, IdentityOperator< Real >, ROL::StdLinearOperator< Real >, ROL::NewtonKrylovStep< Real >::PrecondNK, ROL::ProjectedNewtonKrylovStep< Real >::HessianPNK, ROL::HouseholderReflector< Real >, ROL::LinearOperatorSum< Real >, ROL::PrimalDualInteriorPointBlock11< Real >, ROL::NewtonKrylovStep< Real >::HessianNK, ROL::LinearOperatorProduct< Real >, TridiagonalToeplitzOperator< Real >, TridiagonalToeplitzOperator< Real >, ROL::BlockOperator2< Real >, ROL::BlockOperator2Determinant< Real >, ROL::DiagonalOperator< Real >, ROL::BlockOperator< Real >, ROL::LinearOperatorFromConstraint< Real >, ROL::AugmentedSystemOperator< Real >, ROL::DyadicOperator< Real >, ROL::AugmentedSystemPrecOperator< Real >, Identity< Real >, Identity< Real >, ROL::IdentityOperator< Real >, and ROL::NullOperator< Real >.

Referenced by ROL::LinearOperator< Real >::applyAdjoint(), ROL::Lanczos< Real >::iterate(), ROL::ConjugateGradients< Real >::run(), ROL::ConjugateResiduals< Real >::run(), ROL::GMRES< Real >::run(), and ROL::details::MINRES< Real >::run().

template<class Real>
virtual void ROL::LinearOperator< Real >::applyInverse ( Vector< Real > &  Hv,
const Vector< Real > &  v,
Real &  tol 
) const
inlinevirtual

Apply inverse of linear operator.

This function applies the inverse of linear operator to a vector.

Parameters
[out]Hvis the output vector.
[in]vis the input vector.
[in]tolis a tolerance for inexact linear operator application.

Reimplemented in ROL::InteriorPoint::PrimalDualSymmetrizer< Real >, ROL::InteriorPoint::PrimalDualSymmetrizer< Real >, ROL::PrimalDualInteriorPointBlock22< Real >, ROL::BoundFletcher< Real >::AugSystemPrecond, ROL::PrimalDualActiveSetStep< Real >::PrecondPD, ROL::PrimalDualInteriorPointBlock21< Real >, ROL::NullSpaceOperator< Real >, ROL::PrimalDualInteriorPointBlock12< Real >, ROL::Secant< Real >, ROL::Fletcher< Real >::AugSystemPrecond, ROL::RangeSpaceOperator< Real >, ROL::ProjectedNewtonKrylovStep< Real >::PrecondPNK, ROL::StdLinearOperator< Real >, ROL::PrimalDualInteriorPointBlock11< Real >, ROL::HouseholderReflector< Real >, TridiagonalToeplitzOperator< Real >, TridiagonalToeplitzOperator< Real >, ROL::NewtonKrylovStep< Real >::PrecondNK, ROL::LinearOperatorSum< Real >, ROL::LinearOperatorProduct< Real >, ROL::BlockOperator2< Real >, ROL::BlockOperator2Determinant< Real >, ROL::AugmentedSystemOperator< Real >, ROL::DiagonalOperator< Real >, ROL::LinearOperatorFromConstraint< Real >, ROL::AugmentedSystemPrecOperator< Real >, ROL::DyadicOperator< Real >, ROL::IdentityOperator< Real >, and ROL::NullOperator< Real >.

Definition at line 101 of file ROL_LinearOperator.hpp.

References ROL::Vector< Real >::set().

Referenced by ROL::LinearOperator< Real >::applyAdjointInverse(), ROL::ConjugateGradients< Real >::run(), ROL::ConjugateResiduals< Real >::run(), and ROL::GMRES< Real >::run().

template<class Real>
virtual void ROL::LinearOperator< Real >::applyAdjoint ( Vector< Real > &  Hv,
const Vector< Real > &  v,
Real &  tol 
) const
inlinevirtual

Apply adjoint of linear operator.

This function applies the adjoint of a linear operator to a vector. Default behavior assumes operator is self-adjoint.

Parameters
[out]Hvis the output vector.
[in]vis the input vector.
[in]tolis a tolerance for inexact linear operator application.

Reimplemented in ROL::NullSpaceOperator< Real >, ROL::RangeSpaceOperator< Real >, ROL::StdLinearOperator< Real >, ROL::AugmentedSystemOperator< Real >, and ROL::AugmentedSystemPrecOperator< Real >.

Definition at line 114 of file ROL_LinearOperator.hpp.

References ROL::LinearOperator< Real >::apply().

template<class Real>
virtual void ROL::LinearOperator< Real >::applyAdjointInverse ( Vector< Real > &  Hv,
const Vector< Real > &  v,
Real &  tol 
) const
inlinevirtual

Apply adjoint of the inverse linear operator.

This function applies the adjoint of the inverse linear operator to a vector. Default behavior assumes operator is self-adjoint.

Parameters
[out]Hvis the output vector.
[in]vis the input vector.
[in]tolis a tolerance for inexact linear operator application.

Reimplemented in ROL::StdLinearOperator< Real >, ROL::NullSpaceOperator< Real >, ROL::RangeSpaceOperator< Real >, ROL::AugmentedSystemOperator< Real >, and ROL::AugmentedSystemPrecOperator< Real >.

Definition at line 127 of file ROL_LinearOperator.hpp.

References ROL::LinearOperator< Real >::applyInverse().


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