33 #ifndef __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__VANKA__
34 #define __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__VANKA__
49 template<
typename Matrix_type,
typename Vector_type>
50 bool Vanka_step(
const Matrix_type &A, Vector_type &x,
const Vector_type &b,
number relax)
58 for(
size_t i=0; i < x.size(); i++)
63 for(
size_t i=0; i < x.size(); i++)
79 for(
typename Matrix_type::const_row_iterator it = A.begin_row(i); it != A.end_row(i) ; ++it){
80 blockind[blocksize] = it.index();
84 mat.resize(blocksize,blocksize);
86 localx.resize(blocksize);
87 for (
size_t j=0;j<blocksize;j++){
89 for (
size_t k=0;k<blocksize;k++){
90 mat.
subassign(j,k,A(blockind[j],blockind[k]));
94 for(
typename Matrix_type::const_row_iterator it = A.begin_row(blockind[j]); it != A.end_row(blockind[j]) ; ++it){
95 MatMultAdd(sj, 1.0, sj, -1.0, it.value(), x[it.index()]);
101 for (
size_t j=0;j<blocksize;j++){
102 x[blockind[j]] += relax*localx[j];
123 template<
typename Matrix_type,
typename Vector_type>
135 for(
size_t i=0; i < x.size(); i++)
140 for(
size_t i=0; i < x.size(); i++)
156 for(
typename Matrix_type::const_row_iterator it = A.begin_row(i); it != A.end_row(i) ; ++it){
157 if (it.index()==i)
continue;
158 blockind[blocksize] = it.index();
159 s[blocksize] = b[blockind[blocksize]];
160 for(
typename Matrix_type::const_row_iterator rowit = A.begin_row(blockind[blocksize]); rowit != A.end_row(blockind[blocksize]) ; ++rowit){
161 if ((rowit.index()==blockind[blocksize])||(rowit.index()==i))
continue;
163 MatMultAdd(
s[blocksize], 1.0,
s[blocksize], -1.0, rowit.value(), x[rowit.index()]);
169 block_type a_ii = A(i,i);
172 for (
size_t j=0;j<blocksize;j++){
173 block_type a_q = A(i,blockind[j]);
174 block_type a_jj = A(blockind[j],blockind[j]);
179 a_ii-=a_q*A(blockind[j],i);
184 for (
size_t j=0;j<blocksize;j++){
186 MatMultAdd(
s[j], 1.0,
s[j], -1.0, A(blockind[j],i), x[i]);
188 InverseMatMult(x[blockind[j]], relax, A(blockind[j],blockind[j]),
s[j]);
196 template <
typename TAlgebra>
229 newInst->set_damp(this->
damping());
248 virtual const char*
name()
const {
return "Vanka";}
259 std::vector<IndexLayout::Element> vIndex;
275 dhelp.resize(d.size()); dhelp = d;
306 template <
typename TAlgebra>
339 newInst->set_damp(this->
damping());
358 virtual const char*
name()
const {
return "DiagVanka";}
369 std::vector<IndexLayout::Element> vIndex;
385 dhelp.resize(d.size()); dhelp = d;
Definition: smart_pointer.h:108
virtual void set_debug(SmartPtr< IDebugWriter< algebra_type > > spDebugWriter)
set debug writer
Definition: debug_writer.h:384
void write_debug(const matrix_type &mat, const char *filename)
write debug output for a matrix (if debug writer set)
Definition: debug_writer.h:399
TAlgebra::vector_type vector_type
type of vector
Definition: debug_writer.h:360
SmartPtr< IDebugWriter< algebra_type > > debug_writer()
returns the debug writer
Definition: debug_writer.h:391
TAlgebra::matrix_type matrix_type
type of matrix
Definition: debug_writer.h:363
Definition: densematrix.h:57
void subassign(size_t r, size_t c, const TStorage2 &t)
Definition: densematrix.h:215
Definition: densevector.h:101
Diagvanka Preconditioner, description see above diagvanka_step function.
Definition: vanka.h:308
number m_relax
Definition: vanka.h:351
virtual ~DiagVanka()
Destructor.
Definition: vanka.h:346
IPreconditioner< TAlgebra > base_type
Base type.
Definition: vanka.h:323
TAlgebra::matrix_type matrix_type
Matrix type.
Definition: vanka.h:317
virtual bool preprocess(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp)
Preprocess routine.
Definition: vanka.h:361
DiagVanka()
default constructor
Definition: vanka.h:332
virtual bool step(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp, vector_type &c, const vector_type &d)
computes a new correction c = B*d
Definition: vanka.h:377
TAlgebra algebra_type
Algebra type.
Definition: vanka.h:311
IPreconditioner< TAlgebra >::matrix_operator_type matrix_operator_type
Matrix Operator type.
Definition: vanka.h:320
virtual SmartPtr< ILinearIterator< vector_type > > clone()
Clone.
Definition: vanka.h:335
void set_relax(number omega)
set relaxation parameter
Definition: vanka.h:351
matrix_type m_A
Definition: vanka.h:411
virtual bool supports_parallel() const
returns if parallel solving is supported
Definition: vanka.h:343
TAlgebra::vector_type vector_type
Vector type.
Definition: vanka.h:314
virtual const char * name() const
Name of preconditioner.
Definition: vanka.h:358
virtual bool postprocess()
Postprocess routine.
Definition: vanka.h:407
SmartPtr< IDamping< TAlgebra::vector_type, TAlgebra::vector_type > > damping()
returns the scaling
Definition: linear_iterator.h:173
describes a linear iterator that is based on a matrix operator
Definition: preconditioner.h:103
Definition: matrix_operator.h:49
Vanka Preconditioner.
Definition: vanka.h:198
Vanka()
default constructor
Definition: vanka.h:222
TAlgebra::matrix_type matrix_type
Matrix type.
Definition: vanka.h:207
TAlgebra algebra_type
Algebra type.
Definition: vanka.h:201
matrix_type m_A
Definition: vanka.h:300
virtual ~Vanka()
Destructor.
Definition: vanka.h:234
virtual bool supports_parallel() const
returns if parallel solving is supported
Definition: vanka.h:237
virtual const char * name() const
Name of preconditioner.
Definition: vanka.h:248
TAlgebra::vector_type vector_type
Vector type.
Definition: vanka.h:204
IPreconditioner< TAlgebra >::matrix_operator_type matrix_operator_type
Matrix Operator type.
Definition: vanka.h:210
void set_relax(number omega)
set relaxation parameter
Definition: vanka.h:241
virtual SmartPtr< ILinearIterator< vector_type > > clone()
Clone.
Definition: vanka.h:225
number m_relax
Definition: vanka.h:241
IPreconditioner< TAlgebra > base_type
Base type.
Definition: vanka.h:213
virtual bool step(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp, vector_type &c, const vector_type &d)
computes a new correction c = B*d
Definition: vanka.h:267
virtual bool postprocess()
Postprocess routine.
Definition: vanka.h:296
virtual bool preprocess(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp)
Preprocess routine.
Definition: vanka.h:251
@ PST_UNIQUE
Definition: parallel_storage_type.h:70
void SetDirichletRow(TSparseMatrix &A, size_t i, size_t alpha)
Definition: sparsematrix_util.h:796
int NumProcs()
returns the number of processes
Definition: pcl_base.cpp:91
void CollectUniqueElements(std::vector< typename TLayout::Element > &elemsOut, const TLayout &layout)
writes all elements in the interfaces into the resulting vector. avoids doubles.
Definition: pcl_layout_util.h:142
#define UG_THROW(msg)
Definition: error.h:57
double number
Definition: types.h:124
static const size_t MAXBLOCKSIZE
Definition: vanka.h:47
bool MakeConsistent(const ParallelMatrix< matrix_type > &_mat, ParallelMatrix< matrix_type > &newMat)
Definition: parallel_matrix_overlap_impl.h:438
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 Diag_Vanka_step(const Matrix_type &A, Vector_type &x, const Vector_type &b, number relax)
Definition: vanka.h:124
bool Vanka_step(const Matrix_type &A, Vector_type &x, const Vector_type &b, number relax)
Definition: vanka.h:50
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