FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
pcr.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<
31 typename Matrix_,
32 typename Filter_>
33 class PCR :
34 public PreconditionedIterativeSolver<typename Matrix_::VectorTypeR>
35 {
36 public:
37 typedef Matrix_ MatrixType;
38 typedef Filter_ FilterType;
39 typedef typename MatrixType::VectorTypeR VectorType;
40 typedef typename MatrixType::DataType DataType;
42
44
45 protected:
47 const MatrixType& _system_matrix;
49 const FilterType& _system_filter;
51 VectorType _vec_r, _vec_s, _vec_p, _vec_q, _vec_t;
52
53 public:
66 explicit PCR(const MatrixType& matrix, const FilterType& filter,
67 std::shared_ptr<PrecondType> precond = nullptr) :
68 BaseClass("PCR", precond),
69 _system_matrix(matrix),
70 _system_filter(filter)
71 {
72 // set communicator by system matrix
73 this->_set_comm_by_matrix(matrix);
74 }
75
95 explicit PCR(const String& section_name, const PropertyMap* section,
96 const MatrixType& matrix, const FilterType& filter,
97 std::shared_ptr<PrecondType> precond = nullptr) :
98 BaseClass("PCR", section_name, section, precond),
99 _system_matrix(matrix),
100 _system_filter(filter)
101 {
102 // set communicator by system matrix
103 this->_set_comm_by_matrix(matrix);
104 }
105
106 virtual String name() const override
107 {
108 return "PCR";
109 }
110
111 virtual void init_symbolic() override
112 {
114 // create five temporary vectors
115 _vec_p = this->_system_matrix.create_vector_r();
116 _vec_q = this->_system_matrix.create_vector_r();
117 _vec_r = this->_system_matrix.create_vector_r();
118 _vec_s = this->_system_matrix.create_vector_r();
119 _vec_t = this->_system_matrix.create_vector_r();
120 }
121
122 virtual void done_symbolic() override
123 {
124 this->_vec_t.clear();
125 this->_vec_s.clear();
126 this->_vec_r.clear();
127 this->_vec_q.clear();
128 this->_vec_p.clear();
130 }
131
132 virtual Status apply(VectorType& vec_cor, const VectorType& vec_def) override
133 {
134 // save defect
135 this->_vec_r.copy(vec_def);
136
137 // clear solution vector
138 vec_cor.format();
139
140 // apply solver
141 this->_status = _apply_intern(vec_cor);
142
143 // plot summary
144 this->plot_summary();
145
146 // return status
147 return this->_status;
148 }
149
150 virtual Status correct(VectorType& vec_sol, const VectorType& vec_rhs) override
151 {
152 // compute defect
153 this->_system_matrix.apply(this->_vec_r, vec_sol, vec_rhs, -DataType(1));
154 this->_system_filter.filter_def(this->_vec_r);
155
156 // apply solver
157 this->_status = _apply_intern(vec_sol);
158
159 // plot summary
160 this->plot_summary();
161
162 // return status
163 return this->_status;
164 }
165
166 protected:
167 virtual Status _apply_intern(VectorType& vec_sol)
168 {
169 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->name()));
170 const MatrixType& matrix(this->_system_matrix);
171 const FilterType& filter(this->_system_filter);
172 VectorType& vec_p(this->_vec_p);
173 VectorType& vec_q(this->_vec_q);
174 VectorType& vec_r(this->_vec_r);
175 VectorType& vec_s(this->_vec_s);
176 // Note: y and z are temporary vectors whose usage does
177 // not overlap, so we use the same vector for them
178 VectorType& vec_y(this->_vec_t);
179 VectorType& vec_z(this->_vec_t);
180
181 // Note:
182 // In the algorithm below, the following relations hold:
183 // q[k] = A * p[k]
184 // s[k] = M^{-1} * r[k]
185 // y[k] = A * s[k]
186 // z[k] = M^{-1} * q[k]
187
188 // set initial defect:
189 // r[0] := b - A*x[0]
190 Status status = this->_set_initial_defect(vec_r, vec_sol);
191 if(status != Status::progress)
192 {
193 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
194 return status;
195 }
196
197 // apply preconditioner to defect vector:
198 // s[0] := M^{-1} * r[0]
199 if(!this->_apply_precond(vec_s, vec_r, filter))
200 {
201 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
202 return Status::aborted;
203 }
204
205 // compute initial p and q
206 // p[0] := s[0]
207 vec_p.copy(vec_s);
208
209 // q[0] := A*p[0]
210 matrix.apply(vec_q, vec_p);
211 filter.filter_def(vec_q);
212
213 // compute initial gamma:
214 // gamma[0] := < p[0], q[0] >
215 DataType gamma = vec_p.dot(vec_q);
216
217 // start iterating
218 while(status == Status::progress)
219 {
220 IterationStats stat(*this);
221
222 // z[k] := M^{-1} * q[k]
223 if(!this->_apply_precond(vec_z, vec_q, filter))
224 {
225 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
226 return Status::aborted;
227 }
228
229 // compute alpha
230 // alpha[k] := gamma[k] / < z[k], q[k] >
231 DataType alpha = gamma / vec_z.dot(vec_q);
232
233 // update solution vector
234 // x[k+1] := x[k] + alpha[k] * p[k]
235 vec_sol.axpy(vec_p, alpha);
236
237 // update real residual vector
238 // r[k+1] := r[k] - alpha[k] * q[k]
239 vec_r.axpy(vec_q, -alpha);
240
241 // compute defect norm and check for convergence
242 status = this->_set_new_defect(vec_r, vec_sol);
243 if(status != Status::progress)
244 {
245 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
246 return status;
247 }
248
249 // update preconditioned residual vector
250 // s[k+1] := s[k] - alpha[k] * z[k]
251 vec_s.axpy(vec_z, -alpha);
252
253 // y[k+1] := A*s[k+1]
254 matrix.apply(vec_y, vec_s);
255 filter.filter_def(vec_y);
256
257 // compute new gamma:
258 // gamma[k+1] := < s[k+1], y[k+1] >
259 DataType gamma2 = gamma;
260 gamma = vec_s.dot(vec_y);
261
262 // compute beta:
263 // beta[k] := gamma[k+1] / gamma[k]
264 DataType beta = gamma / gamma2;
265
266 // update direction vectors:
267 // p[k+1] := s[k+1] + beta[k] * p[k]
268 vec_p.scale(vec_p, beta);
269 vec_p.axpy(vec_s);
270
271 // q[k+1] := y[k+1] + beta[k] * q[k]
272 vec_q.scale(vec_q, beta);
273 vec_q.axpy(vec_y);
274 }
275
276 // we should never reach this point...
277 {
278 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::undefined, this->get_num_iter()));
279 return Status::undefined;
280 }
281 }
282 }; // class PCR<...>
283
299 template<typename Matrix_, typename Filter_>
300 inline std::shared_ptr<PCR<Matrix_, Filter_>> new_pcr(
301 const Matrix_& matrix, const Filter_& filter,
302 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
303 {
304 return std::make_shared<PCR<Matrix_, Filter_>>(matrix, filter, precond);
305 }
306
328 template<typename Matrix_, typename Filter_>
329 inline std::shared_ptr<PCR<Matrix_, Filter_>> new_pcr(
330 const String& section_name, const PropertyMap* section,
331 const Matrix_& matrix, const Filter_& filter,
332 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
333 {
334 return std::make_shared<PCR<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
335 }
336 } // namespace Solver
337} // namespace FEAT
A class organizing a tree of key-value pairs.
Helper class for iteration statistics collection.
Definition: base.hpp:392
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-Residual solver implementation
Definition: pcr.hpp:35
virtual void done_symbolic() override
Symbolic finalization method.
Definition: pcr.hpp:122
VectorType _vec_r
vectors
Definition: pcr.hpp:51
virtual String name() const override
Returns a descriptive string.
Definition: pcr.hpp:106
PCR(const String &section_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor using a PropertyMap.
Definition: pcr.hpp:95
virtual void init_symbolic() override
Symbolic initialization method.
Definition: pcr.hpp:111
virtual Status _apply_intern(VectorType &vec_sol)
Definition: pcr.hpp:167
const MatrixType & _system_matrix
the matrix for the solver
Definition: pcr.hpp:47
PCR(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
Definition: pcr.hpp:66
const FilterType & _system_filter
the filter for the solver
Definition: pcr.hpp:49
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< PCR< Matrix_, Filter_ > > new_pcr(const Matrix_ &matrix, const Filter_ &filter, std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond=nullptr)
Creates a new PCR solver object.
Definition: pcr.hpp:300
FEAT namespace.
Definition: adjactor.hpp:12