Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Amesos_Umfpack.h
Go to the documentation of this file.
1 
2 // @HEADER
3 // ***********************************************************************
4 //
5 // Amesos: Direct Sparse Solver Package
6 // Copyright (2004) Sandia Corporation
7 //
8 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9 // license for use of this work by or on behalf of the U.S. Government.
10 //
11 // This library is free software; you can redistribute it and/or modify
12 // it under the terms of the GNU Lesser General Public License as
13 // published by the Free Software Foundation; either version 2.1 of the
14 // License, or (at your option) any later version.
15 //
16 // This library is distributed in the hope that it will be useful, but
17 // WITHOUT ANY WARRANTY; without even the implied warranty of
18 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 // Lesser General Public License for more details.
20 //
21 // You should have received a copy of the GNU Lesser General Public
22 // License along with this library; if not, write to the Free Software
23 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
24 // USA
25 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
26 //
27 // ***********************************************************************
28 // @HEADER
29 
40 #ifndef AMESOS_UMFPACK_H
41 #define AMESOS_UMFPACK_H
42 
43 #include "Amesos_ConfigDefs.h"
44 #include "Amesos_BaseSolver.h"
45 #include "Amesos_BaseSolver.h"
46 #include "Amesos_NoCopiable.h"
47 #include "Amesos_Utils.h"
48 #include "Amesos_Time.h"
49 #include "Amesos_Status.h"
50 #include "Amesos_Control.h"
51 #include "Epetra_CrsMatrix.h"
52 #include "Epetra_LinearProblem.h"
53 #include "Epetra_Time.h"
54 #include "Epetra_Import.h"
55 #include "Teuchos_RCP.hpp"
56 #ifdef EPETRA_MPI
57 #include "Epetra_MpiComm.h"
58 #else
59 #include "Epetra_Comm.h"
60 #endif
61 
63 
71  private Amesos_Time,
72  private Amesos_NoCopiable,
73  private Amesos_Utils,
74  private Amesos_Control,
75  private Amesos_Status
76 {
77 public:
78 
80 
88  Amesos_Umfpack( const Epetra_LinearProblem& LinearProblem );
89 
91 
93  ~Amesos_Umfpack(void);
94 
96 
97 
99 
100  int NumericFactorization();
101 
102  int Solve();
103 
105 
106 
107  const Epetra_LinearProblem *GetProblem() const { return(Problem_); };
108 
110 
113  bool MatrixShapeOK() const ;
114 
115  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
116 
117  bool UseTranspose() const {return(UseTranspose_);};
118 
119  const Epetra_Comm & Comm() const {return(GetProblem()->GetOperator()->Comm());};
120 
122  /* Rcond is an estimate of the reciprocal of the condition number of the
123  matrix at the time of the most recent call to NumericFactorization()
124  Rcond = min(abs(diag))/max(abs(diag)) see Umfpack documentatoin
125  for details.
126  */
127  double GetRcond() const ;
128 
130 
133 
136 
138  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
139 
141  void PrintTiming() const;
142 
144  void PrintStatus() const;
145 
147  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming(TimingParameterList); }
148 
149 private:
150 
152 
153 
156  {
157  return(dynamic_cast<Epetra_RowMatrix *>(Problem_->GetOperator()));
158  }
159 
161  int ConvertToSerial(const bool FirstTime);
162 
163  /*
164  ConvertToUmfpackCRS - Convert matirx to form expected by Umfpack: Ai, Ap, Aval
165  Preconditions:
166  numentries_, NumGloalElements_ and SerialMatrix_ must be set.
167  Postconditions:
168  Ai, Ap, and Aval are resized and populated with a compresses row storage
169  version of the input matrix A.
170  */
171  int ConvertToUmfpackCRS();
172 
173  /*
174  PerformSymbolicFactorization - Call Umfpack to perform symbolic factorization
175  Preconditions:
176  IsLocal must be set to 1 if the input matrix is entirely stored on process 0
177  Ap, Ai and Aval are a compressed row storage version of the input matrix A.
178  Postconditions:
179  Symbolic points to an UMFPACK internal opaque object containing the
180  symbolic factorization and accompanying information.
181  SymbolicFactorizationOK_ = true;
182  Note: All action is performed on process 0
183  */
184 
186 
187  /*
188  PerformNumericFactorization - Call Umfpack to perform numeric factorization
189  Preconditions:
190  IsLocal must be set
191  Ap, Ai and Aval are a compressed row storage version of the input matrix A.
192  Symbolic must be set
193  Postconditions:
194  Numeric points to an UMFPACK internal opaque object containing the
195  numeric factorization and accompanying information.
196  NumericFactorizationOK_ = true;
197  Note: All action is performed on process 0
198  */
200 
201  inline const Epetra_Import& Importer() const
202  {
203  return(*(ImportToSerial_.get()));
204  }
205 
206  inline const Epetra_Map& SerialMap() const
207  {
208  return(*(SerialMap_.get()));
209  }
210 
211  inline const Epetra_CrsMatrix& SerialCrsMatrix() const
212  {
213  return(*(SerialCrsMatrixA_.get()));
214  }
215 
217  {
218  return(*(SerialCrsMatrixA_.get()));
219  }
220 
221  // @}
222 
224  void *Symbolic;
226  void *Numeric;
227 
229  std::vector <int> Ap;
230  std::vector <int> Ai;
231  std::vector <double> Aval;
232 
234  int IsLocal_;
239 
243  /* If IsLocal==1 - Points to the original matrix
244  * If IsLocal==0 - Points to SerialCrsMatrixA
245  */
247 
249 
255  mutable double Rcond_;
256  // True if Rcond_ is the same on all processes
257  mutable bool RcondValidOnAllProcs_;
260 
264 
265 }; // class Amesos_Umfpack
266 #endif /* AMESOS_UMFPACK_H */
std::vector< int > Ap
Ap, Ai, Aval form the compressed row storage used by Umfpack.
int numentries_
Number of non-zero entries in Problem_-&gt;GetOperator()
int NumSymbolicFact_
Number of symbolic factorization phases.
Definition: Amesos_Status.h:67
Amesos_Control: Container for some control variables.
void GetTiming(Teuchos::ParameterList &TimingParameterList) const
Extracts timing information from the current solver and places it in the parameter list...
int PerformSymbolicFactorization()
int MtxConvTime_
Quick access pointers to internal timer data.
int SetUseTranspose(bool UseTranspose_in)
If set true, X will be set to the solution of AT X = B (not A X = B)
const Epetra_LinearProblem * Problem_
Pointer to the linear problem to solve.
void * Numeric
Umfpack internal opaque object.
const Epetra_CrsMatrix & SerialCrsMatrix() const
int NumSolve() const
Returns the number of solves performed by this object.
Teuchos::RCP< Epetra_CrsMatrix > SerialCrsMatrixA_
T * get() const
int NumNumericFact_
Number of numeric factorization phases.
Definition: Amesos_Status.h:69
~Amesos_Umfpack(void)
Amesos_Umfpack Destructor.
int ConvertToSerial(const bool FirstTime)
Converts matrix to a serial Epetra_CrsMatrix.
void PrintStatus() const
Prints information about the factorization and solution phases.
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
const Epetra_LinearProblem * GetProblem() const
Returns the Epetra_LinearProblem.
int NumSolve_
Number of solves.
Definition: Amesos_Status.h:71
int NumGlobalElements_
Number of rows and columns in the Problem_-&gt;GetOperator()
int SymbolicFactorization()
Performs SymbolicFactorization on the matrix A.
void GetTiming(Teuchos::ParameterList &list) const
Load up the current timing information into the parameter list.
Definition: Amesos_Time.h:124
Epetra_RowMatrix * Matrix()
Returns a pointer to the linear system matrix.
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:20
Teuchos::RCP< Epetra_Import > ImportToSerial_
Importer from distributed to serial (all rows on process 0).
Amesos_Time: Container for timing information.
Definition: Amesos_Time.h:50
Epetra_RowMatrix * SerialMatrix_
Points to a Serial Copy of A.
const Epetra_Map & SerialMap() const
std::vector< double > Aval
Amesos_Umfpack(const Epetra_LinearProblem &LinearProblem)
Amesos_Umfpack Constructor.
std::vector< int > Ai
bool MatrixShapeOK() const
Returns true if UMFPACK can handle this matrix shape.
double Rcond_
Reciprocal condition number estimate.
int NumSymbolicFact() const
Returns the number of symbolic factorizations performed by this object.
int Solve()
Solves A X = B (or AT x = B)
int PerformNumericFactorization()
int NumericFactorization()
Performs NumericFactorization on the matrix A.
bool UseTranspose_
If true, solve the problem with the transpose.
double GetRcond() const
Returns an estimate of the reciprocal of the condition number.
const Epetra_Import & Importer() const
Epetra_Operator * GetOperator() const
int NumNumericFact() const
Returns the number of numeric factorizations performed by this object.
int IsLocal_
1 if Problem_-&gt;GetOperator() is stored entirely on process 0
void PrintTiming() const
Prints timing information.
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
bool RcondValidOnAllProcs_
int SetParameters(Teuchos::ParameterList &ParameterList)
Updates internal variables.
Class Amesos_Umfpack: An object-oriented wrapper for UMFPACK.
Amesos_NoCopiable: Simple class to prevent the usage of copy constructor and operator =...
Epetra_CrsMatrix & SerialCrsMatrix()
bool UseTranspose() const
Returns the current UseTranspose setting.
Teuchos::RCP< Epetra_Map > SerialMap_
Points to a Serial Map (unused if IsLocal == 1 )
Amesos_Utils: Collections of basic utilities.
Definition: Amesos_Utils.h:19
void * Symbolic
Umfpack internal opaque object.