FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
rgcr.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 {
32 template<
33 typename Matrix_,
34 typename Filter_>
35 class RGCR :
36 public PreconditionedIterativeSolver<typename Matrix_::VectorTypeR>
37 {
38 public:
39 typedef Matrix_ MatrixType;
40 typedef Filter_ FilterType;
41 typedef typename MatrixType::VectorTypeR VectorType;
42 typedef typename MatrixType::DataType DataType;
44
46
47 protected:
49 const MatrixType& _system_matrix;
51 const FilterType& _system_filter;
53 VectorType _vec_r, _vec_p_hat, _vec_q_hat;
55 std::vector<VectorType> p_list, q_list;
56
57 public:
70 explicit RGCR(const MatrixType& matrix, const FilterType& filter,
71 std::shared_ptr<PrecondType> precond = nullptr) :
72 BaseClass("RGCR", precond),
73 _system_matrix(matrix),
74 _system_filter(filter)
75 {
76 // set communicator by system matrix
77 this->_set_comm_by_matrix(matrix);
78 }
79
99 explicit RGCR(const String& section_name, const PropertyMap* section,
100 const MatrixType& matrix, const FilterType& filter,
101 std::shared_ptr<PrecondType> precond = nullptr) :
102 BaseClass("RGCR", section_name, section, precond),
103 _system_matrix(matrix),
104 _system_filter(filter)
105 {
106 // set communicator by system matrix
107 this->_set_comm_by_matrix(matrix);
108 }
109
110 virtual String name() const override
111 {
112 return "RGCR";
113 }
114
115 virtual void init_symbolic() override
116 {
118 _vec_r = this->_system_matrix.create_vector_r();
119 _vec_p_hat = this->_system_matrix.create_vector_r();
120 _vec_q_hat = this->_system_matrix.create_vector_r();
121 }
122
123 virtual void done_symbolic() override
124 {
125 this->_vec_q_hat.clear();
126 this->_vec_p_hat.clear();
127 this->_vec_r.clear();
128 q_list.clear();
129 p_list.clear();
131 }
132
133 virtual Status apply(VectorType& vec_cor, const VectorType& vec_def) override
134 {
135 // save defect
136 this->_vec_r.copy(vec_def);
137 //this->_system_filter.filter_def(this->_vec_r);
138
139 // clear solution vector
140 vec_cor.format();
141
142 // apply solver
143 this->_status = _apply_intern(vec_cor);
144
145 // plot summary
146 this->plot_summary();
147
148 p_list.resize(p_list.size() / 4);
149 q_list.resize(q_list.size() / 4);
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 p_list.resize(p_list.size() / 4);
168 q_list.resize(q_list.size() / 4);
169
170 // return status
171 return this->_status;
172 }
173
174 protected:
175 virtual Status _apply_intern(VectorType& vec_sol)
176 {
177 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->name()));
178
179 const MatrixType& matrix(this->_system_matrix);
180 const FilterType& filter(this->_system_filter);
181 VectorType& vec_r(this->_vec_r);
182 VectorType& vec_p_hat(this->_vec_p_hat);
183 VectorType& vec_q_hat(this->_vec_q_hat);
184
185 // set initial defect:
186 // r[0] := b - A*x[0]
187 Status status = this->_set_initial_defect(vec_r, vec_sol);
188 if(status != Status::progress)
189 {
190 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
191 return status;
192 }
193
194 // start iterating
195 while(status == Status::progress)
196 {
197 IterationStats stat(*this);
198
199 // no stored descent vectors left, calculate new one
200 if (this->_num_iter >= p_list.size())
201 {
202 // apply preconditioner to defect vector
203 if(!this->_apply_precond(vec_p_hat, vec_r, filter))
204 {
205 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
206 return Status::aborted;
207 }
208
209 //vec_p_hat.copy(vec_r);
210 matrix.apply(vec_q_hat, vec_p_hat);
211 filter.filter_def(vec_q_hat);
212 for (Index j(0) ; j < this->_num_iter ; ++j)
213 {
214 auto beta = vec_q_hat.dot(q_list.at(j));
215 vec_q_hat.axpy(q_list.at(j), -beta);
216 vec_p_hat.axpy(p_list.at(j), -beta);
217 }
218 const auto gamma = DataType(1) / vec_q_hat.norm2();
219 vec_q_hat.scale(vec_q_hat, gamma);
220 q_list.push_back(vec_q_hat.clone());
221 vec_p_hat.scale(vec_p_hat, gamma);
222 p_list.push_back(vec_p_hat.clone());
223 }
224
225 auto alpha = vec_r.dot(q_list.at(this->_num_iter));
226 vec_sol.axpy(p_list.at(this->_num_iter), alpha);
227 vec_r.axpy(q_list.at(this->_num_iter), -alpha);
228 // compute defect norm
229 status = this->_set_new_defect(vec_r, vec_sol);
230 if(status != Status::progress)
231 {
232 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
233 return status;
234 }
235 }
236
237 // we should never reach this point...
238 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::undefined, this->get_num_iter()));
239 return Status::undefined;
240 }
241 }; // class RGCR<...>
242
258 template<typename Matrix_, typename Filter_>
259 inline std::shared_ptr<RGCR<Matrix_, Filter_>> new_rgcr(
260 const Matrix_& matrix, const Filter_& filter,
261 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
262 {
263 return std::make_shared<RGCR<Matrix_, Filter_>>(matrix, filter, precond);
264 }
265
287 template<typename Matrix_, typename Filter_>
288 inline std::shared_ptr<RGCR<Matrix_, Filter_>> new_rgcr(
289 const String& section_name, const PropertyMap* section,
290 const Matrix_& matrix, const Filter_& filter,
291 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
292 {
293 return std::make_shared<RGCR<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
294 }
295 } // namespace Solver
296} // 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
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
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
(Preconditioned) Recycling Generalized Conjugate Residual solver implementation
Definition: rgcr.hpp:37
const FilterType & _system_filter
the filter for the solver
Definition: rgcr.hpp:51
RGCR(const String &section_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor using a PropertyMap.
Definition: rgcr.hpp:99
VectorType _vec_r
temporary vectors
Definition: rgcr.hpp:53
const MatrixType & _system_matrix
the matrix for the solver
Definition: rgcr.hpp:49
virtual void done_symbolic() override
Symbolic finalization method.
Definition: rgcr.hpp:123
std::vector< VectorType > p_list
descent vectors for later recycling
Definition: rgcr.hpp:55
RGCR(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
Definition: rgcr.hpp:70
virtual String name() const override
Returns a descriptive string.
Definition: rgcr.hpp:110
virtual void init_symbolic() override
Symbolic initialization method.
Definition: rgcr.hpp:115
Polymorphic solver interface.
Definition: base.hpp:183
String class implementation.
Definition: string.hpp:47
std::shared_ptr< RGCR< Matrix_, Filter_ > > new_rgcr(const Matrix_ &matrix, const Filter_ &filter, std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond=nullptr)
Creates a new RGCR solver object.
Definition: rgcr.hpp:259
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)
FEAT namespace.
Definition: adjactor.hpp:12
std::uint64_t Index
Index data type.