FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
pcgnr.hpp
1// FEAT3: Finite Element Analysis Toolbox, Version 3
2// Copyright (C) 2010 by Stefan Turek & the FEAT group
3// FEAT3 is released under the GNU General Public License version 3,
4// see the file 'copyright.txt' in the top level directory for details.
5
6#pragma once
7
8// includes, FEAT
9#include <kernel/solver/iterative.hpp>
10
11namespace FEAT
12{
13 namespace Solver
14 {
41 template<
42 typename Matrix_,
43 typename Filter_>
44 class PCGNR :
45 public IterativeSolver<typename Matrix_::VectorTypeR>
46 {
47 public:
48 typedef Matrix_ MatrixType;
49 typedef Filter_ FilterType;
50 typedef typename MatrixType::VectorTypeR VectorType;
51 typedef typename MatrixType::DataType DataType;
53
55
56 protected:
58 const MatrixType& _system_matrix;
60 const FilterType& _system_filter;
62 std::shared_ptr<PrecondType> _precond_l;
64 std::shared_ptr<PrecondType> _precond_r;
66 MatrixType _transp_matrix;
68 VectorType _vec_r, _vec_p, _vec_q, _vec_s, _vec_t;
69
70 public:
89 explicit PCGNR(const MatrixType& matrix, const FilterType& filter,
90 std::shared_ptr<PrecondType> precond_l = nullptr,
91 std::shared_ptr<PrecondType> precond_r = nullptr) :
92 BaseClass("PCGNR"),
93 _system_matrix(matrix),
94 _system_filter(filter),
95 _precond_l(precond_l),
96 _precond_r(precond_r)
97 {
98 // set communicator by system matrix
99 this->_set_comm_by_matrix(matrix);
100 }
101
120 explicit PCGNR(const String& section_name, const PropertyMap* section,
121 const MatrixType& matrix, const FilterType& filter,
122 std::shared_ptr<PrecondType> precond_l = nullptr,
123 std::shared_ptr<PrecondType> precond_r = nullptr) :
124 BaseClass("PCGNR", section_name, section),
125 _system_matrix(matrix),
126 _system_filter(filter),
127 _precond_l(precond_l),
128 _precond_r(precond_r)
129 {
130 // set communicator by system matrix
131 this->_set_comm_by_matrix(matrix);
132 }
133
134 virtual String name() const override
135 {
136 return "PCGNR";
137 }
138
139 virtual void init_symbolic() override
140 {
142
143 // initialize preconditioners
144 if(_precond_l)
145 _precond_l->init_symbolic();
146 if((_precond_r) && (_precond_r != _precond_l))
147 _precond_r->init_symbolic();
148
149 // create five temporary vectors
150 _vec_p = this->_system_matrix.create_vector_r();
151 _vec_q = this->_system_matrix.create_vector_r();
152 _vec_r = this->_system_matrix.create_vector_r();
153 _vec_s = this->_system_matrix.create_vector_r();
154 _vec_t = this->_system_matrix.create_vector_r();
155 }
156
157 virtual void done_symbolic() override
158 {
159 this->_vec_t.clear();
160 this->_vec_s.clear();
161 this->_vec_r.clear();
162 this->_vec_q.clear();
163 this->_vec_p.clear();
164
165 // release preconditioners
166 if((_precond_r) && (_precond_r != _precond_l))
167 _precond_r->done_symbolic();
168 if(_precond_l)
169 _precond_l->done_symbolic();
170
172 }
173
174 virtual void init_numeric() override
175 {
177
178 // initialize preconditioners
179 if(_precond_l)
180 _precond_l->init_numeric();
181 if((_precond_r) && (_precond_r != _precond_l))
182 _precond_r->init_numeric();
183
184 // transpose matrix
185 this->_transp_matrix = this->_system_matrix.transpose();
186 }
187
188 virtual void done_numeric() override
189 {
190 // release transposed matrix
191 this->_transp_matrix.clear();
192
193 // release preconditioners
194 if((_precond_r) && (_precond_r != _precond_l))
195 _precond_r->done_numeric();
196 if(_precond_l)
197 _precond_l->done_numeric();
198
200 }
201
202 virtual Status apply(VectorType& vec_cor, const VectorType& vec_def) override
203 {
204 // save defect
205 this->_vec_r.copy(vec_def);
206
207 // clear solution vector
208 vec_cor.format();
209
210 // apply solver
211 this->_status = _apply_intern(vec_cor);
212
213 // plot summary
214 this->plot_summary();
215
216 // return status
217 return this->_status;
218 }
219
220 virtual Status correct(VectorType& vec_sol, const VectorType& vec_rhs) override
221 {
222 // compute defect
223 this->_system_matrix.apply(this->_vec_r, vec_sol, vec_rhs, -DataType(1));
224 this->_system_filter.filter_def(this->_vec_r);
225
226 // apply solver
227 this->_status = _apply_intern(vec_sol);
228
229 // plot summary
230 this->plot_summary();
231
232 // return status
233 return this->_status;
234 }
235
236 protected:
237 bool _apply_precond_l(VectorType& vec_cor, const VectorType& vec_def)
238 {
239 if(_precond_l)
240 {
241 Statistics::add_solver_expression(std::make_shared<ExpressionCallPrecond>(this->name(), this->_precond_l->name()));
242 return status_success(_precond_l->apply(vec_cor, vec_def));
243 }
244 vec_cor.copy(vec_def);
245 this->_system_filter.filter_cor(vec_cor);
246 return true;
247 }
248
249 bool _apply_precond_r(VectorType& vec_cor, const VectorType& vec_def)
250 {
251 if(_precond_r)
252 {
253 Statistics::add_solver_expression(std::make_shared<ExpressionCallPrecond>(this->name(), this->_precond_r->name()));
254 return status_success(_precond_r->apply(vec_cor, vec_def));
255 }
256 vec_cor.copy(vec_def);
257 this->_system_filter.filter_cor(vec_cor);
258 return true;
259 }
260
261 virtual Status _apply_intern(VectorType& vec_sol)
262 {
263 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->name()));
264
265 const MatrixType& matrix(this->_system_matrix);
266 const MatrixType& transp(this->_transp_matrix);
267 const FilterType& filter(this->_system_filter);
268
269 VectorType& vec_p(this->_vec_p);
270 VectorType& vec_q(this->_vec_q);
271 VectorType& vec_r(this->_vec_r);
272 VectorType& vec_s(this->_vec_s);
273 VectorType& vec_t(this->_vec_t);
274 VectorType& vec_y(this->_vec_s); // same as s
275 VectorType& vec_z(this->_vec_t); // same as t
276
277 // set initial defect
278 // r[0] := b - A*x[0]
279 Status status = this->_set_initial_defect(vec_r, vec_sol);
280 if(status != Status::progress)
281 {
282 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
283 return status;
284 }
285
286 // apply left preconditioner to defect vector
287 // p[0] := M_L^{-1} * r[0]
288 if(!this->_apply_precond_l(vec_p, vec_r))
289 {
290 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
291 return Status::aborted;
292 }
293
294 // apply transposed matrix
295 // s[0] := A^T * p[0]
296 transp.apply(vec_s, vec_p);
297
298 // apply right preconditioner
299 // q[0] := M_R^{-1} * s[0]
300 if(!this->_apply_precond_r(vec_q, vec_s))
301 {
302 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
303 return Status::aborted;
304 }
305
306 // compute initial gamma
307 // gamma[0] := < s[0], q[0] >
308 DataType gamma = vec_s.dot(vec_q);
309
310 // start iterating
311 while(status == Status::progress)
312 {
313 IterationStats stat(*this);
314
315 // y[k] := A * q[k]
316 matrix.apply(vec_y, vec_q);
317 filter.filter_def(vec_y);
318
319 // z[k] := M_L^{-1} * y[k]
320 if(!this->_apply_precond_l(vec_z, vec_y))
321 {
322 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
323 return Status::aborted;
324 }
325
326 // compute alpha
327 // alpha[k] := gamma[k] / < y[k], z[k] >
328 DataType alpha = gamma / vec_y.dot(vec_z);
329
330 // update solution vector
331 // x[k+1] := x[k] + alpha[k] * q[k]
332 vec_sol.axpy(vec_q, alpha);
333
334 // update defect vector
335 // r[k+1] := r[k] - alpha[k] * y[k]
336 vec_r.axpy(vec_y, -alpha);
337
338 // compute defect norm and check for convergence
339 status = this->_set_new_defect(vec_r, vec_sol);
340 if(status != Status::progress)
341 {
342 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
343 return status;
344 }
345
346 // p[k+1] := p[k] - alpha[k] * z[k]
347 vec_p.axpy(vec_z, -alpha);
348
349 // s[k+1] := A^T * p[k+1]
350 transp.apply(vec_s, vec_p);
351 filter.filter_def(vec_s);
352
353 // t[k+1] := M_R^{-1} * s[k+1]
354 if(!this->_apply_precond_r(vec_t, vec_s))
355 {
356 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
357 return Status::aborted;
358 }
359
360 // compute new gamma
361 // gamma[k+1] := < s[k+1], t[k+1] >
362 DataType gamma2 = gamma;
363 gamma = vec_s.dot(vec_t);
364
365 // compute beta
366 // beta[k] := gamma[k+1] / gamma[k]
367 DataType beta = gamma / gamma2;
368
369 // update direction vector
370 // q[k+1] := t[k+1] + beta * q[k]
371 vec_q.scale(vec_q, beta);
372 vec_q.axpy(vec_t);
373 }
374
375 // we should never reach this point...
376 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::undefined, this->get_num_iter()));
377 return Status::undefined;
378 }
379 }; // class PCGNR<...>
380
399 template<typename Matrix_, typename Filter_>
400 inline std::shared_ptr<PCGNR<Matrix_, Filter_>> new_pcgnr(
401 const Matrix_& matrix, const Filter_& filter,
402 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond_l = nullptr,
403 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond_r = nullptr)
404 {
405 return std::make_shared<PCGNR<Matrix_, Filter_>>(matrix, filter, precond_l, precond_r);
406 }
407
432 template<typename Matrix_, typename Filter_>
433 inline std::shared_ptr<PCGNR<Matrix_, Filter_>> new_pcgnr(
434 const String& section_name, const PropertyMap* section,
435 const Matrix_& matrix, const Filter_& filter,
436 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond_l = nullptr,
437 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond_r = nullptr)
438 {
439 return std::make_shared<PCGNR<Matrix_, Filter_>>(section_name, section, matrix, filter, precond_l, precond_r);
440 }
441 } // namespace Solver
442} // namespace FEAT
A class organizing a tree of key-value pairs.
Helper class for iteration statistics collection.
Definition: base.hpp:392
Abstract base-class for iterative solvers.
Definition: iterative.hpp:198
Index get_num_iter() const
Returns number of performed iterations.
Definition: iterative.hpp:462
Status _status
current status of the solver
Definition: iterative.hpp:213
void _set_comm_by_matrix(const Matrix_ &matrix)
Sets the communicator for the solver from a matrix.
Definition: iterative.hpp:680
virtual Status _set_new_defect(const VectorType &vec_def, const VectorType &vec_sol)
Internal function: sets the new (next) defect vector.
Definition: iterative.hpp:924
virtual void plot_summary() const
Plot a summary of the last solver run.
Definition: iterative.hpp:627
virtual Status _set_initial_defect(const VectorType &vec_def, const VectorType &vec_sol)
Internal function: sets the initial defect vector.
Definition: iterative.hpp:802
(Preconditioned) Conjugate-Gradient on Normal Equations solver implementation
Definition: pcgnr.hpp:46
virtual void done_symbolic() override
Symbolic finalization method.
Definition: pcgnr.hpp:157
PCGNR(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond_l=nullptr, std::shared_ptr< PrecondType > precond_r=nullptr)
Constructor.
Definition: pcgnr.hpp:89
virtual void init_symbolic() override
Symbolic initialization method.
Definition: pcgnr.hpp:139
MatrixType _transp_matrix
the transposed system matrix
Definition: pcgnr.hpp:66
VectorType _vec_r
temporary vectors
Definition: pcgnr.hpp:68
virtual void init_numeric() override
Numeric initialization method.
Definition: pcgnr.hpp:174
const FilterType & _system_filter
the filter for the solver
Definition: pcgnr.hpp:60
virtual void done_numeric() override
Numeric finalization method.
Definition: pcgnr.hpp:188
std::shared_ptr< PrecondType > _precond_l
left preconditioner
Definition: pcgnr.hpp:62
const MatrixType & _system_matrix
the matrix for the solver
Definition: pcgnr.hpp:58
virtual String name() const override
Returns a descriptive string.
Definition: pcgnr.hpp:134
PCGNR(const String &section_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond_l=nullptr, std::shared_ptr< PrecondType > precond_r=nullptr)
Constructor using a PropertyMap.
Definition: pcgnr.hpp:120
std::shared_ptr< PrecondType > _precond_r
right preconditioner
Definition: pcgnr.hpp:64
virtual Status _apply_intern(VectorType &vec_sol)
Definition: pcgnr.hpp:261
Polymorphic solver interface.
Definition: base.hpp:183
virtual void init_symbolic()
Symbolic initialization method.
Definition: base.hpp:227
virtual void init_numeric()
Numeric initialization method.
Definition: base.hpp:237
virtual void done_symbolic()
Symbolic finalization method.
Definition: base.hpp:255
virtual void done_numeric()
Numeric finalization method.
Definition: base.hpp:246
String class implementation.
Definition: string.hpp:47
bool status_success(Status status)
Status success check function.
Definition: base.hpp:108
Status
Solver status return codes enumeration.
Definition: base.hpp:47
@ progress
continue iteration (internal use only)
@ undefined
undefined status
@ aborted
premature abort (solver aborted due to internal errors or preconditioner failure)
std::shared_ptr< PCGNR< Matrix_, Filter_ > > new_pcgnr(const Matrix_ &matrix, const Filter_ &filter, std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond_l=nullptr, std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond_r=nullptr)
Creates a new PCGNR solver object.
Definition: pcgnr.hpp:400
FEAT namespace.
Definition: adjactor.hpp:12