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::DefaultScaledAdjointLinearOp< Scalar > Class Template Reference

Concrete decorator LinearOpBase subclass that wraps a LinearOpBase object and adds on an extra scaling factor and/or a transpose enum. More...

#include <Thyra_DefaultScaledAdjointLinearOp_decl.hpp>

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

Related Functions

(Note that these are not member functions.)

template<class Scalar >
RCP< LinearOpBase< Scalar > > nonconstScale (const Scalar &scalar, const RCP< LinearOpBase< Scalar > > &Op, const std::string &label="")
 Build an implicit non-const scaled linear operator. More...
 
template<class Scalar >
RCP< const LinearOpBase< Scalar > > scale (const Scalar &scalar, const RCP< const LinearOpBase< Scalar > > &Op, const std::string &label="")
 Build an implicit const scaled linear operator. More...
 
template<class Scalar >
RCP< LinearOpBase< Scalar > > nonconstAdjoint (const RCP< LinearOpBase< Scalar > > &Op, const std::string &label="")
 Build an implicit non-const adjoined linear operator. More...
 
template<class Scalar >
RCP< const LinearOpBase< Scalar > > adjoint (const RCP< const LinearOpBase< Scalar > > &Op, const std::string &label="")
 Build an implicit const adjoined linear operator. More...
 
template<class Scalar >
RCP< LinearOpBase< Scalar > > nonconstTranspose (const RCP< LinearOpBase< Scalar > > &Op, const std::string &label="")
 Build an implicit non-const transposed linear operator. More...
 
template<class Scalar >
RCP< const LinearOpBase< Scalar > > transpose (const RCP< const LinearOpBase< Scalar > > &Op, const std::string &label="")
 Build an implicit const transposed linear operator. More...
 
template<class Scalar >
RCP< LinearOpBase< Scalar > > nonconstScaleAndAdjoint (const Scalar &scalar, const EOpTransp &transp, const RCP< LinearOpBase< Scalar > > &Op, const std::string &label="")
 Build an implicit non-const scaled and/or adjoined (transposed) linear operator. More...
 
template<class Scalar >
RCP< const LinearOpBase< Scalar > > scaleAndAdjoint (const Scalar &scalar, const EOpTransp &transp, const RCP< const LinearOpBase< Scalar > > &Op, const std::string &label="")
 Build an implicit const scaled and/or adjoined (transposed) linear operator. More...
 

Constructors/initializers/accessors

 DefaultScaledAdjointLinearOp ()
 Constructs to uninitialized. More...
 
 DefaultScaledAdjointLinearOp (const Scalar &scalar, const EOpTransp &transp, const RCP< LinearOpBase< Scalar > > &Op)
 Calls initialize(). More...
 
 DefaultScaledAdjointLinearOp (const Scalar &scalar, const EOpTransp &transp, const RCP< const LinearOpBase< Scalar > > &Op)
 Calls initialize(). More...
 
void initialize (const Scalar &scalar, const EOpTransp &transp, const RCP< LinearOpBase< Scalar > > &Op)
 Initialize with an operator with by defining adjoint (transpose) and scaling arguments. More...
 
void initialize (const Scalar &scalar, const EOpTransp &transp, const RCP< const LinearOpBase< Scalar > > &Op)
 Initialize with an operator with by defining adjoint (transpose) and scaling arguments. More...
 
RCP< LinearOpBase< Scalar > > getNonconstOp ()
 Return the non-const linear operator passed into initialize(). More...
 
RCP< const LinearOpBase< Scalar > > getOp () const
 Return the const linear operator passed into initialize(). More...
 
void uninitialize ()
 Set to uninitialized and (optionally) extract the objects passed into initialize(). More...
 

Overridden from Teuchos::Describable

std::string description () const
 Outputs DefaultScaledAdjointLinearOp<Scalar>{this->getOrigOp().description()) along with the dimensions. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
 Prints out the original operator as well as all of the scalings and transpositions in the order that they occurred. More...
 

Overridden from LinearOpBase

RCP< const VectorSpaceBase
< Scalar > > 
range () const
 Return the range space of the logical linear operator. More...
 
RCP< const VectorSpaceBase
< Scalar > > 
domain () const
 Return the domain space of the logical linear operator. More...
 
RCP< const LinearOpBase< Scalar > > clone () const
 
bool opSupportedImpl (EOpTransp M_trans) const
 Return if the operation is supported on the logical linear operator. More...
 
void applyImpl (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
 Apply the linear operator (or its transpose) to a multi-vector : Y = alpha*op(M)*X + beta*Y. More...
 

Overridden from ScaledAdointLinearOpBase

Scalar overallScalar () const
 
EOpTransp overallTransp () const
 
RCP< LinearOpBase< Scalar > > getNonconstOrigOp ()
 
RCP< const LinearOpBase< Scalar > > getOrigOp () const
 

Detailed Description

template<class Scalar>
class Thyra::DefaultScaledAdjointLinearOp< Scalar >

Concrete decorator LinearOpBase subclass that wraps a LinearOpBase object and adds on an extra scaling factor and/or a transpose enum.

This class represents a scaled, adjointed (transposed) linear operator M of the form:

  M = scalar * op(Op)

where Op is another LinearOpBase object, scalar is a Scalar, and the operation op(Op) is specified by a EOpTransp and is given as op(Op) = Op (NOTRANS), or op(Op) = Op^T (TRANS), or op(Op) = Op^H (CONJTRANS). Of course the operator M is not constructed explicitly but instead just applies the decorated operator Op by modifying the apply() function that calls Op.apply().

This subclass is designed to allow the efficient handling of multiple implicit scalings and/or adjoints (transposes) and allow these implicit transformations to be reversed. A sequence of scalings/adjoints from some original LinearOpBase object origOp is shown as:

 M = scalar_n * op_n( ... scalar_2 * op_2( scalar_1 * op_1( origOp ) ) ... )

   =>

 M  = overallScalar * overall_op(origOp)

where overallScalar and overall_op(...) are the cumulative transformations given as:

 overallScalar = scalar_n * ... * scalar_2 * ... * scalar_1

 overall_op(origOp) = op_n( ... op_2( op_1( ... op_n( origOp ) ... ) ) )

Each individual transformation pair (scalar_i,op_i(...)) is specified with arguments Scalar scalar and EOpTransp transp. The overall scaling is returned using this->overallScalar(), the overall adjoint enum is returned using this->overallTransp(), and the original linear operator is returned from this->getOrigOp().

The operator most recently wrapped is returned from this->getOp(). The individual scalings and transformations are not exposed from this interface but can be viewed by calling this->description() with a verbosity level of ???. The arguments passed into the constructor DefaultScaledAdjointLinearOp() or initialize() can always be extracted using uninitialize().

This subclass keeps track of all of the individual scalings scalar_i and adjoining operations op_i(...) so that the description() function can print this out and also so that the operations can be reversed.

The copy constructor and assignment operators are declared private since some thought needs to be put into what they should mean.

Note: This class does not maintain any specific cached information about the original operator so it is safe, in some cases, to manipulate the original operator and still keep this intact and automatically updated for any changes that are made.

Definition at line 126 of file Thyra_DefaultScaledAdjointLinearOp_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Thyra::DefaultScaledAdjointLinearOp< Scalar >::DefaultScaledAdjointLinearOp ( )
inline

Constructs to uninitialized.

Postconditions:

Definition at line 581 of file Thyra_DefaultScaledAdjointLinearOp_decl.hpp.

template<class Scalar>
Thyra::DefaultScaledAdjointLinearOp< Scalar >::DefaultScaledAdjointLinearOp ( const Scalar &  scalar,
const EOpTransp transp,
const RCP< LinearOpBase< Scalar > > &  Op 
)
inline

Calls initialize().

Note, instead of calling this constructor directly consider using the non-member functions described belos which create dynamically allocated RCP-wrapped objects.

Definition at line 590 of file Thyra_DefaultScaledAdjointLinearOp_decl.hpp.

template<class Scalar>
Thyra::DefaultScaledAdjointLinearOp< Scalar >::DefaultScaledAdjointLinearOp ( const Scalar &  scalar,
const EOpTransp transp,
const RCP< const LinearOpBase< Scalar > > &  Op 
)
inline

Calls initialize().

Note, instead of calling this constructor directly consider using the non-member functions described belos which create dynamically allocated RCP-wrapped objects.

Definition at line 604 of file Thyra_DefaultScaledAdjointLinearOp_decl.hpp.

Member Function Documentation

template<class Scalar>
void Thyra::DefaultScaledAdjointLinearOp< Scalar >::initialize ( const Scalar &  scalar,
const EOpTransp transp,
const RCP< LinearOpBase< Scalar > > &  Op 
)

Initialize with an operator with by defining adjoint (transpose) and scaling arguments.

Parameters
scalar[in] Scalar argument defining scalar_0 (see introduction).
transp[in] Transpose argument defining op_0(...) (see introduction).
Op[in] Smart pointer to linear operator (persisting relationship).

Preconditions:

  • Op.get() != NULL

Postconditions:

  • ToDo: Fill these in!!!!

Definition at line 59 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar>
void Thyra::DefaultScaledAdjointLinearOp< Scalar >::initialize ( const Scalar &  scalar,
const EOpTransp transp,
const RCP< const LinearOpBase< Scalar > > &  Op 
)

Initialize with an operator with by defining adjoint (transpose) and scaling arguments.

Parameters
scalar[in] Scalar argument defining scalar_0 (see introduction).
transp[in] Transpose argument defining op_0(...) (see introduction).
Op[in] Smart pointer to linear operator (persisting relationship).

Preconditions:

  • Op.get() != NULL

Postconditions:

  • ToDo: Fill these in!!!!

Definition at line 70 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar >
Teuchos::RCP< LinearOpBase< Scalar > > Thyra::DefaultScaledAdjointLinearOp< Scalar >::getNonconstOp ( )

Return the non-const linear operator passed into initialize().

Definition at line 82 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar >
Teuchos::RCP< const LinearOpBase< Scalar > > Thyra::DefaultScaledAdjointLinearOp< Scalar >::getOp ( ) const

Return the const linear operator passed into initialize().

Definition at line 90 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar >
void Thyra::DefaultScaledAdjointLinearOp< Scalar >::uninitialize ( )

Set to uninitialized and (optionally) extract the objects passed into initialize().

Postconditions:

Definition at line 97 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar >
std::string Thyra::DefaultScaledAdjointLinearOp< Scalar >::description ( ) const
virtual

Outputs DefaultScaledAdjointLinearOp<Scalar>{this->getOrigOp().description()) along with the dimensions.

Reimplemented from Teuchos::Describable.

Definition at line 112 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar >
void Thyra::DefaultScaledAdjointLinearOp< Scalar >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel 
) const
virtual

Prints out the original operator as well as all of the scalings and transpositions in the order that they occurred.

This function outputs different levels of detail based on the value passed in for verbLevel:

ToDo: Finish documentation!

Reimplemented from Teuchos::Describable.

Definition at line 124 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar >
Teuchos::RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultScaledAdjointLinearOp< Scalar >::range ( ) const
virtual

Return the range space of the logical linear operator.

Simply returns:

return ( this->overallTransp()==NOTRANS ? this->getOrigOp()->range() : this->getOrigOp()->domain() );

Implements Thyra::LinearOpBase< Scalar >.

Definition at line 181 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar >
Teuchos::RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultScaledAdjointLinearOp< Scalar >::domain ( ) const
virtual

Return the domain space of the logical linear operator.

Simply returns:

return ( this->overallTransp()==NOTRANS ? this->getOrigOp()->domain() : this->getOrigOp()->range() );

Implements Thyra::LinearOpBase< Scalar >.

Definition at line 191 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

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

Reimplemented from Thyra::LinearOpBase< Scalar >.

Definition at line 201 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar >
Scalar Thyra::DefaultScaledAdjointLinearOp< Scalar >::overallScalar ( ) const
virtual
template<class Scalar >
EOpTransp Thyra::DefaultScaledAdjointLinearOp< Scalar >::overallTransp ( ) const
virtual
template<class Scalar >
Teuchos::RCP< LinearOpBase< Scalar > > Thyra::DefaultScaledAdjointLinearOp< Scalar >::getNonconstOrigOp ( )
virtual
template<class Scalar >
Teuchos::RCP< const LinearOpBase< Scalar > > Thyra::DefaultScaledAdjointLinearOp< Scalar >::getOrigOp ( ) const
virtual
template<class Scalar >
bool Thyra::DefaultScaledAdjointLinearOp< Scalar >::opSupportedImpl ( EOpTransp  M_trans) const
protectedvirtual

Return if the operation is supported on the logical linear operator.

Simply returns:

return this->getOrigOp()->opSupported(trans_trans(this->overallTransp(),M_trans));

Implements Thyra::LinearOpBase< Scalar >.

Definition at line 247 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

template<class Scalar>
void Thyra::DefaultScaledAdjointLinearOp< Scalar >::applyImpl ( const EOpTransp  M_trans,
const MultiVectorBase< Scalar > &  X,
const Ptr< MultiVectorBase< Scalar > > &  Y,
const Scalar  alpha,
const Scalar  beta 
) const
protectedvirtual

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

Simply calls:

this->getOrigOp()->apply(trans_trans(M_trans,this->overallTransp()),X,Y,(this->overallScalar()*alpha),beta)

Implements Thyra::LinearOpBase< Scalar >.

Definition at line 256 of file Thyra_DefaultScaledAdjointLinearOp_def.hpp.

Friends And Related Function Documentation

template<class Scalar >
RCP< LinearOpBase< Scalar > > nonconstScale ( const Scalar &  scalar,
const RCP< LinearOpBase< Scalar > > &  Op,
const std::string &  label = "" 
)
related

Build an implicit non-const scaled linear operator.

Returns Teuchos::rcp(new DefaultScaledAdjointLinearOp<Scalar>(scalar,NOTRANS,Op).

Preconditions:

  • Op.get()!=NULL

Postconditions:

  • return.get()!=NULL
template<class Scalar >
RCP< const LinearOpBase< Scalar > > scale ( const Scalar &  scalar,
const RCP< const LinearOpBase< Scalar > > &  Op,
const std::string &  label = "" 
)
related

Build an implicit const scaled linear operator.

Returns Teuchos::rcp(new DefaultScaledAdjointLinearOp<Scalar>(scalar,NOTRANS,Op).

Preconditions:

  • Op.get()!=NULL

Postconditions:

  • return.get()!=NULL
template<class Scalar >
RCP< LinearOpBase< Scalar > > nonconstAdjoint ( const RCP< LinearOpBase< Scalar > > &  Op,
const std::string &  label = "" 
)
related

Build an implicit non-const adjoined linear operator.

Returns Teuchos::rcp(new DefaultScaledAdjointLinearOp<Scalar>(Teuchos::ScalarTraits<Scalar>::one(),CONJTRANS,Op).

Preconditions:

  • Op.get()!=NULL

Postconditions:

  • return.get()!=NULL
template<class Scalar >
RCP< const LinearOpBase< Scalar > > adjoint ( const RCP< const LinearOpBase< Scalar > > &  Op,
const std::string &  label = "" 
)
related

Build an implicit const adjoined linear operator.

Returns Teuchos::rcp(new DefaultScaledAdjointLinearOp<Scalar>(Teuchos::ScalarTraits<Scalar>::one(),CONJTRANS,Op).

Preconditions:

  • Op.get()!=NULL

Postconditions:

  • return.get()!=NULL
template<class Scalar >
RCP< LinearOpBase< Scalar > > nonconstTranspose ( const RCP< LinearOpBase< Scalar > > &  Op,
const std::string &  label = "" 
)
related

Build an implicit non-const transposed linear operator.

Returns Teuchos::rcp(new DefaultScaledAdjointLinearOp<Scalar>(Teuchos::ScalarTraits<Scalar>::one(),TRANS,Op).

Preconditions:

  • Op.get()!=NULL

Postconditions:

  • return.get()!=NULL
template<class Scalar >
RCP< const LinearOpBase< Scalar > > transpose ( const RCP< const LinearOpBase< Scalar > > &  Op,
const std::string &  label = "" 
)
related

Build an implicit const transposed linear operator.

Returns Teuchos::rcp(new DefaultScaledAdjointLinearOp<Scalar>(Teuchos::ScalarTraits<Scalar>::one(),TRANS,Op).

Preconditions:

  • Op.get()!=NULL

Postconditions:

  • return.get()!=NULL
template<class Scalar >
RCP< LinearOpBase< Scalar > > nonconstScaleAndAdjoint ( const Scalar &  scalar,
const EOpTransp transp,
const RCP< LinearOpBase< Scalar > > &  Op,
const std::string &  label = "" 
)
related

Build an implicit non-const scaled and/or adjoined (transposed) linear operator.

Returns Teuchos::rcp(new DefaultScaledAdjointLinearOp<Scalar>(scale,transp,Op).

Preconditions:

  • Op.get()!=NULL

Postconditions:

  • return.get()!=NULL
template<class Scalar >
RCP< const LinearOpBase< Scalar > > scaleAndAdjoint ( const Scalar &  scalar,
const EOpTransp transp,
const RCP< const LinearOpBase< Scalar > > &  Op,
const std::string &  label = "" 
)
related

Build an implicit const scaled and/or adjoined (transposed) linear operator.

Returns Teuchos::rcp(new DefaultScaledAdjointLinearOp<Scalar>(scale,transp,Op).

Preconditions:

  • Op.get()!=NULL

Postconditions:

  • return.get()!=NULL

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