FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
psd.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 PSD :
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_q, _vec_z;
52
53 public:
66 explicit PSD(const MatrixType& matrix, const FilterType& filter,
67 std::shared_ptr<PrecondType> precond = nullptr) :
68 BaseClass("PSD", 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
94 explicit PSD(const String& section_name, const PropertyMap* section,
95 const MatrixType& matrix, const FilterType& filter,
96 std::shared_ptr<PrecondType> precond = nullptr) :
97 BaseClass("PSD", section_name, section, precond),
98 _system_matrix(matrix),
99 _system_filter(filter)
100 {
101 // set communicator by system matrix
102 this->_set_comm_by_matrix(matrix);
103 }
104
105 virtual String name() const override
106 {
107 return "PSD";
108 }
109
110 virtual void init_symbolic() override
111 {
113 // create three temporary vectors
114 _vec_r = this->_system_matrix.create_vector_r();
115 _vec_q = this->_system_matrix.create_vector_r();
116 _vec_z = this->_system_matrix.create_vector_r();
117 }
118
119 virtual void done_symbolic() override
120 {
121 this->_vec_z.clear();
122 this->_vec_q.clear();
123 this->_vec_r.clear();
125 }
126
127 virtual Status apply(VectorType& vec_cor, const VectorType& vec_def) override
128 {
129 // save defect
130 this->_vec_r.copy(vec_def);
131
132 // clear solution vector
133 vec_cor.format();
134
135 // apply solver
136 this->_status = _apply_intern(vec_cor);
137
138 // plot summary
139 this->plot_summary();
140
141 // return status
142 return this->_status;
143 }
144
145 virtual Status correct(VectorType& vec_sol, const VectorType& vec_rhs) override
146 {
147 // compute defect
148 this->_system_matrix.apply(this->_vec_r, vec_sol, vec_rhs, -DataType(1));
149 this->_system_filter.filter_def(this->_vec_r);
150
151 // apply solver
152 this->_status = _apply_intern(vec_sol);
153
154 // plot summary
155 this->plot_summary();
156
157 // return status
158 return this->_status;
159 }
160
161 protected:
162 virtual Status _apply_intern(VectorType& vec_sol)
163 {
164 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->name()));
165
166 const MatrixType& matrix(this->_system_matrix);
167 const FilterType& filter(this->_system_filter);
168 VectorType& vec_r(this->_vec_r);
169 VectorType& vec_q(this->_vec_q);
170 VectorType& vec_z(this->_vec_z);
171
172 // set initial defect:
173 // r[0] := b - A*x[0]
174 Status status = this->_set_initial_defect(vec_r, vec_sol);
175 if(status != Status::progress)
176 {
177 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
178 return status;
179 }
180
181 // start iterating
182 while(status == Status::progress)
183 {
184 IterationStats stat(*this);
185
186 // apply preconditioner to defect vector
187 // z[k] := M^{-1} * r[k]
188 if(!this->_apply_precond(vec_z, vec_r, filter))
189 {
190 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
191 return Status::aborted;
192 }
193
194 // q[k] := A*z[k]
195 matrix.apply(vec_q, vec_z);
196 filter.filter_def(vec_q);
197
198 // alpha[k] := < z[k], r[k] > / < z[k], q[k] >
199 DataType alpha = vec_z.dot(vec_r) / vec_z.dot(vec_q);
200
201 // update solution vector:
202 // x[k+1] := x[k] + alpha[k] * z[k]
203 vec_sol.axpy(vec_z, alpha);
204
205 // update defect vector:
206 // r[k+1] := r[k] - alpha[k] * q[k]
207 vec_r.axpy(vec_q, -alpha);
208
209 // compute defect norm
210 status = this->_set_new_defect(vec_r, vec_sol);
211 if(status != Status::progress)
212 {
213 stat.destroy();
214 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
215 return status;
216 }
217 }
218
219 // we should never reach this point...
220 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::undefined, this->get_num_iter()));
221 return Status::undefined;
222 }
223 }; // class PSD<...>
224
240 template<typename Matrix_, typename Filter_>
241 inline std::shared_ptr<PSD<Matrix_, Filter_>> new_psd(
242 const Matrix_& matrix, const Filter_& filter,
243 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
244 {
245 return std::make_shared<PSD<Matrix_, Filter_>>(matrix, filter, precond);
246 }
247
269 template<typename Matrix_, typename Filter_>
270 inline std::shared_ptr<PSD<Matrix_, Filter_>> new_psd(
271 const String& section_name, const PropertyMap* section,
272 const Matrix_& matrix, const Filter_& filter,
273 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
274 {
275 return std::make_shared<PSD<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
276 }
277 } // namespace Solver
278} // 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) Steepest-Descent solver implementation
Definition: psd.hpp:35
const FilterType & _system_filter
the filter for the solver
Definition: psd.hpp:49
virtual void init_symbolic() override
Symbolic initialization method.
Definition: psd.hpp:110
virtual String name() const override
Returns a descriptive string.
Definition: psd.hpp:105
PSD(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
Definition: psd.hpp:66
const MatrixType & _system_matrix
the matrix for the solver
Definition: psd.hpp:47
virtual void done_symbolic() override
Symbolic finalization method.
Definition: psd.hpp:119
PSD(const String &section_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor using a PropertyMap.
Definition: psd.hpp:94
VectorType _vec_r
temporary vectors
Definition: psd.hpp:51
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< PSD< Matrix_, Filter_ > > new_psd(const Matrix_ &matrix, const Filter_ &filter, std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond=nullptr)
Creates a new PSD solver object.
Definition: psd.hpp:241
FEAT namespace.
Definition: adjactor.hpp:12