33 #ifndef __H__UG__LIB_ALGEBRA__OPERATOR__PRECONDITIONER__PROJECTED_GAUSS_SEIDEL__PROJ_GAUSS_SEIDEL_IMPL__
34 #define __H__UG__LIB_ALGEBRA__OPERATOR__PRECONDITIONER__PROJECTED_GAUSS_SEIDEL__PROJ_GAUSS_SEIDEL_IMPL__
41 template<
typename Matrix_type,
typename Vector_type>
43 const size_t i,
const number relaxFactor)
47 for(
typename Matrix_type::const_row_iterator it = A.begin_row(i);
48 it != A.end_row(i) && it.index() < i; ++it)
56 template<
typename Matrix_type,
typename Vector_type>
58 const size_t i,
const number relaxFactor)
62 typename Matrix_type::const_row_iterator diag = A.get_connection(i, i);
63 typename Matrix_type::const_row_iterator it = diag; ++it;
65 for(; it != A.end_row(i); ++it)
74 template <
typename TDomain,
typename TAlgebra>
79 for(
size_t i = 0; i < c.size(); i++)
84 this->project_correction(c[i], i);
88 template <
typename TDomain,
typename TAlgebra>
93 if(c.size() == 0)
return;
94 size_t i = c.size()-1;
100 this->project_correction(c[i], i);
105 template <
typename TDomain,
typename TAlgebra>
110 for(
size_t i = 0; i < c.size(); i++)
117 MatMult(c[i], 1.0, A(i, i), c[i]);
124 this->project_correction(c[i], i);
TAlgebra::vector_type vector_type
type of vector
Definition: debug_writer.h:360
TAlgebra::matrix_type matrix_type
type of matrix
Definition: debug_writer.h:363
virtual void step(const matrix_type &mat, vector_type &c, const vector_type &d, const number relax)
computes a new correction c = B*d and projects on the underlying constraint
Definition: proj_gauss_seidel_impl.h:91
virtual void step(const matrix_type &mat, vector_type &c, const vector_type &d, const number relax)
computes a new correction c = B*d and projects on the underlying constraint
Definition: proj_gauss_seidel_impl.h:77
virtual void step(const matrix_type &mat, vector_type &c, const vector_type &d, const number relax)
computes a new correction c = B*d and projects on the underlying constraint
Definition: proj_gauss_seidel_impl.h:108
double number
Definition: types.h:124
void backward_gs_step(Vector_type &c, const Matrix_type &A, const Vector_type &d, const size_t i, const number relaxFactor)
Definition: proj_gauss_seidel_impl.h:57
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
void forward_gs_step(Vector_type &c, const Matrix_type &A, const Vector_type &d, const size_t i, const number relaxFactor)
commmon GaussSeidel-step-calls for a single index 'i'
Definition: proj_gauss_seidel_impl.h:42
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