Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | Protected Attributes | List of all members
Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver > Class Template Reference

Amesos2 solver adapter for MP::Vector scalar type. More...

#include <Amesos2_Solver_MP_Vector.hpp>

Inherits Solver< Tpetra::CrsMatrix< Sacado::MP::Vector< Storage >, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::MultiVector< Sacado::MP::Vector< Storage >, LocalOrdinal, GlobalOrdinal, Node > >.

Public Types

typedef Sacado::MP::Vector
< Storage
Scalar
 
typedef Tpetra::CrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node
Matrix
 
typedef Tpetra::MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node
Vector
 
typedef Scalar::value_type BaseScalar
 
typedef Tpetra::Map
< LocalOrdinal, GlobalOrdinal,
Node
Map
 
typedef Tpetra::CrsGraph
< LocalOrdinal, GlobalOrdinal,
Node
FlatGraph
 
typedef Tpetra::CrsMatrix
< BaseScalar, LocalOrdinal,
GlobalOrdinal, Node
FlatMatrix
 
typedef Tpetra::MultiVector
< BaseScalar, LocalOrdinal,
GlobalOrdinal, Node
FlatVector
 
typedef ConcreteSolver
< FlatMatrix, FlatVector
FlatConcreteSolver
 
typedef Solver< FlatMatrix,
FlatVector
FlatSolver
 
typedef Solver< Matrix, Vectorsolver_type
 
typedef solver_type::type type
 

Public Member Functions

 MPVectorSolverAdapter (const Teuchos::RCP< const Matrix > &A_, const Teuchos::RCP< Vector > &X_, const Teuchos::RCP< const Vector > &B_)
 Constructor. More...
 

Protected Attributes

Teuchos::RCP< const MatrixA
 
Teuchos::RCP< VectorX
 
Teuchos::RCP< const VectorB
 
Teuchos::RCP< const Mapflat_X_map
 
Teuchos::RCP< const Mapflat_B_map
 
Teuchos::RCP< const FlatGraphflat_graph
 
Teuchos::RCP< const FlatMatrixflat_A
 
Teuchos::RCP< FlatVectorflat_X
 
Teuchos::RCP< const FlatVectorflat_B
 
Teuchos::RCP< FlatSolverflat_solver
 

Mathematical Functions

virtual typepreOrdering (void)
 Pre-orders the matrix. More...
 
virtual typesymbolicFactorization (void)
 Performs symbolic factorization on the matrix. More...
 
virtual typenumericFactorization (void)
 Performs numeric factorization on the matrix. More...
 
virtual void solve (void)
 Solves $ A X = B$ (or $ A^T X = B$ ) More...
 
virtual void solve (const Teuchos::Ptr< Vector > XX, const Teuchos::Ptr< const Vector > BB) const
 Solve $ A X = B$ using the given XX and BB (multi)vectors. More...
 
virtual void solve (Vector *XX, const Vector *BB) const
 Solve $ A X = B$ using the given XX and BB (multi)vectors. More...
 

Parameter Methods

virtual typesetParameters (const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
 Set/update internal variables and solver options. More...
 
virtual Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters (void) const
 Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept. More...
 

Accessor Methods

virtual void setA (const Teuchos::RCP< const Matrix > a, EPhase keep_phase=CLEAN)
 Sets the matrix A of this solver. More...
 
virtual void setA (const Matrix *a, EPhase keep_phase=CLEAN)
 Sets the matrix A of this solver. More...
 
virtual bool matrixShapeOK (void)
 Returns true if the solver can handle the matrix shape. More...
 
virtual void setX (const Teuchos::RCP< Vector > x)
 Sets the LHS vector X. More...
 
virtual void setX (Vector *x)
 Sets the LHS vector X using a raw pointer. More...
 
virtual const Teuchos::RCP
< Vector
getX (void)
 Returns the vector that is the LHS of the linear system. More...
 
virtual VectorgetXRaw (void)
 Returns a raw pointer to the LHS of the linear system. More...
 
virtual void setB (const Teuchos::RCP< const Vector > b)
 Sets the RHS vector B. More...
 
virtual void setB (const Vector *b)
 Sets the RHS vector B using a raw pointer. More...
 
virtual const Teuchos::RCP
< const Vector
getB (void)
 Returns the vector that is the RHS of the linear system. More...
 
virtual const VectorgetBRaw (void)
 Returns a raw pointer to the RHS of the linear system. More...
 
virtual Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm (void) const
 Returns a pointer to the Teuchos::Comm communicator with this matrix. More...
 
virtual Status & getStatus () const
 Returns a reference to this solver's internal status object. More...
 
virtual std::string name (void) const
 Return the name of this solver. More...
 

Methods implementing Describable

virtual std::string description (void) const
 Returns a short description of this Solver. More...
 
virtual void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 

Performance and Timing

virtual void printTiming (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Prints timing information about the current solver. More...
 
virtual void getTiming (Teuchos::ParameterList &timingParameterList) const
 Extracts timing information from the current solver. More...
 

Detailed Description

template<class Storage, class LocalOrdinal, class GlobalOrdinal, class Node, template< class, class > class ConcreteSolver>
class Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >

Amesos2 solver adapter for MP::Vector scalar type.

This adapter enables Amesos2 solvers to work with Tpetra matrices and vectors of the Sacado::MP::Vector scalar type by "flattening" these matrices and vectors into ones with a standard (e.g., double) scalar type.

Definition at line 46 of file Amesos2_Solver_MP_Vector.hpp.

Member Typedef Documentation

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Sacado::MP::Vector<Storage> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::Scalar

Definition at line 59 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::Matrix

Definition at line 60 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::Vector

Definition at line 61 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Scalar::value_type Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::BaseScalar

Definition at line 63 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::Map

Definition at line 64 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatGraph

Definition at line 65 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Tpetra::CrsMatrix<BaseScalar,LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatMatrix

Definition at line 66 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Tpetra::MultiVector<BaseScalar,LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatVector

Definition at line 67 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef ConcreteSolver<FlatMatrix,FlatVector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatConcreteSolver

Definition at line 68 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Solver<FlatMatrix,FlatVector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatSolver

Definition at line 69 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef Solver<Matrix,Vector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::solver_type

Definition at line 71 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
typedef solver_type::type Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::type

Definition at line 72 of file Amesos2_Solver_MP_Vector.hpp.

Constructor & Destructor Documentation

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::MPVectorSolverAdapter ( const Teuchos::RCP< const Matrix > &  A_,
const Teuchos::RCP< Vector > &  X_,
const Teuchos::RCP< const Vector > &  B_ 
)
inline

Constructor.

Definition at line 75 of file Amesos2_Solver_MP_Vector.hpp.

Member Function Documentation

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual type& Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::preOrdering ( void  )
inlinevirtual

Pre-orders the matrix.

Uses the default solver option, unless a solver-specific pre-ordering parameter is given.

See Also
setParameters

Definition at line 105 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual type& Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::symbolicFactorization ( void  )
inlinevirtual

Performs symbolic factorization on the matrix.

Precondition
  • The matrix A must not be null

Definition at line 116 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual type& Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::numericFactorization ( void  )
inlinevirtual

Performs numeric factorization on the matrix.

numericFactorization checks first that symbolicFactorization has successfully been called, and if not, calls it before continuing.

Precondition
  • The matrix A must not be null
Postcondition
  • The factors L and U of A are computed

Definition at line 134 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::solve ( void  )
inlinevirtual

Solves $ A X = B$ (or $ A^T X = B$ )

solve checks first that numericFactorization has successfully been called, and if not, calls it before continuing.

Precondition
  • The (multi)vectors X and B must not be null
Postcondition
  • The (multi)vector X (given at construction time) contains the solution to the system.

Definition at line 152 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::solve ( const Teuchos::Ptr< Vector XX,
const Teuchos::Ptr< const Vector BB 
) const
inlinevirtual

Solve $ A X = B$ using the given XX and BB (multi)vectors.

This overload of solve uses the given XX and BB (multi)vectors when solving. These XX and BB (multi)vectors are used in place of any X and B that were given upon construction of the Amesos2 solver instance. XX and BB are used only for this solve.

If a permanent change of X and B are required, see the setX() and setB() methods.

Postcondition
  • The (multi)vector XX contains the solution to the system
  • The (multi)vectors X and B given at construction time (if any) are unchanged.

Definition at line 174 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::solve ( Vector XX,
const Vector BB 
) const
inlinevirtual

Solve $ A X = B$ using the given XX and BB (multi)vectors.

This overload of solve uses the given XX and BB (multi)vectors when solving. These XX and BB (multi)vectors are used in place of any X and B that were given upon construction of the Amesos2 solver instance. XX and BB are used only for this solve.

If a permanent change of X and B are required, see the setX() and setB() methods.

Postcondition
  • The (multi)vector XX contains the solution to the system
  • The (multi)vectors X and B given at construction time (if any) are unchanged.

Definition at line 199 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual type& Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setParameters ( const Teuchos::RCP< Teuchos::ParameterList > &  parameterList)
inlinevirtual

Set/update internal variables and solver options.

Expects that parameterList be named "Amesos2". That list may contain Amesos2-specific parameters. In addition, it may contain sublist for solver-specific parameters. These sublists should be named according to what is returned by the name() function (i.e. The solver's name when enabling for Amesos2 during configuration).

See each solver interface directly for a list of the supported parameters for that solver.

Definition at line 224 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual Teuchos::RCP<const Teuchos::ParameterList> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::getValidParameters ( void  ) const
inlinevirtual

Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept.

Definition at line 236 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setA ( const Teuchos::RCP< const Matrix a,
EPhase  keep_phase = CLEAN 
)
inlinevirtual

Sets the matrix A of this solver.

Parameters
[in]aAn RCP to a matrix will will be used for future computation steps
[in]keep_phaseThis parameter tells the solver what state it should keep. For example, you may want to replace the matrix but keep the symbolic factorization because you know the structure of the new matrix is the same as the structure of the old matrix. In this case you would pass Amesos2::SYMBFACT as this parameter.

The default value for the second parameter is Amesos2::CLEAN, which means that the internal state of the solver will be completely reset. It will be as if no previous computational steps were performed.

Definition at line 266 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setA ( const Matrix a,
EPhase  keep_phase = CLEAN 
)
inlinevirtual

Sets the matrix A of this solver.

Parameters
[in]aAn raw C pointer to a matrix will will be used for future computation steps.
[in]keep_phaseThis parameter tells the solver what state it should keep. For example, you may want to replace the matrix but keep the symbolic factorization because you know the structure of the new matrix is the same as the structure of the old matrix. In this case you would pass Amesos2::SYMBFACT as this parameter.

The default value for the second parameter is Amesos2::CLEAN, which means that the internal state of the solver will be completely reset. It will be as if no previous computational steps were performed.

Definition at line 306 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual bool Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::matrixShapeOK ( void  )
inlinevirtual

Returns true if the solver can handle the matrix shape.

Definition at line 312 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setX ( const Teuchos::RCP< Vector x)
inlinevirtual

Sets the LHS vector X.

Definition at line 318 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setX ( Vector x)
inlinevirtual

Sets the LHS vector X using a raw pointer.

Definition at line 329 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual const Teuchos::RCP<Vector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::getX ( void  )
inlinevirtual

Returns the vector that is the LHS of the linear system.

Definition at line 343 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual Vector* Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::getXRaw ( void  )
inlinevirtual

Returns a raw pointer to the LHS of the linear system.

Definition at line 349 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setB ( const Teuchos::RCP< const Vector b)
inlinevirtual

Sets the RHS vector B.

Definition at line 355 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setB ( const Vector b)
inlinevirtual

Sets the RHS vector B using a raw pointer.

Definition at line 366 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual const Teuchos::RCP<const Vector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::getB ( void  )
inlinevirtual

Returns the vector that is the RHS of the linear system.

Definition at line 380 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual const Vector* Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::getBRaw ( void  )
inlinevirtual

Returns a raw pointer to the RHS of the linear system.

Definition at line 386 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual Teuchos::RCP<const Teuchos::Comm<int> > Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::getComm ( void  ) const
inlinevirtual

Returns a pointer to the Teuchos::Comm communicator with this matrix.

Definition at line 392 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual Status& Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::getStatus ( ) const
inlinevirtual

Returns a reference to this solver's internal status object.

Definition at line 398 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual std::string Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::name ( void  ) const
inlinevirtual

Return the name of this solver.

Definition at line 404 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual std::string Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::description ( void  ) const
inlinevirtual

Returns a short description of this Solver.

Definition at line 416 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
inlinevirtual

Prints the status information about the current solver with some level of verbosity.

Definition at line 423 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::printTiming ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
inlinevirtual

Prints timing information about the current solver.

Definition at line 436 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
virtual void Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::getTiming ( Teuchos::ParameterList timingParameterList) const
inlinevirtual

Extracts timing information from the current solver.

Results are placed in the parameter list timingParameterList

Parameters
timingParameterListAccepts timing information from the current solver

Definition at line 450 of file Amesos2_Solver_MP_Vector.hpp.

Member Data Documentation

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<const Matrix> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::A
protected

Definition at line 458 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<Vector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::X
protected

Definition at line 459 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<const Vector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::B
protected

Definition at line 460 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<const Map> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::flat_X_map
protected

Definition at line 461 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<const Map> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::flat_B_map
protected

Definition at line 461 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<const FlatGraph> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::flat_graph
protected

Definition at line 462 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<const FlatMatrix> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::flat_A
protected

Definition at line 463 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<FlatVector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::flat_X
protected

Definition at line 464 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<const FlatVector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::flat_B
protected

Definition at line 465 of file Amesos2_Solver_MP_Vector.hpp.

template<class Storage , class LocalOrdinal , class GlobalOrdinal , class Node , template< class, class > class ConcreteSolver>
Teuchos::RCP<FlatSolver> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::flat_solver
protected

Definition at line 466 of file Amesos2_Solver_MP_Vector.hpp.


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