33 #ifndef __H__UG__CPU_ALGEBRA__CORE_SMOOTHERS__
34 #define __H__UG__CPU_ALGEBRA__CORE_SMOOTHERS__
52 template<
typename Matrix_type,
typename Vector_type>
53 bool gs_step_LL(
const Matrix_type &A, Vector_type &x,
const Vector_type &b)
60 for(
size_t i=0; i < x.size(); i++)
64 for(
typename Matrix_type::const_row_iterator it = A.begin_row(i); it != A.end_row(i) && it.index() < i; ++it)
84 template<
typename Matrix_type,
typename Vector_type>
85 bool gs_step_UR(
const Matrix_type &A, Vector_type &x,
const Vector_type &b)
91 if(x.size() == 0)
return true;
92 size_t i = x.size()-1;
96 typename Matrix_type::const_row_iterator diag = A.get_connection(i, i);
98 typename Matrix_type::const_row_iterator it = diag; ++it;
99 for(; it != A.end_row(i); ++it)
119 template<
typename Matrix_type,
typename Vector_type>
120 bool sgs_step(
const Matrix_type &A, Vector_type &x,
const Vector_type &b)
129 for(
size_t i = 0; i<x.size(); i++)
130 MatMult(x[i], 1.0, A(i, i), x[i]);
147 template<
typename Matrix_type,
typename Vector_type>
148 bool diag_step(
const Matrix_type& A, Vector_type& x,
const Vector_type& b,
number damp)
151 UG_ASSERT(x.size() == b.size() && x.size() == A.num_rows(), x <<
", " << b <<
" and " << A <<
" need to have same size.");
153 for(
size_t i=0; i < x.size(); i++)
void gs_step_UR(const Matrix_type &A, Vector_type &c, const Vector_type &d, const number relaxFactor)
Performs a backward gauss-seidel-step, that is, solve on the upper right of A. Using gs_step_UR withi...
Definition: core_smoothers.h:150
void gs_step_LL(const Matrix_type &A, Vector_type &c, const Vector_type &d, const number relaxFactor)
Gauss-Seidel-Iterations.
Definition: core_smoothers.h:106
void diag_step(const Matrix_type &A, Vector_type &c, const Vector_type &d, number damp)
Performs a jacobi-step.
Definition: core_smoothers.h:218
void sgs_step(const Matrix_type &A, Vector_type &c, const Vector_type &d, const number relaxFactor)
Performs a symmetric gauss-seidel step. Using sgs_step within a preconditioner-scheme leads to the fa...
Definition: core_smoothers.h:189
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
double number
Definition: types.h:124
bool MatMultAdd(vector_t &dest, const number &alpha1, const vector_t &v1, const number &beta1, const matrix_t &A1, const vector_t &w1)
calculates dest = alpha1*v1 + beta1 * A1 *w1;
Definition: operations_mat.h:68
bool MatMult(vector_t &dest, const number &beta1, const matrix_t &A1, const vector_t &w1)
calculates dest = beta1 * A1;
Definition: operations_mat.h:59
bool InverseMatMult(number &dest, const double &beta, const TMat &mat, const TVec &vec)
you can implement this function with GetInverse and MatMult
T value_type
Definition: sparsematrix_interface.h:2