ITL_Solver.hpp 15.6 KB
Newer Older
1
2
3
/** \file ITL_Solver.h */

#pragma once
4
5
6
#ifdef HAVE_CONFIG_H
# include "config.h"
#endif
7
8

// MTL4 includes
9
10
11
12
13
14
15
16
17
18
19
#include <boost/numeric/itl/itl.hpp>
// #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>
20
21
22
23
24
25
26
27
28
29
30
31

// more solvers defined in AMDiS
#include "itl/minres.hpp"
#include "itl/gcr.hpp"
#include "itl/fgmres.hpp"
#include "itl/fgmres_householder.hpp"
#include "itl/gmres2.hpp"
#include "itl/gmres_householder.hpp"
#include "itl/preonly.hpp"

#include <dune/amdis/Initfile.hpp>
#include <dune/amdis/linear_algebra/mtl/LinearSolver.hpp>
32
#include <dune/amdis/linear_algebra/mtl/KrylovRunner.hpp>
33
#include <dune/amdis/linear_algebra/mtl/UmfpackRunner.hpp>
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58

namespace AMDiS
{
  /**
   * \ingroup Solver
   *
   * \brief
   * Creator for MTL4 itl-solvers.
   *
   * One of the following solvers can be chosen:
   * - @ref CGSolver "cg" (conjugate gradient method)
   * - @ref CGSSolver "cgs" (squared conjugate gradient method)
   * - @ref BiCGSolver "bicg" (biconjugate gradient method)
   * - @ref BiCGStabSolver "bicgstab" (stabilized BiCG method)
   * - @ref BiCGStab2Solver "bicgstab2" (stabilized BiCG(l) method with l=2)
   * - @ref QMRSolver "qmr" (Quasi-Minimal Residual method)
   * - @ref TFQMRSolver "tfqmr" (Transposed-Free Quasi-Minimal Residual method)
   * - @ref BiCGStabEllSolver "bicgstab_ell" (stabilized BiCG(l) method)
   * - @ref GMResSolver "gmres" (generalized minimal residual method)
   * - @ref IDRsSolver "idr_s" (Induced Dimension Reduction method)
   * - @ref MinResSolver "minres" (minimal residual method)
   * - @ref GcrSolver "gcr" (generalized conjugate residual method)
   * - @ref FGMResSolver "fgmres" (flexible GMRes method)
   * - @ref PreOnly "preonly" (solver that implements pure preconditioning applied to the rhs)
   */
59
  
60
61
62
63
64
65
66
67
68
69
70
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
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
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
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
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
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
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::CGSolver
   * \brief ITL_Solver <\ref cg_solver_type> implementation of conjugate gradient
   * method \implements ITL_Solver
   *
   * Solves a linear system \f$ Ax=b \f$ by the conjugate gradient method (CG)
   * and can be used for symmetric positive definite system matrices.
   * Right preconditioner is ignored.
   */

  class cg_solver_type
  {
  public:
    cg_solver_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      return itl::cg(A, x, b, l, r, iter);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::CGSSolver
   * \brief ITL_Solver <\ref cgs_solver_type> implementation of squared conjugate
   * gradient method \implements ITL_Solver
   *
   * Solves a linear system \f$ Ax=b \f$ by the squared conjugate gradient method
   * (CGS). Right preconditioner is ignored.
   */

  class cgs_solver_type
  {
  public:
    cgs_solver_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const&, I& iter)
    {
      return itl::cgs(A, x, b, l, iter);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::BiCGSolver
   * \brief ITL_Solver <\ref bicg_solver_type> implementation of bi-conjugate
   * gradient method \implements ITL_Solver
   *
   * Solves a linear system \f$ Ax=b \f$ by a BiCG method and can be used for
   * system matrices.
   */

  class bicg_solver_type
  {
  public:
    bicg_solver_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const&, I& iter)
    {
      return itl::bicg(A, x, b, l, iter);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::BiCGStabSolver
   * \brief ITL_Solver <\ref bicgstab_type> implementation of stabilized
   * bi-conjugate gradient method \implements ITL_Solver
   *
   * Solves a linear system \f$ Ax=b \f$ by a stabilized BiCG method and can be
   * used for system matrices.
   */

  class bicgstab_type
  {
  public:
    bicgstab_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const&, I& iter)
    {
      return itl::bicgstab(A, x, b, l, iter);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::BiCGStab2Solver
   * \brief ITL_Solver <\ref bicgstab2_type> implementation of BiCGStab(l) method
   * with l=2 \implements ITL_Solver
   *
   * Solves a linear system \f$ Ax=b \f$ by a stabilized BiCG(2) method and can
   * be used for system matrices.
   */

  class bicgstab2_type
  {
  public:
    bicgstab2_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const&, I& iter)
    {
      return itl::bicgstab_2(A, x, b, l, iter);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::QMRSolver
   * \brief ITL_Solver <\ref qmr_solver_type> implementation of Quasi-Minimal
   * Residual method \implements ITL_Solver
   *
   * Solves a linear system \f$ Ax=b \f$ by the Quasi-Minimal Residual method (QMR).
   */

  class qmr_solver_type
  {
  public:
    qmr_solver_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      return itl::qmr(A, x, b, l, r, iter);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::TFQMRSolver
   * \brief ITL_Solver <\ref tfqmr_solver_type> implementation of Transposed-Free
   * Quasi-Minimal Residual method \implements ITL_Solver
   *
   * Solves a linear system by the Transposed-Free Quasi-Minimal Residual method
   * (TFQMR). Does not use preconditioning currently.
   */

  class tfqmr_solver_type
  {
  public:
    tfqmr_solver_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      return itl::tfqmr(A, x, b, l, r, iter);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::BiCGStabEllSolver
   * \brief ITL_Solver <\ref bicgstab_ell_type> implementation of stabilized
   * BiCG(ell) method \implements ITL_Solver
   *
   * Solves a linear system by a stabilized BiCG(ell) method and can be used for
   * system matrices. The parameter ell [3] can be specified.
   */

  class bicgstab_ell_type
  {
    int ell;
  public:
    bicgstab_ell_type(std::string name) : ell(3)
    {
      Parameters::get(name + "->ell", ell);
    }
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      return itl::bicgstab_ell(A, x, b, l, r, iter, ell);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::GMResSolver
   * \brief ITL_Solver <\ref gmres_type> implementation of generalized minimal
   * residual method \implements ITL_Solver
   *
   * Solves a linear system by the GMRES method.
   * The parameter restart [30] is the maximal number of orthogonalized vectors.
   * The method is not preconditioned
   */

  class gmres_type
  {

  public:
    gmres_type(std::string name) : restart(30), ortho(GRAM_SCHMIDT)
    {
      Parameters::get(name + "->restart", restart);
      Parameters::get(name + "->orthogonalization", ortho);
    }
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      switch (ortho)
      {
      default:
      case GRAM_SCHMIDT:
        return itl::gmres2(A, x, b, l, r, iter, restart);
        break;
      case HOUSEHOLDER:
        return itl::gmres_householder(A, x, b, l, iter, restart);
        break;
      }
    }
    
  private:
    static constexpr int GRAM_SCHMIDT = 1;
    static constexpr int HOUSEHOLDER = 2;
    
    int restart;
    int ortho;
  };


  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::IDRsSolver
   * \brief ITL_Solver <\ref idr_s_type> implementation of Induced Dimension
   * Reduction method \implements ITL_Solver
   *
   * Solves a linear system by an Induced Dimension Reduction method and can be
   * used for system matrices.  The parameter s [30] can be specified.
   *
   * 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 idr_s_type
  {
    int s;
  public:
    idr_s_type(std::string name) : s(30)
    {
      Parameters::get(name + "->s", s);
    }
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      return itl::idr_s(A, x, b, l, r, iter, s);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::MinResSolver
   * \brief ITL_Solver <\ref minres_solver_type> implementation of minimal
   * residual method \implements ITL_Solver
   *
   * Solves a linear system by the Minres method. Can be used for symmetric
   * indefinite systems.
   */

  class minres_solver_type
  {
  public:
    minres_solver_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      return itl::minres(A, x, b, l, r, iter);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::GcrSolver
   * \brief ITL_Solver <\ref gcr_type> implementation of generalized conjugate
   * residual method \implements ITL_Solver
   *
   * Solves a linear system by the GCR method - generalized conjugate residual
   * method. The parameter restart [30] is the maximal number of orthogonalized
   * vectors.
   */

  class gcr_type
  {
    int restart;

  public:
    gcr_type(std::string name) : restart(30)
    {
      Parameters::get(name + "->restart", restart);
    }
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      return itl::gcr(A, x, b, l, r, iter, restart);
    }
  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::FGMResSolver
   * \brief ITL_Solver <\ref fgmres_type> implementation of flexible GMRes method
   * \implements ITL_Solver
   *
   * Solves a linear system by the FGMRES method.
   * The parameter restart [30] is the maximal number of orthogonalized vectors.
   * See reference "A Flexible Inner-Outer Preconditiones GMRES Algorithm",
   * Youcef Saad, (1993)
   */

  class fgmres_type
  {
  public:
    fgmres_type(std::string name) : restart(30), orthogonalization(GRAM_SCHMIDT)
    {
      Parameters::get(name + "->restart", restart);
      Parameters::get(name + "->orthogonalization", orthogonalization);
    }
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& r, I& iter)
    {
      switch (orthogonalization)
      {
      default:
      case GRAM_SCHMIDT:
        return itl::fgmres(A, x, b, l, r, iter, restart);
        break;
      case HOUSEHOLDER:
        return itl::fgmres_householder(A, x, b, r, iter, restart);
        break;
      }
    }
    
  private:
    static constexpr int GRAM_SCHMIDT = 1;
    static constexpr int HOUSEHOLDER = 2;
    
    int restart;
    int orthogonalization;

  };

  // ===========================================================================

  /**
   * \ingroup Solver
   * \class AMDiS::PreOnly
   * \brief ITL_Solver <\ref preonly_type> implementation of preconditioner as
   * \implements ITL_Solver
   *
   * Solves a linear system by applying a preconditioner only.
   */
  class preonly_type
  {
  public:
    preonly_type(std::string /*name*/) {}
    template <class LinOp, class X, class B, class L, class R, class I>
    int operator()(LinOp const& A, X& x, B const& b, L const& l, R const& /*r*/, I& iter)
    {
      return itl::preonly(A, x, b, l, iter);
    }
  };

  // ===========================================================================

  
  
448
449
450
451
452
453
454
455
456
  template <class T, class Param, class Vector>
  struct SolverCreator< mtl::compressed2D<T, Param>, Vector >
  {
    using Matrix = mtl::compressed2D<T, Param>;
    using SolverBase = LinearSolverInterface<Matrix, Vector, Vector>;
    
    template <class ITLSolver>
    using Solver = LinearSolver<Matrix, Vector, 
                      KrylovRunner<ITLSolver, Matrix, BaseVector<Vector>>>;
457
    
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
    /// Instantiate a new linear solver
    static std::unique_ptr<SolverBase> create(std::string name, std::string prefix)
    { 
      if (name == "cg")
        return std::make_unique<Solver<cg_solver_type>>(prefix);
      else if (name == "cgs")
        return std::make_unique<Solver<cgs_solver_type>>(prefix);
      else if (name == "bicg")
        return std::make_unique<Solver<bicg_solver_type>>(prefix);
      else if (name == "bicgstab" || name == "bcgs")
        return std::make_unique<Solver<bicgstab_type>>(prefix);
      else if (name == "bicgstab2")
        return std::make_unique<Solver<bicgstab2_type>>(prefix);
      else if (name == "qmr")
        return std::make_unique<Solver<qmr_solver_type>>(prefix);
      else if (name == "tfqmr")
        return std::make_unique<Solver<tfqmr_solver_type>>(prefix);
      else if (name == "bicgstab_ell" || name == "bcgsl")
        return std::make_unique<Solver<bicgstab_ell_type>>(prefix);
      else if (name == "gmres")
        return std::make_unique<Solver<gmres_type>>(prefix);
      else if (name == "fgmres")
        return std::make_unique<Solver<fgmres_type>>(prefix);
      else if (name == "minres")
        return std::make_unique<Solver<minres_solver_type>>(prefix);
      else if (name == "idrs")
        return std::make_unique<Solver<idr_s_type>>(prefix);
      else if (name == "gcr")
        return std::make_unique<Solver<gcr_type>>(prefix);
      else if (name == "preonly")
        return std::make_unique<Solver<preonly_type>>(prefix);
489
490
491
492
#ifdef HAVE_UMFPACK
      else if (name == "umfpack" || name == "direct")
        return std::make_unique< LinearSolver<Matrix, Vector, UmfpackRunner<Matrix, BaseVector<Vector>>> >(prefix);
#endif
493
494
495
496
      else
        AMDIS_ERROR_EXIT("Unknown Solver-name!");
    }
  };
497
498
  
} // end namespace AMDiS