Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Amesos_Paraklete.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_PARAKLETE_H
40 #define AMESOS_PARAKLETE_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 #ifdef HAVE_AMESOS_EPETRAEXT
67 #endif
68 
69 
71 
90 // Amesos_Paraklete_Pimpl contains a pointer to structures defined in
91 // paraklete.h. This prevents Amesos_Paraklete.h
92 // from having to include paraklete.h.
93 //
94 // Doxygen does not handle forward class references well.
95 #ifndef DOXYGEN_SHOULD_SKIP_THIS
97 class Amesos_StandardIndex ;
98 #endif
99 
101  private Amesos_Time,
102  private Amesos_NoCopiable,
103  private Amesos_Utils,
104  private Amesos_Control,
105  private Amesos_Status {
106 
107 public:
108 
110 
118  Amesos_Paraklete(const Epetra_LinearProblem& LinearProblem );
119 
121  ~Amesos_Paraklete(void);
122 
124 
125 
126  int SymbolicFactorization() ;
127 
128  int NumericFactorization() ;
129 
130  int Solve();
131 
133 
134 
136  const Epetra_LinearProblem *GetProblem() const { return(Problem_); };
137 
139 
142  bool MatrixShapeOK() const ;
143 
145 
149  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
150 
151  bool UseTranspose() const {return(UseTranspose_);};
152 
153  const Epetra_Comm & Comm() const {return(GetProblem()->GetOperator()->Comm());};
154 
156 
159 
162 
164  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
165 
167  void PrintTiming() const;
168 
170  void PrintStatus() const;
171 
173  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming(TimingParameterList); }
174 
175 private:
176 
178 
179 
180  /*
181  CreateLocalMatrixAndExporters - Prepare to convert matrix and vectors to serial
182  Preconditions:
183  Problem_ must be set
184 
185  Postconditions:
186  UseDataInPlace_ is set to 1 if the input matrix can be used in place, i.e.
187  1) is entirely stored on process 0
188  2) range map and domain map are same as the row map
189  The following are only set if (! UseDataInPlace_ )"
190  SerialMap_
191  ImportToSerial_
192  SerialCrsMatrixA_
193 
194  SerialMatrix_
195  */
197  /*
198  ExportToSerial
199  Preconditions:
200  UseDataInPlace_ must be set
201  ImportToSerial and SerialCrsMatrixA_ must be set if UseDataInPlace_ != 1
202  Postconditions
203  SerialMatrix_ points to a serial version of the matrix
204  */
205  int ExportToSerial() ;
206  /*
207  ConvertToParakleteCRS - Convert matrix to form expected by Paraklete: Ai, Ap, Aval
208  Preconditions:
209  numentries_, RowMatrixA_, ImportToSerial_, StdIndexMatrix_, Reindex_
210  Postconditions:
211  SerialCrsMatrixA_
212  */
213  int ConvertToParakleteCRS(bool firsttime);
214 
215  /*
216  PerformSymbolicFactorization - Call Paraklete to perform symbolic factorization
217  Preconditions:
218  UseDataInPlace_ must be set to 1 if the input matrix is entirely stored on process 0
219  Ap, Ai and Aval point to a compressed row storage version of the input matrix A.
220  Postconditions:
221  Symbolic points to an PARAKLETE internal opaque object containing the
222  symbolic factorization and accompanying information.
223  SymbolicFactorizationOK_ = true;
224  Note: All action is performed on process 0
225  */
226 
228 
229  /*
230  PerformNumericFactorization - Call Paraklete to perform numeric factorization
231  Preconditions:
232  UseDataInPlace_ must be set
233  Ap, Ai and Aval point to a compressed row storage version of the input matrix A.
234  Symbolic must be set
235  Postconditions:
236  Numeric points to an PARAKLETE internal opaque object containing the
237  numeric factorization and accompanying information.
238  NumericFactorizationOK_ = true;
239  Note: All action is performed on process 0
240  */
242 
243  // @}
244 
245  bool IamInGroup_; // True if this process is involved in the computation. Set by SymbolicFactorization
246 
248  //
249  // PrivateParakleteData_ contains pointers to data needed by paraklete whose
250  // data structures are defined by paraklete.h
251  //
256  MPI_Comm ParakleteComm_;
257 
262  std::vector <long> Ap;
263  std::vector <long> Ai;
264  std::vector <double> VecAval;
265  double* Aval;
266 
273 
278  //
279  // transposer_ transposes a CrsMatrix
280  // Created in CreateLocalMatrixAndExporters
281  // Used in ExportToSerial()
282  //
283 #ifdef HAVE_AMESOS_EPETRAEXT
285 #endif
296 
298  // serially, StorageOptimized, the LHS and RHS are assumed to be available
299  // when SymbolicFactorization is called and not to change (address or number
300  // of vectors) thereafter.
301  bool TrustMe_;
305  double *SerialXBvalues_ ;
306  double *SerialBvalues_ ;
313 
318 
320  std::vector<int> ColIndicesV_;
322  std::vector<double> RowValuesV_;
327 
331 
332 }; // class Amesos_Paraklete
333 
334 #endif /* AMESOS_PARAKLETE_H */
int NumSymbolicFact_
Number of symbolic factorization phases.
Definition: Amesos_Status.h:73
int MtxConvTime_
Quick access pointers to internal timing information.
Amesos_Control: Container for some control variables.
int SymbolicFactorization()
Performs SymbolicFactorization on the matrix A.
Teuchos::RCP< Epetra_Import > ImportToSerial_
Importer to process 0.
~Amesos_Paraklete(void)
Amesos_Paraklete Destructor.
Amesos_Paraklete: A serial, unblocked code ideal for getting started and for very sparse matrices...
Epetra_MultiVector * SerialB_
Serial versions of the LHS and RHS (may point to the original vector if serial)
int NumSymbolicFact() const
Returns the number of symbolic factorizations performed by this object.
std::vector< double > RowValuesV_
Only used for RowMatrices to extract copies.
Epetra_RowMatrix * StdIndexMatrix_
Points to a Contiguous Copy of A.
int numentries_
Number of non-zero entries in Problem_-&gt;GetOperator()
Teuchos::RCP< Epetra_MultiVector > SerialBextract_
void GetTiming(Teuchos::ParameterList &TimingParameterList) const
Extracts timing information from the current solver and places it in the parameter list...
Amesos_Paraklete(const Epetra_LinearProblem &LinearProblem)
Amesos_Paraklete Constructor.
std::vector< int > ColIndicesV_
Only used for RowMatrices to extract copies.
Epetra_MultiVector * StdIndexRangeVector_
int NumNumericFact_
Number of numeric factorization phases.
Definition: Amesos_Status.h:75
bool UseTranspose_
If true, the transpose of A is used.
Epetra_MultiVector * SerialX_
int NumSolve() const
Returns the number of solves performed by this object.
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
Teuchos::RCP< Epetra_CrsMatrix > SerialCrsMatrixA_
Points to a Serial Copy of A (unused if UseDataInPlace_==1)
int NumSolve_
Number of solves.
Definition: Amesos_Status.h:77
Epetra_MultiVector * StdIndexDomainVector_
Teuchos::RCP< Epetra_Import > ImportDomainToSerial_
int SetUseTranspose(bool UseTranspose_in)
SetUseTranpose()
Teuchos::RCP< Epetra_Map > SerialMap_
Points to a Serial Map (unused if UseDataInPlace_ == 1 )
int NumNumericFact() const
Returns the number of numeric factorizations performed by this object.
int NumGlobalElements_
Number of rows and columns in the Problem_-&gt;GetOperator()
const Epetra_LinearProblem * GetProblem() const
Get a pointer to the Problem.
Epetra_RowMatrix * RowMatrixA_
Operator converted to a RowMatrix.
void GetTiming(Teuchos::ParameterList &list) const
Load up the current timing information into the parameter list.
Definition: Amesos_Time.h:130
void PrintStatus() const
Prints information about the factorization and solution phases.
Teuchos::RCP< Amesos_Paraklete_Pimpl > PrivateParakleteData_
Epetra_CrsMatrix * CrsMatrixA_
Operator converted to a CrsMatrix.
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:26
int NumVectors_
Number of vectors in RHS and LHS.
Teuchos::RCP< Amesos_StandardIndex > StdIndexDomain_
int ConvertToParakleteCRS(bool firsttime)
Amesos_Time: Container for timing information.
Definition: Amesos_Time.h:56
int Solve()
Solves A X = B (or AT x = B)
void PrintTiming() const
Prints timing information.
std::vector< long > Ai
int NumericFactorization()
Performs NumericFactorization on the matrix A.
int SetParameters(Teuchos::ParameterList &ParameterList)
Updates internal variables.
Teuchos::RCP< Amesos_StandardIndex > StdIndexRange_
bool MatrixShapeOK() const
Returns true if PARAKLETE can handle this matrix shape.
int UseDataInPlace_
1 if Problem_-&gt;GetOperator() is stored entirely on process 0
Epetra_RowMatrix * SerialMatrix_
Points to a Serial Copy of A.
Teuchos::RCP< Epetra_MultiVector > SerialXextract_
Serial versions of the LHS and RHS (if necessary)
bool TrustMe_
If true, no checks are made and the matrix is assume to be distributed.
const Epetra_LinearProblem * Problem_
Pointer to the linear system problem.
Teuchos::RCP< Epetra_Import > ImportRangeToSerial_
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
double * SerialXBvalues_
Pointer to the actual values in the serial version of X and B.
Amesos_NoCopiable: Simple class to prevent the usage of copy constructor and operator =...
std::vector< long > Ap
Ap, Ai, Aval form the compressed row storage used by Paraklete Ai and Aval can point directly into a ...
std::vector< double > VecAval
Amesos_Utils: Collections of basic utilities.
Definition: Amesos_Utils.h:25
bool UseTranspose() const
Returns the current UseTranspose setting.
Teuchos::RCP< Amesos_StandardIndex > StdIndex_