9#include <kernel/solver/iterative.hpp> 
   42      typedef Matrix_ MatrixType;
 
   43      typedef Filter_ FilterType;
 
   44      typedef typename MatrixType::VectorTypeR VectorType;
 
   45      typedef typename MatrixType::DataType DataType;
 
   56      VectorType 
_vec_r, _vec_u, _vec_w, _vec_z, _vec_q, _vec_s, _vec_p, _vec_m, _vec_n;
 
   71      explicit PipePCG(
const MatrixType& matrix, 
const FilterType& filter,
 
   72        std::shared_ptr<PrecondType> precond = 
nullptr) :
 
  101        const MatrixType& matrix, 
const FilterType& filter,
 
  102        std::shared_ptr<PrecondType> precond = 
nullptr) :
 
  103        BaseClass(
"PipePCG", section_name, section, precond),
 
  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();
 
  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();
 
  145      virtual Status apply(VectorType& vec_cor, 
const VectorType& vec_def)
 override 
  148        this->_vec_r.copy(vec_def);
 
  163      virtual Status correct(VectorType& vec_sol, 
const VectorType& vec_rhs)
 override 
  166        this->_system_matrix.apply(this->_vec_r, vec_sol, vec_rhs, -DataType(1));
 
  167        this->_system_filter.filter_def(this->_vec_r);
 
  183        Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->
name()));
 
  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);
 
  197        DataType gamma, gamma_old, delta, beta, alpha;
 
  199        gamma_old = DataType(0);
 
  208          Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->
name(), status, this->
get_num_iter()));
 
  219        matrix.apply(vec_w, vec_u);
 
  220        filter.filter_def(vec_w);
 
  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);
 
  240          matrix.apply(vec_n, vec_m);
 
  241          filter.filter_def(vec_n);
 
  243          gamma = dot_gamma.wait();
 
  244          delta = dot_delta.wait();
 
  250            Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->
name(), status, this->
get_num_iter()));
 
  256            alpha = gamma / delta;
 
  264            beta = gamma / gamma_old;
 
  265            alpha = gamma / (delta - beta / alpha * gamma);
 
  266            vec_z.scale(vec_z, beta);
 
  268            vec_q.scale(vec_q, beta);
 
  270            vec_p.scale(vec_p, beta);
 
  272            vec_s.scale(vec_s, beta);
 
  276          vec_sol.axpy(vec_p, alpha);
 
  290            vec_u.axpy(vec_q, -alpha);
 
  291            vec_w.axpy(vec_z, -alpha);
 
  292            vec_r.axpy(vec_s, -alpha);
 
  320    template<
typename Matrix_, 
typename Filter_>
 
  322      const Matrix_& matrix, 
const Filter_& filter,
 
  325      return std::make_shared<PipePCG<Matrix_, Filter_>>(matrix, filter, precond);
 
  349    template<
typename Matrix_, 
typename Filter_>
 
  352      const Matrix_& matrix, 
const Filter_& filter,
 
  355      return std::make_shared<PipePCG<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
 
A class organizing a tree of key-value pairs.
Helper class for iteration statistics collection.
void destroy()
destroy the objects contents (and generate Statistics::expression) before the actual destructor call
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 _update_defect(const DataType def_cur_norm)
Internal function: sets the new (next) defect norm.
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.
(Preconditioned) pipelined Conjugate-Gradient solver implementation from Ghysels and Vnaroose
const FilterType & _system_filter
the filter for the solver
virtual String name() const override
Returns a descriptive string.
virtual void init_symbolic() override
Symbolic initialization method.
PipePCG(const String §ion_name, const PropertyMap *section, const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor using a PropertyMap.
PipePCG(const MatrixType &matrix, const FilterType &filter, std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
virtual Status _apply_intern(VectorType &vec_sol)
const MatrixType & _system_matrix
the matrix for the solver
virtual void done_symbolic() override
Symbolic finalization method.
VectorType _vec_r
temporary vectors
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.
Polymorphic solver interface.
String class implementation.
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::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.