Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Amesos_Mumps.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_MUMPS_H
30 #define AMESOS_MUMPS_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 class Epetra_Import;
39 class Epetra_RowMatrix;
40 class Epetra_MultiVector;
41 #include "Epetra_Import.h"
42 #include "Epetra_CrsMatrix.h"
43 #include "Epetra_Map.h"
47 class Amesos_EpetraInterface;
48 
49 #include "Amesos_ConfigDefs.h"
50 #include "Amesos_BaseSolver.h"
51 #include "Amesos_NoCopiable.h"
52 #include "Amesos_Utils.h"
53 #include "Amesos_Time.h"
54 #include "Amesos_Status.h"
55 #include "Amesos_Control.h"
56 #include "Epetra_LinearProblem.h"
57 #ifdef EPETRA_MPI
58 #include "Epetra_MpiComm.h"
59 #else
60 #include "Epetra_Comm.h"
61 #endif
62 #include "Teuchos_RCP.hpp"
63 #include <map>
64 using namespace Teuchos;
65 
67 
114 extern "C" {
115 #include "dmumps_c.h"
116 }
117 
119  private Amesos_Time,
120  private Amesos_NoCopiable,
121  private Amesos_Utils,
122  private Amesos_Control,
123  private Amesos_Status {
124 
125 public:
126 
128 
132  Amesos_Mumps(const Epetra_LinearProblem& LinearProblem);
133 
135 
137  ~Amesos_Mumps(void);
139 
141 
142  int SymbolicFactorization() ;
143 
144  int NumericFactorization() ;
145 
146  int Solve();
147 
149  void Destroy();
150 
151  int SetUseTranspose(bool UseTranspose_in)
152  {
153  UseTranspose_ = UseTranspose_in;
154  if (UseTranspose_in)
155  return (-1);
156  return (0);
157  };
158 
159  bool UseTranspose() const {return(UseTranspose_);};
160 
161  int SetParameters( Teuchos::ParameterList &ParameterList );
162 
164 
166 
169 
172 
174  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
175 
177 
179 
181 
185  void PrintTiming() const;
186 
188 
191  void PrintStatus() const;
192 
194  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming(TimingParameterList); }
195 
197 
199 
200 
201 #if 0
202 
211  int ComputeSchurComplement(bool flag,
212  int NumSchurComplementRows, int * SchurComplementRows);
213 
215 
220  Epetra_CrsMatrix * GetCrsSchurComplement();
221 
223 
228  Epetra_SerialDenseMatrix * GetDenseSchurComplement();
229 #endif
230 
232 
239  int SetPrecscaling(double * ColSca, double * RowSca )
240  {
241  ColSca_ = ColSca;
242  RowSca_ = RowSca;
243  return 0;
244  }
245 
247 
252  int SetRowScaling(double * RowSca )
253  {
254  RowSca_ = RowSca;
255  return 0;
256  }
257 
259 
264  int SetColScaling(double * ColSca )
265  {
266  ColSca_ = ColSca;
267  return 0;
268  }
269 
271 
275  int SetOrdering(int * PermIn)
276  {
277  PermIn_ = PermIn;
278  return 0;
279  }
280 
282 
286  double * GetRINFO() ;
287 
289 
291  int * GetINFO() ;
292 
294 
298  double * GetRINFOG() ;
299 
301 
304  int * GetINFOG() ;
305 
307  void SetICNTL(int pos, int value);
308 
310  void SetCNTL(int pos, double value);
311 
313 
314  bool MatrixShapeOK() const
315  {
316  bool OK = true;
317 
318  if ( GetProblem()->GetOperator()->OperatorRangeMap().NumGlobalPoints() !=
319  GetProblem()->GetOperator()->OperatorDomainMap().NumGlobalPoints() ) OK = false;
320  return OK;
321 }
322 
323 
325  const Epetra_Comm & Comm() const {return(GetProblem()->GetOperator()->Comm());};
326 
328  const Epetra_LinearProblem * GetProblem() const { return(Problem_); };
329 
330 protected:
331 
333  Epetra_RowMatrix& Matrix();
334 
335  const Epetra_RowMatrix& Matrix() const;
336 
338  Epetra_Map& RedistrMap();
339 
341  Epetra_Import& RedistrImporter();
342 
344  Epetra_RowMatrix& RedistrMatrix(const bool ImportMatrix = false);
345 
347  Epetra_Map& SerialMap();
348 
350  Epetra_Import& SerialImporter();
351 
353  int ConvertToTriplet(const bool OnlyValues);
354 
356  int CheckError();
357 
359  void CheckParameters();
360 
361  void SetICNTLandCNTL();
362 
367 
368 
369  bool NoDestroy_ ; // Set true to prevent memory freeing
370 
372  std::vector <int> Row;
374  std::vector<int> Col;
376  std::vector<double> Val;
377 
380 
383 
385  int MtxConvTime_, MtxRedistTime_, VecRedistTime_;
386  int SymFactTime_, NumFactTime_, SolveTime_;
387 
389  double * RowSca_, * ColSca_;
390 
392  int * PermIn_;
393 
398 
403 
404 #ifdef HAVE_MPI
405  MPI_Comm MUMPSComm_;
407 #endif
408 
411 
422 
423  DMUMPS_STRUC_C MDS;
424 
425  std::map<int, int> ICNTL;
426  std::map<int, double> CNTL;
427 }; // class Amesos_Mumps
428 
429 #endif /* AMESOS_MUMPS_H */
RCP< Epetra_Import > SerialImporter_
Importer from Matrix.OperatorDomainMap() to SerialMap_.
Definition: Amesos_Mumps.h:421
int NumSymbolicFact_
Number of symbolic factorization phases.
Definition: Amesos_Status.h:73
Amesos_Control: Container for some control variables.
bool MatrixShapeOK() const
Returns true if the solver can handle this matrix shape.
Definition: Amesos_Mumps.h:314
bool UseTranspose() const
Returns the current UseTranspose setting.
Definition: Amesos_Mumps.h:159
std::map< int, int > ICNTL
Definition: Amesos_Mumps.h:425
RCP< Epetra_CrsMatrix > CrsSchurComplement_
Pointer to the Schur complement, as CrsMatrix.
Definition: Amesos_Mumps.h:400
RCP< Epetra_Map > SerialMap_
Map with all elements on process 0 (for solution and rhs).
Definition: Amesos_Mumps.h:419
int NumSchurComplementRows_
Number of rows in the Schur complement (if required)
Definition: Amesos_Mumps.h:395
Amesos_Mumps: An object-oriented wrapper for the double precision version of MUMPS.
Definition: Amesos_Mumps.h:118
std::map< int, double > CNTL
Definition: Amesos_Mumps.h:426
int * PermIn_
PermIn for MUMPS.
Definition: Amesos_Mumps.h:392
int * SchurComplementRows_
Rows for the Schur complement (if required)
Definition: Amesos_Mumps.h:397
int NumNumericFact() const
Returns the number of numeric factorizations performed by this object.
Definition: Amesos_Mumps.h:171
RCP< Epetra_Map > RedistrMap_
Redistributed matrix.
Definition: Amesos_Mumps.h:413
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this matrix.
Definition: Amesos_Mumps.h:325
const Epetra_LinearProblem * Problem_
Pointer to the linear problem to be solved.
Definition: Amesos_Mumps.h:410
std::vector< int > Col
column indices of nonzero elements
Definition: Amesos_Mumps.h:374
int NumNumericFact_
Number of numeric factorization phases.
Definition: Amesos_Status.h:75
void GetTiming(Teuchos::ParameterList &TimingParameterList) const
Extracts timing information from the current solver and places it in the parameter list...
Definition: Amesos_Mumps.h:194
bool IsConvertToTripletOK_
true if matrix has already been converted to COO format
Definition: Amesos_Mumps.h:364
int NumSymbolicFact() const
Returns the number of symbolic factorizations performed by this object.
Definition: Amesos_Mumps.h:168
int NumSolve_
Number of solves.
Definition: Amesos_Status.h:77
int NumSolve() const
Returns the number of solves performed by this object.
Definition: Amesos_Mumps.h:174
int SetPrecscaling(double *ColSca, double *RowSca)
Set prescaling.
Definition: Amesos_Mumps.h:239
int SetUseTranspose(bool UseTranspose_in)
If set true, X will be set to the solution of AT X = B (not A X = B)
Definition: Amesos_Mumps.h:151
int SetColScaling(double *ColSca)
Set column scaling.
Definition: Amesos_Mumps.h:264
void GetTiming(Teuchos::ParameterList &list) const
Load up the current timing information into the parameter list.
Definition: Amesos_Time.h:130
int MaxProcs_
Maximum number of processors in the MUMPS&#39; communicator.
Definition: Amesos_Mumps.h:379
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:26
Amesos_Time: Container for timing information.
Definition: Amesos_Time.h:56
DMUMPS_STRUC_C MDS
Definition: Amesos_Mumps.h:423
RCP< Epetra_CrsMatrix > RedistrMatrix_
Redistributed matrix (only if MaxProcs_ &gt; 1).
Definition: Amesos_Mumps.h:417
int SetOrdering(int *PermIn)
Sets ordering.
Definition: Amesos_Mumps.h:275
std::vector< int > Row
row indices of nonzero elements
Definition: Amesos_Mumps.h:372
RCP< Epetra_SerialDenseMatrix > DenseSchurComplement_
Pointer to the Schur complement,as DenseMatrix.
Definition: Amesos_Mumps.h:402
std::vector< double > Val
values of nonzero elements
Definition: Amesos_Mumps.h:376
const Epetra_LinearProblem * GetProblem() const
Gets a pointer to the Epetra_LinearProblem.
Definition: Amesos_Mumps.h:328
bool IsComputeSchurComplementOK_
true if the Schur complement has been computed (need to free memory)
Definition: Amesos_Mumps.h:366
bool CheckError(const std::string SolverType, const std::string Descriptor, const Epetra_RowMatrix &A, const Epetra_MultiVector &x, const Epetra_MultiVector &b, const Epetra_MultiVector &x_exact)
int Solve(int, TYPE *, TYPE *, TYPE *)
int SetRowScaling(double *RowSca)
Set row scaling.
Definition: Amesos_Mumps.h:252
RCP< Epetra_Import > RedistrImporter_
Redistributed importer (from Matrix().RowMatrixRowMap() to RedistrMatrix().RowMatrixRowMap()).
Definition: Amesos_Mumps.h:415
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
double * RowSca_
Row and column scaling.
Definition: Amesos_Mumps.h:389
Amesos_NoCopiable: Simple class to prevent the usage of copy constructor and operator =...
bool UseTranspose_
If true, solve the problem with AT.
Definition: Amesos_Mumps.h:382
Amesos_Utils: Collections of basic utilities.
Definition: Amesos_Utils.h:25