Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Amesos_Superlu.h
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Amesos: Direct Sparse Solver Package
5 // Copyright (2004) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // This library is free software; you can redistribute it and/or modify
11 // it under the terms of the GNU Lesser General Public License as
12 // published by the Free Software Foundation; either version 2.1 of the
13 // License, or (at your option) any later version.
14 //
15 // This library is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 // Lesser General Public License for more details.
19 //
20 // You should have received a copy of the GNU Lesser General Public
21 // License along with this library; if not, write to the Free Software
22 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
23 // USA
24 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
25 //
26 // ***********************************************************************
27 // @HEADER
28 
29 //
30 // Coding tasks:
31 // 1) Create the dense matrices in Solve() DONE
32 // 2) Make the call to dgsvx() in Factor() dONE
33 // 3) Make the call to dgsvx() in Solve() DONE
34 //
35 // Later coding tasks:
36 // 0) Factor called twice
37 // 1) Refactor()
38 // 2) Parameter list
39 // 3) Transpose - DONE
40 // 4) Destructor - In particular, need to call the SuperLU_FREE routines.
41 // 5) Coments - especially in Amesos_Superlu.h
42 //
43 // SymbolicFactorization() performs no action other than making sure that Factorization
44 // s performed
45 // NumericFactorization() performs the factorization but no solve (right hand side is
46 // is set to 0 vectors) reuses factors only if ReuseFactorization_ is set.
47 // If FactorizationOK_() && ReuseSymbolic_
48 // ReFactor()
49 // else
50 // Factor()
51 //
52 // Solve()
53 //
54 // Factor() does everything from scratch:
55 // Redistributes the data if necessary
56 // Deletes any data structures left over from the previous call to Factor()
57 // Copies the data into the format that SuperLU wants it
58 // Calls dgssvx to factor the matrix with factor set to true
59 // ReFactor()
60 // Redistributes the data if necessary
61 // - Attempting to check to make sure that the non-zero structure is unchanged
62 // Copies the data into the format that SuperLU already has it
63 // FIRST PASS - assert( false )
64 //
65 
66 #ifndef AMESOS_SUPERLU_H
67 #define AMESOS_SUPERLU_H
68 
69 #include "Amesos_ConfigDefs.h"
70 #include "Amesos_BaseSolver.h"
71 #include "Amesos_NoCopiable.h"
72 #include "Amesos_Utils.h"
73 #include "Amesos_Time.h"
74 #include "Amesos_Status.h"
75 #include "Amesos_Control.h"
76 #include "Teuchos_RCP.hpp"
77 
78 class SLUData;
79 class Epetra_Comm;
80 class Epetra_CrsMatrix;
82 
84 
93  private Amesos_Time,
94  private Amesos_NoCopiable,
95  private Amesos_Utils,
96  private Amesos_Control,
97  private Amesos_Status {
98 
99 public:
100 
102 
110  Amesos_Superlu(const Epetra_LinearProblem& LinearProblem );
111 
113  ~Amesos_Superlu();
114 
116 
117 
118  int SymbolicFactorization();
119 
120  int NumericFactorization();
121 
122  int Solve();
123 
125 
126 
127  const Epetra_LinearProblem *GetProblem() const { return(Problem_); };
128 
129  bool MatrixShapeOK() const ;
130 
131  int SetUseTranspose (bool useTheTranspose) {
132  UseTranspose_ = useTheTranspose; return(0);
133  }
134 
135  bool UseTranspose() const {return(UseTranspose_);};
136 
137  const Epetra_Comm& Comm() const {return(GetProblem()->GetOperator()->Comm());};
138 
140 
143 
146 
148  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
149 
151  void PrintTiming() const;
152 
154  void PrintStatus() const;
155 
157  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming(TimingParameterList); }
158 
159 private:
160 
162 
163 
165  // Note: this method is delicate!
166  const Epetra_Map& SerialMap() const
167  {
168  return(*(SerialMap_.get()));
169  }
170 
172  // Note: this method is delicate!
174  {
175  return(*(ImportToSerial_.get()));
176  }
177 
179  int Factor();
181  int ReFactor();
182 
184  int ConvertToSerial();
185 
187  // Note: All action is performed on process 0
189 
191 
194  std::vector<double> berr_;
195  std::vector<double> ferr_;
196  std::vector<int> perm_r_;
197  std::vector<int> perm_c_;
198  std::vector<int> etree_;
199  std::vector<double> R_;
200  std::vector<double> C_;
201  char equed_;
202  // no idea of the following.
203  double* DummyArray;
204 
206  std::vector <int> Ap_;
208  std::vector <int> Ai_;
210  std::vector <double> Aval_;
212  long long NumGlobalRows_;
222  int iam_;
238 
239 }; // End of class Amesos_Superlu
240 #endif /* AMESOS_SUPERLU_H */
int NumSymbolicFact_
Number of symbolic factorization phases.
Definition: Amesos_Status.h:67
void PrintStatus() const
Prints status information.
Amesos_Control: Container for some control variables.
std::vector< int > etree_
Teuchos::RCP< Epetra_Map > SerialMap_
Contains a map with all elements assigned to processor 0.
std::vector< double > Aval_
int SetParameters(Teuchos::ParameterList &ParameterList)
Updates internal variables.
Epetra_RowMatrix * SerialMatrix_
For parallel runs, stores the matrix defined on SerialMap_.
int Solve()
Solves A X = B (or AT x = B)
std::vector< int > perm_c_
std::vector< int > Ap_
stores the matrix in SuperLU format.
std::vector< double > C_
const Epetra_Map & SerialMap() const
Returns a reference to the serial map.
T * get() const
const Epetra_Import & ImportToSerial() const
Returns a reference to the importer.
int Factor()
Factors the matrix, no previous factorization available.
int PerformNumericFactorization()
PerformNumericFactorization - Call Superlu to perform numeric factorization.
int MtxConvTime_
Quick access pointer to internal timing data.
bool UseTranspose_
If true, solve the linear system with the transpose of the matrix.
Amesos_Superlu: Amesos interface to Xioye Li&#39;s SuperLU 3.0 serial code.
long long NumGlobalNonzeros_
Global number of nonzeros in the matrix.
int NumNumericFact_
Number of numeric factorization phases.
Definition: Amesos_Status.h:69
int ReFactor()
Re-factors the matrix.
bool UseTranspose() const
Returns the current UseTranspose setting.
int NumSolve_
Number of solves.
Definition: Amesos_Status.h:71
void GetTiming(Teuchos::ParameterList &TimingParameterList) const
Extracts timing information from the current solver and places it in the parameter list...
bool MatrixShapeOK() const
Returns true if the solver can handle this matrix shape.
Amesos_Superlu(const Epetra_LinearProblem &LinearProblem)
Amesos_Superlu Constructor.
const Epetra_LinearProblem * GetProblem() const
Returns the Epetra_LinearProblem.
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
int NumSolve() const
Returns the number of solves performed by this object.
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:20
std::vector< double > ferr_
Amesos_Time: Container for timing information.
Definition: Amesos_Time.h:50
int SetUseTranspose(bool useTheTranspose)
If set true, X will be set to the solution of AT X = B (not A X = B)
Epetra_RowMatrix * RowMatrixA_
Pointer to the linear system matrix.
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
double * DummyArray
stores the matrix in SuperLU format.
int ConvertToSerial()
Sets up the matrix on processor 0.
int NumNumericFact() const
Returns the number of numeric factorizations performed by this object.
std::vector< double > R_
Teuchos::RCP< Epetra_CrsMatrix > SerialCrsMatrixA_
Contains a matrix with all rows assigned to processor 0.
Teuchos::RCP< Epetra_Import > ImportToSerial_
Importer from distributed to SerialMap_.
SLUData * data_
Main structure for SuperLU.
~Amesos_Superlu()
Amesos_Superlu Destructor.
std::vector< double > berr_
std::vector< int > perm_r_
bool FactorizationOK_
If true, the factorization has been successfully computed.
int iam_
Process number (i.e. Comm().MyPID()
void PrintTiming() const
Prints timing information.
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
long long NumGlobalRows_
Global size of the matrix.
Amesos_NoCopiable: Simple class to prevent the usage of copy constructor and operator =...
std::vector< int > Ai_
stores the matrix in SuperLU format.
const Epetra_LinearProblem * Problem_
Pointer to the user&#39;s defined linear problem.
Amesos_Utils: Collections of basic utilities.
Definition: Amesos_Utils.h:19
int NumericFactorization()
Performs NumericFactorization on the matrix A.
int NumSymbolicFact() const
Returns the number of symbolic factorizations performed by this object.