Dune Core Modules (2.3.1)

amg.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_AMG_AMG_HH
5#define DUNE_AMG_AMG_HH
6
7#include <memory>
12#include <dune/istl/solvers.hh>
14#include <dune/istl/superlu.hh>
15#include <dune/istl/umfpack.hh>
19
20namespace Dune
21{
22 namespace Amg
23 {
41 template<class M, class X, class S, class P, class K, class A>
42 class KAMG;
43
44 template<class T>
45 class KAmgTwoGrid;
46
54 template<class M, class X, class S, class PI=SequentialInformation,
55 class A=std::allocator<X> >
56 class AMG : public Preconditioner<X,X>
57 {
58 template<class M1, class X1, class S1, class P1, class K1, class A1>
59 friend class KAMG;
60
61 friend class KAmgTwoGrid<AMG>;
62
63 public:
65 typedef M Operator;
77
79 typedef X Domain;
81 typedef X Range;
89 typedef S Smoother;
90
93
94 enum {
96 category = S::category
97 };
98
115 AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
116 const SmootherArgs& smootherArgs, std::size_t gamma,
117 std::size_t preSmoothingSteps,
118 std::size_t postSmoothingSteps,
119 bool additive=false) DUNE_DEPRECATED;
120
130 AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
131 const SmootherArgs& smootherArgs, const Parameters& parms);
132
152 template<class C>
153 AMG(const Operator& fineOperator, const C& criterion,
154 const SmootherArgs& smootherArgs, std::size_t gamma,
155 std::size_t preSmoothingSteps,
156 std::size_t postSmoothingSteps,
157 bool additive=false,
159
171 template<class C>
172 AMG(const Operator& fineOperator, const C& criterion,
173 const SmootherArgs& smootherArgs=SmootherArgs(),
175
179 AMG(const AMG& amg);
180
181 ~AMG();
182
184 void pre(Domain& x, Range& b);
185
187 void apply(Domain& v, const Range& d);
188
190 void post(Domain& x);
191
196 template<class A1>
197 void getCoarsestAggregateNumbers(std::vector<std::size_t,A1>& cont);
198
199 std::size_t levels();
200
201 std::size_t maxlevels();
202
212 {
213 matrices_->recalculateGalerkin(NegateSet<typename PI::OwnerSet>());
214 }
215
221
222 private:
229 template<class C>
230 void createHierarchies(C& criterion, Operator& matrix,
231 const PI& pinfo);
238 struct LevelContext
239 {
240 typedef Smoother SmootherType;
256 typename OperatorHierarchy::RedistributeInfoList::const_iterator redist;
260 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates;
276 std::size_t level;
277 };
278
279
284 void mgc(LevelContext& levelContext);
285
286 void additiveMgc();
287
294 void moveToFineLevel(LevelContext& levelContext,bool processedFineLevel);
295
300 bool moveToCoarseLevel(LevelContext& levelContext);
301
306 void initIteratorsWithFineLevel(LevelContext& levelContext);
307
311 SmootherArgs smootherArgs_;
317 Hierarchy<Range,A>* rhs_;
321 Hierarchy<Domain,A>* update_;
325 typedef typename ScalarProductChooser::ScalarProduct ScalarProduct;
328 ScalarProductPointer scalarProduct_;
330 std::size_t gamma_;
332 std::size_t preSteps_;
334 std::size_t postSteps_;
335 bool buildHierarchy_;
336 bool additive;
337 bool coarsesolverconverged;
338 Dune::shared_ptr<Smoother> coarseSmoother_;
340 std::size_t verbosity_;
341 };
342
343 template<class M, class X, class S, class PI, class A>
344 inline AMG<M,X,S,PI,A>::AMG(const AMG& amg)
345 : matrices_(amg.matrices_), smootherArgs_(amg.smootherArgs_),
346 smoothers_(amg.smoothers_), solver_(amg.solver_),
347 rhs_(), lhs_(), update_(),
348 scalarProduct_(amg.scalarProduct_), gamma_(amg.gamma_),
349 preSteps_(amg.preSteps_), postSteps_(amg.postSteps_),
350 buildHierarchy_(amg.buildHierarchy_),
351 additive(amg.additive), coarsesolverconverged(amg.coarsesolverconverged),
352 coarseSmoother_(amg.coarseSmoother_), verbosity_(amg.verbosity_)
353 {
354 if(amg.rhs_)
355 rhs_=new Hierarchy<Range,A>(*amg.rhs_);
356 if(amg.lhs_)
357 lhs_=new Hierarchy<Domain,A>(*amg.lhs_);
358 if(amg.update_)
359 update_=new Hierarchy<Domain,A>(*amg.update_);
360 }
361
362 template<class M, class X, class S, class PI, class A>
364 const SmootherArgs& smootherArgs,
365 std::size_t gamma, std::size_t preSmoothingSteps,
366 std::size_t postSmoothingSteps, bool additive_)
367 : matrices_(&matrices), smootherArgs_(smootherArgs),
368 smoothers_(new Hierarchy<Smoother,A>), solver_(&coarseSolver),
369 rhs_(), lhs_(), update_(), scalarProduct_(0),
370 gamma_(gamma), preSteps_(preSmoothingSteps), postSteps_(postSmoothingSteps), buildHierarchy_(false),
371 additive(additive_), coarsesolverconverged(true),
372 coarseSmoother_(), verbosity_(2)
373 {
374 assert(matrices_->isBuilt());
375
376 // build the necessary smoother hierarchies
377 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
378 }
379
380 template<class M, class X, class S, class PI, class A>
382 const SmootherArgs& smootherArgs,
383 const Parameters& parms)
384 : matrices_(&matrices), smootherArgs_(smootherArgs),
385 smoothers_(new Hierarchy<Smoother,A>), solver_(&coarseSolver),
386 rhs_(), lhs_(), update_(), scalarProduct_(0),
387 gamma_(parms.getGamma()), preSteps_(parms.getNoPreSmoothSteps()),
388 postSteps_(parms.getNoPostSmoothSteps()), buildHierarchy_(false),
389 additive(parms.getAdditive()), coarsesolverconverged(true),
390 coarseSmoother_(), verbosity_(parms.debugLevel())
391 {
392 assert(matrices_->isBuilt());
393
394 // build the necessary smoother hierarchies
395 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
396 }
397
398 template<class M, class X, class S, class PI, class A>
399 template<class C>
401 const C& criterion,
402 const SmootherArgs& smootherArgs,
403 std::size_t gamma, std::size_t preSmoothingSteps,
404 std::size_t postSmoothingSteps,
405 bool additive_,
406 const PI& pinfo)
407 : smootherArgs_(smootherArgs),
408 smoothers_(new Hierarchy<Smoother,A>), solver_(),
409 rhs_(), lhs_(), update_(), scalarProduct_(), gamma_(gamma),
410 preSteps_(preSmoothingSteps), postSteps_(postSmoothingSteps), buildHierarchy_(true),
411 additive(additive_), coarsesolverconverged(true),
412 coarseSmoother_(), verbosity_(criterion.debugLevel())
413 {
414 dune_static_assert(static_cast<int>(M::category)==static_cast<int>(S::category),
415 "Matrix and Solver must match in terms of category!");
416 // TODO: reestablish compile time checks.
417 //dune_static_assert(static_cast<int>(PI::category)==static_cast<int>(S::category),
418 // "Matrix and Solver must match in terms of category!");
419 createHierarchies(criterion, const_cast<Operator&>(matrix), pinfo);
420 }
421
422 template<class M, class X, class S, class PI, class A>
423 template<class C>
425 const C& criterion,
426 const SmootherArgs& smootherArgs,
427 const PI& pinfo)
428 : smootherArgs_(smootherArgs),
429 smoothers_(new Hierarchy<Smoother,A>), solver_(),
430 rhs_(), lhs_(), update_(), scalarProduct_(),
431 gamma_(criterion.getGamma()), preSteps_(criterion.getNoPreSmoothSteps()),
432 postSteps_(criterion.getNoPostSmoothSteps()), buildHierarchy_(true),
433 additive(criterion.getAdditive()), coarsesolverconverged(true),
434 coarseSmoother_(), verbosity_(criterion.debugLevel())
435 {
436 dune_static_assert(static_cast<int>(M::category)==static_cast<int>(S::category),
437 "Matrix and Solver must match in terms of category!");
438 // TODO: reestablish compile time checks.
439 //dune_static_assert(static_cast<int>(PI::category)==static_cast<int>(S::category),
440 // "Matrix and Solver must match in terms of category!");
441 createHierarchies(criterion, const_cast<Operator&>(matrix), pinfo);
442 }
443
444
445 template<class M, class X, class S, class PI, class A>
447 {
448 if(buildHierarchy_) {
449 if(solver_)
450 solver_.reset();
451 if(coarseSmoother_)
452 coarseSmoother_.reset();
453 }
454 if(lhs_)
455 delete lhs_;
456 lhs_=nullptr;
457 if(update_)
458 delete update_;
459 update_=nullptr;
460 if(rhs_)
461 delete rhs_;
462 rhs_=nullptr;
463 }
464
465 template<class M, class X, class S, class PI, class A>
466 template<class C>
467 void AMG<M,X,S,PI,A>::createHierarchies(C& criterion, Operator& matrix,
468 const PI& pinfo)
469 {
470 Timer watch;
471 matrices_.reset(new OperatorHierarchy(matrix, pinfo));
472
473 matrices_->template build<NegateSet<typename PI::OwnerSet> >(criterion);
474
475 // build the necessary smoother hierarchies
476 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
477
478 if(buildHierarchy_ && matrices_->levels()==matrices_->maxlevels()) {
479 // We have the carsest level. Create the coarse Solver
480 SmootherArgs sargs(smootherArgs_);
481 sargs.iterations = 1;
482
484 cargs.setArgs(sargs);
485 if(matrices_->redistributeInformation().back().isSetup()) {
486 // Solve on the redistributed partitioning
487 cargs.setMatrix(matrices_->matrices().coarsest().getRedistributed().getmat());
488 cargs.setComm(matrices_->parallelInformation().coarsest().getRedistributed());
489 }else{
490 cargs.setMatrix(matrices_->matrices().coarsest()->getmat());
491 cargs.setComm(*matrices_->parallelInformation().coarsest());
492 }
493
494 coarseSmoother_.reset(ConstructionTraits<Smoother>::construct(cargs));
495 scalarProduct_.reset(ScalarProductChooser::construct(cargs.getComm()));
496
497#if HAVE_SUPERLU || HAVE_UMFPACK
498#if HAVE_UMFPACK
499#define DIRECTSOLVER UMFPack
500#else
501#define DIRECTSOLVER SuperLU
502#endif
503 // Use superlu if we are purely sequential or with only one processor on the coarsest level.
504 if(is_same<ParallelInformation,SequentialInformation>::value // sequential mode
505 || matrices_->parallelInformation().coarsest()->communicator().size()==1 //parallel mode and only one processor
506 || (matrices_->parallelInformation().coarsest().isRedistributed()
507 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()==1
508 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()>0)) { // redistribute and 1 proc
509 if(matrices_->parallelInformation().coarsest().isRedistributed())
510 {
511 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0)
512 // We are still participating on this level
513 solver_.reset(new DIRECTSOLVER<typename M::matrix_type>(matrices_->matrices().coarsest().getRedistributed().getmat(), false, false));
514 else
515 solver_.reset();
516 }else
517 solver_.reset(new DIRECTSOLVER<typename M::matrix_type>(matrices_->matrices().coarsest()->getmat(), false, false));
518 if(verbosity_>0 && matrices_->parallelInformation().coarsest()->communicator().rank()==0)
519 std::cout<< "Using a direct coarse solver (" << static_cast< DIRECTSOLVER<typename M::matrix_type>* >(solver_.get())->name() << ")" << std::endl;
520 }else
521#undef DIRECTSOLVER
522#endif
523 {
524 if(matrices_->parallelInformation().coarsest().isRedistributed())
525 {
526 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0)
527 // We are still participating on this level
528 solver_.reset(new BiCGSTABSolver<X>(const_cast<M&>(matrices_->matrices().coarsest().getRedistributed()),
529 *scalarProduct_,
530 *coarseSmoother_, 1E-2, 1000, 0));
531 else
532 solver_.reset();
533 }else
534 solver_.reset(new BiCGSTABSolver<X>(const_cast<M&>(*matrices_->matrices().coarsest()),
535 *scalarProduct_,
536 *coarseSmoother_, 1E-2, 1000, 0));
537 }
538 }
539
540 if(verbosity_>0 && matrices_->parallelInformation().finest()->communicator().rank()==0)
541 std::cout<<"Building hierarchy of "<<matrices_->maxlevels()<<" levels "
542 <<"(inclusive coarse solver) took "<<watch.elapsed()<<" seconds."<<std::endl;
543 }
544
545
546 template<class M, class X, class S, class PI, class A>
548 {
549 // Detect Matrix rows where all offdiagonal entries are
550 // zero and set x such that A_dd*x_d=b_d
551 // Thus users can be more careless when setting up their linear
552 // systems.
553 typedef typename M::matrix_type Matrix;
554 typedef typename Matrix::ConstRowIterator RowIter;
555 typedef typename Matrix::ConstColIterator ColIter;
556 typedef typename Matrix::block_type Block;
557 Block zero;
558 zero=typename Matrix::field_type();
559
560 const Matrix& mat=matrices_->matrices().finest()->getmat();
561 for(RowIter row=mat.begin(); row!=mat.end(); ++row) {
562 bool isDirichlet = true;
563 bool hasDiagonal = false;
564 Block diagonal;
565 for(ColIter col=row->begin(); col!=row->end(); ++col) {
566 if(row.index()==col.index()) {
567 diagonal = *col;
568 hasDiagonal = false;
569 }else{
570 if(*col!=zero)
571 isDirichlet = false;
572 }
573 }
574 if(isDirichlet && hasDiagonal)
575 diagonal.solve(x[row.index()], b[row.index()]);
576 }
577
578 if(smoothers_->levels()>0)
579 smoothers_->finest()->pre(x,b);
580 else
581 // No smoother to make x consistent! Do it by hand
582 matrices_->parallelInformation().coarsest()->copyOwnerToAll(x,x);
583 Range* copy = new Range(b);
584 if(rhs_)
585 delete rhs_;
586 rhs_ = new Hierarchy<Range,A>(copy);
587 Domain* dcopy = new Domain(x);
588 if(lhs_)
589 delete lhs_;
590 lhs_ = new Hierarchy<Domain,A>(dcopy);
591 dcopy = new Domain(x);
592 if(update_)
593 delete update_;
594 update_ = new Hierarchy<Domain,A>(dcopy);
595 matrices_->coarsenVector(*rhs_);
596 matrices_->coarsenVector(*lhs_);
597 matrices_->coarsenVector(*update_);
598
599 // Preprocess all smoothers
600 typedef typename Hierarchy<Smoother,A>::Iterator Iterator;
601 typedef typename Hierarchy<Range,A>::Iterator RIterator;
602 typedef typename Hierarchy<Domain,A>::Iterator DIterator;
603 Iterator coarsest = smoothers_->coarsest();
604 Iterator smoother = smoothers_->finest();
605 RIterator rhs = rhs_->finest();
606 DIterator lhs = lhs_->finest();
607 if(smoothers_->levels()>0) {
608
609 assert(lhs_->levels()==rhs_->levels());
610 assert(smoothers_->levels()==lhs_->levels() || matrices_->levels()==matrices_->maxlevels());
611 assert(smoothers_->levels()+1==lhs_->levels() || matrices_->levels()<matrices_->maxlevels());
612
613 if(smoother!=coarsest)
614 for(++smoother, ++lhs, ++rhs; smoother != coarsest; ++smoother, ++lhs, ++rhs)
615 smoother->pre(*lhs,*rhs);
616 smoother->pre(*lhs,*rhs);
617 }
618
619
620 // The preconditioner might change x and b. So we have to
621 // copy the changes to the original vectors.
622 x = *lhs_->finest();
623 b = *rhs_->finest();
624
625 }
626 template<class M, class X, class S, class PI, class A>
627 std::size_t AMG<M,X,S,PI,A>::levels()
628 {
629 return matrices_->levels();
630 }
631 template<class M, class X, class S, class PI, class A>
632 std::size_t AMG<M,X,S,PI,A>::maxlevels()
633 {
634 return matrices_->maxlevels();
635 }
636
638 template<class M, class X, class S, class PI, class A>
640 {
641 LevelContext levelContext;
642
643 if(additive) {
644 *(rhs_->finest())=d;
645 additiveMgc();
646 v=*lhs_->finest();
647 }else{
648 // Init all iterators for the current level
649 initIteratorsWithFineLevel(levelContext);
650
651
652 *levelContext.lhs = v;
653 *levelContext.rhs = d;
654 *levelContext.update=0;
655 levelContext.level=0;
656
657 mgc(levelContext);
658
659 if(postSteps_==0||matrices_->maxlevels()==1)
660 levelContext.pinfo->copyOwnerToAll(*levelContext.update, *levelContext.update);
661
662 v=*levelContext.update;
663 }
664
665 }
666
667 template<class M, class X, class S, class PI, class A>
668 void AMG<M,X,S,PI,A>::initIteratorsWithFineLevel(LevelContext& levelContext)
669 {
670 levelContext.smoother = smoothers_->finest();
671 levelContext.matrix = matrices_->matrices().finest();
672 levelContext.pinfo = matrices_->parallelInformation().finest();
673 levelContext.redist =
674 matrices_->redistributeInformation().begin();
675 levelContext.aggregates = matrices_->aggregatesMaps().begin();
676 levelContext.lhs = lhs_->finest();
677 levelContext.update = update_->finest();
678 levelContext.rhs = rhs_->finest();
679 }
680
681 template<class M, class X, class S, class PI, class A>
682 bool AMG<M,X,S,PI,A>
683 ::moveToCoarseLevel(LevelContext& levelContext)
684 {
685
686 bool processNextLevel=true;
687
688 if(levelContext.redist->isSetup()) {
689 levelContext.redist->redistribute(static_cast<const Range&>(*levelContext.rhs),
690 levelContext.rhs.getRedistributed());
691 processNextLevel = levelContext.rhs.getRedistributed().size()>0;
692 if(processNextLevel) {
693 //restrict defect to coarse level right hand side.
694 typename Hierarchy<Range,A>::Iterator fineRhs = levelContext.rhs++;
695 ++levelContext.pinfo;
696 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
697 ::restrictVector(*(*levelContext.aggregates), *levelContext.rhs,
698 static_cast<const Range&>(fineRhs.getRedistributed()),
699 *levelContext.pinfo);
700 }
701 }else{
702 //restrict defect to coarse level right hand side.
703 typename Hierarchy<Range,A>::Iterator fineRhs = levelContext.rhs++;
704 ++levelContext.pinfo;
705 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
706 ::restrictVector(*(*levelContext.aggregates),
707 *levelContext.rhs, static_cast<const Range&>(*fineRhs),
708 *levelContext.pinfo);
709 }
710
711 if(processNextLevel) {
712 // prepare coarse system
713 ++levelContext.lhs;
714 ++levelContext.update;
715 ++levelContext.matrix;
716 ++levelContext.level;
717 ++levelContext.redist;
718
719 if(levelContext.matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()) {
720 // next level is not the globally coarsest one
721 ++levelContext.smoother;
722 ++levelContext.aggregates;
723 }
724 // prepare the update on the next level
725 *levelContext.update=0;
726 }
727 return processNextLevel;
728 }
729
730 template<class M, class X, class S, class PI, class A>
731 void AMG<M,X,S,PI,A>
732 ::moveToFineLevel(LevelContext& levelContext, bool processNextLevel)
733 {
734 if(processNextLevel) {
735 if(levelContext.matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()) {
736 // previous level is not the globally coarsest one
737 --levelContext.smoother;
738 --levelContext.aggregates;
739 }
740 --levelContext.redist;
741 --levelContext.level;
742 //prolongate and add the correction (update is in coarse left hand side)
743 --levelContext.matrix;
744
745 //typename Hierarchy<Domain,A>::Iterator coarseLhs = lhs--;
746 --levelContext.lhs;
747 --levelContext.pinfo;
748 }
749 if(levelContext.redist->isSetup()) {
750 // Need to redistribute during prolongateVector
751 levelContext.lhs.getRedistributed()=0;
752 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
753 ::prolongateVector(*(*levelContext.aggregates), *levelContext.update, *levelContext.lhs,
754 levelContext.lhs.getRedistributed(),
755 matrices_->getProlongationDampingFactor(),
756 *levelContext.pinfo, *levelContext.redist);
757 }else{
758 *levelContext.lhs=0;
759 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
760 ::prolongateVector(*(*levelContext.aggregates), *levelContext.update, *levelContext.lhs,
761 matrices_->getProlongationDampingFactor(),
762 *levelContext.pinfo);
763 }
764
765
766 if(processNextLevel) {
767 --levelContext.update;
768 --levelContext.rhs;
769 }
770
771 *levelContext.update += *levelContext.lhs;
772 }
773
774 template<class M, class X, class S, class PI, class A>
776 {
777 return IsDirectSolver< CoarseSolver>::value;
778 }
779
780 template<class M, class X, class S, class PI, class A>
781 void AMG<M,X,S,PI,A>::mgc(LevelContext& levelContext){
782 if(levelContext.matrix == matrices_->matrices().coarsest() && levels()==maxlevels()) {
783 // Solve directly
785 res.converged=true; // If we do not compute this flag will not get updated
786 if(levelContext.redist->isSetup()) {
787 levelContext.redist->redistribute(*levelContext.rhs, levelContext.rhs.getRedistributed());
788 if(levelContext.rhs.getRedistributed().size()>0) {
789 // We are still participating in the computation
790 levelContext.pinfo.getRedistributed().copyOwnerToAll(levelContext.rhs.getRedistributed(),
791 levelContext.rhs.getRedistributed());
792 solver_->apply(levelContext.update.getRedistributed(),
793 levelContext.rhs.getRedistributed(), res);
794 }
795 levelContext.redist->redistributeBackward(*levelContext.update, levelContext.update.getRedistributed());
796 levelContext.pinfo->copyOwnerToAll(*levelContext.update, *levelContext.update);
797 }else{
798 levelContext.pinfo->copyOwnerToAll(*levelContext.rhs, *levelContext.rhs);
799 solver_->apply(*levelContext.update, *levelContext.rhs, res);
800 }
801
802 if (!res.converged)
803 coarsesolverconverged = false;
804 }else{
805 // presmoothing
806 presmooth(levelContext, preSteps_);
807
808#ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
809 bool processNextLevel = moveToCoarseLevel(levelContext);
810
811 if(processNextLevel) {
812 // next level
813 for(std::size_t i=0; i<gamma_; i++)
814 mgc(levelContext);
815 }
816
817 moveToFineLevel(levelContext, processNextLevel);
818#else
819 *lhs=0;
820#endif
821
822 if(levelContext.matrix == matrices_->matrices().finest()) {
823 coarsesolverconverged = matrices_->parallelInformation().finest()->communicator().prod(coarsesolverconverged);
824 if(!coarsesolverconverged)
825 DUNE_THROW(MathError, "Coarse solver did not converge");
826 }
827 // postsmoothing
828 postsmooth(levelContext, postSteps_);
829
830 }
831 }
832
833 template<class M, class X, class S, class PI, class A>
834 void AMG<M,X,S,PI,A>::additiveMgc(){
835
836 // restrict residual to all levels
837 typename ParallelInformationHierarchy::Iterator pinfo=matrices_->parallelInformation().finest();
838 typename Hierarchy<Range,A>::Iterator rhs=rhs_->finest();
839 typename Hierarchy<Domain,A>::Iterator lhs = lhs_->finest();
840 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates=matrices_->aggregatesMaps().begin();
841
842 for(typename Hierarchy<Range,A>::Iterator fineRhs=rhs++; fineRhs != rhs_->coarsest(); fineRhs=rhs++, ++aggregates) {
843 ++pinfo;
844 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
845 ::restrictVector(*(*aggregates), *rhs, static_cast<const Range&>(*fineRhs), *pinfo);
846 }
847
848 // pinfo is invalid, set to coarsest level
849 //pinfo = matrices_->parallelInformation().coarsest
850 // calculate correction for all levels
851 lhs = lhs_->finest();
852 typename Hierarchy<Smoother,A>::Iterator smoother = smoothers_->finest();
853
854 for(rhs=rhs_->finest(); rhs != rhs_->coarsest(); ++lhs, ++rhs, ++smoother) {
855 // presmoothing
856 *lhs=0;
857 smoother->apply(*lhs, *rhs);
858 }
859
860 // Coarse level solve
861#ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
862 InverseOperatorResult res;
863 pinfo->copyOwnerToAll(*rhs, *rhs);
864 solver_->apply(*lhs, *rhs, res);
865
866 if(!res.converged)
867 DUNE_THROW(MathError, "Coarse solver did not converge");
868#else
869 *lhs=0;
870#endif
871 // Prologate and add up corrections from all levels
872 --pinfo;
873 --aggregates;
874
875 for(typename Hierarchy<Domain,A>::Iterator coarseLhs = lhs--; coarseLhs != lhs_->finest(); coarseLhs = lhs--, --aggregates, --pinfo) {
876 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
877 ::prolongateVector(*(*aggregates), *coarseLhs, *lhs, 1.0, *pinfo);
878 }
879 }
880
881
883 template<class M, class X, class S, class PI, class A>
885 {
887 // Postprocess all smoothers
888 typedef typename Hierarchy<Smoother,A>::Iterator Iterator;
889 typedef typename Hierarchy<Domain,A>::Iterator DIterator;
890 Iterator coarsest = smoothers_->coarsest();
891 Iterator smoother = smoothers_->finest();
892 DIterator lhs = lhs_->finest();
893 if(smoothers_->levels()>0) {
894 if(smoother != coarsest || matrices_->levels()<matrices_->maxlevels())
895 smoother->post(*lhs);
896 if(smoother!=coarsest)
897 for(++smoother, ++lhs; smoother != coarsest; ++smoother, ++lhs)
898 smoother->post(*lhs);
899 smoother->post(*lhs);
900 }
901 //delete &(*lhs_->finest());
902 delete lhs_;
903 lhs_=nullptr;
904 //delete &(*update_->finest());
905 delete update_;
906 update_=nullptr;
907 //delete &(*rhs_->finest());
908 delete rhs_;
909 rhs_=nullptr;
910 }
911
912 template<class M, class X, class S, class PI, class A>
913 template<class A1>
914 void AMG<M,X,S,PI,A>::getCoarsestAggregateNumbers(std::vector<std::size_t,A1>& cont)
915 {
916 matrices_->getCoarsestAggregatesOnFinest(cont);
917 }
918
919 } // end namespace Amg
920} // end namespace Dune
921
922#endif
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:57
LevelIterator< Hierarchy< Range, A >, Range > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:258
LevelIterator< const Hierarchy< MatrixOperator, Allocator >, const MatrixOperator > ConstIterator
Type of the const iterator.
Definition: hierarchy.hh:261
an algebraic multigrid method using a Krylov-cycle.
Definition: kamg.hh:134
Two grid operator for AMG with Krylov cycle.
Definition: kamg.hh:31
The hierarchies build by the coarsening process.
Definition: hierarchy.hh:318
All parameters for AMG.
Definition: parameters.hh:391
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
T::field_type field_type
Export the type representing the underlying field.
Definition: matrix.hh:29
row_type::const_iterator ConstColIterator
Const iterator for the entries of each row.
Definition: matrix.hh:53
T block_type
Export the type representing the components.
Definition: matrix.hh:32
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:95
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:26
ConstIterator class for sequential access.
Definition: vbvector.hh:647
A reference counting smart pointer.
Definition: shared_ptr.hh:64
A few common exception classes.
#define dune_static_assert(COND, MSG)
Helper template so that compilation fails if condition is not true.
Definition: static_assert.hh:79
#define DUNE_DEPRECATED
Mark some entity as deprecated.
Definition: deprecated.hh:84
#define DUNE_THROW(E, m)
Definition: exceptions.hh:244
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: amg.hh:547
Hierarchy< Domain, A >::Iterator update
The iterator over the updates.
Definition: amg.hh:268
Hierarchy< Range, A >::Iterator rhs
The iterator over the right hand sided.
Definition: amg.hh:272
bool usesDirectCoarseLevelSolver() const
Check whether the coarse solver used is a direct solver.
Definition: amg.hh:775
AMG(const OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, std::size_t gamma, std::size_t preSmoothingSteps, std::size_t postSmoothingSteps, bool additive=false) DUNE_DEPRECATED
Construct a new amg with a specific coarse solver.
Definition: amg.hh:363
X Domain
The domain type.
Definition: amg.hh:79
ParallelInformationHierarchy::Iterator pinfo
The iterator over the parallel information.
Definition: amg.hh:252
OperatorHierarchy::AggregatesMapList::const_iterator aggregates
The iterator over the aggregates maps.
Definition: amg.hh:260
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:92
S Smoother
The type of the smoother.
Definition: amg.hh:89
Hierarchy< Smoother, A >::Iterator smoother
The iterator over the smoothers.
Definition: amg.hh:244
M Operator
The matrix operator type.
Definition: amg.hh:65
OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator matrix
The iterator over the matrices.
Definition: amg.hh:248
OperatorHierarchy::RedistributeInfoList::const_iterator redist
The iterator over the redistribution information.
Definition: amg.hh:256
X Range
The range type.
Definition: amg.hh:81
void presmooth(LevelContext &levelContext, size_t steps)
Apply pre smoothing on the current level.
Definition: smoother.hh:408
void getCoarsestAggregateNumbers(std::vector< std::size_t, A1 > &cont)
Get the aggregate number of each unknown on the coarsest level.
Definition: amg.hh:914
Hierarchy< Domain, A >::Iterator lhs
The iterator over the left hand side.
Definition: amg.hh:264
const void * Arguments
A type holding all the arguments needed to call the constructor.
Definition: construction.hh:45
void recalculateHierarchy()
Recalculate the matrix hierarchy.
Definition: amg.hh:211
Iterator coarsest()
Get an iterator positioned at the coarsest level.
Definition: hierarchy.hh:1386
OperatorHierarchy::ParallelInformationHierarchy ParallelInformationHierarchy
The parallal data distribution hierarchy type.
Definition: amg.hh:76
InverseOperator< X, X > CoarseSolver
the type of the coarse solver.
Definition: amg.hh:83
void post(Domain &x)
Clean up.
Definition: amg.hh:884
void postsmooth(LevelContext &levelContext, size_t steps)
Apply post smoothing on the current level.
Definition: smoother.hh:430
static T * construct(Arguments &args)
Construct an object with the specified arguments.
Definition: construction.hh:53
std::size_t level
The level index.
Definition: amg.hh:276
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: amg.hh:639
MatrixHierarchy< M, ParallelInformation, A > OperatorHierarchy
The operator hierarchy type.
Definition: amg.hh:74
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:72
@ category
The solver category.
Definition: amg.hh:96
Provides a classes representing the hierarchies in AMG.
Dune namespace.
Definition: alignment.hh:14
STL namespace.
Define base class for scalar product and norm.
Classes for the generic construction and application of the smoothers.
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Traits class for getting the attribute class of a smoother.
Definition: smoother.hh:64
Statistics about the application of an inverse operator.
Definition: solver.hh:32
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
Choose the approriate scalar product for a solver category.
Definition: scalarproducts.hh:76
Classes for using SuperLU with ISTL matrices.
Prolongation and restriction for amg.
Traits for type conversions and type information.
Classes for using UMFPack with ISTL matrices.
#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 (Jul 15, 22:36, 2024)