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 #if defined(Amesos_SHOW_DEPRECATED_WARNINGS)
44 #ifdef __GNUC__
45 #warning "The Amesos package is deprecated"
46 #endif
47 #endif
48 
49 #include "Amesos_ConfigDefs.h"
50 #include "Amesos_BaseSolver.h"
51 #include "Amesos_BaseSolver.h"
52 #include "Amesos_NoCopiable.h"
53 #include "Amesos_Utils.h"
54 #include "Amesos_Time.h"
55 #include "Amesos_Status.h"
56 #include "Amesos_Control.h"
57 #include "Epetra_CrsMatrix.h"
58 #include "Epetra_LinearProblem.h"
59 #include "Epetra_Time.h"
60 #include "Epetra_Import.h"
61 #include "Teuchos_RCP.hpp"
62 #ifdef EPETRA_MPI
63 #include "Epetra_MpiComm.h"
64 #else
65 #include "Epetra_Comm.h"
66 #endif
67 
69 
77  private Amesos_Time,
78  private Amesos_NoCopiable,
79  private Amesos_Utils,
80  private Amesos_Control,
81  private Amesos_Status
82 {
83 public:
84 
86 
94  Amesos_Umfpack( const Epetra_LinearProblem& LinearProblem );
95 
97 
99  ~Amesos_Umfpack(void);
100 
102 
103 
104  int SymbolicFactorization();
105 
106  int NumericFactorization();
107 
108  int Solve();
109 
111 
112 
113  const Epetra_LinearProblem *GetProblem() const { return(Problem_); };
114 
116 
119  bool MatrixShapeOK() const ;
120 
121  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
122 
123  bool UseTranspose() const {return(UseTranspose_);};
124 
125  const Epetra_Comm & Comm() const {return(GetProblem()->GetOperator()->Comm());};
126 
128  /* Rcond is an estimate of the reciprocal of the condition number of the
129  matrix at the time of the most recent call to NumericFactorization()
130  Rcond = min(abs(diag))/max(abs(diag)) see Umfpack documentatoin
131  for details.
132  */
133  double GetRcond() const ;
134 
136 
139 
142 
144  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
145 
147  void PrintTiming() const;
148 
150  void PrintStatus() const;
151 
153  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming(TimingParameterList); }
154 
155 private:
156 
158 
159 
162  {
163  return(dynamic_cast<Epetra_RowMatrix *>(Problem_->GetOperator()));
164  }
165 
167  int ConvertToSerial(const bool FirstTime);
168 
169  /*
170  ConvertToUmfpackCRS - Convert matirx to form expected by Umfpack: Ai, Ap, Aval
171  Preconditions:
172  numentries_, NumGloalElements_ and SerialMatrix_ must be set.
173  Postconditions:
174  Ai, Ap, and Aval are resized and populated with a compresses row storage
175  version of the input matrix A.
176  */
177  int ConvertToUmfpackCRS();
178 
179  /*
180  PerformSymbolicFactorization - Call Umfpack to perform symbolic factorization
181  Preconditions:
182  IsLocal must be set to 1 if the input matrix is entirely stored on process 0
183  Ap, Ai and Aval are a compressed row storage version of the input matrix A.
184  Postconditions:
185  Symbolic points to an UMFPACK internal opaque object containing the
186  symbolic factorization and accompanying information.
187  SymbolicFactorizationOK_ = true;
188  Note: All action is performed on process 0
189  */
190 
192 
193  /*
194  PerformNumericFactorization - Call Umfpack to perform numeric factorization
195  Preconditions:
196  IsLocal must be set
197  Ap, Ai and Aval are a compressed row storage version of the input matrix A.
198  Symbolic must be set
199  Postconditions:
200  Numeric points to an UMFPACK internal opaque object containing the
201  numeric factorization and accompanying information.
202  NumericFactorizationOK_ = true;
203  Note: All action is performed on process 0
204  */
206 
207  inline const Epetra_Import& Importer() const
208  {
209  return(*(ImportToSerial_.get()));
210  }
211 
212  inline const Epetra_Map& SerialMap() const
213  {
214  return(*(SerialMap_.get()));
215  }
216 
217  inline const Epetra_CrsMatrix& SerialCrsMatrix() const
218  {
219  return(*(SerialCrsMatrixA_.get()));
220  }
221 
223  {
224  return(*(SerialCrsMatrixA_.get()));
225  }
226 
227  // @}
228 
230  void *Symbolic;
232  void *Numeric;
233 
235  std::vector <int> Ap;
236  std::vector <int> Ai;
237  std::vector <double> Aval;
238 
240  int IsLocal_;
245 
249  /* If IsLocal==1 - Points to the original matrix
250  * If IsLocal==0 - Points to SerialCrsMatrixA
251  */
253 
255 
261  mutable double Rcond_;
262  // True if Rcond_ is the same on all processes
263  mutable bool RcondValidOnAllProcs_;
266 
270 
271 }; // class Amesos_Umfpack
272 #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:73
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:75
~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:77
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:130
Epetra_RowMatrix * Matrix()
Returns a pointer to the linear system matrix.
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:26
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:56
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:25
void * Symbolic
Umfpack internal opaque object.