Dune Core Modules (2.3.1)

hierarchy.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// $Id$
4#ifndef DUNE_AMGHIERARCHY_HH
5#define DUNE_AMGHIERARCHY_HH
6
7#include <list>
8#include <memory>
9#include <limits>
10#include <algorithm>
11#include "aggregates.hh"
12#include "graph.hh"
13#include "galerkin.hh"
14#include "renumberer.hh"
15#include "graphcreator.hh"
17#include <dune/common/unused.hh>
18#include <dune/common/timer.hh>
19#include <dune/common/tuples.hh>
21#include <dune/istl/bvector.hh>
32
33namespace Dune
34{
35 namespace Amg
36 {
48 enum {
56 MAX_PROCESSES = 72000
57 };
58
66 template<typename T, typename A=std::allocator<T> >
68 {
69 public:
73 typedef T MemberType;
74
75 template<typename T1, typename T2>
76 class LevelIterator;
77
78 private:
82 struct Element
83 {
84 friend class LevelIterator<Hierarchy<T,A>, T>;
85 friend class LevelIterator<const Hierarchy<T,A>, const T>;
86
88 Element* coarser_;
89
91 Element* finer_;
92
94 MemberType* element_;
95
97 MemberType* redistributed_;
98 };
99 public:
100 // enum{
101 // /**
102 // * @brief If true only the method addCoarser will be usable
103 // * otherwise only the method addFiner will be usable.
104 // */
105 // coarsen = b
106 // };
107
111 typedef typename A::template rebind<Element>::other Allocator;
112
113 typedef typename ConstructionTraits<T>::Arguments Arguments;
114
120
128
133
137 Hierarchy(const Hierarchy& other);
142 void addCoarser(Arguments& args);
143
144 void addRedistributedOnCoarsest(Arguments& args);
145
150 void addFiner(Arguments& args);
151
158 template<class C, class T1>
160 : public BidirectionalIteratorFacade<LevelIterator<C,T1>,T1,T1&>
161 {
162 friend class LevelIterator<typename remove_const<C>::type,
163 typename remove_const<T1>::type >;
164 friend class LevelIterator<const typename remove_const<C>::type,
165 const typename remove_const<T1>::type >;
166
167 public:
170 : element_(0)
171 {}
172
173 LevelIterator(Element* element)
174 : element_(element)
175 {}
176
178 LevelIterator(const LevelIterator<typename remove_const<C>::type,
179 typename remove_const<T1>::type>& other)
180 : element_(other.element_)
181 {}
182
184 LevelIterator(const LevelIterator<const typename remove_const<C>::type,
185 const typename remove_const<T1>::type>& other)
186 : element_(other.element_)
187 {}
188
192 bool equals(const LevelIterator<typename remove_const<C>::type,
193 typename remove_const<T1>::type>& other) const
194 {
195 return element_ == other.element_;
196 }
197
201 bool equals(const LevelIterator<const typename remove_const<C>::type,
202 const typename remove_const<T1>::type>& other) const
203 {
204 return element_ == other.element_;
205 }
206
208 T1& dereference() const
209 {
210 return *(element_->element_);
211 }
212
215 {
216 element_ = element_->coarser_;
217 }
218
221 {
222 element_ = element_->finer_;
223 }
224
229 bool isRedistributed() const
230 {
231 return element_->redistributed_;
232 }
233
239 {
240 assert(element_->redistributed_);
241 return *element_->redistributed_;
242 }
243 void addRedistributed(T1* t)
244 {
245 element_->redistributed_ = t;
246 }
247
248 void deleteRedistributed()
249 {
250 element_->redistributed_ = nullptr;
251 }
252
253 private:
254 Element* element_;
255 };
256
259
262
268
274
275
281
287
292 std::size_t levels() const;
293
296
297 private:
299 Element* finest_;
301 Element* coarsest_;
303 Element* nonAllocated_;
305 Allocator allocator_;
307 int levels_;
308 };
309
316 template<class M, class PI, class A=std::allocator<M> >
318 {
319 public:
321 typedef M MatrixOperator;
322
324 typedef typename MatrixOperator::matrix_type Matrix;
325
328
330 typedef A Allocator;
331
334
337
340
342 typedef typename Allocator::template rebind<AggregatesMap*>::other AAllocator;
343
345 typedef std::list<AggregatesMap*,AAllocator> AggregatesMapList;
346
348 typedef RedistributeInformation<ParallelInformation> RedistributeInfoType;
349
351 typedef typename Allocator::template rebind<RedistributeInfoType>::other RILAllocator;
352
354 typedef std::list<RedistributeInfoType,RILAllocator> RedistributeInfoList;
355
361 MatrixHierarchy(const MatrixOperator& fineMatrix,
363
364
366
372 template<typename O, typename T>
373 void build(const T& criterion);
374
382 template<class F>
383 void recalculateGalerkin(const F& copyFlags);
384
389 template<class V, class TA>
390 void coarsenVector(Hierarchy<BlockVector<V,TA> >& hierarchy) const;
391
397 template<class S, class TA>
398 void coarsenSmoother(Hierarchy<S,TA>& smoothers,
399 const typename SmootherTraits<S>::Arguments& args) const;
400
405 std::size_t levels() const;
406
411 std::size_t maxlevels() const;
412
413 bool hasCoarsest() const;
414
419 bool isBuilt() const;
420
425 const ParallelMatrixHierarchy& matrices() const;
426
432
437 const AggregatesMapList& aggregatesMaps() const;
438
445
446
447 typename MatrixOperator::field_type getProlongationDampingFactor() const
448 {
449 return prolongDamp_;
450 }
451
462 void getCoarsestAggregatesOnFinest(std::vector<std::size_t>& data) const;
463
464 private:
465 typedef typename ConstructionTraits<MatrixOperator>::Arguments MatrixArgs;
466 typedef typename ConstructionTraits<ParallelInformation>::Arguments CommunicationArgs;
468 AggregatesMapList aggregatesMaps_;
470 RedistributeInfoList redistributes_;
472 ParallelMatrixHierarchy matrices_;
474 ParallelInformationHierarchy parallelInformation_;
475
477 bool built_;
478
480 int maxlevels_;
481
482 typename MatrixOperator::field_type prolongDamp_;
483
487 template<class Matrix, bool print>
488 struct MatrixStats
489 {
490
494 static void stats(const Matrix& matrix)
495 {
496 DUNE_UNUSED_PARAMETER(matrix);
497 }
498 };
499
500 template<class Matrix>
501 struct MatrixStats<Matrix,true>
502 {
503 struct calc
504 {
505 typedef typename Matrix::size_type size_type;
506 typedef typename Matrix::row_type matrix_row;
507
508 calc()
509 {
510 min=std::numeric_limits<size_type>::max();
511 max=0;
512 sum=0;
513 }
514
515 void operator()(const matrix_row& row)
516 {
517 min=std::min(min, row.size());
518 max=std::max(max, row.size());
519 sum += row.size();
520 }
521
522 size_type min;
523 size_type max;
524 size_type sum;
525 };
529 static void stats(const Matrix& matrix)
530 {
531 calc c= for_each(matrix.begin(), matrix.end(), calc());
532 dinfo<<"Matrix row: min="<<c.min<<" max="<<c.max
533 <<" average="<<static_cast<double>(c.sum)/matrix.N()
534 <<std::endl;
535 }
536 };
537 };
538
542 template<class T>
543 class CoarsenCriterion : public T
544 {
545 public:
551
562 CoarsenCriterion(int maxLevel=100, int coarsenTarget=1000, double minCoarsenRate=1.2,
563 double prolongDamp=1.6, AccumulationMode accumulate=successiveAccu)
564 : AggregationCriterion(Dune::Amg::Parameters(maxLevel, coarsenTarget, minCoarsenRate, prolongDamp, accumulate))
565 {}
566
568 : AggregationCriterion(parms)
569 {}
570
571 };
572
573 template<typename M, typename C1>
574 bool repartitionAndDistributeMatrix(const M& origMatrix, M& newMatrix,
575 SequentialInformation& origSequentialInformationomm,
576 SequentialInformation*& newComm,
577 RedistributeInformation<SequentialInformation>& ri,
578 int nparts, C1& criterion)
579 {
580 DUNE_UNUSED_PARAMETER(origMatrix);
581 DUNE_UNUSED_PARAMETER(newMatrix);
582 DUNE_UNUSED_PARAMETER(origSequentialInformationomm);
583 DUNE_UNUSED_PARAMETER(newComm);
585 DUNE_UNUSED_PARAMETER(nparts);
586 DUNE_UNUSED_PARAMETER(criterion);
587 DUNE_THROW(NotImplemented, "Redistribution does not make sense in sequential code!");
588 }
589
590
591 template<typename M, typename C, typename C1>
592 bool repartitionAndDistributeMatrix(const M& origMatrix, M& newMatrix, C& origComm, C*& newComm,
593 RedistributeInformation<C>& ri,
594 int nparts, C1& criterion)
595 {
596 Timer time;
597#ifdef AMG_REPART_ON_COMM_GRAPH
598 // Done not repartition the matrix graph, but a graph of the communication scheme.
599 bool existentOnRedist=Dune::commGraphRepartition(origMatrix, origComm, nparts, newComm,
600 ri.getInterface(),
601 criterion.debugLevel()>1);
602
603#else
604 typedef Dune::Amg::MatrixGraph<const M> MatrixGraph;
605 typedef Dune::Amg::PropertiesGraph<MatrixGraph,
606 VertexProperties,
607 EdgeProperties,
608 IdentityMap,
609 IdentityMap> PropertiesGraph;
610 MatrixGraph graph(origMatrix);
611 PropertiesGraph pgraph(graph);
612 buildDependency(pgraph, origMatrix, criterion, false);
613
614#ifdef DEBUG_REPART
615 if(origComm.communicator().rank()==0)
616 std::cout<<"Original matrix"<<std::endl;
617 origComm.communicator().barrier();
618 printGlobalSparseMatrix(origMatrix, origComm, std::cout);
619#endif
620 bool existentOnRedist=Dune::graphRepartition(pgraph, origComm, nparts,
621 newComm, ri.getInterface(),
622 criterion.debugLevel()>1);
623#endif // if else AMG_REPART
624
625 if(origComm.communicator().rank()==0 && criterion.debugLevel()>1)
626 std::cout<<"Repartitioning took "<<time.elapsed()<<" seconds."<<std::endl;
627
628 ri.setSetup();
629
630#ifdef DEBUG_REPART
631 ri.checkInterface(origComm.indexSet(), newComm->indexSet(), origComm.communicator());
632#endif
633
634 redistributeMatrix(const_cast<M&>(origMatrix), newMatrix, origComm, *newComm, ri);
635
636#ifdef DEBUG_REPART
637 if(origComm.communicator().rank()==0)
638 std::cout<<"Original matrix"<<std::endl;
639 origComm.communicator().barrier();
640 if(newComm->communicator().size()>0)
641 printGlobalSparseMatrix(newMatrix, *newComm, std::cout);
642 origComm.communicator().barrier();
643#endif
644
645 if(origComm.communicator().rank()==0 && criterion.debugLevel()>1)
646 std::cout<<"Redistributing matrix took "<<time.elapsed()<<" seconds."<<std::endl;
647 return existentOnRedist;
648
649 }
650
651 template<typename M>
652 bool repartitionAndDistributeMatrix(M& origMatrix, M& newMatrix,
653 SequentialInformation& origComm,
654 SequentialInformation& newComm,
655 RedistributeInformation<SequentialInformation>& ri)
656 {
657 return true;
658 }
659
660 template<class M, class IS, class A>
662 const ParallelInformation& pinfo)
663 : matrices_(const_cast<MatrixOperator&>(fineOperator)),
664 parallelInformation_(const_cast<ParallelInformation&>(pinfo))
665 {
666 dune_static_assert((static_cast<int>(MatrixOperator::category) ==
667 static_cast<int>(SolverCategory::sequential) ||
668 static_cast<int>(MatrixOperator::category) ==
669 static_cast<int>(SolverCategory::overlapping) ||
670 static_cast<int>(MatrixOperator::category) ==
671 static_cast<int>(SolverCategory::nonoverlapping)),
672 "MatrixOperator must be of category sequential or overlapping or nonoverlapping");
673 if (static_cast<int>(MatrixOperator::category) != static_cast<int>(pinfo.getSolverCategory()))
674 DUNE_THROW(ISTLError, "MatrixOperator and ParallelInformation must belong to the same category!");
675
676 }
677
678 template<class M, class IS, class A>
679 template<typename O, typename T>
680 void MatrixHierarchy<M,IS,A>::build(const T& criterion)
681 {
682 prolongDamp_ = criterion.getProlongationDampingFactor();
683 typedef O OverlapFlags;
684 typedef typename ParallelMatrixHierarchy::Iterator MatIterator;
685 typedef typename ParallelInformationHierarchy::Iterator PInfoIterator;
686
687 static const int noints=(Dune::Amg::MAX_PROCESSES/4096>0) ? (Dune::Amg::MAX_PROCESSES/4096) : 1;
688
689 typedef bigunsignedint<sizeof(int)*8*noints> BIGINT;
690 GalerkinProduct<ParallelInformation> productBuilder;
691 MatIterator mlevel = matrices_.finest();
692 MatrixStats<typename M::matrix_type,MINIMAL_DEBUG_LEVEL<=INFO_DEBUG_LEVEL>::stats(mlevel->getmat());
693
694 PInfoIterator infoLevel = parallelInformation_.finest();
695 BIGINT finenonzeros=countNonZeros(mlevel->getmat());
696 finenonzeros = infoLevel->communicator().sum(finenonzeros);
697 BIGINT allnonzeros = finenonzeros;
698
699
700 int level = 0;
701 int rank = 0;
702
703 BIGINT unknowns = mlevel->getmat().N();
704
705 unknowns = infoLevel->communicator().sum(unknowns);
706 double dunknowns=unknowns.todouble();
707 infoLevel->buildGlobalLookup(mlevel->getmat().N());
708 redistributes_.push_back(RedistributeInfoType());
709
710 for(; level < criterion.maxLevel(); ++level, ++mlevel) {
711 assert(matrices_.levels()==redistributes_.size());
712 rank = infoLevel->communicator().rank();
713 if(rank==0 && criterion.debugLevel()>1)
714 std::cout<<"Level "<<level<<" has "<<dunknowns<<" unknowns, "<<dunknowns/infoLevel->communicator().size()
715 <<" unknowns per proc (procs="<<infoLevel->communicator().size()<<")"<<std::endl;
716
717 MatrixOperator* matrix=&(*mlevel);
718 ParallelInformation* info =&(*infoLevel);
719
720 if((
721#if HAVE_PARMETIS
722 criterion.accumulate()==successiveAccu
723#else
724 false
725#endif
726 || (criterion.accumulate()==atOnceAccu
727 && dunknowns < 30*infoLevel->communicator().size()))
728 && infoLevel->communicator().size()>1 &&
729 dunknowns/infoLevel->communicator().size() <= criterion.coarsenTarget())
730 {
731 // accumulate to fewer processors
732 Matrix* redistMat= new Matrix();
733 ParallelInformation* redistComm=0;
734 std::size_t nodomains = (std::size_t)std::ceil(dunknowns/(criterion.minAggregateSize()
735 *criterion.coarsenTarget()));
736 if( nodomains<=criterion.minAggregateSize()/2 ||
737 dunknowns <= criterion.coarsenTarget() )
738 nodomains=1;
739
740 bool existentOnNextLevel =
741 repartitionAndDistributeMatrix(mlevel->getmat(), *redistMat, *infoLevel,
742 redistComm, redistributes_.back(), nodomains,
743 criterion);
744 BIGINT unknowns = redistMat->N();
745 unknowns = infoLevel->communicator().sum(unknowns);
746 dunknowns= unknowns.todouble();
747 if(redistComm->communicator().rank()==0 && criterion.debugLevel()>1)
748 std::cout<<"Level "<<level<<" (redistributed) has "<<dunknowns<<" unknowns, "<<dunknowns/redistComm->communicator().size()
749 <<" unknowns per proc (procs="<<redistComm->communicator().size()<<")"<<std::endl;
750 MatrixArgs args(*redistMat, *redistComm);
751 mlevel.addRedistributed(ConstructionTraits<MatrixOperator>::construct(args));
752 assert(mlevel.isRedistributed());
753 infoLevel.addRedistributed(redistComm);
754 infoLevel->freeGlobalLookup();
755
756 if(!existentOnNextLevel)
757 // We do not hold any data on the redistributed partitioning
758 break;
759
760 // Work on the redistributed Matrix from now on
761 matrix = &(mlevel.getRedistributed());
762 info = &(infoLevel.getRedistributed());
763 info->buildGlobalLookup(matrix->getmat().N());
764 }
765
766 rank = info->communicator().rank();
767 if(dunknowns <= criterion.coarsenTarget())
768 // No further coarsening needed
769 break;
770
771 typedef PropertiesGraphCreator<MatrixOperator> GraphCreator;
772 typedef typename GraphCreator::PropertiesGraph PropertiesGraph;
773 typedef typename GraphCreator::GraphTuple GraphTuple;
774
775 typedef typename PropertiesGraph::VertexDescriptor Vertex;
776
777 std::vector<bool> excluded(matrix->getmat().N(), false);
778
779 GraphTuple graphs = GraphCreator::create(*matrix, excluded, *info, OverlapFlags());
780
781 AggregatesMap* aggregatesMap=new AggregatesMap(get<1>(graphs)->maxVertex()+1);
782
783 aggregatesMaps_.push_back(aggregatesMap);
784
785 Timer watch;
786 watch.reset();
787 int noAggregates, isoAggregates, oneAggregates, skippedAggregates;
788
789 tie(noAggregates, isoAggregates, oneAggregates, skippedAggregates) =
790 aggregatesMap->buildAggregates(matrix->getmat(), *(get<1>(graphs)), criterion, level==0);
791
792 if(rank==0 && criterion.debugLevel()>2)
793 std::cout<<" Have built "<<noAggregates<<" aggregates totally ("<<isoAggregates<<" isolated aggregates, "<<
794 oneAggregates<<" aggregates of one vertex, and skipped "<<
795 skippedAggregates<<" aggregates)."<<std::endl;
796#ifdef TEST_AGGLO
797 {
798 // calculate size of local matrix in the distributed direction
799 int start, end, overlapStart, overlapEnd;
800 int procs=info->communicator().rank();
801 int n = UNKNOWNS/procs; // number of unknowns per process
802 int bigger = UNKNOWNS%procs; // number of process with n+1 unknows
803
804 // Compute owner region
805 if(rank<bigger) {
806 start = rank*(n+1);
807 end = (rank+1)*(n+1);
808 }else{
809 start = bigger + rank * n;
810 end = bigger + (rank + 1) * n;
811 }
812
813 // Compute overlap region
814 if(start>0)
815 overlapStart = start - 1;
816 else
817 overlapStart = start;
818
819 if(end<UNKNOWNS)
820 overlapEnd = end + 1;
821 else
822 overlapEnd = end;
823
824 assert((UNKNOWNS)*(overlapEnd-overlapStart)==aggregatesMap->noVertices());
825 for(int j=0; j< UNKNOWNS; ++j)
826 for(int i=0; i < UNKNOWNS; ++i)
827 {
828 if(i>=overlapStart && i<overlapEnd)
829 {
830 int no = (j/2)*((UNKNOWNS)/2)+i/2;
831 (*aggregatesMap)[j*(overlapEnd-overlapStart)+i-overlapStart]=no;
832 }
833 }
834 }
835#endif
836 if(criterion.debugLevel()>1 && info->communicator().rank()==0)
837 std::cout<<"aggregating finished."<<std::endl;
838
839 BIGINT gnoAggregates=noAggregates;
840 gnoAggregates = info->communicator().sum(gnoAggregates);
841 double dgnoAggregates = gnoAggregates.todouble();
842#ifdef TEST_AGGLO
843 BIGINT gnoAggregates=((UNKNOWNS)/2)*((UNKNOWNS)/2);
844#endif
845
846 if(criterion.debugLevel()>2 && rank==0)
847 std::cout << "Building "<<dgnoAggregates<<" aggregates took "<<watch.elapsed()<<" seconds."<<std::endl;
848
849 if(dgnoAggregates==0 || dunknowns/dgnoAggregates<criterion.minCoarsenRate())
850 {
851 if(rank==0)
852 {
853 if(dgnoAggregates>0)
854 std::cerr << "Stopped coarsening because of rate breakdown "<<dunknowns<<"/"<<dgnoAggregates
855 <<"="<<dunknowns/dgnoAggregates<<"<"
856 <<criterion.minCoarsenRate()<<std::endl;
857 else
858 std::cerr<< "Could not build any aggregates. Probably no connected nodes."<<std::endl;
859 }
860 aggregatesMap->free();
861 delete aggregatesMap;
862 aggregatesMaps_.pop_back();
863
864 if(criterion.accumulate() && mlevel.isRedistributed() && info->communicator().size()>1) {
865 // coarse level matrix was already redistributed, but to more than 1 process
866 // Therefore need to delete the redistribution. Further down it will
867 // then be redistributed to 1 process
868 delete &(mlevel.getRedistributed().getmat());
869 mlevel.deleteRedistributed();
870 delete &(infoLevel.getRedistributed());
871 infoLevel.deleteRedistributed();
872 redistributes_.back().resetSetup();
873 }
874
875 break;
876 }
877 unknowns = noAggregates;
878 dunknowns = dgnoAggregates;
879
880 CommunicationArgs commargs(info->communicator(),info->getSolverCategory());
881 parallelInformation_.addCoarser(commargs);
882
883 ++infoLevel; // parallel information on coarse level
884
886 get(VertexVisitedTag(), *(get<1>(graphs)));
887
888 watch.reset();
889 int aggregates = IndicesCoarsener<ParallelInformation,OverlapFlags>
890 ::coarsen(*info,
891 *(get<1>(graphs)),
892 visitedMap,
893 *aggregatesMap,
894 *infoLevel,
895 noAggregates);
896 GraphCreator::free(graphs);
897
898 if(criterion.debugLevel()>2) {
899 if(rank==0)
900 std::cout<<"Coarsening of index sets took "<<watch.elapsed()<<" seconds."<<std::endl;
901 }
902
903 watch.reset();
904
905 infoLevel->buildGlobalLookup(aggregates);
906 AggregatesPublisher<Vertex,OverlapFlags,ParallelInformation>::publish(*aggregatesMap,
907 *info,
908 infoLevel->globalLookup());
909
910
911 if(criterion.debugLevel()>2) {
912 if(rank==0)
913 std::cout<<"Communicating global aggregate numbers took "<<watch.elapsed()<<" seconds."<<std::endl;
914 }
915
916 watch.reset();
917 std::vector<bool>& visited=excluded;
918
919 typedef std::vector<bool>::iterator Iterator;
921 Iterator end = visited.end();
922 for(Iterator iter= visited.begin(); iter != end; ++iter)
923 *iter=false;
924
925 VisitedMap2 visitedMap2(visited.begin(), Dune::IdentityMap());
926
927 typename MatrixOperator::matrix_type* coarseMatrix;
928
929 coarseMatrix = productBuilder.build(*(get<0>(graphs)), visitedMap2,
930 *info,
931 *aggregatesMap,
932 aggregates,
933 OverlapFlags());
934 dverb<<"Building of sparsity pattern took "<<watch.elapsed()<<std::endl;
935 watch.reset();
936 info->freeGlobalLookup();
937
938 delete get<0>(graphs);
939 productBuilder.calculate(matrix->getmat(), *aggregatesMap, *coarseMatrix, *infoLevel, OverlapFlags());
940
941 if(criterion.debugLevel()>2) {
942 if(rank==0)
943 std::cout<<"Calculation entries of Galerkin product took "<<watch.elapsed()<<" seconds."<<std::endl;
944 }
945
946 BIGINT nonzeros = countNonZeros(*coarseMatrix);
947 allnonzeros = allnonzeros + infoLevel->communicator().sum(nonzeros);
948 MatrixArgs args(*coarseMatrix, *infoLevel);
949
950 matrices_.addCoarser(args);
951 redistributes_.push_back(RedistributeInfoType());
952 } // end level loop
953
954
955 infoLevel->freeGlobalLookup();
956
957 built_=true;
958 AggregatesMap* aggregatesMap=new AggregatesMap(0);
959 aggregatesMaps_.push_back(aggregatesMap);
960
961 if(criterion.debugLevel()>0) {
962 if(level==criterion.maxLevel()) {
963 BIGINT unknowns = mlevel->getmat().N();
964 unknowns = infoLevel->communicator().sum(unknowns);
965 double dunknowns = unknowns.todouble();
966 if(rank==0 && criterion.debugLevel()>1) {
967 std::cout<<"Level "<<level<<" has "<<dunknowns<<" unknowns, "<<dunknowns/infoLevel->communicator().size()
968 <<" unknowns per proc (procs="<<infoLevel->communicator().size()<<")"<<std::endl;
969 }
970 }
971 }
972
973 if(criterion.accumulate() && !redistributes_.back().isSetup() &&
974 infoLevel->communicator().size()>1) {
975#if HAVE_MPI && !HAVE_PARMETIS
976 if(criterion.accumulate()==successiveAccu &&
977 infoLevel->communicator().rank()==0)
978 std::cerr<<"Successive accumulation of data on coarse levels only works with ParMETIS installed."
979 <<" Fell back to accumulation to one domain on coarsest level"<<std::endl;
980#endif
981
982 // accumulate to fewer processors
983 Matrix* redistMat= new Matrix();
984 ParallelInformation* redistComm=0;
985 int nodomains = 1;
986
987 repartitionAndDistributeMatrix(mlevel->getmat(), *redistMat, *infoLevel,
988 redistComm, redistributes_.back(), nodomains,criterion);
989 MatrixArgs args(*redistMat, *redistComm);
990 BIGINT unknowns = redistMat->N();
991 unknowns = infoLevel->communicator().sum(unknowns);
992
993 if(redistComm->communicator().rank()==0 && criterion.debugLevel()>1) {
994 double dunknowns= unknowns.todouble();
995 std::cout<<"Level "<<level<<" redistributed has "<<dunknowns<<" unknowns, "<<dunknowns/redistComm->communicator().size()
996 <<" unknowns per proc (procs="<<redistComm->communicator().size()<<")"<<std::endl;
997 }
998 mlevel.addRedistributed(ConstructionTraits<MatrixOperator>::construct(args));
999 infoLevel.addRedistributed(redistComm);
1000 infoLevel->freeGlobalLookup();
1001 }
1002
1003 int levels = matrices_.levels();
1004 maxlevels_ = parallelInformation_.finest()->communicator().max(levels);
1005 assert(matrices_.levels()==redistributes_.size());
1006 if(hasCoarsest() && rank==0 && criterion.debugLevel()>1)
1007 std::cout<<"operator complexity: "<<allnonzeros.todouble()/finenonzeros.todouble()<<std::endl;
1008
1009 }
1010
1011 template<class M, class IS, class A>
1014 {
1015 return matrices_;
1016 }
1017
1018 template<class M, class IS, class A>
1021 {
1022 return parallelInformation_;
1023 }
1024
1025 template<class M, class IS, class A>
1026 void MatrixHierarchy<M,IS,A>::getCoarsestAggregatesOnFinest(std::vector<std::size_t>& data) const
1027 {
1028 int levels=aggregatesMaps().size();
1029 int maxlevels=parallelInformation_.finest()->communicator().max(levels);
1030 std::size_t size=(*(aggregatesMaps().begin()))->noVertices();
1031 // We need an auxiliary vector for the consecutive prolongation.
1032 std::vector<std::size_t> tmp;
1033 std::vector<std::size_t> *coarse, *fine;
1034
1035 // make sure the allocated space suffices.
1036 tmp.reserve(size);
1037 data.reserve(size);
1038
1039 // Correctly assign coarse and fine for the first prolongation such that
1040 // we end up in data in the end.
1041 if(levels%2==0) {
1042 coarse=&tmp;
1043 fine=&data;
1044 }else{
1045 coarse=&data;
1046 fine=&tmp;
1047 }
1048
1049 // Number the unknowns on the coarsest level consecutively for each process.
1050 if(levels==maxlevels) {
1051 const AggregatesMap& map = *(*(++aggregatesMaps().rbegin()));
1052 std::size_t m=0;
1053
1054 for(typename AggregatesMap::const_iterator iter = map.begin(); iter != map.end(); ++iter)
1055 if(*iter< AggregatesMap::ISOLATED)
1056 m=std::max(*iter,m);
1057
1058 coarse->resize(m+1);
1059 std::size_t i=0;
1060 srand((unsigned)std::clock());
1061 std::set<size_t> used;
1062 for(typename std::vector<std::size_t>::iterator iter=coarse->begin(); iter != coarse->end();
1063 ++iter, ++i)
1064 {
1065 std::pair<std::set<std::size_t>::iterator,bool> ibpair
1066 = used.insert(static_cast<std::size_t>((((double)rand())/(RAND_MAX+1.0)))*coarse->size());
1067
1068 while(!ibpair.second)
1069 ibpair = used.insert(static_cast<std::size_t>((((double)rand())/(RAND_MAX+1.0))*coarse->size()));
1070 *iter=*(ibpair.first);
1071 }
1072 }
1073
1074 typename ParallelInformationHierarchy::Iterator pinfo = parallelInformation().coarsest();
1075 --pinfo;
1076
1077 // Now consecutively project the numbers to the finest level.
1078 for(typename AggregatesMapList::const_reverse_iterator aggregates=++aggregatesMaps().rbegin();
1079 aggregates != aggregatesMaps().rend(); ++aggregates,--levels) {
1080
1081 fine->resize((*aggregates)->noVertices());
1082 fine->assign(fine->size(), 0);
1083 Transfer<typename AggregatesMap::AggregateDescriptor, std::vector<std::size_t>, ParallelInformation>
1084 ::prolongateVector(*(*aggregates), *coarse, *fine, static_cast<std::size_t>(1), *pinfo);
1085 --pinfo;
1086 std::swap(coarse, fine);
1087 }
1088
1089 // Assertion to check that we really projected to data on the last step.
1090 assert(coarse==&data);
1091 }
1092
1093 template<class M, class IS, class A>
1096 {
1097 return aggregatesMaps_;
1098 }
1099 template<class M, class IS, class A>
1102 {
1103 return redistributes_;
1104 }
1105
1106 template<class M, class IS, class A>
1108 {
1109 typedef typename AggregatesMapList::reverse_iterator AggregatesMapIterator;
1110 typedef typename ParallelMatrixHierarchy::Iterator Iterator;
1111 typedef typename ParallelInformationHierarchy::Iterator InfoIterator;
1112
1113 AggregatesMapIterator amap = aggregatesMaps_.rbegin();
1114 InfoIterator info = parallelInformation_.coarsest();
1115 for(Iterator level=matrices_.coarsest(), finest=matrices_.finest(); level != finest; --level, --info, ++amap) {
1116 (*amap)->free();
1117 delete *amap;
1118 delete &level->getmat();
1119 if(level.isRedistributed())
1120 delete &(level.getRedistributed().getmat());
1121 }
1122 delete *amap;
1123 }
1124
1125 template<class M, class IS, class A>
1126 template<class V, class TA>
1128 {
1129 assert(hierarchy.levels()==1);
1130 typedef typename ParallelMatrixHierarchy::ConstIterator Iterator;
1131 typedef typename RedistributeInfoList::const_iterator RIter;
1132 RIter redist = redistributes_.begin();
1133
1134 Iterator matrix = matrices_.finest(), coarsest = matrices_.coarsest();
1135 int level=0;
1136 if(redist->isSetup())
1137 hierarchy.addRedistributedOnCoarsest(matrix.getRedistributed().getmat().N());
1138 Dune::dvverb<<"Level "<<level<<" has "<<matrices_.finest()->getmat().N()<<" unknowns!"<<std::endl;
1139
1140 while(matrix != coarsest) {
1141 ++matrix; ++level; ++redist;
1142 Dune::dvverb<<"Level "<<level<<" has "<<matrix->getmat().N()<<" unknowns!"<<std::endl;
1143
1144 hierarchy.addCoarser(matrix->getmat().N());
1145 if(redist->isSetup())
1146 hierarchy.addRedistributedOnCoarsest(matrix.getRedistributed().getmat().N());
1147
1148 }
1149
1150 }
1151
1152 template<class M, class IS, class A>
1153 template<class S, class TA>
1155 const typename SmootherTraits<S>::Arguments& sargs) const
1156 {
1157 assert(smoothers.levels()==0);
1158 typedef typename ParallelMatrixHierarchy::ConstIterator MatrixIterator;
1159 typedef typename ParallelInformationHierarchy::ConstIterator PinfoIterator;
1160 typedef typename AggregatesMapList::const_iterator AggregatesIterator;
1161
1162 typename ConstructionTraits<S>::Arguments cargs;
1163 cargs.setArgs(sargs);
1164 PinfoIterator pinfo = parallelInformation_.finest();
1165 AggregatesIterator aggregates = aggregatesMaps_.begin();
1166 int level=0;
1167 for(MatrixIterator matrix = matrices_.finest(), coarsest = matrices_.coarsest();
1168 matrix != coarsest; ++matrix, ++pinfo, ++aggregates, ++level) {
1169 cargs.setMatrix(matrix->getmat(), **aggregates);
1170 cargs.setComm(*pinfo);
1171 smoothers.addCoarser(cargs);
1172 }
1173 if(maxlevels()>levels()) {
1174 // This is not the globally coarsest level and therefore smoothing is needed
1175 cargs.setMatrix(matrices_.coarsest()->getmat(), **aggregates);
1176 cargs.setComm(*pinfo);
1177 smoothers.addCoarser(cargs);
1178 ++level;
1179 }
1180 }
1181
1182 template<class M, class IS, class A>
1183 template<class F>
1185 {
1186 typedef typename AggregatesMapList::iterator AggregatesMapIterator;
1187 typedef typename ParallelMatrixHierarchy::Iterator Iterator;
1188 typedef typename ParallelInformationHierarchy::Iterator InfoIterator;
1189
1190 AggregatesMapIterator amap = aggregatesMaps_.begin();
1191 BaseGalerkinProduct productBuilder;
1192 InfoIterator info = parallelInformation_.finest();
1193 typename RedistributeInfoList::iterator riIter = redistributes_.begin();
1194 Iterator level = matrices_.finest(), coarsest=matrices_.coarsest();
1195 if(level.isRedistributed()) {
1196 info->buildGlobalLookup(info->indexSet().size());
1197 redistributeMatrixEntries(const_cast<Matrix&>(level->getmat()),
1198 const_cast<Matrix&>(level.getRedistributed().getmat()),
1199 *info,info.getRedistributed(), *riIter);
1200 info->freeGlobalLookup();
1201 }
1202
1203 for(; level!=coarsest; ++amap) {
1204 const Matrix& fine = (level.isRedistributed() ? level.getRedistributed() : *level).getmat();
1205 ++level;
1206 ++info;
1207 ++riIter;
1208 productBuilder.calculate(fine, *(*amap), const_cast<Matrix&>(level->getmat()), *info, copyFlags);
1209 if(level.isRedistributed()) {
1210 info->buildGlobalLookup(info->indexSet().size());
1211 redistributeMatrixEntries(const_cast<Matrix&>(level->getmat()),
1212 const_cast<Matrix&>(level.getRedistributed().getmat()), *info,
1213 info.getRedistributed(), *riIter);
1214 info->freeGlobalLookup();
1215 }
1216 }
1217 }
1218
1219 template<class M, class IS, class A>
1221 {
1222 return matrices_.levels();
1223 }
1224
1225 template<class M, class IS, class A>
1227 {
1228 return maxlevels_;
1229 }
1230
1231 template<class M, class IS, class A>
1233 {
1234 return levels()==maxlevels() &&
1235 (!matrices_.coarsest().isRedistributed() ||matrices_.coarsest()->getmat().N()>0);
1236 }
1237
1238 template<class M, class IS, class A>
1240 {
1241 return built_;
1242 }
1243
1244 template<class T, class A>
1246 : finest_(0), coarsest_(0), nonAllocated_(0), allocator_(), levels_(0)
1247 {}
1248
1249 template<class T, class A>
1251 : allocator_()
1252 {
1253 finest_ = allocator_.allocate(1,0);
1254 finest_->element_ = &first;
1255 finest_->redistributed_ = nullptr;
1256 nonAllocated_ = finest_;
1257 coarsest_ = finest_;
1258 coarsest_->coarser_ = coarsest_->finer_ = nullptr;
1259 levels_ = 1;
1260 }
1261
1262 template<class T, class A>
1264 : allocator_()
1265 {
1266 finest_ = allocator_.allocate(1,0);
1267 finest_->element_ = first;
1268 finest_->redistributed_ = nullptr;
1269 nonAllocated_ = nullptr;
1270 coarsest_ = finest_;
1271 coarsest_->coarser_ = coarsest_->finer_ = nullptr;
1272 levels_ = 1;
1273 }
1274 template<class T, class A>
1276 {
1277 while(coarsest_) {
1278 Element* current = coarsest_;
1279 coarsest_ = coarsest_->finer_;
1280 if(current != nonAllocated_) {
1281 if(current->redistributed_)
1282 ConstructionTraits<T>::deconstruct(current->redistributed_);
1283 ConstructionTraits<T>::deconstruct(current->element_);
1284 }
1285 allocator_.deallocate(current, 1);
1286 current=nullptr;
1287 //coarsest_->coarser_ = nullptr;
1288 }
1289 }
1290
1291 template<class T, class A>
1293 : nonAllocated_(), allocator_(other.allocator_),
1294 levels_(other.levels_)
1295 {
1296 if(!other.finest_)
1297 {
1298 finest_=coarsest_=nonAllocated_=nullptr;
1299 return;
1300 }
1301 finest_=allocator_.allocate(1,0);
1302 Element* finer_ = nullptr;
1303 Element* current_ = finest_;
1304 Element* otherCurrent_ = other.finest_;
1305
1306 while(otherCurrent_)
1307 {
1308 T* t=new T(*(otherCurrent_->element_));
1309 current_->element_=t;
1310 current_->finer_=finer_;
1311 if(otherCurrent_->redistributed_)
1312 current_->redistributed_ = new T(*otherCurrent_->redistributed_);
1313 else
1314 current_->redistributed_= nullptr;
1315 finer_=current_;
1316 if(otherCurrent_->coarser_)
1317 {
1318 current_->coarser_=allocator_.allocate(1,0);
1319 current_=current_->coarser_;
1320 }else
1321 current_->coarser_=nullptr;
1322 otherCurrent_=otherCurrent_->coarser_;
1323 }
1324 coarsest_=current_;
1325 }
1326
1327 template<class T, class A>
1328 std::size_t Hierarchy<T,A>::levels() const
1329 {
1330 return levels_;
1331 }
1332
1333 template<class T, class A>
1334 void Hierarchy<T,A>::addRedistributedOnCoarsest(Arguments& args)
1335 {
1336 coarsest_->redistributed_ = ConstructionTraits<MemberType>::construct(args);
1337 }
1338
1339 template<class T, class A>
1340 void Hierarchy<T,A>::addCoarser(Arguments& args)
1341 {
1342 if(!coarsest_) {
1343 assert(!finest_);
1344 coarsest_ = allocator_.allocate(1,0);
1345 coarsest_->element_ = ConstructionTraits<MemberType>::construct(args);
1346 finest_ = coarsest_;
1347 coarsest_->finer_ = nullptr;
1348 }else{
1349 coarsest_->coarser_ = allocator_.allocate(1,0);
1350 coarsest_->coarser_->finer_ = coarsest_;
1351 coarsest_ = coarsest_->coarser_;
1352 coarsest_->element_ = ConstructionTraits<MemberType>::construct(args);
1353 }
1354 coarsest_->redistributed_ = nullptr;
1355 coarsest_->coarser_=nullptr;
1356 ++levels_;
1357 }
1358
1359
1360 template<class T, class A>
1361 void Hierarchy<T,A>::addFiner(Arguments& args)
1362 {
1363 if(!finest_) {
1364 assert(!coarsest_);
1365 finest_ = allocator_.allocate(1,0);
1366 finest_->element = ConstructionTraits<T>::construct(args);
1367 coarsest_ = finest_;
1368 coarsest_->coarser_ = coarsest_->finer_ = nullptr;
1369 }else{
1370 finest_->finer_ = allocator_.allocate(1,0);
1371 finest_->finer_->coarser_ = finest_;
1372 finest_ = finest_->finer_;
1373 finest_->finer = nullptr;
1374 finest_->element = ConstructionTraits<T>::construct(args);
1375 }
1376 ++levels_;
1377 }
1378
1379 template<class T, class A>
1381 {
1382 return Iterator(finest_);
1383 }
1384
1385 template<class T, class A>
1387 {
1388 return Iterator(coarsest_);
1389 }
1390
1391 template<class T, class A>
1393 {
1394 return ConstIterator(finest_);
1395 }
1396
1397 template<class T, class A>
1399 {
1400 return ConstIterator(coarsest_);
1401 }
1403 } // namespace Amg
1404} // namespace Dune
1405
1406#endif
Provides classes for the Coloring process of AMG.
Portable very large unsigned integers.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Class providing information about the mapping of the vertices onto aggregates.
Definition: aggregates.hh:498
Base class of all aggregation criterions.
Definition: aggregates.hh:46
The criterion describing the stop criteria for the coarsening process.
Definition: hierarchy.hh:544
T AggregationCriterion
The criterion for tagging connections as strong and nodes as isolated. This might be e....
Definition: hierarchy.hh:550
CoarsenCriterion(int maxLevel=100, int coarsenTarget=1000, double minCoarsenRate=1.2, double prolongDamp=1.6, AccumulationMode accumulate=successiveAccu)
Constructor.
Definition: hierarchy.hh:562
Traits class for generically constructing non default constructable types.
Definition: construction.hh:39
Iterator over the levels in the hierarchy.
Definition: hierarchy.hh:161
LevelIterator(const LevelIterator< const typename remove_const< C >::type, const typename remove_const< T1 >::type > &other)
Copy constructor.
Definition: hierarchy.hh:184
T1 & dereference() const
Dereference the iterator.
Definition: hierarchy.hh:208
bool isRedistributed() const
Check whether there was a redistribution at the current level.
Definition: hierarchy.hh:229
void increment()
Move to the next coarser level.
Definition: hierarchy.hh:214
void decrement()
Move to the next fine level.
Definition: hierarchy.hh:220
LevelIterator(const LevelIterator< typename remove_const< C >::type, typename remove_const< T1 >::type > &other)
Copy constructor.
Definition: hierarchy.hh:178
bool equals(const LevelIterator< typename remove_const< C >::type, typename remove_const< T1 >::type > &other) const
Equality check.
Definition: hierarchy.hh:192
T1 & getRedistributed() const
Get the redistributed container.
Definition: hierarchy.hh:238
LevelIterator()
Constructor.
Definition: hierarchy.hh:169
bool equals(const LevelIterator< const typename remove_const< C >::type, const typename remove_const< T1 >::type > &other) const
Equality check.
Definition: hierarchy.hh:201
A hierarchy of coantainers (e.g. matrices or vectors)
Definition: hierarchy.hh:68
T MemberType
The type of the container we store.
Definition: hierarchy.hh:73
LevelIterator< Hierarchy< T, A >, T > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:258
LevelIterator< const Hierarchy< T, A >, const T > ConstIterator
Type of the const iterator.
Definition: hierarchy.hh:261
A::template rebind< Element >::other Allocator
The allocator to use for the list elements.
Definition: hierarchy.hh:111
The hierarchies build by the coarsening process.
Definition: hierarchy.hh:318
Dune::Amg::Hierarchy< ParallelInformation, Allocator > ParallelInformationHierarchy
The type of the parallel informarion hierarchy.
Definition: hierarchy.hh:339
std::list< AggregatesMap *, AAllocator > AggregatesMapList
The type of the aggregates maps list.
Definition: hierarchy.hh:345
PI ParallelInformation
The type of the index set.
Definition: hierarchy.hh:327
Dune::Amg::Hierarchy< MatrixOperator, Allocator > ParallelMatrixHierarchy
The type of the parallel matrix hierarchy.
Definition: hierarchy.hh:336
A Allocator
The allocator to use.
Definition: hierarchy.hh:330
RedistributeInformation< ParallelInformation > RedistributeInfoType
The type of the redistribute information.
Definition: hierarchy.hh:348
std::list< RedistributeInfoType, RILAllocator > RedistributeInfoList
The type of the list of redistribute information.
Definition: hierarchy.hh:354
Dune::Amg::AggregatesMap< typename MatrixGraph< Matrix >::VertexDescriptor > AggregatesMap
The type of the aggregates map we use.
Definition: hierarchy.hh:333
Allocator::template rebind< AggregatesMap * >::other AAllocator
Allocator for pointers.
Definition: hierarchy.hh:342
MatrixOperator::matrix_type Matrix
The type of the matrix.
Definition: hierarchy.hh:324
Allocator::template rebind< RedistributeInfoType >::other RILAllocator
Allocator for RedistributeInfoType.
Definition: hierarchy.hh:351
M MatrixOperator
The type of the matrix operator.
Definition: hierarchy.hh:321
All parameters for AMG.
Definition: parameters.hh:391
Attaches properties to the edges and vertices of a graph.
Definition: graph.hh:977
Graph::VertexDescriptor VertexDescriptor
The vertex descriptor.
Definition: graph.hh:987
Facade class for stl conformant bidirectional iterators.
Definition: iteratorfacades.hh:273
A vector of blocks with memory management.
Definition: bvector.hh:254
derive error class from the base class in common
Definition: istlexception.hh:16
Adapter to turn a random access iterator into a property map.
Definition: propertymap.hh:109
Enables iteration over all entities of a given codimension and level of a grid. See also the document...
Definition: leveliterator.hh:31
A generic dynamic dense matrix.
Definition: matrix.hh:25
VariableBlockVector< T, A >::window_type row_type
The type implementing a matrix row.
Definition: matrix.hh:38
A::size_type size_type
Type for indices and sizes.
Definition: matrix.hh:41
A simple stop watch.
Definition: timer.hh:52
void reset()
Reset timer while keeping the running/stopped state.
Definition: timer.hh:66
double elapsed() const
Get elapsed user-time from last reset until now/last stop in seconds.
Definition: timer.hh:86
Portable very large unsigned integers.
Definition: bigunsignedint.hh:42
Helper classes for the construction of classes without empty constructor.
Provides classes for initializing the link attributes of a matrix graph.
Provides a class for building the galerkin product based on a aggregation scheme.
Provdes class for identifying aggregates globally.
Provides classes for building the matrix graph.
#define dune_static_assert(COND, MSG)
Helper template so that compilation fails if condition is not true.
Definition: static_assert.hh:79
#define DUNE_THROW(E, m)
Definition: exceptions.hh:244
const AggregatesMapList & aggregatesMaps() const
Get the hierarchy of the mappings of the nodes onto aggregates.
Definition: hierarchy.hh:1095
bool isBuilt() const
Whether the hierarchy was built.
Definition: hierarchy.hh:1239
Hierarchy(const Hierarchy &other)
Copy constructor.
Definition: hierarchy.hh:1292
tuple< int, int, int, int > buildAggregates(const M &matrix, G &graph, const C &criterion, bool finestLevel)
Build the aggregates.
std::size_t levels() const
Get the number of levels in the hierarchy.
Definition: hierarchy.hh:1328
std::size_t levels() const
Get the number of levels in the hierarchy.
Definition: hierarchy.hh:1220
ConstIterator coarsest() const
Get an iterator positioned at the coarsest level.
Definition: hierarchy.hh:1398
void addCoarser(Arguments &args)
Add an element on a coarser level.
Definition: hierarchy.hh:1340
const RedistributeInfoList & redistributeInformation() const
Get the hierachy of the information about redistributions,.
Definition: hierarchy.hh:1101
void coarsenVector(Hierarchy< BlockVector< V, TA > > &hierarchy) const
Coarsen the vector hierarchy according to the matrix hierarchy.
Definition: hierarchy.hh:1127
const ParallelInformationHierarchy & parallelInformation() const
Get the hierarchy of the parallel data distribution information.
Definition: hierarchy.hh:1020
void addFiner(Arguments &args)
Add an element on a finer level.
Definition: hierarchy.hh:1361
const ParallelMatrixHierarchy & matrices() const
Get the matrix hierarchy.
Definition: hierarchy.hh:1013
std::size_t maxlevels() const
Get the max number of levels in the hierarchy of processors.
Definition: hierarchy.hh:1226
static const V ISOLATED
Identifier of isolated vertices.
Definition: aggregates.hh:509
void recalculateGalerkin(const F &copyFlags)
Recalculate the galerkin products.
Definition: hierarchy.hh:1184
const void * Arguments
A type holding all the arguments needed to call the constructor.
Definition: construction.hh:45
std::size_t noVertices() const
Get the number of vertices.
Iterator coarsest()
Get an iterator positioned at the coarsest level.
Definition: hierarchy.hh:1386
Hierarchy(MemberType &first)
Construct a new hierarchy.
Definition: hierarchy.hh:1250
ConstIterator finest() const
Get an iterator positioned at the finest level.
Definition: hierarchy.hh:1392
Hierarchy(MemberType *first)
Construct a new hierarchy.
Definition: hierarchy.hh:1263
Hierarchy()
Construct a new empty hierarchy.
Definition: hierarchy.hh:1245
~Hierarchy()
Destructor.
Definition: hierarchy.hh:1275
AccumulationMode
Identifiers for the different accumulation modes.
Definition: parameters.hh:230
static T * construct(Arguments &args)
Construct an object with the specified arguments.
Definition: construction.hh:53
Iterator finest()
Get an iterator positioned at the finest level.
Definition: hierarchy.hh:1380
void build(const T &criterion)
Build the matrix hierarchy using aggregation.
Definition: hierarchy.hh:680
MatrixHierarchy(const MatrixOperator &fineMatrix, const ParallelInformation &pinfo=ParallelInformation())
Constructor.
Definition: hierarchy.hh:661
void free()
Free the allocated memory.
static void deconstruct(T *t)
Destroys an object.
Definition: construction.hh:62
void coarsenSmoother(Hierarchy< S, TA > &smoothers, const typename SmootherTraits< S >::Arguments &args) const
Coarsen the smoother hierarchy according to the matrix hierarchy.
Definition: hierarchy.hh:1154
void buildDependency(G &graph, const typename C::Matrix &matrix, C criterion, bool finestLevel)
Build the dependency of the matrix graph.
void getCoarsestAggregatesOnFinest(std::vector< std::size_t > &data) const
Get the mapping of fine level unknowns to coarse level aggregates.
Definition: hierarchy.hh:1026
@ MAX_PROCESSES
Hard limit for the number of processes allowed.
Definition: hierarchy.hh:56
@ atOnceAccu
Accumulate data to on process at once.
Definition: parameters.hh:242
@ successiveAccu
Successively accumulate to fewer processes.
Definition: parameters.hh:246
int countNonZeros(const M &matrix)
Get the number of nonzero fields in the matrix.
Definition: matrixutils.hh:155
DVVerbType dvverb(std::cout)
stream for very verbose output.
Definition: stdstreams.hh:94
DInfoType dinfo(std::cout)
Stream for informative output.
Definition: stdstreams.hh:139
DVerbType dverb(std::cout)
Singleton of verbose debug stream.
Definition: stdstreams.hh:115
Provides a map between global and local indices.
Provides a class for building the index set and remote indices on the coarse level.
Functionality for redistributing a sparse matrix.
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignment.hh:14
void redistributeMatrix(M &origMatrix, M &newMatrix, C &origComm, C &newComm, RedistributeInformation< C > &ri)
Redistribute a matrix according to given domain decompositions.
Definition: matrixredistribute.hh:838
bool graphRepartition(const G &graph, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm, int nparts, Dune::OwnerOverlapCopyCommunication< T1, T2 > *&outcomm, RedistributeInterface &redistInf, bool verbose=false)
execute a graph repartition for a giving graph and indexset.
Definition: repartition.hh:1233
Classes for the generic construction and application of the smoothers.
Standard Dune debug streams.
Traits class for getting the attribute class of a smoother.
Definition: smoother.hh:64
Tag idnetifying the visited property of a vertex.
Definition: properties.hh:27
Get the N-th element of a tuple.
Definition: tuples.hh:458
A property map that applies the identity function to integers.
Definition: propertymap.hh:294
Selector for the property map type.
Definition: propertymap.hh:321
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:22
@ nonoverlapping
Category for on overlapping solvers.
Definition: solvercategory.hh:24
@ overlapping
Category for ovelapping solvers.
Definition: solvercategory.hh:26
Removes a const qualifier while preserving others.
Definition: typetraits.hh:176
A simple timing class.
Prolongation and restriction for amg.
Fallback implementation of the std::tuple class.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)