FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
pipepcg.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 {
35 template<
36 typename Matrix_,
37 typename Filter_>
38 class PipePCG :
39 public PreconditionedIterativeSolver<typename Matrix_::VectorTypeR>
40 {
41 public:
42 typedef Matrix_ MatrixType;
43 typedef Filter_ FilterType;
44 typedef typename MatrixType::VectorTypeR VectorType;
45 typedef typename MatrixType::DataType DataType;
47
49
50 protected:
52 const MatrixType& _system_matrix;
54 const FilterType& _system_filter;
56 VectorType _vec_r, _vec_u, _vec_w, _vec_z, _vec_q, _vec_s, _vec_p, _vec_m, _vec_n;
57
58 public:
71 explicit PipePCG(const MatrixType& matrix, const FilterType& filter,
72 std::shared_ptr<PrecondType> precond = nullptr) :
73 BaseClass("PipePCG", precond),
74 _system_matrix(matrix),
75 _system_filter(filter)
76 {
77 // set communicator by system matrix
78 this->_set_comm_by_matrix(matrix);
79 }
80
100 explicit PipePCG(const String& section_name, const PropertyMap* section,
101 const MatrixType& matrix, const FilterType& filter,
102 std::shared_ptr<PrecondType> precond = nullptr) :
103 BaseClass("PipePCG", section_name, section, precond),
104 _system_matrix(matrix),
105 _system_filter(filter)
106 {
107 // set communicator by system matrix
108 this->_set_comm_by_matrix(matrix);
109 }
110
111 virtual String name() const override
112 {
113 return "PipePCG";
114 }
115
116 virtual void init_symbolic() override
117 {
119 // create temporary vectors
120 _vec_r = this->_system_matrix.create_vector_r();
121 _vec_u = this->_system_matrix.create_vector_r();
122 _vec_w = this->_system_matrix.create_vector_r();
123 _vec_z = this->_system_matrix.create_vector_r();
124 _vec_q = this->_system_matrix.create_vector_r();
125 _vec_s = this->_system_matrix.create_vector_r();
126 _vec_p = this->_system_matrix.create_vector_r();
127 _vec_m = this->_system_matrix.create_vector_r();
128 _vec_n = this->_system_matrix.create_vector_r();
129 }
130
131 virtual void done_symbolic() override
132 {
133 this->_vec_r.clear();
134 this->_vec_u.clear();
135 this->_vec_w.clear();
136 this->_vec_z.clear();
137 this->_vec_q.clear();
138 this->_vec_s.clear();
139 this->_vec_p.clear();
140 this->_vec_m.clear();
141 this->_vec_n.clear();
143 }
144
145 virtual Status apply(VectorType& vec_cor, const VectorType& vec_def) override
146 {
147 // save defect
148 this->_vec_r.copy(vec_def);
149
150 // clear solution vector
151 vec_cor.format();
152
153 // apply solver
154 this->_status = _apply_intern(vec_cor);
155
156 // plot summary
157 this->plot_summary();
158
159 // return status
160 return this->_status;
161 }
162
163 virtual Status correct(VectorType& vec_sol, const VectorType& vec_rhs) override
164 {
165 // compute defect
166 this->_system_matrix.apply(this->_vec_r, vec_sol, vec_rhs, -DataType(1));
167 this->_system_filter.filter_def(this->_vec_r);
168
169 // apply solver
170 this->_status = _apply_intern(vec_sol);
171
172 // plot summary
173 this->plot_summary();
174
175 // return status
176 return this->_status;
177 }
178
179 protected:
180 virtual Status _apply_intern(VectorType& vec_sol)
181 {
182 IterationStats pre_iter(*this);
183 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->name()));
184
185 const MatrixType& matrix(this->_system_matrix);
186 const FilterType& filter(this->_system_filter);
187 VectorType& vec_r(this->_vec_r);
188 VectorType& vec_u(this->_vec_u);
189 VectorType& vec_w(this->_vec_w);
190 VectorType& vec_z(this->_vec_z);
191 VectorType& vec_q(this->_vec_q);
192 VectorType& vec_s(this->_vec_s);
193 VectorType& vec_p(this->_vec_p);
194 VectorType& vec_m(this->_vec_m);
195 VectorType& vec_n(this->_vec_n);
196
197 DataType gamma, gamma_old, delta, beta, alpha;
198 gamma = DataType(0);
199 gamma_old = DataType(0);
200 delta = DataType(0);
201 beta = DataType(0);
202 alpha = DataType(0);
203
204 Status status = this->_set_initial_defect(vec_r, vec_sol);
205 if(status != Status::progress)
206 {
207 pre_iter.destroy();
208 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
209 return status;
210 }
211
212 if(!this->_apply_precond(vec_u, vec_r, filter))
213 {
214 pre_iter.destroy();
215 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
216 return Status::aborted;
217 }
218
219 matrix.apply(vec_w, vec_u);
220 filter.filter_def(vec_w);
221
222 pre_iter.destroy();
223
224 // start iterating
225 while(status == Status::progress)
226 {
227 IterationStats stat(*this);
228
229 auto norm_def_cur = vec_r.norm2_async();
230 auto dot_gamma = vec_r.dot_async(vec_u);
231 auto dot_delta = vec_w.dot_async(vec_u);
232
233 if(!this->_apply_precond(vec_m, vec_w, filter))
234 {
235 stat.destroy();
236 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
237 return Status::aborted;
238 }
239
240 matrix.apply(vec_n, vec_m);
241 filter.filter_def(vec_n);
242
243 gamma = dot_gamma.wait();
244 delta = dot_delta.wait();
245
246 status = this->_update_defect(norm_def_cur.wait());
247 if(status != Status::progress)
248 {
249 stat.destroy();
250 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
251 return status;
252 }
253
254 if (this->_num_iter == 1) // num_iter has already been increased to 1 by previous _update_defect call
255 {
256 alpha = gamma / delta;
257 vec_z.copy(vec_n);
258 vec_q.copy(vec_m);
259 vec_p.copy(vec_u);
260 vec_s.copy(vec_w);
261 }
262 else
263 {
264 beta = gamma / gamma_old;
265 alpha = gamma / (delta - beta / alpha * gamma);
266 vec_z.scale(vec_z, beta);
267 vec_z.axpy(vec_n);
268 vec_q.scale(vec_q, beta);
269 vec_q.axpy(vec_m);
270 vec_p.scale(vec_p, beta);
271 vec_p.axpy(vec_u);
272 vec_s.scale(vec_s, beta);
273 vec_s.axpy(vec_w);
274 }
275
276 vec_sol.axpy(vec_p, alpha);
277
278 /* Stabilization, recalculate real values
279 if(this->_num_iter > 40 && this->_num_iter % 50 == 0)
280 {
281 matrix.apply(vec_r, vec_sol, DataType(-1));
282 filter.filter_def(vec_r);
283 if(!this->_apply_precond(vec_u, vec_r, filter))
284 return Status::aborted;
285 matrix.apply(vec_w, vec_u);
286 filter.filter_def(vec_w);
287 }
288 else */
289 {
290 vec_u.axpy(vec_q, -alpha);
291 vec_w.axpy(vec_z, -alpha);
292 vec_r.axpy(vec_s, -alpha);
293 }
294
295 gamma_old = gamma;
296 }
297
298 // we should never reach this point...
299 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::undefined, this->get_num_iter()));
300 return Status::undefined;
301 }
302
303 }; // class PipePCG<...>
304
320 template<typename Matrix_, typename Filter_>
321 inline std::shared_ptr<PipePCG<Matrix_, Filter_>> new_pipepcg(
322 const Matrix_& matrix, const Filter_& filter,
323 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
324 {
325 return std::make_shared<PipePCG<Matrix_, Filter_>>(matrix, filter, precond);
326 }
327
349 template<typename Matrix_, typename Filter_>
350 inline std::shared_ptr<PipePCG<Matrix_, Filter_>> new_pipepcg(
351 const String& section_name, const PropertyMap* section,
352 const Matrix_& matrix, const Filter_& filter,
353 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
354 {
355 return std::make_shared<PipePCG<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
356 }
357 } // namespace Solver
358} // 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 _update_defect(const DataType def_cur_norm)
Internal function: sets the new (next) defect norm.
Definition: iterative.hpp:970
Index _num_iter
number of performed iterations
Definition: iterative.hpp:231
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) pipelined Conjugate-Gradient solver implementation from Ghysels and Vnaroose
Definition: pipepcg.hpp:40
const FilterType & _system_filter
the filter for the solver
Definition: pipepcg.hpp:54
virtual String name() const override
Returns a descriptive string.
Definition: pipepcg.hpp:111
virtual void init_symbolic() override
Symbolic initialization method.
Definition: pipepcg.hpp:116
PipePCG(const String &section_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor using a PropertyMap.
Definition: pipepcg.hpp:100
PipePCG(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
Definition: pipepcg.hpp:71
virtual Status _apply_intern(VectorType &vec_sol)
Definition: pipepcg.hpp:180
const MatrixType & _system_matrix
the matrix for the solver
Definition: pipepcg.hpp:52
virtual void done_symbolic() override
Symbolic finalization method.
Definition: pipepcg.hpp:131
VectorType _vec_r
temporary vectors
Definition: pipepcg.hpp:56
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:46
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< PipePCG< Matrix_, Filter_ > > new_pipepcg(const Matrix_ &matrix, const Filter_ &filter, std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond=nullptr)
Creates a new PipePCG solver object.
Definition: pipepcg.hpp:321
FEAT namespace.
Definition: adjactor.hpp:12