9#include <kernel/solver/iterative.hpp>
39 typedef Matrix_ MatrixType;
40 typedef Filter_ FilterType;
41 typedef typename MatrixType::VectorTypeR VectorType;
42 typedef typename MatrixType::DataType DataType;
53 VectorType
_vec_r, _vec_p_hat, _vec_q_hat;
55 std::vector<VectorType>
p_list, q_list;
70 explicit RGCR(
const MatrixType& matrix,
const FilterType& filter,
71 std::shared_ptr<PrecondType> precond =
nullptr) :
100 const MatrixType& matrix,
const FilterType& filter,
101 std::shared_ptr<PrecondType> precond =
nullptr) :
102 BaseClass(
"RGCR", section_name, section, precond),
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();
125 this->_vec_q_hat.clear();
126 this->_vec_p_hat.clear();
127 this->_vec_r.clear();
133 virtual Status apply(VectorType& vec_cor,
const VectorType& vec_def)
override
136 this->_vec_r.copy(vec_def);
143 this->
_status = _apply_intern(vec_cor);
148 p_list.resize(
p_list.size() / 4);
149 q_list.resize(q_list.size() / 4);
155 virtual Status correct(VectorType& vec_sol,
const VectorType& vec_rhs)
override
158 this->_system_matrix.apply(this->_vec_r, vec_sol, vec_rhs, -DataType(1));
159 this->_system_filter.filter_def(this->_vec_r);
162 this->
_status = _apply_intern(vec_sol);
167 p_list.resize(
p_list.size() / 4);
168 q_list.resize(q_list.size() / 4);
175 virtual Status _apply_intern(VectorType& vec_sol)
177 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->
name()));
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);
190 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->
name(), status, this->
get_num_iter()));
197 IterationStats stat(*
this);
210 matrix.apply(vec_q_hat, vec_p_hat);
211 filter.filter_def(vec_q_hat);
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);
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());
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);
232 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->
name(), status, this->
get_num_iter()));
258 template<
typename Matrix_,
typename Filter_>
259 inline std::shared_ptr<RGCR<Matrix_, Filter_>>
new_rgcr(
260 const Matrix_& matrix,
const Filter_& filter,
263 return std::make_shared<RGCR<Matrix_, Filter_>>(matrix, filter, precond);
287 template<
typename Matrix_,
typename Filter_>
288 inline std::shared_ptr<RGCR<Matrix_, Filter_>>
new_rgcr(
290 const Matrix_& matrix,
const Filter_& filter,
293 return std::make_shared<RGCR<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
A class organizing a tree of key-value pairs.
Index get_num_iter() const
Returns number of performed iterations.
Status _status
current status of the solver
void _set_comm_by_matrix(const Matrix_ &matrix)
Sets the communicator for the solver from a matrix.
virtual Status _set_new_defect(const VectorType &vec_def, const VectorType &vec_sol)
Internal function: sets the new (next) defect vector.
Index _num_iter
number of performed iterations
virtual void plot_summary() const
Plot a summary of the last solver run.
virtual Status _set_initial_defect(const VectorType &vec_def, const VectorType &vec_sol)
Internal function: sets the initial defect vector.
Abstract base-class for preconditioned iterative solvers.
virtual void init_symbolic() override
Symbolic initialization method.
bool _apply_precond(VectorType &vec_cor, const VectorType &vec_def, const Filter_ &filter)
Applies the preconditioner onto a defect vector.
virtual void done_symbolic() override
Symbolic finalization method.
(Preconditioned) Recycling Generalized Conjugate Residual solver implementation
const FilterType & _system_filter
the filter for the solver
RGCR(const String §ion_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor using a PropertyMap.
VectorType _vec_r
temporary vectors
const MatrixType & _system_matrix
the matrix for the solver
virtual void done_symbolic() override
Symbolic finalization method.
std::vector< VectorType > p_list
descent vectors for later recycling
RGCR(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
virtual String name() const override
Returns a descriptive string.
virtual void init_symbolic() override
Symbolic initialization method.
Polymorphic solver interface.
String class implementation.
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.
Status
Solver status return codes enumeration.
@ progress
continue iteration (internal use only)
@ undefined
undefined status
@ aborted
premature abort (solver aborted due to internal errors or preconditioner failure)
std::uint64_t Index
Index data type.