52 #ifndef AMESOS2_SOLVERCORE_DECL_HPP
53 #define AMESOS2_SOLVERCORE_DECL_HPP
57 #include <Teuchos_RCP.hpp>
58 #include <Teuchos_ParameterList.hpp>
60 #include "Amesos2_Solver.hpp"
61 #include "Amesos2_MatrixTraits.hpp"
62 #include "Amesos2_MatrixAdapter_decl.hpp"
63 #include "Amesos2_MultiVecAdapter_decl.hpp"
102 template <
template <
class,
class>
class ConcreteSolver,
112 typedef ConcreteSolver<Matrix,Vector> solver_type;
113 typedef Matrix matrix_type;
114 typedef Vector vector_type;
135 Teuchos::RCP<Vector> X,
136 Teuchos::RCP<const Vector> B );
237 void solve(
const Teuchos::Ptr<Vector> X,
const Teuchos::Ptr<const Vector> B)
const;
240 void solve(Vector* X,
const Vector* B)
const;
257 void setA(
const Teuchos::RCP<const Matrix> a,
EPhase keep_phase = CLEAN );
259 void setA(
const Matrix* a,
EPhase keep_phase = CLEAN ){
setA(Teuchos::rcp(a), keep_phase); }
303 const Teuchos::RCP<Teuchos::ParameterList> & parameterList );
340 return Teuchos::null;
351 return Teuchos::null;
362 Teuchos::RCP<const Teuchos::Comm<int> >
getComm()
const
381 void describe(Teuchos::FancyOStream &out,
382 const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const;
402 Teuchos::FancyOStream &out,
403 const Teuchos::EVerbosityLevel verbLevel)
const;
414 void getTiming(Teuchos::ParameterList& timingParameterList)
const;
424 std::string
name()
const;
454 Teuchos::RCP<const MatrixAdapter<Matrix> >
matrixA_;
516 #endif // AMESOS2_SOLVERCORE_DECL_HPP
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶meterList)
Set or update internal variables and solver options.
Definition: Amesos2_SolverCore_decl.hpp:327
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers...
Definition: Amesos2_SolverCore_decl.hpp:105
~SolverCore()
Destructor.
Definition: Amesos2_SolverCore_def.hpp:91
const Teuchos::RCP< const Vector > getB()
Returns the vector that is the RHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:273
void setA(const Teuchos::RCP< const Matrix > a, EPhase keep_phase=CLEAN)
Sets the matrix A of this solver.
Definition: Amesos2_SolverCore_def.hpp:239
Holds internal status data about the owning Amesos2 solver.
Definition: Amesos2_Status.hpp:73
EPhase
Used to indicate a phase in the direct solution.
Definition: Amesos2_TypeDecl.hpp:65
global_size_type columnIndexBase_
Index base of column map of matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:487
global_size_type globalNumCols_
Number of global columns in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:478
void printTiming(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Prints timing information about the current solver.
Definition: Amesos2_SolverCore_def.hpp:417
bool root_
If true, then this is the root processor.
Definition: Amesos2_SolverCore_decl.hpp:506
const Teuchos::RCP< Vector > getX()
Returns the vector that is the LHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:265
global_size_type globalNumRows_
Number of global rows in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:475
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Definition: Amesos2_SolverCore_def.hpp:358
std::string name() const
Return the name of this solver.
Definition: Amesos2_SolverCore_def.hpp:509
Control control_
Parameters for solving.
Definition: Amesos2_SolverCore_decl.hpp:494
bool matrixShapeOK()
Returns true if the solver can handle this matrix shape.
Definition: Amesos2_SolverCore_def.hpp:225
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept.
Definition: Amesos2_SolverCore_def.hpp:314
super_type & preOrdering()
Pre-orders the matrix A for minimal fill-in.
Definition: Amesos2_SolverCore_def.hpp:99
std::string description() const
Returns a short description of this Solver.
Definition: Amesos2_SolverCore_def.hpp:348
Status & getStatus() const
Returns a reference to this solver's internal status object.
Definition: Amesos2_SolverCore_decl.hpp:369
int nprocs_
Number of process images in the matrix communicator.
Definition: Amesos2_SolverCore_decl.hpp:509
SolverCore(Teuchos::RCP< const Matrix > A, Teuchos::RCP< Vector > X, Teuchos::RCP< const Vector > B)
Initialize a Solver instance.
Definition: Amesos2_SolverCore_def.hpp:66
void setX(Vector *x)
Sets the LHS vector X using a raw pointer.
Definition: Amesos2_SolverCore_decl.hpp:263
Container class for control variables.
void solve()
Solves (or )
Definition: Amesos2_SolverCore_def.hpp:163
super_type & operator=(const solver_type *rhs)
Do not allow copying of this Solver by assignment.
super_type & numericFactorization()
Performs numeric factorization on the matrix A.
Definition: Amesos2_SolverCore_def.hpp:140
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns a pointer to the Teuchos::Comm communicator with this operator.
Definition: Amesos2_SolverCore_decl.hpp:362
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
This is an empty stub.
Definition: Amesos2_SolverCore_decl.hpp:349
void setB(const Teuchos::RCP< const Vector > b)
Sets the RHS vector B.
Definition: Amesos2_SolverCore_decl.hpp:269
bool matrix_loaded_
Definition: Amesos2_SolverCore_decl.hpp:461
A Matrix adapter interface for Amesos2.
Definition: Amesos2_MatrixAdapter_decl.hpp:76
void loadA(EPhase current_phase)
Refresh this solver's internal data about A.
Definition: Amesos2_SolverCore_def.hpp:517
Interface to Amesos2 solver objects.
Definition: Amesos2_Solver_decl.hpp:78
void getTiming(Teuchos::ParameterList &timingParameterList) const
Extracts timing information from the current solver.
Definition: Amesos2_SolverCore_def.hpp:499
global_size_type globalNumNonZeros_
Number of global non-zero values in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:481
void setNnzLU(size_t nnz)
Set the number of non-zero values in the and factors.
Definition: Amesos2_SolverCore_decl.hpp:451
void setX(const Teuchos::RCP< Vector > x)
Sets the LHS vector X.
Definition: Amesos2_SolverCore_decl.hpp:261
void setB(const Vector *b)
Sets the RHS vector B using a raw pointer.
Definition: Amesos2_SolverCore_decl.hpp:271
super_type & setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶meterList)
Set/update internal variables and solver options.
Definition: Amesos2_SolverCore_def.hpp:282
Enum and other types declarations for Amesos2.
int rank_
The MPI rank of this image.
Definition: Amesos2_SolverCore_decl.hpp:503
Vector * getXRaw()
Returns a raw pointer to the LHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:267
Teuchos::RCP< const Vector > multiVecB_
The RHS vector/multi-vector.
Definition: Amesos2_SolverCore_decl.hpp:472
Container class for Timers used with the Amesos2::Solver class.
Timers timers_
Various timing statistics.
Definition: Amesos2_SolverCore_decl.hpp:497
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
This is a empty stub.
Definition: Amesos2_SolverCore_decl.hpp:338
size_t lu_nnz_
The number of non-zeros in the factors.
Definition: Amesos2_Status.hpp:149
super_type & symbolicFactorization()
Performs symbolic factorization on the matrix A.
Definition: Amesos2_SolverCore_def.hpp:117
Container class for status variables.
Teuchos::RCP< Vector > multiVecX_
The LHS vector/multi-vector.
Definition: Amesos2_SolverCore_decl.hpp:465
Teuchos::RCP< const MatrixAdapter< Matrix > > matrixA_
The LHS operator.
Definition: Amesos2_SolverCore_decl.hpp:454
const Vector * getBRaw()
Returns a raw pointer to the RHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:275
global_size_type rowIndexBase_
Index base of rowmap of matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:484
Status status_
Holds status information about a solver.
Definition: Amesos2_SolverCore_decl.hpp:491
void setA(const Matrix *a, EPhase keep_phase=CLEAN)
Sets the matrix A of this solver.
Definition: Amesos2_SolverCore_decl.hpp:259