Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
component_gauss_seidel.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3 * Authors: Christian Wehner, Andreas Vogel
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33/*
34 * Remark: Base on Vanka idea and implementation
35 */
36
37#ifndef __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__COMPONENT_GAUSS_SEIDEL__
38#define __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__COMPONENT_GAUSS_SEIDEL__
39
41
42#include <vector>
43#include <algorithm>
44
45#ifdef UG_PARALLEL
46 #include "pcl/pcl_util.h"
49#endif
50
51namespace ug{
52
54template <typename TDomain, typename TAlgebra>
55class ComponentGaussSeidel : public IPreconditioner<TAlgebra>
56{
57 public:
60 typedef typename GF::element_type Element;
61 typedef typename GF::side_type Side;
62 static const int dim = TDomain::dim;
63
66 typedef TAlgebra algebra_type;
67 typedef typename TAlgebra::vector_type vector_type;
68 typedef typename TAlgebra::matrix_type matrix_type;
70
71 protected:
76
77 public:
79 ComponentGaussSeidel(const std::vector<std::string>& vFullRowCmp)
80 : m_bInit(false), m_relax(1), m_alpha(1.0), m_beta(1.0), m_bWeighted(false), m_vFullRowCmp(vFullRowCmp) {}
81
83 ComponentGaussSeidel(number relax, const std::vector<std::string>& vFullRowCmp)
84 : m_bInit(false), m_relax(relax), m_alpha(1.0), m_beta(1.0), m_bWeighted(false), m_vFullRowCmp(vFullRowCmp){}
85
87 ComponentGaussSeidel(number relax, const std::vector<std::string>& vFullRowCmp,
88 const std::vector<int>& vSmooth, const std::vector<number>& vDamp)
89 : m_bInit(false), m_relax(relax), m_alpha(1.0), m_beta(1.0), m_bWeighted(false), m_vFullRowCmp(vFullRowCmp),
90 m_vGroupObj(vSmooth), m_vDamp(vDamp) {};
91
94
96 virtual bool supports_parallel() const {return true;}
97
98 void set_alpha(number alpha)
99 {m_alpha = alpha;}
100
101 void set_beta(number beta)
102 {m_beta = beta;}
103
104 void set_weights(bool b)
105 {m_bWeighted = b;}
106
107 protected:
109 virtual const char* name() const {return "ComponentGaussSeidel";}
110
113
114 struct DimCache{
115 std::vector<std::vector<DoFIndex> > vvDoFIndex;
116 std::vector<DenseMatrix< VariableArray2<number> > > vBlockInverse;
117 };
118
119 void apply_blocks(const matrix_type& A, GF& c,
120 const vector_type& d, number relax,
121 const DimCache& dimCache,
122 bool bReverse);
123
124 void apply_blocks_weighted(const matrix_type& A, GF& c,
125 const vector_type& d, number relax,
126 const DimCache& dimCache,
127 bool bReverse);
128
129 void extract_blocks(const matrix_type& A, const GF& c);
130
133
135 virtual bool postprocess() {return true;}
136
137 protected:
138#ifdef UG_PARALLEL
140#endif
141
144
149
152
154 std::vector<std::string> m_vFullRowCmp;
155
157 std::vector<int> m_vGroupObj;
158 std::vector<number> m_vDamp;
159
162};
163
164template <typename TDomain, typename TAlgebra>
167{
168#ifdef UG_PARALLEL
169 if(pcl::NumProcs() > 1)
170 {
171 // copy original matrix
172 MakeConsistent(*pOp, m_A);
173 // set zero on slaves
174 std::vector<IndexLayout::Element> vIndex;
175 CollectUniqueElements(vIndex, m_A.layouts()->slave());
176 SetDirichletRow(m_A, vIndex);
177 }
178#endif
179
180 m_bInit = false;
181
182 return true;
183}
184
185template <typename TDomain, typename TAlgebra>
187apply_blocks(const matrix_type& A, GF& c,
188 const vector_type& d, number relax,
189 const DimCache& dimCache,
190 bool bReverse)
191{
192// memory for local algebra
195
196// get caches
197 const std::vector<std::vector<DoFIndex> >& vvDoFIndex = dimCache.vvDoFIndex;
198 const std::vector<DenseMatrix< VariableArray2<number> > >& vBlockInverse = dimCache.vBlockInverse;
199
200// loop blocks
201 for(size_t b = 0; b < vvDoFIndex.size(); ++b)
202 {
203 size_t block = (!bReverse) ? b : (vvDoFIndex.size()-1 - b);
204
205 // get storage
206 const std::vector<DoFIndex>& vDoFIndex = vvDoFIndex[block];
207 const DenseMatrix<VariableArray2<number> >& BlockInv = vBlockInverse[block];
208 const size_t numIndex = vDoFIndex.size();
209
210 // compute s[j] := d[j] - sum_k A(j,k)*c[k]
211 // note: the loop over k is the whole matrix row (not only selected indices)
212 typedef typename matrix_type::const_row_iterator const_row_iterator;
213 const static int blockSize = TAlgebra::blockSize;
214 s.resize(numIndex);
215 for (size_t j = 0; j<numIndex; j++){
216 typename vector_type::value_type sj = d[vDoFIndex[j][0]];
217 for(const_row_iterator it = A.begin_row(vDoFIndex[j][0]);
218 it != A.end_row(vDoFIndex[j][0]) ; ++it){
219 MatMultAdd(sj, 1.0, sj, -1.0, it.value(), c[it.index()]);
220 };
221 s.subassign(j*blockSize,sj);
222 };
223
224 // solve block
225 x.resize(numIndex);
226 MatMult(x, relax, BlockInv, s);
227
228 // add to global correction
229 for (size_t j=0;j<numIndex;j++)
230 DoFRef(c, vDoFIndex[j]) += x[j];
231 }
232}
233
234
235template <typename TDomain, typename TAlgebra>
238 const vector_type& d, number relax,
239 const DimCache& dimCache,
240 bool bReverse)
241{
242// memory for local algebra
245
246// get caches
247 const std::vector<std::vector<DoFIndex> >& vvDoFIndex = dimCache.vvDoFIndex;
248 const std::vector<DenseMatrix< VariableArray2<number> > >& vBlockInverse = dimCache.vBlockInverse;
249
250// loop blocks
251 for(size_t b = 0; b < vvDoFIndex.size(); ++b)
252 {
253 size_t block = (!bReverse) ? b : (vvDoFIndex.size()-1 - b);
254
255 // get storage
256 const std::vector<DoFIndex>& vDoFIndex = vvDoFIndex[block];
257 const DenseMatrix<VariableArray2<number> >& BlockInv = vBlockInverse[block];
258 const size_t numIndex = vDoFIndex.size();
259
260 // compute s[i] := d[i] - sum_k A(i,k)*c[k]
261 // note: the loop over k is the whole matrix row (not only selected indices)
262 typedef typename matrix_type::const_row_iterator const_row_iterator;
263 const static int blockSize = TAlgebra::blockSize;
264 s.resize(numIndex);
265 for (size_t i = 0; i<numIndex; i++)
266 {
267 typename vector_type::value_type si = d[vDoFIndex[i][0]];
268 for(const_row_iterator it = A.begin_row(vDoFIndex[i][0]);
269 it != A.end_row(vDoFIndex[i][0]) ; ++it)
270 {
271 MatMultAdd(si, 1.0, si, -1.0, it.value(), c[it.index()]);
272 }
273 s.subassign(i*blockSize, si);
274 }
275
276 // restrict
277 for (size_t i = 0; i<numIndex; ++i)
278 {
279 number wi = sqrt(DoFRef(*m_weight, vDoFIndex[i]));
280 UG_ASSERT(wi!=0.0, "Huhh: Zero weights found!" << i << "dofIndex" << vDoFIndex[i])
281 s[i] /= wi;
282 }
283
284 // solve block
285 x.resize(numIndex);
286 MatMult(x, relax, BlockInv, s);
287
288 // add to global correction
289 for (size_t i=0; i<numIndex;i++)
290 {
291 number wi = sqrt(DoFRef(*m_weight, vDoFIndex[i]));
292 DoFRef(c, vDoFIndex[i]) += x[i]/wi;
293 }
294 }
295}
296
297// extract_by_grouping(std::vector<std::vector<DoFIndex> >& vvDoFIndex,
298template<typename TGroupObj, typename GF>
299void extract_by_grouping(std::vector<std::vector<DoFIndex> >& vvDoFIndex,
300 const GF& c,
301 const std::vector<size_t>& vFullRowCmp,
302 const std::vector<size_t>& vRemainCmp)
303{
304
305 typedef typename GF::element_type Element;
306
307// memory for local algebra
308 std::vector<DoFIndex> vFullRowDoFIndex;
309 std::vector<Element*> vElem;
310
311// clear indices
312// \todo: improve this by precomputing size
313 vvDoFIndex.clear();
314
315// loop all grouping objects
316 typedef typename GF::template traits<TGroupObj>::const_iterator GroupObjIter;
317 for(GroupObjIter iter = c.template begin<TGroupObj>();
318 iter != c.template end<TGroupObj>(); ++iter)
319 {
320 // get grouping obj
321 TGroupObj* groupObj = *iter;
322
323 // get all dof indices on obj associated to cmps
324 vFullRowDoFIndex.clear();
325 for(size_t f = 0; f < vFullRowCmp.size(); ++f)
326 c.inner_dof_indices(groupObj, vFullRowCmp[f], vFullRowDoFIndex, false);
327
328 // check if equation present
329 if(vFullRowDoFIndex.empty())
330 UG_THROW("extract_by_grouping: Should not happen.");
331
332 // get all algebraic indices on element
333 if(TGroupObj::dim <= VERTEX){
334 std::vector<Edge*> vSub;
335 c.collect_associated(vSub, groupObj);
336 for(size_t i = 0; i < vSub.size(); ++i)
337 for(size_t f = 0; f < vFullRowCmp.size(); ++f)
338 c.inner_dof_indices(vSub[i], vFullRowCmp[f], vFullRowDoFIndex, false);
339 }
340 if(TGroupObj::dim <= EDGE){
341 std::vector<Face*> vSub;
342 c.collect_associated(vSub, groupObj);
343 for(size_t i = 0; i < vSub.size(); ++i)
344 for(size_t f = 0; f < vFullRowCmp.size(); ++f)
345 c.inner_dof_indices(vSub[i], vFullRowCmp[f], vFullRowDoFIndex, false);
346 }
347
348 // collect elems associated to grouping object
349 c.collect_associated(vElem, groupObj);
350
351 // get all algebraic indices on element
352 std::vector<DoFIndex> vDoFIndex;
353 for(size_t i = 0; i < vElem.size(); ++i)
354 for(size_t f = 0; f < vRemainCmp.size(); ++f)
355 c.dof_indices(vElem[i], vRemainCmp[f], vDoFIndex, false, false);
356
357 // check for dublicates
358 if(vElem.size() > 1){
359 std::sort(vDoFIndex.begin(), vDoFIndex.end());
360 vDoFIndex.erase(std::unique(vDoFIndex.begin(), vDoFIndex.end()), vDoFIndex.end());
361 }
362
363 // concat all indices
364 vDoFIndex.insert(vDoFIndex.end(), vFullRowDoFIndex.begin(), vFullRowDoFIndex.end());
365
366 // add
367 vvDoFIndex.push_back(vDoFIndex);
368 }
369}
370
371template <typename TDomain, typename TAlgebra>
373extract_blocks(const matrix_type& A, const GF& c)
374{
375
377 c.approx_space()->dof_distribution_info();
378
379 // tokenize passed functions
380 if(m_vFullRowCmp.empty())
381 UG_THROW("ComponentGaussSeidel: No components set.")
382
383 // get ids of selected functions
384 std::vector<size_t> vFullRowCmp;
385 for(size_t i = 0; i < m_vFullRowCmp.size(); ++i)
386 vFullRowCmp.push_back(ddinfo->fct_id_by_name(m_vFullRowCmp[i].c_str()));
387
388 // compute remaining functions
389 std::vector<size_t> vRemainCmp;
390 for(size_t f = 0; f < ddinfo->num_fct(); ++f)
391 if(std::find(vFullRowCmp.begin(), vFullRowCmp.end(), f) == vFullRowCmp.end())
392 vRemainCmp.push_back(f);
393
395 if(m_vGroupObj.empty()){
396 for(int d = VERTEX; d <= VOLUME; ++d){
397 bool bCarryDoFs = false;
398 for(size_t f = 0; f < vFullRowCmp.size(); ++f)
399 if(ddinfo->max_fct_dofs(vFullRowCmp[f], d) > 0)
400 bCarryDoFs = true;
401 if(bCarryDoFs)
402 m_vGroupObj.push_back(d);
403 }
404 }
405
406// extract for each dim-grouping objects
407 for(int d = VERTEX; d <= VOLUME; ++d)
408 {
409 // only extract if needed
410 if(std::find(m_vGroupObj.begin(), m_vGroupObj.end(), d) == m_vGroupObj.end())
411 continue;
412
413 // only extract if carrying dofs
414 bool bCarryDoFs = false;
415 for(size_t f = 0; f < vFullRowCmp.size(); ++f)
416 if(ddinfo->max_fct_dofs(vFullRowCmp[f], d) > 0)
417 bCarryDoFs = true;
418 if(!bCarryDoFs)
419 continue;
420
421 // extract
422 std::vector<std::vector<DoFIndex> >& vvDoFIndex = m_vDimCache[d].vvDoFIndex;
423 switch(d){
424 case VERTEX: extract_by_grouping<Vertex, GF>(vvDoFIndex, c, vFullRowCmp, vRemainCmp); break;
425 case EDGE: extract_by_grouping<Edge, GF>(vvDoFIndex, c, vFullRowCmp, vRemainCmp); break;
426 case FACE: extract_by_grouping<Face, GF>(vvDoFIndex, c, vFullRowCmp, vRemainCmp); break;
427 case VOLUME: extract_by_grouping<Volume, GF>(vvDoFIndex, c, vFullRowCmp, vRemainCmp); break;
428 default: UG_THROW("wrong dim");
429 }
430 }
431
432 // compute weights v(i) (= number of groups j that i belongs to)
433 if (m_bWeighted)
434 {
435 m_weight = c.clone_without_values();
436 m_weight->set(0.0);
437 for(int d = VERTEX; d <= VOLUME; ++d)
438 {
439 std::vector<std::vector<DoFIndex> >& vvDoFIndex = m_vDimCache[d].vvDoFIndex;
440
441 for(size_t j = 0; j < vvDoFIndex.size(); ++j)
442 {
443
444 std::vector<DoFIndex>& vDoFIndex = vvDoFIndex[j];
445 const size_t numIndex = vDoFIndex.size();
446
447 // count number of connections for this row
448 DoFRef(*m_weight, vDoFIndex[numIndex-1]) = 1.0;
449 for (size_t i =0; i<numIndex-1; ++i)
450 DoFRef(*m_weight, vDoFIndex[i]) += 1.0;
451 }
452
453 }
454 }
455
456 // extract local matrices
457 for(int d = VERTEX; d <= VOLUME; ++d)
458 {
459 std::vector<std::vector<DoFIndex> >& vvDoFIndex = m_vDimCache[d].vvDoFIndex;
460 std::vector<DenseMatrix<VariableArray2<number> > >& vBlockInverse = m_vDimCache[d].vBlockInverse;
461
462 vBlockInverse.resize(vvDoFIndex.size());
463 for(size_t b = 0; b < vvDoFIndex.size(); ++b)
464 {
465 // get storage
466 std::vector<DoFIndex>& vDoFIndex = vvDoFIndex[b];
467 DenseMatrix<VariableArray2<number> >& BlockInv = vBlockInverse[b];
468
469 // get number of indices on patch
470 const size_t numIndex = vDoFIndex.size();
471 // std::cerr << "locSize" << numIndex << std::endl;
472
473 // fill local block matrix
474 BlockInv.resize(numIndex, numIndex);
475 BlockInv = 0.0;
476
477 // copy matrix rows (only including cols of selected indices)
478 for (size_t i = 0; i < numIndex; i++)
479 {
480 // Diag (A)
481 BlockInv(i,i) = m_alpha*DoFRef(A, vDoFIndex[i], vDoFIndex[i]);
482 BlockInv(numIndex-1,i) = DoFRef(A, vDoFIndex[numIndex-1], vDoFIndex[i]);
483 BlockInv(i,numIndex-1) = DoFRef(A, vDoFIndex[i], vDoFIndex[numIndex-1]);
484
485 //for (size_t k = 0; k < numIndex; k++)
486 // BlockInv(i,k) = DoFRef(A, vDoFIndex[i], vDoFIndex[k]);
487 }
488
489
490 // schur complement correction
491
492
493 number schur = 0.0;
494 if (m_weight.invalid())
495 {
496 for (size_t i = 0; i < numIndex-1; i++)
497 {
498 schur += BlockInv(numIndex-1,i)/BlockInv(i,i)*BlockInv(i,numIndex-1);
499 }
500 // std::cerr << "unweighted:" << schur << " " << BlockInv(numIndex-1, numIndex-1)<< m_beta<< "=>";
501
502 BlockInv(numIndex-1, numIndex-1) = schur + (BlockInv(numIndex-1, numIndex-1)-schur)/m_beta;
503
504 //std::cerr << BlockInv(numIndex-1, numIndex-1) << std::endl;
505
506 } else {
507 for (size_t i = 0; i < numIndex-1; i++)
508 {
509 number wi2 = DoFRef(*m_weight, vDoFIndex[i]);
510 schur += wi2*BlockInv(numIndex-1,i)/BlockInv(i,i)*BlockInv(i,numIndex-1);
511 }
512 //std::cerr << "weighted:"<< schur << " " << BlockInv(numIndex-1, numIndex-1)<< m_beta<< "=>";
513
514 // BlockInv(numIndex-1, numIndex-1) = - schur - (BlockInv(numIndex-1, numIndex-1)-schur)/m_beta; // worked with alpha=1, beta=-2.0
515 BlockInv(numIndex-1, numIndex-1) = schur + (BlockInv(numIndex-1, numIndex-1) - schur)/m_beta;
516
517 //std::cerr << BlockInv(numIndex-1, numIndex-1) << std::endl;
518 }
519
520
521 //BlockInv(numIndex-1, numIndex-1) = - 2.0*schur;
522
523
524
525
526 // get inverse
527 if(!Invert(BlockInv)){
528 std::stringstream ss;
529 ss << d <<"dim - Elem-Mat "<<b<<" singular. size: "<<numIndex<<"\n";
530 for (size_t j = 0; j < numIndex; j++)
531 ss << vDoFIndex[j][0] << ", ";
532 ss << "\n";
533
534 BlockInv = 0.0;
535 for (size_t j = 0; j < numIndex; j++)
536 for (size_t k = 0; k < numIndex; k++)
537 BlockInv(j,k) = DoFRef(A, vDoFIndex[j], vDoFIndex[k]);
538
539 ss << BlockInv;
540
541 UG_THROW(ss.str());
542 }
543 }
544 }
545
546 m_bInit = true;
547}
548
549template <typename TDomain, typename TAlgebra>
552{
553// check that grid funtion passed
555 = dynamic_cast<GridFunction<TDomain, TAlgebra>*>(&c);
556 if(pC == NULL)
557 UG_THROW("ComponentGaussSeidel: expects correction to be a GridFunction.");
558
559 const vector_type* pD = &d;
560 const matrix_type* pMat = pOp.get();
561#ifdef UG_PARALLEL
563 if(pcl::NumProcs() > 1){
564 // make defect unique
565 spDtmp = d.clone();
566 spDtmp->change_storage_type(PST_UNIQUE);
567 pD = spDtmp.get();
568 pMat = &m_A;
569 }
570#endif
571
572// check if initialized
573 if(!m_bInit)
574 extract_blocks(*pMat, *pC);
575
576// set correction to zero
577 pC->set(0.0);
578
579// loop Grouping objects
580 for(size_t i = 0; i < m_vGroupObj.size(); ++i)
581 {
582 // get its dimension
583 const int d = m_vGroupObj[i];
584
585 // get relax param
586 number damp = m_relax;
587 if(d < (int)m_vDamp.size()) damp *= m_vDamp[d];
588
589 // apply
590 if (m_weight.invalid()) {
591 apply_blocks(*pMat, *pC, *pD, damp, m_vDimCache[d], false);
592 apply_blocks(*pMat, *pC, *pD, damp, m_vDimCache[d], true);
593 } else {
594 apply_blocks_weighted(*pMat, *pC, *pD, damp, m_vDimCache[d], false);
595 apply_blocks_weighted(*pMat, *pC, *pD, damp, m_vDimCache[d], true);
596 }
597
598 }
599
600#ifdef UG_PARALLEL
601 pC->set_storage_type(PST_UNIQUE);
602#endif
603
604 return true;
605}
606
607
608template <typename TDomain, typename TAlgebra>
611{
613 newInst(new ComponentGaussSeidel<TDomain, TAlgebra>(m_relax, m_vFullRowCmp, m_vGroupObj, m_vDamp));
614 newInst->set_debug(debug_writer());
615 newInst->set_damp(this->damping());
616 newInst->set_alpha(this->m_alpha);
617 newInst->set_beta(this->m_beta);
618 newInst->set_weights(this->m_bWeighted);
619 return newInst;
620}
621
622
623} // end namespace ug
624
625#endif /* __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__COMPONENT_GAUSS_SEIDEL__ */
parameterString s
Definition smart_pointer.h:296
Definition smart_pointer.h:108
T * get()
returns encapsulated pointer
Definition smart_pointer.h:197
ComponentGaussSeidel Preconditioner.
Definition component_gauss_seidel.h:56
virtual const char * name() const
Name of preconditioner.
Definition component_gauss_seidel.h:109
ComponentGaussSeidel(const std::vector< std::string > &vFullRowCmp)
default constructor
Definition component_gauss_seidel.h:79
virtual bool postprocess()
Postprocess routine.
Definition component_gauss_seidel.h:135
ComponentGaussSeidel(number relax, const std::vector< std::string > &vFullRowCmp, const std::vector< int > &vSmooth, const std::vector< number > &vDamp)
constructor setting relaxation and type
Definition component_gauss_seidel.h:87
static const int dim
Definition component_gauss_seidel.h:62
void set_beta(number beta)
Definition component_gauss_seidel.h:101
IPreconditioner< TAlgebra > base_type
Definition component_gauss_seidel.h:72
number m_relax
relaxing parameter
Definition component_gauss_seidel.h:146
std::vector< std::string > m_vFullRowCmp
components, whose matrix row must be fulfilled completely on gs-block
Definition component_gauss_seidel.h:154
TAlgebra::matrix_type matrix_type
Definition component_gauss_seidel.h:68
virtual bool preprocess(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp)
Preprocess routine.
Definition component_gauss_seidel.h:166
SmartPtr< GF > m_weight
Definition component_gauss_seidel.h:151
TAlgebra::vector_type vector_type
Definition component_gauss_seidel.h:67
GridFunction< TDomain, TAlgebra > GF
World dimension.
Definition component_gauss_seidel.h:59
number m_beta
Definition component_gauss_seidel.h:148
bool m_bWeighted
Definition component_gauss_seidel.h:150
void extract_blocks(const matrix_type &A, const GF &c)
Definition component_gauss_seidel.h:373
std::vector< int > m_vGroupObj
smooth order
Definition component_gauss_seidel.h:157
std::vector< number > m_vDamp
Definition component_gauss_seidel.h:158
bool m_bInit
init flag
Definition component_gauss_seidel.h:143
GF::element_type Element
Definition component_gauss_seidel.h:60
matrix_type m_A
Definition component_gauss_seidel.h:139
GF::side_type Side
Definition component_gauss_seidel.h:61
void set_weights(bool b)
Definition component_gauss_seidel.h:104
number m_alpha
Definition component_gauss_seidel.h:147
void apply_blocks_weighted(const matrix_type &A, GF &c, const vector_type &d, number relax, const DimCache &dimCache, bool bReverse)
Definition component_gauss_seidel.h:237
virtual bool step(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp, vector_type &c, const vector_type &d)
step method
Definition component_gauss_seidel.h:551
virtual SmartPtr< ILinearIterator< vector_type > > clone()
Clone.
Definition component_gauss_seidel.h:610
DimCache m_vDimCache[VOLUME+1]
caching storage
Definition component_gauss_seidel.h:161
TAlgebra algebra_type
Definition component_gauss_seidel.h:66
ComponentGaussSeidel(number relax, const std::vector< std::string > &vFullRowCmp)
constructor setting relaxation and type
Definition component_gauss_seidel.h:83
void apply_blocks(const matrix_type &A, GF &c, const vector_type &d, number relax, const DimCache &dimCache, bool bReverse)
Definition component_gauss_seidel.h:187
void set_alpha(number alpha)
Definition component_gauss_seidel.h:98
virtual bool supports_parallel() const
returns if parallel solving is supported
Definition component_gauss_seidel.h:96
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
SmartPtr< IDebugWriter< algebra_type > > debug_writer()
returns the debug writer
Definition debug_writer.h:391
Definition densematrix.h:57
Definition densevector.h:101
represents numerical solutions on a grid using an algebraic vector
Definition grid_function.h:121
SmartPtr< ApproximationSpace< TDomain > > approx_space()
returns approx space
Definition grid_function.h:348
SmartPtr< this_type > clone_without_values() const
clone excluding values
Definition grid_function.h:219
dim_traits< dim >::grid_base_object element_type
Definition grid_function.h:159
dim_traits< dim-1 >::grid_base_object side_type
Definition grid_function.h:163
describes a linear iterator that is based on a matrix operator
Definition preconditioner.h:103
Definition matrix_operator.h:49
@ 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
#define UG_ASSERT(expr, msg)
Definition assert.h:70
#define UG_THROW(msg)
Definition error.h:57
double number
Definition types.h:124
the ug namespace
bool MakeConsistent(const ParallelMatrix< matrix_type > &_mat, ParallelMatrix< matrix_type > &newMat)
Definition parallel_matrix_overlap_impl.h:438
bool Invert(T &m)
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 extract_by_grouping(std::vector< std::vector< DoFIndex > > &vvDoFIndex, const GF &c, const std::vector< size_t > &vFullRowCmp, const std::vector< size_t > &vRemainCmp)
Definition component_gauss_seidel.h:299
@ VOLUME
Definition grid_base_objects.h:63
@ VERTEX
Definition grid_base_objects.h:60
@ EDGE
Definition grid_base_objects.h:61
@ FACE
Definition grid_base_objects.h:62
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
number & DoFRef(TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
Definition multi_index.h:276
Definition component_gauss_seidel.h:114
std::vector< DenseMatrix< VariableArray2< number > > > vBlockInverse
Definition component_gauss_seidel.h:116
std::vector< std::vector< DoFIndex > > vvDoFIndex
Definition component_gauss_seidel.h:115