Thyra
Version of the Day
|
Concrete composite LinearOpBase
subclass that creates single linear operator object out of a set of constituent LinearOpBase
blocks.
More...
#include <Thyra_DefaultBlockedLinearOp_decl.hpp>
Related Functions | |
(Note that these are not member functions.) | |
template<class Scalar > | |
RCP< DefaultBlockedLinearOp < Scalar > > | defaultBlockedLinearOp () |
Nonmember default constructor. More... | |
template<class Scalar > | |
Teuchos::RCP< const LinearOpBase< Scalar > > | block1x1 (const Teuchos::RCP< const LinearOpBase< Scalar > > &A00, const std::string &label="") |
Form an implicit block 1x1 linear operator [ A00 ] . More... | |
template<class Scalar > | |
Teuchos::RCP< const LinearOpBase< Scalar > > | block1x2 (const Teuchos::RCP< const LinearOpBase< Scalar > > &A00, const Teuchos::RCP< const LinearOpBase< Scalar > > &A01, const std::string &label="") |
Form an implicit block 1x2 linear operator [ A00, A01 ] . More... | |
template<class Scalar > | |
Teuchos::RCP< const LinearOpBase< Scalar > > | block2x1 (const Teuchos::RCP< const LinearOpBase< Scalar > > &A00, const Teuchos::RCP< const LinearOpBase< Scalar > > &A10, const std::string &label="") |
Form an implicit block 2x1 linear operator [ A00; A10 ] . More... | |
template<class Scalar > | |
Teuchos::RCP< const LinearOpBase< Scalar > > | block2x2 (const Teuchos::RCP< const LinearOpBase< Scalar > > &A00, const Teuchos::RCP< const LinearOpBase< Scalar > > &A01, const Teuchos::RCP< const LinearOpBase< Scalar > > &A10, const Teuchos::RCP< const LinearOpBase< Scalar > > &A11, const std::string &label="") |
Form an implicit block 2x2 linear operator [ A00, A01; A10, A11 ] . More... | |
template<class Scalar > | |
Teuchos::RCP< LinearOpBase < Scalar > > | nonconstBlock1x1 (const Teuchos::RCP< LinearOpBase< Scalar > > &A00, const std::string &label="") |
Form an implicit block 1x1 linear operator [ A00 ] . More... | |
template<class Scalar > | |
Teuchos::RCP< LinearOpBase < Scalar > > | nonconstBlock1x2 (const Teuchos::RCP< LinearOpBase< Scalar > > &A00, const Teuchos::RCP< LinearOpBase< Scalar > > &A01, const std::string &label="") |
Form an implicit block 1x2 linear operator [ A00, A01 ] . More... | |
template<class Scalar > | |
Teuchos::RCP< LinearOpBase < Scalar > > | nonconstBlock2x1 (const Teuchos::RCP< LinearOpBase< Scalar > > &A00, const Teuchos::RCP< LinearOpBase< Scalar > > &A10, const std::string &label="") |
Form an implicit block 2x1 linear operator [ A00; A10 ] . More... | |
template<class Scalar > | |
Teuchos::RCP< LinearOpBase < Scalar > > | nonconstBlock2x2 (const Teuchos::RCP< LinearOpBase< Scalar > > &A00, const Teuchos::RCP< LinearOpBase< Scalar > > &A01, const Teuchos::RCP< LinearOpBase< Scalar > > &A10, const Teuchos::RCP< LinearOpBase< Scalar > > &A11, const std::string &label="") |
Form an implicit block 2x2 linear operator [ A00, A01; A10, A11 ] . 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 | |
DefaultBlockedLinearOp () | |
Overridden from PhysicallyBlockedLinearOpBase | |
void | beginBlockFill () |
void | beginBlockFill (const int numRowBlocks, const int numColBlocks) |
void | beginBlockFill (const Teuchos::RCP< const ProductVectorSpaceBase< Scalar > > &productRange, const Teuchos::RCP< const ProductVectorSpaceBase< Scalar > > &productDomain) |
bool | blockFillIsActive () const |
bool | acceptsBlock (const int i, const int j) const |
void | setNonconstBlock (const int i, const int j, const Teuchos::RCP< LinearOpBase< Scalar > > &block) |
void | setBlock (const int i, const int j, const Teuchos::RCP< const LinearOpBase< Scalar > > &block) |
void | endBlockFill () |
void | uninitialize () |
Overridden from BlockedLinearOpBase | |
Teuchos::RCP< const ProductVectorSpaceBase< Scalar > > | productRange () const |
Teuchos::RCP< const ProductVectorSpaceBase< Scalar > > | productDomain () const |
bool | blockExists (const int i, const int j) const |
bool | blockIsConst (const int i, const int j) const |
Teuchos::RCP< LinearOpBase < Scalar > > | getNonconstBlock (const int i, const int j) |
Teuchos::RCP< const LinearOpBase< Scalar > > | getBlock (const int i, const int j) const |
Overridden from LinearOpBase | |
Teuchos::RCP< const VectorSpaceBase< Scalar > > | range () const |
Teuchos::RCP< const VectorSpaceBase< Scalar > > | domain () const |
Teuchos::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 DefaultBlockedLinearOp 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... | |
Overridden from RowStatLinearOpBase | |
virtual bool | rowStatIsSupportedImpl (const RowStatLinearOpBaseUtils::ERowStat rowStat) const |
virtual void | getRowStatImpl (const RowStatLinearOpBaseUtils::ERowStat rowStat, const Teuchos::Ptr< VectorBase< Scalar > > &rowStatVec) const |
Overridden from ScaledLinearOpBase | |
virtual bool | supportsScaleLeftImpl () const |
virtual bool | supportsScaleRightImpl () const |
virtual void | scaleLeftImpl (const VectorBase< Scalar > &row_scaling) |
virtual void | scaleRightImpl (const VectorBase< Scalar > &col_scaling) |
Additional Inherited Members | |
Public Member Functions inherited from Thyra::LinearOpBase< Scalar > | |
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... | |
Public Member Functions inherited from Thyra::RowStatLinearOpBase< Scalar > | |
bool | rowStatIsSupported (const RowStatLinearOpBaseUtils::ERowStat rowStat) const |
Determine if a given row stat is supported. More... | |
void | getRowStat (const RowStatLinearOpBaseUtils::ERowStat rowStat, const Ptr< VectorBase< Scalar > > &rowStatVec) const |
Get some statistics about a supported row. More... | |
Public Member Functions inherited from Thyra::ScaledLinearOpBase< Scalar > | |
bool | supportsScaleLeft () const |
Determines if this objects supports left scaling. More... | |
bool | supportsScaleRight () const |
Determines if this objects supports right scaling. More... | |
void | scaleLeft (const VectorBase< Scalar > &row_scaling) |
Left scales operator with diagonal scaling operator. More... | |
void | scaleRight (const VectorBase< Scalar > &col_scaling) |
Right scales operator with diagonal scaling operator. More... | |
Protected Member Functions inherited from Thyra::LinearOpBase< Scalar > | |
Protected Member Functions inherited from Thyra::RowStatLinearOpBase< Scalar > | |
Protected Member Functions inherited from Thyra::ScaledLinearOpBase< Scalar > |
Concrete composite LinearOpBase
subclass that creates single linear operator object out of a set of constituent LinearOpBase
blocks.
This class represents a blocked linear operator M
of the form:
M = [ Op[0,0], Op[0,1], ... , Op[0,N]; Op[1,0], Op[1,1], ... , Op[1,N]; . . . Op[M,0], Op[M,1], ... , Op[M,N]; ]
where Op[]
is a logical 2D array of LinearOpBase
objects and M=this->productRange()->getNumBlocks()
and N=this->productDomain()->getNumBlocks()
. Of course the operator M
is not constructed explicitly but instead just applies the constituent linear operators with each set of blocks.
ToDo: Finish Documentation!
Definition at line 54 of file Thyra_DefaultBlockedLinearOp_decl.hpp.
Thyra::DefaultBlockedLinearOp< Scalar >::DefaultBlockedLinearOp | ( | ) |
Definition at line 31 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 40 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 49 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 60 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 80 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 87 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 98 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 108 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 118 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 194 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 213 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 221 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 228 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 239 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 254 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 267 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 283 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 291 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Reimplemented from Thyra::LinearOpBase< Scalar >.
Definition at line 299 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
virtual |
Prints just the name DefaultBlockedLinearOp
along with the overall dimensions and the number of constituent operators.
Reimplemented from Teuchos::Describable.
Definition at line 309 of file Thyra_DefaultBlockedLinearOp_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 323 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
protectedvirtual |
Returns true
only if all constituent operators support M_trans
.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 381 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 397 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::RowStatLinearOpBase< Scalar >.
Definition at line 472 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::RowStatLinearOpBase< Scalar >.
Definition at line 537 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::ScaledLinearOpBase< Scalar >.
Definition at line 627 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::ScaledLinearOpBase< Scalar >.
Definition at line 663 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::ScaledLinearOpBase< Scalar >.
Definition at line 698 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::ScaledLinearOpBase< Scalar >.
Definition at line 751 of file Thyra_DefaultBlockedLinearOp_def.hpp.
|
related |
Nonmember default constructor.
|
related |
Form an implicit block 1x1 linear operator [ A00 ]
.
|
related |
Form an implicit block 1x2 linear operator [ A00, A01 ]
.
|
related |
Form an implicit block 2x1 linear operator [ A00; A10 ]
.
|
related |
Form an implicit block 2x2 linear operator [ A00, A01; A10, A11 ]
.
|
related |
Form an implicit block 1x1 linear operator [ A00 ]
.
|
related |
Form an implicit block 1x2 linear operator [ A00, A01 ]
.
|
related |
Form an implicit block 2x1 linear operator [ A00; A10 ]
.
|
related |
Form an implicit block 2x2 linear operator [ A00, A01; A10, A11 ]
.