Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Amesos2_SolverCore_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Amesos2: Templated Direct Sparse Solver Package
6 // Copyright 2011 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ***********************************************************************
41 //
42 // @HEADER
43 
52 #ifndef AMESOS2_SOLVERCORE_DECL_HPP
53 #define AMESOS2_SOLVERCORE_DECL_HPP
54 
55 #include <string>
56 
57 #include <Teuchos_RCP.hpp>
58 #include <Teuchos_ParameterList.hpp>
59 
60 #include "Amesos2_Solver.hpp"
61 #include "Amesos2_MatrixTraits.hpp"
62 #include "Amesos2_MatrixAdapter_decl.hpp"
63 #include "Amesos2_MultiVecAdapter_decl.hpp"
64 #include "Amesos2_TypeDecl.hpp"
65 
66 #include "Amesos2_Control.hpp"
67 #include "Amesos2_Status.hpp"
68 #include "Amesos2_Timers.hpp"
69 
70 namespace Amesos2 {
71 
72 
73  /* This is the base class to be used in a *statically* polymorphic
74  * way. E.g. for the Superlu solver:
75  *
76  * In Amesos2_Superlu.hpp:
77  * class Superlu : SolverCore<Superlu> { ... }
78  *
79  * Each concrete solver will implement several private sub-functions
80  * that will be called within the common code for each function.
81  */
82 
102  template <template <class,class> class ConcreteSolver,
103  class Matrix,
104  class Vector >
105  class SolverCore : public Amesos2::Solver<Matrix,Vector>
106  {
107  public:
108 
109  // Grant public access to contained types
112  typedef ConcreteSolver<Matrix,Vector> solver_type;
113  typedef Matrix matrix_type;
114  typedef Vector vector_type;
115  typedef typename MatrixAdapter<matrix_type>::scalar_t scalar_type;
116  typedef typename MatrixAdapter<matrix_type>::local_ordinal_t local_ordinal_type;
117  typedef typename MatrixAdapter<matrix_type>::global_ordinal_t global_ordinal_type;
118  typedef typename MatrixAdapter<matrix_type>::global_size_t global_size_type;
119  typedef typename MatrixAdapter<matrix_type>::node_t node_type;
120 
122 
123 
134  SolverCore( Teuchos::RCP<const Matrix> A,
135  Teuchos::RCP<Vector> X,
136  Teuchos::RCP<const Vector> B );
137 
138 
140  ~SolverCore( );
141 
142 
144  SolverCore(const solver_type& rhs);
145 
146 
148  super_type& operator=(const solver_type* rhs);
149 
151 
153 
154 
166  super_type& preOrdering() override;
167 
168 
186  super_type& symbolicFactorization() override;
187 
188 
212  super_type& numericFactorization() override;
213 
214 
234  void solve() override;
235  void solve(const Teuchos::Ptr<Vector> X, const Teuchos::Ptr<const Vector> B) const override;
236  void solve(Vector* X, const Vector* B) const override;
237 
238  int solve_ir(const Teuchos::Ptr<Vector> X, const Teuchos::Ptr<const Vector> B,
239  const int maxNumIters, const bool verbose) const;
240  int solve_ir(Vector* X, const Vector* B,
241  const int maxNumIters, const bool verbose) const;
242  int solve_ir(const int maxNumIters, const bool verbose);
243 
245 
246 
256  bool matrixShapeOK() override;
257 
258  void setA( const Teuchos::RCP<const Matrix> a, EPhase keep_phase = CLEAN ) override;
259 
260  void setA( const Matrix* a, EPhase keep_phase = CLEAN ) override { setA(Teuchos::rcp(a), keep_phase); }
261 
262  void setX(const Teuchos::RCP<Vector> x) override { multiVecX_ = x; }
263 
264  void setX(Vector* x) override { multiVecX_ = Teuchos::rcp(x); }
265 
266  const Teuchos::RCP<Vector> getX() override { return( multiVecX_ ); }
267 
268  Vector* getXRaw() override { return multiVecX_.getRawPtr(); }
269 
270  void setB(const Teuchos::RCP<const Vector> b) override { multiVecB_ = b; }
271 
272  void setB(const Vector* b) override { multiVecB_ = Teuchos::rcp(b); }
273 
274  const Teuchos::RCP<const Vector> getB() override { return( multiVecB_ ); }
275 
276  const Vector* getBRaw() override { return multiVecB_.getRawPtr(); }
277 
278 
280 
281 
303  super_type& setParameters(
304  const Teuchos::RCP<Teuchos::ParameterList> & parameterList ) override;
305 
306 
316  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const override;
317 
318 
328  void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & parameterList)
329  {
330  setParameters(parameterList);
331  }
332 
333 
339  Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList()
340  {
341  return Teuchos::null;
342  }
343 
344 
350  Teuchos::RCP<Teuchos::ParameterList> unsetParameterList()
351  {
352  return Teuchos::null;
353  }
354 
355 
357 
358 
360 
361 
363  Teuchos::RCP<const Teuchos::Comm<int> > getComm() const override
364  {
365  return matrixA_->getComm();
366  }
367 
368 
370  inline Status& getStatus() const override { return( status_ ); }
371 
372 
374 
375 
377  std::string description() const override;
378 
379 
382  void describe(Teuchos::FancyOStream &out,
383  const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override;
384 
385 
402  void printTiming(
403  Teuchos::FancyOStream &out,
404  const Teuchos::EVerbosityLevel verbLevel) const override;
405 
406 
415  void getTiming(Teuchos::ParameterList& timingParameterList) const override;
416 
417 
425  std::string name() const override;
426 
427  private:
428 
442  void loadA(EPhase current_phase);
443 
444  protected:
445 
452  void setNnzLU(size_t nnz){ status_.lu_nnz_ = nnz; }
453 
455  Teuchos::RCP<const MatrixAdapter<Matrix> > matrixA_;
456 
463 
464 
466  Teuchos::RCP<Vector> multiVecX_;
467 
473  Teuchos::RCP<const Vector> multiVecB_;
474 
476  global_size_type globalNumRows_;
477 
479  global_size_type globalNumCols_;
480 
482  global_size_type globalNumNonZeros_;
483 
485  global_size_type rowIndexBase_;
486 
488  global_size_type columnIndexBase_;
489 
491  mutable Status status_;
492 
494  Control control_;
495 
497  mutable Timers timers_;
498 
499 
500  /* Useful MPI vars */
501 
503  int rank_;
504 
506  bool root_;
507 
509  int nprocs_;
510 
511  }; // End class Amesos2::SolverCore
512 
513 
514 } // end namespace Amesos2
515 
516 #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:558
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
Set or update internal variables and solver options.
Definition: Amesos2_SolverCore_decl.hpp:328
void solve() override
Solves (or )
Definition: Amesos2_SolverCore_def.hpp:173
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers...
Definition: Amesos2_SolverCore_decl.hpp:105
Vector * getXRaw() override
Returns a raw pointer to the LHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:268
~SolverCore()
Destructor.
Definition: Amesos2_SolverCore_def.hpp:95
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:488
global_size_type globalNumCols_
Number of global columns in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:479
const Teuchos::RCP< const Vector > getB() override
Returns the vector that is the RHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:274
super_type & setParameters(const Teuchos::RCP< Teuchos::ParameterList > &parameterList) override
Set/update internal variables and solver options.
Definition: Amesos2_SolverCore_def.hpp:526
bool root_
If true, then this is the root processor.
Definition: Amesos2_SolverCore_decl.hpp:506
global_size_type globalNumRows_
Number of global rows in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:476
std::string description() const override
Returns a short description of this Solver.
Definition: Amesos2_SolverCore_def.hpp:596
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:483
super_type & preOrdering() override
Pre-orders the matrix A for minimal fill-in.
Definition: Amesos2_SolverCore_def.hpp:103
Control control_
Parameters for solving.
Definition: Amesos2_SolverCore_decl.hpp:494
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:363
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Definition: Amesos2_SolverCore_def.hpp:606
void printTiming(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
Prints timing information about the current solver.
Definition: Amesos2_SolverCore_def.hpp:667
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:70
void setA(const Matrix *a, EPhase keep_phase=CLEAN) override
Sets the matrix A of this solver.
Definition: Amesos2_SolverCore_decl.hpp:260
bool matrixShapeOK() override
Returns true if the solver can handle this matrix shape.
Definition: Amesos2_SolverCore_def.hpp:469
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:123
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
This is an empty stub.
Definition: Amesos2_SolverCore_decl.hpp:350
bool matrix_loaded_
Definition: Amesos2_SolverCore_decl.hpp:462
A Matrix adapter interface for Amesos2.
Definition: Amesos2_MatrixAdapter_decl.hpp:76
void setB(const Teuchos::RCP< const Vector > b) override
Sets the RHS vector B.
Definition: Amesos2_SolverCore_decl.hpp:270
void loadA(EPhase current_phase)
Refresh this solver&#39;s internal data about A.
Definition: Amesos2_SolverCore_def.hpp:767
std::string name() const override
Return the name of this solver.
Definition: Amesos2_SolverCore_def.hpp:759
Interface to Amesos2 solver objects.
Definition: Amesos2_Solver_decl.hpp:78
void setB(const Vector *b) override
Sets the RHS vector B using a raw pointer.
Definition: Amesos2_SolverCore_decl.hpp:272
void setX(Vector *x) override
Sets the LHS vector X using a raw pointer.
Definition: Amesos2_SolverCore_decl.hpp:264
void getTiming(Teuchos::ParameterList &timingParameterList) const override
Extracts timing information from the current solver.
Definition: Amesos2_SolverCore_def.hpp:749
const Vector * getBRaw() override
Returns a raw pointer to the RHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:276
global_size_type globalNumNonZeros_
Number of global non-zero values in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:482
void setNnzLU(size_t nnz)
Set the number of non-zero values in the and factors.
Definition: Amesos2_SolverCore_decl.hpp:452
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:262
int rank_
The MPI rank of this image.
Definition: Amesos2_SolverCore_decl.hpp:503
super_type & numericFactorization() override
Performs numeric factorization on the matrix A.
Definition: Amesos2_SolverCore_def.hpp:148
Teuchos::RCP< const Vector > multiVecB_
The RHS vector/multi-vector.
Definition: Amesos2_SolverCore_decl.hpp:473
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:339
size_t lu_nnz_
The number of non-zeros in the factors.
Definition: Amesos2_Status.hpp:149
Container class for status variables.
Teuchos::RCP< Vector > multiVecX_
The LHS vector/multi-vector.
Definition: Amesos2_SolverCore_decl.hpp:466
Teuchos::RCP< const MatrixAdapter< Matrix > > matrixA_
The LHS operator.
Definition: Amesos2_SolverCore_decl.hpp:455
global_size_type rowIndexBase_
Index base of rowmap of matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:485
Status & getStatus() const override
Returns a reference to this solver&#39;s internal status object.
Definition: Amesos2_SolverCore_decl.hpp:370
const Teuchos::RCP< Vector > getX() override
Returns the vector that is the LHS of the linear system.
Definition: Amesos2_SolverCore_decl.hpp:266
Status status_
Holds status information about a solver.
Definition: Amesos2_SolverCore_decl.hpp:491