dune-istl  2.2.1
pardiso.hh
Go to the documentation of this file.
1 #ifndef DUNE_PARDISO_HH
2 #define DUNE_PARDISO_HH
3 
6 /* Change this, if your Fortran compiler does not append underscores. */
7 /* e.g. the AIX compiler: #define F77_FUNC(func) func */
8 
9 #ifdef AIX
10 #define F77_FUNC(func) func
11 #else
12 #define F77_FUNC(func) func ## _
13 #endif
14 
15 
16 #ifdef HAVE_PARDISO
17 /* PARDISO prototype. */
18 extern "C" int F77_FUNC(pardisoinit)
19  (void *, int *, int *);
20 
21 extern "C" int F77_FUNC(pardiso)
22  (void *, int *, int *, int *, int *, int *,
23  double *, int *, int *, int *, int *, int *,
24  int *, double *, double *, int *);
25 #endif
26 
27 namespace Dune {
28 
29 
34  template<class M, class X, class Y>
35  class SeqPardiso : public Preconditioner<X,Y> {
36  public:
38  typedef M matrix_type;
40  typedef X domain_type;
42  typedef Y range_type;
44  typedef typename X::field_type field_type;
45 
46  typedef typename M::RowIterator RowIterator;
47  typedef typename M::ColIterator ColIterator;
48 
49  // define the category
50  enum {
53  };
54 
60  SeqPardiso (const M& A)
61  : A_(A)
62  {
63 #ifdef HAVE_PARDISO
64 
65  mtype_ = 11;
66  nrhs_ = 1;
67  num_procs_ = 1;
68  maxfct_ = 1;
69  mnum_ = 1;
70  msglvl_ = 0;
71  error_ = 0;
72 
73  n_ = A_.rowdim();
74  int nnz = 0;
75  RowIterator endi = A_.end();
76  for (RowIterator i = A_.begin(); i != endi; ++i)
77  {
78  if (A_.rowdim(i.index()) != 1)
79  DUNE_THROW(NotImplemented, "SeqPardiso: row blocksize != 1.");
80  ColIterator endj = (*i).end();
81  for (ColIterator j = (*i).begin(); j != endj; ++j) {
82  if (A_.coldim(j.index()) != 1)
83  DUNE_THROW(NotImplemented, "SeqPardiso: column blocksize != 1.");
84  nnz++;
85  }
86  }
87 
88  std::cout << "dimension = " << n_ << ", number of nonzeros = " << nnz << std::endl;
89 
90  a_ = new double[nnz];
91  ia_ = new int[n_+1];
92  ja_ = new int[nnz];
93 
94  int count = 0;
95  for (RowIterator i = A_.begin(); i != endi; ++i)
96  {
97  ia_[i.index()] = count+1;
98  ColIterator endj = (*i).end();
99  for (ColIterator j = (*i).begin(); j != endj; ++j) {
100  a_[count] = *j;
101  ja_[count] = j.index()+1;
102 
103  count++;
104  }
105  }
106  ia_[n_] = count+1;
107 
108  F77_FUNC(pardisoinit) (pt_, &mtype_, iparm_);
109 
110  int phase = 11;
111  int idum;
112  double ddum;
113  iparm_[2] = num_procs_;
114 
115  F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
116  &n_, a_, ia_, ja_, &idum, &nrhs_,
117  iparm_, &msglvl_, &ddum, &ddum, &error_);
118 
119  if (error_ != 0)
120  DUNE_THROW(MathError, "Constructor SeqPardiso: Factorization failed. Error code " << error_);
121 
122  std::cout << "Constructor SeqPardiso: Factorization completed." << std::endl;
123 
124 #else
125  DUNE_THROW(NotImplemented, "no Pardiso library available, reconfigure with correct --with-pardiso options");
126 #endif
127  }
128 
134  virtual void pre (X& x, Y& b) {}
135 
141  virtual void apply (X& v, const Y& d)
142  {
143 #ifdef HAVE_PARDISO
144  int phase = 33;
145 
146  iparm_[7] = 1; /* Max numbers of iterative refinement steps. */
147  int idum;
148 
149  double x[n_];
150  double b[n_];
151  for (int i = 0; i < n_; i++) {
152  x[i] = v[i];
153  b[i] = d[i];
154  }
155 
156  F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
157  &n_, a_, ia_, ja_, &idum, &nrhs_,
158  iparm_, &msglvl_, b, x, &error_);
159 
160  if (error_ != 0)
161  DUNE_THROW(MathError, "SeqPardiso.apply: Backsolve failed. Error code " << error_);
162 
163  for (int i = 0; i < n_; i++)
164  v[i] = x[i];
165 
166  std::cout << "SeqPardiso: Backsolve completed." << std::endl;
167 #endif
168  }
169 
175  virtual void post (X& x) {}
176 
178  {
179 #ifdef HAVE_PARDISO
180  int phase = -1; // Release internal memory.
181  int idum;
182  double ddum;
183 
184  F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
185  &n_, &ddum, ia_, ja_, &idum, &nrhs_,
186  iparm_, &msglvl_, &ddum, &ddum, &error_);
187  delete[] a_;
188  delete[] ia_;
189  delete[] ja_;
190 #endif
191  }
192 
193  private:
194  M A_;
195  int n_;
196  double *a_;
197  int *ia_;
198  int *ja_;
199  int mtype_;
200  int nrhs_;
201  void *pt_[64];
202  int iparm_[64];
203  int maxfct_;
204  int mnum_;
205  int msglvl_;
206  int error_;
207  int num_procs_;
208  };
209 
210  template<class M, class X, class Y>
211  struct IsDirectSolver<SeqPardiso<M,X,Y> >
212  {
213  enum{ value=true};
214  };
215 
216 
217 } // end namespace Dune
218 #endif
219