Loading [MathJax]/extensions/tex2jax.js

Dune Core Modules (unstable)

aggregates.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_AMG_AGGREGATES_HH
6#define DUNE_AMG_AGGREGATES_HH
7
8
9#include "parameters.hh"
10#include "graph.hh"
11#include "properties.hh"
12#include "combinedfunctor.hh"
13
14#include <dune/common/timer.hh>
17#include <dune/common/sllist.hh>
21
22#include <utility>
23#include <set>
24#include <algorithm>
25#include <complex>
26#include <limits>
27#include <ostream>
28#include <tuple>
29#include <cmath>
30
31namespace Dune
32{
33 namespace Amg
34 {
35
49 template<class T>
50 class AggregationCriterion : public T
51 {
52
53 public:
58
69 : T()
70 {}
71
73 : T(parms)
74 {}
84 void setDefaultValuesIsotropic(std::size_t dim, std::size_t diameter=2)
85 {
86 this->setMaxDistance(diameter-1);
87 std::size_t csize=1;
88
89 for(; dim>0; dim--) {
90 csize*=diameter;
91 this->setMaxDistance(this->maxDistance()+diameter-1);
92 }
93 this->setMinAggregateSize(csize);
94 this->setMaxAggregateSize(static_cast<std::size_t>(csize*1.5));
95 }
96
107 void setDefaultValuesAnisotropic(std::size_t dim,std::size_t diameter=2)
108 {
109 setDefaultValuesIsotropic(dim, diameter);
110 this->setMaxDistance(this->maxDistance()+dim-1);
111 }
112 };
113
114 template<class T>
115 std::ostream& operator<<(std::ostream& os, const AggregationCriterion<T>& criterion)
116 {
117 os<<"{ maxdistance="<<criterion.maxDistance()<<" minAggregateSize="
118 <<criterion.minAggregateSize()<< " maxAggregateSize="<<criterion.maxAggregateSize()
119 <<" connectivity="<<criterion.maxConnectivity()<<" debugLevel="<<criterion.debugLevel()<<"}";
120 return os;
121 }
122
134 template<class M, class N>
136 {
137 public:
141 typedef M Matrix;
142
146 typedef N Norm;
147
151 typedef typename Matrix::row_type Row;
152
157
158 void init(const Matrix* matrix);
159
160 void initRow(const Row& row, int index);
161
162 void examine(const ColIter& col);
163
164 template<class G>
165 void examine(G& graph, const typename G::EdgeIterator& edge, const ColIter& col);
166
167 bool isIsolated();
168
169
171 : Parameters(parms)
172 {}
174 : Parameters()
175 {}
176
177 protected:
182 typedef typename FieldTraits<field_type>::real_type real_type;
183 real_type maxValue_;
187 int row_;
189 real_type diagonal_;
190 std::vector<real_type> vals_;
191 typename std::vector<real_type>::iterator valIter_;
192
193 };
194
195
196 template<class M, class N>
197 inline void SymmetricMatrixDependency<M,N>::init(const Matrix* matrix)
198 {
199 matrix_ = matrix;
200 }
201
202 template<class M, class N>
203 inline void SymmetricMatrixDependency<M,N>::initRow(const Row& row, int index)
204 {
205 using std::min;
206 vals_.assign(row.size(), 0.0);
207 assert(vals_.size()==row.size());
208 valIter_=vals_.begin();
209
211 diagonal_=norm_(row[index]);
212 row_ = index;
213 }
214
215 template<class M, class N>
216 inline void SymmetricMatrixDependency<M,N>::examine(const ColIter& col)
217 {
218 using std::max;
219 // skip positive offdiagonals if norm preserves sign of them.
220 real_type eij = norm_(*col);
221 if(!N::is_sign_preserving || eij<0) // || eji<0)
222 {
223 *valIter_ = eij/diagonal_*eij/norm_(matrix_->operator[](col.index())[col.index()]);
224 maxValue_ = max(maxValue_, *valIter_);
225 }else
226 *valIter_ =0;
227 ++valIter_;
228 }
229
230 template<class M, class N>
231 template<class G>
232 inline void SymmetricMatrixDependency<M,N>::examine(G&, const typename G::EdgeIterator& edge, const ColIter&)
233 {
234 if(*valIter_ > alpha() * maxValue_) {
235 edge.properties().setDepends();
236 edge.properties().setInfluences();
237 }
238 ++valIter_;
239 }
240
241 template<class M, class N>
242 inline bool SymmetricMatrixDependency<M,N>::isIsolated()
243 {
244 if(diagonal_==0)
245 DUNE_THROW(Dune::ISTLError, "No diagonal entry for row "<<row_<<".");
246 valIter_=vals_.begin();
247 return maxValue_ < beta();
248 }
249
253 template<class M, class N>
254 class Dependency : public Parameters
255 {
256 public:
260 typedef M Matrix;
261
265 typedef N Norm;
266
270 typedef typename Matrix::row_type Row;
271
276
277 void init(const Matrix* matrix);
278
279 void initRow(const Row& row, int index);
280
281 void examine(const ColIter& col);
282
283 template<class G>
284 void examine(G& graph, const typename G::EdgeIterator& edge, const ColIter& col);
285
286 bool isIsolated();
287
288 Dependency(const Parameters& parms)
289 : Parameters(parms)
290 {}
291
292 Dependency()
293 : Parameters()
294 {}
295
296 protected:
301 typedef typename FieldTraits<field_type>::real_type real_type;
302 real_type maxValue_;
306 int row_;
308 real_type diagonal_;
309 };
310
314 template<class M, class N>
316 {
317 public:
321 typedef M Matrix;
322
326 typedef N Norm;
327
331 typedef typename Matrix::row_type Row;
332
337
338 void init(const Matrix* matrix);
339
340 void initRow(const Row& row, int index);
341
342 void examine(const ColIter& col);
343
344 template<class G>
345 void examine(G& graph, const typename G::EdgeIterator& edge, const ColIter& col);
346
347 bool isIsolated();
348
349
350 SymmetricDependency(const Parameters& parms)
351 : Parameters(parms)
352 {}
354 : Parameters()
355 {}
356
357 protected:
362 typedef typename FieldTraits<field_type>::real_type real_type;
363 real_type maxValue_;
367 int row_;
369 real_type diagonal_;
370 private:
371 void initRow(const Row& row, int index, const std::true_type&);
372 void initRow(const Row& row, int index, const std::false_type&);
373 };
374
379 template<int N>
381 {
382 public:
383 enum { /* @brief We preserve the sign.*/
384 is_sign_preserving = true
385 };
386
391 template<class M>
392 typename FieldTraits<typename M::field_type>::real_type operator()(const M& m,
393 [[maybe_unused]] typename std::enable_if_t<!Dune::IsNumber<M>::value>* sfinae = nullptr) const
394 {
395 typedef typename M::field_type field_type;
396 typedef typename FieldTraits<field_type>::real_type real_type;
397 static_assert( std::is_convertible<field_type, real_type >::value,
398 "use of diagonal norm in AMG not implemented for complex field_type");
399 return m[N][N];
400 // possible implementation for complex types: return signed_abs(m[N][N]);
401 }
402
407 template<class M>
408 auto operator()(const M& m,
409 typename std::enable_if_t<Dune::IsNumber<M>::value>* sfinae = nullptr) const
410 {
411 typedef typename FieldTraits<M>::real_type real_type;
412 static_assert( std::is_convertible<M, real_type >::value,
413 "use of diagonal norm in AMG not implemented for complex field_type");
414 return m;
415 // possible implementation for complex types: return signed_abs(m[N][N]);
416 }
417
418 private:
419
421 template<typename T>
422 static T signed_abs(const T & v)
423 {
424 return v;
425 }
426
428 template<typename T>
429 static T signed_abs(const std::complex<T> & v)
430 {
431 // return sign * abs_value
432 // in case of complex numbers this extends to using the csgn function to determine the sign
433 return csgn(v) * std::abs(v);
434 }
435
437 template<typename T>
438 static T csgn(const T & v)
439 {
440 return (T(0) < v) - (v < T(0));
441 }
442
444 template<typename T>
445 static T csgn(std::complex<T> a)
446 {
447 return csgn(a.real())+(a.real() == 0.0)*csgn(a.imag());
448 }
449
450 };
451
456 class FirstDiagonal : public Diagonal<0>
457 {};
458
464 struct RowSum
465 {
466
467 enum { /* @brief We preserve the sign.*/
468 is_sign_preserving = false
469 };
474 template<class M>
475 auto operator()(const M& m) const
476 {
477 using std::abs;
478 if constexpr(Dune::IsNumber<M>::value)
479 return abs(m);
480 else
481 return m.infinity_norm();
482 }
483 };
484
485 struct FrobeniusNorm
486 {
487
488 enum { /* @brief We preserve the sign.*/
489 is_sign_preserving = false
490 };
495 template<class M>
496 typename FieldTraits<typename M::field_type>::real_type operator()(const M& m) const
497 {
498 return m.frobenius_norm();
499 }
500 };
501 struct AlwaysOneNorm
502 {
503
504 enum { /* @brief We preserve the sign.*/
505 is_sign_preserving = false
506 };
511 template<class M>
512 typename FieldTraits<typename M::field_type>::real_type operator()(const M& m) const
513 {
514 return 1;
515 }
516 };
523 template<class M, class Norm>
524 class SymmetricCriterion : public AggregationCriterion<SymmetricDependency<M,Norm> >
525 {
526 public:
527 SymmetricCriterion(const Parameters& parms)
529 {}
531 {}
532 };
533
534
543 template<class M, class Norm>
544 class UnSymmetricCriterion : public AggregationCriterion<Dependency<M,Norm> >
545 {
546 public:
547 UnSymmetricCriterion(const Parameters& parms)
549 {}
551 {}
552 };
553 // forward declaration
554 template<class G> class Aggregator;
555
556
564 template<class V>
566 {
567 public:
568
572 static const V UNAGGREGATED;
573
577 static const V ISOLATED;
582
587
593
599
604 {
605 public:
606 template<class EdgeIterator>
607 void operator()([[maybe_unused]] const EdgeIterator& edge) const
608 {}
609 };
610
611
616
623
628
640 template<class M, class G, class C>
641 std::tuple<int,int,int,int> buildAggregates(const M& matrix, G& graph, const C& criterion,
642 bool finestLevel);
643
661 template<bool reset, class G, class F, class VM>
662 std::size_t breadthFirstSearch(const VertexDescriptor& start,
663 const AggregateDescriptor& aggregate,
664 const G& graph,
665 F& aggregateVisitor,
666 VM& visitedMap) const;
667
691 template<bool remove, bool reset, class G, class L, class F1, class F2, class VM>
692 std::size_t breadthFirstSearch(const VertexDescriptor& start,
693 const AggregateDescriptor& aggregate,
694 const G& graph, L& visited, F1& aggregateVisitor,
695 F2& nonAggregateVisitor,
696 VM& visitedMap) const;
697
703 void allocate(std::size_t noVertices);
704
708 std::size_t noVertices() const;
709
713 void free();
714
721
728
729 typedef const AggregateDescriptor* const_iterator;
730
731 const_iterator begin() const
732 {
733 return aggregates_;
734 }
735
736 const_iterator end() const
737 {
738 return aggregates_+noVertices();
739 }
740
741 typedef AggregateDescriptor* iterator;
742
743 iterator begin()
744 {
745 return aggregates_;
746 }
747
748 iterator end()
749 {
750 return aggregates_+noVertices();
751 }
752 private:
754 AggregatesMap(const AggregatesMap<V>&) = delete;
756 AggregatesMap<V>& operator=(const AggregatesMap<V>&) = delete;
757
761 AggregateDescriptor* aggregates_;
762
766 std::size_t noVertices_;
767 };
768
772 template<class G, class C>
773 void buildDependency(G& graph,
774 const typename C::Matrix& matrix,
775 C criterion,
776 bool finestLevel);
777
782 template<class G, class S>
784 {
785
786 public:
787
788 /***
789 * @brief The type of the matrix graph we work with.
790 */
791 typedef G MatrixGraph;
796
802
807 typedef S VertexSet;
808
810 typedef typename VertexSet::const_iterator const_iterator;
811
815 typedef std::size_t* SphereMap;
816
825 Aggregate(MatrixGraph& graph, AggregatesMap<Vertex>& aggregates,
826 VertexSet& connectivity, std::vector<Vertex>& front_);
827
828 void invalidate()
829 {
830 --id_;
831 }
832
840
844 void seed(const Vertex& vertex);
845
849 void add(const Vertex& vertex);
850
851 void add(std::vector<Vertex>& vertex);
855 void clear();
856
860 typename VertexSet::size_type size();
864 typename VertexSet::size_type connectSize();
865
869 int id();
870
873
876
877 private:
881 VertexSet vertices_;
882
887 int id_;
888
892 MatrixGraph& graph_;
893
897 AggregatesMap<Vertex>& aggregates_;
898
902 VertexSet& connected_;
903
907 std::vector<Vertex>& front_;
908 };
909
913 template<class G>
915 {
916 public:
917
921 typedef G MatrixGraph;
922
927
930
935
940
957 template<class M, class C>
958 std::tuple<int,int,int,int> build(const M& m, G& graph,
959 AggregatesMap<Vertex>& aggregates, const C& c,
960 bool finestLevel);
961 private:
967
972
976 typedef std::set<Vertex,std::less<Vertex>,Allocator> VertexSet;
977
981 typedef std::size_t* SphereMap;
982
986 MatrixGraph* graph_;
987
992
996 std::vector<Vertex> front_;
997
1001 VertexSet connected_;
1002
1006 int size_;
1007
1011 class Stack
1012 {
1013 public:
1014 static const Vertex NullEntry;
1015
1016 Stack(const MatrixGraph& graph,
1017 const Aggregator<G>& aggregatesBuilder,
1018 const AggregatesMap<Vertex>& aggregates);
1019 ~Stack();
1020 Vertex pop();
1021 private:
1022 enum { N = 1300000 };
1023
1025 const MatrixGraph& graph_;
1027 const Aggregator<G>& aggregatesBuilder_;
1029 const AggregatesMap<Vertex>& aggregates_;
1031 int size_;
1032 Vertex maxSize_;
1034 typename MatrixGraph::ConstVertexIterator begin_;
1036
1038 Vertex* vals_;
1039
1040 };
1041
1042 friend class Stack;
1043
1054 template<class V>
1055 void visitAggregateNeighbours(const Vertex& vertex, const AggregateDescriptor& aggregate,
1056 const AggregatesMap<Vertex>& aggregates,
1057 V& visitor) const;
1058
1063 template<class V>
1064 class AggregateVisitor
1065 {
1066 public:
1070 typedef V Visitor;
1079 Visitor& visitor);
1080
1087 void operator()(const typename MatrixGraph::ConstEdgeIterator& edge);
1088
1089 private:
1091 const AggregatesMap<Vertex>& aggregates_;
1093 AggregateDescriptor aggregate_;
1095 Visitor* visitor_;
1096 };
1097
1101 class Counter
1102 {
1103 public:
1107 int value();
1108
1109 protected:
1114
1115 private:
1116 int count_;
1117 };
1118
1119
1124 class FrontNeighbourCounter : public Counter
1125 {
1126 public:
1132
1133 void operator()(const typename MatrixGraph::ConstEdgeIterator& edge);
1134
1135 private:
1136 const MatrixGraph& graph_;
1137 };
1138
1143 int noFrontNeighbours(const Vertex& vertex) const;
1144
1148 class TwoWayCounter : public Counter
1149 {
1150 public:
1151 void operator()(const typename MatrixGraph::ConstEdgeIterator& edge);
1152 };
1153
1165 int twoWayConnections(const Vertex&, const AggregateDescriptor& aggregate,
1166 const AggregatesMap<Vertex>& aggregates) const;
1167
1171 class OneWayCounter : public Counter
1172 {
1173 public:
1174 void operator()(const typename MatrixGraph::ConstEdgeIterator& edge);
1175 };
1176
1188 int oneWayConnections(const Vertex&, const AggregateDescriptor& aggregate,
1189 const AggregatesMap<Vertex>& aggregates) const;
1190
1197 class ConnectivityCounter : public Counter
1198 {
1199 public:
1206 ConnectivityCounter(const VertexSet& connected, const AggregatesMap<Vertex>& aggregates);
1207
1208 void operator()(const typename MatrixGraph::ConstEdgeIterator& edge);
1209
1210 private:
1212 const VertexSet& connected_;
1214 const AggregatesMap<Vertex>& aggregates_;
1215
1216 };
1217
1229 double connectivity(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates) const;
1237 bool connected(const Vertex& vertex, const AggregateDescriptor& aggregate,
1238 const AggregatesMap<Vertex>& aggregates) const;
1239
1247 bool connected(const Vertex& vertex, const SLList<AggregateDescriptor>& aggregateList,
1248 const AggregatesMap<Vertex>& aggregates) const;
1249
1257 class DependencyCounter : public Counter
1258 {
1259 public:
1264
1265 void operator()(const typename MatrixGraph::ConstEdgeIterator& edge);
1266 };
1267
1274 class FrontMarker
1275 {
1276 public:
1283 FrontMarker(std::vector<Vertex>& front, MatrixGraph& graph);
1284
1285 void operator()(const typename MatrixGraph::ConstEdgeIterator& edge);
1286
1287 private:
1289 std::vector<Vertex>& front_;
1291 MatrixGraph& graph_;
1292 };
1293
1297 void unmarkFront();
1298
1313 int unusedNeighbours(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates) const;
1314
1328 std::pair<int,int> neighbours(const Vertex& vertex,
1329 const AggregateDescriptor& aggregate,
1330 const AggregatesMap<Vertex>& aggregates) const;
1347 int aggregateNeighbours(const Vertex& vertex, const AggregateDescriptor& aggregate, const AggregatesMap<Vertex>& aggregates) const;
1348
1356 bool admissible(const Vertex& vertex, const AggregateDescriptor& aggregate, const AggregatesMap<Vertex>& aggregates) const;
1357
1365 std::size_t distance(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates);
1366
1375 Vertex mergeNeighbour(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates) const;
1376
1385 void nonisoNeighbourAggregate(const Vertex& vertex,
1386 const AggregatesMap<Vertex>& aggregates,
1387 SLList<Vertex>& neighbours) const;
1388
1396 template<class C>
1397 void growAggregate(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates, const C& c);
1398 template<class C>
1399 void growIsolatedAggregate(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates, const C& c);
1400 };
1401
1402#ifndef DOXYGEN
1403
1404 template<class M, class N>
1405 inline void SymmetricDependency<M,N>::init(const Matrix* matrix)
1406 {
1407 matrix_ = matrix;
1408 }
1409
1410 template<class M, class N>
1411 inline void SymmetricDependency<M,N>::initRow(const Row& row, int index)
1412 {
1413 initRow(row, index, std::is_convertible<field_type, real_type>());
1414 }
1415
1416 template<class M, class N>
1417 inline void SymmetricDependency<M,N>::initRow(const Row& row, int index, const std::false_type&)
1418 {
1419 DUNE_THROW(InvalidStateException, "field_type needs to convertible to real_type");
1420 }
1421
1422 template<class M, class N>
1423 inline void SymmetricDependency<M,N>::initRow([[maybe_unused]] const Row& row, int index, const std::true_type&)
1424 {
1425 using std::min;
1427 row_ = index;
1428 diagonal_ = norm_(matrix_->operator[](row_)[row_]);
1429 }
1430
1431 template<class M, class N>
1432 inline void SymmetricDependency<M,N>::examine(const ColIter& col)
1433 {
1434 using std::max;
1435 real_type eij = norm_(*col);
1436 typename Matrix::ConstColIterator opposite_entry =
1437 matrix_->operator[](col.index()).find(row_);
1438 if ( opposite_entry == matrix_->operator[](col.index()).end() )
1439 {
1440 // Consider this a weak connection we disregard.
1441 return;
1442 }
1443 real_type eji = norm_(*opposite_entry);
1444
1445 // skip positive offdiagonals if norm preserves sign of them.
1446 if(!N::is_sign_preserving || eij<0 || eji<0)
1447 maxValue_ = max(maxValue_,
1448 eij /diagonal_ * eji/
1449 norm_(matrix_->operator[](col.index())[col.index()]));
1450 }
1451
1452 template<class M, class N>
1453 template<class G>
1454 inline void SymmetricDependency<M,N>::examine(G& graph, const typename G::EdgeIterator& edge, const ColIter& col)
1455 {
1456 real_type eij = norm_(*col);
1457 typename Matrix::ConstColIterator opposite_entry =
1458 matrix_->operator[](col.index()).find(row_);
1459
1460 if ( opposite_entry == matrix_->operator[](col.index()).end() )
1461 {
1462 // Consider this as a weak connection we disregard.
1463 return;
1464 }
1465 real_type eji = norm_(*opposite_entry);
1466 // skip positive offdiagonals if norm preserves sign of them.
1467 if(!N::is_sign_preserving || (eij<0 || eji<0))
1468 if(eji / norm_(matrix_->operator[](edge.target())[edge.target()]) *
1469 eij/ diagonal_ > alpha() * maxValue_) {
1470 edge.properties().setDepends();
1471 edge.properties().setInfluences();
1472 typename G::EdgeProperties& other = graph.getEdgeProperties(edge.target(), edge.source());
1473 other.setInfluences();
1474 other.setDepends();
1475 }
1476 }
1477
1478 template<class M, class N>
1479 inline bool SymmetricDependency<M,N>::isIsolated()
1480 {
1481 return maxValue_ < beta();
1482 }
1483
1484
1485 template<class M, class N>
1486 inline void Dependency<M,N>::init(const Matrix* matrix)
1487 {
1488 matrix_ = matrix;
1489 }
1490
1491 template<class M, class N>
1492 inline void Dependency<M,N>::initRow([[maybe_unused]] const Row& row, int index)
1493 {
1494 using std::min;
1496 row_ = index;
1497 diagonal_ = norm_(matrix_->operator[](row_)[row_]);
1498 }
1499
1500 template<class M, class N>
1501 inline void Dependency<M,N>::examine(const ColIter& col)
1502 {
1503 using std::max;
1504 maxValue_ = max(maxValue_, -norm_(*col));
1505 }
1506
1507 template<class M, class N>
1508 template<class G>
1509 inline void Dependency<M,N>::examine(G& graph, const typename G::EdgeIterator& edge, const ColIter& col)
1510 {
1511 if(-norm_(*col) >= maxValue_ * alpha()) {
1512 edge.properties().setDepends();
1513 typedef typename G::EdgeDescriptor ED;
1514 ED e= graph.findEdge(edge.target(), edge.source());
1516 {
1517 typename G::EdgeProperties& other = graph.getEdgeProperties(e);
1518 other.setInfluences();
1519 }
1520 }
1521 }
1522
1523 template<class M, class N>
1524 inline bool Dependency<M,N>::isIsolated()
1525 {
1526 return maxValue_ < beta() * diagonal_;
1527 }
1528
1529 template<class G,class S>
1530 Aggregate<G,S>::Aggregate(MatrixGraph& graph, AggregatesMap<Vertex>& aggregates,
1531 VertexSet& connected, std::vector<Vertex>& front)
1532 : vertices_(), id_(-1), graph_(graph), aggregates_(aggregates),
1533 connected_(connected), front_(front)
1534 {}
1535
1536 template<class G,class S>
1538 {
1539 /*
1540 vertices_.push_back(vertex);
1541 typedef typename VertexList::const_iterator iterator;
1542 iterator begin = vertices_.begin();
1543 iterator end = vertices_.end();*/
1544 throw "Not yet implemented";
1545
1546 // while(begin!=end){
1547 //for();
1548 // }
1549
1550 }
1551
1552 template<class G,class S>
1553 inline void Aggregate<G,S>::seed(const Vertex& vertex)
1554 {
1555 dvverb<<"Connected cleared"<<std::endl;
1556 connected_.clear();
1557 vertices_.clear();
1558 connected_.insert(vertex);
1559 dvverb << " Inserting "<<vertex<<" size="<<connected_.size();
1560 ++id_ ;
1561 add(vertex);
1562 }
1563
1564
1565 template<class G,class S>
1566 inline void Aggregate<G,S>::add(const Vertex& vertex)
1567 {
1568 vertices_.insert(vertex);
1569 aggregates_[vertex]=id_;
1570 if(front_.size())
1571 front_.erase(std::lower_bound(front_.begin(), front_.end(), vertex));
1572
1573
1574 typedef typename MatrixGraph::ConstEdgeIterator iterator;
1575 const iterator end = graph_.endEdges(vertex);
1576 for(iterator edge = graph_.beginEdges(vertex); edge != end; ++edge) {
1577 dvverb << " Inserting "<<aggregates_[edge.target()];
1578 connected_.insert(aggregates_[edge.target()]);
1579 dvverb <<" size="<<connected_.size();
1580 if(aggregates_[edge.target()]==AggregatesMap<Vertex>::UNAGGREGATED &&
1581 !graph_.getVertexProperties(edge.target()).front())
1582 {
1583 front_.push_back(edge.target());
1584 graph_.getVertexProperties(edge.target()).setFront();
1585 }
1586 }
1587 dvverb <<std::endl;
1588 std::sort(front_.begin(), front_.end());
1589 }
1590
1591 template<class G,class S>
1592 inline void Aggregate<G,S>::add(std::vector<Vertex>& vertices)
1593 {
1594#ifndef NDEBUG
1595 std::size_t oldsize = vertices_.size();
1596#endif
1597 typedef typename std::vector<Vertex>::iterator Iterator;
1598
1599 typedef typename VertexSet::iterator SIterator;
1600
1601 SIterator pos=vertices_.begin();
1602 std::vector<Vertex> newFront;
1603 newFront.reserve(front_.capacity());
1604
1605 std::set_difference(front_.begin(), front_.end(), vertices.begin(), vertices.end(),
1606 std::back_inserter(newFront));
1607 front_=newFront;
1608
1609 for(Iterator vertex=vertices.begin(); vertex != vertices.end(); ++vertex)
1610 {
1611 pos=vertices_.insert(pos,*vertex);
1612 vertices_.insert(*vertex);
1613 graph_.getVertexProperties(*vertex).resetFront(); // Not a front node any more.
1614 aggregates_[*vertex]=id_;
1615
1616 typedef typename MatrixGraph::ConstEdgeIterator iterator;
1617 const iterator end = graph_.endEdges(*vertex);
1618 for(iterator edge = graph_.beginEdges(*vertex); edge != end; ++edge) {
1619 dvverb << " Inserting "<<aggregates_[edge.target()];
1620 connected_.insert(aggregates_[edge.target()]);
1621 if(aggregates_[edge.target()]==AggregatesMap<Vertex>::UNAGGREGATED &&
1622 !graph_.getVertexProperties(edge.target()).front())
1623 {
1624 front_.push_back(edge.target());
1625 graph_.getVertexProperties(edge.target()).setFront();
1626 }
1627 dvverb <<" size="<<connected_.size();
1628 }
1629 dvverb <<std::endl;
1630 }
1631 std::sort(front_.begin(), front_.end());
1632 assert(oldsize+vertices.size()==vertices_.size());
1633 }
1634 template<class G,class S>
1635 inline void Aggregate<G,S>::clear()
1636 {
1637 vertices_.clear();
1638 connected_.clear();
1639 id_=-1;
1640 }
1641
1642 template<class G,class S>
1643 inline typename Aggregate<G,S>::VertexSet::size_type
1645 {
1646 return vertices_.size();
1647 }
1648
1649 template<class G,class S>
1650 inline typename Aggregate<G,S>::VertexSet::size_type
1652 {
1653 return connected_.size();
1654 }
1655
1656 template<class G,class S>
1657 inline int Aggregate<G,S>::id()
1658 {
1659 return id_;
1660 }
1661
1662 template<class G,class S>
1664 {
1665 return vertices_.begin();
1666 }
1667
1668 template<class G,class S>
1670 {
1671 return vertices_.end();
1672 }
1673
1674 template<class V>
1676
1677 template<class V>
1679
1680 template<class V>
1682 : aggregates_(0)
1683 {}
1684
1685 template<class V>
1687 {
1688 if(aggregates_!=0)
1689 delete[] aggregates_;
1690 }
1691
1692
1693 template<class V>
1694 inline AggregatesMap<V>::AggregatesMap(std::size_t noVertices)
1695 {
1696 allocate(noVertices);
1697 }
1698
1699 template<class V>
1700 inline std::size_t AggregatesMap<V>::noVertices() const
1701 {
1702 return noVertices_;
1703 }
1704
1705 template<class V>
1706 inline void AggregatesMap<V>::allocate(std::size_t noVertices)
1707 {
1708 aggregates_ = new AggregateDescriptor[noVertices];
1709 noVertices_ = noVertices;
1710
1711 for(std::size_t i=0; i < noVertices; i++)
1712 aggregates_[i]=UNAGGREGATED;
1713 }
1714
1715 template<class V>
1716 inline void AggregatesMap<V>::free()
1717 {
1718 assert(aggregates_ != 0);
1719 delete[] aggregates_;
1720 aggregates_=0;
1721 }
1722
1723 template<class V>
1725 AggregatesMap<V>::operator[](const VertexDescriptor& v)
1726 {
1727 return aggregates_[v];
1728 }
1729
1730 template<class V>
1731 inline const typename AggregatesMap<V>::AggregateDescriptor&
1732 AggregatesMap<V>::operator[](const VertexDescriptor& v) const
1733 {
1734 return aggregates_[v];
1735 }
1736
1737 template<class V>
1738 template<bool reset, class G, class F,class VM>
1739 inline std::size_t AggregatesMap<V>::breadthFirstSearch(const V& start,
1740 const AggregateDescriptor& aggregate,
1741 const G& graph, F& aggregateVisitor,
1742 VM& visitedMap) const
1743 {
1744 VertexList vlist;
1745
1746 DummyEdgeVisitor dummy;
1747 return breadthFirstSearch<true,reset>(start, aggregate, graph, vlist, aggregateVisitor, dummy, visitedMap);
1748 }
1749
1750 template<class V>
1751 template<bool remove, bool reset, class G, class L, class F1, class F2, class VM>
1752 std::size_t AggregatesMap<V>::breadthFirstSearch(const V& start,
1753 const AggregateDescriptor& aggregate,
1754 const G& graph,
1755 L& visited,
1756 F1& aggregateVisitor,
1757 F2& nonAggregateVisitor,
1758 VM& visitedMap) const
1759 {
1760 typedef typename L::const_iterator ListIterator;
1761 int visitedSpheres = 0;
1762
1763 visited.push_back(start);
1764 put(visitedMap, start, true);
1765
1766 ListIterator current = visited.begin();
1767 ListIterator end = visited.end();
1768 std::size_t i=0, size=visited.size();
1769
1770 // visit the neighbours of all vertices of the
1771 // current sphere.
1772 while(current != end) {
1773
1774 for(; i<size; ++current, ++i) {
1775 typedef typename G::ConstEdgeIterator EdgeIterator;
1776 const EdgeIterator endEdge = graph.endEdges(*current);
1777
1778 for(EdgeIterator edge = graph.beginEdges(*current);
1779 edge != endEdge; ++edge) {
1780
1781 if(aggregates_[edge.target()]==aggregate) {
1782 if(!get(visitedMap, edge.target())) {
1783 put(visitedMap, edge.target(), true);
1784 visited.push_back(edge.target());
1785 aggregateVisitor(edge);
1786 }
1787 }else
1788 nonAggregateVisitor(edge);
1789 }
1790 }
1791 end = visited.end();
1792 size = visited.size();
1793 if(current != end)
1794 visitedSpheres++;
1795 }
1796
1797 if(reset)
1798 for(current = visited.begin(); current != end; ++current)
1799 put(visitedMap, *current, false);
1800
1801
1802 if(remove)
1803 visited.clear();
1804
1805 return visitedSpheres;
1806 }
1807
1808 template<class G>
1810 : graph_(0), aggregate_(0), front_(), connected_(), size_(-1)
1811 {}
1812
1813 template<class G>
1815 {
1816 size_=-1;
1817 }
1818
1819 template<class G, class C>
1820 void buildDependency(G& graph,
1821 const typename C::Matrix& matrix,
1822 C criterion, bool firstlevel)
1823 {
1824 // assert(graph.isBuilt());
1825 typedef typename C::Matrix Matrix;
1826 typedef typename G::VertexIterator VertexIterator;
1827
1828 criterion.init(&matrix);
1829
1830 for(VertexIterator vertex = graph.begin(); vertex != graph.end(); ++vertex) {
1831 typedef typename Matrix::row_type Row;
1832
1833 const Row& row = matrix[*vertex];
1834
1835 // Tell the criterion what row we will examine now
1836 // This might for example be used for calculating the
1837 // maximum offdiagonal value
1838 criterion.initRow(row, *vertex);
1839
1840 // On a first path all columns are examined. After this
1841 // the calculator should know whether the vertex is isolated.
1842 typedef typename Matrix::ConstColIterator ColIterator;
1843 ColIterator end = row.end();
1844 typename FieldTraits<typename Matrix::field_type>::real_type absoffdiag=0.;
1845
1846 using std::max;
1847 if(firstlevel) {
1848 for(ColIterator col = row.begin(); col != end; ++col)
1849 if(col.index()!=*vertex) {
1850 criterion.examine(col);
1851 absoffdiag = max(absoffdiag, Impl::asMatrix(*col).frobenius_norm());
1852 }
1853
1854 if(absoffdiag==0)
1855 vertex.properties().setExcludedBorder();
1856 }
1857 else
1858 for(ColIterator col = row.begin(); col != end; ++col)
1859 if(col.index()!=*vertex)
1860 criterion.examine(col);
1861
1862 // reset the vertex properties
1863 //vertex.properties().reset();
1864
1865 // Check whether the vertex is isolated.
1866 if(criterion.isIsolated()) {
1867 //std::cout<<"ISOLATED: "<<*vertex<<std::endl;
1868 vertex.properties().setIsolated();
1869 }else{
1870 // Examine all the edges beginning at this vertex.
1871 auto eEnd = vertex.end();
1872 auto col = matrix[*vertex].begin();
1873
1874 for(auto edge = vertex.begin(); edge!= eEnd; ++edge, ++col) {
1875 // Move to the right column.
1876 while(col.index()!=edge.target())
1877 ++col;
1878 criterion.examine(graph, edge, col);
1879 }
1880 }
1881
1882 }
1883 }
1884
1885
1886 template<class G>
1887 template<class V>
1888 inline Aggregator<G>::AggregateVisitor<V>::AggregateVisitor(const AggregatesMap<Vertex>& aggregates,
1889 const AggregateDescriptor& aggregate, V& visitor)
1890 : aggregates_(aggregates), aggregate_(aggregate), visitor_(&visitor)
1891 {}
1892
1893 template<class G>
1894 template<class V>
1895 inline void Aggregator<G>::AggregateVisitor<V>::operator()(const typename MatrixGraph::ConstEdgeIterator& edge)
1896 {
1897 if(aggregates_[edge.target()]==aggregate_)
1898 visitor_->operator()(edge);
1899 }
1900
1901 template<class G>
1902 template<class V>
1903 inline void Aggregator<G>::visitAggregateNeighbours(const Vertex& vertex,
1904 const AggregateDescriptor& aggregate,
1905 const AggregatesMap<Vertex>& aggregates,
1906 V& visitor) const
1907 {
1908 // Only evaluates for edge pointing to the aggregate
1909 AggregateVisitor<V> v(aggregates, aggregate, visitor);
1910 visitNeighbours(*graph_, vertex, v);
1911 }
1912
1913
1914 template<class G>
1915 inline Aggregator<G>::Counter::Counter()
1916 : count_(0)
1917 {}
1918
1919 template<class G>
1920 inline void Aggregator<G>::Counter::increment()
1921 {
1922 ++count_;
1923 }
1924
1925 template<class G>
1926 inline void Aggregator<G>::Counter::decrement()
1927 {
1928 --count_;
1929 }
1930 template<class G>
1931 inline int Aggregator<G>::Counter::value()
1932 {
1933 return count_;
1934 }
1935
1936 template<class G>
1937 inline void Aggregator<G>::TwoWayCounter::operator()(const typename MatrixGraph::ConstEdgeIterator& edge)
1938 {
1939 if(edge.properties().isTwoWay())
1940 Counter::increment();
1941 }
1942
1943 template<class G>
1944 int Aggregator<G>::twoWayConnections(const Vertex& vertex, const AggregateDescriptor& aggregate,
1945 const AggregatesMap<Vertex>& aggregates) const
1946 {
1947 TwoWayCounter counter;
1948 visitAggregateNeighbours(vertex, aggregate, aggregates, counter);
1949 return counter.value();
1950 }
1951
1952 template<class G>
1953 int Aggregator<G>::oneWayConnections(const Vertex& vertex, const AggregateDescriptor& aggregate,
1954 const AggregatesMap<Vertex>& aggregates) const
1955 {
1956 OneWayCounter counter;
1957 visitAggregateNeighbours(vertex, aggregate, aggregates, counter);
1958 return counter.value();
1959 }
1960
1961 template<class G>
1962 inline void Aggregator<G>::OneWayCounter::operator()(const typename MatrixGraph::ConstEdgeIterator& edge)
1963 {
1964 if(edge.properties().isOneWay())
1965 Counter::increment();
1966 }
1967
1968 template<class G>
1969 inline Aggregator<G>::ConnectivityCounter::ConnectivityCounter(const VertexSet& connected,
1970 const AggregatesMap<Vertex>& aggregates)
1971 : Counter(), connected_(connected), aggregates_(aggregates)
1972 {}
1973
1974
1975 template<class G>
1976 inline void Aggregator<G>::ConnectivityCounter::operator()(const typename MatrixGraph::ConstEdgeIterator& edge)
1977 {
1978 if(connected_.find(aggregates_[edge.target()]) == connected_.end() || aggregates_[edge.target()]==AggregatesMap<Vertex>::UNAGGREGATED)
1979 // Would be a new connection
1980 Counter::increment();
1981 else{
1982 Counter::increment();
1983 Counter::increment();
1984 }
1985 }
1986
1987 template<class G>
1988 inline double Aggregator<G>::connectivity(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates) const
1989 {
1990 ConnectivityCounter counter(connected_, aggregates);
1991 double noNeighbours=visitNeighbours(*graph_, vertex, counter);
1992 return (double)counter.value()/noNeighbours;
1993 }
1994
1995 template<class G>
1996 inline Aggregator<G>::DependencyCounter::DependencyCounter()
1997 : Counter()
1998 {}
1999
2000 template<class G>
2001 inline void Aggregator<G>::DependencyCounter::operator()(const typename MatrixGraph::ConstEdgeIterator& edge)
2002 {
2003 if(edge.properties().depends())
2004 Counter::increment();
2005 if(edge.properties().influences())
2006 Counter::increment();
2007 }
2008
2009 template<class G>
2010 int Aggregator<G>::unusedNeighbours(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates) const
2011 {
2012 return aggregateNeighbours(vertex, AggregatesMap<Vertex>::UNAGGREGATED, aggregates);
2013 }
2014
2015 template<class G>
2016 std::pair<int,int> Aggregator<G>::neighbours(const Vertex& vertex,
2017 const AggregateDescriptor& aggregate,
2018 const AggregatesMap<Vertex>& aggregates) const
2019 {
2020 DependencyCounter unused, aggregated;
2021 typedef AggregateVisitor<DependencyCounter> CounterT;
2022 typedef std::tuple<CounterT,CounterT> CounterTuple;
2023 CombinedFunctor<CounterTuple> visitors(CounterTuple(CounterT(aggregates, AggregatesMap<Vertex>::UNAGGREGATED, unused), CounterT(aggregates, aggregate, aggregated)));
2024 visitNeighbours(*graph_, vertex, visitors);
2025 return std::make_pair(unused.value(), aggregated.value());
2026 }
2027
2028
2029 template<class G>
2030 int Aggregator<G>::aggregateNeighbours(const Vertex& vertex, const AggregateDescriptor& aggregate, const AggregatesMap<Vertex>& aggregates) const
2031 {
2032 DependencyCounter counter;
2033 visitAggregateNeighbours(vertex, aggregate, aggregates, counter);
2034 return counter.value();
2035 }
2036
2037 template<class G>
2038 std::size_t Aggregator<G>::distance(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates)
2039 {
2040 return 0;
2041 typename PropertyMapTypeSelector<VertexVisitedTag,G>::Type visitedMap = get(VertexVisitedTag(), *graph_);
2042 VertexList vlist;
2043 typename AggregatesMap<Vertex>::DummyEdgeVisitor dummy;
2044 return aggregates.template breadthFirstSearch<true,true>(vertex,
2045 aggregate_->id(), *graph_,
2046 vlist, dummy, dummy, visitedMap);
2047 }
2048
2049 template<class G>
2050 inline Aggregator<G>::FrontMarker::FrontMarker(std::vector<Vertex>& front, MatrixGraph& graph)
2051 : front_(front), graph_(graph)
2052 {}
2053
2054 template<class G>
2055 inline void Aggregator<G>::FrontMarker::operator()(const typename MatrixGraph::ConstEdgeIterator& edge)
2056 {
2057 Vertex target = edge.target();
2058
2059 if(!graph_.getVertexProperties(target).front()) {
2060 front_.push_back(target);
2061 graph_.getVertexProperties(target).setFront();
2062 }
2063 }
2064
2065 template<class G>
2066 inline bool Aggregator<G>::admissible(const Vertex& vertex, const AggregateDescriptor& aggregate, const AggregatesMap<Vertex>& aggregates) const
2067 {
2068 // Todo
2069 Dune::dvverb<<" Admissible not yet implemented!"<<std::endl;
2070 return true;
2071 //Situation 1: front node depends on two nodes. Then these
2072 // have to be strongly connected to each other
2073
2074 // Iterate over all all neighbours of front node
2075 typedef typename MatrixGraph::ConstEdgeIterator Iterator;
2076 Iterator vend = graph_->endEdges(vertex);
2077 for(Iterator edge = graph_->beginEdges(vertex); edge != vend; ++edge) {
2078 // if(edge.properties().depends() && !edge.properties().influences()
2079 if(edge.properties().isStrong()
2080 && aggregates[edge.target()]==aggregate)
2081 {
2082 // Search for another link to the aggregate
2083 Iterator edge1 = edge;
2084 for(++edge1; edge1 != vend; ++edge1) {
2085 //if(edge1.properties().depends() && !edge1.properties().influences()
2086 if(edge1.properties().isStrong()
2087 && aggregates[edge.target()]==aggregate)
2088 {
2089 //Search for an edge connecting the two vertices that is
2090 //strong
2091 bool found=false;
2092 Iterator v2end = graph_->endEdges(edge.target());
2093 for(Iterator edge2 = graph_->beginEdges(edge.target()); edge2 != v2end; ++edge2) {
2094 if(edge2.target()==edge1.target() &&
2095 edge2.properties().isStrong()) {
2096 found =true;
2097 break;
2098 }
2099 }
2100 if(found)
2101 {
2102 return true;
2103 }
2104 }
2105 }
2106 }
2107 }
2108
2109 // Situation 2: cluster node depends on front node and other cluster node
2111 vend = graph_->endEdges(vertex);
2112 for(Iterator edge = graph_->beginEdges(vertex); edge != vend; ++edge) {
2113 //if(!edge.properties().depends() && edge.properties().influences()
2114 if(edge.properties().isStrong()
2115 && aggregates[edge.target()]==aggregate)
2116 {
2117 // Search for a link from target that stays within the aggregate
2118 Iterator v1end = graph_->endEdges(edge.target());
2119
2120 for(Iterator edge1=graph_->beginEdges(edge.target()); edge1 != v1end; ++edge1) {
2121 //if(edge1.properties().depends() && !edge1.properties().influences()
2122 if(edge1.properties().isStrong()
2123 && aggregates[edge1.target()]==aggregate)
2124 {
2125 bool found=false;
2126 // Check if front node is also connected to this one
2127 Iterator v2end = graph_->endEdges(vertex);
2128 for(Iterator edge2 = graph_->beginEdges(vertex); edge2 != v2end; ++edge2) {
2129 if(edge2.target()==edge1.target()) {
2130 if(edge2.properties().isStrong())
2131 found=true;
2132 break;
2133 }
2134 }
2135 if(found)
2136 {
2137 return true;
2138 }
2139 }
2140 }
2141 }
2142 }
2143 return false;
2144 }
2145
2146 template<class G>
2147 void Aggregator<G>::unmarkFront()
2148 {
2149 typedef typename std::vector<Vertex>::const_iterator Iterator;
2150
2151 for(Iterator vertex=front_.begin(); vertex != front_.end(); ++vertex)
2152 graph_->getVertexProperties(*vertex).resetFront();
2153
2154 front_.clear();
2155 }
2156
2157 template<class G>
2158 inline void
2159 Aggregator<G>::nonisoNeighbourAggregate(const Vertex& vertex,
2160 const AggregatesMap<Vertex>& aggregates,
2161 SLList<Vertex>& neighbours) const
2162 {
2163 typedef typename MatrixGraph::ConstEdgeIterator Iterator;
2164 Iterator end=graph_->beginEdges(vertex);
2165 neighbours.clear();
2166
2167 for(Iterator edge=graph_->beginEdges(vertex); edge!=end; ++edge)
2168 {
2169 if(aggregates[edge.target()]!=AggregatesMap<Vertex>::UNAGGREGATED && graph_->getVertexProperties(edge.target()).isolated())
2170 neighbours.push_back(aggregates[edge.target()]);
2171 }
2172 }
2173
2174 template<class G>
2175 inline typename G::VertexDescriptor Aggregator<G>::mergeNeighbour(const Vertex& vertex, const AggregatesMap<Vertex>& aggregates) const
2176 {
2177 typedef typename MatrixGraph::ConstEdgeIterator Iterator;
2178
2179 Iterator end = graph_->endEdges(vertex);
2180 for(Iterator edge = graph_->beginEdges(vertex); edge != end; ++edge) {
2181 if(aggregates[edge.target()] != AggregatesMap<Vertex>::UNAGGREGATED &&
2182 graph_->getVertexProperties(edge.target()).isolated() == graph_->getVertexProperties(edge.source()).isolated()) {
2183 if( graph_->getVertexProperties(vertex).isolated() ||
2184 ((edge.properties().depends() || edge.properties().influences())
2185 && admissible(vertex, aggregates[edge.target()], aggregates)))
2186 return edge.target();
2187 }
2188 }
2190 }
2191
2192 template<class G>
2193 Aggregator<G>::FrontNeighbourCounter::FrontNeighbourCounter(const MatrixGraph& graph)
2194 : Counter(), graph_(graph)
2195 {}
2196
2197 template<class G>
2198 void Aggregator<G>::FrontNeighbourCounter::operator()(const typename MatrixGraph::ConstEdgeIterator& edge)
2199 {
2200 if(graph_.getVertexProperties(edge.target()).front())
2201 Counter::increment();
2202 }
2203
2204 template<class G>
2205 int Aggregator<G>::noFrontNeighbours(const Vertex& vertex) const
2206 {
2207 FrontNeighbourCounter counter(*graph_);
2208 visitNeighbours(*graph_, vertex, counter);
2209 return counter.value();
2210 }
2211 template<class G>
2212 inline bool Aggregator<G>::connected(const Vertex& vertex,
2213 const AggregateDescriptor& aggregate,
2214 const AggregatesMap<Vertex>& aggregates) const
2215 {
2216 typedef typename G::ConstEdgeIterator iterator;
2217 const iterator end = graph_->endEdges(vertex);
2218 for(iterator edge = graph_->beginEdges(vertex); edge != end; ++edge)
2219 if(aggregates[edge.target()]==aggregate)
2220 return true;
2221 return false;
2222 }
2223 template<class G>
2224 inline bool Aggregator<G>::connected(const Vertex& vertex,
2225 const SLList<AggregateDescriptor>& aggregateList,
2226 const AggregatesMap<Vertex>& aggregates) const
2227 {
2228 typedef typename SLList<AggregateDescriptor>::const_iterator Iter;
2229 for(Iter i=aggregateList.begin(); i!=aggregateList.end(); ++i)
2230 if(connected(vertex, *i, aggregates))
2231 return true;
2232 return false;
2233 }
2234
2235 template<class G>
2236 template<class C>
2237 void Aggregator<G>::growIsolatedAggregate(const Vertex& seed, const AggregatesMap<Vertex>& aggregates, const C& c)
2238 {
2239 SLList<Vertex> connectedAggregates;
2240 nonisoNeighbourAggregate(seed, aggregates,connectedAggregates);
2241
2242 while(aggregate_->size()< c.minAggregateSize() && aggregate_->connectSize() < c.maxConnectivity()) {
2243 double maxCon=-1;
2244 std::size_t maxFrontNeighbours=0;
2245
2247
2248 typedef typename std::vector<Vertex>::const_iterator Iterator;
2249
2250 for(Iterator vertex = front_.begin(); vertex != front_.end(); ++vertex) {
2251 if(distance(*vertex, aggregates)>c.maxDistance())
2252 continue; // distance of proposes aggregate too big
2253
2254 if(connectedAggregates.size()>0) {
2255 // there is already a neighbour cluster
2256 // front node must be connected to same neighbour cluster
2257
2258 if(!connected(*vertex, connectedAggregates, aggregates))
2259 continue;
2260 }
2261
2262 double con = connectivity(*vertex, aggregates);
2263
2264 if(con == maxCon) {
2265 std::size_t frontNeighbours = noFrontNeighbours(*vertex);
2266
2267 if(frontNeighbours >= maxFrontNeighbours) {
2268 maxFrontNeighbours = frontNeighbours;
2269 candidate = *vertex;
2270 }
2271 }else if(con > maxCon) {
2272 maxCon = con;
2273 maxFrontNeighbours = noFrontNeighbours(*vertex);
2274 candidate = *vertex;
2275 }
2276 }
2277
2279 break;
2280
2281 aggregate_->add(candidate);
2282 }
2283 }
2284
2285 template<class G>
2286 template<class C>
2287 void Aggregator<G>::growAggregate(const Vertex& seed, const AggregatesMap<Vertex>& aggregates, const C& c)
2288 {
2289 using std::min;
2290
2291 std::size_t distance_ =0;
2292 while(aggregate_->size() < c.minAggregateSize()&& distance_<c.maxDistance()) {
2293 int maxTwoCons=0, maxOneCons=0, maxNeighbours=-1;
2294 double maxCon=-1;
2295
2296 std::vector<Vertex> candidates;
2297 candidates.reserve(30);
2298
2299 typedef typename std::vector<Vertex>::const_iterator Iterator;
2300
2301 for(Iterator vertex = front_.begin(); vertex != front_.end(); ++vertex) {
2302 // Only nonisolated nodes are considered
2303 if(graph_->getVertexProperties(*vertex).isolated())
2304 continue;
2305
2306 int twoWayCons = twoWayConnections(*vertex, aggregate_->id(), aggregates);
2307
2308 /* The case of two way connections. */
2309 if( maxTwoCons == twoWayCons && twoWayCons > 0) {
2310 double con = connectivity(*vertex, aggregates);
2311
2312 if(con == maxCon) {
2313 int neighbours = noFrontNeighbours(*vertex);
2314
2315 if(neighbours > maxNeighbours) {
2316 maxNeighbours = neighbours;
2317 candidates.clear();
2318 candidates.push_back(*vertex);
2319 }else{
2320 candidates.push_back(*vertex);
2321 }
2322 }else if( con > maxCon) {
2323 maxCon = con;
2324 maxNeighbours = noFrontNeighbours(*vertex);
2325 candidates.clear();
2326 candidates.push_back(*vertex);
2327 }
2328 }else if(twoWayCons > maxTwoCons) {
2329 maxTwoCons = twoWayCons;
2330 maxCon = connectivity(*vertex, aggregates);
2331 maxNeighbours = noFrontNeighbours(*vertex);
2332 candidates.clear();
2333 candidates.push_back(*vertex);
2334
2335 // two way connections precede
2336 maxOneCons = std::numeric_limits<int>::max();
2337 }
2338
2339 if(twoWayCons > 0)
2340 {
2341 continue; // THis is a two-way node, skip tests for one way nodes
2342 }
2343
2344 /* The one way case */
2345 int oneWayCons = oneWayConnections(*vertex, aggregate_->id(), aggregates);
2346
2347 if(oneWayCons==0)
2348 continue; // No strong connections, skip the tests.
2349
2350 if(!admissible(*vertex, aggregate_->id(), aggregates))
2351 continue;
2352
2353 if( maxOneCons == oneWayCons && oneWayCons > 0) {
2354 double con = connectivity(*vertex, aggregates);
2355
2356 if(con == maxCon) {
2357 int neighbours = noFrontNeighbours(*vertex);
2358
2359 if(neighbours > maxNeighbours) {
2360 maxNeighbours = neighbours;
2361 candidates.clear();
2362 candidates.push_back(*vertex);
2363 }else{
2364 if(neighbours==maxNeighbours)
2365 {
2366 candidates.push_back(*vertex);
2367 }
2368 }
2369 }else if( con > maxCon) {
2370 maxCon = con;
2371 maxNeighbours = noFrontNeighbours(*vertex);
2372 candidates.clear();
2373 candidates.push_back(*vertex);
2374 }
2375 }else if(oneWayCons > maxOneCons) {
2376 maxOneCons = oneWayCons;
2377 maxCon = connectivity(*vertex, aggregates);
2378 maxNeighbours = noFrontNeighbours(*vertex);
2379 candidates.clear();
2380 candidates.push_back(*vertex);
2381 }
2382 }
2383
2384
2385 if(!candidates.size())
2386 break; // No more candidates found
2387 distance_=distance(seed, aggregates);
2388 candidates.resize(min(candidates.size(), c.maxAggregateSize()-
2389 aggregate_->size()));
2390 aggregate_->add(candidates);
2391 }
2392 }
2393
2394 template<typename V>
2395 template<typename M, typename G, typename C>
2396 std::tuple<int,int,int,int> AggregatesMap<V>::buildAggregates(const M& matrix, G& graph, const C& criterion,
2397 bool finestLevel)
2398 {
2399 Aggregator<G> aggregator;
2400 return aggregator.build(matrix, graph, *this, criterion, finestLevel);
2401 }
2402
2403 template<class G>
2404 template<class M, class C>
2405 std::tuple<int,int,int,int> Aggregator<G>::build(const M& m, G& graph, AggregatesMap<Vertex>& aggregates, const C& c,
2406 bool finestLevel)
2407 {
2408 using std::max;
2409 using std::min;
2410 // Stack for fast vertex access
2411 Stack stack_(graph, *this, aggregates);
2412
2413 graph_ = &graph;
2414
2415 aggregate_ = new Aggregate<G,VertexSet>(graph, aggregates, connected_, front_);
2416
2417 Timer watch;
2418 watch.reset();
2419
2420 buildDependency(graph, m, c, finestLevel);
2421
2422 dverb<<"Build dependency took "<< watch.elapsed()<<" seconds."<<std::endl;
2423 int noAggregates, conAggregates, isoAggregates, oneAggregates;
2424 std::size_t maxA=0, minA=1000000, avg=0;
2425 int skippedAggregates;
2426 noAggregates = conAggregates = isoAggregates = oneAggregates =
2427 skippedAggregates = 0;
2428
2429 while(true) {
2430 Vertex seed = stack_.pop();
2431
2432 if(seed == Stack::NullEntry)
2433 // No more unaggregated vertices. We are finished!
2434 break;
2435
2436 // Debugging output
2437 if((noAggregates+1)%10000 == 0)
2438 Dune::dverb<<"c";
2439 unmarkFront();
2440
2441 if(graph.getVertexProperties(seed).excludedBorder()) {
2442 aggregates[seed]=AggregatesMap<Vertex>::ISOLATED;
2443 ++skippedAggregates;
2444 continue;
2445 }
2446
2447 if(graph.getVertexProperties(seed).isolated()) {
2448 if(c.skipIsolated()) {
2449 // isolated vertices are not aggregated but skipped on the coarser levels.
2450 aggregates[seed]=AggregatesMap<Vertex>::ISOLATED;
2451 ++skippedAggregates;
2452 // skip rest as no agglomeration is done.
2453 continue;
2454 }else{
2455 aggregate_->seed(seed);
2456 growIsolatedAggregate(seed, aggregates, c);
2457 }
2458 }else{
2459 aggregate_->seed(seed);
2460 growAggregate(seed, aggregates, c);
2461 }
2462
2463 /* The rounding step. */
2464 while(!(graph.getVertexProperties(seed).isolated()) && aggregate_->size() < c.maxAggregateSize()) {
2465
2466 std::vector<Vertex> candidates;
2467 candidates.reserve(30);
2468
2469 typedef typename std::vector<Vertex>::const_iterator Iterator;
2470
2471 for(Iterator vertex = front_.begin(); vertex != front_.end(); ++vertex) {
2472
2473 if(graph.getVertexProperties(*vertex).isolated())
2474 continue; // No isolated nodes here
2475
2476 if(twoWayConnections( *vertex, aggregate_->id(), aggregates) == 0 &&
2477 (oneWayConnections( *vertex, aggregate_->id(), aggregates) == 0 ||
2478 !admissible( *vertex, aggregate_->id(), aggregates) ))
2479 continue;
2480
2481 std::pair<int,int> neighbourPair=neighbours(*vertex, aggregate_->id(),
2482 aggregates);
2483
2484 //if(aggregateNeighbours(*vertex, aggregate_->id(), aggregates) <= unusedNeighbours(*vertex, aggregates))
2485 // continue;
2486
2487 if(neighbourPair.first >= neighbourPair.second)
2488 continue;
2489
2490 if(distance(*vertex, aggregates) > c.maxDistance())
2491 continue; // Distance too far
2492 candidates.push_back(*vertex);
2493 break;
2494 }
2495
2496 if(!candidates.size()) break; // no more candidates found.
2497
2498 candidates.resize(min(candidates.size(), c.maxAggregateSize()-
2499 aggregate_->size()));
2500 aggregate_->add(candidates);
2501
2502 }
2503
2504 // try to merge aggregates consisting of only one nonisolated vertex with other aggregates
2505 if(aggregate_->size()==1 && c.maxAggregateSize()>1) {
2506 if(!graph.getVertexProperties(seed).isolated()) {
2507 Vertex mergedNeighbour = mergeNeighbour(seed, aggregates);
2508
2509 if(mergedNeighbour != AggregatesMap<Vertex>::UNAGGREGATED) {
2510 // assign vertex to the neighbouring cluster
2511 aggregates[seed] = aggregates[mergedNeighbour];
2512 aggregate_->invalidate();
2513 }else{
2514 ++avg;
2515 minA=min(minA,static_cast<std::size_t>(1));
2516 maxA=max(maxA,static_cast<std::size_t>(1));
2517 ++oneAggregates;
2518 ++conAggregates;
2519 }
2520 }else{
2521 ++avg;
2522 minA=min(minA,static_cast<std::size_t>(1));
2523 maxA=max(maxA,static_cast<std::size_t>(1));
2524 ++oneAggregates;
2525 ++isoAggregates;
2526 }
2527 ++avg;
2528 }else{
2529 avg+=aggregate_->size();
2530 minA=min(minA,aggregate_->size());
2531 maxA=max(maxA,aggregate_->size());
2532 if(graph.getVertexProperties(seed).isolated())
2533 ++isoAggregates;
2534 else
2535 ++conAggregates;
2536 }
2537
2538 }
2539
2540 Dune::dinfo<<"connected aggregates: "<<conAggregates;
2541 Dune::dinfo<<" isolated aggregates: "<<isoAggregates;
2542 if(conAggregates+isoAggregates>0)
2543 Dune::dinfo<<" one node aggregates: "<<oneAggregates<<" min size="
2544 <<minA<<" max size="<<maxA
2545 <<" avg="<<avg/(conAggregates+isoAggregates)<<std::endl;
2546
2547 delete aggregate_;
2548 return std::make_tuple(conAggregates+isoAggregates,isoAggregates,
2549 oneAggregates,skippedAggregates);
2550 }
2551
2552
2553 template<class G>
2554 Aggregator<G>::Stack::Stack(const MatrixGraph& graph, const Aggregator<G>& aggregatesBuilder,
2555 const AggregatesMap<Vertex>& aggregates)
2556 : graph_(graph), aggregatesBuilder_(aggregatesBuilder), aggregates_(aggregates), begin_(graph.begin()), end_(graph.end())
2557 {
2558 //vals_ = new Vertex[N];
2559 }
2560
2561 template<class G>
2562 Aggregator<G>::Stack::~Stack()
2563 {
2564 //Dune::dverb << "Max stack size was "<<maxSize_<<" filled="<<filled_<<std::endl;
2565 //delete[] vals_;
2566 }
2567
2568 template<class G>
2569 const typename Aggregator<G>::Vertex Aggregator<G>::Stack::NullEntry
2571
2572 template<class G>
2573 inline typename G::VertexDescriptor Aggregator<G>::Stack::pop()
2574 {
2575 for(; begin_!=end_ && aggregates_[*begin_] != AggregatesMap<Vertex>::UNAGGREGATED; ++begin_) ;
2576
2577 if(begin_!=end_)
2578 {
2579 typename G::VertexDescriptor current=*begin_;
2580 ++begin_;
2581 return current;
2582 }else
2583 return NullEntry;
2584 }
2585
2586#endif // DOXYGEN
2587
2588 template<class V>
2589 void printAggregates2d(const AggregatesMap<V>& aggregates, int n, int m, std::ostream& os)
2590 {
2591 using std::max;
2592
2593 std::ios_base::fmtflags oldOpts=os.flags();
2594
2595 os.setf(std::ios_base::right, std::ios_base::adjustfield);
2596
2597 V maxVal=0;
2598 int width=1;
2599
2600 for(int i=0; i< n*m; i++)
2601 maxVal=max(maxVal, aggregates[i]);
2602
2603 for(int i=10; i < 1000000; i*=10)
2604 if(maxVal/i>0)
2605 width++;
2606 else
2607 break;
2608
2609 for(int j=0, entry=0; j < m; j++) {
2610 for(int i=0; i<n; i++, entry++) {
2611 os.width(width);
2612 os<<aggregates[entry]<<" ";
2613 }
2614
2615 os<<std::endl;
2616 }
2617 os<<std::endl;
2618 os.flags(oldOpts);
2619 }
2620
2621
2622 } // namespace Amg
2623
2624} // namespace Dune
2625
2626
2627#endif
A class for temporarily storing the vertices of an aggregate in.
Definition: aggregates.hh:784
A Dummy visitor that does nothing for each visited edge.
Definition: aggregates.hh:604
Class providing information about the mapping of the vertices onto aggregates.
Definition: aggregates.hh:566
Base class of all aggregation criterions.
Definition: aggregates.hh:51
Class for building the aggregates.
Definition: aggregates.hh:915
Dependency policy for symmetric matrices.
Definition: aggregates.hh:255
Norm that uses only the [N][N] entry of the block to determine couplings.
Definition: aggregates.hh:381
Norm that uses only the [0][0] entry of the block to determine couplings.
Definition: aggregates.hh:457
Iterator over all edges starting from a vertex.
Definition: graph.hh:95
The vertex iterator type of the graph.
Definition: graph.hh:209
The (undirected) graph of a matrix.
Definition: graph.hh:51
M::size_type VertexDescriptor
The vertex descriptor.
Definition: graph.hh:73
EdgeIteratorT< const MatrixGraph< Matrix > > ConstEdgeIterator
The constant edge iterator type.
Definition: graph.hh:298
All parameters for AMG.
Definition: parameters.hh:416
Criterion taking advantage of symmetric matrices.
Definition: aggregates.hh:525
Dependency policy for symmetric matrices.
Definition: aggregates.hh:316
Dependency policy for symmetric matrices.
Definition: aggregates.hh:136
Criterion suitable for unsymmetric matrices.
Definition: aggregates.hh:545
derive error class from the base class in common
Definition: istlexception.hh:19
A generic dynamic dense matrix.
Definition: matrix.hh:561
typename Imp::BlockTraits< T >::field_type field_type
Export the type representing the underlying field.
Definition: matrix.hh:565
row_type::const_iterator ConstColIterator
Const iterator for the entries of each row.
Definition: matrix.hh:589
MatrixImp::DenseMatrixBase< T, A >::window_type row_type
The type implementing a matrix row.
Definition: matrix.hh:574
An allocator managing a pool of objects for reuse.
Definition: poolallocator.hh:223
A single linked list.
Definition: sllist.hh:44
Type traits to determine the type of reals (when working with complex numbers)
Provides classes for building the matrix graph.
SLListConstIterator< T, A > const_iterator
The constant iterator of the list.
Definition: sllist.hh:74
#define DUNE_THROW(E,...)
Definition: exceptions.hh:312
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:492
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:485
constexpr auto min
Function object that returns the smaller of the given values.
Definition: hybridutilities.hh:507
Matrix::ConstColIterator ColIter
Constant column iterator of the matrix.
Definition: aggregates.hh:275
Matrix::ConstColIterator ColIter
Constant column iterator of the matrix.
Definition: aggregates.hh:156
std::size_t breadthFirstSearch(const VertexDescriptor &start, const AggregateDescriptor &aggregate, const G &graph, L &visited, F1 &aggregateVisitor, F2 &nonAggregateVisitor, VM &visitedMap) const
Breadth first search within an aggregate.
PoolAllocator< VertexDescriptor, 100 > Allocator
The allocator we use for our lists and the set.
Definition: aggregates.hh:592
int id()
Get the id identifying the aggregate.
Norm norm_
The functor for calculating the norm.
Definition: aggregates.hh:304
MatrixGraph::VertexDescriptor Vertex
The vertex identifier.
Definition: aggregates.hh:926
AggregationCriterion()
Constructor.
Definition: aggregates.hh:68
const Matrix * matrix_
The matrix we work on.
Definition: aggregates.hh:359
auto operator()(const M &m, typename std::enable_if_t< Dune::IsNumber< M >::value > *sfinae=nullptr) const
Compute the norm of a scalar.
Definition: aggregates.hh:408
M Matrix
The matrix type we build the dependency of.
Definition: aggregates.hh:260
G MatrixGraph
The matrix graph type used.
Definition: aggregates.hh:921
Norm norm_
The functor for calculating the norm.
Definition: aggregates.hh:365
M Matrix
The matrix type we build the dependency of.
Definition: aggregates.hh:321
real_type diagonal_
The norm of the current diagonal.
Definition: aggregates.hh:189
N Norm
The norm to use for examining the matrix entries.
Definition: aggregates.hh:265
int row_
index of the currently evaluated row.
Definition: aggregates.hh:187
std::tuple< int, int, int, int > build(const M &m, G &graph, AggregatesMap< Vertex > &aggregates, const C &c, bool finestLevel)
Build the aggregates.
FrontNeighbourCounter(const MatrixGraph &front)
Constructor.
Matrix::row_type Row
Constant Row iterator of the matrix.
Definition: aggregates.hh:331
const Matrix * matrix_
The matrix we work on.
Definition: aggregates.hh:298
const AggregateDescriptor & operator[](const VertexDescriptor &v) const
Get the aggregate a vertex belongs to.
AggregateVisitor(const AggregatesMap< Vertex > &aggregates, const AggregateDescriptor &aggregate, Visitor &visitor)
Constructor.
Matrix::ConstColIterator ColIter
Constant column iterator of the matrix.
Definition: aggregates.hh:336
~AggregatesMap()
Destructor.
Matrix::field_type field_type
The current max value.
Definition: aggregates.hh:181
void decrement()
Decrement counter.
Aggregate(MatrixGraph &graph, AggregatesMap< Vertex > &aggregates, VertexSet &connectivity, std::vector< Vertex > &front_)
Constructor.
V Visitor
The type of the adapted visitor.
Definition: aggregates.hh:1070
std::size_t * SphereMap
Type of the mapping of aggregate members onto distance spheres.
Definition: aggregates.hh:815
Matrix::row_type Row
Constant Row iterator of the matrix.
Definition: aggregates.hh:270
VertexSet::size_type connectSize()
Get the number of connections to other aggregates.
N Norm
The norm to use for examining the matrix entries.
Definition: aggregates.hh:326
Norm norm_
The functor for calculating the norm.
Definition: aggregates.hh:185
VertexSet::const_iterator const_iterator
Const iterator over a vertex list.
Definition: aggregates.hh:810
MatrixGraph::VertexDescriptor AggregateDescriptor
The type of the aggregate descriptor.
Definition: aggregates.hh:929
real_type diagonal_
The norm of the current diagonal.
Definition: aggregates.hh:369
~Aggregator()
Destructor.
void add(const Vertex &vertex)
Add a vertex to the aggregate.
T DependencyPolicy
The policy for calculating the dependency graph.
Definition: aggregates.hh:57
auto operator()(const M &m) const
compute the norm of a matrix.
Definition: aggregates.hh:475
Aggregator()
Constructor.
void operator()(const typename MatrixGraph::ConstEdgeIterator &edge)
Examine an edge.
FrontMarker(std::vector< Vertex > &front, MatrixGraph &graph)
Constructor.
FieldTraits< typenameM::field_type >::real_type operator()(const M &m) const
compute the norm of a matrix.
Definition: aggregates.hh:512
int visitNeighbours(const G &graph, const typename G::VertexDescriptor &vertex, V &visitor)
Visit all neighbour vertices of a vertex in a graph.
void setDefaultValuesIsotropic(std::size_t dim, std::size_t diameter=2)
Sets reasonable default values for an isotropic problem.
Definition: aggregates.hh:84
V AggregateDescriptor
The aggregate descriptor type.
Definition: aggregates.hh:586
static const V ISOLATED
Identifier of isolated vertices.
Definition: aggregates.hh:577
int row_
index of the currently evaluated row.
Definition: aggregates.hh:367
real_type diagonal_
The norm of the current diagonal.
Definition: aggregates.hh:308
std::size_t noVertices() const
Get the number of vertices.
void setDefaultValuesAnisotropic(std::size_t dim, std::size_t diameter=2)
Sets reasonable default values for an aisotropic problem.
Definition: aggregates.hh:107
AggregatesMap(std::size_t noVertices)
Constructs with allocating memory.
Matrix::field_type field_type
The current max value.
Definition: aggregates.hh:361
AggregateDescriptor & operator[](const VertexDescriptor &v)
Get the aggregate a vertex belongs to.
AggregatesMap()
Constructs without allocating memory.
int value()
Access the current count.
SLList< VertexDescriptor, Allocator > VertexList
The type of a single linked list of vertex descriptors.
Definition: aggregates.hh:598
ConnectivityCounter(const VertexSet &connected, const AggregatesMap< Vertex > &aggregates)
Constructor.
VertexSet::size_type size()
Get the size of the aggregate.
const_iterator end() const
get an iterator over the vertices of the aggregate.
int row_
index of the currently evaluated row.
Definition: aggregates.hh:306
M Matrix
The matrix type we build the dependency of.
Definition: aggregates.hh:141
const Matrix * matrix_
The matrix we work on.
Definition: aggregates.hh:179
S VertexSet
The type of a single linked list of vertex descriptors.
Definition: aggregates.hh:807
FieldTraits< typenameM::field_type >::real_type operator()(const M &m) const
compute the norm of a matrix.
Definition: aggregates.hh:496
static const V UNAGGREGATED
Identifier of not yet aggregated vertices.
Definition: aggregates.hh:572
std::size_t breadthFirstSearch(const VertexDescriptor &start, const AggregateDescriptor &aggregate, const G &graph, F &aggregateVisitor, VM &visitedMap) const
Breadth first search within an aggregate.
Matrix::field_type field_type
The current max value.
Definition: aggregates.hh:300
void allocate(std::size_t noVertices)
Allocate memory for holding the information.
N Norm
The norm to use for examining the matrix entries.
Definition: aggregates.hh:146
void reconstruct(const Vertex &vertex)
Reconstruct the aggregat from an seed node.
const_iterator begin() const
get an iterator over the vertices of the aggregate.
FieldTraits< typenameM::field_type >::real_type operator()(const M &m, typename std::enable_if_t<!Dune::IsNumber< M >::value > *sfinae=nullptr) const
compute the norm of a matrix.
Definition: aggregates.hh:392
MatrixGraph::VertexDescriptor Vertex
The vertex descriptor type.
Definition: aggregates.hh:795
void seed(const Vertex &vertex)
Initialize the aggregate with one vertex.
void clear()
Clear the aggregate.
void free()
Free the allocated memory.
void increment()
Increment counter.
void buildDependency(G &graph, const typename C::Matrix &matrix, C criterion, bool finestLevel)
Build the dependency of the matrix graph.
V VertexDescriptor
The vertex descriptor type.
Definition: aggregates.hh:581
std::tuple< int, int, int, int > buildAggregates(const M &matrix, G &graph, const C &criterion, bool finestLevel)
Build the aggregates.
Matrix::row_type Row
Constant Row iterator of the matrix.
Definition: aggregates.hh:151
PoolAllocator< Vertex, 100 > Allocator
The allocator we use for our lists and the set.
Definition: aggregates.hh:801
DVVerbType dvverb(std::cout)
stream for very verbose output.
Definition: stdstreams.hh:96
DInfoType dinfo(std::cout)
Stream for informative output.
Definition: stdstreams.hh:141
DVerbType dverb(std::cout)
Singleton of verbose debug stream.
Definition: stdstreams.hh:117
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
constexpr std::integral_constant< T, I0 > front(std::integer_sequence< T, I0, II... >)
Return the first entry of the sequence.
Definition: integersequence.hh:39
constexpr auto get(std::integer_sequence< T, II... >, std::integral_constant< std::size_t, pos >={})
Return the entry at position pos of the given sequence.
Definition: integersequence.hh:22
Parameter classes for customizing AMG.
An stl-compliant pool allocator.
Provides classes for handling internal properties in a graph.
Implements a scalar matrix view wrapper around an existing scalar.
Implements a singly linked list together with the necessary iterators.
Standard Dune debug streams.
Functor using the row sum (infinity) norm to determine strong couplings.
Definition: aggregates.hh:465
Whether this type acts as a scalar in the context of (hierarchically blocked) containers.
Definition: typetraits.hh:194
A simple timing class.
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Mar 12, 23:28, 2025)