Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Amesos2_SolverCore_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Amesos2: Templated Direct Sparse Solver Package
4 //
5 // Copyright 2011 NTESS and the Amesos2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
18 #ifndef AMESOS2_SOLVERCORE_DECL_HPP
19 #define AMESOS2_SOLVERCORE_DECL_HPP
20 
21 #include <string>
22 
23 #include <Teuchos_RCP.hpp>
24 #include <Teuchos_ParameterList.hpp>
25 
26 #include "Amesos2_Solver.hpp"
27 #include "Amesos2_MatrixTraits.hpp"
28 #include "Amesos2_MatrixAdapter_decl.hpp"
29 #include "Amesos2_MultiVecAdapter_decl.hpp"
30 #include "Amesos2_TypeDecl.hpp"
31 
32 #include "Amesos2_Control.hpp"
33 #include "Amesos2_Status.hpp"
34 #include "Amesos2_Timers.hpp"
35 
36 namespace Amesos2 {
37 
38 
39  /* This is the base class to be used in a *statically* polymorphic
40  * way. E.g. for the Superlu solver:
41  *
42  * In Amesos2_Superlu.hpp:
43  * class Superlu : SolverCore<Superlu> { ... }
44  *
45  * Each concrete solver will implement several private sub-functions
46  * that will be called within the common code for each function.
47  */
48 
68  template <template <class,class> class ConcreteSolver,
69  class Matrix,
70  class Vector >
71  class SolverCore : public Amesos2::Solver<Matrix,Vector>
72  {
73  public:
74 
75  // Grant public access to contained types
78  typedef ConcreteSolver<Matrix,Vector> solver_type;
79  typedef Matrix matrix_type;
80  typedef Vector vector_type;
81  typedef typename MatrixAdapter<matrix_type>::scalar_t scalar_type;
82  typedef typename MatrixAdapter<matrix_type>::local_ordinal_t local_ordinal_type;
83  typedef typename MatrixAdapter<matrix_type>::global_ordinal_t global_ordinal_type;
84  typedef typename MatrixAdapter<matrix_type>::global_size_t global_size_type;
85  typedef typename MatrixAdapter<matrix_type>::node_t node_type;
86 
88 
89 
100  SolverCore( Teuchos::RCP<const Matrix> A,
101  Teuchos::RCP<Vector> X,
102  Teuchos::RCP<const Vector> B );
103 
104 
106  ~SolverCore( );
107 
108 
110  SolverCore(const solver_type& rhs);
111 
112 
114  super_type& operator=(const solver_type* rhs);
115 
117 
119 
120 
132  super_type& preOrdering() override;
133 
134 
152  super_type& symbolicFactorization() override;
153 
154 
178  super_type& numericFactorization() override;
179 
180 
200  void solve() override;
201  void solve(const Teuchos::Ptr<Vector> X, const Teuchos::Ptr<const Vector> B) const override;
202  void solve(Vector* X, const Vector* B) const override;
203 
204  int solve_ir(const Teuchos::Ptr<Vector> X, const Teuchos::Ptr<const Vector> B,
205  const int maxNumIters, const bool verbose) const;
206  int solve_ir(Vector* X, const Vector* B,
207  const int maxNumIters, const bool verbose) const;
208  int solve_ir(const int maxNumIters, const bool verbose);
209 
211 
212 
222  bool matrixShapeOK() override;
223 
224  void setA( const Teuchos::RCP<const Matrix> a, EPhase keep_phase = CLEAN ) override;
225 
226  void setA( const Matrix* a, EPhase keep_phase = CLEAN ) override { setA(Teuchos::rcp(a), keep_phase); }
227 
228  void setX(const Teuchos::RCP<Vector> x) override { multiVecX_ = x; }
229 
230  void setX(Vector* x) override { multiVecX_ = Teuchos::rcp(x); }
231 
232  const Teuchos::RCP<Vector> getX() override { return( multiVecX_ ); }
233 
234  Vector* getXRaw() override { return multiVecX_.getRawPtr(); }
235 
236  void setB(const Teuchos::RCP<const Vector> b) override { multiVecB_ = b; }
237 
238  void setB(const Vector* b) override { multiVecB_ = Teuchos::rcp(b); }
239 
240  const Teuchos::RCP<const Vector> getB() override { return( multiVecB_ ); }
241 
242  const Vector* getBRaw() override { return multiVecB_.getRawPtr(); }
243 
244 
246 
247 
269  super_type& setParameters(
270  const Teuchos::RCP<Teuchos::ParameterList> & parameterList ) override;
271 
272 
282  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const override;
283 
284 
294  void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & parameterList)
295  {
296  setParameters(parameterList);
297  }
298 
299 
305  Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList()
306  {
307  return Teuchos::null;
308  }
309 
310 
316  Teuchos::RCP<Teuchos::ParameterList> unsetParameterList()
317  {
318  return Teuchos::null;
319  }
320 
321 
323 
324 
326 
327 
329  Teuchos::RCP<const Teuchos::Comm<int> > getComm() const override
330  {
331  return matrixA_->getComm();
332  }
333 
334 
336  inline Status& getStatus() const override { return( status_ ); }
337 
338 
340 
341 
343  std::string description() const override;
344 
345 
348  void describe(Teuchos::FancyOStream &out,
349  const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override;
350 
351 
368  void printTiming(
369  Teuchos::FancyOStream &out,
370  const Teuchos::EVerbosityLevel verbLevel) const override;
371 
372 
381  void getTiming(Teuchos::ParameterList& timingParameterList) const override;
382 
383 
391  std::string name() const override;
392 
393  private:
394 
408  void loadA(EPhase current_phase);
409 
410  protected:
411 
418  void setNnzLU(size_t nnz){ status_.lu_nnz_ = nnz; }
419 
421  Teuchos::RCP<const MatrixAdapter<Matrix> > matrixA_;
422 
429 
430 
432  Teuchos::RCP<Vector> multiVecX_;
433 
439  Teuchos::RCP<const Vector> multiVecB_;
440 
442  global_size_type globalNumRows_;
443 
445  global_size_type globalNumCols_;
446 
448  global_size_type globalNumNonZeros_;
449 
451  global_size_type rowIndexBase_;
452 
454  global_size_type columnIndexBase_;
455 
457  mutable Status status_;
458 
460  Control control_;
461 
463  mutable Timers timers_;
464 
465 
466  /* Useful MPI vars */
467 
469  int rank_;
470 
472  bool root_;
473 
475  int nprocs_;
476 
477  }; // End class Amesos2::SolverCore
478 
479 
480 } // end namespace Amesos2
481 
482 #endif // AMESOS2_SOLVERCORE_DECL_HPP
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Return a const parameter list of all of the valid parameters that this-&gt;setParameterList(...) will accept.
Definition: Amesos2_SolverCore_def.hpp:524
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
Set or update internal variables and solver options.
Definition: Amesos2_SolverCore_decl.hpp:294
void solve() override
Solves (or )
Definition: Amesos2_SolverCore_def.hpp:139
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers...
Definition: Amesos2_SolverCore_decl.hpp:71
Vector * getXRaw() override
Returns a raw pointer to the LHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:234
~SolverCore()
Destructor.
Definition: Amesos2_SolverCore_def.hpp:61
Holds internal status data about the owning Amesos2 solver.
Definition: Amesos2_Status.hpp:39
EPhase
Used to indicate a phase in the direct solution.
Definition: Amesos2_TypeDecl.hpp:31
global_size_type columnIndexBase_
Index base of column map of matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:454
global_size_type globalNumCols_
Number of global columns in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:445
const Teuchos::RCP< const Vector > getB() override
Returns the vector that is the RHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:240
super_type & setParameters(const Teuchos::RCP< Teuchos::ParameterList > &parameterList) override
Set/update internal variables and solver options.
Definition: Amesos2_SolverCore_def.hpp:492
bool root_
If true, then this is the root processor.
Definition: Amesos2_SolverCore_decl.hpp:472
global_size_type globalNumRows_
Number of global rows in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:442
std::string description() const override
Returns a short description of this Solver.
Definition: Amesos2_SolverCore_def.hpp:562
void setA(const Teuchos::RCP< const Matrix > a, EPhase keep_phase=CLEAN) override
Sets the matrix A of this solver.
Definition: Amesos2_SolverCore_def.hpp:449
super_type & preOrdering() override
Pre-orders the matrix A for minimal fill-in.
Definition: Amesos2_SolverCore_def.hpp:69
Control control_
Parameters for solving.
Definition: Amesos2_SolverCore_decl.hpp:460
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const override
Returns a pointer to the Teuchos::Comm communicator with this operator.
Definition: Amesos2_SolverCore_decl.hpp:329
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Definition: Amesos2_SolverCore_def.hpp:572
void printTiming(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
Prints timing information about the current solver.
Definition: Amesos2_SolverCore_def.hpp:633
int nprocs_
Number of process images in the matrix communicator.
Definition: Amesos2_SolverCore_decl.hpp:475
SolverCore(Teuchos::RCP< const Matrix > A, Teuchos::RCP< Vector > X, Teuchos::RCP< const Vector > B)
Initialize a Solver instance.
Definition: Amesos2_SolverCore_def.hpp:36
void setA(const Matrix *a, EPhase keep_phase=CLEAN) override
Sets the matrix A of this solver.
Definition: Amesos2_SolverCore_decl.hpp:226
bool matrixShapeOK() override
Returns true if the solver can handle this matrix shape.
Definition: Amesos2_SolverCore_def.hpp:435
Container class for control variables.
super_type & operator=(const solver_type *rhs)
Do not allow copying of this Solver by assignment.
super_type & symbolicFactorization() override
Performs symbolic factorization on the matrix A.
Definition: Amesos2_SolverCore_def.hpp:89
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
This is an empty stub.
Definition: Amesos2_SolverCore_decl.hpp:316
bool matrix_loaded_
Definition: Amesos2_SolverCore_decl.hpp:428
A Matrix adapter interface for Amesos2.
Definition: Amesos2_MatrixAdapter_decl.hpp:42
void setB(const Teuchos::RCP< const Vector > b) override
Sets the RHS vector B.
Definition: Amesos2_SolverCore_decl.hpp:236
void loadA(EPhase current_phase)
Refresh this solver&#39;s internal data about A.
Definition: Amesos2_SolverCore_def.hpp:733
std::string name() const override
Return the name of this solver.
Definition: Amesos2_SolverCore_def.hpp:725
Interface to Amesos2 solver objects.
Definition: Amesos2_Solver_decl.hpp:44
void setB(const Vector *b) override
Sets the RHS vector B using a raw pointer.
Definition: Amesos2_SolverCore_decl.hpp:238
void setX(Vector *x) override
Sets the LHS vector X using a raw pointer.
Definition: Amesos2_SolverCore_decl.hpp:230
void getTiming(Teuchos::ParameterList &timingParameterList) const override
Extracts timing information from the current solver.
Definition: Amesos2_SolverCore_def.hpp:715
const Vector * getBRaw() override
Returns a raw pointer to the RHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:242
global_size_type globalNumNonZeros_
Number of global non-zero values in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:448
void setNnzLU(size_t nnz)
Set the number of non-zero values in the and factors.
Definition: Amesos2_SolverCore_decl.hpp:418
Enum and other types declarations for Amesos2.
void setX(const Teuchos::RCP< Vector > x) override
Sets the LHS vector X.
Definition: Amesos2_SolverCore_decl.hpp:228
int rank_
The MPI rank of this image.
Definition: Amesos2_SolverCore_decl.hpp:469
super_type & numericFactorization() override
Performs numeric factorization on the matrix A.
Definition: Amesos2_SolverCore_def.hpp:114
Teuchos::RCP< const Vector > multiVecB_
The RHS vector/multi-vector.
Definition: Amesos2_SolverCore_decl.hpp:439
Container class for Timers used with the Amesos2::Solver class.
Timers timers_
Various timing statistics.
Definition: Amesos2_SolverCore_decl.hpp:463
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
This is a empty stub.
Definition: Amesos2_SolverCore_decl.hpp:305
size_t lu_nnz_
The number of non-zeros in the factors.
Definition: Amesos2_Status.hpp:115
Container class for status variables.
Teuchos::RCP< Vector > multiVecX_
The LHS vector/multi-vector.
Definition: Amesos2_SolverCore_decl.hpp:432
Teuchos::RCP< const MatrixAdapter< Matrix > > matrixA_
The LHS operator.
Definition: Amesos2_SolverCore_decl.hpp:421
global_size_type rowIndexBase_
Index base of rowmap of matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:451
Status & getStatus() const override
Returns a reference to this solver&#39;s internal status object.
Definition: Amesos2_SolverCore_decl.hpp:336
const Teuchos::RCP< Vector > getX() override
Returns the vector that is the LHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:232
Status status_
Holds status information about a solver.
Definition: Amesos2_SolverCore_decl.hpp:457