Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Amesos_Pardiso.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_PARDISO_H
30 #define AMESOS_PARDISO_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_LinearProblem.h"
40 #include "Epetra_Time.h"
41 #include "Epetra_Map.h"
42 #include "Epetra_Import.h"
43 #include "Epetra_Comm.h"
44 #include "Epetra_RowMatrix.h"
45 #include "Epetra_CrsMatrix.h"
47 #include "Teuchos_RCP.hpp"
48 
50 
58  private Amesos_Time,
59  private Amesos_NoCopiable,
60  private Amesos_Utils,
61  private Amesos_Control,
62  private Amesos_Status {
63 
64 public:
65 
67  Amesos_Pardiso(const Epetra_LinearProblem& LinearProblem );
69 
73 
75 
77  int SymbolicFactorization() ;
78 
80  int NumericFactorization() ;
81 
83  int Solve();
85 
87 
89  const Epetra_LinearProblem* GetProblem() const { return(Problem_); }
90 
92 
95  bool MatrixShapeOK() const;
96 
98 
103 
105  bool UseTranspose() const { return(UseTranspose_); }
106 
108  const Epetra_Comm& Comm() const { return(GetProblem()->GetOperator()->Comm()); }
109 
112 
115 
118 
120  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
121 
123  void PrintTiming() const;
124 
126  void PrintStatus() const;
127 
129  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming(TimingParameterList); }
130 
132 
133 private:
134 
135  int CheckError(const int error) const;
136 
137  inline const Epetra_Map& Map() const
138  {
139  return(Matrix_->RowMatrixRowMap());
140  }
141 
142  inline const Epetra_RowMatrix& Matrix() const
143  {
144  return(*Matrix_);
145  }
146 
148  {
149  return(*(SerialMap_.get()));
150  }
151 
153  {
154  return(*(SerialMatrix_.get()));
155  }
156 
158  {
159  return(*(SerialCrsMatrix_.get()));
160  }
161 
163  {
164  return(*(Importer_.get()));
165  }
166 
167  int ConvertToSerial();
168  int ConvertToPardiso();
171 
176 
177  const Epetra_Map* Map_;
179 
184 
188 
189  // Data for PARDISO
190  std::vector<double> aa_;
191  std::vector<int> ia_;
192  std::vector<int> ja_;
193 
195  int mtype_;
196  void* pt_[64];
197 
198  int iparm_[64];
199  double dparm_[64];
200  int maxfct_; // Maximal number of factors with idential nonzero pattern (always 1)
201  int mnum_;
202  int msglvl_;
203  int nrhs_; // Number of RHS
204 
206 
207 }; // class Amesos_Pardiso
208 #endif
const Epetra_LinearProblem * Problem_
Pointer to the linear system problem.
int NumSymbolicFact_
Number of symbolic factorization phases.
Definition: Amesos_Status.h:67
Teuchos::RCP< Epetra_RowMatrix > SerialMatrix_
const Epetra_LinearProblem * GetProblem() const
Get a pointer to the Problem.
Amesos_Control: Container for some control variables.
const Epetra_RowMatrix & Matrix() const
Epetra_Import & Importer()
virtual const Epetra_Map & RowMatrixRowMap() const =0
int PerformNumericFactorization()
int NumNumericFact() const
Returns the number of numeric factorizations performed by this object.
Amesos_Pardiso(const Epetra_LinearProblem &LinearProblem)
Constructor.
void PrintStatus() const
Prints information about the factorization and solution phases.
int SetUseTranspose(bool UseTranspose)
SetUseTranpose()
Epetra_CrsMatrix & SerialCrsMatrix()
void PrintTiming() const
Prints timing information.
T * get() const
Teuchos::RCP< Epetra_CrsMatrix > SerialCrsMatrix_
int SymbolicFactorization()
Performs SymbolicFactorization on the matrix A.
int NumSolve() const
Returns the number of solves performed by this object.
Amesos_Pardiso: Interface to the PARDISO package.
int nrhs_
Output level.
void GetTiming(Teuchos::ParameterList &TimingParameterList) const
Extracts timing information from the current solver and places it in the parameter list...
int NumNumericFact_
Number of numeric factorization phases.
Definition: Amesos_Status.h:69
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this matrix.
int msglvl_
Actual matrix for solution phase (always 1)
int NumSymbolicFact() const
Returns the number of symbolic factorizations performed by this object.
int NumSolve_
Number of solves.
Definition: Amesos_Status.h:71
std::vector< int > ia_
std::vector< double > aa_
int CheckError(const int error) const
const Epetra_Map * Map_
bool UseTranspose() const
Returns the current UseTranspose setting.
int SetParameters(Teuchos::ParameterList &ParameterList)
Set parameters from the input parameters list, returns 0 if successful.
const Epetra_RowMatrix * Matrix_
void GetTiming(Teuchos::ParameterList &list) const
Load up the current timing information into the parameter list.
Definition: Amesos_Time.h:124
const Epetra_Map & Map() const
int MtxConvTime_
Quick access pointers to the internal timing data.
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:20
int NumericFactorization()
Performs NumericFactorization on the matrix A.
void * pt_[64]
Amesos_Time: Container for timing information.
Definition: Amesos_Time.h:50
int PerformSymbolicFactorization()
~Amesos_Pardiso()
Destructor.
Epetra_RowMatrix & SerialMatrix()
Teuchos::RCP< Epetra_Map > SerialMap_
std::vector< int > ja_
Teuchos::RCP< Epetra_Import > Importer_
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
int Solve()
Solves A X = B (or AT X = B)
bool UseTranspose_
If true, the transpose of A is used.
double dparm_[64]
Amesos_NoCopiable: Simple class to prevent the usage of copy constructor and operator =...
bool MatrixShapeOK() const
Returns true if PARDISO can handle this matrix shape.
Teuchos::ParameterList param_
Epetra_Map & SerialMap()
Amesos_Utils: Collections of basic utilities.
Definition: Amesos_Utils.h:19