Stokhos  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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>

Inheritance diagram for Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >:
Inheritance graph
[legend]
Collaboration diagram for Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >:
Collaboration graph
[legend]

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, Vector > solver_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.
 
Mathematical Functions
virtual type & preOrdering (void)
 Pre-orders the matrix. More...
 
virtual type & symbolicFactorization (void)
 Performs symbolic factorization on the matrix. More...
 
virtual type & numericFactorization (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 type & setParameters (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.
 
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.
 
virtual void setX (const Teuchos::RCP< Vector > x)
 Sets the LHS vector X.
 
virtual void setX (Vector *x)
 Sets the LHS vector X using a raw pointer.
 
virtual const Teuchos::RCP
< Vector > 
getX (void)
 Returns the vector that is the LHS of the linear system.
 
virtual Vector * getXRaw (void)
 Returns a raw pointer to the LHS of the linear system.
 
virtual void setB (const Teuchos::RCP< const Vector > b)
 Sets the RHS vector B.
 
virtual void setB (const Vector *b)
 Sets the RHS vector B using a raw pointer.
 
virtual const Teuchos::RCP
< const Vector > 
getB (void)
 Returns the vector that is the RHS of the linear system.
 
virtual const Vector * getBRaw (void)
 Returns a raw pointer to the RHS of the linear system.
 
virtual Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm (void) const
 Returns a pointer to the Teuchos::Comm communicator with this matrix.
 
virtual Status & getStatus () const
 Returns a reference to this solver's internal status object.
 
virtual std::string name (void) const
 Return the name of this solver.
 
Methods implementing Describable
virtual std::string description (void) const
 Returns a short description of this Solver.
 
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.
 
virtual void getTiming (Teuchos::ParameterList &timingParameterList) const
 Extracts timing information from the current solver. More...
 

Protected Attributes

Teuchos::RCP< const Matrix > A
 
Teuchos::RCP< Vector > X
 
Teuchos::RCP< const Vector > B
 
Teuchos::RCP< const Map > flat_X_map
 
Teuchos::RCP< const Map > flat_B_map
 
Teuchos::RCP< const FlatGraph > flat_graph
 
Teuchos::RCP< const FlatMatrix > flat_A
 
Teuchos::RCP< FlatVector > flat_X
 
Teuchos::RCP< const FlatVector > flat_B
 
Teuchos::RCP< FlatSolver > flat_solver
 

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.

Member Function Documentation

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.

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
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
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
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.

Referenced by Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setA().

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.

References Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::setA().

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.

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.
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.
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.
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

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