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 
167 
168 
187 
188 
213 
214 
234  void solve();
235 
236 
237  void solve(const Teuchos::Ptr<Vector> X, const Teuchos::Ptr<const Vector> B) const;
238 
239 
240  void solve(Vector* X, const Vector* B) const;
241 
242 
244 
245 
255  bool matrixShapeOK();
256 
257  void setA( const Teuchos::RCP<const Matrix> a, EPhase keep_phase = CLEAN );
258 
259  void setA( const Matrix* a, EPhase keep_phase = CLEAN ){ setA(Teuchos::rcp(a), keep_phase); }
260 
261  void setX(const Teuchos::RCP<Vector> x){ multiVecX_ = x; }
262 
263  void setX(Vector* x){ multiVecX_ = Teuchos::rcp(x); }
264 
265  const Teuchos::RCP<Vector> getX(){ return( multiVecX_ ); }
266 
267  Vector* getXRaw(){ return multiVecX_.getRawPtr(); }
268 
269  void setB(const Teuchos::RCP<const Vector> b){ multiVecB_ = b; }
270 
271  void setB(const Vector* b){ multiVecB_ = Teuchos::rcp(b); }
272 
273  const Teuchos::RCP<const Vector> getB(){ return( multiVecB_ ); }
274 
275  const Vector* getBRaw(){ return multiVecB_.getRawPtr(); }
276 
277 
279 
280 
302  super_type& setParameters(
303  const Teuchos::RCP<Teuchos::ParameterList> & parameterList );
304 
305 
315  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
316 
317 
327  void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & parameterList)
328  {
329  setParameters(parameterList);
330  }
331 
332 
338  Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList()
339  {
340  return Teuchos::null;
341  }
342 
343 
349  Teuchos::RCP<Teuchos::ParameterList> unsetParameterList()
350  {
351  return Teuchos::null;
352  }
353 
354 
356 
357 
359 
360 
362  Teuchos::RCP<const Teuchos::Comm<int> > getComm() const
363  {
364  return matrixA_->getComm();
365  }
366 
367 
369  inline Status& getStatus() const { return( status_ ); }
370 
371 
373 
374 
376  std::string description() const;
377 
378 
381  void describe(Teuchos::FancyOStream &out,
382  const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const;
383 
384 
401  void printTiming(
402  Teuchos::FancyOStream &out,
403  const Teuchos::EVerbosityLevel verbLevel) const;
404 
405 
414  void getTiming(Teuchos::ParameterList& timingParameterList) const;
415 
416 
424  std::string name() const;
425 
426  private:
427 
441  void loadA(EPhase current_phase);
442 
443  protected:
444 
451  void setNnzLU(size_t nnz){ status_.lu_nnz_ = nnz; }
452 
454  Teuchos::RCP<const MatrixAdapter<Matrix> > matrixA_;
455 
462 
463 
465  Teuchos::RCP<Vector> multiVecX_;
466 
472  Teuchos::RCP<const Vector> multiVecB_;
473 
475  global_size_type globalNumRows_;
476 
478  global_size_type globalNumCols_;
479 
481  global_size_type globalNumNonZeros_;
482 
484  global_size_type rowIndexBase_;
485 
487  global_size_type columnIndexBase_;
488 
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
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
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-&gt;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&#39;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&#39;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 > &parameterList)
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