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
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
// ============================================================================
// == ==
// == AMDiS - Adaptive multidimensional simulations ==
// == ==
// ============================================================================
// == ==
// == TU Dresden ==
// == ==
// == Institut für Wissenschaftliches Rechnen ==
// == Zellescher Weg 12-14 ==
// == 01069 Dresden ==
// == germany ==
// == ==
// ============================================================================
// == ==
// == https://gforge.zih.tu-dresden.de/projects/amdis/ ==
// == ==
// ============================================================================
/** \file PITL_Solver.h */
#ifndef AMDIS_PITL_SOLVER_H
#define AMDIS_PITL_SOLVER_H
#include "ITL_Solver.h"
namespace AMDiS {
namespace Parallel {
typedef mtl::matrix::distributed< MTLTypes::MTLMatrix > MTLMatrix;
typedef mtl::vector::distributed< MTLTypes::MTLVector > MTLVector;
typedef CreatorInterface< ITL_BasePreconditioner< MTLVector > > ParallelPreconditionCreator;
typedef ITL_BasePreconditioner< MTLVector > ParallelPreconditioner;
template< typename ITLSolver >
class PITL_Solver : public AMDiS::MTL4Solver< MTLMatrix, MTLVector, ITL_OEMSolver_runner< ITLSolver > > {
public:
PITL_Solver(std::string name):
AMDiS::MTL4Solver< MTLMatrix, MTLVector, ITL_OEMSolver_runner< ITLSolver > >(name) {}
int solveSystem(const SolverMatrix<Matrix<DOFMatrix*> >& A,
SystemVector& x,
SystemVector& b,
ParallelMapper& m)
{
return AMDiS::MTL4Solver< MTLMatrix, MTLVector, ITL_OEMSolver_runner< ITLSolver > >::solve(A,x,b,m);
}
/// Creator class used in the OEMSolverMap.
class Creator : public OEMSolverCreator
{
public:
virtual ~Creator() {}
/// Returns a new CGSolver object.
OEMSolver* create()
{
return new PITL_Solver<ITLSolver>(this->name);
}
};
};
template< typename ITLSolver >
class PITL_Solver_para : public AMDiS::MTL4Solver< MTLMatrix, MTLVector, ITL_OEMSolver_para_runner< ITLSolver > > {
public:
PITL_Solver_para(std::string name):
AMDiS::MTL4Solver< MTLMatrix, MTLVector, ITL_OEMSolver_para_runner< ITLSolver > >(name) {}
/// Creator class used in the OEMSolverMap.
class Creator : public OEMSolverCreator
{
public:
virtual ~Creator() {}
/// Returns a new CGSolver object.
OEMSolver* create()
{
return new PITL_Solver_para<ITLSolver>(this->name);
}
};
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by the conjugate gradient method (CG) and can be used for
* symmetric positive definite system matrices.
* Right preconditioner is ignored.
*/
class CGSolver : public PITL_Solver<cg_solver_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
CGSolver(std::string name) : PITL_Solver<cg_solver_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by the squared conjugate gradient method (CGS).
* Right preconditioner is ignored.
*/
class CGSSolver : public PITL_Solver<cgs_solver_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
CGSSolver(std::string name) : PITL_Solver<cgs_solver_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by a stabilized BiCG method and can be used for
* system matrices.
*/
class BiCGSolver : public PITL_Solver<bicg_solver_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
BiCGSolver(std::string name) : PITL_Solver<bicg_solver_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by a stabilized BiCG method and can be used for
* system matrices.
*/
class BiCGStabSolver : public PITL_Solver<bicgstab_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
BiCGStabSolver(std::string name) : PITL_Solver<bicgstab_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by a stabilized BiCG method and can be used for
* system matrices.
*/
class BiCGStab2Solver : public PITL_Solver<bicgstab2_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
BiCGStab2Solver(std::string name) : PITL_Solver<bicgstab2_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by the Quasi-Minimal Residual method (QMR).
*/
class QMRSolver : public PITL_Solver<qmr_solver_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
QMRSolver(std::string name) : PITL_Solver<qmr_solver_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by the Transposed-Free Quasi-Minimal Residual method (TFQMR).
* Does not use preconditioning currently.
*/
class TFQMRSolver : public PITL_Solver<tfqmr_solver_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
TFQMRSolver(std::string name) : PITL_Solver<tfqmr_solver_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by a stabilized BiCG(l) method and can be used for
* system matrices.
*/
class BiCGStabEllSolver : public PITL_Solver_para<bicgstab_ell_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
BiCGStabEllSolver(std::string name) : PITL_Solver_para<bicgstab_ell_type>(name) {}
};
/**
* \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 PITL_Solver_para<gmres_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
GMResSolver(std::string name) : PITL_Solver_para<gmres_type>(name) {}
};
/**
* \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 PITL_Solver_para<idr_s_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
IDRsSolver(std::string name) : PITL_Solver_para<idr_s_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Solves a linear system by the Minres method. Can be used for symmetric
* indefinite systems.
*/
class MinResSolver : public PITL_Solver<minres_solver_type>
{
public:
/// The constructor reads required parameters and sets solvers \ref name.
MinResSolver(std::string name) : PITL_Solver<minres_solver_type>(name) {}
};
/**
* \ingroup Solver
*
* \brief
* Diagonal preconditioner.
*/
class DiagonalPreconditioner
: public ITL_Preconditioner<itl::pc::diagonal<MTLMatrix >, MTLVector, MTLMatrix >
{
typedef ITL_Preconditioner<itl::pc::diagonal<MTLMatrix >, MTLVector, MTLMatrix > base;
public:
DiagonalPreconditioner(const MTLMatrix& A) : base(A) {}
};
/**
* \ingroup Solver
*
* \brief
* Identity preconditioner. Behaves like no preconditioning.
*/
class IdentityPreconditioner
: public ITL_Preconditioner<itl::pc::identity<MTLMatrix >, MTLVector, MTLMatrix >
{
typedef ITL_Preconditioner<itl::pc::identity<MTLMatrix >, MTLVector, MTLMatrix > base;
public:
IdentityPreconditioner(const MTLMatrix& 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<MTLMatrix >,MTLVector, MTLMatrix >
{};
class ICPreconditioner
: public ITL_Preconditioner< itl::pc::ic_0<MTLMatrix >, MTLVector, MTLMatrix >
{};
}
template< >
struct Collection< Parallel::MTLMatrix > {
typedef mtl::matrix::inserter< Parallel::MTLMatrix , mtl::update_plus< Parallel::MTLMatrix::value_type > > Inserter;
};
template< >
struct Collection< Parallel::MTLVector > {
typedef mtl::vector::inserter< Parallel::MTLVector , mtl::update_plus< Parallel::MTLVector::value_type > > Inserter;
typedef Parallel::MTLMatrix PreconditionMatrix;
};
} // namespace AMDiS
#endif // AMDIS_PITL_SOLVER