Stokhos Package Browser (Single Doxygen Collection)
Version of the Day
|
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, 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. 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 |
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 (or ) More... | |
virtual void | solve (const Teuchos::Ptr< Vector > XX, const Teuchos::Ptr< const Vector > BB) const |
Solve using the given XX and BB (multi)vectors. More... | |
virtual void | solve (Vector *XX, const Vector *BB) const |
Solve using the given XX and BB (multi)vectors. More... | |
Parameter Methods | |
virtual type & | setParameters (const Teuchos::RCP< Teuchos::ParameterList > ¶meterList) |
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 Vector * | getXRaw (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 Vector * | getBRaw (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... | |
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 90 of file Amesos2_Solver_MP_Vector.hpp.
typedef Sacado::MP::Vector<Storage> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::Scalar |
Definition at line 103 of file Amesos2_Solver_MP_Vector.hpp.
typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::Matrix |
Definition at line 104 of file Amesos2_Solver_MP_Vector.hpp.
typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::Vector |
Definition at line 105 of file Amesos2_Solver_MP_Vector.hpp.
typedef Scalar::value_type Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::BaseScalar |
Definition at line 107 of file Amesos2_Solver_MP_Vector.hpp.
typedef Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::Map |
Definition at line 108 of file Amesos2_Solver_MP_Vector.hpp.
typedef Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatGraph |
Definition at line 109 of file Amesos2_Solver_MP_Vector.hpp.
typedef Tpetra::CrsMatrix<BaseScalar,LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatMatrix |
Definition at line 110 of file Amesos2_Solver_MP_Vector.hpp.
typedef Tpetra::MultiVector<BaseScalar,LocalOrdinal,GlobalOrdinal,Node> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatVector |
Definition at line 111 of file Amesos2_Solver_MP_Vector.hpp.
typedef ConcreteSolver<FlatMatrix,FlatVector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatConcreteSolver |
Definition at line 112 of file Amesos2_Solver_MP_Vector.hpp.
typedef Solver<FlatMatrix,FlatVector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::FlatSolver |
Definition at line 113 of file Amesos2_Solver_MP_Vector.hpp.
typedef Solver<Matrix,Vector> Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::solver_type |
Definition at line 115 of file Amesos2_Solver_MP_Vector.hpp.
typedef solver_type::type Amesos2::MPVectorSolverAdapter< Storage, LocalOrdinal, GlobalOrdinal, Node, ConcreteSolver >::type |
Definition at line 116 of file Amesos2_Solver_MP_Vector.hpp.
|
inline |
Constructor.
Definition at line 119 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Pre-orders the matrix.
Uses the default solver option, unless a solver-specific pre-ordering parameter is given.
Definition at line 149 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Performs symbolic factorization on the matrix.
null
Definition at line 160 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Performs numeric factorization on the matrix.
numericFactorization checks first that symbolicFactorization has successfully been called, and if not, calls it before continuing.
null
Definition at line 178 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Solves (or )
solve checks first that numericFactorization has successfully been called, and if not, calls it before continuing.
X
and B
must not be null
X
(given at construction time) contains the solution to the system. Definition at line 196 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Solve 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.
XX
contains the solution to the systemX
and B
given at construction time (if any) are unchanged. Definition at line 218 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Solve 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.
XX
contains the solution to the systemX
and B
given at construction time (if any) are unchanged. Definition at line 243 of file Amesos2_Solver_MP_Vector.hpp.
|
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 268 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept.
Definition at line 280 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Sets the matrix A of this solver.
[in] | a | An RCP to a matrix will will be used for future computation steps |
[in] | keep_phase | This 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 310 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Sets the matrix A of this solver.
[in] | a | An raw C pointer to a matrix will will be used for future computation steps. |
[in] | keep_phase | This 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 350 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Returns true
if the solver can handle the matrix shape.
Definition at line 356 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Sets the LHS vector X.
Definition at line 362 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Sets the LHS vector X using a raw pointer.
Definition at line 373 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Returns the vector that is the LHS of the linear system.
Definition at line 387 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Returns a raw pointer to the LHS of the linear system.
Definition at line 393 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Sets the RHS vector B.
Definition at line 399 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Sets the RHS vector B using a raw pointer.
Definition at line 410 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Returns the vector that is the RHS of the linear system.
Definition at line 424 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Returns a raw pointer to the RHS of the linear system.
Definition at line 430 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Returns a pointer to the Teuchos::Comm communicator with this matrix.
Definition at line 436 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Returns a reference to this solver's internal status object.
Definition at line 442 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Return the name of this solver.
Definition at line 448 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Returns a short description of this Solver.
Definition at line 460 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Prints the status information about the current solver with some level of verbosity.
Definition at line 467 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Prints timing information about the current solver.
Definition at line 480 of file Amesos2_Solver_MP_Vector.hpp.
|
inlinevirtual |
Extracts timing information from the current solver.
Results are placed in the parameter list timingParameterList
timingParameterList | Accepts timing information from the current solver |
Definition at line 494 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 502 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 503 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 504 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 505 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 505 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 506 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 507 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 508 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 509 of file Amesos2_Solver_MP_Vector.hpp.
|
protected |
Definition at line 510 of file Amesos2_Solver_MP_Vector.hpp.