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 #if defined(Amesos_SHOW_DEPRECATED_WARNINGS)
33 #ifdef __GNUC__
34 #warning "The Amesos package is deprecated"
35 #endif
36 #endif
37 
38 #include "Amesos_ConfigDefs.h"
39 #include "Amesos_BaseSolver.h"
40 #include "Amesos_NoCopiable.h"
41 #include "Amesos_Utils.h"
42 #include "Amesos_Time.h"
43 #include "Amesos_Status.h"
44 #include "Amesos_Control.h"
45 #include "Epetra_Comm.h"
46 #include "Epetra_Map.h"
50 #include "Epetra_CrsMatrix.h"
51 #include "Epetra_Import.h"
52 #include "Epetra_Export.h"
53 class Epetra_RowMatrix;
55 #include "Teuchos_RCP.hpp"
56 
58 
72  private Amesos_Time,
73  private Amesos_NoCopiable,
74  private Amesos_Utils,
75  private Amesos_Control,
76  private Amesos_Status
77 {
78 public:
79 
81 
89  Amesos_Lapack(const Epetra_LinearProblem& LinearProblem );
90 
92 
94  ~Amesos_Lapack(void);
95 
97 
98 
99  int SymbolicFactorization() ;
100 
101  int NumericFactorization() ;
102 
103  int Solve();
104 
106 
107 
108  const Epetra_LinearProblem *GetProblem() const { return(Problem_); };
109 
110  bool MatrixShapeOK() const;
111 
112  int SetUseTranspose(bool UseTranspose_in) {
113  UseTranspose_ = UseTranspose_in;
114  return(0);
115  }
116 
117  bool UseTranspose() const {return(UseTranspose_);};
118 
119  const Epetra_Comm & Comm() const {
120  return(GetProblem()->GetOperator()->Comm());
121  }
122 
124 
127 
129 
130 
133 
135 
145  int GEEV(Epetra_Vector& Er, Epetra_Vector& Ei);
146 
149 
152 
154  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
155 
157  void PrintTiming() const;
158 
160  void PrintStatus() const;
161 
163  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming(TimingParameterList); }
164 
166 
167 protected:
168 
170  const Epetra_RowMatrix* Matrix() const
171  {
172  return(Problem_->GetMatrix());
173  }
174 
176 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
177  inline int NumGlobalRows() const
178  {
179  return(Matrix()->NumGlobalRows());
180  }
181 #endif
182 
183  inline long long NumGlobalRows64() const
184  {
185  return(Matrix()->NumGlobalRows64());
186  }
187 
189  inline int NumMyRows() const
190  {
191  return(Matrix()->NumMyRows());
192  }
193 
195  inline const Epetra_Map& SerialMap()
196  {
197  return(*(SerialMap_.get()));
198  }
199 
202  {
203  return(*(SerialMatrix_.get()));
204  }
205 
207  {
208  return(*(SerialCrsMatrix_.get()));
209  }
210 
213  {
214  return(*(MatrixImporter_.get()));
215  }
216 
219  {
220  return(*(RhsExporter_.get()));
221  }
222 
225  {
226  return(*(SolutionImporter_.get()));
227  }
228 
230 
233  const Epetra_MultiVector& B);
234 
237  const Epetra_MultiVector& B);
238 
240  int DistributedToSerial();
241 
243  int SerialToDense();
244 
246  int DenseToFactored();
247 
254 
263 
268 
271 
272  long long NumGlobalRows_;
274 
276 
277 }; // End of class Amesos_Lapack
278 #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:73
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:75
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:77
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:130
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:26
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:56
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:71
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:25
const Epetra_Import & MatrixImporter()
Returns a reference to the matrix importer (from row map to serial map).