00001
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 namespace Dune
00014 {
00015 namespace Amg
00016 {
00037 template<class M, class X, class S, class PI=SequentialInformation,
00038 class A=std::allocator<X> >
00039 class AMG : public Preconditioner<X,X>
00040 {
00041 public:
00043 typedef M Operator;
00050 typedef PI ParallelInformation;
00052 typedef MatrixHierarchy<M, ParallelInformation, A> OperatorHierarchy;
00054 typedef typename OperatorHierarchy::ParallelInformationHierarchy ParallelInformationHierarchy;
00055
00057 typedef X Domain;
00059 typedef X Range;
00061 typedef InverseOperator<X,X> CoarseSolver;
00067 typedef S Smoother;
00068
00070 typedef typename SmootherTraits<Smoother>::Arguments SmootherArgs;
00071
00072 enum {
00074 category = S::category
00075 };
00076
00088 AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
00089 const SmootherArgs& smootherArgs, std::size_t gamma,
00090 std::size_t preSmoothingSteps,
00091 std::size_t postSmoothingSteps, bool additive=false);
00092
00107 template<class C>
00108 AMG(const Operator& fineOperator, const C& criterion,
00109 const SmootherArgs& smootherArgs, std::size_t gamma=1,
00110 std::size_t preSmoothingSteps=2, std::size_t postSmoothingSteps=2,
00111 bool additive=false, const ParallelInformation& pinfo=ParallelInformation());
00112
00113 ~AMG();
00114
00116 void pre(Domain& x, Range& b);
00117
00119 void apply(Domain& v, const Range& d);
00120
00122 void post(Domain& x);
00123
00124 private:
00126 void mgc(typename Hierarchy<Smoother,A>::Iterator& smoother,
00127 typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator& matrix,
00128 typename ParallelInformationHierarchy::Iterator& pinfo,
00129 typename OperatorHierarchy::AggregatesMapList::const_iterator& aggregates,
00130 typename Hierarchy<Domain,A>::Iterator& lhs,
00131 typename Hierarchy<Domain,A>::Iterator& update,
00132 typename Hierarchy<Range,A>::Iterator& rhs);
00133
00134 void additiveMgc();
00135
00136
00137
00139 const OperatorHierarchy* matrices_;
00141 SmootherArgs smootherArgs_;
00143 Hierarchy<Smoother,A> smoothers_;
00145 CoarseSolver* solver_;
00147 Hierarchy<Range,A>* rhs_;
00149 Hierarchy<Domain,A>* lhs_;
00151 Hierarchy<Domain,A>* update_;
00153 typedef ScalarProductChooser<X,PI,M::category> ScalarProductChooser;
00155 typedef typename ScalarProductChooser::ScalarProduct ScalarProduct;
00157 ScalarProduct* scalarProduct_;
00159 std::size_t gamma_;
00161 std::size_t preSteps_;
00163 std::size_t postSteps_;
00164 std::size_t level;
00165 bool buildHierarchy_;
00166 bool additive;
00167 Smoother *coarseSmoother_;
00168 };
00169
00170 template<class M, class X, class S, class P, class A>
00171 AMG<M,X,S,P,A>::AMG(const OperatorHierarchy& matrices, CoarseSolver& coarseSolver,
00172 const SmootherArgs& smootherArgs,
00173 std::size_t gamma, std::size_t preSmoothingSteps,
00174 std::size_t postSmoothingSteps, bool additive_)
00175 : matrices_(&matrices), smootherArgs_(smootherArgs),
00176 smoothers_(), solver_(&coarseSolver), scalarProduct_(0),
00177 gamma_(gamma), preSteps_(preSmoothingSteps), postSteps_(postSmoothingSteps), buildHierarchy_(false),
00178 additive(additive_)
00179 {
00180 assert(matrices_->isBuilt());
00181
00182
00183
00184 matrices_->coarsenSmoother(smoothers_, smootherArgs_);
00185
00186 }
00187
00188 template<class M, class X, class S, class P, class A>
00189 template<class C>
00190 AMG<M,X,S,P,A>::AMG(const Operator& matrix,
00191 const C& criterion,
00192 const SmootherArgs& smootherArgs,
00193 std::size_t gamma, std::size_t preSmoothingSteps,
00194 std::size_t postSmoothingSteps,
00195 bool additive_,
00196 const P& pinfo)
00197 : smootherArgs_(smootherArgs),
00198 smoothers_(), scalarProduct_(0), gamma_(gamma),
00199 preSteps_(preSmoothingSteps), postSteps_(postSmoothingSteps), buildHierarchy_(true),
00200 additive(additive_)
00201 {
00202 IsTrue<static_cast<int>(M::category)==static_cast<int>(S::category)>::yes();
00203 IsTrue<static_cast<int>(P::category)==static_cast<int>(S::category)>::yes();
00204 OperatorHierarchy* matrices = new OperatorHierarchy(const_cast<Operator&>(matrix), pinfo);
00205
00206 matrices->template build<typename P::CopyFlags>(criterion);
00207
00208 matrices_ = matrices;
00209
00210 matrices_->coarsenSmoother(smoothers_, smootherArgs_);
00211
00212 }
00213
00214 template<class M, class X, class S, class P, class A>
00215 AMG<M,X,S,P,A>::~AMG()
00216 {
00217 if(buildHierarchy_){
00218 delete matrices_;
00219 }
00220 if(scalarProduct_)
00221 delete scalarProduct_;
00222 }
00223
00224 template<class M, class X, class S, class P, class A>
00225 void AMG<M,X,S,P,A>::pre(Domain& x, Range& b)
00226 {
00227 Range* copy = new Range(b);
00228 rhs_ = new Hierarchy<Range,A>(*copy);
00229 Domain* dcopy = new Domain(x);
00230 lhs_ = new Hierarchy<Domain,A>(*dcopy);
00231 dcopy = new Domain(x);
00232 update_ = new Hierarchy<Domain,A>(*dcopy);
00233 matrices_->coarsenVector(*rhs_);
00234 matrices_->coarsenVector(*lhs_);
00235 matrices_->coarsenVector(*update_);
00236
00237
00238 typedef typename Hierarchy<Smoother,A>::Iterator Iterator;
00239 typedef typename Hierarchy<Range,A>::Iterator RIterator;
00240 typedef typename Hierarchy<Domain,A>::Iterator DIterator;
00241 Iterator coarsest = smoothers_.coarsest();
00242 Iterator smoother = smoothers_.finest();
00243 RIterator rhs = rhs_->finest();
00244 DIterator lhs = lhs_->finest();
00245
00246 if(rhs!=rhs_->coarsest()){
00247 for(; smoother != coarsest; ++smoother, ++lhs, ++rhs)
00248 smoother->pre(*lhs,*rhs);
00249 smoother->pre(*lhs,*rhs);
00250 }
00251
00252
00253
00254 x = *lhs_->finest();
00255 b = *rhs_->finest();
00256
00257 if(buildHierarchy_){
00258
00259 SmootherArgs sargs(smootherArgs_);
00260 sargs.iterations = 1;
00261
00262 typename ConstructionTraits<Smoother>::Arguments cargs;
00263 cargs.setArgs(sargs);
00264 cargs.setMatrix(matrices_->matrices().coarsest()->getmat());
00265 cargs.setComm(*matrices_->parallelInformation().coarsest());
00266
00267 coarseSmoother_ = ConstructionTraits<Smoother>::construct(cargs);
00268 scalarProduct_ = ScalarProductChooser::construct(*matrices_->parallelInformation().coarsest());
00269
00270 #ifdef HAVE_SUPERLU
00271 std::cout<<"Using superlu"<<std::endl;
00272 solver_ = new SuperLU<typename M::matrix_type>(matrices_->matrices().coarsest()->getmat());
00273 #else
00274 solver_ = new BiCGSTABSolver<X>(const_cast<M&>(*matrices_->matrices().coarsest()),
00275 *scalarProduct_,
00276 *coarseSmoother_, 1E-2, 10000, 0);
00277 #endif
00278 }
00279 }
00280
00282 template<class M, class X, class S, class P, class A>
00283 void AMG<M,X,S,P,A>::apply(Domain& v, const Range& d)
00284 {
00285 if(additive){
00286 *(rhs_->finest())=d;
00287 additiveMgc();
00288 v=*lhs_->finest();
00289 }else{
00290 typename Hierarchy<Smoother,A>::Iterator smoother = smoothers_.finest();
00291 typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator matrix = matrices_->matrices().finest();
00292 typename ParallelInformationHierarchy::Iterator pinfo = matrices_->parallelInformation().finest();
00293 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates = matrices_->aggregatesMaps().begin();
00294 typename Hierarchy<Domain,A>::Iterator lhs = lhs_->finest();
00295 typename Hierarchy<Domain,A>::Iterator update = update_->finest();
00296 typename Hierarchy<Range,A>::Iterator rhs = rhs_->finest();
00297
00298 *lhs = v;
00299 *rhs = d;
00300 *update=0;
00301 level=0;
00302 mgc(smoother, matrix, pinfo, aggregates, lhs, update, rhs);
00303 v=*update;
00304 }
00305
00306 }
00307
00308 template<class M, class X, class S, class P, class A>
00309 void AMG<M,X,S,P,A>::mgc(typename Hierarchy<Smoother,A>::Iterator& smoother,
00310 typename OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator& matrix,
00311 typename ParallelInformationHierarchy::Iterator& pinfo,
00312 typename OperatorHierarchy::AggregatesMapList::const_iterator& aggregates,
00313 typename Hierarchy<Domain,A>::Iterator& lhs,
00314 typename Hierarchy<Domain,A>::Iterator& update,
00315 typename Hierarchy<Range,A>::Iterator& rhs){
00316 if(matrix == matrices_->matrices().coarsest()){
00317
00318 InverseOperatorResult res;
00319 pinfo->copyOwnerToAll(*rhs, *rhs);
00320 solver_->apply(*update, *rhs, res);
00321 if(!res.converged)
00322 DUNE_THROW(MathError, "Coarse solver did not converge");
00323 }else{
00324
00325 *lhs=0;
00326 for(std::size_t i=0; i < preSteps_; ++i)
00327 SmootherApplier<S>::preSmooth(*smoother, *lhs, *rhs);
00328
00329
00330 *update += *lhs;
00331
00332 matrix->applyscaleadd(-1,static_cast<const Domain&>(*lhs), *rhs);
00333
00334 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
00335
00336 typename Hierarchy<Range,A>::Iterator fineRhs = rhs++;
00337 ++pinfo;
00338
00339 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00340 ::restrict(*(*aggregates), *rhs, static_cast<const Range&>(*fineRhs), *pinfo);
00341
00342
00343 ++lhs;
00344 ++update;
00345 ++matrix;
00346 ++level;
00347
00348 if(matrix != matrices_->matrices().coarsest()){
00349 ++smoother;
00350 ++aggregates;
00351 }
00352
00353
00354 *update=0;
00355
00356
00357 mgc(smoother, matrix, pinfo, aggregates, lhs, update, rhs);
00358
00359 if(matrix != matrices_->matrices().coarsest()){
00360 --smoother;
00361 --aggregates;
00362 }
00363 --level;
00364
00365 --matrix;
00366
00367
00368 --lhs;
00369 --pinfo;
00370
00371
00372 *lhs=0;
00373 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00374 ::prolongate(*(*aggregates), *update, *lhs, 1.6, *pinfo);
00375
00376 --update;
00377 --rhs;
00378
00379 *update += *lhs;
00380
00381
00382 #endif
00383
00384
00385 matrix->applyscaleadd(-1,static_cast<const Domain&>(*lhs), *rhs);
00386
00387
00388 *lhs=0;
00389 for(std::size_t i=0; i < postSteps_; ++i)
00390 SmootherApplier<S>::postSmooth(*smoother, *lhs, *rhs);
00391
00392 *update += *lhs;
00393 }
00394 }
00395
00396 template<class M, class X, class S, class P, class A>
00397 void AMG<M,X,S,P,A>::additiveMgc(){
00398
00399
00400 typename ParallelInformationHierarchy::Iterator pinfo=matrices_->parallelInformation().finest();
00401 typename Hierarchy<Range,A>::Iterator rhs=rhs_->finest();
00402 typename Hierarchy<Domain,A>::Iterator lhs = lhs_->finest();
00403 typename OperatorHierarchy::AggregatesMapList::const_iterator aggregates=matrices_->aggregatesMaps().begin();
00404
00405 for(typename Hierarchy<Range,A>::Iterator fineRhs=rhs++; fineRhs != rhs_->coarsest(); fineRhs=rhs++, ++aggregates){
00406 ++pinfo;
00407 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00408 ::restrict(*(*aggregates), *rhs, static_cast<const Range&>(*fineRhs), *pinfo);
00409 }
00410
00411
00412
00413
00414 lhs = lhs_->finest();
00415 typename Hierarchy<Smoother,A>::Iterator smoother = smoothers_.finest();
00416
00417 for(rhs=rhs_->finest(); rhs != rhs_->coarsest(); ++lhs, ++rhs, ++smoother){
00418
00419 *lhs=0;
00420 smoother->apply(*lhs, *rhs);
00421 }
00422
00423
00424 #ifndef DUNE_AMG_NO_COARSEGRIDCORRECTION
00425 InverseOperatorResult res;
00426 pinfo->copyOwnerToAll(*rhs, *rhs);
00427 solver_->apply(*lhs, *rhs, res);
00428
00429 if(!res.converged)
00430 DUNE_THROW(MathError, "Coarse solver did not converge");
00431 #else
00432 *lhs=0;
00433 #endif
00434
00435 --pinfo;
00436 --aggregates;
00437
00438 for(typename Hierarchy<Domain,A>::Iterator coarseLhs = lhs--; coarseLhs != lhs_->finest(); coarseLhs = lhs--, --aggregates, --pinfo){
00439 Transfer<typename OperatorHierarchy::AggregatesMap::AggregateDescriptor,Range,ParallelInformation>
00440 ::prolongate(*(*aggregates), *coarseLhs, *lhs, 1, *pinfo);
00441 }
00442 }
00443
00444
00446 template<class M, class X, class S, class P, class A>
00447 void AMG<M,X,S,P,A>::post(Domain& x)
00448 {
00449 if(buildHierarchy_){
00450 delete solver_;
00451 ConstructionTraits<Smoother>::deconstruct(coarseSmoother_);
00452 }
00453
00454
00455 typedef typename Hierarchy<Smoother,A>::Iterator Iterator;
00456 typedef typename Hierarchy<Range,A>::Iterator RIterator;
00457 typedef typename Hierarchy<Domain,A>::Iterator DIterator;
00458 Iterator coarsest = smoothers_.coarsest();
00459 Iterator smoother = smoothers_.finest();
00460 DIterator lhs = lhs_->finest();
00461
00462 if(lhs!= lhs_->coarsest())
00463 smoother->post(*lhs);
00464
00465 if(smoother != coarsest)
00466 for(++smoother; smoother != coarsest; ++smoother, ++lhs)
00467 smoother->post(*lhs);
00468
00469 delete &(*lhs_->finest());
00470 delete lhs_;
00471 delete &(*update_->finest());
00472 delete update_;
00473 delete &(*rhs_->finest());
00474 delete rhs_;
00475 }
00476 }
00477 }
00478
00479 #endif