FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
pcg.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 {
30 template<typename Matrix_, typename Filter_>
31 class PCG :
32 public PreconditionedIterativeSolver<typename Matrix_::VectorTypeR>
33 {
34 public:
36 typedef Matrix_ MatrixType;
38 typedef Filter_ FilterType;
40 typedef typename MatrixType::VectorTypeR VectorType;
42 typedef typename MatrixType::DataType DataType;
47
48 protected:
59
60 public:
73 explicit PCG(const MatrixType& matrix, const FilterType& filter,
74 std::shared_ptr<PrecondType> precond = nullptr) :
75 BaseClass("PCG", precond),
76 _system_matrix(matrix),
77 _system_filter(filter)
78 {
79 // set communicator by system matrix
80 this->_set_comm_by_matrix(matrix);
81 }
82
102 explicit PCG(const String& section_name, const PropertyMap* section,
103 const MatrixType& matrix, const FilterType& filter, std::shared_ptr<PrecondType> precond = nullptr) :
104 BaseClass("PCG", section_name, section, precond),
105 _system_matrix(matrix),
106 _system_filter(filter)
107 {
108 // set communicator by system matrix
109 this->_set_comm_by_matrix(matrix);
110 }
111
113 virtual String name() const override
114 {
115 return "PCG";
116 }
117
119 virtual void init_symbolic() override
120 {
122 // create three temporary vectors
123 _vec_r = this->_system_matrix.create_vector_r();
124 _vec_p = this->_system_matrix.create_vector_r();
125 _vec_t = this->_system_matrix.create_vector_r();
126 }
127
129 virtual void done_symbolic() override
130 {
131 this->_vec_t.clear();
132 this->_vec_p.clear();
133 this->_vec_r.clear();
135 }
136
138 virtual Status apply(VectorType& vec_cor, const VectorType& vec_def) override
139 {
140 // save defect
141 this->_vec_r.copy(vec_def);
142
143 // clear solution vector
144 vec_cor.format();
145
146 // apply solver
147 this->_status = _apply_intern(vec_cor);
148
149 // plot summary
150 this->plot_summary();
151
152 // return status
153 return this->_status;
154 }
155
157 virtual Status correct(VectorType& vec_sol, const VectorType& vec_rhs) override
158 {
159 // compute defect
160 this->_system_matrix.apply(this->_vec_r, vec_sol, vec_rhs, -DataType(1));
161 this->_system_filter.filter_def(this->_vec_r);
162
163 // apply solver
164 this->_status = _apply_intern(vec_sol);
165
166 // plot summary
167 this->plot_summary();
168
169 // return status
170 return this->_status;
171 }
172
173 protected:
187 {
188 IterationStats pre_iter(*this);
189 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->name()));
190
191 const MatrixType& matrix(this->_system_matrix);
192 const FilterType& filter(this->_system_filter);
193 VectorType& vec_r(this->_vec_r);
194 VectorType& vec_p(this->_vec_p);
195 // Note: q and z are temporary vectors whose usage does
196 // not overlap, so we use the same vector for them
197 VectorType& vec_q(this->_vec_t);
198 VectorType& vec_z(this->_vec_t);
199
200 // Note:
201 // In the algorithm below, the following relations hold:
202 // q[k] = A * p[k]
203 // z[k] = M^{-1} * r[k]
204
205 // set initial defect:
206 // r[0] := b - A*x[0]
207 Status status = this->_set_initial_defect(vec_r, vec_sol);
208 if(status != Status::progress)
209 {
210 pre_iter.destroy();
211 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
212 return status;
213 }
214
215 // apply preconditioner to defect vector
216 // p[0] := M^{-1} * r[0]
217 if(!this->_apply_precond(vec_p, vec_r, filter))
218 {
219 pre_iter.destroy();
220 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
221 return Status::aborted;
222 }
223
224 // compute initial gamma:
225 // gamma[0] := < r[0], p[0] >
226 DataType gamma = vec_r.dot(vec_p);
227
228 pre_iter.destroy();
229
230 // start iterating
231 while(status == Status::progress)
232 {
233 IterationStats stat(*this);
234
235 // q[k] := A*p[k]
236 matrix.apply(vec_q, vec_p);
237 filter.filter_def(vec_q);
238
239 // compute alpha
240 // alpha[k] := gamma[k] / < q[k], p[k] >
241 DataType alpha = gamma / vec_q.dot(vec_p);
242
243 // update solution vector:
244 // x[k+1] := x[k] + alpha[k] * p[k]
245 vec_sol.axpy(vec_p, alpha);
246
247 // update defect vector:
248 // r[k+1] := r[k] - alpha[k] * q[k]
249 vec_r.axpy(vec_q, -alpha);
250
251 // compute defect norm
252 status = this->_set_new_defect(vec_r, vec_sol);
253 if(status != Status::progress)
254 {
255 stat.destroy();
256 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
257 return status;
258 }
259
260 // apply preconditioner
261 // z[k+1] := M^{-1} * r[k+1]
262 if(!this->_apply_precond(vec_z, vec_r, filter))
263 {
264 stat.destroy();
265 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
266 return Status::aborted;
267 }
268
269 // compute new gamma:
270 // gamma[k+1] := < r[k+1] , z[k+1] >
271 DataType gamma2 = gamma;
272 gamma = vec_r.dot(vec_z);
273
274 // compute beta:
275 // beta[k] := gamma[k+1] / gamma[k]
276 DataType beta = gamma / gamma2;
277
278 // update direction vector:
279 // p[k+1] := z[k+1] + beta[k] * p[k]
280 vec_p.scale(vec_p, beta);
281 vec_p.axpy(vec_z);
282 }
283
284 // we should never reach this point...
285 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::undefined, this->get_num_iter()));
286 return Status::undefined;
287 }
288 }; // class PCG<...>
289
305 template<typename Matrix_, typename Filter_>
306 inline std::shared_ptr<PCG<Matrix_, Filter_>> new_pcg(
307 const Matrix_& matrix, const Filter_& filter,
308 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
309 {
310 return std::make_shared<PCG<Matrix_, Filter_>>(matrix, filter, precond);
311 }
312
334 template<typename Matrix_, typename Filter_>
335 inline std::shared_ptr<PCG<Matrix_, Filter_>> new_pcg(
336 const String& section_name, const PropertyMap* section,
337 const Matrix_& matrix, const Filter_& filter,
338 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
339 {
340 return std::make_shared<PCG<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
341 }
342 } // namespace Solver
343} // namespace FEAT
A class organizing a tree of key-value pairs.
Helper class for iteration statistics collection.
Definition: base.hpp:392
void destroy()
destroy the objects contents (and generate Statistics::expression) before the actual destructor call
Definition: base.hpp:464
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 solver implementation
Definition: pcg.hpp:33
virtual Status apply(VectorType &vec_cor, const VectorType &vec_def) override
Solver application method.
Definition: pcg.hpp:138
PCG(const String &section_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor using a PropertyMap.
Definition: pcg.hpp:102
MatrixType::VectorTypeR VectorType
The vector type this solver can be applied to.
Definition: pcg.hpp:40
VectorType _vec_t
Temporary vector, used for e.g. the preconditioned defect.
Definition: pcg.hpp:58
virtual void init_symbolic() override
Symbolic initialization method.
Definition: pcg.hpp:119
MatrixType::DataType DataType
The floating point precision.
Definition: pcg.hpp:42
Filter_ FilterType
The filter for projecting solution, rhs, defect and correction vectors to subspaces.
Definition: pcg.hpp:38
PCG(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
Definition: pcg.hpp:73
PreconditionedIterativeSolver< VectorType > BaseClass
Our base class.
Definition: pcg.hpp:44
VectorType _vec_r
The defect vector.
Definition: pcg.hpp:54
const FilterType & _system_filter
the filter for the solver
Definition: pcg.hpp:52
SolverBase< VectorType > PrecondType
The type of the preconditioner that can be used.
Definition: pcg.hpp:46
const MatrixType & _system_matrix
the matrix for the solver
Definition: pcg.hpp:50
virtual Status _apply_intern(VectorType &vec_sol)
Internal function, applies the solver.
Definition: pcg.hpp:186
virtual void done_symbolic() override
Symbolic finalization method.
Definition: pcg.hpp:129
virtual Status correct(VectorType &vec_sol, const VectorType &vec_rhs) override
Solver correction method.
Definition: pcg.hpp:157
virtual String name() const override
Returns a descriptive string.
Definition: pcg.hpp:113
VectorType _vec_p
The update (or search) direction.
Definition: pcg.hpp:56
Matrix_ MatrixType
The type of matrix this solver can be applied to.
Definition: pcg.hpp:36
Abstract base-class for preconditioned iterative solvers.
Definition: iterative.hpp:1027
virtual void init_symbolic() override
Symbolic initialization method.
Definition: iterative.hpp:1086
bool _apply_precond(VectorType &vec_cor, const VectorType &vec_def, const Filter_ &filter)
Applies the preconditioner onto a defect vector.
Definition: iterative.hpp:1140
virtual void done_symbolic() override
Symbolic finalization method.
Definition: iterative.hpp:1110
Polymorphic solver interface.
Definition: base.hpp:183
String class implementation.
Definition: string.hpp:47
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< PCG< Matrix_, Filter_ > > new_pcg(const Matrix_ &matrix, const Filter_ &filter, std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond=nullptr)
Creates a new PCG solver object.
Definition: pcg.hpp:306
FEAT namespace.
Definition: adjactor.hpp:12