MTLPreconPfc.h 5.08 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
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
/******************************************************************************
 *
 * Extension of AMDiS - Adaptive multidimensional simulations
 *
 * Copyright (C) 2013 Dresden University of Technology. All Rights Reserved.
 * Web: https://fusionforge.zih.tu-dresden.de/projects/amdis
 *
 * Authors: Simon Praetorius et al.
 *
 * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
 * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
 *
 *
 * See also license.opensource.txt in the distribution.
 * 
 ******************************************************************************/
#ifndef MTL_PFC_PRECONDITIONER_H
#define MTL_PFC_PRECONDITIONER_H

#include "AMDiS.h"
#include "BlockPreconditioner.h"

using namespace std;
using namespace AMDiS;    

template< typename Base >
struct PfcSchurMatrix
{
  PfcSchurMatrix() : super(NULL) {};
  PfcSchurMatrix(Base* super) : super(super)
  {
    m = num_rows(super->getM());
    y1.change_dim(m);
    y2.change_dim(m);
    y3.change_dim(m);
  }

  template <typename VectorIn, typename VectorOut, typename Assign>
  void mult(const VectorIn& b, VectorOut& x, Assign) const;

  template <typename VectorIn>
  mtl::vector::mat_cvec_multiplier<PfcSchurMatrix<Base>, VectorIn> operator*(const VectorIn& v) const
  {   return mtl::vector::mat_cvec_multiplier<PfcSchurMatrix<Base>, VectorIn>(*this, v);    }

  int m;
  
private:
  Base* super;
  
  mutable MTLVector y1;
  mutable MTLVector y2;
  mutable MTLVector y3;
};


struct MTLPreconPfc : BlockPreconditioner
{
  typedef BlockPreconditioner super;
 
  class Creator : public PreconditionCreator
  {
  public:
    virtual ~Creator() {}
    
    BasePreconditioner* create() { 
      return new MTLPreconPfc(this->name);
    }
  };  
  
  MTLPreconPfc(std::string name)
  : BlockPreconditioner(),
    name(name), 
    tau(NULL),
    solverM(NULL), solverM2(NULL), solverMpL(NULL), runnerM(NULL), runnerM2(NULL), runnerMpL(NULL),
    Pid(NULL), Pdiag(NULL),
76
    maxIterS(30), restartS(100)
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
  {
    createSubSolver(name + "->subsolver M", solverM, runnerM, "cg");
    createSubSolver(name + "->subsolver M2", solverM2, runnerM2, "cg");
    createSubSolver(name + "->subsolver MpL", solverMpL, runnerMpL, "cg");
    
    Parameters::get(name + "->subsolver S->solver->max iteration", maxIterS);
    Parameters::get(name + "->subsolver S->solver->restart", restartS);
  }
  
  virtual ~MTLPreconPfc()
  {
    delete solverM;
    delete solverMpL;
    
    if (Pid) { delete Pid; Pid = NULL; }
    if (Pdiag) { delete Pdiag; Pdiag = NULL; }
  }
  
  void setData(double* tau_)
  {
    tau = tau_;
  }
  
  virtual void init(const typename super::BlockMatrix& A, const MTLTypes::MTLMatrix& fullMatrix);
  
  virtual void exit()
  {
    runnerM->exit();
    runnerM2->exit();
    runnerMpL->exit();
    
    if (Pid) { delete Pid; Pid = NULL; }
    if (Pdiag) { delete Pdiag; Pdiag = NULL; }
  }
  
  virtual void solve(const MTLVector& b, MTLVector& x) const;  
    
  template <typename VectorX, typename VectorB>
  void solveM(VectorX& x, const VectorB& b, int nIter = -1) const
  {
    runnerM->solve(getM(), x, b);
  }
  template <typename VectorX, typename VectorB>
  void solveM2(VectorX& x, const VectorB& b) const
  {
    runnerM2->solve(getM(), x, b);
  }
  
  template <typename VectorX, typename VectorB>
  void solveMpK(VectorX& x, const VectorB& b) const
  {
    runnerMpL->solve(MpL, x, b);
  }
  
  template <typename SchurMatrix, typename VectorX, typename VectorB>
  void solveS(const SchurMatrix& S, VectorX& x, const VectorB& b) const
  { 
    itl::basic_iteration<double> iter(b, maxIterS, 0, 0);
    x = 0.0;
136
    itl::fgmres_householder(S, x, b, *Pdiag, iter, restartS);	// S*x = b
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
  }
  
  const MTLMatrix& getM() const { return A->getSubMatrix(2,2); }
  const MTLMatrix& getL() const { return A->getSubMatrix(2,1); }
  double getTau() const { return *tau; }
  
protected: 
  MTLMatrix MpL;
  PfcSchurMatrix<MTLPreconPfc> S;
   
  mutable MTLVector y0;
  mutable MTLVector y1;
  mutable MTLVector tmp;
  
  LinearSolver* solverM;
  LinearSolver* solverM2;
  LinearSolver* solverMpL;
  
  MTL4Runner<MTLTypes::MTLMatrix, MTLTypes::MTLVector>* runnerM;
  MTL4Runner<MTLTypes::MTLMatrix, MTLTypes::MTLVector>* runnerM2;
  MTL4Runner<MTLTypes::MTLMatrix, MTLTypes::MTLVector>* runnerMpL;
  
  itl::pc::identity<MTLTypes::MTLMatrix, double> *Pid;
  itl::pc::diagonal<MTLTypes::MTLMatrix, double> *Pdiag;
  
  std::string name;
  
  double* tau;
  
  int maxIterS;
  int restartS;
};

// some necessary functions that must be available for matrix-types

inline std::size_t size(const PfcSchurMatrix<MTLPreconPfc>& A) { return A.m * A.m; }
inline std::size_t num_rows(const PfcSchurMatrix<MTLPreconPfc>& A) { return A.m; }
inline std::size_t num_cols(const PfcSchurMatrix<MTLPreconPfc>& A) { return A.m; }

namespace mtl {
  template <>
  struct Collection<PfcSchurMatrix<MTLPreconPfc> >
  {
    typedef double value_type;
    typedef int    size_type;
  };

  namespace ashape {
    template <> struct ashape_aux<PfcSchurMatrix<MTLPreconPfc> > 
    {       typedef nonscal type;    };
  }
}

#include "MTLPreconPfc.hh"

#endif // MTL_PFC_PRECONDITIONER_H