ITL_Solver.h 12.7 KB
Newer Older
1
2
3
4
5
6
// ============================================================================
// ==                                                                        ==
// == AMDiS - Adaptive multidimensional simulations                          ==
// ==                                                                        ==
// ============================================================================
// ==                                                                        ==
7
// ==  TU Dresden                                                            ==
8
// ==                                                                        ==
9
10
11
// ==  Institut fr Wissenschaftliches Rechnen                               ==
// ==  Zellescher Weg 12-14                                                  ==
// ==  01069 Dresden                                                         ==
12
13
14
15
// ==  germany                                                               ==
// ==                                                                        ==
// ============================================================================
// ==                                                                        ==
16
// ==  https://gforge.zih.tu-dresden.de/projects/amdis/                      ==
17
18
19
20
21
22
23
24
25
// ==                                                                        ==
// ============================================================================

/** \file ITL_Preconditioner.h */

#ifndef AMDIS_ITL_SOLVER_H
#define AMDIS_ITL_SOLVER_H

#include "ITL_OEMSolver.h"
26
27
28
29
30
31
32
33
34
35
#include <boost/numeric/itl/krylov/bicg.hpp>
#include <boost/numeric/itl/krylov/bicgstab_2.hpp>
#include <boost/numeric/itl/krylov/bicgstab_ell.hpp>
#include <boost/numeric/itl/krylov/bicgstab.hpp>
#include <boost/numeric/itl/krylov/cg.hpp>
#include <boost/numeric/itl/krylov/cgs.hpp>
#include <boost/numeric/itl/krylov/gmres.hpp>
#include <boost/numeric/itl/krylov/idr_s.hpp>
#include <boost/numeric/itl/krylov/qmr.hpp>
#include <boost/numeric/itl/krylov/tfqmr.hpp>
36
#include <boost/numeric/itl/krylov/minres.hpp>
37

38
39
40
41
42
43
44
45
46
47
48

namespace AMDiS {

  // ============================================================================
  // ===== class CGSolver =======================================================
  // ============================================================================
  
  /// Helper class to establish a type from a function
  class cg_solver_type
  {
  public:
49
50
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R& r, I& iter)
51
52
53
54
55
56
57
58
    { return itl::cg(A, x, b, l, r, iter); }
  };


  /**
   * \ingroup Solver
   * 
   * \brief
59
   * Solves a linear system by the conjugate gradient method (CG) and can be used for
60
   * symmetric positive definite system matrices.
61
   * Right preconditioner is ignored.
62
63
64
65
66
67
68
69
70
   */
  class CGSolver : public ITL_OEMSolver<cg_solver_type>
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
    CGSolver(std::string name) : ITL_OEMSolver<cg_solver_type>(name) {}
  };


71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
  // ============================================================================
  // ===== class CGSSolver =======================================================
  // ============================================================================
  
  /// Helper class to establish a type from a function
  class cgs_solver_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R&, I& iter)
    { return itl::cgs(A, x, b, l, iter); }
  };


  /**
   * \ingroup Solver
   * 
   * \brief
   * Solves a linear system by the squared conjugate gradient method (CGS).
   * Right preconditioner is ignored.
   */
  class CGSSolver : public ITL_OEMSolver<cgs_solver_type>
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
    CGSSolver(std::string name) : ITL_OEMSolver<cgs_solver_type>(name) {}
  };


100
101
102
103
104
105
106
107
  // ============================================================================
  // ===== class BiCGSolver =====================================================
  // ============================================================================

  /// Helper class to establish a type from a function
  class bicg_solver_type
  {
  public:
108
109
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R&, I& iter)
110
111
112
113
114
115
116
117
118
119
120
121
122
123
    { return itl::bicg(A, x, b, l, iter); }
  };

  /**  
   * \ingroup Solver
   *
   * \brief
   * Solves a linear system by a stabilized BiCG method and can be used for 
   * system matrices.
   */
  class BiCGSolver : public ITL_OEMSolver<bicg_solver_type>
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
124
    BiCGSolver(std::string name) : ITL_OEMSolver<bicg_solver_type>(name) {}
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
  };


  // ============================================================================
  // ===== class BiCGStab =======================================================
  // ============================================================================

  /// Helper class to establish a type from a function
  class bicgstab_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R&, I& iter)
    { return itl::bicgstab(A, x, b, l, iter); }
  };


  /**  
   * \ingroup Solver
   *
   * \brief
   * Solves a linear system by a stabilized BiCG method and can be used for 
   * system matrices.
   */
149
  class BiCGStabSolver : public ITL_OEMSolver<bicgstab_type>
150
151
152
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
153
    BiCGStabSolver(std::string name) : ITL_OEMSolver<bicgstab_type>(name) {}
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
  };


  // ============================================================================
  // ===== class BiCGStab2 ======================================================
  // ============================================================================

  /// Helper class to establish a type from a function
  class bicgstab2_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R&, I& iter)
    { return itl::bicgstab_2(A, x, b, l, iter); }
  };


  /**  
   * \ingroup Solver
   *
   * \brief
   * Solves a linear system by a stabilized BiCG method and can be used for 
   * system matrices.
   */
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
  class BiCGStab2Solver : public ITL_OEMSolver<bicgstab2_type>
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
    BiCGStab2Solver(std::string name) : ITL_OEMSolver<bicgstab2_type>(name) {}
  };


  // ============================================================================
  // ===== class QMRSolver =======================================================
  // ============================================================================
  
  /// Helper class to establish a type from a function
  class qmr_solver_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R& r, I& iter)
    { return itl::qmr(A, x, b, l, r, iter); }
  };


  /**
   * \ingroup Solver
   * 
   * \brief
   * Solves a linear system by the Quasi-Minimal Residual method (QMR).
   */
  class QMRSolver : public ITL_OEMSolver<qmr_solver_type>
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
    QMRSolver(std::string name) : ITL_OEMSolver<qmr_solver_type>(name) {}
  };


  // ============================================================================
  // ===== class TFQMRSolver =======================================================
  // ============================================================================
  
  /// Helper class to establish a type from a function
  class tfqmr_solver_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
223
224
225
226
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R& r, I& iter)
    { 
      return itl::tfqmr(A, x, b, l, r, iter); 
    }
227
228
229
230
231
232
233
234
235
236
237
  };


  /**
   * \ingroup Solver
   * 
   * \brief
   * Solves a linear system by the Transposed-Free Quasi-Minimal Residual method (TFQMR).
   * Does not use preconditioning currently.
   */
  class TFQMRSolver : public ITL_OEMSolver<tfqmr_solver_type>
238
239
240
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
241
    TFQMRSolver(std::string name) : ITL_OEMSolver<tfqmr_solver_type>(name) {}
242
243
244
  };


245

246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
  // ============================================================================
  // ===== class BiCGStabEll ======================================================
  // ============================================================================

  /// Helper class to establish a type from a function
  class bicgstab_ell_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R& r, I& iter, int ell)
    { return itl::bicgstab_ell(A, x, b, l, r, iter, ell); }
  };

  /**  
   * \ingroup Solver
   *
   * \brief
263
   * Solves a linear system by a stabilized BiCG(l) method and can be used for 
264
265
   * system matrices.
   */
266
  class BiCGStabEllSolver : public ITL_OEMSolver_para<bicgstab_ell_type>
267
268
269
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
270
    BiCGStabEllSolver(std::string name) : ITL_OEMSolver_para<bicgstab_ell_type>(name) {}
271
272
  };

273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
  // ============================================================================
  // ===== class GMRES ======================================================
  // ============================================================================

  /// Helper class to establish a type from a function
  class gmres_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L&, const R&, I& iter, int ell)
    { return itl::gmres(A, x, b, iter, ell); }
  };

  /**  
   * \ingroup Solver
   *
   * \brief
   * Solves a linear system by the GMRES method.
   * The parameter ell is the maximal number of orthogonalized vectors.
   * The method is not preconditioned 
   */
  class GMResSolver : public ITL_OEMSolver_para<gmres_type>
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
    GMResSolver(std::string name) : ITL_OEMSolver_para<gmres_type>(name) {}
  };


  // ============================================================================
  // ===== class IDRs ======================================================
  // ============================================================================

  /// Helper class to establish a type from a function
  class idr_s_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R& r, I& iter, int ell)
    { return itl::idr_s(A, x, b, l, r, iter, ell); }
  };

  /**  
   * \ingroup Solver
   *
   * \brief
   * Solves a linear system by an Induced Dimension Reduction method and can be used for 
   * system matrices.
   * The parameter s can be specified as ell.
   * Peter Sonneveld and Martin B. van Gijzen, IDR(s): a family of simple and fast algorithms for solving large nonsymmetric linear systems. 
   * SIAM J. Sci. Comput. Vol. 31, No. 2, pp. 1035-1062 (2008). (copyright SIAM)
   * 
   */
  class IDRsSolver : public ITL_OEMSolver_para<idr_s_type>
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
    IDRsSolver(std::string name) : ITL_OEMSolver_para<idr_s_type>(name) {}
  };


334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
  // ============================================================================
  // ===== class Minres  ========================================================
  // ============================================================================
  
  /// Helper class to establish a type from a function
  class minres_solver_type
  {
  public:
    template < class LinOp, class X, class B, class L, class R, class I >
    int operator()(const LinOp& A, X& x, const B& b, const L& l, const R& r, I& iter)
    { 
      return itl::minres(A, x, b, l, r, iter); 
    }
  };


  /**
   * \ingroup Solver
   * 
   * \brief
   * Solves a linear system by the Minres method. Can be used for symmetric 
   * indefinite systems.
   */
  class MinResSolver : public ITL_OEMSolver<minres_solver_type>
  {
  public:
    /// The constructor reads required parameters and sets solvers \ref name.
    MinResSolver(std::string name) : ITL_OEMSolver<minres_solver_type>(name) {}
  };


365
366
367
368
369

} // namespace AMDiS

#endif // AMDIS_ITL_SOLVER