Thyra
Version of the Day
|
Concrete composite LinearOpWithSolveBase
subclass that creates single upper or lower block triangular LOWSB object out of a set of LOWSB objects along the diagonal with LOB objects off diagonal.
More...
#include <Thyra_DefaultBlockedTriangularLinearOpWithSolve_decl.hpp>
Related Functions | |
(Note that these are not member functions.) | |
template<class Scalar > | |
RCP < DefaultBlockedTriangularLinearOpWithSolve < Scalar > > | defaultBlockedTriangularLinearOpWithSolve () |
Nonmember constructor. More... | |
Related Functions inherited from Thyra::LinearOpWithSolveBase< Scalar > | |
template<class Scalar > | |
bool | solveSupports (const LinearOpWithSolveBase< Scalar > &A, const EOpTransp transp) |
Call solveSupports() as a non-member function. More... | |
template<class Scalar > | |
bool | solveSupports (const LinearOpWithSolveBase< Scalar > &A, const EOpTransp transp, const Ptr< const SolveCriteria< Scalar > > solveCriteria) |
Call solveSupports() as a non-member function. More... | |
template<class Scalar > | |
bool | solveSupportsSolveMeasureType (const LinearOpWithSolveBase< Scalar > &A, const EOpTransp transp, const SolveMeasureType &solveMeasureType) |
Call solveSupportsSolveMeasureType() as a non-member function. More... | |
template<class Scalar > | |
SolveStatus< Scalar > | solve (const LinearOpWithSolveBase< Scalar > &A, const EOpTransp A_trans, const MultiVectorBase< Scalar > &B, const Ptr< MultiVectorBase< Scalar > > &X, const Ptr< const SolveCriteria< Scalar > > solveCriteria=Teuchos::null) |
Call solve() as a non-member function. More... | |
template<class Scalar > | |
void | assertSolveSupports (const LinearOpWithSolveBase< Scalar > &lows, const EOpTransp M_trans, const Ptr< const SolveCriteria< Scalar > > solveCriteria=Teuchos::null) |
Assert that a LOWSB object supports a particular solve type. 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 | |
DefaultBlockedTriangularLinearOpWithSolve () | |
void | setNonconstBlocks (const RCP< PhysicallyBlockedLinearOpBase< Scalar > > &blocks) |
void | setBlocks (const RCP< const PhysicallyBlockedLinearOpBase< Scalar > > &blocks) |
RCP < PhysicallyBlockedLinearOpBase < Scalar > > | getNonconstBlocks () |
RCP< const PhysicallyBlockedLinearOpBase < Scalar > > | getBlocks () |
Overridden from PhysicallyBlockedLinearOpWithSolveBase | |
bool | acceptsLOWSBlock (const int i, const int j) const |
void | setNonconstLOWSBlock (const int i, const int j, const RCP< LinearOpWithSolveBase< Scalar > > &block) |
void | setLOWSBlock (const int i, const int j, const RCP< const LinearOpWithSolveBase< Scalar > > &block) |
Overridden from PhysicallyBlockedLinearOpBase | |
void | beginBlockFill () |
void | beginBlockFill (const int numRowBlocks, const int numColBlocks) |
void | beginBlockFill (const RCP< const ProductVectorSpaceBase< Scalar > > &productRange, const 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 RCP< LinearOpBase< Scalar > > &block) |
void | setBlock (const int i, const int j, const RCP< const LinearOpBase< Scalar > > &block) |
void | endBlockFill () |
void | uninitialize () |
Overridden from BlockedLinearOpWithSolveBase | |
RCP< LinearOpWithSolveBase < Scalar > > | getNonconstLOWSBlock (const int i, const int j) |
RCP< const LinearOpWithSolveBase< Scalar > > | getLOWSBlock (const int i, const int j) const |
Overridden from BlockedLinearOpBase | |
RCP< const ProductVectorSpaceBase< Scalar > > | productRange () const |
RCP< const ProductVectorSpaceBase< Scalar > > | productDomain () const |
bool | blockExists (const int i, const int j) const |
bool | blockIsConst (const int i, const int j) const |
RCP< LinearOpBase< Scalar > > | getNonconstBlock (const int i, const int j) |
RCP< const LinearOpBase< Scalar > > | getBlock (const int i, const int j) const |
Overridden from LinearOpBase | |
RCP< const VectorSpaceBase < Scalar > > | range () const |
RCP< const VectorSpaceBase < Scalar > > | domain () const |
RCP< const LinearOpBase< Scalar > > | clone () const |
bool | opSupportedImpl (EOpTransp M_trans) const |
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 DefaultBlockedTriangularLinearOpWithSolve 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 LinearOpWithSolveBase | |
bool | solveSupportsImpl (EOpTransp M_trans) const |
bool | solveSupportsSolveMeasureTypeImpl (EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const |
SolveStatus< Scalar > | solveImpl (const EOpTransp transp, const MultiVectorBase< Scalar > &B, const Ptr< MultiVectorBase< Scalar > > &X, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const |
Additional Inherited Members | |
Public Member Functions inherited from Thyra::LinearOpWithSolveBase< Scalar > | |
bool | solveSupports (EOpTransp transp) const |
Return if solve() supports the argument transp . More... | |
bool | solveSupports (EOpTransp transp, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const |
Return if solve() supports a given transpose and solve criteria specification. More... | |
bool | solveSupportsSolveMeasureType (EOpTransp transp, const SolveMeasureType &solveMeasureType) const |
Return if solve() supports the given the solve measure type. More... | |
SolveStatus< Scalar > | solve (const EOpTransp A_trans, const MultiVectorBase< Scalar > &B, const Ptr< MultiVectorBase< Scalar > > &X, const Ptr< const SolveCriteria< Scalar > > solveCriteria=Teuchos::null) const |
Request the solution of a block linear system. More... | |
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... | |
Protected Member Functions inherited from Thyra::LinearOpWithSolveBase< Scalar > | |
virtual bool | solveSupportsNewImpl (EOpTransp, const Ptr< const SolveCriteria< Scalar > >) const |
Virtual implementation of solveSupports() . More... | |
Protected Member Functions inherited from Thyra::LinearOpBase< Scalar > |
Concrete composite LinearOpWithSolveBase
subclass that creates single upper or lower block triangular LOWSB object out of a set of LOWSB objects along the diagonal with LOB objects off diagonal.
This subclass implements a strictly block upper or lower triangular LOWSB object. With LOWSB objects on the diagonal, the block system can be solved by solving for each of the diagonal blocks and moving lower (or upper) blocks to the RHS on each successive solve.
For example, consider the lower block triangular linear operator:
[ M(0,0) ] M = [ M(1,0) M(1,1) ] [ M(2,0) M(2,1) M(2,2) ]
A linear system of the form:
M * x = b => [ M(0,0) ] [ x(0) ] [ b(0) ] [ M(1,0) M(1,1) ] * [ x(1) ] = [ b(1) ] [ M(2,0) M(2,1) M(2,2) ] [ x(2) ] [ b(2) ]
is solved as:
x(0) = inv(M(0,0)) * b(0) x(1) = inv(M(1,1)) * ( b(1) - M(1,0)*x(0) ) x(2) = inv(M(2,2)) * ( b(2) - M(2,0)*x(0) - M(2,1)*x(1) )
The same approach can be used for block upper triangular linear operators as well of course.
See the class DefaultBlockedTriangularLinearOpWithSolveFactory
for an example of how one of these objects can be created from any PhysicallyBlockedLinearOpBase
object and compatible LinearWithSolveBase
objects.
ToDo: Finish Documentation!
Definition at line 112 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_decl.hpp.
Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >::DefaultBlockedTriangularLinearOpWithSolve | ( | ) |
Definition at line 62 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
void Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >::setNonconstBlocks | ( | const RCP< PhysicallyBlockedLinearOpBase< Scalar > > & | blocks | ) |
Definition at line 68 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
void Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >::setBlocks | ( | const RCP< const PhysicallyBlockedLinearOpBase< Scalar > > & | blocks | ) |
Definition at line 78 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
RCP< PhysicallyBlockedLinearOpBase< Scalar > > Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >::getNonconstBlocks | ( | ) |
Definition at line 89 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
RCP< const PhysicallyBlockedLinearOpBase< Scalar > > Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >::getBlocks | ( | ) |
Definition at line 97 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpWithSolveBase< Scalar >.
Definition at line 107 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpWithSolveBase< Scalar >.
Definition at line 117 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpWithSolveBase< Scalar >.
Definition at line 127 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 140 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 148 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 168 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 188 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 195 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 206 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 217 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 228 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::PhysicallyBlockedLinearOpBase< Scalar >.
Definition at line 253 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpWithSolveBase< Scalar >.
Definition at line 268 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpWithSolveBase< Scalar >.
Definition at line 282 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 299 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 307 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 314 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 327 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 339 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::BlockedLinearOpBase< Scalar >.
Definition at line 353 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 370 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 378 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Reimplemented from Thyra::LinearOpBase< Scalar >.
Definition at line 386 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
virtual |
Prints just the name DefaultBlockedTriangularLinearOpWithSolve
along with the overall dimensions and the number of constituent operators.
Reimplemented from Teuchos::Describable.
Definition at line 397 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_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 410 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 428 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 446 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::LinearOpWithSolveBase< Scalar >.
Definition at line 496 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::LinearOpWithSolveBase< Scalar >.
Definition at line 514 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpWithSolveBase< Scalar >.
Definition at line 537 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_def.hpp.
|
related |
Nonmember constructor.
Definition at line 346 of file Thyra_DefaultBlockedTriangularLinearOpWithSolve_decl.hpp.