Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Related Functions | List of all members
Thyra::LinearOpBase< Scalar > Class Template Referenceabstract

Base class for all linear operators. More...

#include <Thyra_LinearOpBase_decl.hpp>

Inheritance diagram for Thyra::LinearOpBase< Scalar >:
Inheritance graph
[legend]

Related Functions

(Note that these are not member functions.)

template<class Scalar >
bool isFullyUninitialized (const LinearOpBase< Scalar > &M)
 Determines if a linear operator is in the "Fully Uninitialized" state or not. More...
 
template<class Scalar >
bool isPartiallyInitialized (const LinearOpBase< Scalar > &M)
 Determines if a linear operator is in the "Partially Initialized" state or not. More...
 
template<class Scalar >
bool isFullyInitialized (const LinearOpBase< Scalar > &M)
 Determines if a linear operator is in the "Fully Initialized" state or not. More...
 
template<class Scalar >
bool opSupported (const LinearOpBase< Scalar > &M, EOpTransp M_trans)
 Determines if an operation is supported for a single scalar type. More...
 
template<class Scalar >
void apply (const LinearOpBase< Scalar > &M, const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha=static_cast< Scalar >(1.0), const Scalar beta=static_cast< Scalar >(0.0))
 Non-member function call for M.apply(...). More...
 
void apply (const LinearOpBase< double > &M, const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha=1.0, const double beta=0.0)
 Calls apply<double>(...). More...
 

Public interface functions

virtual RCP< const
VectorSpaceBase< Scalar > > 
range () const =0
 Return a smart pointer for the range space for this operator. More...
 
virtual RCP< const
VectorSpaceBase< Scalar > > 
domain () const =0
 Return a smart pointer for the domain space for this operator. More...
 
bool opSupported (EOpTransp M_trans) const
 Return if the M_trans operation of apply() is supported or not. More...
 
void apply (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
 Apply the linear operator to a multi-vector : Y = alpha*op(M)*X + beta*Y. More...
 
virtual RCP< const
LinearOpBase< Scalar > > 
clone () const
 Clone the linear operator object (if supported). More...
 

Protected virtual functions to be overridden by subclasses.

virtual bool opSupportedImpl (EOpTransp M_trans) const =0
 Override in subclass. More...
 
virtual void applyImpl (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const =0
 Override in subclass. More...
 

Detailed Description

template<class Scalar>
class Thyra::LinearOpBase< Scalar >

Base class for all linear operators.

Outline

Introduction

A linear operator can optionally perform following operations

through the apply() function where Y and X are MultiVectorBase objects. The reason for the exact form of the above operations is that there are direct BLAS and equivalent versions of these operations and performing a sum-into multiplication is more efficient in general.

Range and domain spaces

A linear operator has vector spaces associated with it for the vectors x and y that lie in the domain and the range spaces of the non-transposed linear operator y = M*x. These spaces are returned by domain() and range().

Scalar products and the adjoint relation

Note that the vector spaces returned from domain() and range() may have specialized implementations of the scalar product $<u,w>$ (i.e. $<u,w> \neq u^H w$ in general). As a result, the operator and adjoint operator must obey the defined scalar products. Specifically, for any two vectors $w\in\mathcal{D}$ (in the domain space of A) and $u\in\mathcal{R}$ (in the range space of A), the adjoint operation must obey the adjoint property

\[ <u,A v>_{\mathcal{R}} =\!= <A^H u, v>_{\mathcal{D}} \]

where $<.,.>_{\mathcal{R}}$ is the scalar product defined by this->range()->scalarProd() and $<.,.>_{\mathcal{D}}$ is the scalar product defined by this->domain()->scalarProd(). This property of the adjoint can be checked numerically, if adjoints are supported, using the testing utility class LinearOpTester.

Aliasing policy

It is strictly forbidden to alias the input/output object Y with the input object X in apply(). Allowing aliasing would greatly complicate the development of concrete subclasses.

Optional support for specific types of operator applications

This interface does not require that a linear operator implementation support all of the different types of operator applications defined in the introduction above. If a LinearOpBase object can not support a particular type of operator application, then this is determined by the functions opSupported().

Testing LinearOpBase objects

The concrete class LinearOpTester provides a full featured set of tests for any LinearOpBase object. This testing class can check if the operator is truly "linear", and/or if the adjoint relationship holds, and/or if an operator is symmetric. All of the tests are controlled by the client, can be turned on and off, and pass/failure is determined by tolerances that the client can specify. In addition, this testing class can also check if two linear operators are approximately equal.

Initialization states

A LinearOpBase object has three different states of initialization. These three initailziation states, a description of their definition, and non-member helper functions that return these states are given below:

These three different states of initialization allow for the simplification of the implementation of many different types of use cases.

Notes for subclass developers

There are only foure functions that a concrete subclass is required to override: domain(), range() opSupportedImpl(), and applyImpl(). Note that the functions domain() and range() should simply return VectorSpaceBase objects for subclasses that are already defined for the vectors that the linear operator interacts with through the function apply(). The function opSupportedImpl() just returns what operations are supported and is therefore trivial to implement. Therefore, given that appropriate VectorSpaceBase and MultiVectorBase (and/or VectorBase) subclasses exist, the only real work involved in implementing a LinearOpBase subclass is in defining a single function applyImpl().

If possible, the subclass should also override the clone() function which allows clients to create copies of a LinearOpBase object. This functionality is useful in some circumstances. However, this functionality is not required and the default clone() implementation returns a null smart pointer object.

Definition at line 191 of file Thyra_LinearOpBase_decl.hpp.

Member Function Documentation

template<class Scalar>
virtual RCP< const VectorSpaceBase<Scalar> > Thyra::LinearOpBase< Scalar >::range ( ) const
pure virtual

Return a smart pointer for the range space for this operator.

Note that a return value of is_null(returnVal) is a flag that *this is not fully initialized.

If nonnull(returnVal), it is required that the object referenced by *returnVal must have lifetime that extends past the lifetime of the returned smart pointer object. However, the object referenced by *returnVal may change if *this modified so this reference should not be maintained for too long.

New Behavior! It is required that the VectorSpaceBase object embedded in return must be valid past the lifetime of *this linear operator object.

Implemented in Thyra::EpetraLinearOp, Thyra::DefaultScaledAdjointLinearOp< Scalar >, Thyra::DefaultInverseLinearOp< Scalar >, Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >, Thyra::DefaultDiagonalLinearOp< Scalar >, Thyra::DefaultMultipliedLinearOp< Scalar >, ExampleTridiagSpmdLinearOp< Scalar >, Thyra::DefaultAddedLinearOp< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::DefaultBlockedLinearOp< Scalar >, Thyra::DefaultProductMultiVector< Scalar >, ExampleTridiagSerialLinearOp< Scalar >, Thyra::DefaultZeroLinearOp< Scalar >, Thyra::DefaultIdentityLinearOp< Scalar >, Thyra::VectorDefaultBase< Scalar >, Thyra::TpetraLinearOp< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultSerialDenseLinearOpWithSolve< Scalar >, Thyra::DelayedLinearOpWithSolve< Scalar >, Thyra::DefaultMultiVectorLinearOpWithSolve< Scalar >, Thyra::MultiVectorAdapterBase< Scalar >, and Thyra::DefaultAdjointLinearOpWithSolve< Scalar >.

template<class Scalar>
virtual RCP< const VectorSpaceBase<Scalar> > Thyra::LinearOpBase< Scalar >::domain ( ) const
pure virtual

Return a smart pointer for the domain space for this operator.

Note that a return value of is_null(returnVal) is a flag that *this is not fully initialized.

If nonnull(returnVal), it is required that the object referenced by *returnVal must have lifetime that extends past the lifetime of the returned smart pointer object. However, the object referenced by *returnVal may change if *this modified so this reference should not be maintained for too long.

New Behavior! It is required that the VectorSpaceBase object embedded in return must be valid past the lifetime of *this linear operator object.

Implemented in Thyra::EpetraLinearOp, Thyra::DefaultScaledAdjointLinearOp< Scalar >, Thyra::DefaultInverseLinearOp< Scalar >, Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >, Thyra::DefaultDiagonalLinearOp< Scalar >, Thyra::DefaultMultipliedLinearOp< Scalar >, ExampleTridiagSpmdLinearOp< Scalar >, Thyra::DefaultAddedLinearOp< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::DefaultBlockedLinearOp< Scalar >, Thyra::DefaultProductMultiVector< Scalar >, ExampleTridiagSerialLinearOp< Scalar >, Thyra::DefaultZeroLinearOp< Scalar >, Thyra::VectorDefaultBase< Scalar >, Thyra::DefaultIdentityLinearOp< Scalar >, Thyra::TpetraLinearOp< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultSerialDenseLinearOpWithSolve< Scalar >, Thyra::DelayedLinearOpWithSolve< Scalar >, Thyra::DefaultMultiVectorLinearOpWithSolve< Scalar >, Thyra::MultiVectorAdapterBase< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Thyra::DefaultAdjointLinearOpWithSolve< Scalar >.

template<class Scalar>
bool Thyra::LinearOpBase< Scalar >::opSupported ( EOpTransp  M_trans) const
inline

Return if the M_trans operation of apply() is supported or not.

Preconditions:

  • isPartiallyInitialized(*this)

Note that an operator must support at least one of the values of ETrans (i.e. the transposed or the non-transposed operations must be supported, both can not be unsupported)

Definition at line 242 of file Thyra_LinearOpBase_decl.hpp.

template<class Scalar>
void Thyra::LinearOpBase< Scalar >::apply ( const EOpTransp  M_trans,
const MultiVectorBase< Scalar > &  X,
const Ptr< MultiVectorBase< Scalar > > &  Y,
const Scalar  alpha,
const Scalar  beta 
) const
inline

Apply the linear operator to a multi-vector : Y = alpha*op(M)*X + beta*Y.

Parameters
M_trans[in] Determines whether the operator is applied or the adjoint for op(M).
X[in] The right hand side multi-vector.
Y[in/out] The target multi-vector being transformed. When beta==0.0, this multi-vector can have uninitialized elements.
alpha[in] Scalar multiplying M, where M==*this. The default value of alpha is 1.0
beta[in] The multiplier for the target multi-vector Y. The default value of beta is 0.0.

Preconditions:

Postconditions:

  • Is it not obvious? After the function returns the multi-vector Y is transformed as indicated above.

Definition at line 292 of file Thyra_LinearOpBase_decl.hpp.

template<class Scalar >
RCP< const LinearOpBase< Scalar > > Thyra::LinearOpBase< Scalar >::clone ( ) const
virtual

Clone the linear operator object (if supported).

The primary purpose for this function is to allow a client to capture the current state of a linear operator object and be guaranteed that some other client will not alter its behavior. A smart implementation will use reference counting and lazy evaluation internally and will not actually copy any large amount of data unless it has to.

The default implementation returns is_null(returnVal) which is allowable. A linear operator object is not required to return a non-NULL value but many good matrix-based linear operator implementations will.

Reimplemented in Thyra::MultiVectorBase< Scalar >, Thyra::EpetraLinearOp, Thyra::DefaultScaledAdjointLinearOp< Scalar >, Thyra::DefaultInverseLinearOp< Scalar >, Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >, Thyra::DefaultDiagonalLinearOp< Scalar >, Thyra::DefaultMultipliedLinearOp< Scalar >, Thyra::DefaultAddedLinearOp< Scalar >, Thyra::DefaultBlockedLinearOp< Scalar >, Thyra::DefaultZeroLinearOp< Scalar >, Thyra::DefaultIdentityLinearOp< Scalar >, Thyra::DelayedLinearOpWithSolve< Scalar >, and Thyra::DefaultMultiVectorLinearOpWithSolve< Scalar >.

Definition at line 58 of file Thyra_LinearOpBase_def.hpp.

template<class Scalar>
virtual bool Thyra::LinearOpBase< Scalar >::opSupportedImpl ( EOpTransp  M_trans) const
protectedpure virtual
template<class Scalar>
virtual void Thyra::LinearOpBase< Scalar >::applyImpl ( const EOpTransp  M_trans,
const MultiVectorBase< Scalar > &  X,
const Ptr< MultiVectorBase< Scalar > > &  Y,
const Scalar  alpha,
const Scalar  beta 
) const
protectedpure virtual

Friends And Related Function Documentation

template<class Scalar >
bool isFullyUninitialized ( const LinearOpBase< Scalar > &  M)
related

Determines if a linear operator is in the "Fully Uninitialized" state or not.

template<class Scalar >
bool isPartiallyInitialized ( const LinearOpBase< Scalar > &  M)
related

Determines if a linear operator is in the "Partially Initialized" state or not.

template<class Scalar >
bool isFullyInitialized ( const LinearOpBase< Scalar > &  M)
related

Determines if a linear operator is in the "Fully Initialized" state or not.

template<class Scalar >
bool opSupported ( const LinearOpBase< Scalar > &  M,
EOpTransp  M_trans 
)
related

Determines if an operation is supported for a single scalar type.

template<class Scalar >
void apply ( const LinearOpBase< Scalar > &  M,
const EOpTransp  M_trans,
const MultiVectorBase< Scalar > &  X,
const Ptr< MultiVectorBase< Scalar > > &  Y,
const Scalar  alpha = static_cast< Scalar >(1.0),
const Scalar  beta = static_cast< Scalar >(0.0) 
)
related

Non-member function call for M.apply(...).

template<class Scalar>
void apply ( const LinearOpBase< double > &  M,
const EOpTransp  M_trans,
const MultiVectorBase< double > &  X,
const Ptr< MultiVectorBase< double > > &  Y,
const double  alpha = 1.0,
const double  beta = 0.0 
)
related

Calls apply<double>(...).

Non-tempalted double inlined non-member helper function.


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