19#include <dune/istl/solverregistry.hh>
47 template<
class M,
class X,
class S,
class P,
class K,
class A>
63 template<
class M,
class X,
class S,
class PI=SequentialInformation,
64 class A=std::allocator<X> >
67 template<
class M1,
class X1,
class S1,
class P1,
class K1,
class A1>
210 std::size_t levels();
212 std::size_t maxlevels();
242 void createCriterionAndHierarchies(std::shared_ptr<const Operator> matrixptr,
243 const PI& pinfo,
const Norm&,
245 std::true_type compiles = std::true_type());
247 void createCriterionAndHierarchies(std::shared_ptr<const Operator> matrixptr,
248 const PI& pinfo,
const Norm&,
256 void createHierarchies(C& criterion, std::shared_ptr<const Operator> matrixptr,
265 void createHierarchies(C& criterion,
266 const std::shared_ptr<const Operator>& matrixptr,
292 typename OperatorHierarchy::RedistributeInfoList::const_iterator
redist;
296 typename OperatorHierarchy::AggregatesMapList::const_iterator
aggregates;
320 void mgc(LevelContext& levelContext);
330 void moveToFineLevel(LevelContext& levelContext,
bool processedFineLevel);
336 bool moveToCoarseLevel(LevelContext& levelContext);
342 void initIteratorsWithFineLevel(LevelContext& levelContext);
345 std::shared_ptr<OperatorHierarchy> matrices_;
349 std::shared_ptr<Hierarchy<Smoother,A> > smoothers_;
351 std::shared_ptr<CoarseSolver> solver_;
353 std::shared_ptr<Hierarchy<Range,A>> rhs_;
355 std::shared_ptr<Hierarchy<Domain,A>> lhs_;
357 std::shared_ptr<Hierarchy<Domain,A>> update_;
361 std::shared_ptr<ScalarProduct> scalarProduct_;
365 std::size_t preSteps_;
367 std::size_t postSteps_;
368 bool buildHierarchy_;
370 bool coarsesolverconverged;
371 std::shared_ptr<Smoother> coarseSmoother_;
375 std::size_t verbosity_;
379 std::string operator()(
const std::string& str)
381 std::stringstream retval;
382 std::ostream_iterator<char> out(retval);
383 std::transform(str.begin(), str.end(), out,
385 return std::tolower(c, std::locale::classic());
392 template<
class M,
class X,
class S,
class PI,
class A>
394 : matrices_(amg.matrices_), smootherArgs_(amg.smootherArgs_),
395 smoothers_(amg.smoothers_), solver_(amg.solver_),
396 rhs_(), lhs_(), update_(),
397 scalarProduct_(amg.scalarProduct_), gamma_(amg.gamma_),
398 preSteps_(amg.preSteps_), postSteps_(amg.postSteps_),
399 buildHierarchy_(amg.buildHierarchy_),
400 additive(amg.additive), coarsesolverconverged(amg.coarsesolverconverged),
401 coarseSmoother_(amg.coarseSmoother_),
402 category_(amg.category_),
403 verbosity_(amg.verbosity_)
406 template<
class M,
class X,
class S,
class PI,
class A>
412 rhs_(), lhs_(), update_(), scalarProduct_(0),
413 gamma_(parms.getGamma()), preSteps_(parms.getNoPreSmoothSteps()),
414 postSteps_(parms.getNoPostSmoothSteps()), buildHierarchy_(false),
415 additive(parms.getAdditive()), coarsesolverconverged(true),
419 verbosity_(parms.debugLevel())
421 assert(matrices_->isBuilt());
424 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
427 template<
class M,
class X,
class S,
class PI,
class A>
433 : smootherArgs_(smootherArgs),
435 rhs_(), lhs_(), update_(), scalarProduct_(),
436 gamma_(criterion.getGamma()), preSteps_(criterion.getNoPreSmoothSteps()),
437 postSteps_(criterion.getNoPostSmoothSteps()), buildHierarchy_(true),
438 additive(criterion.getAdditive()), coarsesolverconverged(true),
441 verbosity_(criterion.debugLevel())
444 DUNE_THROW(InvalidSolverCategory,
"Matrix and Communication must have the same SolverCategory!");
449 createHierarchies(criterion, matrixptr, pinfo);
452 template<
class M,
class X,
class S,
class PI,
class A>
457 solver_(), rhs_(), lhs_(), update_(), scalarProduct_(), buildHierarchy_(true),
458 coarsesolverconverged(true), coarseSmoother_(),
462 if (configuration.
hasKey (
"smootherIterations"))
463 smootherArgs_.iterations = configuration.
get<
int>(
"smootherIterations");
465 if (configuration.
hasKey (
"smootherRelaxation"))
466 smootherArgs_.relaxationFactor = configuration.
get<
typename SmootherArgs::RelaxationFactor>(
"smootherRelaxation");
468 auto normName = ToLower()(configuration.
get(
"strengthMeasure",
"diagonal"));
469 auto index = configuration.
get<
int>(
"diagonalRowIndex", 0);
471 if ( normName ==
"diagonal")
474 using real_type =
typename FieldTraits<field_type>::real_type;
475 std::is_convertible<field_type, real_type> compiles;
480 createCriterionAndHierarchies(matrixptr, pinfo,
Diagonal<0>(), configuration, compiles);
483 createCriterionAndHierarchies(matrixptr, pinfo,
Diagonal<1>(), configuration, compiles);
486 createCriterionAndHierarchies(matrixptr, pinfo,
Diagonal<2>(), configuration, compiles);
489 createCriterionAndHierarchies(matrixptr, pinfo,
Diagonal<3>(), configuration, compiles);
492 createCriterionAndHierarchies(matrixptr, pinfo,
Diagonal<4>(), configuration, compiles);
498 else if (normName ==
"rowsum")
499 createCriterionAndHierarchies(matrixptr, pinfo,
RowSum(), configuration);
500 else if (normName ==
"frobenius")
501 createCriterionAndHierarchies(matrixptr, pinfo, FrobeniusNorm(), configuration);
502 else if (normName ==
"one")
503 createCriterionAndHierarchies(matrixptr, pinfo, AlwaysOneNorm(), configuration);
508 template<
class M,
class X,
class S,
class PI,
class A>
513 << className<typename M::field_type>() <<
") as it is lacking a conversion to"
514 <<
className<
typename FieldTraits<typename M::field_type>::real_type>() <<
".");
517 template<
class M,
class X,
class S,
class PI,
class A>
519 void AMG<M,X,S,PI,A>::createCriterionAndHierarchies(std::shared_ptr<const Operator> matrixptr,
const PI& pinfo,
const Norm&,
const ParameterTree& configuration, std::true_type)
521 if (configuration.
get<
bool>(
"criterionSymmetric",
true))
526 createHierarchies(criterion, matrixptr, pinfo, configuration);
533 createHierarchies(criterion, matrixptr, pinfo, configuration);
537 template<
class M,
class X,
class S,
class PI,
class A>
539 void AMG<M,X,S,PI,A>::createHierarchies(C& criterion, std::shared_ptr<const Operator> matrixptr,
const PI& pinfo,
const ParameterTree& configuration)
541 if (configuration.hasKey (
"maxLevel"))
542 criterion.setMaxLevel(configuration.get<
int>(
"maxLevel"));
544 if (configuration.hasKey (
"minCoarseningRate"))
545 criterion.setMinCoarsenRate(configuration.get<
int>(
"minCoarseningRate"));
547 if (configuration.hasKey (
"coarsenTarget"))
548 criterion.setCoarsenTarget (configuration.get<
int>(
"coarsenTarget"));
550 if (configuration.hasKey (
"accumulationMode"))
552 std::string mode = ToLower()(configuration.get<std::string>(
"accumulationMode"));
555 else if ( mode ==
"atonce" )
557 else if ( mode ==
"successive")
560 DUNE_THROW(InvalidSolverFactoryConfiguration,
"Parameter accumulationMode does not allow value "
564 if (configuration.hasKey (
"prolongationDampingFactor"))
565 criterion.setProlongationDampingFactor (configuration.get<
double>(
"prolongationDampingFactor"));
567 if (configuration.hasKey(
"defaultAggregationSizeMode"))
569 auto mode = ToLower()(configuration.get<std::string>(
"defaultAggregationSizeMode"));
570 auto dim = configuration.get<std::size_t>(
"defaultAggregationDimension");
571 std::size_t maxDistance = 2;
572 if (configuration.hasKey(
"MaxAggregateDistance"))
573 maxDistance = configuration.get<std::size_t>(
"maxAggregateDistance");
574 if (mode ==
"isotropic")
575 criterion.setDefaultValuesIsotropic(dim, maxDistance);
576 else if(mode ==
"anisotropic")
577 criterion.setDefaultValuesAnisotropic(dim, maxDistance);
579 DUNE_THROW(InvalidSolverFactoryConfiguration,
"Parameter accumulationMode does not allow value "
583 if (configuration.hasKey(
"maxAggregateDistance"))
584 criterion.setMaxDistance(configuration.get<std::size_t>(
"maxAggregateDistance"));
586 if (configuration.hasKey(
"minAggregateSize"))
587 criterion.setMinAggregateSize(configuration.get<std::size_t>(
"minAggregateSize"));
589 if (configuration.hasKey(
"maxAggregateSize"))
590 criterion.setMaxAggregateSize(configuration.get<std::size_t>(
"maxAggregateSize"));
592 if (configuration.hasKey(
"maxAggregateConnectivity"))
593 criterion.setMaxConnectivity(configuration.get<std::size_t>(
"maxAggregateConnectivity"));
595 if (configuration.hasKey (
"alpha"))
596 criterion.setAlpha (configuration.get<
double> (
"alpha"));
598 if (configuration.hasKey (
"beta"))
599 criterion.setBeta (configuration.get<
double> (
"beta"));
601 if (configuration.hasKey (
"gamma"))
602 criterion.setGamma (configuration.get<std::size_t> (
"gamma"));
603 gamma_ = criterion.getGamma();
605 if (configuration.hasKey (
"additive"))
606 criterion.setAdditive (configuration.get<
bool>(
"additive"));
607 additive = criterion.getAdditive();
609 if (configuration.hasKey (
"preSteps"))
610 criterion.setNoPreSmoothSteps (configuration.get<std::size_t> (
"preSteps"));
611 preSteps_ = criterion.getNoPreSmoothSteps ();
613 if (configuration.hasKey (
"postSteps"))
614 criterion.setNoPostSmoothSteps (configuration.get<std::size_t> (
"postSteps"));
615 postSteps_ = criterion.getNoPostSmoothSteps ();
617 verbosity_ = configuration.get(
"verbosity", 0);
618 criterion.setDebugLevel (verbosity_);
620 createHierarchies(criterion, matrixptr, pinfo);
623 template <
class Matrix,
625 struct DirectSolverSelector
628 enum SolverType { umfpack, superlu,
none };
630 static constexpr SolverType solver =
631#if DISABLE_AMG_DIRECTSOLVER
633#elif HAVE_SUITESPARSE_UMFPACK
634 UMFPackMethodChooser< field_type > :: valid ? umfpack :
none ;
641 template <
class M, SolverType>
644 typedef InverseOperator<Vector,Vector> type;
645 static type* create(
const M& mat,
bool verbose,
bool reusevector )
647 DUNE_THROW(NotImplemented,
"DirectSolver not selected");
650 static std::string name () {
return "None"; }
652#if HAVE_SUITESPARSE_UMFPACK
654 struct Solver< M, umfpack >
656 typedef UMFPack< M > type;
657 static type* create(
const M& mat,
bool verbose,
bool reusevector )
659 return new type(mat, verbose, reusevector );
661 static std::string name () {
return "UMFPack"; }
666 struct Solver< M, superlu >
668 typedef SuperLU< M > type;
669 static type* create(
const M& mat,
bool verbose,
bool reusevector )
671 return new type(mat, verbose, reusevector );
673 static std::string name () {
return "SuperLU"; }
678 typedef Solver< Matrix, solver > SelectedSolver ;
679 typedef typename SelectedSolver :: type DirectSolver;
680 static constexpr bool isDirectSolver = solver !=
none;
681 static std::string name() {
return SelectedSolver :: name (); }
682 static DirectSolver* create(
const Matrix& mat,
bool verbose,
bool reusevector )
684 return SelectedSolver :: create( mat, verbose, reusevector );
688 template<
class M,
class X,
class S,
class PI,
class A>
690 void AMG<M,X,S,PI,A>::createHierarchies(C& criterion,
691 const std::shared_ptr<const Operator>& matrixptr,
695 matrices_ = std::make_shared<OperatorHierarchy>(
696 std::const_pointer_cast<Operator>(matrixptr),
699 matrices_->template build<NegateSet<typename PI::OwnerSet> >(criterion);
702 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
708 if(buildHierarchy_ && matrices_->levels()==matrices_->maxlevels()
709 && ( ! matrices_->redistributeInformation().back().isSetup() ||
710 matrices_->parallelInformation().coarsest().getRedistributed().communicator().size() ) )
713 SmootherArgs sargs(smootherArgs_);
714 sargs.iterations = 1;
717 cargs.setArgs(sargs);
718 if(matrices_->redistributeInformation().back().isSetup()) {
720 cargs.setMatrix(matrices_->matrices().coarsest().getRedistributed().getmat());
721 cargs.setComm(matrices_->parallelInformation().coarsest().getRedistributed());
723 cargs.setMatrix(matrices_->matrices().coarsest()->getmat());
724 cargs.setComm(*matrices_->parallelInformation().coarsest());
728 scalarProduct_ = createScalarProduct<X>(cargs.getComm(),category());
730 typedef DirectSolverSelector< typename M::matrix_type, X > SolverSelector;
733 if( SolverSelector::isDirectSolver &&
734 (std::is_same<ParallelInformation,SequentialInformation>::value
735 || matrices_->parallelInformation().coarsest()->communicator().size()==1
736 || (matrices_->parallelInformation().coarsest().isRedistributed()
737 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()==1
738 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()>0) )
741 if(matrices_->parallelInformation().coarsest().isRedistributed())
743 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0)
746 solver_.reset(SolverSelector::create(matrices_->matrices().coarsest().getRedistributed().getmat(),
false,
false));
753 solver_.reset(SolverSelector::create(matrices_->matrices().coarsest()->getmat(),
false,
false));
755 if(verbosity_>0 && matrices_->parallelInformation().coarsest()->communicator().rank()==0)
756 std::cout<<
"Using a direct coarse solver (" << SolverSelector::name() <<
")" << std::endl;
760 if(matrices_->parallelInformation().coarsest().isRedistributed())
762 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0)
767 solver_.reset(
new BiCGSTABSolver<X>(
const_cast<M&
>(matrices_->matrices().coarsest().getRedistributed()),
769 *coarseSmoother_, 1E-2, 1000, 0));
774 solver_.reset(
new BiCGSTABSolver<X>(
const_cast<M&
>(*matrices_->matrices().coarsest()),
776 *coarseSmoother_, 1E-2, 1000, 0));
795 if(verbosity_>0 && matrices_->parallelInformation().finest()->communicator().rank()==0)
796 std::cout<<
"Building hierarchy of "<<matrices_->maxlevels()<<
" levels "
797 <<
"(including coarse solver) took "<<watch.elapsed()<<
" seconds."<<std::endl;
801 template<
class M,
class X,
class S,
class PI,
class A>
808 typedef typename M::matrix_type
Matrix;
815 const Matrix& mat=matrices_->matrices().finest()->getmat();
816 for(RowIter row=mat.begin(); row!=mat.end(); ++row) {
817 bool isDirichlet =
true;
818 bool hasDiagonal =
false;
820 for(ColIter col=row->begin(); col!=row->end(); ++col) {
821 if(row.index()==col.index()) {
829 if(isDirichlet && hasDiagonal)
831 auto&& xEntry = Impl::asVector(x[row.index()]);
832 auto&& bEntry = Impl::asVector(b[row.index()]);
833 Impl::asMatrix(diagonal).solve(xEntry, bEntry);
837 if(smoothers_->levels()>0)
838 smoothers_->finest()->pre(x,b);
841 matrices_->parallelInformation().coarsest()->copyOwnerToAll(x,x);
842 rhs_ = std::make_shared<Hierarchy<Range,A>>(std::make_shared<Range>(b));
843 lhs_ = std::make_shared<Hierarchy<Domain,A>>(std::make_shared<Domain>(x));
844 update_ = std::make_shared<Hierarchy<Domain,A>>(std::make_shared<Domain>(x));
845 matrices_->coarsenVector(*rhs_);
846 matrices_->coarsenVector(*lhs_);
847 matrices_->coarsenVector(*update_);
853 Iterator coarsest = smoothers_->
coarsest();
854 Iterator smoother = smoothers_->finest();
855 RIterator rhs = rhs_->finest();
856 DIterator lhs = lhs_->finest();
857 if(smoothers_->levels()>1) {
859 assert(lhs_->levels()==rhs_->levels());
860 assert(smoothers_->levels()==lhs_->levels() || matrices_->levels()==matrices_->maxlevels());
861 assert(smoothers_->levels()+1==lhs_->levels() || matrices_->levels()<matrices_->maxlevels());
863 if(smoother!=coarsest)
864 for(++smoother, ++lhs, ++rhs; smoother != coarsest; ++smoother, ++lhs, ++rhs)
865 smoother->pre(*lhs,*rhs);
866 smoother->pre(*lhs,*rhs);
876 template<
class M,
class X,
class S,
class PI,
class A>
879 return matrices_->levels();
881 template<
class M,
class X,
class S,
class PI,
class A>
882 std::size_t AMG<M,X,S,PI,A>::maxlevels()
884 return matrices_->maxlevels();
888 template<
class M,
class X,
class S,
class PI,
class A>
891 LevelContext levelContext;
899 initIteratorsWithFineLevel(levelContext);
902 *levelContext.lhs = v;
903 *levelContext.rhs = d;
904 *levelContext.update=0;
905 levelContext.level=0;
909 if(postSteps_==0||matrices_->maxlevels()==1)
910 levelContext.pinfo->copyOwnerToAll(*levelContext.update, *levelContext.update);
912 v=*levelContext.update;
917 template<
class M,
class X,
class S,
class PI,
class A>
920 levelContext.smoother = smoothers_->finest();
921 levelContext.matrix = matrices_->matrices().finest();
922 levelContext.pinfo = matrices_->parallelInformation().finest();
923 levelContext.redist =
924 matrices_->redistributeInformation().begin();
925 levelContext.aggregates = matrices_->aggregatesMaps().begin();
926 levelContext.lhs = lhs_->finest();
927 levelContext.update = update_->finest();
928 levelContext.rhs = rhs_->finest();
931 template<
class M,
class X,
class S,
class PI,
class A>
933 ::moveToCoarseLevel(LevelContext& levelContext)
936 bool processNextLevel=
true;
938 if(levelContext.redist->isSetup()) {
939 levelContext.redist->redistribute(
static_cast<const Range&
>(*levelContext.rhs),
940 levelContext.rhs.getRedistributed());
941 processNextLevel = levelContext.rhs.getRedistributed().size()>0;
942 if(processNextLevel) {
945 ++levelContext.pinfo;
946 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
947 ::restrictVector(*(*levelContext.aggregates), *levelContext.rhs,
948 static_cast<const Range&
>(fineRhs.getRedistributed()),
949 *levelContext.pinfo);
954 ++levelContext.pinfo;
955 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
956 ::restrictVector(*(*levelContext.aggregates),
957 *levelContext.rhs,
static_cast<const Range&
>(*fineRhs),
958 *levelContext.pinfo);
961 if(processNextLevel) {
964 ++levelContext.update;
965 ++levelContext.matrix;
966 ++levelContext.level;
967 ++levelContext.redist;
969 if(levelContext.matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()) {
971 ++levelContext.smoother;
972 ++levelContext.aggregates;
975 *levelContext.update=0;
977 return processNextLevel;
980 template<
class M,
class X,
class S,
class PI,
class A>
982 ::moveToFineLevel(LevelContext& levelContext,
bool processNextLevel)
984 if(processNextLevel) {
985 if(levelContext.matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()) {
987 --levelContext.smoother;
988 --levelContext.aggregates;
990 --levelContext.redist;
991 --levelContext.level;
993 --levelContext.matrix;
997 --levelContext.pinfo;
999 if(levelContext.redist->isSetup()) {
1001 levelContext.lhs.getRedistributed()=0;
1002 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
1003 ::prolongateVector(*(*levelContext.aggregates), *levelContext.update, *levelContext.lhs,
1004 levelContext.lhs.getRedistributed(),
1005 matrices_->getProlongationDampingFactor(),
1006 *levelContext.pinfo, *levelContext.redist);
1008 *levelContext.lhs=0;
1009 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
1010 ::prolongateVector(*(*levelContext.aggregates), *levelContext.update, *levelContext.lhs,
1011 matrices_->getProlongationDampingFactor(),
1012 *levelContext.pinfo);
1016 if(processNextLevel) {
1017 --levelContext.update;
1021 *levelContext.update += *levelContext.lhs;
1024 template<
class M,
class X,
class S,
class PI,
class A>
1027 return IsDirectSolver< CoarseSolver>::value;
1030 template<
class M,
class X,
class S,
class PI,
class A>
1032 if(levelContext.matrix == matrices_->matrices().coarsest() && levels()==maxlevels()) {
1036 if(levelContext.redist->isSetup()) {
1037 levelContext.redist->redistribute(*levelContext.rhs, levelContext.rhs.getRedistributed());
1038 if(levelContext.rhs.getRedistributed().size()>0) {
1040 levelContext.pinfo.getRedistributed().copyOwnerToAll(levelContext.rhs.getRedistributed(),
1041 levelContext.rhs.getRedistributed());
1042 solver_->apply(levelContext.update.getRedistributed(),
1043 levelContext.rhs.getRedistributed(), res);
1045 levelContext.redist->redistributeBackward(*levelContext.update, levelContext.update.getRedistributed());
1046 levelContext.pinfo->copyOwnerToAll(*levelContext.update, *levelContext.update);
1048 levelContext.pinfo->copyOwnerToAll(*levelContext.rhs, *levelContext.rhs);
1049 solver_->apply(*levelContext.update, *levelContext.rhs, res);
1053 coarsesolverconverged =
false;
1058#ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
1059 bool processNextLevel = moveToCoarseLevel(levelContext);
1061 if(processNextLevel) {
1063 for(std::size_t i=0; i<gamma_; i++){
1065 if (levelContext.matrix == matrices_->matrices().coarsest() && levels()==maxlevels())
1068 levelContext.matrix->applyscaleadd(-1., *levelContext.lhs, *levelContext.rhs);
1073 moveToFineLevel(levelContext, processNextLevel);
1078 if(levelContext.matrix == matrices_->matrices().finest()) {
1079 coarsesolverconverged = matrices_->parallelInformation().finest()->communicator().prod(coarsesolverconverged);
1080 if(!coarsesolverconverged)
1081 DUNE_THROW(MathError,
"Coarse solver did not converge");
1089 template<
class M,
class X,
class S,
class PI,
class A>
1090 void AMG<M,X,S,PI,A>::additiveMgc(){
1093 typename ParallelInformationHierarchy::Iterator pinfo=matrices_->parallelInformation().finest();
1096 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates=matrices_->aggregatesMaps().begin();
1100 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
1101 ::restrictVector(*(*aggregates), *rhs,
static_cast<const Range&
>(*fineRhs), *pinfo);
1107 lhs = lhs_->finest();
1110 for(rhs=rhs_->finest(); rhs != rhs_->coarsest(); ++lhs, ++rhs, ++smoother) {
1113 smoother->apply(*lhs, *rhs);
1117#ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
1118 InverseOperatorResult res;
1119 pinfo->copyOwnerToAll(*rhs, *rhs);
1120 solver_->apply(*lhs, *rhs, res);
1123 DUNE_THROW(MathError,
"Coarse solver did not converge");
1132 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
1133 ::prolongateVector(*(*aggregates), *coarseLhs, *lhs, 1.0, *pinfo);
1139 template<
class M,
class X,
class S,
class PI,
class A>
1145 Iterator coarsest = smoothers_->
coarsest();
1146 Iterator smoother = smoothers_->finest();
1147 DIterator lhs = lhs_->finest();
1148 if(smoothers_->levels()>0) {
1149 if(smoother != coarsest || matrices_->levels()<matrices_->maxlevels())
1150 smoother->post(*lhs);
1151 if(smoother!=coarsest)
1152 for(++smoother, ++lhs; smoother != coarsest; ++smoother, ++lhs)
1153 smoother->post(*lhs);
1154 smoother->post(*lhs);
1161 template<
class M,
class X,
class S,
class PI,
class A>
1165 matrices_->getCoarsestAggregatesOnFinest(cont);
1171 template<
class>
struct isValidMatrix : std::false_type{};
1172 template<
class T,
int n,
int m,
class A>
struct isValidMatrix<BCRSMatrix<FieldMatrix<T,n,m>, A>> : std::true_type{};
1175 std::shared_ptr<Dune::Preconditioner<typename OP::element_type::domain_type, typename OP::element_type::range_type> >
1176 makeAMG(
const OP& op,
const std::string& smoother,
const Dune::ParameterTree& config)
const
1181 template<
class M,
class X,
class Y>
1182 std::shared_ptr<Dune::Preconditioner<X,Y> >
1183 makeAMG(
const std::shared_ptr<MatrixAdapter<M,X,Y>>& op,
const std::string& smoother,
1186 using OP = MatrixAdapter<M,X,Y>;
1188 if(smoother ==
"ssor")
1189 return std::make_shared<Amg::AMG<OP, X, SeqSSOR<M,X,Y>>>(op, config);
1190 if(smoother ==
"sor")
1191 return std::make_shared<Amg::AMG<OP, X, SeqSOR<M,X,Y>>>(op, config);
1192 if(smoother ==
"jac")
1193 return std::make_shared<Amg::AMG<OP, X, SeqJac<M,X,Y>>>(op, config);
1194 if(smoother ==
"gs")
1195 return std::make_shared<Amg::AMG<OP, X, SeqGS<M,X,Y>>>(op, config);
1196 if(smoother ==
"ilu")
1197 return std::make_shared<Amg::AMG<OP, X, SeqILU<M,X,Y>>>(op, config);
1202 template<
class M,
class X,
class Y,
class C>
1203 std::shared_ptr<Dune::Preconditioner<X,Y> >
1204 makeAMG(
const std::shared_ptr<OverlappingSchwarzOperator<M,X,Y,C>>& op,
const std::string& smoother,
1207 using OP = OverlappingSchwarzOperator<M,X,Y,C>;
1209 auto cop = std::static_pointer_cast<const OP>(op);
1211 if(smoother ==
"ssor")
1212 return std::make_shared<Amg::AMG<OP, X, BlockPreconditioner<X,Y,C,SeqSSOR<M,X,Y>>,C>>(cop, config, op->getCommunication());
1213 if(smoother ==
"sor")
1214 return std::make_shared<Amg::AMG<OP, X, BlockPreconditioner<X,Y,C,SeqSOR<M,X,Y>>,C>>(cop, config, op->getCommunication());
1215 if(smoother ==
"jac")
1216 return std::make_shared<Amg::AMG<OP, X, BlockPreconditioner<X,Y,C,SeqJac<M,X,Y>>,C>>(cop, config, op->getCommunication());
1217 if(smoother ==
"gs")
1218 return std::make_shared<Amg::AMG<OP, X, BlockPreconditioner<X,Y,C,SeqGS<M,X,Y>>,C>>(cop, config, op->getCommunication());
1219 if(smoother ==
"ilu")
1220 return std::make_shared<Amg::AMG<OP, X, BlockPreconditioner<X,Y,C,SeqILU<M,X,Y>>,C>>(cop, config, op->getCommunication());
1225 template<
class M,
class X,
class Y,
class C>
1226 std::shared_ptr<Dune::Preconditioner<X,Y> >
1227 makeAMG(
const std::shared_ptr<NonoverlappingSchwarzOperator<M,X,Y,C>>& op,
const std::string& smoother,
1230 using OP = NonoverlappingSchwarzOperator<M,X,Y,C>;
1232 if(smoother ==
"ssor")
1233 return std::make_shared<Amg::AMG<OP, X, NonoverlappingBlockPreconditioner<C,SeqSSOR<M,X,Y>>,C>>(op, config, op->getCommunication());
1234 if(smoother ==
"sor")
1235 return std::make_shared<Amg::AMG<OP, X, NonoverlappingBlockPreconditioner<C,SeqSOR<M,X,Y>>,C>>(op, config, op->getCommunication());
1236 if(smoother ==
"jac")
1237 return std::make_shared<Amg::AMG<OP, X, NonoverlappingBlockPreconditioner<C,SeqJac<M,X,Y>>,C>>(op, config, op->getCommunication());
1238 if(smoother ==
"gs")
1239 return std::make_shared<Amg::AMG<OP, X, NonoverlappingBlockPreconditioner<C,SeqGS<M,X,Y>>,C>>(op, config, op->getCommunication());
1240 if(smoother ==
"ilu")
1241 return std::make_shared<Amg::AMG<OP, X, NonoverlappingBlockPreconditioner<C,SeqILU<M,X,Y>>,C>>(op, config, op->getCommunication());
1246 template<
typename OpTraits,
typename OP>
1248 typename OpTraits::range_type>>
1249 operator() (OpTraits opTraits,
const std::shared_ptr<OP>& op,
const Dune::ParameterTree& config,
1250 std::enable_if_t<isValidMatrix<typename OpTraits::matrix_type>::value,
int> = 0)
const
1252 using field_type =
typename OpTraits::matrix_type::field_type;
1253 using real_type =
typename FieldTraits<field_type>::real_type;
1254 if (!std::is_convertible<field_type, real_type>())
1255 DUNE_THROW(UnsupportedType,
"AMG needs field_type(" <<
1256 className<field_type>() <<
1257 ") to be convertible to its real_type (" <<
1258 className<real_type>() <<
1260 std::string smoother = config.
get(
"smoother",
"ssor");
1266 return makeAMG(op, smoother, config);
1269 template<
typename OpTraits,
typename OP>
1271 typename OpTraits::range_type>>
1272 operator() (OpTraits opTraits,
const std::shared_ptr<OP>& op,
const Dune::ParameterTree& config,
1273 std::enable_if_t<!isValidMatrix<typename OpTraits::matrix_type>::value,
int> = 0)
const
1275 DUNE_THROW(UnsupportedType,
"AMG needs a FieldMatrix as Matrix block_type");
1279 DUNE_REGISTER_PRECONDITIONER(
"amg", AMGCreator());
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:66
The criterion describing the stop criteria for the coarsening process.
Definition: matrixhierarchy.hh:283
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:220
LevelIterator< const Hierarchy< MatrixOperator, Allocator >, const MatrixOperator > ConstIterator
Type of the const iterator.
Definition: hierarchy.hh:223
an algebraic multigrid method using a Krylov-cycle.
Definition: kamg.hh:140
Two grid operator for AMG with Krylov cycle.
Definition: kamg.hh:33
The hierarchies build by the coarsening process.
Definition: matrixhierarchy.hh:61
All parameters for AMG.
Definition: parameters.hh:416
Criterion taking advantage of symmetric matrices.
Definition: aggregates.hh:519
Criterion suitable for unsymmetric matrices.
Definition: aggregates.hh:539
Base class for Dune-Exceptions.
Definition: exceptions.hh:96
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:281
ConstIterator class for sequential access.
Definition: matrix.hh:404
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
T block_type
Export the type representing the components.
Definition: matrix.hh:568
The negation of a set. An item is contained in the set if and only if it is not contained in the nega...
Definition: enumset.hh:96
Default exception for dummy implementations.
Definition: exceptions.hh:263
Hierarchical structure of string parameters.
Definition: parametertree.hh:37
std::string get(const std::string &key, const std::string &defaultValue) const
get value as string
Definition: parametertree.cc:188
bool hasKey(const std::string &key) const
test for key
Definition: parametertree.cc:51
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:33
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioner.hh:40
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:52
A few common exception classes.
A hierarchical structure of string parameters.
Define base class for scalar product and norm.
Classes for using UMFPack with ISTL matrices.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr GeometryType none(unsigned int dim)
Returns a GeometryType representing a singular of dimension dim.
Definition: type.hh:471
AMG(const AMG &amg)
Copy constructor.
Definition: amg.hh:393
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: amg.hh:802
Hierarchy< Domain, A >::Iterator update
The iterator over the updates.
Definition: amg.hh:304
Hierarchy< Range, A >::Iterator rhs
The iterator over the right hand sided.
Definition: amg.hh:308
bool usesDirectCoarseLevelSolver() const
Check whether the coarse solver used is a direct solver.
Definition: amg.hh:1025
X Domain
The domain type.
Definition: amg.hh:88
AMG(OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, const Parameters &parms)
Construct a new amg with a specific coarse solver.
Definition: amg.hh:407
AMG(std::shared_ptr< const Operator > fineOperator, const ParameterTree &configuration, const ParallelInformation &pinfo=ParallelInformation())
Constructor an AMG via ParameterTree.
Definition: amg.hh:453
ParallelInformationHierarchy::Iterator pinfo
The iterator over the parallel information.
Definition: amg.hh:288
OperatorHierarchy::AggregatesMapList::const_iterator aggregates
The iterator over the aggregates maps.
Definition: amg.hh:296
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:101
S Smoother
The type of the smoother.
Definition: amg.hh:98
Hierarchy< Smoother, A >::Iterator smoother
The iterator over the smoothers.
Definition: amg.hh:280
M Operator
The matrix operator type.
Definition: amg.hh:74
OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator matrix
The iterator over the matrices.
Definition: amg.hh:284
OperatorHierarchy::RedistributeInfoList::const_iterator redist
The iterator over the redistribution information.
Definition: amg.hh:292
X Range
The range type.
Definition: amg.hh:90
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:406
void getCoarsestAggregateNumbers(std::vector< std::size_t, A1 > &cont)
Get the aggregate number of each unknown on the coarsest level.
Definition: amg.hh:1163
Hierarchy< Domain, A >::Iterator lhs
The iterator over the left hand side.
Definition: amg.hh:300
const void * Arguments
A type holding all the arguments needed to call the constructor.
Definition: construction.hh:44
static std::shared_ptr< T > construct(Arguments &args)
Construct an object with the specified arguments.
Definition: construction.hh:52
void recalculateHierarchy()
Recalculate the matrix hierarchy.
Definition: amg.hh:222
Iterator coarsest()
Get an iterator positioned at the coarsest level.
Definition: hierarchy.hh:387
OperatorHierarchy::ParallelInformationHierarchy ParallelInformationHierarchy
The parallal data distribution hierarchy type.
Definition: amg.hh:85
InverseOperator< X, X > CoarseSolver
the type of the coarse solver.
Definition: amg.hh:92
void post(Domain &x)
Clean up.
Definition: amg.hh:1140
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:428
std::size_t level
The level index.
Definition: amg.hh:312
AMG(const Operator &fineOperator, const C &criterion, const SmootherArgs &smootherArgs=SmootherArgs(), const ParallelInformation &pinfo=ParallelInformation())
Construct an AMG with an inexact coarse solver based on the smoother.
Definition: amg.hh:429
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: amg.hh:889
MatrixHierarchy< M, ParallelInformation, A > OperatorHierarchy
The operator hierarchy type.
Definition: amg.hh:83
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: amg.hh:195
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:81
@ atOnceAccu
Accumulate data to one process at once.
Definition: parameters.hh:243
@ noAccu
No data accumulution.
Definition: parameters.hh:237
@ successiveAccu
Successively accumulate to fewer processes.
Definition: parameters.hh:247
Provides a classes representing the hierarchies in AMG.
Dune namespace.
Definition: alignedallocator.hh:13
std::shared_ptr< T > stackobject_to_shared_ptr(T &t)
Create a shared_ptr for a stack-allocated object.
Definition: shared_ptr.hh:72
std::string className()
Provide the demangled class name of a type T as a string.
Definition: classname.hh:47
Implements a scalar matrix view wrapper around an existing scalar.
Implements a scalar vector view wrapper around an existing scalar.
Classes for the generic construction and application of the smoothers.
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Functor using the row sum (infinity) norm to determine strong couplings.
Definition: aggregates.hh:463
Traits class for getting the attribute class of a smoother.
Definition: smoother.hh:66
Statistics about the application of an inverse operator.
Definition: solver.hh:50
bool converged
True if convergence criterion has been met.
Definition: solver.hh:75
Categories for the solvers.
Definition: solvercategory.hh:22
Category
Definition: solvercategory.hh:23
static Category category(const OP &op, decltype(op.category()) *=nullptr)
Helperfunction to extract the solver category either from an enum, or from the newly introduced virtu...
Definition: solvercategory.hh:34
Classes for using SuperLU with ISTL matrices.
Prolongation and restriction for amg.
Traits for type conversions and type information.