DUNE-FEM (unstable)

gmres.hh
1#ifndef DUNE_FEM_GMRES_HH
2#define DUNE_FEM_GMRES_HH
3
4#include <cmath>
5#include <cassert>
6#include <iostream>
7
8#include <utility>
9
10#include <dune/fem/solver/parameter.hh>
11
12namespace Dune
13{
14namespace Fem
15{
16namespace LinearSolver
17{
18
19 namespace detail {
20 template <class T>
21 class Matrix
22 {
23 public:
24 Matrix(int n, int m) : n_(n), m_(m), data_( n*m, T(0) ) {}
25 Matrix(int n) : Matrix(n,n) {}
26
27 // element access
28 T& operator()(int i, int j)
29 {
30 assert(i>=0 && i<n_ && j>=0 && j<m_);
31 return data_[i*m_ + j];
32 }
33
34 // const element access
35 T operator()(int i, int j) const
36 {
37 assert(i>=0 && i<n_ && j>=0 && j<m_);
38 return data_[i*m_ + j];
39 }
40
41 // conversion operators
42 operator T*() { return data_.data(); }
43 operator const T *() const { return data_.data(); }
44
45 protected:
46 const int n_, m_;
47 std::vector< T > data_;
48 };
49 }
50
52 template <class FieldType>
53 FieldType scalarProduct( const int dim, const FieldType *x, const FieldType* y )
54 {
55 FieldType scp = 0;
56 for( int i=0; i<dim; ++i )
57 {
58 scp += x[ i ] * y[ i ];
59 }
60 return scp;
61 }
62
63 // computes y = beta y + alpha op(A) x
64 template <class Communication, class FieldType, class DiscreteFunction>
65 void gemv(const Communication& comm,
66 const int m, // j+1
67 std::vector< DiscreteFunction >& v,
68 const DiscreteFunction& vjp,
69 FieldType *y // global_dot
70 )
71 {
72 for(int l=0; l<m; ++l)
73 {
74 y[ l ] = 0;
75 }
76
77 const auto& auxiliaryDofs = vjp.space().auxiliaryDofs();
78 const auto& vj = vjp.dofVector();
79
80 auto scp = [&m, &y, &vj, &v] (const size_t dof)
81 {
82 for(int l=0; l<m; ++l)
83 {
84 y[ l ] += (vj[ dof ] * v[ l ].dofVector()[ dof ]);
85 }
86 };
87 // see dune/fem/space/common/auiliarydofs.hh
88 forEachPrimaryDof( auxiliaryDofs, scp );
89
90 // communicate sum
91 comm.sum( y, m );
92 }
93
95 template<class FieldType>
96 void rotate( const int dim,
97 FieldType* x, FieldType* y,
98 const FieldType c, const FieldType s)
99 {
100 int i = dim;
101 while (i--)
102 {
103 const FieldType _x = *x;
104 const FieldType _y = *y;
105 *x = c*_x + s*_y;
106 *y = c*_y - s*_x;
107 ++x;
108 ++y;
109 }
110 }
111
112 // Saad, Youcef; Schultz, Martin H.
113 // GMRES: A generalized minimal residual algorithm for solving nonsymmetric
114 // linear systems. (English)
115 // [J] SIAM J. Sci. Stat. Comput. 7, 856-869 (1986). [ISSN 0196-5204]
116 template <class Operator, class Preconditioner, class DiscreteFunction>
117 inline int gmres( Operator& op, Preconditioner* preconditioner,
118 std::vector< DiscreteFunction >& v,
119 DiscreteFunction& u,
120 const DiscreteFunction& b,
121 const int m, // gmres inner iterations
122 const double tolerance,
123 const int maxIterations,
124 const int toleranceCriteria,
125 std::ostream* os = nullptr )
126 {
127 typedef typename DiscreteFunction :: RangeFieldType FieldType;
128
129 const auto& comm = u.space().gridPart().comm();
130
131 detail::Matrix< FieldType > H( m+1, m ); // \in \R^{m+1 \times m}
132 std::vector< FieldType > g_( 6*m, 0.0 );
133
134 FieldType* g = g_.data();
135 FieldType* s = g + (m+1);
136 FieldType* c = s + m;
137 FieldType* y = c + m;
138
139 DiscreteFunction& v0 = v[ 0 ];
140
141 std::vector< FieldType > global_dot( m+1, FieldType(0) );
142
143 double _tolerance = tolerance;
144 if (toleranceCriteria == ToleranceCriteria::relative)
145 {
146 global_dot[ 0 ] = b.scalarProductDofs( b );
147 _tolerance *= std::sqrt(global_dot[0]);
148 }
149
150 int iterations = 0;
151 while (true)
152 {
153 // start
154 op(u, v0);
155
156 v0 -= b ;
157
158 // scalarProduct( dim, v, v );
159 // contains 1 allreduce (global zum)
160 global_dot[ 0 ] = v0.scalarProductDofs( v0 );
161
162 //comm.allreduce(1, local_dot, global_dot, MPI_SUM);
163 FieldType res = std::sqrt(global_dot[0]);
164
165 if (toleranceCriteria == ToleranceCriteria::residualReduction && iterations==0)
166 {
167 _tolerance *= res;
168 }
169
170 assert(res == res);
171 if (os)
172 {
173 (*os) << "Fem::GMRES outer iteration : " << res << std::endl;
174 }
175
176 // make sure this also works if reisual is 0
177 if (res <= _tolerance*(1+1e-15)) break;
178
179 g[0] = -res;
180 for(int i=1; i<=m; i++) g[i] = 0.0;
181
182 // cblas_dscal(dim, 1.0/res, v, 1);
183 v0 *= (1.0/res);
184
185 //scale( dim, 1.0/res, v );
186
187 // iterate
188 for(int j=0; j<m; j++)
189 {
190 DiscreteFunction& vj = v[ j ];
191 DiscreteFunction& vjp = v[ j + 1 ];
192
193 // apply the linear operator (perhaps in combination with the
194 // preconditioner)
195 if (preconditioner)
196 {
197 DiscreteFunction& z = v[ m+1 ];
198 (*preconditioner)(vj, z );
199 op( z, vjp);
200 }
201 else
202 {
203 op(vj, vjp);
204 }
205
206 // contains 1 allreduce (global sum)
207 gemv(comm, j+1, v, vjp, global_dot.data());
208
209 for(int i=0; i<=j; i++) H(i,j) = global_dot[i];
210
211 // assuming beta == 1.0
212 for(int l=0; l<j+1; ++l)
213 {
214 vjp.axpy( -global_dot[l], v[l] );
215 }
216
217 // 1 allreduce (global sum)
218 global_dot[ 0 ] = vjp.scalarProductDofs( vjp );
219
220 H(j+1,j) = std::sqrt(global_dot[0]);
221 // scale(dim, 1.0/H(j+1,j), vjp );
222 vjp *= 1.0/H(j+1,j);
223
224 // perform Givens rotation
225 for(int i=0; i<j; i++)
226 {
227 rotate(1, &H(i+1,j), &H(i,j), c[i], s[i]);
228 }
229
230 const FieldType h_j_j = H(j,j);
231 const FieldType h_jp_j = H(j+1,j);
232 const FieldType norm = std::sqrt(h_j_j*h_j_j + h_jp_j*h_jp_j);
233 c[j] = h_j_j / norm;
234 s[j] = -h_jp_j / norm;
235 rotate(1, &H(j+1,j), &H(j,j), c[j], s[j]);
236 rotate(1, &g[j+1], &g[j], c[j], s[j]);
237
238 assert(g[j+1] == g[j+1]);
239 if ( os )
240 {
241 (*os) << "Fem::GMRES it: " << iterations << " : " << std::abs(g[j+1]) << std::endl;
242 }
243
244 ++iterations;
245 if (std::abs(g[j+1]) < _tolerance
246 || iterations >= maxIterations ) break;
247 }
248
249 //
250 // form the approximate solution
251 //
252
253 int last = iterations%m;
254 if (last == 0) last = m;
255
256 // compute y via backsubstitution
257 for(int i=last-1; i>=0; --i)
258 {
259 const FieldType dot = scalarProduct( last-(i+1), &H(i,i)+1, &y[i+1] );
260 y[i] = (g[i] - dot)/ H(i,i);
261 }
262
263 // update the approx. solution
264 if (preconditioner)
265 {
266 // u += M^{-1} (v[0], ..., v[last-1]) y
267 DiscreteFunction& u_tmp = v[ m ]; // we don't need this vector anymore
268 DiscreteFunction& z = v[ m+1 ];
269 u_tmp.clear();
270
271 // u += (v[0], ..., v[last-1]) y
272 for(int i=0; i<last; ++i)
273 {
274 u_tmp.axpy( y[ i ], v[ i ] );
275 }
276
277 (*preconditioner)(u_tmp, z);
278 u += z;
279 }
280 else{
281 // u += (v[0], ..., v[last-1]) y
282 for(int i=0; i<last; ++i)
283 {
284 u.axpy( y[ i ], v[ i ] );
285 }
286 }
287
288 if (std::abs(g[last]) < _tolerance) break;
289 }
290
291 // output
292 if ( os ) {
293 (*os) << "Fem::GMRES: number of iterations: "
294 << iterations
295 << std::endl;
296 }
297
298 return (iterations < maxIterations) ? iterations : -iterations;
299 }
300
301
302} // end namespace Solver
303
304} // end namespace Fem
305
306} // end namespace Dune
307
308#endif
auto dot(const A &a, const B &b) -> typename std::enable_if< IsNumber< A >::value &&!IsVector< A >::value &&!std::is_same< typename FieldTraits< A >::field_type, typename FieldTraits< A >::real_type > ::value, decltype(conj(a) *b)>::type
computes the dot product for fundamental data types according to Petsc's VectDot function: dot(a,...
Definition: dotproduct.hh:42
static void forEachPrimaryDof(const AuxiliaryDofs &auxiliaryDofs, F &&f)
Apply action encoded in Functor f to all primary dofs.
Definition: auxiliarydofs.hh:303
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Sep 4, 22:38, 2025)