Thyra
Version of the Day
|
Base class for all linear operators that can support a high-level solve operation. More...
#include <Thyra_LinearOpWithSolveBase_decl.hpp>
Related Functions | |
(Note that these are not member functions.) | |
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... | |
Public interface funtions. | |
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... | |
Protected virtual functions to be overridden by subclasses. | |
virtual bool | solveSupportsImpl (EOpTransp transp) const |
Virtual implementation for solveSupports(). More... | |
virtual bool | solveSupportsNewImpl (EOpTransp, const Ptr< const SolveCriteria< Scalar > >) const |
Virtual implementation of solveSupports() . More... | |
virtual bool | solveSupportsSolveMeasureTypeImpl (EOpTransp transp, const SolveMeasureType &solveMeasureType) const |
Virtual implementation for solveSupportsSolveMeasureType(). More... | |
virtual SolveStatus< Scalar > | solveImpl (const EOpTransp transp, const MultiVectorBase< Scalar > &B, const Ptr< MultiVectorBase< Scalar > > &X, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const =0 |
Virtual implementation for solve(). More... | |
Base class for all linear operators that can support a high-level solve operation.
This interface supports linear operators (with potentially different range and domain scalar types) that can also support a forward solve operation (using solve()
) of the form:
and/or a transpose solve operation (using A_trans==TRANS
) of the form:
and/or an adjoint solve operation (using A_trans==CONJTRANS()
) of the form:
where is *this
linear operator, is an appropriate RHS multi-vector with columns, and is a LHS multi-vector with columns that is computed by this interface. Note that if the underlying operator has real-valued entries then the transpose and the adjoint are the same.
Note that this interface does not assume that the linear operator itself is nonsingular or invertible in the classic sense (i.e. an inverse operator may not exist).
Let signify either the forward operator , the transpose operator or the adjoint operator . What this interface assumes is that for any appropriately selected consistent multi-vector RHS B
that a solve of will yield an approximate solution LHS multi-vector X
such that A X == B
. Note that this interface does not assume that a solution can be computed for any random RHS multi-vector . Solutions for any random RHS can on be expected for relatively well conditioned non-singular operators.
Note: It is recommended that clients use the non-member helper functions defined here rather than call these member functions directly as they support a number of other simpler use cases.
This interface potentially allows clients to specify a relative tolerance on either the relative residual norm or the relative norm of the solution error and can target different solution criteria to different blocks of linear systems. This interface tries to allow for mathematically rigorous solution tolerances that are not based only any implementation-dependent features like the number of iterations of some solver algorithm. This interface, however, allows *this
operator to exclude support certain types of solve measures (see the functions solveSupportsSolveMeasureType()
and solveTransposeSupportsSolveMeasureType()
). Also, this interface assumes that all implementations can support a "default" solve criteria that is determined internally to *this
.
This interface is meant to support direct and iterative linear solvers as well as combinations of the two in a variety of configurations. Because of the almost infinite number of types of linear solver configurations possible, this interface tries not to specify any particular solver-specific solution control options. The one exception is a maximum number of iterations which is totally implementation defined. These types of control options are better specified in lower lever implementations and should be kept out of an interface such as this.
Let solveCriteria
be a SolveCriteria
object setup by a client to be passed into a solve operation. This object can be set up in a variety of ways to support several different use cases which are described below:
Unspecified (Default) solution criteria type [ solveCriteria.solveMeasureType.useDefault()==true
]: In this mode, criteria for the solution tolerance is determined internally by *this
object. Usually, it would be assumed that *this
would solve the linear systems to a sufficient tolerance for the given ANA client. This is the mode that many ANAs not designed to control inexactness would work with. In this mode, the solution criteria will be determined in the end by the application or the user in an "appropriate" manner. In this case, the value of solveCriteria.requestedTol
is ignored and no meaningful solve status can be returned to the client.
Relative residual tolerance solution criteria [ solveCriteria.solveMeasureType.numerator==SOLVE_MEASURE_NORM_RESIDUAL && solveCriteria.solveMeasureType.denominator==SOLVE_MEASURE_NORM_RHS && solveCriteria.requestedTol!=SolveCriteria::unspecifiedTolerance()
]: In this mode, the solution algorithm would be requested to solve the block system to the relative residual tolerance of
for , where solveCriteria.requestedTol
and where the norm is given by the natural norm defined by the range space of and computed from Thyra::norm()
. Many linear solvers should be able to monitor this tolerance and be able to achieve it, or on failure be able to report the actual tolerance achieved.
ToDo: Add examples of other types of solve measures when they are needed!
After the solve()
and solveTranspose()
functions return, the client can optionally get back a solution status for each block of linear systems for of block solve criteria. Specifically, for each block of linear systems
whose solution criteria is specified by a SolveCriteria
object, a SolveStatus
object can optionally be returned that lets the client know the status of the linear solve.
A note about direct solvers is in order. The "inexact" solve features of this interface are primarily designed to support "loose" solve tolerances that exploit the properties of iterative linear solvers. With that said, any decent direct solver can assume that it has met the convergence criteria as requested by the client but does not have to return an estimate of the actual tolerance achieved.
If solveStatus
is a SolveStatus
object returned for the above block linear system the the following return status are significant:
Converged [ solveStatus.solveStatus==SOLVE_STATUS_CONVERGED
]: This status is returned by the linear solver if the solution criteria was likely achieved (within some acceptable cushion cased by round-off etc.). This should almost always be the return value for a direct linear solver. The maximum actual tolerance achieved may or may not be returned in the field solveStatus.achievedTol
. The two sub-cases are:
Known tolerance [ solveStatus.achievedTol!=SolveStatus::unknownTolerance()
] : The linear solver knows the approximate order-of-magnitude estimate of the maximum tolerance achieved. An order-of-magnitude (or so) estimate of the achieved tolerance would likely be known by any iterative linear solver when solveCriteria.solveMeasureType.numerator==SOLVE_MEASURE_NORM_RESIDUAL && solveCriteria.solveMeasureType.denominator==SOLVE_MEASURE_NORM_RHS
. Most direct linear solvers will not return a known tolerance except for the case where solveCriteria.solveMeasureType.numerator==SOLVE_MEASURE_NORM_RESIDUAL && solveCriteria.solveMeasureType.denominator==SOLVE_MEASURE_NORM_RHS
and where iterative refinement is used.
Unknown tolerance [ solveStatus.achievedTol==SolveStatus::unknownTolerance()
] : The linear solver does not know the tolerance that was achieved but the achieved tolerance should be very close to the requested tolerance. This would be the most likely return status for a direct linear solver.
Unconverged [ solveStatus.solveStatus==SOLVE_STATUS_UNCONVERGED
]: The linear solver was most likely not able to achieve the requested tolerance. A direct linear solver should almost never return this status except for in extreme cases (e.g. highly ill conditioned matrix and a tight requested tolerance). The linear solver may not be able to return the actual tolerance achieved and the same two cases as for the <it>unconverged</it> case are possible and the two subdcases are:
Known tolerance [ solveStatus.achievedTol!===SolveStatus::unknownTolerance()0
] : The linear solver knows the approximate order-of-magnitude estimate of the maximum tolerance achieved.
Unknown tolerance [ solveStatus.achievedTol==SolveStatus::unknownTolerance()
] : The linear solver does not know the tolerance that was achieved but the achieved tolerance is most likely significantly greater than the requested tolerance.
Unknown status [ solveStatus.solveStatus==SOLVE_STATUS_UNKNOWN
]: The linear solver does not know if the solution status was achieved or not. In this case, the value of solveStatus.achievedTol==SolveStatus::unknownTolerance()
may or many not be returned. This may be the return value when there is no reasonable way that the linear solver algorithm can know now to compute or estimate the requested tolerance. This will also always be the return status when solveCriteria.solveMeasureType.useDefault()==true
since the client would have no way to interpret this tolerance. The value of solveStatus.achievedTol!=SolveStatus::unknownTolerance()
in this case should only be returned when solveCriteria.solveMeasureType.useDefault()==true
and therefore the client would have no way to interpret this tolerance as a residual or an solution norm.
The implementation of the function accumulateSolveStatus()
defines how to accumulate the individual solve status for each RHS in a block into the overall solve status for an entire block returned in the SolveStatus
object from teh solve()
.
This interface supports a variety of use cases where where described, more or less, in the above sections. Here, we give specific examples for a number of important use cases and show how to use the non-member helper functions defined here.
ToDo: Finish documentation!
This interface assumes, by default, that subclasses will only support the forward solve operation. See LinearOpBase
for what other virtual functions must be overridden to completely define a concrete subclass.
Definition at line 276 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
inline |
Return if solve()
supports the argument transp
.
The default implementation returns true
for non-transposed, non-conjugate solves..
Definition at line 348 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
inline |
Return if solve()
supports a given transpose and solve criteria specification.
Definition at line 355 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
inline |
Return if solve()
supports the given the solve measure type.
The default implementation returns true
for solveMeasureType.useDefault()==true
.
Definition at line 365 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
inline |
Request the solution of a block linear system.
A_trans | [in] Determines if the elements are non-conjugate non-transposed (NONTRANS ) or conjugate transposed (CONJTRANS ). |
B | [in] The RHS multi-vector with m = B.domain()->dim() columns. |
X | [in/out] The LHS multi-vector with with m = X->domain()->dim() columns. On input, contains the initial guess for the solution (only significant for iterative solvers) and on output contains an estimate of the solution. |
solveCriteria | [in] Gives the desired solution criteria for linear systems. A value of solveCriteria==null means to use the default solve criteria. |
Preconditions:
this->solveSupports(transp, solveCriteria)==true
nonnull(X)==true
Postconditions:
If any progress on solving the linear systems could be achieved, then the function will return normally with the solve status. However, if no progress could be made in solving the linear systems, then an exception of type CatastrophicSolveFailure
will be thrown. If the function returns normally, then the following postconditions apply.
See the above introduction for a more complete description of how this function behaves and the meaning of the argument solveCriteria
and the return value.
Definition at line 418 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
protectedvirtual |
Virtual implementation for solveSupports().
Reimplemented in Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >, Thyra::DelayedLinearOpWithSolve< Scalar >, Thyra::DefaultSerialDenseLinearOpWithSolve< Scalar >, Thyra::DefaultMultiVectorLinearOpWithSolve< Scalar >, Thyra::DefaultDiagonalLinearOpWithSolve< Scalar >, and Thyra::DefaultAdjointLinearOpWithSolve< Scalar >.
Definition at line 24 of file Thyra_LinearOpWithSolveBase_def.hpp.
|
inlineprotectedvirtual |
Virtual implementation of solveSupports()
.
Definition at line 437 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
protectedvirtual |
Virtual implementation for solveSupportsSolveMeasureType().
Reimplemented in Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >, Thyra::DelayedLinearOpWithSolve< Scalar >, Thyra::DefaultSerialDenseLinearOpWithSolve< Scalar >, Thyra::DefaultMultiVectorLinearOpWithSolve< Scalar >, Thyra::DefaultDiagonalLinearOpWithSolve< Scalar >, and Thyra::DefaultAdjointLinearOpWithSolve< Scalar >.
Definition at line 33 of file Thyra_LinearOpWithSolveBase_def.hpp.
|
protectedpure virtual |
Virtual implementation for solve().
Implemented in Thyra::DefaultBlockedTriangularLinearOpWithSolve< Scalar >, Thyra::DelayedLinearOpWithSolve< Scalar >, Thyra::DefaultSerialDenseLinearOpWithSolve< Scalar >, Thyra::DefaultMultiVectorLinearOpWithSolve< Scalar >, Thyra::DefaultDiagonalLinearOpWithSolve< Scalar >, and Thyra::DefaultAdjointLinearOpWithSolve< Scalar >.
|
related |
Call solveSupports()
as a non-member function.
Definition at line 476 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
related |
Call solveSupports()
as a non-member function.
Definition at line 488 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
related |
Call solveSupportsSolveMeasureType()
as a non-member function.
Definition at line 505 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
related |
Call solve()
as a non-member function.
Definition at line 521 of file Thyra_LinearOpWithSolveBase_decl.hpp.
|
related |
Assert that a LOWSB object supports a particular solve type.
This function will throw an excetion with a very good error message if the requested solve type is not supported.