Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Amesos_Klu.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 
39 #ifndef AMESOS_KLU_H
40 #define AMESOS_KLU_H
41 
42 #if defined(Amesos_SHOW_DEPRECATED_WARNINGS)
43 #ifdef __GNUC__
44 #warning "The Amesos package is deprecated"
45 #endif
46 #endif
47 
48 #include "Amesos_ConfigDefs.h"
49 #include "Amesos_BaseSolver.h"
50 #include "Amesos_NoCopiable.h"
51 #include "Amesos_Utils.h"
52 #include "Amesos_Time.h"
53 #include "Amesos_Status.h"
54 #include "Amesos_Control.h"
55 #include "Epetra_LinearProblem.h"
56 #include "Epetra_Time.h"
57 #include "Epetra_Import.h"
58 #ifdef EPETRA_MPI
59 #include "Epetra_MpiComm.h"
60 #else
61 #include "Epetra_Comm.h"
62 #endif
63 #include "Epetra_CrsGraph.h"
64 #include "Epetra_CrsMatrix.h"
65 
66 // class EpetraExt::MultiVector_Reindex ;
67 // class EpetraExt::CrsMatrix_Reindex ;
69 
107 // Amesos_Klu_Pimpl contains a pointer to two structures defined in
108 // klu.h: trilinos_klu_symbolic and trilinos_klu_numeric. This prevents Amesos_Klu.h
109 // from having to include klu.h.
110 //
111 // Doxygen does not handle forward class references well.
112 #ifndef DOXYGEN_SHOULD_SKIP_THIS
113 class Amesos_Klu_Pimpl ;
114 class Amesos_StandardIndex ;
115 #endif
116 
118  private Amesos_Time,
119  private Amesos_NoCopiable,
120  private Amesos_Utils,
121  private Amesos_Control,
122  private Amesos_Status {
123 
124 public:
125 
127 
135  Amesos_Klu(const Epetra_LinearProblem& LinearProblem );
136 
138  ~Amesos_Klu(void);
139 
141 
142 
143  int SymbolicFactorization() ;
144 
145  int NumericFactorization() ;
146 
147  int Solve();
148 
150 
151 
153  const Epetra_LinearProblem *GetProblem() const { return(Problem_); };
154 
156 
159  bool MatrixShapeOK() const ;
160 
162 
166  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
167 
168  bool UseTranspose() const {return(UseTranspose_);};
169 
170  const Epetra_Comm & Comm() const {return(GetProblem()->GetOperator()->Comm());};
171 
173 
176 
179 
181  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
182 
184  void PrintTiming() const;
185 
187  void PrintStatus() const;
188 
190  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming( TimingParameterList ); }
191 
192 private:
193 
195 
196 
197  /*
198  CreateLocalMatrixAndExporters - Prepare to convert matrix and vectors to serial
199  Preconditions:
200  Problem_ must be set
201 
202  Postconditions:
203  UseDataInPlace_ is set to 1 if the input matrix can be used in place, i.e.
204  1) is entirely stored on process 0
205  2) range map and domain map are same as the row map
206  The following are only set if (! UseDataInPlace_ )"
207  SerialMap_
208  ImportToSerial_
209  SerialCrsMatrixA_
210 
211  SerialMatrix_
212  */
214  /*
215  ExportToSerial
216  Preconditions:
217  UseDataInPlace_ must be set
218  ImportToSerial and SerialCrsMatrixA_ must be set if UseDataInPlace_ != 1
219  Postconditions
220  SerialMatrix_ points to a serial version of the matrix
221  */
222  int ExportToSerial() ;
223  /*
224  ConvertToKluCRS - Convert matrix to form expected by Klu: Ai, Ap, Aval
225  Preconditions:
226  numentries_, RowMatrixA_, ImportToSerial_, StdIndexMatrix_, Reindex_
227  Postconditions:
228  SerialCrsMatrixA_
229  */
230  int ConvertToKluCRS(bool firsttime);
231 
232  /*
233  PerformSymbolicFactorization - Call Klu to perform symbolic factorization
234  Preconditions:
235  UseDataInPlace_ must be set to 1 if the input matrix is entirely stored on process 0
236  Ap, Ai and Aval point to a compressed row storage version of the input matrix A.
237  Postconditions:
238  Symbolic points to an KLU internal opaque object containing the
239  symbolic factorization and accompanying information.
240  SymbolicFactorizationOK_ = true;
241  Note: All action is performed on process 0
242 
243  Returns non-zero if the symbolic factorization failed
244  */
245 
247 
248  /*
249  PerformNumericFactorization - Call Klu to perform numeric factorization
250  Preconditions:
251  UseDataInPlace_ must be set
252  Ap, Ai and Aval point to a compressed row storage version of the input matrix A.
253  Symbolic must be set
254  Postconditions:
255  Numeric points to an KLU internal opaque object containing the
256  numeric factorization and accompanying information.
257  NumericFactorizationOK_ = true;
258  Note: All action is performed on process 0
259  */
261 
262  // @}
263 
265 
266 #ifdef Bug_8212
267  int *lose_this_;
268 #endif
269  //
270  // PrivateKluData_ contains pointers to data needed by klu whose
271  // data structures are defined by klu.h
272  //
277 
282  std::vector <int> Ap;
283  std::vector <int> VecAi;
284  std::vector <double> VecAval;
285  double* Aval;
286  int *Ai;
287 
291  long long numentries_;
294 
299 #if 0
305  Teuchos::RCP<Epetra_Map> ContiguousMap_;
306 #endif
315 
317  // serially, StorageOptimized, the LHS and RHS are assumed to be available
318  // when SymbolicFactorization is called and not to change (address or number
319  // of vectors) thereafter.
320  bool TrustMe_;
324  double *SerialXBvalues_ ;
325  double *SerialBvalues_ ;
332 
337 
339  std::vector<int> ColIndicesV_;
341  std::vector<double> RowValuesV_;
349 
350 }; // class Amesos_Klu
351 
352 #endif /* AMESOS_KLU_H */
Teuchos::RCP< Epetra_Import > ImportRangeToSerial_
Definition: Amesos_Klu.h:344
int NumSymbolicFact_
Number of symbolic factorization phases.
Definition: Amesos_Status.h:73
long long numentries_
Number of non-zero entries in Problem_-&gt;GetOperator()
Definition: Amesos_Klu.h:291
int PerformNumericFactorization()
Definition: Amesos_Klu.cpp:516
Amesos_Klu: A serial, unblocked code ideal for getting started and for very sparse matrices...
Definition: Amesos_Klu.h:117
Amesos_Control: Container for some control variables.
bool TrustMe_
If true, no checks are made and the matrix is assume to be distributed.
Definition: Amesos_Klu.h:320
std::vector< int > VecAi
Definition: Amesos_Klu.h:283
std::vector< int > Ap
Ap, Ai, Aval form the compressed row storage used by Klu Ai and Aval can point directly into a matrix...
Definition: Amesos_Klu.h:282
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
Definition: Amesos_Klu.h:170
Teuchos::RCP< Amesos_StandardIndex > StdIndexRange_
Definition: Amesos_Klu.h:275
Epetra_RowMatrix * RowMatrixA_
Operator converted to a RowMatrix.
Definition: Amesos_Klu.h:296
long long NumGlobalElements_
Number of rows and columns in the Problem_-&gt;GetOperator()
Definition: Amesos_Klu.h:293
int SolveTime_
Definition: Amesos_Klu.h:348
int Solve()
Solves A X = B (or AT x = B)
Definition: Amesos_Klu.cpp:729
int * Ai
Definition: Amesos_Klu.h:286
const Epetra_LinearProblem * GetProblem() const
Get a pointer to the Problem.
Definition: Amesos_Klu.h:153
Teuchos::RCP< Epetra_Map > SerialMap_
Points to a Serial Map (unused if UseDataInPlace_ == 1 )
Definition: Amesos_Klu.h:308
bool MatrixShapeOK() const
Returns true if KLU can handle this matrix shape.
Definition: Amesos_Klu.cpp:603
int ExportToSerial()
Definition: Amesos_Klu.cpp:119
Teuchos::RCP< Epetra_CrsMatrix > SerialCrsMatrixA_
Points to a Serial Copy of A (unused if UseDataInPlace_==1)
Definition: Amesos_Klu.h:310
int SerialXlda_
Definition: Amesos_Klu.h:264
bool UseTranspose() const
Returns the current UseTranspose setting.
Definition: Amesos_Klu.h:168
Teuchos::RCP< Epetra_Import > ImportDomainToSerial_
Definition: Amesos_Klu.h:345
Epetra_CrsMatrix * CrsMatrixA_
Operator converted to a CrsMatrix.
Definition: Amesos_Klu.h:298
int MtxRedistTime_
Quick access ids for the individual timings.
Definition: Amesos_Klu.h:347
Teuchos::RCP< Epetra_MultiVector > SerialX_
Definition: Amesos_Klu.h:328
int NumSolve() const
Returns the number of solves performed by this object.
Definition: Amesos_Klu.h:181
Teuchos::RCP< Amesos_StandardIndex > StdIndexDomain_
Definition: Amesos_Klu.h:276
int OverheadTime_
Definition: Amesos_Klu.h:348
Epetra_RowMatrix * SerialMatrix_
Points to a Serial Copy of A.
Definition: Amesos_Klu.h:314
double * Aval
Definition: Amesos_Klu.h:285
int NumNumericFact_
Number of numeric factorization phases.
Definition: Amesos_Status.h:75
double * SerialXBvalues_
Pointer to the actual values in the serial version of X and B.
Definition: Amesos_Klu.h:324
std::vector< double > RowValuesV_
Only used for RowMatrices to extract copies.
Definition: Amesos_Klu.h:341
int CreateLocalMatrixAndExporters()
Definition: Amesos_Klu.cpp:196
int NumSolve_
Number of solves.
Definition: Amesos_Status.h:77
int NumSymbolicFact() const
Returns the number of symbolic factorizations performed by this object.
Definition: Amesos_Klu.h:175
int NumNumericFact() const
Returns the number of numeric factorizations performed by this object.
Definition: Amesos_Klu.h:178
std::vector< int > ColIndicesV_
Only used for RowMatrices to extract copies.
Definition: Amesos_Klu.h:339
Teuchos::RCP< Epetra_Import > ImportToSerial_
Importer to process 0.
Definition: Amesos_Klu.h:343
int SymFactTime_
Definition: Amesos_Klu.h:348
void GetTiming(Teuchos::ParameterList &list) const
Load up the current timing information into the parameter list.
Definition: Amesos_Time.h:130
int NumFactTime_
Definition: Amesos_Klu.h:348
int ConvertToKluCRS(bool firsttime)
Definition: Amesos_Klu.cpp:337
int PerformSymbolicFactorization()
Definition: Amesos_Klu.cpp:478
Epetra_RowMatrix * StdIndexMatrix_
Points to a Contiguous Copy of A.
Definition: Amesos_Klu.h:312
int NumericFactorization()
Performs NumericFactorization on the matrix A.
Definition: Amesos_Klu.cpp:681
bool UseTranspose_
If true, the transpose of A is used.
Definition: Amesos_Klu.h:334
~Amesos_Klu(void)
Amesos_Klu Destructor.
Definition: Amesos_Klu.cpp:111
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:26
Amesos_Time: Container for timing information.
Definition: Amesos_Time.h:56
int NumVectors_
Number of vectors in RHS and LHS.
Definition: Amesos_Klu.h:322
int SetUseTranspose(bool UseTranspose_in)
SetUseTranpose(true) is more efficient in Amesos_Klu.
Definition: Amesos_Klu.h:166
Teuchos::RCP< Amesos_Klu_Pimpl > PrivateKluData_
Definition: Amesos_Klu.h:273
int SymbolicFactorization()
Performs SymbolicFactorization on the matrix A.
Definition: Amesos_Klu.cpp:617
void GetTiming(Teuchos::ParameterList &TimingParameterList) const
Extracts timing information and places in parameter list.
Definition: Amesos_Klu.h:190
const Epetra_LinearProblem * Problem_
Pointer to the linear system problem.
Definition: Amesos_Klu.h:336
int UseDataInPlace_
1 if Problem_-&gt;GetOperator() is stored entirely on process 0
Definition: Amesos_Klu.h:289
int MtxConvTime_
Definition: Amesos_Klu.h:347
Teuchos::RCP< Epetra_MultiVector > SerialBextract_
Definition: Amesos_Klu.h:331
void PrintStatus() const
Prints information about the factorization and solution phases.
Definition: Amesos_Klu.cpp:926
double * SerialBvalues_
Definition: Amesos_Klu.h:325
Teuchos::RCP< Epetra_MultiVector > SerialB_
Serial versions of the LHS and RHS (may point to the original vector if serial)
Definition: Amesos_Klu.h:327
Amesos_Klu(const Epetra_LinearProblem &LinearProblem)
Amesos_Klu Constructor.
Definition: Amesos_Klu.cpp:89
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
std::vector< double > VecAval
Definition: Amesos_Klu.h:284
int VecRedistTime_
Definition: Amesos_Klu.h:347
Teuchos::RCP< Epetra_MultiVector > SerialXextract_
Serial versions of the LHS and RHS (if necessary)
Definition: Amesos_Klu.h:330
Teuchos::RCP< Amesos_StandardIndex > StdIndex_
Definition: Amesos_Klu.h:274
Amesos_NoCopiable: Simple class to prevent the usage of copy constructor and operator =...
int SetParameters(Teuchos::ParameterList &ParameterList)
Updates internal variables.
Definition: Amesos_Klu.cpp:449
void PrintTiming() const
Prints timing information.
Definition: Amesos_Klu.cpp:949
Amesos_Utils: Collections of basic utilities.
Definition: Amesos_Utils.h:25