FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
pmr.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 PMR :
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_s, _vec_q, _vec_z;
57
58 public:
71 explicit PMR(const MatrixType& matrix, const FilterType& filter,
72 std::shared_ptr<PrecondType> precond = nullptr) :
73 BaseClass("PMR", 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
102 explicit PMR(const String& section_name, const PropertyMap* section,
103 const MatrixType& matrix, const FilterType& filter,
104 std::shared_ptr<PrecondType> precond = nullptr) :
105 BaseClass("PMR", section_name, section, precond),
106 _system_matrix(matrix),
107 _system_filter(filter)
108 {
109 // set communicator by system matrix
110 this->_set_comm_by_matrix(matrix);
111 }
112
113 virtual String name() const override
114 {
115 return "PMR";
116 }
117
118 virtual void init_symbolic() override
119 {
121
122 _vec_r = this->_system_matrix.create_vector_r();
123 _vec_s = this->_system_matrix.create_vector_r();
124 _vec_q = this->_system_matrix.create_vector_r();
125 _vec_z = this->_system_matrix.create_vector_r();
126 }
127
128 virtual void done_symbolic() override
129 {
130 this->_vec_z.clear();
131 this->_vec_q.clear();
132 this->_vec_s.clear();
133 this->_vec_r.clear();
135 }
136
137 virtual Status apply(VectorType& vec_cor, const VectorType& vec_def) override
138 {
139 // save defect
140 this->_vec_r.copy(vec_def);
141
142 // clear solution vector
143 vec_cor.format();
144
145 // apply solver
146 this->_status = _apply_intern(vec_cor);
147
148 // plot summary
149 this->plot_summary();
150
151 // return status
152 return this->_status;
153 }
154
155 virtual Status correct(VectorType& vec_sol, const VectorType& vec_rhs) override
156 {
157 // compute defect
158 this->_system_matrix.apply(this->_vec_r, vec_sol, vec_rhs, -DataType(1));
159 this->_system_filter.filter_def(this->_vec_r);
160
161 // apply solver
162 this->_status = _apply_intern(vec_sol);
163
164 // plot summary
165 this->plot_summary();
166
167 // return status
168 return this->_status;
169 }
170
171 protected:
172 virtual Status _apply_intern(VectorType& vec_sol)
173 {
174 IterationStats pre_iter(*this);
175 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->name()));
176
177 const MatrixType& matrix(this->_system_matrix);
178 const FilterType& filter(this->_system_filter);
179 VectorType& vec_r(this->_vec_r);
180 VectorType& vec_s(this->_vec_s);
181 VectorType& vec_q(this->_vec_q);
182 VectorType& vec_z(this->_vec_z);
183
184 // set initial defect:
185 // r[0] := b - A*x[0]
186 Status status = this->_set_initial_defect(vec_r, vec_sol);
187 if(status != Status::progress)
188 {
189 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
190 return status;
191 }
192
193 // apply preconditioner to defect vector
194 // s[0] := M^{-1} * r[0]
195 if(!this->_apply_precond(vec_s, vec_r, filter))
196 {
197 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
198 return Status::aborted;
199 }
200 pre_iter.destroy();
201
202 // start iterating
203 while(status == Status::progress)
204 {
205 IterationStats stat(*this);
206
207 // q[k] := A*s[k]
208 matrix.apply(vec_q, vec_s);
209 filter.filter_def(vec_q);
210
211 // apply preconditioner to defect vector
212 // z[k] := M^{-1} * q[k]
213 if(!this->_apply_precond(vec_z, vec_q, filter))
214 {
215 stat.destroy();
216 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
217 return Status::aborted;
218 }
219
220 // alpha[k] := < q[k], s[k] > / < z[k], q[k] >
221 DataType alpha = vec_q.dot(vec_s) / vec_z.dot(vec_q);
222
223 // update solution vector:
224 // x[k+1] := x[k] + alpha[k] * s[k]
225 vec_sol.axpy(vec_s, alpha);
226
227 // update defect vector:
228 // r[k+1] := r[k] - alpha[k] * q[k]
229 vec_r.axpy(vec_q, -alpha);
230
231 // compute defect norm
232 status = this->_set_new_defect(vec_r, vec_sol);
233 if(status != Status::progress)
234 {
235 stat.destroy();
236 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
237 return status;
238 }
239
240 // update preconditioned defect:
241 // s[k+1] := s[k] - alpha[k] * z[k]
242 vec_s.axpy(vec_z, -alpha);
243 }
244
245 // we should never reach this point...
246 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::undefined, this->get_num_iter()));
247 return Status::undefined;
248 }
249 }; // class PMR<...>
250
266 template<typename Matrix_, typename Filter_>
267 inline std::shared_ptr<PMR<Matrix_, Filter_>> new_pmr(
268 const Matrix_& matrix, const Filter_& filter,
269 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
270 {
271 return std::make_shared<PMR<Matrix_, Filter_>>(matrix, filter, precond);
272 }
273
295 template<typename Matrix_, typename Filter_>
296 inline std::shared_ptr<PMR<Matrix_, Filter_>> new_pmr(
297 const String& section_name, const PropertyMap* section,
298 const Matrix_& matrix, const Filter_& filter,
299 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
300 {
301 return std::make_shared<PMR<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
302 }
303 } // namespace Solver
304} // namespace FEAT
A class organizing a tree of key-value pairs.
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) Minimal-Residual-Iteration solver implementation
Definition: pmr.hpp:40
VectorType _vec_r
temporary vectors
Definition: pmr.hpp:56
const FilterType & _system_filter
the filter for the solver
Definition: pmr.hpp:54
virtual void done_symbolic() override
Symbolic finalization method.
Definition: pmr.hpp:128
virtual void init_symbolic() override
Symbolic initialization method.
Definition: pmr.hpp:118
PMR(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
Definition: pmr.hpp:71
const MatrixType & _system_matrix
the matrix for the solver
Definition: pmr.hpp:52
PMR(const String &section_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor using a PropertyMap.
Definition: pmr.hpp:102
virtual String name() const override
Returns a descriptive string.
Definition: pmr.hpp:113
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< Descent< Matrix_, Filter_ > > new_pmr(const String &section_name, const PropertyMap *section, const Matrix_ &matrix, const Filter_ &filter, std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond=nullptr)
Creates a new Descent solver object using a PropertyMap.
Definition: descent.hpp:364
FEAT namespace.
Definition: adjactor.hpp:12