dune-istl  2.1.1
novlpschwarz.hh
Go to the documentation of this file.
00001 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
00002 // vi: set ts=4 sw=2 et sts=2:
00003 #ifndef DUNE_NOVLPSCHWARZ_HH
00004 #define DUNE_NOVLPSCHWARZ_HH
00005 
00006 #include <iostream>              // for input/output to shell
00007 #include <fstream>               // for input/output to files
00008 #include <vector>                // STL vector class
00009 #include <sstream>
00010 
00011 #include <cmath>                // Yes, we do some math here
00012 #include <sys/times.h>           // for timing measurements
00013 
00014 #include <dune/common/timer.hh>
00015 
00016 #include"io.hh"
00017 #include"bvector.hh"
00018 #include"vbvector.hh"
00019 #include"bcrsmatrix.hh"
00020 #include"io.hh"
00021 #include"gsetc.hh"
00022 #include"ilu.hh"
00023 #include"operators.hh"
00024 #include"solvers.hh"
00025 #include"preconditioners.hh"
00026 #include"scalarproducts.hh"
00027 #include"owneroverlapcopy.hh"
00028 
00029 namespace Dune {
00030 
00059   template<class M, class X, class Y, class C>
00060   class NonoverlappingSchwarzOperator : public AssembledLinearOperator<M,X,Y>
00061   {
00062   public:
00064     typedef M matrix_type;
00066     typedef X domain_type;
00068     typedef Y range_type;
00070     typedef typename X::field_type field_type;
00072     typedef C communication_type;
00073 
00074     typedef typename C::PIS PIS;
00075     typedef typename C::RI RI;
00076     typedef typename RI::RemoteIndexList RIL;
00077     typedef typename RI::const_iterator RIIterator;
00078     typedef typename RIL::const_iterator RILIterator;
00079     typedef typename M::ConstColIterator ColIterator;
00080     typedef typename M::ConstRowIterator RowIterator;
00081     typedef std::multimap<int,int> MM;
00082     
00083         enum {
00085           category=SolverCategory::nonoverlapping
00086         };
00087 
00095     NonoverlappingSchwarzOperator (const matrix_type& A, const communication_type& com) 
00096       : _A_(A), communication(com), buildcomm(true)
00097     {}
00098 
00100     virtual void apply (const X& x, Y& y) const
00101     {
00102       y = 0;
00103       novlp_op_apply(x,y,1);
00104       communication.addOwnerCopyToOwnerCopy(y,y); 
00105     }
00106 
00108     virtual void applyscaleadd (field_type alpha, const X& x, Y& y) const
00109     {
00110           // only apply communication to alpha*A*x to make it consistent,
00111           // y already has to be consistent. 
00112           Y y1(y);
00113           y = 0;
00114           novlp_op_apply(x,y,alpha);
00115           communication.addOwnerCopyToOwnerCopy(y,y);
00116           y += y1;
00117     }
00118     
00120     virtual const matrix_type& getmat () const
00121     {
00122       return _A_;
00123     }
00124 
00125     void novlp_op_apply (const X& x, Y& y, field_type alpha) const
00126     { 
00127       //get index sets
00128       const PIS& pis=communication.indexSet();
00129       const RI& ri = communication.remoteIndices();
00130       
00131       // at the beginning make a multimap "bordercontribution".
00132       // process has i and j as border dofs but is not the owner
00133       // => only contribute to Ax if i,j is in bordercontribution     
00134       if (buildcomm == true){
00135       
00136         // set up mask vector
00137         if (mask.size()!=static_cast<typename std::vector<double>::size_type>(x.size())){
00138           mask.resize(x.size());
00139           for (typename std::vector<double>::size_type i=0; i<mask.size(); i++) 
00140             mask[i] = 1; 
00141           for (typename PIS::const_iterator i=pis.begin(); i!=pis.end(); ++i)
00142             if (i->local().attribute()==OwnerOverlapCopyAttributeSet::copy)
00143               mask[i->local().local()] = 0;
00144             else if (i->local().attribute()==OwnerOverlapCopyAttributeSet::overlap)
00145               mask[i->local().local()] = 2;
00146         }
00147         
00148         for (MM::iterator iter = bordercontribution.begin();
00149              iter != bordercontribution.end(); ++iter)
00150           bordercontribution.erase(iter);
00151         for (RowIterator i = _A_.begin(); i != _A_.end(); ++i){
00152           if (mask[i.index()] == 0){
00153             std::set<int> neighbours; //processes have i as interior/border dof
00154             int iowner; //process which owns i
00155             for (RIIterator remote = ri.begin(); remote != ri.end(); ++remote){
00156               RIL& ril = *(remote->second.first);
00157               for (RILIterator rindex = ril.begin(); rindex != ril.end(); ++rindex){
00158                 if (rindex->attribute() != OwnerOverlapCopyAttributeSet::overlap){
00159                   if (rindex->localIndexPair().local().local() == i.index()){
00160                     neighbours.insert(remote->first);
00161                     if(rindex->attribute()==OwnerOverlapCopyAttributeSet::owner)
00162                       iowner = remote->first;
00163                   }
00164                 }
00165               }
00166             }
00167             for (ColIterator j = _A_[i.index()].begin(); j != _A_[i.index()].end(); ++j){
00168               if (mask[j.index()] == 0){
00169                 bool flag = true;
00170                 for (RIIterator remote = ri.begin(); remote != ri.end(); ++remote)
00171                   if (neighbours.find(remote->first) != neighbours.end()){  
00172                     RIL& ril = *(remote->second.first); 
00173                     for (RILIterator rindex = ril.begin(); rindex != ril.end(); 
00174                          ++rindex)
00175                       if (rindex->attribute() != OwnerOverlapCopyAttributeSet::overlap 
00176                           && rindex->localIndexPair().local().local()==j.index()){
00177                         if (rindex->attribute() == OwnerOverlapCopyAttributeSet::owner
00178                             || remote->first == iowner
00179                             || remote->first < communication.communicator().rank()){
00180                           flag = false;
00181                           continue;
00182                         }
00183                       }
00184                     if (flag == false)
00185                       continue;
00186                   }
00187                 // don“t contribute to Ax if 
00188                 // 1. the owner of j has i as interior/border dof
00189                 // 2. iowner has j as interior/border dof 
00190                 // 3. there is another process with smaller rank that has i and j 
00191                 // as interor/border dofs
00192                 // if the owner of j does not have i as interior/border dof, 
00193                 // it will not be taken into account
00194                 if (flag==true)
00195                   bordercontribution.insert(std::pair<int,int>(i.index(),j.index()));
00196               }
00197             }
00198           }
00199         }
00200         buildcomm = false;
00201       }
00202       
00203       //compute alpha*A*x nonoverlapping case
00204       for (RowIterator i = _A_.begin(); i != _A_.end(); ++i){ 
00205         if (mask[i.index()] == 0){ 
00206           //dof doesn't belong to process but is border (not ghost)
00207           for (ColIterator j = _A_[i.index()].begin(); j != _A_[i.index()].end(); ++j){
00208             if (mask[j.index()] == 1) //j is owner => then sum entries
00209               (*j).usmv(alpha,x[j.index()],y[i.index()]);
00210             else if (mask[j.index()] == 0){
00211               std::pair<MM::iterator, MM::iterator> itp = 
00212                 bordercontribution.equal_range(i.index());
00213               for (MM::iterator it = itp.first; it != itp.second; ++it) 
00214                 if ((*it).second == (int)j.index())
00215                   (*j).usmv(alpha,x[j.index()],y[i.index()]); 
00216             }
00217           }
00218         }
00219         else if (mask[i.index()] == 1){
00220           for (ColIterator j = _A_[i.index()].begin(); j != _A_[i.index()].end(); ++j)
00221             if (mask[j.index()] != 2)
00222               (*j).usmv(alpha,x[j.index()],y[i.index()]);
00223         }
00224       }
00225     }
00226     
00227   private:
00228     const matrix_type& _A_;
00229     const communication_type& communication;
00230     mutable bool buildcomm;
00231     mutable std::vector<double> mask;
00232     mutable std::multimap<int,int>  bordercontribution;
00233   };
00234   
00246   template<class X, class C>
00247   class NonoverlappingSchwarzScalarProduct : public ScalarProduct<X>
00248   {
00249   public:
00251     typedef X domain_type;
00253     typedef typename X::field_type field_type;
00255     typedef C communication_type;
00256     
00258         enum {category=SolverCategory::nonoverlapping};
00259 
00264         NonoverlappingSchwarzScalarProduct (const communication_type& com)
00265           : communication(com)
00266         {}
00267 
00272         virtual field_type dot (const X& x, const X& y)
00273         {
00274           field_type result;
00275           communication.dot(x,y,result);
00276           return result;
00277         }
00278 
00282         virtual double norm (const X& x)
00283         {
00284           return communication.norm(x);
00285         }
00286     
00289     void make_consistent (X& x) const
00290     {
00291       communication.copyOwnerToAll(x,x);
00292     }
00293     
00294   private:
00295         const communication_type& communication;
00296   };
00297 
00298   template<class X, class C>
00299   struct ScalarProductChooser<X,C,SolverCategory::nonoverlapping>
00300   {
00302     typedef NonoverlappingSchwarzScalarProduct<X,C> ScalarProduct;
00304     typedef C communication_type;
00305 
00306     enum{
00308       solverCategory=SolverCategory::nonoverlapping
00309     };
00310 
00311     static ScalarProduct* construct(const communication_type& comm)
00312     {
00313       return new ScalarProduct(comm);
00314     }
00315   };
00316 
00317   namespace Amg
00318   {
00319     template<class T> class ConstructionTraits;
00320   }
00321 
00331   template<class C, class P>
00332   class NonoverlappingBlockPreconditioner
00333     : public Dune::Preconditioner<typename P::domain_type,typename P::range_type> {
00334     friend class Amg::ConstructionTraits<NonoverlappingBlockPreconditioner<C,P> >;
00335   public:
00337     typedef typename P::domain_type domain_type;
00339     typedef typename P::range_type range_type;
00341     typedef C communication_type;
00342 
00343     // define the category
00344     enum {
00346       category=SolverCategory::nonoverlapping};
00347     
00355     NonoverlappingBlockPreconditioner (P& prec, const communication_type& c)
00356       : preconditioner(prec), communication(c)
00357     {}
00358     
00364     virtual void pre (domain_type& x, range_type& b) 
00365         {
00366           preconditioner.pre(x,b);
00367         }
00368 
00374     virtual void apply (domain_type& v, const range_type& d)
00375     {
00376       // block preconditioner equivalent to WrappedPreconditioner from 
00377       // pdelab/backend/ovlpistsolverbackend.hh, 
00378       // but not to BlockPreconditioner from schwarz.hh
00379       preconditioner.apply(v,d);
00380       communication.addOwnerCopyToOwnerCopy(v,v);
00381     }
00382 
00388     virtual void post (domain_type& x) 
00389         {
00390           preconditioner.post(x);
00391         }
00392 
00393   private:
00395         P& preconditioner;
00396 
00398         const communication_type& communication;
00399   };
00400 
00403 } // end namespace
00404 
00405 #endif