3 #ifndef DUNE_AMG_AMG_HH
4 #define DUNE_AMG_AMG_HH
7 #include <dune/common/exceptions.hh>
16 #include <dune/common/typetraits.hh>
17 #include <dune/common/exceptions.hh>
40 template<
class M,
class X,
class S,
class P,
class K,
class A>
57 class A=std::allocator<X> >
60 template<
class M1,
class X1,
class S1,
class P1,
class K1,
class A1>
120 AMG(
const Operator& fineOperator,
const C& criterion,
167 matrices_->recalculateGalerkin(NegateSet<typename PI::OwnerSet>());
184 void createHierarchies(C& criterion,
Operator& matrix,
210 typename OperatorHierarchy::RedistributeInfoList::const_iterator
redist;
214 typename OperatorHierarchy::AggregatesMapList::const_iterator
aggregates;
238 void mgc(LevelContext& levelContext);
248 void moveToFineLevel(LevelContext& levelContext,
bool processedFineLevel);
254 bool moveToCoarseLevel(LevelContext& levelContext);
260 void initIteratorsWithFineLevel(LevelContext& levelContext);
263 std::shared_ptr<OperatorHierarchy> matrices_;
267 std::shared_ptr<Hierarchy<Smoother,A> > smoothers_;
269 std::shared_ptr<CoarseSolver> solver_;
279 std::shared_ptr<ScalarProduct> scalarProduct_;
283 std::size_t preSteps_;
285 std::size_t postSteps_;
286 bool buildHierarchy_;
288 bool coarsesolverconverged;
289 std::shared_ptr<Smoother> coarseSmoother_;
293 std::size_t verbosity_;
296 template<
class M,
class X,
class S,
class PI,
class A>
298 : matrices_(amg.matrices_), smootherArgs_(amg.smootherArgs_),
299 smoothers_(amg.smoothers_), solver_(amg.solver_),
300 rhs_(), lhs_(), update_(),
301 scalarProduct_(amg.scalarProduct_), gamma_(amg.gamma_),
302 preSteps_(amg.preSteps_), postSteps_(amg.postSteps_),
303 buildHierarchy_(amg.buildHierarchy_),
304 additive(amg.additive), coarsesolverconverged(amg.coarsesolverconverged),
305 coarseSmoother_(amg.coarseSmoother_),
306 category_(amg.category_),
307 verbosity_(amg.verbosity_)
317 template<
class M,
class X,
class S,
class PI,
class A>
321 : matrices_(stackobject_to_shared_ptr(matrices)), smootherArgs_(smootherArgs),
323 rhs_(), lhs_(), update_(), scalarProduct_(0),
324 gamma_(parms.getGamma()), preSteps_(parms.getNoPreSmoothSteps()),
325 postSteps_(parms.getNoPostSmoothSteps()), buildHierarchy_(false),
326 additive(parms.getAdditive()), coarsesolverconverged(true),
330 verbosity_(parms.debugLevel())
332 assert(matrices_->isBuilt());
335 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
338 template<
class M,
class X,
class S,
class PI,
class A>
344 : smootherArgs_(smootherArgs),
346 rhs_(), lhs_(), update_(), scalarProduct_(),
347 gamma_(criterion.getGamma()), preSteps_(criterion.getNoPreSmoothSteps()),
348 postSteps_(criterion.getNoPostSmoothSteps()), buildHierarchy_(true),
349 additive(criterion.getAdditive()), coarsesolverconverged(true),
352 verbosity_(criterion.debugLevel())
354 if(SolverCategory::category(matrix) != SolverCategory::category(pinfo))
359 createHierarchies(criterion,
const_cast<Operator&
>(matrix), pinfo);
363 template<
class M,
class X,
class S,
class PI,
class A>
366 if(buildHierarchy_) {
370 coarseSmoother_.reset();
391 #if DISABLE_AMG_DIRECTSOLVER
393 #elif HAVE_SUITESPARSE_UMFPACK
401 template <
class M, SolverType>
407 DUNE_THROW(NotImplemented,
"DirectSolver not selected");
410 static std::string
name () {
return "None"; }
412 #if HAVE_SUITESPARSE_UMFPACK
414 struct Solver< M, umfpack >
417 static type* create(
const M&
mat,
bool verbose,
bool reusevector )
419 return new type(
mat, verbose, reusevector );
421 static std::string name () {
return "UMFPack"; }
431 return new type(
mat, verbose, reusevector );
433 static std::string
name () {
return "SuperLU"; }
440 static constexpr
bool isDirectSolver = solver != none;
441 static std::string
name() {
return SelectedSolver :: name (); }
444 return SelectedSolver :: create(
mat, verbose, reusevector );
448 template<
class M,
class X,
class S,
class PI,
class A>
450 void AMG<M,X,S,PI,A>::createHierarchies(C& criterion, Operator& matrix,
454 matrices_.reset(
new OperatorHierarchy(matrix, pinfo));
456 matrices_->template build<NegateSet<typename PI::OwnerSet> >(criterion);
459 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
465 if(buildHierarchy_ && matrices_->levels()==matrices_->maxlevels()
466 && ( ! matrices_->redistributeInformation().back().isSetup() ||
467 matrices_->parallelInformation().coarsest().getRedistributed().communicator().size() ) )
470 SmootherArgs sargs(smootherArgs_);
471 sargs.iterations = 1;
473 typename ConstructionTraits<Smoother>::Arguments cargs;
474 cargs.setArgs(sargs);
475 if(matrices_->redistributeInformation().back().isSetup()) {
477 cargs.setMatrix(matrices_->matrices().coarsest().getRedistributed().getmat());
478 cargs.setComm(matrices_->parallelInformation().coarsest().getRedistributed());
480 cargs.setMatrix(matrices_->matrices().coarsest()->getmat());
481 cargs.setComm(*matrices_->parallelInformation().coarsest());
484 coarseSmoother_.reset(ConstructionTraits<Smoother>::construct(cargs));
485 scalarProduct_ = createScalarProduct<X>(cargs.getComm(),category());
487 typedef DirectSolverSelector< typename M::matrix_type, X > SolverSelector;
490 if( SolverSelector::isDirectSolver &&
491 (std::is_same<ParallelInformation,SequentialInformation>::value
492 || matrices_->parallelInformation().coarsest()->communicator().size()==1
493 || (matrices_->parallelInformation().coarsest().isRedistributed()
494 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()==1
495 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()>0) )
498 if(matrices_->parallelInformation().coarsest().isRedistributed())
500 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0)
503 solver_.reset(SolverSelector::create(matrices_->matrices().coarsest().getRedistributed().getmat(),
false,
false));
510 solver_.reset(SolverSelector::create(matrices_->matrices().coarsest()->getmat(),
false,
false));
512 if(verbosity_>0 && matrices_->parallelInformation().coarsest()->communicator().rank()==0)
513 std::cout<<
"Using a direct coarse solver (" << SolverSelector::name() <<
")" << std::endl;
517 if(matrices_->parallelInformation().coarsest().isRedistributed())
519 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0)
524 solver_.reset(
new BiCGSTABSolver<X>(
const_cast<M&
>(matrices_->matrices().coarsest().getRedistributed()),
526 *coarseSmoother_, 1E-2, 1000, 0));
531 solver_.reset(
new BiCGSTABSolver<X>(
const_cast<M&
>(*matrices_->matrices().coarsest()),
533 *coarseSmoother_, 1E-2, 1000, 0));
552 if(verbosity_>0 && matrices_->parallelInformation().finest()->communicator().rank()==0)
553 std::cout<<
"Building hierarchy of "<<matrices_->maxlevels()<<
" levels "
554 <<
"(inclusive coarse solver) took "<<watch.elapsed()<<
" seconds."<<std::endl;
558 template<
class M,
class X,
class S,
class PI,
class A>
565 typedef typename M::matrix_type
Matrix;
572 const Matrix&
mat=matrices_->matrices().finest()->getmat();
573 for(RowIter row=
mat.begin(); row!=
mat.end(); ++row) {
574 bool isDirichlet =
true;
575 bool hasDiagonal =
false;
577 for(ColIter
col=row->begin();
col!=row->end(); ++
col) {
578 if(row.index()==
col.index()) {
586 if(isDirichlet && hasDiagonal)
587 diagonal.solve(x[row.index()], b[row.index()]);
590 if(smoothers_->levels()>0)
591 smoothers_->finest()->pre(x,b);
594 matrices_->parallelInformation().coarsest()->copyOwnerToAll(x,x);
607 matrices_->coarsenVector(*rhs_);
608 matrices_->coarsenVector(*lhs_);
609 matrices_->coarsenVector(*update_);
615 Iterator coarsest = smoothers_->
coarsest();
616 Iterator smoother = smoothers_->finest();
617 RIterator rhs = rhs_->finest();
618 DIterator lhs = lhs_->finest();
619 if(smoothers_->levels()>0) {
621 assert(lhs_->levels()==rhs_->levels());
622 assert(smoothers_->levels()==lhs_->levels() || matrices_->levels()==matrices_->maxlevels());
623 assert(smoothers_->levels()+1==lhs_->levels() || matrices_->levels()<matrices_->maxlevels());
625 if(smoother!=coarsest)
626 for(++smoother, ++lhs, ++rhs; smoother != coarsest; ++smoother, ++lhs, ++rhs)
627 smoother->pre(*lhs,*rhs);
628 smoother->pre(*lhs,*rhs);
638 template<
class M,
class X,
class S,
class PI,
class A>
641 return matrices_->
levels();
643 template<
class M,
class X,
class S,
class PI,
class A>
650 template<
class M,
class X,
class S,
class PI,
class A>
653 LevelContext levelContext;
661 initIteratorsWithFineLevel(levelContext);
664 *levelContext.lhs = v;
665 *levelContext.rhs = d;
666 *levelContext.update=0;
667 levelContext.level=0;
671 if(postSteps_==0||matrices_->maxlevels()==1)
672 levelContext.pinfo->copyOwnerToAll(*levelContext.update, *levelContext.update);
674 v=*levelContext.update;
679 template<
class M,
class X,
class S,
class PI,
class A>
682 levelContext.smoother = smoothers_->finest();
683 levelContext.matrix = matrices_->matrices().finest();
684 levelContext.pinfo = matrices_->parallelInformation().finest();
685 levelContext.redist =
686 matrices_->redistributeInformation().begin();
687 levelContext.aggregates = matrices_->aggregatesMaps().begin();
688 levelContext.lhs = lhs_->finest();
689 levelContext.update = update_->finest();
690 levelContext.rhs = rhs_->finest();
693 template<
class M,
class X,
class S,
class PI,
class A>
695 ::moveToCoarseLevel(LevelContext& levelContext)
698 bool processNextLevel=
true;
700 if(levelContext.redist->isSetup()) {
701 levelContext.redist->redistribute(
static_cast<const Range&
>(*levelContext.rhs),
702 levelContext.rhs.getRedistributed());
703 processNextLevel = levelContext.rhs.getRedistributed().size()>0;
704 if(processNextLevel) {
707 ++levelContext.pinfo;
708 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
709 ::restrictVector(*(*levelContext.aggregates), *levelContext.rhs,
710 static_cast<const Range&
>(fineRhs.getRedistributed()),
711 *levelContext.pinfo);
716 ++levelContext.pinfo;
717 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
718 ::restrictVector(*(*levelContext.aggregates),
719 *levelContext.rhs,
static_cast<const Range&
>(*fineRhs),
720 *levelContext.pinfo);
723 if(processNextLevel) {
726 ++levelContext.update;
727 ++levelContext.matrix;
728 ++levelContext.level;
729 ++levelContext.redist;
731 if(levelContext.matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()) {
733 ++levelContext.smoother;
734 ++levelContext.aggregates;
737 *levelContext.update=0;
739 return processNextLevel;
742 template<
class M,
class X,
class S,
class PI,
class A>
744 ::moveToFineLevel(LevelContext& levelContext,
bool processNextLevel)
746 if(processNextLevel) {
747 if(levelContext.matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()) {
749 --levelContext.smoother;
750 --levelContext.aggregates;
752 --levelContext.redist;
753 --levelContext.level;
755 --levelContext.matrix;
759 --levelContext.pinfo;
761 if(levelContext.redist->isSetup()) {
763 levelContext.lhs.getRedistributed()=0;
764 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
765 ::prolongateVector(*(*levelContext.aggregates), *levelContext.update, *levelContext.lhs,
766 levelContext.lhs.getRedistributed(),
767 matrices_->getProlongationDampingFactor(),
768 *levelContext.pinfo, *levelContext.redist);
771 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
772 ::prolongateVector(*(*levelContext.aggregates), *levelContext.update, *levelContext.lhs,
773 matrices_->getProlongationDampingFactor(),
774 *levelContext.pinfo);
778 if(processNextLevel) {
779 --levelContext.update;
783 *levelContext.update += *levelContext.lhs;
786 template<
class M,
class X,
class S,
class PI,
class A>
792 template<
class M,
class X,
class S,
class PI,
class A>
794 if(levelContext.matrix == matrices_->matrices().coarsest() && levels()==maxlevels()) {
798 if(levelContext.redist->isSetup()) {
799 levelContext.redist->redistribute(*levelContext.rhs, levelContext.rhs.getRedistributed());
800 if(levelContext.rhs.getRedistributed().size()>0) {
802 levelContext.pinfo.getRedistributed().copyOwnerToAll(levelContext.rhs.getRedistributed(),
803 levelContext.rhs.getRedistributed());
804 solver_->apply(levelContext.update.getRedistributed(),
805 levelContext.rhs.getRedistributed(), res);
807 levelContext.redist->redistributeBackward(*levelContext.update, levelContext.update.getRedistributed());
808 levelContext.pinfo->copyOwnerToAll(*levelContext.update, *levelContext.update);
810 levelContext.pinfo->copyOwnerToAll(*levelContext.rhs, *levelContext.rhs);
811 solver_->apply(*levelContext.update, *levelContext.rhs, res);
815 coarsesolverconverged =
false;
820 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
821 bool processNextLevel = moveToCoarseLevel(levelContext);
823 if(processNextLevel) {
825 for(std::size_t i=0; i<gamma_; i++)
829 moveToFineLevel(levelContext, processNextLevel);
834 if(levelContext.matrix == matrices_->matrices().finest()) {
835 coarsesolverconverged = matrices_->parallelInformation().finest()->communicator().prod(coarsesolverconverged);
836 if(!coarsesolverconverged)
837 DUNE_THROW(MathError,
"Coarse solver did not converge");
845 template<
class M,
class X,
class S,
class PI,
class A>
846 void AMG<M,X,S,PI,A>::additiveMgc(){
849 typename ParallelInformationHierarchy::Iterator pinfo=matrices_->parallelInformation().finest();
852 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates=matrices_->aggregatesMaps().begin();
856 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
857 ::restrictVector(*(*aggregates), *rhs,
static_cast<const Range&
>(*fineRhs), *pinfo);
863 lhs = lhs_->finest();
866 for(rhs=rhs_->finest(); rhs != rhs_->coarsest(); ++lhs, ++rhs, ++smoother) {
869 smoother->apply(*lhs, *rhs);
873 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
874 InverseOperatorResult res;
875 pinfo->copyOwnerToAll(*rhs, *rhs);
876 solver_->apply(*lhs, *rhs, res);
879 DUNE_THROW(MathError,
"Coarse solver did not converge");
888 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
889 ::prolongateVector(*(*aggregates), *coarseLhs, *lhs, 1.0, *pinfo);
895 template<
class M,
class X,
class S,
class PI,
class A>
898 DUNE_UNUSED_PARAMETER(x);
902 Iterator coarsest = smoothers_->
coarsest();
903 Iterator smoother = smoothers_->finest();
904 DIterator lhs = lhs_->finest();
905 if(smoothers_->levels()>0) {
906 if(smoother != coarsest || matrices_->levels()<matrices_->maxlevels())
907 smoother->post(*lhs);
908 if(smoother!=coarsest)
909 for(++smoother, ++lhs; smoother != coarsest; ++smoother, ++lhs)
910 smoother->post(*lhs);
911 smoother->post(*lhs);
924 template<
class M,
class X,
class S,
class PI,
class A>
928 matrices_->getCoarsestAggregatesOnFinest(cont);