Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
sequential_subspace_correction.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2018: G-CSC, Goethe University Frankfurt
3 * Authors: Arne Naegel
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: Based on implementation of element/component Seidel
35 */
36
37#ifndef __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__SEQUENTIAL_SSC__
38#define __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__SEQUENTIAL_SSC__
39
41
42#include <vector>
43#include <algorithm>
44
45#ifdef UG_PARALLEL
46 #include "pcl/pcl_util.h"
51#endif
52
53namespace ug{
54
55
57template <typename TDomain, typename TAlgebra, typename TObject>
59{
60public:
62 typedef TAlgebra algebra_type;
63
65 typedef typename TAlgebra::vector_type vector_type;
66
68 typedef typename TAlgebra::matrix_type matrix_type;
69
74
77
79 virtual bool preprocess(const vector_type &c)
80 { return true; }
81
82
84 virtual void init(TObject*, const vector_type &) = 0;
85
87 virtual void extract_matrix(const matrix_type &A) = 0;
88
90 virtual void extract_rhs(const vector_type &d) = 0;
91
93 virtual void extract_rhs(const vector_type &d, const matrix_type &A, const vector_type &c) = 0;
94
96 virtual void update_solution(vector_type &u, double omega=1.0) = 0;
97
98 virtual size_t size() { return 0; }
99};
100
101
102
104template <typename TDomain, typename TAlgebra, typename TObject>
105class LocalIndexSubspace : public ILocalSubspace<TDomain, TAlgebra, TObject>
106{
107public:
109 typedef TAlgebra algebra_type;
110
112 typedef typename TAlgebra::vector_type vector_type;
113
115 typedef typename TAlgebra::matrix_type matrix_type;
116
118 typedef std::vector<size_t> index_vector;
121
124
126 /*virtual bool preprocess(const vector_type &c, obj_iterator_type &objBegin, obj_iterator_type &objEnd)
127 {
128 objBegin = c.template begin<TGroupObj>();
129 objEnd = c.template end<TGroupObj>();
130 return true;
131 }
132*/
134 virtual void init(TObject*, const vector_type &) = 0;
135
137 virtual void extract_matrix(const matrix_type &A)
138 {
139 typedef typename TAlgebra::matrix_type::const_row_iterator const_row_iterator;
140 const static int blockSize = TAlgebra::blockSize;
141 const size_t numIndex = this->size();
142
143 // fill local block matrix
144 bool bFound;
145 m_Aloc.resize(numIndex, numIndex);
146 m_Aloc = 0.0;
147 for (size_t j = 0; j<numIndex; j++)
148 {
149 for (size_t k=0;k<numIndex;k++)
150 {
151 const_row_iterator it = A.get_connection(m_vInd[j],m_vInd[k], bFound);
152 if(bFound){
153 m_Aloc.subassign(j*blockSize,k*blockSize,it.value());
154 }
155 }
156 }
157 }
158
159
161 virtual void extract_rhs(const vector_type &d)
162 {
163 const size_t numIndex = this->size();
164 const static int blockSize = TAlgebra::blockSize;
165
166 // compute s[j] := d[j] - sum_k A(j,k)*c[k]
167 // note: the loop over k is the whole matrix row (not only selected indices)
168 m_dloc.resize(numIndex);
169 for (size_t j = 0; j<numIndex; j++)
170 {
171 typename TAlgebra::vector_type::value_type sumj = d[m_vInd[j]];
172 m_dloc.subassign(j*blockSize,sumj);
173 }
174 }
175
177 virtual void extract_rhs(const vector_type &d, const matrix_type &A, const vector_type &c)
178 {
179 typedef typename TAlgebra::matrix_type::const_row_iterator const_row_iterator;
180 const static int blockSize = TAlgebra::blockSize;
181 const size_t numIndex = this->size();
182
183 // compute s[j] := d[j] - sum_k A(j,k)*c[k]
184 // note: the loop over k is the whole matrix row (not only selected indices)
185 m_dloc.resize(numIndex);
186 for (size_t j = 0; j<numIndex; j++)
187 {
188 typename TAlgebra::vector_type::value_type sumj = d[m_vInd[j]];
189 for(const_row_iterator it = A.begin_row(m_vInd[j]); it != A.end_row(m_vInd[j]) ; ++it)
190 {
191 MatMultAdd(sumj, 1.0, sumj, -1.0, it.value(), c[it.index()]);
192 };
193
194 m_dloc.subassign(j*blockSize,sumj);
195
196 };
197 };
198
200 virtual void update_solution(vector_type &u, double omega=1.0)
201 {
202 const size_t numIndex = this->size();
203
204 // solve block
205 m_uloc.resize(numIndex);
207 for (size_t j=0;j<numIndex;j++)
208 {
209 u[m_vInd[j]] += omega*m_uloc[j];
210 }
211
212 }
213
214
215 virtual size_t size() { return m_vInd.size(); }
216
217protected:
218
221
226
227};
228
229
231template <typename TDomain, typename TAlgebra, typename TObject>
232class LocalDoFSubspace : public ILocalSubspace<TDomain, TAlgebra, TObject>
233{
234public:
236 typedef TAlgebra algebra_type;
237
239 typedef typename TAlgebra::vector_type vector_type;
240
242 typedef typename TAlgebra::matrix_type matrix_type;
243
245 typedef std::vector<DoFIndex> index_vector;
248
251
252 bool check(void *obj) const
253 { return (dynamic_cast<TObject*> (obj) != NULL); }
254
256 virtual void init(TObject*, const vector_type &) = 0;
257
259 virtual void extract_matrix(const matrix_type &A)
260 {
261 //typedef typename TAlgebra::matrix_type::const_row_iterator const_row_iterator;
262 // const static int blockSize = TAlgebra::blockSize;
263 const size_t numIndex = this->size();
264
265 // fill local block matrix
266 m_Aloc.resize(numIndex, numIndex);
267 m_Aloc = 0.0;
268 for (size_t j = 0; j<numIndex; ++j)
269 {
270 for (size_t k=0; k<numIndex; ++k)
271 {
272 m_Aloc.entry(j,k) = DoFRef(A,m_vInd[j],m_vInd[k]);
273 }
274 }
275
276 }
277
278
280 virtual void extract_rhs(const vector_type &d)
281 {
282 const size_t numIndex = this->size();
283 m_dloc.resize(numIndex);
284
285 for (size_t j = 0; j<numIndex; j++)
286 {
287 m_dloc[j] = DoFRef(d, m_vInd[j]);
288 }
289 }
290
292 virtual void extract_rhs(const vector_type &d, const matrix_type &A, const vector_type &c)
293 {
294 typedef typename TAlgebra::matrix_type::const_row_iterator const_row_iterator;
295 const static int blockSize = TAlgebra::blockSize;
296 const size_t numIndex = this->size();
297
298 // compute s[j] := d[j] - sum_k A(j,k)*c[k]
299 // note: the loop over k is the whole matrix row (not only selected indices)
300 // (code taken from component Gauss-Seidel)
301 m_dloc.resize(numIndex);
302 for (size_t j = 0; j<numIndex; j++)
303 {
304 typename TAlgebra::vector_type::value_type sumdj = d[m_vInd[j][0]];
305 for(const_row_iterator it = A.begin_row(m_vInd[j][0]); it != A.end_row(m_vInd[j][0]) ; ++it)
306 { MatMultAdd(sumdj, 1.0, sumdj, -1.0, it.value(), c[it.index()]); }
307
308 // TODO: CHECK!!!
309 m_dloc.subassign(j*blockSize, sumdj);
310 };
311 };
312
313
314
316 virtual void update_solution(vector_type &u, double omega=1.0)
317 {
318 const size_t numIndex = this->size();
319
320 // solve block
321 m_uloc.resize(numIndex);
323 for (size_t j=0;j<numIndex;j++)
324 {
325 DoFRef(u, m_vInd[j]) += omega*m_uloc[j];
326 }
327
328 }
329
330
331 virtual size_t size() { return m_vInd.size(); }
332
333protected:
335 std::vector<DoFIndex> m_vInd;
336
337 // Memory for local algebra
341
342};
343
344
345
347template <typename TDomain, typename TAlgebra>
348class VertexBasedSubspace : public LocalIndexSubspace<TDomain, TAlgebra,Vertex>
349{
350protected:
352
353public:
355 typedef TAlgebra algebra_type;
356
358 typedef typename TAlgebra::vector_type vector_type;
359
361 typedef typename TAlgebra::matrix_type matrix_type;
362
365
368
371
372 bool check(void *obj) const
373 { return true; /*return (dynamic_cast<TObject*> (obj) != NULL);*/ }
374
376 void init(void *obj, const vector_type &cvec)
377 {
378 UG_ASSERT(check(obj), "HUHH: Expecting vertex!");
379 Vertex *groupObj = reinterpret_cast<Vertex*> (obj);
380
381 // We will modify index list of base class.
382 std::vector<size_t> &vInd = base_type::m_vInd;
383 vInd.clear();
384
385 // Union of elements (associated with grouping object).
386 typedef GridFunction<TDomain, TAlgebra> TGridFunction;
387 typedef typename GridFunction<TDomain, TAlgebra>::element_type TElement;
388 const TGridFunction *c = dynamic_cast<const TGridFunction *> (&cvec); // Need a grid function here!
389 std::vector<TElement*> vElem;
390 c->collect_associated(vElem, groupObj);
391
392 // Union of algebraic indices.
393 for(size_t i = 0; i < vElem.size(); ++i)
394 {
395 c->algebra_indices(vElem[i], vInd, false);
396 }
397
398 // Remove dublicates.
399 if(vElem.size() > 1)
400 {
401 std::sort(vInd.begin(), vInd.end());
402 vInd.erase(std::unique(vInd.begin(), vInd.end()), vInd.end());
403 }
404 }
405
406};
407
408
410template <typename TDomain, typename TAlgebra>
411class VertexCenteredVankaSubspace : public LocalDoFSubspace<TDomain, TAlgebra, Vertex>
412{
413public:
415 typedef TAlgebra algebra_type;
416
418 typedef typename TAlgebra::vector_type vector_type;
419
421 typedef typename TAlgebra::matrix_type matrix_type;
422
424
427
429 VertexCenteredVankaSubspace<TDomain, TAlgebra>(const std::vector<std::string> &vVtxCmp, const std::vector<std::string> &vElemCmp)
430 : m_vVtxCmp(vVtxCmp), m_vElemCmp(vElemCmp) {}
431
434
436protected:
438 {
440 c->approx_space()->dof_distribution_info();
441 UG_COND_THROW(ddinfo.invalid(), "Requiring valid ddinfo!");
442
443 // Vertex functions.
444 m_vVtxFct.reserve(m_vVtxCmp.size());
445 m_vVtxFct.clear();
446 for(size_t i = 0; i < m_vVtxCmp.size(); ++i)
447 m_vVtxFct.push_back(ddinfo->fct_id_by_name(m_vVtxCmp[i].c_str()));
448
449 // Element functions.
450 m_vElemFct.reserve(m_vElemCmp.size());
451 m_vElemFct.clear();
452 for(size_t i = 0; i < m_vElemCmp.size(); ++i)
453 m_vElemFct.push_back(ddinfo->fct_id_by_name(m_vElemCmp[i].c_str()));
454 }
455public:
457 bool preprocess(const vector_type &cvec)
458 {
459 typedef GridFunction<TDomain, TAlgebra> TGridFunction;
460 const TGridFunction *c = dynamic_cast<const TGridFunction *> (&cvec); // Need a grid function here!
461
462 UG_COND_THROW(c==NULL, "Requiring a grid function here!");
463 extract_ids(c);
464
465 return true;
466 }
467
469 void init(Vertex *groupObj, const vector_type &cvec)
470 {
471 // We will modify index list of base class.
473 vInd.clear();
474
475 // Union of elements (associated with grouping object).
476 typedef GridFunction<TDomain, TAlgebra> TGridFunction;
477 typedef typename GridFunction<TDomain, TAlgebra>::element_type TElement;
478 const TGridFunction *c = dynamic_cast<const TGridFunction *> (&cvec); // Need a grid function here!
479 std::vector<TElement*> vElem;
480 c->collect_associated(vElem, groupObj);
481
482 // Collect associated indices.
483 for(size_t i = 0; i < vElem.size(); ++i)
484 {
485 for(size_t f = 0; f < m_vElemFct.size(); ++f)
486 c->dof_indices(vElem[i], m_vElemFct[f], vInd, false, false);
487 }
488
489 // Collect vertex indices.
490 for(size_t f = 0; f < m_vVtxFct.size(); ++f)
491 c->dof_indices(groupObj, m_vVtxFct[f], vInd, false, false);
492
493 // Remove dublicates.
494 if(vElem.size() > 1)
495 {
496 std::sort(vInd.begin(), vInd.end());
497 vInd.erase(std::unique(vInd.begin(), vInd.end()), vInd.end());
498 }
499
500 /*
501 UG_DLOG("VertexBasedVankaSubspace [for " << groupObj << "]:")
502 for (typename base_type::index_vector::iterator it = vInd.begin() ; it != vInd.end(); ++it)
503 { UG_DLOG(*it <<","); }
504
505 UG_DLOG(std::endl);
506 */
507 }
508
509
510
511protected:
512 std::vector<std::string> m_vVtxCmp; // primary (vertex) components
513 std::vector<std::string> m_vElemCmp; // secondary (element) components
514
515 std::vector<size_t> m_vVtxFct; // mapping to fct IDs
516 std::vector<size_t> m_vElemFct; // mapping to fct IDs
517
518};
519
520
521/*
523template <typename TDomain, typename TAlgebra>
524class LocalSubspaceCorrection :
525 public ILinearIterator<typename TAlgebra::vector_type>,
526 public DebugWritingObject<TAlgebra>
527{
528public:
530 typedef TAlgebra algebra_type;
531
533 typedef typename TAlgebra::vector_type vector_type;
534
536 typedef typename TAlgebra::matrix_type matrix_type;
537
538
539 typedef MatrixOperator<matrix_type, vector_type> TLinearOperator;
540
541 virtual ~LocalSubspaceCorrection<TDomain, TAlgebra>() {}
542
544 virtual const char* name() const = 0;
545
547 bool supports_parallel() {return false};
548
550 bool init(SmartPtr<TLinearOperator> J, const vector_type& u) {
551
552 };
553
555 virtual bool init(SmartPtr<TLinearOperator> L) = 0;
556
558 virtual bool apply(vector_type& c, const vector_type& d)
559 {
560
561 }
562
564 virtual bool apply_update_defect(vector_type& c, vector_type& d)
565 {
566
567 }
568
569};
570
571
572*/
573
574
576template<typename TGroupObj, typename TDomain, typename TAlgebra>
577void ParallelSubspaceCorrectionLoop(const typename TAlgebra::matrix_type& A,
579 const typename TAlgebra::vector_type& d,
580 number omega_relax,
582 typename GridFunction<TDomain, TAlgebra>::template traits<TGroupObj>::const_iterator objIterBegin,
583 typename GridFunction<TDomain, TAlgebra>::template traits<TGroupObj>::const_iterator objIterEnd
584)
585{
586 // Loop over all grouping objects.
587 typedef typename GridFunction<TDomain, TAlgebra>::template traits<TGroupObj>::const_iterator GroupObjIter;
588 for(GroupObjIter iter = objIterBegin; iter != objIterEnd; ++iter)
589 {
590 // Apply subspace correction (w.r.t. obj.)
591 TGroupObj* groupObj = *iter;
592 subspace.init(groupObj, c);
593 subspace.extract_matrix(A);
594 subspace.extract_rhs(d); // w/o updates in c => parallel
595 subspace.update_solution(c, omega_relax);
596 }
597}
598
600template<typename TGroupObj, typename TGroupObjIter, typename TDomain, typename TAlgebra>
601void SequentialSubspaceCorrectionLoop(const typename TAlgebra::matrix_type& A,
603 const typename TAlgebra::vector_type& d,
604 number omega_relax,
606 TGroupObjIter objIterBegin, TGroupObjIter objIterEnd)
607{
608 // Loop over all grouping objects.
609 size_t count=0;
610 for(TGroupObjIter iter = objIterBegin; iter != objIterEnd; ++iter)
611 {
612 // Apply subspace correction (w.r.t. obj.)
613 TGroupObj* groupObj = *iter;
614 if (groupObj == NULL)
615 {
616 UG_LOG("WARNING: Found empty object!: "<< count);
617 }
618 else
619 {
620 subspace.init(groupObj, c);
621 subspace.extract_matrix(A);
622 subspace.extract_rhs(d,A,c); // w/ updates for c => sequential
623 subspace.update_solution(c, omega_relax);
624 }
625 count++;
626 }
627}
628
630template<typename TGroupObj, typename TDomain, typename TAlgebra>
631void SequentialSubspaceCorrectionLoopBackward(const typename TAlgebra::matrix_type& A,
633 const typename TAlgebra::vector_type& d,
634 number omega_relax,
636 typename GridFunction<TDomain, TAlgebra>::template traits<TGroupObj>::const_iterator objIterBegin,
637 typename GridFunction<TDomain, TAlgebra>::template traits<TGroupObj>::const_iterator objIterEnd
638)
639{
640 // Loop over all grouping objects.
641 typedef typename GridFunction<TDomain, TAlgebra>::template traits<TGroupObj>::const_iterator GroupObjIter;
642
643 int ncount = 0;
644 for(GroupObjIter iter = objIterBegin; iter != objIterEnd; ++iter)
645 { ncount++; }
646
647 // std::cerr << "Found " << ncount << "objects!" << std::endl;
648
649 std::vector<TGroupObj*> objects(ncount);
650 for(GroupObjIter iter = objIterBegin; iter != objIterEnd; ++iter)
651 {
652 // std::cerr << "Found1 " << *iter << std::endl;
653 objects.push_back(*iter);
654 }
655
656 // std::cerr << "Created vector !" << std::endl;
657
658 for(typename std::vector<TGroupObj*>::reverse_iterator riter = objects.rbegin(); riter != objects.rend(); ++riter)
659 {
660 // Apply subspace correction (w.r.t. obj.)
661 TGroupObj* groupObj = *riter;
662 if (groupObj == NULL) continue;
663
664 //std::cerr << "Found2 " << groupObj << std::endl;
665
666 subspace.init(groupObj, c);
667 subspace.extract_matrix(A);
668 subspace.extract_rhs(d,A,c); // w/ updates for c => sequential
669 subspace.update_solution(c, omega_relax);
670 }
671}
672
673template <class TDomain>
675 public:
676 customLexLess<TDomain>(const typename TDomain::position_accessor_type& aaPos): _aaPos(aaPos){}
677 bool operator()(Vertex* a, Vertex *b) const
678 { return _aaPos[a] < _aaPos[b]; }
679 protected:
680 const typename TDomain::position_accessor_type& _aaPos;
681};
682
683
685template <typename TDomain, typename TAlgebra>
687
688public:
690 typedef TAlgebra algebra_type;
691
693 typedef typename TAlgebra::vector_type vector_type;
694
696 typedef typename TAlgebra::matrix_type matrix_type;
697
700
701
703
705 virtual void init(const vector_type &x) = 0;
706
708 static const int dim = TDomain::dim;
709
710};
711
712
713template <typename TDomain, typename TAlgebra>
714class FullVertexCover: public ISpaceDecomposition<TDomain, TAlgebra>
715{
716public:
717
719
721
722 static const int dim = TDomain::dim;
723
724 typedef typename TGridFunction::element_type TElement;
725
727
730
731 virtual void init(const typename base_type::vector_type &x)
732 {
733 TGridFunction* m_pSol = dynamic_cast<TGridFunction*>(&x);
734 UG_ASSERT(m_pSol !=NULL, "Error: Cast failed!");
735 }
736
737 template <typename TElem>
738 typename TGridFunction::template traits<TElem>::const_iterator
739 begin() const
740 { return m_pSol->template begin<Vertex>();}
741
742 template <typename TElem>
743 typename TGridFunction::template traits<TElem>::const_iterator
744 end() const
745 { return m_pSol->template end<Vertex>(); }
746protected:
748
749};
750
751
753template <typename TDomain, typename TAlgebra>
755{
756public:
758 typedef TAlgebra algebra_type;
759
761 typedef typename TAlgebra::vector_type vector_type;
762
764 typedef typename TAlgebra::matrix_type matrix_type;
765
768
770 typedef ILocalSubspace <TDomain, TAlgebra, Vertex> TVertexSubspace;
771
774
775protected:
779
781
782public:
785
787 SequentialSubspaceCorrection(number relax) : m_relax(relax), m_type("vertex") {};
788
789
792 {
795 newInst->set_debug(debug_writer());
796 newInst->set_damp(this->damping());
797 newInst->set_relax(m_relax);
798 newInst->set_type(m_type);
799 newInst->set_vertex_subspace(m_spVertexSubspace);
800 return newInst;
801 }
802
806
808 virtual bool supports_parallel() const {return true;}
809
811 void set_relax(number omega){ m_relax=omega; }
812
814 void set_type(const std::string& type){ m_type=type; }
815
818 { m_spVertexSubspace = spVertexSubspace; }
819
820
821protected:
823 virtual const char* name() const
824 {return "SequentialSubspaceCorrection";}
825
828 {
829 PROFILE_BEGIN_GROUP(SSC_preprocess, "algebra ssc");
830
831 // Creating overlap 1 matrix and vectors.
832 matrix_type *pA=NULL;
833#ifdef UG_PARALLEL
834 if(pcl::NumProcs() > 1)
835 {
836 UG_ASSERT(0, "SequentialSubspaceCorrection not implemented in parallel. Need to think about this");
837 // We should work with element overlap 1 here!
838 m_A = *pOp;
840 m_oD.set_layouts(m_A.layouts());
841 m_oC.set_layouts(m_A.layouts());
842 m_oD.resize(m_A.num_rows(), false);
843 m_oC.resize(m_A.num_rows(), false);
844 pA = &m_A;
845 }
846 else
847#endif
848 pA = &(*pOp);
849
850
851
852 // Checking.
853 THROW_IF_NOT_EQUAL(pA->num_rows(), pA->num_cols());
854 // UG_COND_THROW(CheckDiagonalInvertible(*pA) == false, name() << ": A has noninvertible diagonal");
855
856 return true;
857 }
858
860 {
861 PROFILE_BEGIN_GROUP(SSC_step, "algebra ssc");
862
863 typedef GridFunction<TDomain, TAlgebra> TGridFunction;
864 TGridFunction* pC = dynamic_cast<TGridFunction*>(&c);
865
866 UG_COND_THROW(pC == NULL, "SequentialSubspaceCorrection expects correction to be a GridFunction.");
867
868 //typedef typename GridFunction<TDomain, TAlgebra>::element_type Element;
869 //typedef typename GridFunction<TDomain, TAlgebra>::side_type Side;
870
871
872 // Set all vector entries to zero.
873 pC->set(0.0);
874
875#ifdef UG_PARALLEL
876 if(pcl::NumProcs() > 1){
877 UG_ASSERT(0, "SequentialSubspaceCorrection not implemented in parallel. Need to think about this");
878
879 // The key problem is that we are iterating over grid functions, not vectors !!!
880 }
881 else
882#endif
883 {
884 matrix_type &A=*pOp;
885 if (m_type == "vertex")
886 {
887 // Prepare
888 typedef typename GridFunction<TDomain, TAlgebra>::template traits<Vertex>::const_iterator TVertexIter;
889 TVertexIter objIterBegin = pC->template begin<Vertex>();
890 TVertexIter objIterEnd = pC->template end<Vertex>();
891
892 m_spVertexSubspace->preprocess(*pC);
893
894 if (true)
895 {
896 int ncount = 0;
897 for(TVertexIter iter = objIterBegin; iter != objIterEnd; ++iter)
898 { ncount++; }
899
900 // std::cerr << "Found " << ncount << "objects!" << std::endl;
901
902 // Create bidirectional container.
903 typedef std::vector<Vertex*> TObjVector;
904 TObjVector objects(ncount);
905 objects.clear();
906 for(TVertexIter iter = objIterBegin; iter != objIterEnd; ++iter)
907 {
908 // std::cerr << "Found1 " << *iter << std::endl;
909 objects.push_back(*iter);
910 }
911
912 // Forward loop.
913 SequentialSubspaceCorrectionLoop<Vertex,TObjVector::const_iterator, TDomain,TAlgebra>
914 (A, *pC, d, m_relax, *m_spVertexSubspace, objects.begin(), objects.end());
915
916 // Backward loop.
917 SequentialSubspaceCorrectionLoop<Vertex,TObjVector::const_reverse_iterator, TDomain,TAlgebra>
918 (A, *pC, d, m_relax, *m_spVertexSubspace, objects.rbegin(), objects.rend());
919
920 } // true
921
922
923 } else if(m_type == "element") {
924
925 /* typedef typename GridFunction<TDomain, TAlgebra>::template traits<Element>::const_iterator TElemIter;
926 SequentialSubspaceCorrectionLoop<Element,TElemIter,TDomain,TAlgebra>
927 (A, *pC, d, m_relax, *m_spElemSubspace, pC->template begin<Element>(), pC->template end<Element>());
928*/
929 }
930 else UG_THROW("SequentialSubspaceCorrectionStep: wrong patch type '"<<m_type<<"'."
931 " Options: vertex.")
932#ifdef UG_PARALLEL
933 pC->set_storage_type(PST_CONSISTENT);
934#endif
935
936 return true;
937 }
938 return false;
939 }
940
942 virtual bool postprocess() {return true;}
943
944
945public:
947protected:
949 std::string m_type;
950
951
953
954#ifdef UG_PARALLEL
955
958
962#endif
963
964
965};
966
967} // end namespace ug
968
969#endif /* __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__ELEMENT_GAUSS_SEIDEL__ */
Definition smart_pointer.h:296
bool invalid() const
returns true if the pointer is invalid, false if not.
Definition smart_pointer.h:420
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
SmartPtr< IDebugWriter< algebra_type > > debug_writer()
returns the debug writer
Definition debug_writer.h:391
Definition densematrix.h:57
void subassign(size_t r, size_t c, const TStorage2 &t)
Definition densematrix.h:215
const value_type & entry(size_type r, size_type c) const
Definition densematrix.h:200
Definition densevector.h:101
void subassign(size_t i, const value_type &t)
Definition densevector.h:290
Definition sequential_subspace_correction.h:715
TGridFunction * m_pSol
Definition sequential_subspace_correction.h:747
TGridFunction::template traits< TElem >::const_iterator begin() const
Definition sequential_subspace_correction.h:739
TGridFunction::element_type TElement
Definition sequential_subspace_correction.h:724
ISpaceDecomposition< TDomain, TAlgebra > base_type
Definition sequential_subspace_correction.h:718
static const int dim
Definition sequential_subspace_correction.h:722
TGridFunction::template traits< TElem >::const_iterator end() const
Definition sequential_subspace_correction.h:744
GridFunction< TDomain, TAlgebra > TGridFunction
Definition sequential_subspace_correction.h:720
virtual void init(const typename base_type::vector_type &x)
Definition sequential_subspace_correction.h:731
represents numerical solutions on a grid using an algebraic vector
Definition grid_function.h:121
void collect_associated(std::vector< TBaseElem * > &vAssElem, TElem *elem, bool clearContainer=true) const
returns the adjacend elements
Definition grid_function.h:290
dim_traits< dim >::grid_base_object element_type
Definition grid_function.h:159
SmartPtr< IDamping< TAlgebra::vector_type, TAlgebra::vector_type > > damping()
returns the scaling
Definition linear_iterator.h:173
Abstract definition for subspace V_k.
Definition sequential_subspace_correction.h:59
virtual void update_solution(vector_type &u, double omega=1.0)=0
u = u + omega*ck.
virtual void init(TObject *, const vector_type &)=0
Extract local data (based on group obj).
virtual void extract_rhs(const vector_type &d)=0
Extract rhs (on local index set) for parallel subspace correction.
virtual void extract_rhs(const vector_type &d, const matrix_type &A, const vector_type &c)=0
Extract rhs (on local index set) for sequential subspace correction.
DenseVector< VariableArray1< number > > vector_type_local
Definition sequential_subspace_correction.h:72
GridFunction< TDomain, TAlgebra > grid_function_type
Grid function type.
Definition sequential_subspace_correction.h:71
virtual bool preprocess(const vector_type &c)
Called once.
Definition sequential_subspace_correction.h:79
TAlgebra::vector_type vector_type
Vector type.
Definition sequential_subspace_correction.h:65
DenseMatrix< VariableArray2< number > > matrix_type_local
Definition sequential_subspace_correction.h:73
virtual void extract_matrix(const matrix_type &A)=0
Extract matrix (on local index set).
virtual size_t size()
Definition sequential_subspace_correction.h:98
TAlgebra::matrix_type matrix_type
Matrix type.
Definition sequential_subspace_correction.h:68
TAlgebra algebra_type
Algebra type.
Definition sequential_subspace_correction.h:62
describes a linear iterator that is based on a matrix operator
Definition preconditioner.h:103
Abstract.
Definition sequential_subspace_correction.h:686
virtual void init(const vector_type &x)=0
Initialize decomposition.
TAlgebra::matrix_type matrix_type
Matrix type.
Definition sequential_subspace_correction.h:696
TAlgebra algebra_type
Algebra type.
Definition sequential_subspace_correction.h:690
TAlgebra::vector_type vector_type
Vector type.
Definition sequential_subspace_correction.h:693
GridFunction< TDomain, TAlgebra > grid_function_type
Grid function type.
Definition sequential_subspace_correction.h:699
static const int dim
World Dimension.
Definition sequential_subspace_correction.h:708
Abstract definition for subspace V_k (based on DoFIndex)
Definition sequential_subspace_correction.h:233
matrix_type_local m_Aloc
Definition sequential_subspace_correction.h:340
std::vector< DoFIndex > m_vInd
Algebraic indices.
Definition sequential_subspace_correction.h:335
DenseMatrix< VariableArray2< number > > matrix_type_local
Definition sequential_subspace_correction.h:247
std::vector< DoFIndex > index_vector
Local data.
Definition sequential_subspace_correction.h:245
bool check(void *obj) const
Definition sequential_subspace_correction.h:252
TAlgebra::matrix_type matrix_type
Matrix type.
Definition sequential_subspace_correction.h:242
vector_type_local m_uloc
Definition sequential_subspace_correction.h:338
DenseVector< VariableArray1< number > > vector_type_local
Definition sequential_subspace_correction.h:246
virtual void extract_rhs(const vector_type &d, const matrix_type &A, const vector_type &c)
Extract rhs (on local index set) for sequential subspace correction.
Definition sequential_subspace_correction.h:292
virtual void update_solution(vector_type &u, double omega=1.0)
u = u + ck
Definition sequential_subspace_correction.h:316
vector_type_local m_dloc
Definition sequential_subspace_correction.h:339
TAlgebra::vector_type vector_type
Vector type.
Definition sequential_subspace_correction.h:239
virtual size_t size()
Definition sequential_subspace_correction.h:331
virtual void extract_rhs(const vector_type &d)
Extract rhs (on local index set) for parallel subspace correction.
Definition sequential_subspace_correction.h:280
virtual void extract_matrix(const matrix_type &A)
Extract matrix (on local index set)
Definition sequential_subspace_correction.h:259
TAlgebra algebra_type
Algebra type.
Definition sequential_subspace_correction.h:236
virtual void init(TObject *, const vector_type &)=0
Extract local data (based on group obj)
Concrete definition of subspace V_k (based on size_t)
Definition sequential_subspace_correction.h:106
std::vector< size_t > index_vector
Local data.
Definition sequential_subspace_correction.h:118
DenseMatrix< VariableArray2< number > > matrix_type_local
Definition sequential_subspace_correction.h:120
matrix_type_local m_Aloc
Definition sequential_subspace_correction.h:225
virtual void update_solution(vector_type &u, double omega=1.0)
u = u + ck
Definition sequential_subspace_correction.h:200
virtual void extract_matrix(const matrix_type &A)
Extract matrix (on local index set)
Definition sequential_subspace_correction.h:137
vector_type_local m_uloc
Memory for local algebra.
Definition sequential_subspace_correction.h:223
index_vector m_vInd
Algebraic indices.
Definition sequential_subspace_correction.h:220
TAlgebra::matrix_type matrix_type
Matrix type.
Definition sequential_subspace_correction.h:115
DenseVector< VariableArray1< number > > vector_type_local
Definition sequential_subspace_correction.h:119
virtual void extract_rhs(const vector_type &d)
Extract rhs (on local index set) for parallel subspace correction.
Definition sequential_subspace_correction.h:161
TAlgebra::vector_type vector_type
Vector type.
Definition sequential_subspace_correction.h:112
virtual void init(TObject *, const vector_type &)=0
Called once.
virtual size_t size()
Definition sequential_subspace_correction.h:215
vector_type_local m_dloc
Definition sequential_subspace_correction.h:224
virtual void extract_rhs(const vector_type &d, const matrix_type &A, const vector_type &c)
Extract rhs (on local index set) for sequential subspace correction.
Definition sequential_subspace_correction.h:177
TAlgebra algebra_type
Algebra type.
Definition sequential_subspace_correction.h:109
Sequential subspace correction preconditioner.
Definition sequential_subspace_correction.h:755
TAlgebra algebra_type
Algebra type.
Definition sequential_subspace_correction.h:758
TAlgebra::vector_type vector_type
Vector type.
Definition sequential_subspace_correction.h:761
vector_type m_oD
for overlaps only
Definition sequential_subspace_correction.h:960
virtual const char * name() const
Name of preconditioner.
Definition sequential_subspace_correction.h:823
TAlgebra::matrix_type matrix_type
Matrix type.
Definition sequential_subspace_correction.h:764
std::string m_type
Definition sequential_subspace_correction.h:949
vector_type m_oC
Definition sequential_subspace_correction.h:961
number m_relax
Definition sequential_subspace_correction.h:948
GridFunction< TDomain, TAlgebra > grid_function_type
Definition sequential_subspace_correction.h:780
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 sequential_subspace_correction.h:859
matrix_type m_A
matrix with overlap
Definition sequential_subspace_correction.h:957
virtual bool supports_parallel() const
returns if parallel solving is supported
Definition sequential_subspace_correction.h:808
void set_vertex_subspace(SmartPtr< TVertexSubspace > spVertexSubspace)
set subspace
Definition sequential_subspace_correction.h:817
SequentialSubspaceCorrection(number relax)
constructor setting relaxation
Definition sequential_subspace_correction.h:787
virtual bool preprocess(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp)
Preprocess routine.
Definition sequential_subspace_correction.h:827
virtual bool postprocess()
Postprocess routine.
Definition sequential_subspace_correction.h:942
virtual ~SequentialSubspaceCorrection()
Destructor.
Definition sequential_subspace_correction.h:804
ILocalSubspace< TDomain, TAlgebra, Vertex > TVertexSubspace
Subspace for vertices.
Definition sequential_subspace_correction.h:770
void set_type(const std::string &type)
set type
Definition sequential_subspace_correction.h:814
virtual SmartPtr< ILinearIterator< vector_type > > clone()
Clone.
Definition sequential_subspace_correction.h:791
IPreconditioner< TAlgebra >::matrix_operator_type matrix_operator_type
Matrix Operator type.
Definition sequential_subspace_correction.h:767
SmartPtr< TVertexSubspace > m_spVertexSubspace
Definition sequential_subspace_correction.h:952
IPreconditioner< TAlgebra > base_type
Base type.
Definition sequential_subspace_correction.h:773
void set_relax(number omega)
set relaxation parameter
Definition sequential_subspace_correction.h:811
SequentialSubspaceCorrection()
default constructor
Definition sequential_subspace_correction.h:784
GridFunction< TDomain, TAlgebra >::element_type Element
Definition sequential_subspace_correction.h:946
bool resize(size_type n, bool bCopyValues=true)
Definition variable_array_impl.h:92
bool resize(size_type newRows, size_type newCols, bool bCopyValues=true)
Definition variable_array_impl.h:237
Collects indices on all elements with v \in Vtx(elem)
Definition sequential_subspace_correction.h:349
void init(void *obj, const vector_type &cvec)
Extract indices for local DoFs.
Definition sequential_subspace_correction.h:376
LocalIndexSubspace< TDomain, TAlgebra, Vertex > base_type
Base type.
Definition sequential_subspace_correction.h:364
TAlgebra::matrix_type matrix_type
Matrix type.
Definition sequential_subspace_correction.h:361
Vertex TObject
Definition sequential_subspace_correction.h:351
bool check(void *obj) const
Definition sequential_subspace_correction.h:372
TAlgebra::vector_type vector_type
Vector type.
Definition sequential_subspace_correction.h:358
TAlgebra algebra_type
Algebra type.
Definition sequential_subspace_correction.h:355
Collects indices on all elements with v \in Vtx(elem)
Definition sequential_subspace_correction.h:412
TAlgebra::vector_type vector_type
Vector type.
Definition sequential_subspace_correction.h:418
void extract_ids(const grid_function_type *c)
Extracts function IDs.
Definition sequential_subspace_correction.h:437
TAlgebra algebra_type
Algebra type.
Definition sequential_subspace_correction.h:415
void init(Vertex *groupObj, const vector_type &cvec)
Extract indices for single vertex.
Definition sequential_subspace_correction.h:469
TAlgebra::matrix_type matrix_type
Matrix type.
Definition sequential_subspace_correction.h:421
std::vector< std::string > m_vVtxCmp
Definition sequential_subspace_correction.h:512
GridFunction< TDomain, TAlgebra > grid_function_type
Definition sequential_subspace_correction.h:423
bool preprocess(const vector_type &cvec)
OVERRIDE.
Definition sequential_subspace_correction.h:457
std::vector< size_t > m_vElemFct
Definition sequential_subspace_correction.h:516
LocalDoFSubspace< TDomain, TAlgebra, Vertex > base_type
Base type.
Definition sequential_subspace_correction.h:426
std::vector< size_t > m_vVtxFct
Definition sequential_subspace_correction.h:515
std::vector< std::string > m_vElemCmp
Definition sequential_subspace_correction.h:513
Base-class for all vertex-types.
Definition grid_base_objects.h:231
Definition sequential_subspace_correction.h:674
const TDomain::position_accessor_type & _aaPos
Definition sequential_subspace_correction.h:680
bool operator()(Vertex *a, Vertex *b) const
Definition sequential_subspace_correction.h:677
@ PST_CONSISTENT
Definition parallel_storage_type.h:68
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
#define UG_LOG(msg)
Definition log.h:367
#define THROW_IF_NOT_EQUAL(s1, s2)
Definition error.h:181
#define UG_COND_THROW(cond, msg)
UG_COND_THROW(cond, msg) : performs a UG_THROW(msg) if cond == true.
Definition error.h:61
double number
Definition types.h:124
the ug namespace
void SequentialSubspaceCorrectionLoop(const typename TAlgebra::matrix_type &A, GridFunction< TDomain, TAlgebra > &c, const typename TAlgebra::vector_type &d, number omega_relax, ILocalSubspace< TDomain, TAlgebra, TGroupObj > &subspace, TGroupObjIter objIterBegin, TGroupObjIter objIterEnd)
Abstract loop over TGroupObj (forward)
Definition sequential_subspace_correction.h:601
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 CreateOverlap(TMatrix &matInOut)
Adds coefficients and connections to create an overlapping matrix.
Definition matrix_overlap_impl.h:434
void ParallelSubspaceCorrectionLoop(const typename TAlgebra::matrix_type &A, GridFunction< TDomain, TAlgebra > &c, const typename TAlgebra::vector_type &d, number omega_relax, ILocalSubspace< TDomain, TAlgebra, TGroupObj > &subspace, typename GridFunction< TDomain, TAlgebra >::template traits< TGroupObj >::const_iterator objIterBegin, typename GridFunction< TDomain, TAlgebra >::template traits< TGroupObj >::const_iterator objIterEnd)
Abstract loop.
Definition sequential_subspace_correction.h:577
void SequentialSubspaceCorrectionLoopBackward(const typename TAlgebra::matrix_type &A, GridFunction< TDomain, TAlgebra > &c, const typename TAlgebra::vector_type &d, number omega_relax, ILocalSubspace< TDomain, TAlgebra, TGroupObj > &subspace, typename GridFunction< TDomain, TAlgebra >::template traits< TGroupObj >::const_iterator objIterBegin, typename GridFunction< TDomain, TAlgebra >::template traits< TGroupObj >::const_iterator objIterEnd)
Abstract loop (backward)
Definition sequential_subspace_correction.h:631
bool InverseMatMult(number &dest, const double &beta, const TMat &mat, const TVec &vec)
you can implement this function with GetInverse and MatMult
number & DoFRef(TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
Definition multi_index.h:276
#define PROFILE_BEGIN_GROUP(name, groups)
Definition profiler.h:255