Thyra
Version of the Day
|
Concrete LinearOpBase
subclass that creates an implicit LinearOpBase
object using the inverse action of a LinearOpWithSolveBase
object.
More...
#include <Thyra_DefaultInverseLinearOp_decl.hpp>
Related Functions | |
(Note that these are not member functions.) | |
template<class Scalar > | |
RCP< LinearOpBase< Scalar > > | nonconstInverse (const RCP< LinearOpWithSolveBase< Scalar > > &A, const Ptr< const SolveCriteria< Scalar > > &fwdSolveCriteria=Teuchos::null, const EThrowOnSolveFailure throwOnFwdSolveFailure=THROW_ON_SOLVE_FAILURE, const Ptr< const SolveCriteria< Scalar > > &adjSolveCriteria=Teuchos::null, const EThrowOnSolveFailure throwOnAdjSolveFailure=THROW_ON_SOLVE_FAILURE) |
Form a non-const implicit inverse operator M = inv(A) . More... | |
template<class Scalar > | |
RCP< LinearOpBase< Scalar > > | inverse (const RCP< const LinearOpWithSolveBase< Scalar > > &A, const Ptr< const SolveCriteria< Scalar > > &fwdSolveCriteria=Teuchos::null, const EThrowOnSolveFailure throwOnFwdSolveFailure=THROW_ON_SOLVE_FAILURE, const Ptr< const SolveCriteria< Scalar > > &adjSolveCriteria=Teuchos::null, const EThrowOnSolveFailure throwOnAdjSolveFailure=THROW_ON_SOLVE_FAILURE) |
Form a const implicit inverse operator M = inv(A) . 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 | |
DefaultInverseLinearOp () | |
Constructs to uninitialized (see postconditions for uninitialize() ). More... | |
DefaultInverseLinearOp (const RCP< LinearOpWithSolveBase< Scalar > > &lows, const SolveCriteria< Scalar > *fwdSolveCriteria=NULL, const EThrowOnSolveFailure throwOnFwdSolveFailure=THROW_ON_SOLVE_FAILURE, const SolveCriteria< Scalar > *adjSolveCriteria=NULL, const EThrowOnSolveFailure throwOnAdjSolveFailure=THROW_ON_SOLVE_FAILURE) | |
DefaultInverseLinearOp (const RCP< const LinearOpWithSolveBase< Scalar > > &lows, const SolveCriteria< Scalar > *fwdSolveCriteria=NULL, const EThrowOnSolveFailure throwOnFwdSolveFailure=THROW_ON_SOLVE_FAILURE, const SolveCriteria< Scalar > *adjSolveCriteria=NULL, const EThrowOnSolveFailure throwOnAdjSolveFailure=THROW_ON_SOLVE_FAILURE) | |
void | initialize (const RCP< LinearOpWithSolveBase< Scalar > > &lows, const SolveCriteria< Scalar > *fwdSolveCriteria=NULL, const EThrowOnSolveFailure throwOnFwdSolveFailure=THROW_ON_SOLVE_FAILURE, const SolveCriteria< Scalar > *adjSolveCriteria=NULL, const EThrowOnSolveFailure throwOnAdjSolveFailure=THROW_ON_SOLVE_FAILURE) |
Initialize given a non-const LinearOpWithSolveBase object and an optional . More... | |
void | initialize (const RCP< const LinearOpWithSolveBase< Scalar > > &lows, const SolveCriteria< Scalar > *fwdSolveCriteria=NULL, const EThrowOnSolveFailure throwOnFwdSolveFailure=THROW_ON_SOLVE_FAILURE, const SolveCriteria< Scalar > *adjSolveCriteria=NULL, const EThrowOnSolveFailure throwOnAdjSolveFailure=THROW_ON_SOLVE_FAILURE) |
Initialize given a non-const LinearOpWithSolveBase object and an optional . More... | |
void | uninitialize () |
Set to uninitialized. More... | |
Overridden from InverseLinearOpBase | |
bool | isLowsConst () const |
RCP< LinearOpWithSolveBase < Scalar > > | getNonconstLows () |
RCP< const LinearOpWithSolveBase< Scalar > > | getLows () const |
Overridden from LinearOpBase | |
RCP< const VectorSpaceBase < Scalar > > | range () const |
Returns this->getLows()->domain() if <t>this->getLows().get()!=NULL and returns Teuchos::null otherwise. More... | |
RCP< const VectorSpaceBase < Scalar > > | domain () const |
Returns this->getLows()->range() if <t>this->getLows().get()!=NULL 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 |
void | describe (FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
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... | |
Protected Member Functions inherited from Thyra::LinearOpBase< Scalar > |
Concrete LinearOpBase
subclass that creates an implicit LinearOpBase
object using the inverse action of a LinearOpWithSolveBase
object.
This class represents an implicit inverse linear operator:
M = inv(A)
where A
is any LinearOpWithSolveBase
object. Specifically, the solve(...)
function A
is used to implement this->apply()
and the solveTranspose(...)
function A
is used to implement this->applyTranspose()
.
SolveCriteria
objects can be associated with A
to define the solve criterion for calling the A.solve(...,fwdSolveCriteria)
and A.solveTranspose(...,adjSolveCriteria)
.
Definition at line 54 of file Thyra_DefaultInverseLinearOp_decl.hpp.
Thyra::DefaultInverseLinearOp< Scalar >::DefaultInverseLinearOp | ( | ) |
Constructs to uninitialized (see postconditions for uninitialize()
).
Definition at line 27 of file Thyra_DefaultInverseLinearOp_def.hpp.
Thyra::DefaultInverseLinearOp< Scalar >::DefaultInverseLinearOp | ( | const RCP< LinearOpWithSolveBase< Scalar > > & | lows, |
const SolveCriteria< Scalar > * | fwdSolveCriteria = NULL , |
||
const EThrowOnSolveFailure | throwOnFwdSolveFailure = THROW_ON_SOLVE_FAILURE , |
||
const SolveCriteria< Scalar > * | adjSolveCriteria = NULL , |
||
const EThrowOnSolveFailure | throwOnAdjSolveFailure = THROW_ON_SOLVE_FAILURE |
||
) |
Calls initialize()
.
Definition at line 32 of file Thyra_DefaultInverseLinearOp_def.hpp.
Thyra::DefaultInverseLinearOp< Scalar >::DefaultInverseLinearOp | ( | const RCP< const LinearOpWithSolveBase< Scalar > > & | lows, |
const SolveCriteria< Scalar > * | fwdSolveCriteria = NULL , |
||
const EThrowOnSolveFailure | throwOnFwdSolveFailure = THROW_ON_SOLVE_FAILURE , |
||
const SolveCriteria< Scalar > * | adjSolveCriteria = NULL , |
||
const EThrowOnSolveFailure | throwOnAdjSolveFailure = THROW_ON_SOLVE_FAILURE |
||
) |
Calls initialize()
.
Rather than calling this constructor directly, consider using the non-member helper functions described here.
Definition at line 48 of file Thyra_DefaultInverseLinearOp_def.hpp.
void Thyra::DefaultInverseLinearOp< Scalar >::initialize | ( | const RCP< LinearOpWithSolveBase< Scalar > > & | lows, |
const SolveCriteria< Scalar > * | fwdSolveCriteria = NULL , |
||
const EThrowOnSolveFailure | throwOnFwdSolveFailure = THROW_ON_SOLVE_FAILURE , |
||
const SolveCriteria< Scalar > * | adjSolveCriteria = NULL , |
||
const EThrowOnSolveFailure | throwOnAdjSolveFailure = THROW_ON_SOLVE_FAILURE |
||
) |
Initialize given a non-const LinearOpWithSolveBase
object and an optional .
lows | [in] The LinearOpWithSolveBase object that will solve(...) and/or solveTranspose(...) will be called on. Note that *this may give up non-const views of *lows so that *lows may be changed through clients of this object. |
fwdSolveCriteria | [in] The criteria used to call lows->solve(...) . If fwdSolveCriteria==NULL then the default solve criteria built into *lows |
adjSolveCriteria |
|
Preconditions:
lows.get() != NULL
Postconditions:
Definition at line 64 of file Thyra_DefaultInverseLinearOp_def.hpp.
void Thyra::DefaultInverseLinearOp< Scalar >::initialize | ( | const RCP< const LinearOpWithSolveBase< Scalar > > & | lows, |
const SolveCriteria< Scalar > * | fwdSolveCriteria = NULL , |
||
const EThrowOnSolveFailure | throwOnFwdSolveFailure = THROW_ON_SOLVE_FAILURE , |
||
const SolveCriteria< Scalar > * | adjSolveCriteria = NULL , |
||
const EThrowOnSolveFailure | throwOnAdjSolveFailure = THROW_ON_SOLVE_FAILURE |
||
) |
Initialize given a non-const LinearOpWithSolveBase
object and an optional .
lows | [in] The LinearOpWithSolveBase object that will solve(...) and/or solveTranspose(...) will be called on. Note that *this may give up non-const views of *lows so that *lows may be changed through clients of this object. |
fwdSolveCriteria | [in] The criteria used to call lows->solve(...) . If fwdSolveCriteria==NULL then the default solve criteria built into *lows |
adjSolveCriteria |
|
Preconditions:
lows.get() != NULL
Postconditions:
Definition at line 80 of file Thyra_DefaultInverseLinearOp_def.hpp.
void Thyra::DefaultInverseLinearOp< Scalar >::uninitialize | ( | ) |
Set to uninitialized.
Postconditions:
Definition at line 96 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
virtual |
Implements Thyra::InverseLinearOpBase< Scalar >.
Definition at line 108 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
virtual |
Implements Thyra::InverseLinearOpBase< Scalar >.
Definition at line 116 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
virtual |
Implements Thyra::InverseLinearOpBase< Scalar >.
Definition at line 124 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
virtual |
Returns this->getLows()->domain() if <t>this->getLows().get()!=NULL
and returns Teuchos::null
otherwise.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 135 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
virtual |
Returns this->getLows()->range() if <t>this->getLows().get()!=NULL
and returns Teuchos::null
otherwise.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 144 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
virtual |
Reimplemented from Thyra::LinearOpBase< Scalar >.
Definition at line 153 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
virtual |
Reimplemented from Teuchos::Describable.
Definition at line 163 of file Thyra_DefaultInverseLinearOp_def.hpp.
void Thyra::DefaultInverseLinearOp< Scalar >::describe | ( | FancyOStream & | out, |
const Teuchos::EVerbosityLevel | verbLevel | ||
) | const |
Definition at line 178 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
protectedvirtual |
Returns true
only if all constituent operators support M_trans
.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 223 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 233 of file Thyra_DefaultInverseLinearOp_def.hpp.
|
related |
Form a non-const implicit inverse operator M = inv(A)
.
|
related |
Form a const implicit inverse operator M = inv(A)
.