|
Thyra
Version of the Day
|
Concrete composite LinearOpBase subclass that creates an implicitly multiplied linear operator out of one or more constituent LinearOpBase objects.
More...
#include <Thyra_DefaultMultipliedLinearOp_decl.hpp>

Related Functions | |
(Note that these are not member functions.) | |
| template<class Scalar > | |
| RCP< DefaultMultipliedLinearOp < Scalar > > | defaultMultipliedLinearOp () |
| Nonmember constructor. More... | |
| template<class Scalar > | |
| RCP< DefaultMultipliedLinearOp < Scalar > > | defaultMultipliedLinearOp (const ArrayView< const RCP< LinearOpBase< Scalar > > > &Ops) |
| Nonmember constructor. More... | |
| template<class Scalar > | |
| RCP< DefaultMultipliedLinearOp < Scalar > > | defaultMultipliedLinearOp (const ArrayView< const RCP< const LinearOpBase< Scalar > > > &Ops) |
| Nonmember constructor. More... | |
| template<class Scalar > | |
| RCP< LinearOpBase< Scalar > > | nonconstMultiply (const RCP< LinearOpBase< Scalar > > &A, const RCP< LinearOpBase< Scalar > > &B, const std::string &M_label="") |
Form an implicit multiplication of two linear operators: M = A. More... | |
| template<class Scalar > | |
| RCP< const LinearOpBase< Scalar > > | multiply (const RCP< const LinearOpBase< Scalar > > &A, const RCP< const LinearOpBase< Scalar > > &B, const std::string &M_label="") |
Form an implicit multiplication of two linear operators: M = A. More... | |
| template<class Scalar > | |
| RCP< const LinearOpBase< Scalar > > | multiply (const RCP< const LinearOpBase< Scalar > > &A, const RCP< const LinearOpBase< Scalar > > &B, const RCP< const LinearOpBase< Scalar > > &C, const std::string &M_label="") |
Form an implicit multiplication of three linear operators: M = A * B * C. More... | |
Related Functions inherited from Thyra::LinearOpBase< Scalar > | |
| 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... | |
Constructors/initializers/accessors | |
| DefaultMultipliedLinearOp () | |
| Constructs to uninitialized. More... | |
| void | initialize (const ArrayView< const RCP< LinearOpBase< Scalar > > > &Ops) |
| Initialize given a list of non-const linear operators. More... | |
| void | initialize (const ArrayView< const RCP< const LinearOpBase< Scalar > > > &Ops) |
| Initialize given a list of const linear operators. More... | |
| void | uninitialize () |
| Set to uninitialized. More... | |
Overridden from MultipliedLinearOpBase | |
| int | numOps () const |
| bool | opIsConst (const int k) const |
| RCP< LinearOpBase< Scalar > > | getNonconstOp (const int k) |
| RCP< const LinearOpBase< Scalar > > | getOp (const int k) const |
Overridden from LinearOpBase | |
| RCP< const VectorSpaceBase < Scalar > > | range () const |
Returns this->getOp(0).range() if <t>this->numOps() > 0 and returns Teuchos::null otherwise. More... | |
| RCP< const VectorSpaceBase < Scalar > > | domain () const |
Returns this->getOp(this->numOps()-1).domain() if <t>this->numOps() > 0 and returns Teuchos::null otherwise. More... | |
| RCP< const LinearOpBase< Scalar > > | clone () const |
| bool | opSupportedImpl (EOpTransp M_trans) const |
Returns true only if all constituent operators support M_trans. More... | |
| void | applyImpl (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const |
Overridden from Teuchos::Describable | |
| std::string | description () const |
Prints just the name DefaultMultipliedLinearOp along with the overall dimensions and the number of constituent operators. More... | |
| void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
| Prints the details about the constituent linear operators. More... | |
Concrete composite LinearOpBase subclass that creates an implicitly multiplied linear operator out of one or more constituent LinearOpBase objects.
This class represents a multiplied linear operator M of the form:
M = Op[0] * Op[1] * ... * Op[numOps-1]
where Op[] is an array of numOps LinearOp objects. Of course the operator M is not constructed explicitly but instead just applies the constituent linear operators accordingly using temporaries.
In other words, this class defines apply() as:
y = alpha*M*x + beta*y = alpha * ( Op[0] * ( Op[1] * ( .... ( Op[numOps-1] * x ) ... ) ) ) + beta * y
for the case where M_trans==NOTRANS and as:
y = alpha*M'*x + beta*y = alpha * ( Op[numOps-1]' * ( .... ( Op[1]' * ( Op[0]' * x ) ) ... ) ) + beta * y
for the case where real_trans(M_trans)!=NOTRANS (where the transpose ' either defines TRANS or CONJTRANS).
Constructing a multiplied operator is easy. For example, suppose one wants to construct the multiplied operator D = A * B' * C. To do so one would do:
Rather than calling the constructor directly, consider using the non-member helper functions described here.
Definition at line 120 of file Thyra_DefaultMultipliedLinearOp_decl.hpp.
| Thyra::DefaultMultipliedLinearOp< Scalar >::DefaultMultipliedLinearOp | ( | ) |
Constructs to uninitialized.
Postconditions:
this->numOps()==0 Definition at line 93 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
| void Thyra::DefaultMultipliedLinearOp< Scalar >::initialize | ( | const ArrayView< const RCP< LinearOpBase< Scalar > > > & | Ops | ) |
Initialize given a list of non-const linear operators.
| Ops | [in] Array (length numOps) of constituent linear operators and their aggregated default definitions of the non-transposed operator. |
Preconditions:
numOps > 0 Ops != NULL Ops[k].op().get()!=NULL, for k=0...numOps-1 Postconditions:
this->numOps()==numOps this->getOp(k).op().get()==Ops[k].op().get(), for k=0...numOps-1 Definition at line 98 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
| void Thyra::DefaultMultipliedLinearOp< Scalar >::initialize | ( | const ArrayView< const RCP< const LinearOpBase< Scalar > > > & | Ops | ) |
Initialize given a list of const linear operators.
| Ops | [in] Array (length numOps) of constituent linear operators and their aggregated default definitions of the non-transposed operator. |
Preconditions:
numOps > 0 Ops != NULL Ops[k].op().get()!=NULL, for k=0...numOps-1 Postconditions:
this->numOps()==numOps this->getOp(k).op().get()==Ops[k].op().get(), for k=0...numOps-1 Definition at line 111 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
| void Thyra::DefaultMultipliedLinearOp< Scalar >::uninitialize | ( | ) |
Set to uninitialized.
Postconditions:
this->numOps()==0 Definition at line 124 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Implements Thyra::MultipliedLinearOpBase< Scalar >.
Definition at line 135 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Implements Thyra::MultipliedLinearOpBase< Scalar >.
Definition at line 142 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Implements Thyra::MultipliedLinearOpBase< Scalar >.
Definition at line 153 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Implements Thyra::MultipliedLinearOpBase< Scalar >.
Definition at line 164 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Returns this->getOp(0).range() if <t>this->numOps() > 0 and returns Teuchos::null otherwise.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 178 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Returns this->getOp(this->numOps()-1).domain() if <t>this->numOps() > 0 and returns Teuchos::null otherwise.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 189 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Reimplemented from Thyra::LinearOpBase< Scalar >.
Definition at line 200 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Prints just the name DefaultMultipliedLinearOp along with the overall dimensions and the number of constituent operators.
Reimplemented from Teuchos::Describable.
Definition at line 210 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
virtual |
Prints the details about the constituent linear operators.
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 220 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
protectedvirtual |
Returns true only if all constituent operators support M_trans.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 262 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 273 of file Thyra_DefaultMultipliedLinearOp_def.hpp.
|
related |
Nonmember constructor.
Definition at line 286 of file Thyra_DefaultMultipliedLinearOp_decl.hpp.
|
related |
Nonmember constructor.
Definition at line 298 of file Thyra_DefaultMultipliedLinearOp_decl.hpp.
|
related |
Nonmember constructor.
Definition at line 312 of file Thyra_DefaultMultipliedLinearOp_decl.hpp.
|
related |
Form an implicit multiplication of two linear operators: M = A.
|
related |
Form an implicit multiplication of two linear operators: M = A.
|
related |
Form an implicit multiplication of three linear operators: M = A * B * C.
1.8.5