Dune Core Modules (2.4.1)

pardiso.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_ISTL_PARDISO_HH
4#define DUNE_ISTL_PARDISO_HH
5
8/* Change this, if your Fortran compiler does not append underscores. */
9/* e.g. the AIX compiler: #define F77_FUNC(func) func */
10
11#ifdef AIX
12#define F77_FUNC(func) func
13#else
14#define F77_FUNC(func) func ## _
15#endif
16
17
18#ifdef HAVE_PARDISO
19/* PARDISO prototype. */
20extern "C" int F77_FUNC(pardisoinit)
21 (void *, int *, int *);
22
23extern "C" int F77_FUNC(pardiso)
24 (void *, int *, int *, int *, int *, int *,
25 double *, int *, int *, int *, int *, int *,
26 int *, double *, double *, int *);
27#endif
28
29namespace Dune {
30
31
36 template<class M, class X, class Y>
37 class SeqPardiso : public Preconditioner<X,Y> {
38 public:
40 typedef M matrix_type;
42 typedef X domain_type;
44 typedef Y range_type;
46 typedef typename X::field_type field_type;
47
48 typedef typename M::RowIterator RowIterator;
49 typedef typename M::ColIterator ColIterator;
50
51 // define the category
52 enum {
55 };
56
62 SeqPardiso (const M& A)
63 : A_(A)
64 {
65#ifdef HAVE_PARDISO
66
67 mtype_ = 11;
68 nrhs_ = 1;
69 num_procs_ = 1;
70 maxfct_ = 1;
71 mnum_ = 1;
72 msglvl_ = 0;
73 error_ = 0;
74
75 n_ = A_.rowdim();
76 int nnz = 0;
77 RowIterator endi = A_.end();
78 for (RowIterator i = A_.begin(); i != endi; ++i)
79 {
80 if (A_.rowdim(i.index()) != 1)
81 DUNE_THROW(NotImplemented, "SeqPardiso: row blocksize != 1.");
82 ColIterator endj = (*i).end();
83 for (ColIterator j = (*i).begin(); j != endj; ++j) {
84 if (A_.coldim(j.index()) != 1)
85 DUNE_THROW(NotImplemented, "SeqPardiso: column blocksize != 1.");
86 nnz++;
87 }
88 }
89
90 std::cout << "dimension = " << n_ << ", number of nonzeros = " << nnz << std::endl;
91
92 a_ = new double[nnz];
93 ia_ = new int[n_+1];
94 ja_ = new int[nnz];
95
96 int count = 0;
97 for (RowIterator i = A_.begin(); i != endi; ++i)
98 {
99 ia_[i.index()] = count+1;
100 ColIterator endj = (*i).end();
101 for (ColIterator j = (*i).begin(); j != endj; ++j) {
102 a_[count] = *j;
103 ja_[count] = j.index()+1;
104
105 count++;
106 }
107 }
108 ia_[n_] = count+1;
109
110 F77_FUNC(pardisoinit) (pt_, &mtype_, iparm_);
111
112 int phase = 11;
113 int idum;
114 double ddum;
115 iparm_[2] = num_procs_;
116
117 F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
118 &n_, a_, ia_, ja_, &idum, &nrhs_,
119 iparm_, &msglvl_, &ddum, &ddum, &error_);
120
121 if (error_ != 0)
122 DUNE_THROW(MathError, "Constructor SeqPardiso: Factorization failed. Error code " << error_);
123
124 std::cout << "Constructor SeqPardiso: Factorization completed." << std::endl;
125
126#else
127 DUNE_THROW(NotImplemented, "no Pardiso library available, reconfigure with correct --with-pardiso options");
128#endif
129 }
130
136 virtual void pre (X& x, Y& b) {}
137
143 virtual void apply (X& v, const Y& d)
144 {
145#ifdef HAVE_PARDISO
146 int phase = 33;
147
148 iparm_[7] = 1; /* Max numbers of iterative refinement steps. */
149 int idum;
150
151 double x[n_];
152 double b[n_];
153 for (int i = 0; i < n_; i++) {
154 x[i] = v[i];
155 b[i] = d[i];
156 }
157
158 F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
159 &n_, a_, ia_, ja_, &idum, &nrhs_,
160 iparm_, &msglvl_, b, x, &error_);
161
162 if (error_ != 0)
163 DUNE_THROW(MathError, "SeqPardiso.apply: Backsolve failed. Error code " << error_);
164
165 for (int i = 0; i < n_; i++)
166 v[i] = x[i];
167
168 std::cout << "SeqPardiso: Backsolve completed." << std::endl;
169#endif
170 }
171
177 virtual void post (X& x) {}
178
180 {
181#ifdef HAVE_PARDISO
182 int phase = -1; // Release internal memory.
183 int idum;
184 double ddum;
185
186 F77_FUNC(pardiso) (pt_, &maxfct_, &mnum_, &mtype_, &phase,
187 &n_, &ddum, ia_, ja_, &idum, &nrhs_,
188 iparm_, &msglvl_, &ddum, &ddum, &error_);
189 delete[] a_;
190 delete[] ia_;
191 delete[] ja_;
192#endif
193 }
194
195 private:
196 M A_;
197 int n_;
198 double *a_;
199 int *ia_;
200 int *ja_;
201 int mtype_;
202 int nrhs_;
203 void *pt_[64];
204 int iparm_[64];
205 int maxfct_;
206 int mnum_;
207 int msglvl_;
208 int error_;
209 int num_procs_;
210 };
211
212 template<class M, class X, class Y>
213 struct IsDirectSolver<SeqPardiso<M,X,Y> >
214 {
215 enum { value=true};
216 };
217
218
219} // end namespace Dune
220#endif
Default exception class for mathematical errors.
Definition: exceptions.hh:266
Default exception for dummy implementations.
Definition: exceptions.hh:288
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:26
The sequential Pardiso preconditioner.
Definition: pardiso.hh:37
virtual void post(X &x)
Clean up.
Definition: pardiso.hh:177
@ category
The category the preconditioner is part of.
Definition: pardiso.hh:54
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: pardiso.hh:136
X domain_type
The domain type of the preconditioner.
Definition: pardiso.hh:42
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: pardiso.hh:143
M matrix_type
The matrix type the preconditioner is for.
Definition: pardiso.hh:40
SeqPardiso(const M &A)
Constructor.
Definition: pardiso.hh:62
Y range_type
The range type of the preconditioner.
Definition: pardiso.hh:44
X::field_type field_type
The field type of the preconditioner.
Definition: pardiso.hh:46
#define DUNE_THROW(E, m)
Definition: exceptions.hh:243
Dune namespace.
Definition: alignment.hh:10
Define general preconditioner interface.
Templates characterizing the type of a solver.
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:21
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)