Stratimikos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Private Member Functions | Private Attributes | List of all members
Thyra::AmesosLinearOpWithSolve Class Reference

Concrete LinearOpWithSolveBase subclass that adapts any Amesos_BaseSolver object. More...

#include <Thyra_AmesosLinearOpWithSolve.hpp>

Inherits LinearOpWithSolveBase< double >.

Private Member Functions

void assertInitialized () const
 

Private Attributes

Teuchos::RCP< const
LinearOpBase< double > > 
fwdOp_
 
Teuchos::RCP< const
LinearOpSourceBase< double > > 
fwdOpSrc_
 
Teuchos::RCP
< Epetra_LinearProblem > 
epetraLP_
 
Teuchos::RCP< Amesos_BaseSolver > amesosSolver_
 
EOpTransp amesosSolverTransp_
 
double amesosSolverScalar_
 

Constructors/initializers/accessors

 AmesosLinearOpWithSolve ()
 Construct to uninitialized. More...
 
 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(). More...
 
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. More...
 
Teuchos::RCP< const
LinearOpSourceBase< double > > 
extract_fwdOpSrc ()
 Extract the LinearOpSourceBase<double> object so that it can be modified. More...
 
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. More...
 

Overridden public functions from LinearOpBase

Teuchos::RCP< const
VectorSpaceBase< double > > 
range () const
 
Teuchos::RCP< const
VectorSpaceBase< double > > 
domain () const
 
Teuchos::RCP< const
LinearOpBase< double > > 
clone () const
 

Overridden public functions from Teuchos::Describable

std::string description () const
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
 

Overridden from LinearOpBase

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
 

Overridden from LinearOpWithSolveBase.

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
 

Detailed Description

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.

Constructor & Destructor Documentation

Thyra::AmesosLinearOpWithSolve::AmesosLinearOpWithSolve ( )

Construct to uninitialized.

Definition at line 23 of file Thyra_AmesosLinearOpWithSolve.cpp.

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 
)

Calls this->initialize().

Definition at line 29 of file Thyra_AmesosLinearOpWithSolve.cpp.

Member Function Documentation

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.

Teuchos::RCP< const LinearOpSourceBase< double > > Thyra::AmesosLinearOpWithSolve::extract_fwdOpSrc ( )

Extract the LinearOpSourceBase<double> object so that it can be modified.

Postconditions:

Definition at line 73 of file Thyra_AmesosLinearOpWithSolve.cpp.

Teuchos::RCP< const LinearOpBase< double > > Thyra::AmesosLinearOpWithSolve::get_fwdOp ( ) const
inline

Definition at line 217 of file Thyra_AmesosLinearOpWithSolve.hpp.

Teuchos::RCP< const LinearOpSourceBase< double > > Thyra::AmesosLinearOpWithSolve::get_fwdOpSrc ( ) const
inline

Definition at line 224 of file Thyra_AmesosLinearOpWithSolve.hpp.

Teuchos::RCP< Epetra_LinearProblem > Thyra::AmesosLinearOpWithSolve::get_epetraLP ( ) const
inline

Definition at line 231 of file Thyra_AmesosLinearOpWithSolve.hpp.

Teuchos::RCP< Amesos_BaseSolver > Thyra::AmesosLinearOpWithSolve::get_amesosSolver ( ) const
inline

Definition at line 238 of file Thyra_AmesosLinearOpWithSolve.hpp.

EOpTransp Thyra::AmesosLinearOpWithSolve::get_amesosSolverTransp ( ) const
inline

Definition at line 244 of file Thyra_AmesosLinearOpWithSolve.hpp.

double Thyra::AmesosLinearOpWithSolve::get_amesosSolverScalar ( ) const
inline

Definition at line 250 of file Thyra_AmesosLinearOpWithSolve.hpp.

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 
)

Uninitialize.

Definition at line 82 of file Thyra_AmesosLinearOpWithSolve.cpp.

Teuchos::RCP< const VectorSpaceBase< double > > Thyra::AmesosLinearOpWithSolve::range ( ) const

Definition at line 113 of file Thyra_AmesosLinearOpWithSolve.cpp.

Teuchos::RCP< const VectorSpaceBase< double > > Thyra::AmesosLinearOpWithSolve::domain ( ) const

Definition at line 120 of file Thyra_AmesosLinearOpWithSolve.cpp.

Teuchos::RCP< const LinearOpBase< double > > Thyra::AmesosLinearOpWithSolve::clone ( ) const

Definition at line 127 of file Thyra_AmesosLinearOpWithSolve.cpp.

std::string Thyra::AmesosLinearOpWithSolve::description ( ) const

Definition at line 136 of file Thyra_AmesosLinearOpWithSolve.cpp.

void Thyra::AmesosLinearOpWithSolve::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel 
) const

Definition at line 149 of file Thyra_AmesosLinearOpWithSolve.cpp.

bool Thyra::AmesosLinearOpWithSolve::opSupportedImpl ( EOpTransp  M_trans) const
protectedvirtual

Definition at line 191 of file Thyra_AmesosLinearOpWithSolve.cpp.

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

Definition at line 197 of file Thyra_AmesosLinearOpWithSolve.cpp.

bool Thyra::AmesosLinearOpWithSolve::solveSupportsImpl ( EOpTransp  M_trans) const
protectedvirtual

Definition at line 212 of file Thyra_AmesosLinearOpWithSolve.cpp.

bool Thyra::AmesosLinearOpWithSolve::solveSupportsSolveMeasureTypeImpl ( EOpTransp  M_trans,
const SolveMeasureType &  solveMeasureType 
) const
protectedvirtual

Definition at line 232 of file Thyra_AmesosLinearOpWithSolve.cpp.

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

Definition at line 241 of file Thyra_AmesosLinearOpWithSolve.cpp.

void Thyra::AmesosLinearOpWithSolve::assertInitialized ( ) const
private

Member Data Documentation

Teuchos::RCP<const LinearOpBase<double> > Thyra::AmesosLinearOpWithSolve::fwdOp_
private

Definition at line 201 of file Thyra_AmesosLinearOpWithSolve.hpp.

Teuchos::RCP<const LinearOpSourceBase<double> > Thyra::AmesosLinearOpWithSolve::fwdOpSrc_
private

Definition at line 202 of file Thyra_AmesosLinearOpWithSolve.hpp.

Teuchos::RCP<Epetra_LinearProblem> Thyra::AmesosLinearOpWithSolve::epetraLP_
private

Definition at line 203 of file Thyra_AmesosLinearOpWithSolve.hpp.

Teuchos::RCP<Amesos_BaseSolver> Thyra::AmesosLinearOpWithSolve::amesosSolver_
private

Definition at line 204 of file Thyra_AmesosLinearOpWithSolve.hpp.

EOpTransp Thyra::AmesosLinearOpWithSolve::amesosSolverTransp_
private

Definition at line 205 of file Thyra_AmesosLinearOpWithSolve.hpp.

double Thyra::AmesosLinearOpWithSolve::amesosSolverScalar_
private

Definition at line 206 of file Thyra_AmesosLinearOpWithSolve.hpp.


The documentation for this class was generated from the following files: