|
Thyra Version of the Day
|
Simple concreate subclass of LinearOpWithSolveBase for serial dense matrices implemented using LAPACK.
More...
#include <Thyra_DefaultSerialDenseLinearOpWithSolve_decl.hpp>

Related Symbols | |
(Note that these are not member symbols.) | |
| template<class Scalar > | |
| RCP< DefaultSerialDenseLinearOpWithSolve< Scalar > > | defaultSerialDenseLinearOpWithSolve () |
| Nonmember constructor. | |
| template<class Scalar > | |
| RCP< DefaultSerialDenseLinearOpWithSolve< Scalar > > | defaultSerialDenseLinearOpWithSolve (const RCP< const MultiVectorBase< Scalar > > &M) |
| Nonmember constructor. | |
Related Symbols inherited from Thyra::LinearOpWithSolveBase< Scalar > | |
| template<class Scalar > | |
| bool | solveSupports (const LinearOpWithSolveBase< Scalar > &A, const EOpTransp transp) |
Call solveSupports() as a non-member function. | |
| 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. | |
| template<class Scalar > | |
| bool | solveSupportsSolveMeasureType (const LinearOpWithSolveBase< Scalar > &A, const EOpTransp transp, const SolveMeasureType &solveMeasureType) |
Call solveSupportsSolveMeasureType() as a non-member function. | |
| 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. | |
| 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. | |
| template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node > | |
| SolveStatus< Scalar > | solve (const LinearOpWithSolveBase< Scalar > &A, const EOpTransp A_trans, const RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraB, const RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraX, const Ptr< const SolveCriteria< Scalar > > solveCriteria=Teuchos::null) |
Call solve() as a non-member function. | |
Related Symbols 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. | |
| template<class Scalar > | |
| bool | isPartiallyInitialized (const LinearOpBase< Scalar > &M) |
| Determines if a linear operator is in the "Partially Initialized" state or not. | |
| template<class Scalar > | |
| bool | isFullyInitialized (const LinearOpBase< Scalar > &M) |
| 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) |
| 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)) |
Non-member function call for M.apply(...). | |
| 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>(...). | |
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 |
Constructors/initializers/accessors | |
| DefaultSerialDenseLinearOpWithSolve () | |
| void | initialize (const RCP< const MultiVectorBase< Scalar > > &M) |
| RCP< const LinearOpBase< Scalar > > | getFwdOp () const |
Overridden from LinearOpBase | |
| RCP< const VectorSpaceBase< Scalar > > | range () const |
| RCP< const VectorSpaceBase< Scalar > > | domain () 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 |
Additional Inherited Members | |
Public Member Functions inherited from Thyra::LinearOpWithSolveBase< Scalar > | |
| bool | solveSupports (EOpTransp transp) const |
Return if solve() supports the argument transp. | |
| bool | solveSupports (EOpTransp transp, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const |
Return if solve() supports a given transpose and solve criteria specification. | |
| bool | solveSupportsSolveMeasureType (EOpTransp transp, const SolveMeasureType &solveMeasureType) const |
Return if solve() supports the given the solve measure type. | |
| 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. | |
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. | |
| 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. | |
| virtual RCP< const LinearOpBase< Scalar > > | clone () const |
| Clone the linear operator object (if supported). | |
Protected Member Functions inherited from Thyra::LinearOpWithSolveBase< Scalar > | |
| virtual bool | solveSupportsNewImpl (EOpTransp, const Ptr< const SolveCriteria< Scalar > >) const |
Virtual implementation of solveSupports(). | |
Protected Member Functions inherited from Thyra::LinearOpBase< Scalar > | |
Simple concreate subclass of LinearOpWithSolveBase for serial dense matrices implemented using LAPACK.
This class uses the helper class DetachedMultiVectorView to extract an explicit view of the matrix elements and then uses Teuchos::LAPACK to factor M = L * U and then do back-solves with the factors L and U.
Even through this class accesses explicit matrix entries and is called SerialDense, it is still considered an ANA subclass since it does not have any direct dependance on a specific computing environment or concreate operator/vector/vectorspace implementation.
ToDo: Finish Documentation!
Definition at line 58 of file Thyra_DefaultSerialDenseLinearOpWithSolve_decl.hpp.
| Thyra::DefaultSerialDenseLinearOpWithSolve< Scalar >::DefaultSerialDenseLinearOpWithSolve | ( | ) |
Definition at line 29 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
| void Thyra::DefaultSerialDenseLinearOpWithSolve< Scalar >::initialize | ( | const RCP< const MultiVectorBase< Scalar > > & | M | ) |
Definition at line 34 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
| RCP< const LinearOpBase< Scalar > > Thyra::DefaultSerialDenseLinearOpWithSolve< Scalar >::getFwdOp | ( | ) | const |
Definition at line 50 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 60 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 70 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 85 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 93 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::LinearOpWithSolveBase< Scalar >.
Definition at line 109 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::LinearOpWithSolveBase< Scalar >.
Definition at line 118 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpWithSolveBase< Scalar >.
Definition at line 128 of file Thyra_DefaultSerialDenseLinearOpWithSolve_def.hpp.
|
related |
Nonmember constructor.
Definition at line 162 of file Thyra_DefaultSerialDenseLinearOpWithSolve_decl.hpp.
|
related |
Nonmember constructor.
Definition at line 174 of file Thyra_DefaultSerialDenseLinearOpWithSolve_decl.hpp.