dune-istl
2.1.1
|
00001 // $Id: amg.hh 1444 2011-01-23 17:24:33Z rebecca $ 00002 #ifndef DUNE_AMG_AMG_HH 00003 #define DUNE_AMG_AMG_HH 00004 00005 #include<memory> 00006 #include<dune/common/exceptions.hh> 00007 #include<dune/istl/paamg/smoother.hh> 00008 #include<dune/istl/paamg/transfer.hh> 00009 #include<dune/istl/paamg/hierarchy.hh> 00010 #include<dune/istl/solvers.hh> 00011 #include<dune/istl/scalarproducts.hh> 00012 #include<dune/istl/superlu.hh> 00013 #include<dune/common/typetraits.hh> 00014 #include<dune/common/exceptions.hh> 00015 00016 namespace Dune 00017 { 00018 namespace Amg 00019 { 00037 template<class M, class X, class S, class P, class K, class A> 00038 class KAMG; 00039 00040 template<class T> 00041 class KAmgTwoGrid; 00042 00050 template<class M, class X, class S, class PI=SequentialInformation, 00051 class A=std::allocator<X> > 00052 class AMG : public Preconditioner<X,X> 00053 { 00054 template<class M1, class X1, class S1, class P1, class K1, class A1> 00055 friend class KAMG; 00056 00057 friend class KAmgTwoGrid<AMG>; 00058 00059 public: 00061 typedef M Operator; 00068 typedef PI ParallelInformation; 00070 typedef MatrixHierarchy<M, ParallelInformation, A> OperatorHierarchy; 00072 typedef typename OperatorHierarchy::ParallelInformationHierarchy ParallelInformationHierarchy; 00073 00075 typedef X Domain; 00077 typedef X Range; 00079 typedef InverseOperator<X,X> CoarseSolver; 00085 typedef S Smoother; 00086 00088 typedef typename SmootherTraits<Smoother>::Arguments SmootherArgs; 00089 00090 enum { 00092 category = S::category 00093 }; 00094 00106 AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver, 00107 const SmootherArgs& smootherArgs, std::size_t gamma, 00108 std::size_t preSmoothingSteps, 00109 std::size_t postSmoothingSteps, bool additive=false); 00110 00125 template<class C> 00126 AMG(const Operator& fineOperator, const C& criterion, 00127 const SmootherArgs& smootherArgs, std::size_t gamma=1, 00128 std::size_t preSmoothingSteps=2, std::size_t postSmoothingSteps=2, 00129 bool additive=false, const ParallelInformation& pinfo=ParallelInformation()); 00130 00131 ~AMG(); 00132 00134 void pre(Domain& x, Range& b); 00135 00137 void apply(Domain& v, const Range& d); 00138 00140 void post(Domain& x); 00141 00146 template<class A1> 00147 void getCoarsestAggregateNumbers(std::vector<std::size_t,A1>& cont); 00148 00149 std::size_t levels(); 00150 00151 std::size_t maxlevels(); 00152 00161 void recalculateHierarchy() 00162 { 00163 matrices_->recalculateGalerkin(NegateSet<typename PI::OwnerSet>()); 00164 } 00165 00166 private: 00168 void mgc(); 00169 00170 typename Hierarchy<Smoother,A>::Iterator smoother; 00171 typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator matrix; 00172 typename ParallelInformationHierarchy::Iterator pinfo; 00173 typename OperatorHierarchy::RedistributeInfoList::const_iterator redist; 00174 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates; 00175 typename Hierarchy<Domain,A>::Iterator lhs; 00176 typename Hierarchy<Domain,A>::Iterator update; 00177 typename Hierarchy<Range,A>::Iterator rhs; 00178 00179 void additiveMgc(); 00180 00182 void presmooth(); 00183 00185 void postsmooth(); 00186 00190 void moveToFineLevel(bool processedFineLevel); 00191 00193 bool moveToCoarseLevel(); 00194 00196 void initIteratorsWithFineLevel(); 00197 00199 OperatorHierarchy* matrices_; 00201 SmootherArgs smootherArgs_; 00203 Hierarchy<Smoother,A> smoothers_; 00205 CoarseSolver* solver_; 00207 Hierarchy<Range,A>* rhs_; 00209 Hierarchy<Domain,A>* lhs_; 00211 Hierarchy<Domain,A>* update_; 00213 typedef Dune::ScalarProductChooser<X,PI,M::category> ScalarProductChooser; 00215 typedef typename ScalarProductChooser::ScalarProduct ScalarProduct; 00217 ScalarProduct* scalarProduct_; 00219 std::size_t gamma_; 00221 std::size_t preSteps_; 00223 std::size_t postSteps_; 00224 std::size_t level; 00225 bool buildHierarchy_; 00226 bool additive; 00227 Smoother *coarseSmoother_; 00228 }; 00229 00230 template<class M, class X, class S, class PI, class A> 00231 AMG<M,X,S,PI,A>::AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver, 00232 const SmootherArgs& smootherArgs, 00233 std::size_t gamma, std::size_t preSmoothingSteps, 00234 std::size_t postSmoothingSteps, bool additive_) 00235 : matrices_(&matrices), smootherArgs_(smootherArgs), 00236 smoothers_(), solver_(&coarseSolver), scalarProduct_(0), 00237 gamma_(gamma), preSteps_(preSmoothingSteps), postSteps_(postSmoothingSteps), buildHierarchy_(false), 00238 additive(additive_), coarseSmoother_() 00239 { 00240 assert(matrices_->isBuilt()); 00241 00242 // build the necessary smoother hierarchies 00243 matrices_->coarsenSmoother(smoothers_, smootherArgs_); 00244 } 00245 00246 template<class M, class X, class S, class PI, class A> 00247 template<class C> 00248 AMG<M,X,S,PI,A>::AMG(const Operator& matrix, 00249 const C& criterion, 00250 const SmootherArgs& smootherArgs, 00251 std::size_t gamma, std::size_t preSmoothingSteps, 00252 std::size_t postSmoothingSteps, 00253 bool additive_, 00254 const PI& pinfo) 00255 : smootherArgs_(smootherArgs), 00256 smoothers_(), solver_(), scalarProduct_(0), gamma_(gamma), 00257 preSteps_(preSmoothingSteps), postSteps_(postSmoothingSteps), buildHierarchy_(true), 00258 additive(additive_), coarseSmoother_() 00259 { 00260 dune_static_assert(static_cast<int>(M::category)==static_cast<int>(S::category), 00261 "Matrix and Solver must match in terms of category!"); 00262 // TODO: reestablish compile time checks. 00263 //dune_static_assert(static_cast<int>(PI::category)==static_cast<int>(S::category), 00264 // "Matrix and Solver must match in terms of category!"); 00265 Timer watch; 00266 matrices_ = new OperatorHierarchy(const_cast<Operator&>(matrix), pinfo); 00267 00268 matrices_->template build<NegateSet<typename PI::OwnerSet> >(criterion); 00269 00270 // build the necessary smoother hierarchies 00271 matrices_->coarsenSmoother(smoothers_, smootherArgs_); 00272 00273 if(criterion.debugLevel()>0 && matrices_->parallelInformation().finest()->communicator().rank()==0) 00274 std::cout<<"Building Hierarchy of "<<matrices_->maxlevels()<<" levels took "<<watch.elapsed()<<" seconds."<<std::endl; 00275 } 00276 00277 template<class M, class X, class S, class PI, class A> 00278 AMG<M,X,S,PI,A>::~AMG() 00279 { 00280 if(buildHierarchy_){ 00281 delete matrices_; 00282 } 00283 if(scalarProduct_) 00284 delete scalarProduct_; 00285 } 00286 00287 00288 template<class M, class X, class S, class PI, class A> 00289 void AMG<M,X,S,PI,A>::pre(Domain& x, Range& b) 00290 { 00291 if(smoothers_.levels()>0) 00292 smoothers_.finest()->pre(x,b); 00293 else 00294 // No smoother to make x consistent! Do it by hand 00295 matrices_->parallelInformation().coarsest()->copyOwnerToAll(x,x); 00296 Range* copy = new Range(b); 00297 rhs_ = new Hierarchy<Range,A>(*copy); 00298 Domain* dcopy = new Domain(x); 00299 lhs_ = new Hierarchy<Domain,A>(*dcopy); 00300 dcopy = new Domain(x); 00301 update_ = new Hierarchy<Domain,A>(*dcopy); 00302 matrices_->coarsenVector(*rhs_); 00303 matrices_->coarsenVector(*lhs_); 00304 matrices_->coarsenVector(*update_); 00305 00306 // Preprocess all smoothers 00307 typedef typename Hierarchy<Smoother,A>::Iterator Iterator; 00308 typedef typename Hierarchy<Range,A>::Iterator RIterator; 00309 typedef typename Hierarchy<Domain,A>::Iterator DIterator; 00310 Iterator coarsest = smoothers_.coarsest(); 00311 Iterator smoother = smoothers_.finest(); 00312 RIterator rhs = rhs_->finest(); 00313 DIterator lhs = lhs_->finest(); 00314 if(smoothers_.levels()>0){ 00315 00316 assert(lhs_->levels()==rhs_->levels()); 00317 assert(smoothers_.levels()==lhs_->levels() || matrices_->levels()==matrices_->maxlevels()); 00318 assert(smoothers_.levels()+1==lhs_->levels() || matrices_->levels()<matrices_->maxlevels()); 00319 00320 if(smoother!=coarsest) 00321 for(++smoother, ++lhs, ++rhs; smoother != coarsest; ++smoother, ++lhs, ++rhs) 00322 smoother->pre(*lhs,*rhs); 00323 smoother->pre(*lhs,*rhs); 00324 } 00325 00326 00327 // The preconditioner might change x and b. So we have to 00328 // copy the changes to the original vectors. 00329 x = *lhs_->finest(); 00330 b = *rhs_->finest(); 00331 00332 if(buildHierarchy_ && matrices_->levels()==matrices_->maxlevels()){ 00333 // We have the carsest level. Create the coarse Solver 00334 SmootherArgs sargs(smootherArgs_); 00335 sargs.iterations = 1; 00336 00337 typename ConstructionTraits<Smoother>::Arguments cargs; 00338 cargs.setArgs(sargs); 00339 if(matrices_->redistributeInformation().back().isSetup()){ 00340 // Solve on the redistributed partitioning 00341 cargs.setMatrix(matrices_->matrices().coarsest().getRedistributed().getmat()); 00342 cargs.setComm(matrices_->parallelInformation().coarsest().getRedistributed()); 00343 00344 coarseSmoother_ = ConstructionTraits<Smoother>::construct(cargs); 00345 scalarProduct_ = ScalarProductChooser::construct(matrices_->parallelInformation().coarsest().getRedistributed()); 00346 }else{ 00347 cargs.setMatrix(matrices_->matrices().coarsest()->getmat()); 00348 cargs.setComm(*matrices_->parallelInformation().coarsest()); 00349 00350 coarseSmoother_ = ConstructionTraits<Smoother>::construct(cargs); 00351 scalarProduct_ = ScalarProductChooser::construct(*matrices_->parallelInformation().coarsest()); 00352 } 00353 #if HAVE_SUPERLU 00354 // Use superlu if we are purely sequential or with only one processor on the coarsest level. 00355 if(is_same<ParallelInformation,SequentialInformation>::value // sequential mode 00356 || matrices_->parallelInformation().coarsest()->communicator().size()==1 //parallel mode and only one processor 00357 || (matrices_->parallelInformation().coarsest().isRedistributed() 00358 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()==1 00359 && matrices_->parallelInformation().coarsest().getRedistributed().communicator().size()>0)){ // redistribute and 1 proc 00360 std::cout<<"Using superlu"<<std::endl; 00361 if(matrices_->parallelInformation().coarsest().isRedistributed()) 00362 { 00363 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0) 00364 // We are still participating on this level 00365 solver_ = new SuperLU<typename M::matrix_type>(matrices_->matrices().coarsest().getRedistributed().getmat()); 00366 else 00367 solver_ = 0; 00368 }else 00369 solver_ = new SuperLU<typename M::matrix_type>(matrices_->matrices().coarsest()->getmat()); 00370 }else 00371 #endif 00372 { 00373 if(matrices_->parallelInformation().coarsest().isRedistributed()) 00374 { 00375 if(matrices_->matrices().coarsest().getRedistributed().getmat().N()>0) 00376 // We are still participating on this level 00377 solver_ = new BiCGSTABSolver<X>(const_cast<M&>(matrices_->matrices().coarsest().getRedistributed()), 00378 *scalarProduct_, 00379 *coarseSmoother_, 1E-2, 10000, 0); 00380 else 00381 solver_ = 0; 00382 }else 00383 solver_ = new BiCGSTABSolver<X>(const_cast<M&>(*matrices_->matrices().coarsest()), 00384 *scalarProduct_, 00385 *coarseSmoother_, 1E-2, 1000, 0); 00386 } 00387 } 00388 } 00389 template<class M, class X, class S, class PI, class A> 00390 std::size_t AMG<M,X,S,PI,A>::levels() 00391 { 00392 return matrices_->levels(); 00393 } 00394 template<class M, class X, class S, class PI, class A> 00395 std::size_t AMG<M,X,S,PI,A>::maxlevels() 00396 { 00397 return matrices_->maxlevels(); 00398 } 00399 00401 template<class M, class X, class S, class PI, class A> 00402 void AMG<M,X,S,PI,A>::apply(Domain& v, const Range& d) 00403 { 00404 if(additive){ 00405 *(rhs_->finest())=d; 00406 additiveMgc(); 00407 v=*lhs_->finest(); 00408 }else{ 00409 // Init all iterators for the current level 00410 initIteratorsWithFineLevel(); 00411 00412 00413 *lhs = v; 00414 *rhs = d; 00415 *update=0; 00416 level=0; 00417 00418 mgc(); 00419 00420 if(postSteps_==0||matrices_->maxlevels()==1) 00421 pinfo->copyOwnerToAll(*update, *update); 00422 00423 v=*update; 00424 } 00425 00426 } 00427 00428 template<class M, class X, class S, class PI, class A> 00429 void AMG<M,X,S,PI,A>::initIteratorsWithFineLevel() 00430 { 00431 smoother = smoothers_.finest(); 00432 matrix = matrices_->matrices().finest(); 00433 pinfo = matrices_->parallelInformation().finest(); 00434 redist = 00435 matrices_->redistributeInformation().begin(); 00436 aggregates = matrices_->aggregatesMaps().begin(); 00437 lhs = lhs_->finest(); 00438 update = update_->finest(); 00439 rhs = rhs_->finest(); 00440 } 00441 00442 template<class M, class X, class S, class PI, class A> 00443 bool AMG<M,X,S,PI,A> 00444 ::moveToCoarseLevel() 00445 { 00446 00447 bool processNextLevel=true; 00448 00449 if(redist->isSetup()){ 00450 redist->redistribute(static_cast<const Range&>(*rhs), rhs.getRedistributed()); 00451 processNextLevel =rhs.getRedistributed().size()>0; 00452 if(processNextLevel){ 00453 //restrict defect to coarse level right hand side. 00454 typename Hierarchy<Range,A>::Iterator fineRhs = rhs++; 00455 ++pinfo; 00456 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation> 00457 ::restrict(*(*aggregates), *rhs, static_cast<const Range&>(fineRhs.getRedistributed()), *pinfo); 00458 } 00459 }else{ 00460 //restrict defect to coarse level right hand side. 00461 typename Hierarchy<Range,A>::Iterator fineRhs = rhs++; 00462 ++pinfo; 00463 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation> 00464 ::restrict(*(*aggregates), *rhs, static_cast<const Range&>(*fineRhs), *pinfo); 00465 } 00466 00467 if(processNextLevel){ 00468 // prepare coarse system 00469 ++lhs; 00470 ++update; 00471 ++matrix; 00472 ++level; 00473 ++redist; 00474 00475 if(matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()){ 00476 // next level is not the globally coarsest one 00477 ++smoother; 00478 ++aggregates; 00479 } 00480 // prepare the update on the next level 00481 *update=0; 00482 } 00483 return processNextLevel; 00484 } 00485 00486 template<class M, class X, class S, class PI, class A> 00487 void AMG<M,X,S,PI,A> 00488 ::moveToFineLevel(bool processNextLevel) 00489 { 00490 if(processNextLevel){ 00491 if(matrix != matrices_->matrices().coarsest() || matrices_->levels()<matrices_->maxlevels()){ 00492 // previous level is not the globally coarsest one 00493 --smoother; 00494 --aggregates; 00495 } 00496 --redist; 00497 --level; 00498 //prolongate and add the correction (update is in coarse left hand side) 00499 --matrix; 00500 00501 //typename Hierarchy<Domain,A>::Iterator coarseLhs = lhs--; 00502 --lhs; 00503 --pinfo; 00504 } 00505 if(redist->isSetup()){ 00506 // Need to redistribute during prolongate 00507 lhs.getRedistributed()=0; 00508 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation> 00509 ::prolongate(*(*aggregates), *update, *lhs, lhs.getRedistributed(), matrices_->getProlongationDampingFactor(), 00510 *pinfo, *redist); 00511 }else{ 00512 *lhs=0; 00513 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation> 00514 ::prolongate(*(*aggregates), *update, *lhs, 00515 matrices_->getProlongationDampingFactor(), *pinfo); 00516 } 00517 00518 00519 if(processNextLevel){ 00520 --update; 00521 --rhs; 00522 } 00523 00524 *update += *lhs; 00525 } 00526 00527 00528 template<class M, class X, class S, class PI, class A> 00529 void AMG<M,X,S,PI,A> 00530 ::presmooth() 00531 { 00532 00533 for(std::size_t i=0; i < preSteps_; ++i){ 00534 *lhs=0; 00535 SmootherApplier<S>::preSmooth(*smoother, *lhs, *rhs); 00536 // Accumulate update 00537 *update += *lhs; 00538 00539 // update defect 00540 matrix->applyscaleadd(-1,static_cast<const Domain&>(*lhs), *rhs); 00541 pinfo->project(*rhs); 00542 } 00543 } 00544 00545 template<class M, class X, class S, class PI, class A> 00546 void AMG<M,X,S,PI,A> 00547 ::postsmooth() 00548 { 00549 00550 for(std::size_t i=0; i < postSteps_; ++i){ 00551 // update defect 00552 matrix->applyscaleadd(-1,static_cast<const Domain&>(*lhs), *rhs); 00553 *lhs=0; 00554 pinfo->project(*rhs); 00555 SmootherApplier<S>::postSmooth(*smoother, *lhs, *rhs); 00556 // Accumulate update 00557 *update += *lhs; 00558 } 00559 } 00560 00561 00562 template<class M, class X, class S, class PI, class A> 00563 void AMG<M,X,S,PI,A>::mgc(){ 00564 if(matrix == matrices_->matrices().coarsest() && levels()==maxlevels()){ 00565 // Solve directly 00566 InverseOperatorResult res; 00567 res.converged=true; // If we do not compute this flag will not get updated 00568 if(redist->isSetup()){ 00569 redist->redistribute(*rhs, rhs.getRedistributed()); 00570 if(rhs.getRedistributed().size()>0){ 00571 // We are still participating in the computation 00572 pinfo.getRedistributed().copyOwnerToAll(rhs.getRedistributed(), rhs.getRedistributed()); 00573 solver_->apply(update.getRedistributed(), rhs.getRedistributed(), res); 00574 } 00575 redist->redistributeBackward(*update, update.getRedistributed()); 00576 pinfo->copyOwnerToAll(*update, *update); 00577 }else{ 00578 pinfo->copyOwnerToAll(*rhs, *rhs); 00579 solver_->apply(*update, *rhs, res); 00580 } 00581 00582 if(!res.converged) 00583 DUNE_THROW(MathError, "Coarse solver did not converge"); 00584 }else{ 00585 // presmoothing 00586 presmooth(); 00587 00588 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION 00589 bool processNextLevel = moveToCoarseLevel(); 00590 00591 if(processNextLevel){ 00592 // next level 00593 for(std::size_t i=0; i<gamma_; i++) 00594 mgc(); 00595 } 00596 00597 moveToFineLevel(processNextLevel); 00598 #else 00599 *lhs=0; 00600 #endif 00601 00602 // postsmoothing 00603 postsmooth(); 00604 00605 } 00606 } 00607 00608 template<class M, class X, class S, class PI, class A> 00609 void AMG<M,X,S,PI,A>::additiveMgc(){ 00610 00611 // restrict residual to all levels 00612 typename ParallelInformationHierarchy::Iterator pinfo=matrices_->parallelInformation().finest(); 00613 typename Hierarchy<Range,A>::Iterator rhs=rhs_->finest(); 00614 typename Hierarchy<Domain,A>::Iterator lhs = lhs_->finest(); 00615 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates=matrices_->aggregatesMaps().begin(); 00616 00617 for(typename Hierarchy<Range,A>::Iterator fineRhs=rhs++; fineRhs != rhs_->coarsest(); fineRhs=rhs++, ++aggregates){ 00618 ++pinfo; 00619 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation> 00620 ::restrict(*(*aggregates), *rhs, static_cast<const Range&>(*fineRhs), *pinfo); 00621 } 00622 00623 // pinfo is invalid, set to coarsest level 00624 //pinfo = matrices_->parallelInformation().coarsest 00625 // calculate correction for all levels 00626 lhs = lhs_->finest(); 00627 typename Hierarchy<Smoother,A>::Iterator smoother = smoothers_.finest(); 00628 00629 for(rhs=rhs_->finest(); rhs != rhs_->coarsest(); ++lhs, ++rhs, ++smoother){ 00630 // presmoothing 00631 *lhs=0; 00632 smoother->apply(*lhs, *rhs); 00633 } 00634 00635 // Coarse level solve 00636 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION 00637 InverseOperatorResult res; 00638 pinfo->copyOwnerToAll(*rhs, *rhs); 00639 solver_->apply(*lhs, *rhs, res); 00640 00641 if(!res.converged) 00642 DUNE_THROW(MathError, "Coarse solver did not converge"); 00643 #else 00644 *lhs=0; 00645 #endif 00646 // Prologate and add up corrections from all levels 00647 --pinfo; 00648 --aggregates; 00649 00650 for(typename Hierarchy<Domain,A>::Iterator coarseLhs = lhs--; coarseLhs != lhs_->finest(); coarseLhs = lhs--, --aggregates, --pinfo){ 00651 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation> 00652 ::prolongate(*(*aggregates), *coarseLhs, *lhs, 1, *pinfo); 00653 } 00654 } 00655 00656 00658 template<class M, class X, class S, class PI, class A> 00659 void AMG<M,X,S,PI,A>::post(Domain& x) 00660 { 00661 if(buildHierarchy_){ 00662 if(solver_) 00663 delete solver_; 00664 if(coarseSmoother_) 00665 ConstructionTraits<Smoother>::deconstruct(coarseSmoother_); 00666 } 00667 00668 // Postprocess all smoothers 00669 typedef typename Hierarchy<Smoother,A>::Iterator Iterator; 00670 typedef typename Hierarchy<Range,A>::Iterator RIterator; 00671 typedef typename Hierarchy<Domain,A>::Iterator DIterator; 00672 Iterator coarsest = smoothers_.coarsest(); 00673 Iterator smoother = smoothers_.finest(); 00674 DIterator lhs = lhs_->finest(); 00675 if(smoothers_.levels()>0){ 00676 if(smoother != coarsest || matrices_->levels()<matrices_->maxlevels()) 00677 smoother->post(*lhs); 00678 if(smoother!=coarsest) 00679 for(++smoother, ++lhs; smoother != coarsest; ++smoother, ++lhs) 00680 smoother->post(*lhs); 00681 smoother->post(*lhs); 00682 } 00683 00684 delete &(*lhs_->finest()); 00685 delete lhs_; 00686 delete &(*update_->finest()); 00687 delete update_; 00688 delete &(*rhs_->finest()); 00689 delete rhs_; 00690 } 00691 00692 template<class M, class X, class S, class PI, class A> 00693 template<class A1> 00694 void AMG<M,X,S,PI,A>::getCoarsestAggregateNumbers(std::vector<std::size_t,A1>& cont) 00695 { 00696 matrices_->getCoarsestAggregatesOnFinest(cont); 00697 } 00698 00699 } // end namespace Amg 00700 }// end namespace Dune 00701 00702 #endif