ITL_Preconditioner.h 5.86 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
26
27
28
29
30
31
32
33
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
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
// ==                                                                        ==
// ============================================================================

/** \file ITL_Preconditioner.h */

#ifndef AMDIS_ITL_PRECONDITIONER_H
#define AMDIS_ITL_PRECONDITIONER_H

#include "DOFMatrix.h"
#include "CreatorInterface.h"

#include <boost/numeric/itl/itl.hpp>
#include <boost/numeric/itl/pc/ilu_0.hpp>
#include <boost/numeric/itl/pc/ic_0.hpp>

namespace AMDiS {

  // ============================================================================
  // ===== class ITL_Preconditioner ==============================================
  // ============================================================================

 
  /**
   * \ingroup Solver
   * 
   * \brief Common base class for wrappers to use ITL preconditioners in AMDiS.
   */
  class ITL_BasePreconditioner 
  {
  public:
      typedef DOFMatrix::value_type value_type;

      virtual ~ITL_BasePreconditioner() {}

      virtual mtl::dense_vector<value_type>
      member_solve(const mtl::dense_vector<value_type>& vin) const = 0;

      virtual mtl::dense_vector<value_type>
      member_adjoint_solve(const mtl::dense_vector<value_type>& vin) const = 0;
  };

    inline mtl::dense_vector<DOFMatrix::value_type>
    solve(const ITL_BasePreconditioner& P, const mtl::dense_vector<DOFMatrix::value_type>& vin)
    {
	return P.member_solve(vin);
    }

    inline mtl::dense_vector<DOFMatrix::value_type>
    adjoint_solve(const ITL_BasePreconditioner& P, const mtl::dense_vector<DOFMatrix::value_type>& vin)
    {
	return P.member_adjoint_solve(vin);
    }

  /**
   * \ingroup Solver
   * 
   * \brief Wrapper for using ITL preconditioners in AMDiS.
   */
  template <typename Preconditioner>
  class ITL_Preconditioner : public ITL_BasePreconditioner
  {
    typedef DOFMatrix::value_type value_type;
  public:
    ITL_Preconditioner(const DOFMatrix::base_matrix_type& A) : precon(A) {}
   
    mtl::dense_vector<value_type>
    member_solve(const mtl::dense_vector<value_type>& vin) const
    {
      return solve(precon, vin);
    }

    mtl::dense_vector<value_type>
    member_adjoint_solve(const mtl::dense_vector<value_type>& vin) const
    {
      return adjoint_solve(precon, vin);
    }

    /// Creator class used in the OEMSolverMap.
    class Creator : public CreatorInterface<ITL_BasePreconditioner>
    {
    public:
Thomas Witkowski's avatar
Thomas Witkowski committed
98
      virtual ~Creator() {}
99
100
101
102
103
      
      /** \brief
       * Creates an ITL preconditioner
       */
      ITL_BasePreconditioner *create(const DOFMatrix::base_matrix_type& A) { 
Thomas Witkowski's avatar
Thomas Witkowski committed
104
	return new ITL_Preconditioner<Preconditioner>(A);
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
      }

      ITL_BasePreconditioner *create() {
	ERROR_EXIT("Must not be called! Only defined to avoid abstract function."); return 0;
      }
    };

  private:
      Preconditioner precon;
  };


  // ============================================================================
  // ===== class DiagonalPreconditioner =========================================
  // ============================================================================

  /**
   * \ingroup Solver
   * 
   * \brief
   * Diagonal preconditioner. 
   */
  class DiagonalPreconditioner 
      : public ITL_Preconditioner<itl::pc::diagonal<DOFMatrix::base_matrix_type> >
  {
      typedef ITL_Preconditioner<itl::pc::diagonal<DOFMatrix::base_matrix_type> > base;

      DiagonalPreconditioner(const DOFMatrix::base_matrix_type& A) : base(A) {}
  };

  /**
   * \ingroup Solver
   * 
   * \brief
   * Identity preconditioner. Behaves like no preconditioning.
   */
  class IdentityPreconditioner 
      : public ITL_Preconditioner<itl::pc::identity<DOFMatrix::base_matrix_type> >
  {
      typedef ITL_Preconditioner<itl::pc::identity<DOFMatrix::base_matrix_type> > base;

      IdentityPreconditioner(const DOFMatrix::base_matrix_type& A) : base(A) {}
  };


  // ============================================================================
  // ===== class ILUPreconditioner ==============================================
  // ============================================================================

  /**
   * \ingroup Solver
   * 
   * \brief
   * ILU (Incomplete LU factorization) preconditioner.
   *
   * The preconditioner is used from ITL. It corresponds for instance to "Iterative Methods for
   * Sparce Linear Systems", second edition, Yousef Saad. The preconditioner is
   * described in chapter 10.3 (algorithm 10.4).
   */
  class ILUPreconditioner 
      : public ITL_Preconditioner< itl::pc::ilu_0<DOFMatrix::base_matrix_type> >
  {};

  class ICPreconditioner 
      : public ITL_Preconditioner< itl::pc::ic_0<DOFMatrix::base_matrix_type> >
  {};


} // namespace AMDiS

#endif // AMDIS_ITL_PRECONDITIONER_H