Thyra
Version of the Day
|
Default concrete LinearOpBase
subclass for diagonal linear operators.
More...
#include <Thyra_DefaultDiagonalLinearOp_decl.hpp>
Related Functions | |
(Note that these are not member functions.) | |
template<class Scalar > | |
RCP< const LinearOpBase< Scalar > > | diagonal (const RCP< VectorBase< Scalar > > &diag, const std::string &label="") |
Nonmember constructor function. 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 | |
DefaultDiagonalLinearOp () | |
Constructs to uninitialized. More... | |
DefaultDiagonalLinearOp (const RCP< const VectorSpaceBase< Scalar > > &space) | |
Calls initialize() to construct given a vector space. More... | |
DefaultDiagonalLinearOp (const RCP< VectorBase< Scalar > > &diag) | |
Calls initialize() to construct for a non-const diagonal vector. More... | |
DefaultDiagonalLinearOp (const RCP< const VectorBase< Scalar > > &diag) | |
Calls initialize() to construct for a const diagonal vector. More... | |
void | initialize (const RCP< const VectorSpaceBase< Scalar > > &space) |
Initialize given a vector space which allocates a vector internally. More... | |
void | initialize (const RCP< VectorBase< Scalar > > &diag) |
Initialize given a non-const diagonal vector. More... | |
void | initialize (const RCP< const VectorBase< Scalar > > &diag) |
Initialize given a const diagonal vector. More... | |
void | uninitialize () |
Uninitialize. More... | |
Overridden from DiagonalLinearOpBase | |
bool | isDiagConst () const |
RCP< VectorBase< Scalar > > | getNonconstDiag () |
RCP< const VectorBase< Scalar > > | getDiag () const |
Overridden from LinearOpBase | |
RCP< const VectorSpaceBase < Scalar > > | range () const |
Returns this->getDiag()->space() . More... | |
RCP< const VectorSpaceBase < Scalar > > | domain () const |
Returns this->getDiag()->space() . More... | |
RCP< const LinearOpBase< Scalar > > | clone () const |
Protected functions overridden from LinearOpBase | |
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 |
Default concrete LinearOpBase
subclass for diagonal linear operators.
This class represents a diagonal linear operator M
of the form:
M = diag(diag)
where diag
is a VectorBase
object.
The defined operator implements this->apply()
as follows:
y = alpha*op(M)*x + beta*y => y(i) = alpha*diag(i)*x(i) + beta*y(i), for i = 0 ... n-1
where n = this->domain()->dim()
.
Definition at line 50 of file Thyra_DefaultDiagonalLinearOp_decl.hpp.
Thyra::DefaultDiagonalLinearOp< Scalar >::DefaultDiagonalLinearOp | ( | ) |
Constructs to uninitialized.
Postconditions:
this->getDiag().get()==NULL
Definition at line 28 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
Thyra::DefaultDiagonalLinearOp< Scalar >::DefaultDiagonalLinearOp | ( | const RCP< const VectorSpaceBase< Scalar > > & | space | ) |
Calls initialize()
to construct given a vector space.
Definition at line 33 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
Thyra::DefaultDiagonalLinearOp< Scalar >::DefaultDiagonalLinearOp | ( | const RCP< VectorBase< Scalar > > & | diag | ) |
Calls initialize()
to construct for a non-const diagonal vector.
Definition at line 42 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
Thyra::DefaultDiagonalLinearOp< Scalar >::DefaultDiagonalLinearOp | ( | const RCP< const VectorBase< Scalar > > & | diag | ) |
Calls initialize()
to construct for a const diagonal vector.
Definition at line 51 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
void Thyra::DefaultDiagonalLinearOp< Scalar >::initialize | ( | const RCP< const VectorSpaceBase< Scalar > > & | space | ) |
Initialize given a vector space which allocates a vector internally.
space | [in] Smart pointer to vector space |
Preconditions:
space.get()!=NULL
Postconditions:
this->getNonconstDiag()->space()->isCompatible(*space)==true
this->getDiag()->space()->isCompatible(*space)==true
this->this->domain().get() == space.get()
this->this->range().get() == space.get()
Definition at line 60 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
void Thyra::DefaultDiagonalLinearOp< Scalar >::initialize | ( | const RCP< VectorBase< Scalar > > & | diag | ) |
Initialize given a non-const diagonal vector.
diag | [in] Smart pointer to diagonal vector. |
Preconditions:
diag.get()!=NULL
Postconditions:
this->getNonconstDiag().get()==diag.get()
this->getDiag().get()==diag.get()
this->this->domain().get() == diag->space().get()
this->this->range().get() == diag->space().get()
Definition at line 72 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
void Thyra::DefaultDiagonalLinearOp< Scalar >::initialize | ( | const RCP< const VectorBase< Scalar > > & | diag | ) |
Initialize given a const diagonal vector.
diag | [in] Smart pointer to diagonal vector. |
Preconditions:
diag.get()!=NULL
Postconditions:
this->getNonconstDiag()
with throw an exception if called this->getDiag().get()==diag.get()
this->this->domain().get() == diag->space().get()
this->this->range().get() == diag->space().get()
Definition at line 81 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
void Thyra::DefaultDiagonalLinearOp< Scalar >::uninitialize | ( | ) |
Uninitialize.
Postconditions:
this->getNonconstDiag().get()==NULL
this->getDiag().get()==NULL
this->this->domain().get()==NULL
this->this->range().get()==NULL
Note: If the client wants to hold on to the underlying wrapped diagonal vector then they had better grab it using this->getDiag()
or this->getNonconstDiag()
before calling this function!
Definition at line 90 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
virtual |
Implements Thyra::DiagonalLinearOpBase< Scalar >.
Definition at line 100 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
virtual |
Implements Thyra::DiagonalLinearOpBase< Scalar >.
Definition at line 108 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
virtual |
Implements Thyra::DiagonalLinearOpBase< Scalar >.
Definition at line 116 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
virtual |
Returns this->getDiag()->space()
.
Preconditions:
this->getDiag().get()!=NULL
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 127 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
virtual |
Returns this->getDiag()->space()
.
Preconditions:
this->getDiag().get()!=NULL
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 135 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
virtual |
Reimplemented from Thyra::LinearOpBase< Scalar >.
Definition at line 143 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 156 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 163 of file Thyra_DefaultDiagonalLinearOp_def.hpp.
|
related |
Nonmember constructor function.
Definition at line 219 of file Thyra_DefaultDiagonalLinearOp_decl.hpp.