Concrete LinearOpWithSolveBase subclass that adapts any Amesos_BaseSolver object.
More...
#include <Thyra_AmesosLinearOpWithSolve.hpp>
Inherits LinearOpWithSolveBase< double >.
|
| virtual bool | solveSupportsImpl (EOpTransp M_trans) const |
| |
| virtual bool | solveSupportsSolveMeasureTypeImpl (EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const |
| |
| SolveStatus< double > | solveImpl (const EOpTransp M_trans, const MultiVectorBase< double > &B, const Ptr< MultiVectorBase< double > > &X, const Ptr< const SolveCriteria< double > > solveCriteria) const |
| |
|
| | AmesosLinearOpWithSolve () |
| | Construct to uninitialized.
|
| |
| | AmesosLinearOpWithSolve (const Teuchos::RCP< const LinearOpBase< double > > &fwdOp, const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const Teuchos::RCP< Epetra_LinearProblem > &epetraLP, const Teuchos::RCP< Amesos_BaseSolver > &amesosSolver, const EOpTransp amesosSolverTransp, const double amesosSolverScalar) |
| | Calls this->initialize().
|
| |
| void | initialize (const Teuchos::RCP< const LinearOpBase< double > > &fwdOp, const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const Teuchos::RCP< Epetra_LinearProblem > &epetraLP, const Teuchos::RCP< Amesos_BaseSolver > &amesosSolver, const EOpTransp amesosSolverTransp, const double amesosSolverScalar) |
| | First initialization.
|
| |
| Teuchos::RCP< const LinearOpSourceBase< double > > | extract_fwdOpSrc () |
| | Extract the LinearOpSourceBase<double> object so that it can be modified.
|
| |
| Teuchos::RCP< const LinearOpBase< double > > | get_fwdOp () const |
| |
| Teuchos::RCP< const LinearOpSourceBase< double > > | get_fwdOpSrc () const |
| |
| Teuchos::RCP< Epetra_LinearProblem > | get_epetraLP () const |
| |
| Teuchos::RCP< Amesos_BaseSolver > | get_amesosSolver () const |
| |
| EOpTransp | get_amesosSolverTransp () const |
| |
| double | get_amesosSolverScalar () const |
| |
| void | uninitialize (Teuchos::RCP< const LinearOpBase< double > > *fwdOp=NULL, Teuchos::RCP< const LinearOpSourceBase< double > > *fwdOpSrc=NULL, Teuchos::RCP< Epetra_LinearProblem > *epetraLP=NULL, Teuchos::RCP< Amesos_BaseSolver > *amesosSolver=NULL, EOpTransp *amesosSolverTransp=NULL, double *amesosSolverScalar=NULL) |
| | Uninitialize.
|
| |
|
| Teuchos::RCP< const VectorSpaceBase< double > > | range () const |
| |
| Teuchos::RCP< const VectorSpaceBase< double > > | domain () const |
| |
| Teuchos::RCP< const LinearOpBase< double > > | clone () const |
| |
|
| std::string | description () const |
| |
| void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
| |
|
| virtual bool | opSupportedImpl (EOpTransp M_trans) const |
| |
| virtual void | applyImpl (const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha, const double beta) const |
| |
Concrete LinearOpWithSolveBase subclass that adapts any Amesos_BaseSolver object.
See the LinearOpWithSolveBase interface for a description of how to use objects of this type.
Note: Clients should not generally directly create objects of this type but instead should use AmesosLinearOpWithSolveFactory. Only very sophisticated users should ever directly interact with an object through this subclass interface.
Definition at line 36 of file Thyra_AmesosLinearOpWithSolve.hpp.
◆ AmesosLinearOpWithSolve() [1/2]
| Thyra::AmesosLinearOpWithSolve::AmesosLinearOpWithSolve |
( |
| ) |
|
◆ AmesosLinearOpWithSolve() [2/2]
| Thyra::AmesosLinearOpWithSolve::AmesosLinearOpWithSolve |
( |
const Teuchos::RCP< const LinearOpBase< double > > & |
fwdOp, |
|
|
const Teuchos::RCP< const LinearOpSourceBase< double > > & |
fwdOpSrc, |
|
|
const Teuchos::RCP< Epetra_LinearProblem > & |
epetraLP, |
|
|
const Teuchos::RCP< Amesos_BaseSolver > & |
amesosSolver, |
|
|
const EOpTransp |
amesosSolverTransp, |
|
|
const double |
amesosSolverScalar |
|
) |
| |
◆ initialize()
| void Thyra::AmesosLinearOpWithSolve::initialize |
( |
const Teuchos::RCP< const LinearOpBase< double > > & |
fwdOp, |
|
|
const Teuchos::RCP< const LinearOpSourceBase< double > > & |
fwdOpSrc, |
|
|
const Teuchos::RCP< Epetra_LinearProblem > & |
epetraLP, |
|
|
const Teuchos::RCP< Amesos_BaseSolver > & |
amesosSolver, |
|
|
const EOpTransp |
amesosSolverTransp, |
|
|
const double |
amesosSolverScalar |
|
) |
| |
First initialization.
- Parameters
-
| fwdOp | [in] The forward operator for which the factorization exists. |
| epetraLP | [in] The Epetra_LinearProblem object that was used to create the Amesos_BaseSolver object *amesosSolver. Note that the RHS and the LHS multi-vector pointers in this object will be set and unset here. |
| amesosSolver | [in] Contains the factored, and ready to go, Amesos_BaseSolver object ready to solve linear system. |
| amesosSolverTransp | [in] Determines if the Amesos solver should be used as its transpose or not. |
| amesosSolverScalar | [in] Determines the scaling factor associated with the Amesos solver. The solution to the linear solve is scaled by 1/amesosSolverScalar. |
Preconditions:
-
fwdOp.get()!=NULL
-
epetraLP.get()!=NULL
-
amesosSolver.get()!=NULL
-
*epetraLP->GetOperator() is compatible with *fwdOp
-
epetraLP->GetLHS()==NULL
-
epetraLP->GetRHS()==NULL
-
*amesosSolver contains the factorization of *fwdOp and is ready to solve linear systems!
Postconditions:
Definition at line 43 of file Thyra_AmesosLinearOpWithSolve.cpp.
◆ extract_fwdOpSrc()
| Teuchos::RCP< const LinearOpSourceBase< double > > Thyra::AmesosLinearOpWithSolve::extract_fwdOpSrc |
( |
| ) |
|
◆ get_fwdOp()
| Teuchos::RCP< const LinearOpBase< double > > Thyra::AmesosLinearOpWithSolve::get_fwdOp |
( |
| ) |
const |
|
inline |
◆ get_fwdOpSrc()
| Teuchos::RCP< const LinearOpSourceBase< double > > Thyra::AmesosLinearOpWithSolve::get_fwdOpSrc |
( |
| ) |
const |
|
inline |
◆ get_epetraLP()
◆ get_amesosSolver()
| Teuchos::RCP< Amesos_BaseSolver > Thyra::AmesosLinearOpWithSolve::get_amesosSolver |
( |
| ) |
const |
|
inline |
◆ get_amesosSolverTransp()
| EOpTransp Thyra::AmesosLinearOpWithSolve::get_amesosSolverTransp |
( |
| ) |
const |
|
inline |
◆ get_amesosSolverScalar()
| double Thyra::AmesosLinearOpWithSolve::get_amesosSolverScalar |
( |
| ) |
const |
|
inline |
◆ uninitialize()
| void Thyra::AmesosLinearOpWithSolve::uninitialize |
( |
Teuchos::RCP< const LinearOpBase< double > > * |
fwdOp = NULL, |
|
|
Teuchos::RCP< const LinearOpSourceBase< double > > * |
fwdOpSrc = NULL, |
|
|
Teuchos::RCP< Epetra_LinearProblem > * |
epetraLP = NULL, |
|
|
Teuchos::RCP< Amesos_BaseSolver > * |
amesosSolver = NULL, |
|
|
EOpTransp * |
amesosSolverTransp = NULL, |
|
|
double * |
amesosSolverScalar = NULL |
|
) |
| |
◆ range()
| Teuchos::RCP< const VectorSpaceBase< double > > Thyra::AmesosLinearOpWithSolve::range |
( |
| ) |
const |
◆ domain()
| Teuchos::RCP< const VectorSpaceBase< double > > Thyra::AmesosLinearOpWithSolve::domain |
( |
| ) |
const |
◆ clone()
| Teuchos::RCP< const LinearOpBase< double > > Thyra::AmesosLinearOpWithSolve::clone |
( |
| ) |
const |
◆ description()
| std::string Thyra::AmesosLinearOpWithSolve::description |
( |
| ) |
const |
◆ describe()
| void Thyra::AmesosLinearOpWithSolve::describe |
( |
Teuchos::FancyOStream & |
out, |
|
|
const Teuchos::EVerbosityLevel |
verbLevel |
|
) |
| const |
◆ opSupportedImpl()
| bool Thyra::AmesosLinearOpWithSolve::opSupportedImpl |
( |
EOpTransp |
M_trans | ) |
const |
|
protectedvirtual |
◆ applyImpl()
| void Thyra::AmesosLinearOpWithSolve::applyImpl |
( |
const EOpTransp |
M_trans, |
|
|
const MultiVectorBase< double > & |
X, |
|
|
const Ptr< MultiVectorBase< double > > & |
Y, |
|
|
const double |
alpha, |
|
|
const double |
beta |
|
) |
| const |
|
protectedvirtual |
◆ solveSupportsImpl()
| bool Thyra::AmesosLinearOpWithSolve::solveSupportsImpl |
( |
EOpTransp |
M_trans | ) |
const |
|
protectedvirtual |
◆ solveSupportsSolveMeasureTypeImpl()
| bool Thyra::AmesosLinearOpWithSolve::solveSupportsSolveMeasureTypeImpl |
( |
EOpTransp |
M_trans, |
|
|
const SolveMeasureType & |
solveMeasureType |
|
) |
| const |
|
protectedvirtual |
◆ solveImpl()
| SolveStatus< double > Thyra::AmesosLinearOpWithSolve::solveImpl |
( |
const EOpTransp |
M_trans, |
|
|
const MultiVectorBase< double > & |
B, |
|
|
const Ptr< MultiVectorBase< double > > & |
X, |
|
|
const Ptr< const SolveCriteria< double > > |
solveCriteria |
|
) |
| const |
|
protected |
The documentation for this class was generated from the following files: