Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Amesos_Lapack.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 #ifndef AMESOS_LAPACK_H
30 #define AMESOS_LAPACK_H
31 
32 #include "Amesos_ConfigDefs.h"
33 #include "Amesos_BaseSolver.h"
34 #include "Amesos_NoCopiable.h"
35 #include "Amesos_Utils.h"
36 #include "Amesos_Time.h"
37 #include "Amesos_Status.h"
38 #include "Amesos_Control.h"
39 #include "Epetra_Comm.h"
40 #include "Epetra_Map.h"
44 #include "Epetra_CrsMatrix.h"
45 #include "Epetra_Import.h"
46 #include "Epetra_Export.h"
47 class Epetra_RowMatrix;
49 #include "Teuchos_RCP.hpp"
50 
52 
66  private Amesos_Time,
67  private Amesos_NoCopiable,
68  private Amesos_Utils,
69  private Amesos_Control,
70  private Amesos_Status
71 {
72 public:
73 
75 
83  Amesos_Lapack(const Epetra_LinearProblem& LinearProblem );
84 
86 
88  ~Amesos_Lapack(void);
89 
91 
92 
93  int SymbolicFactorization() ;
94 
95  int NumericFactorization() ;
96 
97  int Solve();
98 
100 
101 
102  const Epetra_LinearProblem *GetProblem() const { return(Problem_); };
103 
104  bool MatrixShapeOK() const;
105 
106  int SetUseTranspose(bool UseTranspose_in) {
107  UseTranspose_ = UseTranspose_in;
108  return(0);
109  }
110 
111  bool UseTranspose() const {return(UseTranspose_);};
112 
113  const Epetra_Comm & Comm() const {
114  return(GetProblem()->GetOperator()->Comm());
115  }
116 
118 
121 
123 
124 
127 
129 
139  int GEEV(Epetra_Vector& Er, Epetra_Vector& Ei);
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 
160 
161 protected:
162 
164  const Epetra_RowMatrix* Matrix() const
165  {
166  return(Problem_->GetMatrix());
167  }
168 
170 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
171  inline int NumGlobalRows() const
172  {
173  return(Matrix()->NumGlobalRows());
174  }
175 #endif
176 
177  inline long long NumGlobalRows64() const
178  {
179  return(Matrix()->NumGlobalRows64());
180  }
181 
183  inline int NumMyRows() const
184  {
185  return(Matrix()->NumMyRows());
186  }
187 
189  inline const Epetra_Map& SerialMap()
190  {
191  return(*(SerialMap_.get()));
192  }
193 
196  {
197  return(*(SerialMatrix_.get()));
198  }
199 
201  {
202  return(*(SerialCrsMatrix_.get()));
203  }
204 
207  {
208  return(*(MatrixImporter_.get()));
209  }
210 
213  {
214  return(*(RhsExporter_.get()));
215  }
216 
219  {
220  return(*(SolutionImporter_.get()));
221  }
222 
224 
227  const Epetra_MultiVector& B);
228 
231  const Epetra_MultiVector& B);
232 
234  int DistributedToSerial();
235 
237  int SerialToDense();
238 
240  int DenseToFactored();
241 
248 
257 
262 
265 
266  long long NumGlobalRows_;
268 
270 
271 }; // End of class Amesos_Lapack
272 #endif /* AMESOS_LAPACK_H */
int NumSolve() const
Returns the number of solves performed by this object.
int NumSymbolicFact_
Number of symbolic factorization phases.
Definition: Amesos_Status.h:67
Amesos_Lapack(const Epetra_LinearProblem &LinearProblem)
Amesos_Lapack Constructor.
Amesos_Control: Container for some control variables.
int SerialToDense()
Converts a serial matrix to dense format.
Epetra_SerialDenseMatrix DenseRHS_
Dense RHS.
Teuchos::RCP< Epetra_CrsMatrix > SerialCrsMatrix_
Epetra_SerialDenseMatrix DenseMatrix_
Dense matrix.
Epetra_SerialDenseSolver DenseSolver_
Linear problem for dense matrix and vectors.
void PrintStatus() const
Print information about the factorization and solution phases.
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const &paramList)
Use this parameter list to read values from.
long long NumGlobalNonzeros_
int SolveDistributed(Epetra_MultiVector &X, const Epetra_MultiVector &B)
Solves the linear system, when more than one process is used.
const Epetra_Import & SolutionImporter()
Returns a reference to the solution importer (to domain map from serial map).
long long NumGlobalRows_
Teuchos::RCP< Teuchos::ParameterList > pl_
T * get() const
int NumGlobalRows() const
Returns the number of global rows, or -1 if Matrix() returns 0.
int NumNumericFact() const
Returns the number of numeric factorizations performed by this object.
int DenseToFactored()
Factors the matrix using LAPACK.
int NumNumericFact_
Number of numeric factorization phases.
Definition: Amesos_Status.h:69
int NumMyRows() const
Returns the number of local rows, or -1 if Matrix() returns 0.
const Epetra_LinearProblem * GetProblem() const
Returns the Epetra_LinearProblem.
Teuchos::RCP< Teuchos::ParameterList > ParameterList_
int NumSolve_
Number of solves.
Definition: Amesos_Status.h:71
int SymbolicFactorization()
Performs SymbolicFactorization on the matrix A.
int GEEV(Epetra_Vector &Er, Epetra_Vector &Ei)
Computes the eigenvalues of the linear system matrix using DGEEV.
Teuchos::RCP< Epetra_Export > RhsExporter_
int SetUseTranspose(bool UseTranspose_in)
If set true, X will be set to the solution of AT X = B (not A X = B)
void GetTiming(Teuchos::ParameterList &list) const
Load up the current timing information into the parameter list.
Definition: Amesos_Time.h:124
bool UseTranspose_
If true, the linear system with the transpose will be solved.
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
int NumSymbolicFact() const
Returns the number of symbolic factorizations performed by this object.
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:20
void GetTiming(Teuchos::ParameterList &TimingParameterList) const
Extracts timing information from the current solver and places it in the parameter list...
Teuchos::RCP< Epetra_RowMatrix > SerialMatrix_
Epetra_SerialDenseMatrix DenseLHS_
Dense LHS.
~Amesos_Lapack(void)
Amesos_Lapack Destructor.
Amesos_Time: Container for timing information.
Definition: Amesos_Time.h:50
const Epetra_Export & RhsExporter()
Returns a reference to the rhs exporter (from range map to serial map).
Amesos_Lapack: an interface to LAPACK.
Definition: Amesos_Lapack.h:65
const Epetra_Map & SerialMap()
Returns a reference to serial map (that with all elements on process 0).
int MtxRedistTime_
Quick access ids for the individual timings.
int SetParameters(Teuchos::ParameterList &ParameterList)
Deprecated - Sets parameters.
Teuchos::RCP< Epetra_Map > SerialMap_
Epetra_RowMatrix * GetMatrix() const
bool MatrixShapeOK() const
Returns true if the solver can handle this matrix shape.
Epetra_RowMatrix & SerialMatrix()
Returns a reference to serial matrix (that with all rows on process 0).
Teuchos::RCP< Epetra_Import > MatrixImporter_
const Epetra_LinearProblem * Problem_
Pointer to the linear problem.
int DistributedToSerial()
Converts a distributed matrix to serial matrix.
Epetra_CrsMatrix & SerialCrsMatrix()
int NumericFactorization()
Performs NumericFactorization on the matrix A.
long long NumGlobalRows64() const
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
Teuchos::RCP< Epetra_Import > SolutionImporter_
const Epetra_RowMatrix * Matrix() const
Returns a pointer to the linear system matrix.
int Solve()
Solves A X = B (or AT x = B)
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
This is an empty stub.
Amesos_NoCopiable: Simple class to prevent the usage of copy constructor and operator =...
int SolveSerial(Epetra_MultiVector &X, const Epetra_MultiVector &B)
Solves the linear system, when only one process is used.
void PrintTiming() const
Print timing information.
bool UseTranspose() const
Returns the current UseTranspose setting.
Amesos_Utils: Collections of basic utilities.
Definition: Amesos_Utils.h:19
const Epetra_Import & MatrixImporter()
Returns a reference to the matrix importer (from row map to serial map).