40 template<
class M,
class X,
class S,
class P,
class K,
class A>
53 template<
class M,
class X,
class S,
class PI=SequentialInformation,
54 class A=std::allocator<X> >
57 template<
class M1,
class X1,
class S1,
class P1,
class K1,
class A1>
149 std::size_t levels();
151 std::size_t maxlevels();
180 void createHierarchies(C& criterion,
Operator& matrix,
206 typename OperatorHierarchy::RedistributeInfoList::const_iterator
redist;
210 typename OperatorHierarchy::AggregatesMapList::const_iterator
aggregates;
234 void mgc(LevelContext& levelContext);
244 void moveToFineLevel(LevelContext& levelContext,
bool processedFineLevel);
250 bool moveToCoarseLevel(LevelContext& levelContext);
256 void initIteratorsWithFineLevel(LevelContext& levelContext);
259 std::shared_ptr<OperatorHierarchy> matrices_;
263 std::shared_ptr<Hierarchy<Smoother,A> > smoothers_;
265 std::shared_ptr<CoarseSolver> solver_;
275 typedef typename ScalarProductChooser::ScalarProduct ScalarProduct;
276 typedef std::shared_ptr<ScalarProduct> ScalarProductPointer;
278 ScalarProductPointer scalarProduct_;
282 std::size_t preSteps_;
284 std::size_t postSteps_;
285 bool buildHierarchy_;
287 bool coarsesolverconverged;
288 std::shared_ptr<Smoother> coarseSmoother_;
290 std::size_t verbosity_;
293 template<
class M,
class X,
class S,
class PI,
class A>
295 : matrices_(amg.matrices_), smootherArgs_(amg.smootherArgs_),
296 smoothers_(amg.smoothers_), solver_(amg.solver_),
297 rhs_(), lhs_(), update_(),
298 scalarProduct_(amg.scalarProduct_), gamma_(amg.gamma_),
299 preSteps_(amg.preSteps_), postSteps_(amg.postSteps_),
300 buildHierarchy_(amg.buildHierarchy_),
301 additive(amg.additive), coarsesolverconverged(amg.coarsesolverconverged),
302 coarseSmoother_(amg.coarseSmoother_), verbosity_(amg.verbosity_)
312 template<
class M,
class X,
class S,
class PI,
class A>
316 : matrices_(&matrices), smootherArgs_(smootherArgs),
318 rhs_(), lhs_(), update_(), scalarProduct_(0),
319 gamma_(parms.getGamma()), preSteps_(parms.getNoPreSmoothSteps()),
320 postSteps_(parms.getNoPostSmoothSteps()), buildHierarchy_(false),
321 additive(parms.getAdditive()), coarsesolverconverged(true),
322 coarseSmoother_(), verbosity_(parms.debugLevel())
324 assert(matrices_->isBuilt());
327 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
330 template<
class M,
class X,
class S,
class PI,
class A>
336 : smootherArgs_(smootherArgs),
338 rhs_(), lhs_(), update_(), scalarProduct_(),
339 gamma_(criterion.getGamma()), preSteps_(criterion.getNoPreSmoothSteps()),
340 postSteps_(criterion.getNoPostSmoothSteps()), buildHierarchy_(true),
341 additive(criterion.getAdditive()), coarsesolverconverged(true),
342 coarseSmoother_(), verbosity_(criterion.debugLevel())
344 static_assert(
static_cast<int>(M::category)==
static_cast<int>(S::category),
345 "Matrix and Solver must match in terms of category!");
349 createHierarchies(criterion,
const_cast<Operator&
>(matrix), pinfo);
353 template<
class M,
class X,
class S,
class PI,
class A>
356 if(buildHierarchy_) {
360 coarseSmoother_.reset();
373 template <
class Matrix,
375 struct DirectSolverSelector
378 enum SolverType { umfpack, superlu, none };
380 static constexpr SolverType solver =
381#if HAVE_SUITESPARSE_UMFPACK
382 UMFPackMethodChooser< field_type > :: valid ? umfpack : none ;
389 template <
class M, SolverType>
392 typedef InverseOperator<Vector,Vector> type;
393 static type* create(
const M& mat,
bool verbose,
bool reusevector )
395 DUNE_THROW(NotImplemented,
"DirectSolver not selected");
398 static std::string name () {
return "None"; }
400#if HAVE_SUITESPARSE_UMFPACK
402 struct Solver< M, umfpack >
404 typedef UMFPack< M > type;
405 static type* create(
const M& mat,
bool verbose,
bool reusevector )
407 return new type(mat, verbose, reusevector );
409 static std::string name () {
return "UMFPack"; }
414 struct Solver< M, superlu >
416 typedef SuperLU< M > type;
417 static type* create(
const M& mat,
bool verbose,
bool reusevector )
419 return new type(mat, verbose, reusevector );
421 static std::string name () {
return "SuperLU"; }
426 typedef Solver< Matrix, solver > SelectedSolver ;
427 typedef typename SelectedSolver :: type DirectSolver;
428 static constexpr bool isDirectSolver = solver != none;
429 static std::string name() {
return SelectedSolver :: name (); }
430 static DirectSolver* create(
const Matrix& mat,
bool verbose,
bool reusevector )
432 return SelectedSolver :: create( mat, verbose, reusevector );
436 template<
class M,
class X,
class S,
class PI,
class A>
438 void AMG<M,X,S,PI,A>::createHierarchies(C& criterion, Operator& matrix,
442 matrices_.reset(
new OperatorHierarchy(matrix, pinfo));
444 matrices_->template build<NegateSet<typename PI::OwnerSet> >(criterion);
447 matrices_->coarsenSmoother(*smoothers_, smootherArgs_);
449 if(buildHierarchy_ && matrices_->levels()==matrices_->maxlevels())
452 SmootherArgs sargs(smootherArgs_);
453 sargs.iterations = 1;
456 cargs.setArgs(sargs);
457 if(matrices_->redistributeInformation().back().isSetup()) {
459 cargs.setMatrix(matrices_->matrices().coarsest().getRedistributed().getmat());
460 cargs.setComm(matrices_->parallelInformation().coarsest().getRedistributed());
462 cargs.setMatrix(matrices_->matrices().coarsest()->getmat());
463 cargs.setComm(*matrices_->parallelInformation().coarsest());
467 scalarProduct_.reset(ScalarProductChooser::construct(cargs.getComm()));
469 typedef DirectSolverSelector< typename M::matrix_type, X > SolverSelector;
472 if( SolverSelector::isDirectSolver &&
473 (std::is_same<ParallelInformation,SequentialInformation>::value
474 || matrices_->parallelInformation().coarsest()->communicator().size()==1
475 || (matrices_->parallelInformation().coarsest().isRedistributed()
476 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()==1
477 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()>0) )
480 if(matrices_->parallelInformation().coarsest().isRedistributed())
482 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0)
485 solver_.reset(SolverSelector::create(matrices_->matrices().coarsest().getRedistributed().getmat(),
false,
false));
492 solver_.reset(SolverSelector::create(matrices_->matrices().coarsest()->getmat(),
false,
false));
494 if(verbosity_>0 && matrices_->parallelInformation().coarsest()->communicator().rank()==0)
495 std::cout<<
"Using a direct coarse solver (" << SolverSelector::name() <<
")" << std::endl;
499 if(matrices_->parallelInformation().coarsest().isRedistributed())
501 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0)
503 solver_.reset(
new BiCGSTABSolver<X>(
const_cast<M&
>(matrices_->matrices().coarsest().getRedistributed()),
505 *coarseSmoother_, 1E-2, 1000, 0));
509 solver_.reset(
new BiCGSTABSolver<X>(
const_cast<M&
>(*matrices_->matrices().coarsest()),
511 *coarseSmoother_, 1E-2, 1000, 0));
515 if(verbosity_>0 && matrices_->parallelInformation().finest()->communicator().rank()==0)
516 std::cout<<
"Building hierarchy of "<<matrices_->maxlevels()<<
" levels "
517 <<
"(inclusive coarse solver) took "<<watch.elapsed()<<
" seconds."<<std::endl;
521 template<
class M,
class X,
class S,
class PI,
class A>
528 typedef typename M::matrix_type
Matrix;
535 const Matrix& mat=matrices_->matrices().finest()->getmat();
536 for(RowIter row=mat.begin(); row!=mat.end(); ++row) {
537 bool isDirichlet =
true;
538 bool hasDiagonal =
false;
540 for(ColIter col=row->begin(); col!=row->end(); ++col) {
541 if(row.index()==col.index()) {
549 if(isDirichlet && hasDiagonal)
550 diagonal.solve(x[row.index()], b[row.index()]);
553 if(smoothers_->levels()>0)
554 smoothers_->finest()->pre(x,b);
557 matrices_->parallelInformation().coarsest()->copyOwnerToAll(x,x);
570 matrices_->coarsenVector(*rhs_);
571 matrices_->coarsenVector(*lhs_);
572 matrices_->coarsenVector(*update_);
578 Iterator coarsest = smoothers_->
coarsest();
579 Iterator smoother = smoothers_->finest();
580 RIterator rhs = rhs_->finest();
581 DIterator lhs = lhs_->finest();
582 if(smoothers_->levels()>0) {
584 assert(lhs_->levels()==rhs_->levels());
585 assert(smoothers_->levels()==lhs_->levels() || matrices_->levels()==matrices_->maxlevels());
586 assert(smoothers_->levels()+1==lhs_->levels() || matrices_->levels()<matrices_->maxlevels());
588 if(smoother!=coarsest)
589 for(++smoother, ++lhs, ++rhs; smoother != coarsest; ++smoother, ++lhs, ++rhs)
590 smoother->pre(*lhs,*rhs);
591 smoother->pre(*lhs,*rhs);
601 template<
class M,
class X,
class S,
class PI,
class A>
604 return matrices_->levels();
606 template<
class M,
class X,
class S,
class PI,
class A>
607 std::size_t AMG<M,X,S,PI,A>::maxlevels()
609 return matrices_->maxlevels();
613 template<
class M,
class X,
class S,
class PI,
class A>
616 LevelContext levelContext;
624 initIteratorsWithFineLevel(levelContext);
627 *levelContext.lhs = v;
628 *levelContext.rhs = d;
629 *levelContext.update=0;
630 levelContext.level=0;
634 if(postSteps_==0||matrices_->maxlevels()==1)
635 levelContext.pinfo->copyOwnerToAll(*levelContext.update, *levelContext.update);
637 v=*levelContext.update;
642 template<
class M,
class X,
class S,
class PI,
class A>
645 levelContext.smoother = smoothers_->finest();
646 levelContext.matrix = matrices_->matrices().finest();
647 levelContext.pinfo = matrices_->parallelInformation().finest();
648 levelContext.redist =
649 matrices_->redistributeInformation().begin();
650 levelContext.aggregates = matrices_->aggregatesMaps().begin();
651 levelContext.lhs = lhs_->finest();
652 levelContext.update = update_->finest();
653 levelContext.rhs = rhs_->finest();
656 template<
class M,
class X,
class S,
class PI,
class A>
658 ::moveToCoarseLevel(LevelContext& levelContext)
661 bool processNextLevel=
true;
663 if(levelContext.redist->isSetup()) {
664 levelContext.redist->redistribute(
static_cast<const Range&
>(*levelContext.rhs),
665 levelContext.rhs.getRedistributed());
666 processNextLevel = levelContext.rhs.getRedistributed().size()>0;
667 if(processNextLevel) {
670 ++levelContext.pinfo;
671 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
672 ::restrictVector(*(*levelContext.aggregates), *levelContext.rhs,
673 static_cast<const Range&
>(fineRhs.getRedistributed()),
674 *levelContext.pinfo);
679 ++levelContext.pinfo;
680 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
681 ::restrictVector(*(*levelContext.aggregates),
682 *levelContext.rhs,
static_cast<const Range&
>(*fineRhs),
683 *levelContext.pinfo);
686 if(processNextLevel) {
689 ++levelContext.update;
690 ++levelContext.matrix;
691 ++levelContext.level;
692 ++levelContext.redist;
694 if(levelContext.matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()) {
696 ++levelContext.smoother;
697 ++levelContext.aggregates;
700 *levelContext.update=0;
702 return processNextLevel;
705 template<
class M,
class X,
class S,
class PI,
class A>
707 ::moveToFineLevel(LevelContext& levelContext,
bool processNextLevel)
709 if(processNextLevel) {
710 if(levelContext.matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()) {
712 --levelContext.smoother;
713 --levelContext.aggregates;
715 --levelContext.redist;
716 --levelContext.level;
718 --levelContext.matrix;
722 --levelContext.pinfo;
724 if(levelContext.redist->isSetup()) {
726 levelContext.lhs.getRedistributed()=0;
727 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
728 ::prolongateVector(*(*levelContext.aggregates), *levelContext.update, *levelContext.lhs,
729 levelContext.lhs.getRedistributed(),
730 matrices_->getProlongationDampingFactor(),
731 *levelContext.pinfo, *levelContext.redist);
734 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
735 ::prolongateVector(*(*levelContext.aggregates), *levelContext.update, *levelContext.lhs,
736 matrices_->getProlongationDampingFactor(),
737 *levelContext.pinfo);
741 if(processNextLevel) {
742 --levelContext.update;
746 *levelContext.update += *levelContext.lhs;
749 template<
class M,
class X,
class S,
class PI,
class A>
752 return IsDirectSolver< CoarseSolver>::value;
755 template<
class M,
class X,
class S,
class PI,
class A>
757 if(levelContext.matrix == matrices_->matrices().coarsest() && levels()==maxlevels()) {
761 if(levelContext.redist->isSetup()) {
762 levelContext.redist->redistribute(*levelContext.rhs, levelContext.rhs.getRedistributed());
763 if(levelContext.rhs.getRedistributed().size()>0) {
765 levelContext.pinfo.getRedistributed().copyOwnerToAll(levelContext.rhs.getRedistributed(),
766 levelContext.rhs.getRedistributed());
767 solver_->apply(levelContext.update.getRedistributed(),
768 levelContext.rhs.getRedistributed(), res);
770 levelContext.redist->redistributeBackward(*levelContext.update, levelContext.update.getRedistributed());
771 levelContext.pinfo->copyOwnerToAll(*levelContext.update, *levelContext.update);
773 levelContext.pinfo->copyOwnerToAll(*levelContext.rhs, *levelContext.rhs);
774 solver_->apply(*levelContext.update, *levelContext.rhs, res);
778 coarsesolverconverged =
false;
783#ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
784 bool processNextLevel = moveToCoarseLevel(levelContext);
786 if(processNextLevel) {
788 for(std::size_t i=0; i<gamma_; i++)
792 moveToFineLevel(levelContext, processNextLevel);
797 if(levelContext.matrix == matrices_->matrices().finest()) {
798 coarsesolverconverged = matrices_->parallelInformation().finest()->communicator().prod(coarsesolverconverged);
799 if(!coarsesolverconverged)
800 DUNE_THROW(MathError,
"Coarse solver did not converge");
808 template<
class M,
class X,
class S,
class PI,
class A>
809 void AMG<M,X,S,PI,A>::additiveMgc(){
812 typename ParallelInformationHierarchy::Iterator pinfo=matrices_->parallelInformation().finest();
815 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates=matrices_->aggregatesMaps().begin();
819 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
820 ::restrictVector(*(*aggregates), *rhs,
static_cast<const Range&
>(*fineRhs), *pinfo);
826 lhs = lhs_->finest();
829 for(rhs=rhs_->finest(); rhs != rhs_->coarsest(); ++lhs, ++rhs, ++smoother) {
832 smoother->apply(*lhs, *rhs);
836#ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
837 InverseOperatorResult res;
838 pinfo->copyOwnerToAll(*rhs, *rhs);
839 solver_->apply(*lhs, *rhs, res);
842 DUNE_THROW(MathError,
"Coarse solver did not converge");
851 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
852 ::prolongateVector(*(*aggregates), *coarseLhs, *lhs, 1.0, *pinfo);
858 template<
class M,
class X,
class S,
class PI,
class A>
865 Iterator coarsest = smoothers_->
coarsest();
866 Iterator smoother = smoothers_->finest();
867 DIterator lhs = lhs_->finest();
868 if(smoothers_->levels()>0) {
869 if(smoother != coarsest || matrices_->levels()<matrices_->maxlevels())
870 smoother->post(*lhs);
871 if(smoother!=coarsest)
872 for(++smoother, ++lhs; smoother != coarsest; ++smoother, ++lhs)
873 smoother->post(*lhs);
874 smoother->post(*lhs);
887 template<
class M,
class X,
class S,
class PI,
class A>
891 matrices_->getCoarsestAggregatesOnFinest(cont);
Parallel algebraic multigrid based on agglomeration.
Definition: amg.hh:56
LevelIterator< Hierarchy< ParallelInformation, Allocator >, ParallelInformation > Iterator
Type of the mutable iterator.
Definition: hierarchy.hh:257
LevelIterator< const Hierarchy< MatrixOperator, Allocator >, const MatrixOperator > ConstIterator
Type of the const iterator.
Definition: hierarchy.hh:260
an algebraic multigrid method using a Krylov-cycle.
Definition: kamg.hh:141
Two grid operator for AMG with Krylov cycle.
Definition: kamg.hh:31
The hierarchies build by the coarsening process.
Definition: hierarchy.hh:317
All parameters for AMG.
Definition: parameters.hh:391
ConstIterator class for sequential access.
Definition: matrix.hh:398
A generic dynamic dense matrix.
Definition: matrix.hh:555
T::field_type field_type
Export the type representing the underlying field.
Definition: matrix.hh:559
row_type::const_iterator ConstColIterator
Const iterator for the entries of each row.
Definition: matrix.hh:583
T block_type
Export the type representing the components.
Definition: matrix.hh:562
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
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
AMG(const AMG &amg)
Copy constructor.
Definition: amg.hh:294
AMG(const OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, const Parameters &parms)
Construct a new amg with a specific coarse solver.
Definition: amg.hh:313
void pre(Domain &x, Range &b)
Prepare the preconditioner.
Definition: amg.hh:522
Hierarchy< Domain, A >::Iterator update
The iterator over the updates.
Definition: amg.hh:218
Hierarchy< Range, A >::Iterator rhs
The iterator over the right hand sided.
Definition: amg.hh:222
bool usesDirectCoarseLevelSolver() const
Check whether the coarse solver used is a direct solver.
Definition: amg.hh:750
X Domain
The domain type.
Definition: amg.hh:78
ParallelInformationHierarchy::Iterator pinfo
The iterator over the parallel information.
Definition: amg.hh:202
OperatorHierarchy::AggregatesMapList::const_iterator aggregates
The iterator over the aggregates maps.
Definition: amg.hh:210
SmootherTraits< Smoother >::Arguments SmootherArgs
The argument type for the construction of the smoother.
Definition: amg.hh:91
S Smoother
The type of the smoother.
Definition: amg.hh:88
Hierarchy< Smoother, A >::Iterator smoother
The iterator over the smoothers.
Definition: amg.hh:194
M Operator
The matrix operator type.
Definition: amg.hh:64
OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator matrix
The iterator over the matrices.
Definition: amg.hh:198
OperatorHierarchy::RedistributeInfoList::const_iterator redist
The iterator over the redistribution information.
Definition: amg.hh:206
X Range
The range type.
Definition: amg.hh:80
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:889
Hierarchy< Domain, A >::Iterator lhs
The iterator over the left hand side.
Definition: amg.hh:214
const void * Arguments
A type holding all the arguments needed to call the constructor.
Definition: construction.hh:44
void recalculateHierarchy()
Recalculate the matrix hierarchy.
Definition: amg.hh:161
Iterator coarsest()
Get an iterator positioned at the coarsest level.
Definition: hierarchy.hh:1385
OperatorHierarchy::ParallelInformationHierarchy ParallelInformationHierarchy
The parallal data distribution hierarchy type.
Definition: amg.hh:75
InverseOperator< X, X > CoarseSolver
the type of the coarse solver.
Definition: amg.hh:82
void post(Domain &x)
Clean up.
Definition: amg.hh:859
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:52
std::size_t level
The level index.
Definition: amg.hh:226
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:332
void apply(Domain &v, const Range &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: amg.hh:614
MatrixHierarchy< M, ParallelInformation, A > OperatorHierarchy
The operator hierarchy type.
Definition: amg.hh:73
PI ParallelInformation
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the...
Definition: amg.hh:71
@ category
The solver category.
Definition: amg.hh:95
Provides a classes representing the hierarchies in AMG.
Dune namespace.
Definition: alignment.hh:11
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:77
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 intentionally unused function parameters with.
Definition: unused.hh:18