overlappingschwarz.hh

Go to the documentation of this file.
00001 #ifndef DUNE_OVERLAPPINGSCHWARZ_HH
00002 #define DUNE_OVERLAPPINGSCHWARZ_HH
00003 #include<cassert>
00004 #include<algorithm>
00005 #include<functional>
00006 #include<vector>
00007 #include<set>
00008 #include<dune/common/sllist.hh>
00009 #include"preconditioners.hh"
00010 #include"superlu.hh"
00011 #include"bvector.hh"
00012 #include"bcrsmatrix.hh"
00013 
00014 namespace Dune
00015 {
00016 
00027 #ifdef HAVE_SUPERLU
00028 
00032   template<class I, class S>
00033   class OverlappingSchwarzInitializer
00034   {
00035   public:
00036     typedef I InitializerList;
00037     typedef typename InitializerList::value_type AtomInitializer;
00038     typedef typename AtomInitializer::Matrix Matrix;
00039     typedef typename Matrix::const_iterator Iter;
00040     typedef typename Matrix::row_type::const_iterator CIter;
00041     
00042     typedef S IndexSet;
00043     typedef typename IndexSet::size_type size_type;
00044     
00045     OverlappingSchwarzInitializer(InitializerList& il,
00046                                   const IndexSet& indices);
00047     
00048     
00049     void addRowNnz(const Iter& row);
00050     
00051     void allocate();
00052         
00053     void countEntries(const Iter& row, const CIter& col) const;
00054 
00055     void calcColstart() const;
00056     
00057     void copyValue(const Iter& row, const CIter& col) const;
00058     
00059     void createMatrix() const;
00060   private:
00061     class IndexMap
00062     {
00063     public:
00064       typedef typename S::size_type size_type;
00065       typedef std::map<size_type,size_type> Map;      
00066       typedef typename Map::iterator iterator;
00067       typedef typename Map::const_iterator const_iterator;
00068       
00069       IndexMap();
00070       
00071       void insert(size_type grow);
00072       
00073       const_iterator find(size_type grow) const;
00074       
00075       iterator find(size_type grow);
00076       
00077       iterator end();
00078       
00079       const_iterator end() const;
00080       
00081     private:
00082       std::map<size_type,size_type> map_;
00083       size_type row;
00084     };
00085     
00086     
00087     typedef typename InitializerList::iterator InitIterator;
00088     typedef typename IndexSet::const_iterator IndexIteratur;
00089     InitializerList* initializers;
00090     const IndexSet *indices;
00091     std::vector<IndexMap> indexMaps;
00092   };
00093 
00097   struct AdditiveSchwarzMode
00098   {};
00099 
00103   struct MultiplicativeSchwarzMode
00104   {};
00105   
00110   struct SymmetricMultiplicativeSchwarzMode
00111   {};  
00112 
00113   namespace
00114   {
00115     template<typename T>
00116     struct AdditiveAdder
00117     {
00118     };
00119     
00120     template<typename T, typename A, int n>
00121     struct AdditiveAdder<BlockVector<FieldVector<T,n>,A> >
00122     {
00123       typedef typename A::size_type size_type;
00124       AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x, const T* t, const T& relax_);
00125       void operator()(const size_type& domain);
00126       void axpy();
00127       
00128     private:
00129       const T* t;
00130       BlockVector<FieldVector<T,n>,A>* v;
00131       BlockVector<FieldVector<T,n>,A>* x;
00132       T relax;
00133       size_type i;
00134     };
00135 
00136     template<typename T>
00137     struct MultiplicativeAdder
00138     {
00139     };
00140     
00141     template<typename T, typename A, int n>
00142     struct MultiplicativeAdder<BlockVector<FieldVector<T,n>,A> >
00143     {
00144       typedef typename A::size_type size_type;
00145       MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x, const T* t, const T& relax_);
00146       void operator()(const size_type& domain);
00147       void axpy();
00148       
00149     private:
00150       const T* t;
00151       BlockVector<FieldVector<T,n>,A>* v;
00152       BlockVector<FieldVector<T,n>,A>* x;
00153       T relax;
00154       size_type i;
00155     };
00156 
00166     template<typename T, class X>
00167     struct AdderSelector
00168     {};
00169     
00170     template<class X>
00171     struct AdderSelector<AdditiveSchwarzMode,X>
00172     {
00173       typedef AdditiveAdder<X> Adder;
00174     };
00175     
00176     template<class X>
00177     struct AdderSelector<MultiplicativeSchwarzMode,X>
00178     {
00179       typedef MultiplicativeAdder<X> Adder;
00180     };    
00181 
00182     template<class X>
00183     struct AdderSelector<SymmetricMultiplicativeSchwarzMode,X>
00184     {
00185       typedef MultiplicativeAdder<X> Adder;
00186     };    
00187 
00199     template<typename T1, typename T2, bool forward>
00200     struct IteratorDirectionSelector
00201     {
00202       typedef T1 solver_vector;
00203       typedef typename solver_vector::iterator solver_iterator;
00204       typedef T2 subdomain_vector;
00205       typedef typename subdomain_vector::const_iterator domain_iterator;
00206     
00207       static solver_iterator begin(solver_vector& sv)
00208       {
00209         return sv.begin();
00210       }
00211       
00212       static solver_iterator end(solver_vector& sv)
00213       {
00214         return sv.end();
00215       }
00216       static domain_iterator begin(const subdomain_vector& sv)
00217       {
00218         return sv.begin();
00219       }
00220       
00221       static domain_iterator end(const subdomain_vector& sv)
00222       {
00223         return sv.end();
00224       }
00225     };
00226 
00227     template<typename T1, typename T2>
00228     struct IteratorDirectionSelector<T1,T2,false>
00229     {
00230       typedef T1 solver_vector;
00231       typedef typename solver_vector::reverse_iterator solver_iterator;
00232       typedef T2 subdomain_vector;
00233       typedef typename subdomain_vector::const_reverse_iterator domain_iterator;
00234     
00235       static solver_iterator begin(solver_vector& sv)
00236       {
00237         return sv.rbegin();
00238       }
00239       
00240       static solver_iterator end(solver_vector& sv)
00241       {
00242         return sv.rend();
00243       }
00244       static domain_iterator begin(const subdomain_vector& sv)
00245       {
00246         return sv.rbegin();
00247       }
00248       
00249       static domain_iterator end(const subdomain_vector& sv)
00250       {
00251         return sv.rend();
00252       }
00253     };
00254     
00263     template<class T>
00264     struct Applier
00265     {
00266       typedef T smoother;
00267       typedef typename smoother::range_type range_type;
00268       
00269       static void apply(smoother& sm, range_type& v, const range_type& b)
00270       {
00271         sm.template apply<true>(v, b);
00272       }
00273     };
00274     
00275     template<class M, class X, bool onTheFly, class TA>
00276     struct Applier<SeqOverlappingSchwarz<M,X,SymmetricMultiplicativeSchwarzMode,onTheFly,TA> >
00277     {
00278       typedef SeqOverlappingSchwarz<M,X,SymmetricMultiplicativeSchwarzMode,onTheFly,TA> smoother;
00279       typedef typename smoother::range_type range_type;
00280       
00281       static void apply(smoother& sm, range_type& v, const range_type& b)
00282       {
00283         sm.template apply<true>(v, b);
00284         sm.template apply<false>(v, b);
00285       }
00286     };
00287   }
00288   
00301   template<class M, class X, class TM=AdditiveSchwarzMode, bool onTheFly=true, class TA=std::allocator<X> >
00302   class SeqOverlappingSchwarz
00303     : public Preconditioner<X,X>
00304   {
00305   public:
00309     typedef M matrix_type;
00310 
00314     typedef X domain_type;
00315 
00319     typedef X range_type;
00320 
00327     typedef TM Mode;
00328 
00332     typedef typename X::field_type field_type;
00333 
00335     typedef typename matrix_type::size_type size_type;
00336 
00338     typedef TA allocator;
00339         
00341     typedef std::set<size_type, std::less<size_type>, 
00342                      typename TA::template rebind<size_type>::other> 
00343       subdomain_type;
00344 
00346     typedef std::vector<subdomain_type, typename TA::template rebind<subdomain_type>::other> subdomain_vector;
00347 
00349     typedef SLList<size_type, typename TA::template rebind<size_type>::other> subdomain_list;
00350       
00352     typedef std::vector<subdomain_list, typename TA::template rebind<subdomain_list>::other > rowtodomain_vector;
00353 
00355     typedef SuperLU<matrix_type> slu;
00356     
00358     typedef std::vector<slu, typename TA::template rebind<slu>::other> slu_vector;
00359 
00360     enum{
00362       category = SolverCategory::sequential
00363         };
00364     
00373     SeqOverlappingSchwarz(const matrix_type& mat, const subdomain_vector& subDomains,
00374                           field_type relaxationFactor=1);
00375 
00382     SeqOverlappingSchwarz(const matrix_type& mat, const rowtodomain_vector& rowToDomain,
00383                           field_type relaxationFactor=1);
00384                           
00390     virtual void pre (X& x, X& b) {}
00391 
00397     virtual void apply (X& v, const X& d);
00398 
00399     template<bool forward>
00400     void apply(X& v, const X& d);
00401     
00407     virtual void post (X& x) {}
00408   
00409   private:
00410     const M& mat;
00411     slu_vector solvers;
00412     subdomain_vector subDomains;
00413     field_type relax;
00414     
00415     typename M::size_type maxlength;
00416 
00417     void initialize(const rowtodomain_vector& rowToDomain, const matrix_type& mat);
00418 
00419     template<typename T>
00420     struct Assigner
00421     {
00422     };
00423     
00424     template<typename T, typename A, int n>
00425     struct Assigner<BlockVector<FieldVector<T,n>,A> >
00426     {
00427       Assigner(const M& mat, T* rhs, const BlockVector<FieldVector<T,n>,A>& b, const BlockVector<FieldVector<T,n>,A>& x);
00428       void operator()(const size_type& domain);
00429     private:
00430       const M* mat;
00431       T* rhs;
00432       const BlockVector<FieldVector<T,n>,A>* b;
00433       const BlockVector<FieldVector<T,n>,A>* x;
00434       size_type i;
00435     };
00436 
00437   };
00438   
00439   
00440   template<class I, class S>
00441   OverlappingSchwarzInitializer<I,S>::OverlappingSchwarzInitializer(InitializerList& il,
00442                                                                     const IndexSet& idx)
00443     : initializers(&il), indices(&idx), indexMaps(il.size())
00444   {
00445   }
00446   
00447   
00448   template<class I, class S>
00449   void OverlappingSchwarzInitializer<I,S>::addRowNnz(const Iter& row)
00450   {
00451     typedef typename IndexSet::value_type::const_iterator iterator;
00452     for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain){
00453       (*initializers)[*domain].addRowNnz(row);
00454       indexMaps[*domain].insert(row.index());
00455     }
00456   }
00457   
00458   template<class I, class S>
00459   void OverlappingSchwarzInitializer<I,S>::allocate()
00460   {
00461     std::for_each(initializers->begin(), initializers->end(),
00462                   std::mem_fun_ref(&AtomInitializer::allocateMatrixStorage));
00463     std::for_each(initializers->begin(), initializers->end(), 
00464                   std::mem_fun_ref(&AtomInitializer::allocateMarker));
00465   }
00466   
00467   template<class I, class S>
00468   void OverlappingSchwarzInitializer<I,S>::countEntries(const Iter& row, const CIter& col) const
00469   {
00470     typedef typename IndexSet::value_type::const_iterator iterator;
00471     for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain){
00472       typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
00473       if(v!= indexMaps[*domain].end()){
00474         (*initializers)[*domain].countEntries(indexMaps[*domain].find(col.index())->second);
00475       }
00476     }
00477   }
00478 
00479   template<class I, class S>
00480   void OverlappingSchwarzInitializer<I,S>::calcColstart() const
00481   {
00482     std::for_each(initializers->begin(), initializers->end(),
00483                   std::mem_fun_ref(&AtomInitializer::calcColstart));
00484   }
00485 
00486   template<class I, class S>
00487   void OverlappingSchwarzInitializer<I,S>::copyValue(const Iter& row, const CIter& col) const
00488   {
00489     typedef typename IndexSet::value_type::const_iterator iterator;
00490     for(iterator domain=(*indices)[row.index()].begin(); domain!= (*indices)[row.index()].end(); ++domain){
00491       typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
00492       if(v!= indexMaps[*domain].end()){
00493         (*initializers)[*domain].copyValue(col, indexMaps[*domain].find(row.index())->second, 
00494                                            indexMaps[*domain].find(col.index())->second);
00495       }
00496     }
00497   }
00498     
00499   template<class I, class S>
00500   void OverlappingSchwarzInitializer<I,S>::createMatrix() const
00501   {
00502     std::for_each(initializers->begin(), initializers->end(),
00503                   std::mem_fun_ref(&AtomInitializer::createMatrix));
00504   }
00505 
00506   template<class I, class S>
00507   OverlappingSchwarzInitializer<I,S>::IndexMap::IndexMap()
00508     : row(0)
00509   {}
00510 
00511   template<class I, class S>
00512   void OverlappingSchwarzInitializer<I,S>::IndexMap::insert(size_type grow)
00513   {
00514     assert(map_.find(grow)==map_.end());
00515     map_.insert(std::make_pair(grow, row++));
00516   }
00517 
00518   template<class I, class S>
00519   typename OverlappingSchwarzInitializer<I,S>::IndexMap::const_iterator 
00520   OverlappingSchwarzInitializer<I,S>::IndexMap::find(size_type grow) const
00521   {
00522     return map_.find(grow);
00523   }
00524   
00525   template<class I, class S>
00526   typename OverlappingSchwarzInitializer<I,S>::IndexMap::iterator 
00527   OverlappingSchwarzInitializer<I,S>::IndexMap::find(size_type grow)
00528   {
00529     return map_.find(grow);
00530   }
00531 
00532   template<class I, class S>
00533   typename OverlappingSchwarzInitializer<I,S>::IndexMap::const_iterator 
00534   OverlappingSchwarzInitializer<I,S>::IndexMap::end() const
00535   {
00536     return map_.end();
00537   }
00538   
00539   template<class I, class S>
00540   typename OverlappingSchwarzInitializer<I,S>::IndexMap::iterator 
00541   OverlappingSchwarzInitializer<I,S>::IndexMap::end()
00542   {
00543     return map_.end();
00544   }
00545 
00546   template<class M, class X, class TM, bool onTheFly, class TA>
00547   SeqOverlappingSchwarz<M,X,TM,onTheFly,TA>::SeqOverlappingSchwarz(const matrix_type& mat_, const rowtodomain_vector& rowToDomain,
00548                         field_type relaxationFactor)
00549     : mat(mat_), relax(relaxationFactor)
00550   {
00551     typedef typename rowtodomain_vector::const_iterator RowDomainIterator;
00552     typedef typename subdomain_list::const_iterator DomainIterator;
00553 #ifdef DUNE_ISTL_WITH_CHECKING
00554     assert(rowToDomain.size()==mat.N());
00555     assert(rowToDomain.size()==mat.M());
00556     
00557     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
00558       assert(iter->size()>0);
00559     
00560 #endif
00561     // calculate the number of domains
00562     size_type domains=0;
00563     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
00564       for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
00565         domains=std::max(domains, *d);
00566     ++domains;
00567     
00568     solvers.resize(domains);
00569     subDomains.resize(domains);
00570     
00571     // initialize subdomains to row mapping from row to subdomain mapping
00572     size_type row=0;
00573     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter, ++row)
00574       for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
00575         subDomains[*d].insert(row);
00576 
00577 #ifdef DUNE_ISTL_WITH_CHECKING
00578     size_type i=0;
00579     typedef typename subdomain_vector::const_iterator iterator;
00580     for(iterator iter=subDomains.begin(); iter != subDomains.end(); ++iter){
00581       typedef typename subdomain_type::const_iterator entry_iterator;
00582       std::cout<<"domain "<<i++<<":";
00583       for(entry_iterator entry = iter->begin(); entry != iter->end(); ++entry){
00584         std::cout<<" "<<*entry;
00585       }
00586       std::cout<<std::endl;
00587     }
00588 #endif
00589     initialize(rowToDomain, mat);
00590   }
00591   
00592   template<class M, class X, class TM, bool onTheFly, class TA>
00593   SeqOverlappingSchwarz<M,X,TM,onTheFly,TA>::SeqOverlappingSchwarz(const matrix_type& mat_,
00594                                                        const subdomain_vector& sd,
00595                                                        field_type relaxationFactor)
00596     :  mat(mat_), solvers(sd.size()), subDomains(sd), relax(relaxationFactor)
00597   {
00598     typedef typename subdomain_vector::const_iterator DomainIterator;
00599 
00600 #ifdef DUNE_ISTL_WITH_CHECKING
00601     size_type i=0;
00602     
00603     for(DomainIterator d=sd.begin(); d != sd.end(); ++d,++i){
00604       //std::cout<<i<<": "<<d->size()<<std::endl;
00605       assert(d->size()>0);
00606       typedef typename DomainIterator::value_type::const_iterator entry_iterator;
00607       std::cout<<"domain "<<i<<":";
00608       for(entry_iterator entry = d->begin(); entry != d->end(); ++entry){
00609         std::cout<<" "<<*entry;
00610       }
00611       std::cout<<std::endl;
00612     }
00613     
00614 #endif
00615     
00616     // Create a row to subdomain mapping
00617     rowtodomain_vector rowToDomain(mat.N());
00618 
00619     size_type domainId=0;
00620     
00621     for(DomainIterator domain=sd.begin(); domain != sd.end(); ++domain, ++domainId){
00622       typedef typename subdomain_type::const_iterator iterator;
00623       for(iterator row=domain->begin(); row != domain->end(); ++row)
00624         rowToDomain[*row].push_back(domainId);
00625     }
00626 
00627     initialize(rowToDomain, mat);
00628   }
00629 
00636   template<class M>
00637   struct SeqOverlappingSchwarzDomainSize {};
00638   
00639   template<typename T, typename A, int n, int m>
00640   struct SeqOverlappingSchwarzDomainSize<BCRSMatrix<FieldMatrix<T,n,m>,A > >
00641   {
00642     template<class Domain>
00643     static int size(const Domain & d)
00644     {
00645       assert(n==m);
00646       return m*d.size();
00647     }
00648   };
00649   
00650   template<class M, class X, class TM, bool onTheFly, class TA>
00651   void SeqOverlappingSchwarz<M,X,TM,onTheFly,TA>::initialize(const rowtodomain_vector& rowToDomain, const matrix_type& mat)
00652   {
00653     typedef typename std::vector<SuperMatrixInitializer<matrix_type> >::iterator InitializerIterator;
00654     typedef typename subdomain_vector::const_iterator DomainIterator;
00655     typedef typename slu_vector::iterator SolverIterator;
00656     
00657     if(onTheFly){
00658       maxlength=0;
00659       for(DomainIterator domain=subDomains.begin();domain!=subDomains.end();++domain)
00660         maxlength=std::max(maxlength, domain->size());
00661       maxlength*=mat[0].begin()->N();
00662     }else{
00663       // initialize the initializers
00664       DomainIterator domain=subDomains.begin();
00665     
00666       // Create the initializers list.
00667       std::vector<SuperMatrixInitializer<matrix_type> > initializers(subDomains.size());    
00668       
00669       SolverIterator solver=solvers.begin();
00670       for(InitializerIterator initializer=initializers.begin(); initializer!=initializers.end();
00671           ++initializer, ++solver, ++domain){
00672         solver->mat.N_=SeqOverlappingSchwarzDomainSize<M>::size(*domain);
00673         solver->mat.M_=SeqOverlappingSchwarzDomainSize<M>::size(*domain);
00674         //solver->setVerbosity(true);
00675         *initializer=SuperMatrixInitializer<matrix_type>(solver->mat);
00676       }
00677 
00678       // Set up the supermatrices according to the subdomains
00679       typedef OverlappingSchwarzInitializer<std::vector<SuperMatrixInitializer<matrix_type> >,
00680         rowtodomain_vector > Initializer;
00681 
00682       Initializer initializer(initializers, rowToDomain);
00683       copyToSuperMatrix(initializer, mat);
00684       if(solvers.size()==1)
00685         assert(solvers[0].mat==mat);
00686     
00687       /*    for(SolverIterator solver=solvers.begin(); solver!=solvers.end(); ++solver)
00688             dPrint_CompCol_Matrix("superlu", &static_cast<SuperMatrix&>(solver->mat)); */
00689       
00690       // Calculate the LU decompositions
00691       std::for_each(solvers.begin(), solvers.end(), std::mem_fun_ref(&slu::decompose));
00692       maxlength=0;
00693       for(SolverIterator solver=solvers.begin(); solver!=solvers.end(); ++solver){
00694         assert(solver->mat.N()==solver->mat.M());
00695         maxlength=std::max(maxlength, solver->mat.N());
00696         //writeCompColMatrixToMatlab(static_cast<SuperLUMatrix<M>&>(solver->mat), std::cout);
00697       }
00698     }
00699   }
00700   
00701   template<class M, class X, class TM, bool onTheFly, class TA>
00702   void SeqOverlappingSchwarz<M,X,TM,onTheFly,TA>::apply(X& x, const X& b)
00703   {
00704     Applier<SeqOverlappingSchwarz>::apply(*this, x, b);
00705   }
00706   
00707   template<class M, class X, class TM, bool onTheFly, class TA>
00708   template<bool forward>
00709   void SeqOverlappingSchwarz<M,X,TM,onTheFly,TA>::apply(X& x, const X& b)
00710   {
00711     typedef typename X::block_type block;
00712     typedef slu_vector solver_vector;
00713     typedef typename IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::solver_iterator iterator;
00714     typedef typename IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::domain_iterator
00715       domain_iterator;
00716     
00717     field_type *lhs=new field_type[maxlength];
00718     field_type *rhs=new field_type[maxlength];
00719     for(size_type i=0; i<maxlength; ++i)
00720       lhs[i]=0;
00721 
00722     
00723     domain_iterator domain=IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::begin(subDomains);
00724     iterator solver = IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::begin(solvers);
00725     X v(x); // temporary for the update
00726     v=0;
00727     
00728     typedef typename AdderSelector<TM,X>::Adder Adder;
00729     for(;domain != IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::end(subDomains); ++domain){
00730       //Copy rhs to C-array for SuperLU
00731       std::for_each(domain->begin(), domain->end(), Assigner<X>(mat, rhs, b, x));
00732       if(onTheFly){
00733         // Create the subdomain solver
00734         slu sdsolver;
00735         sdsolver.setSubMatrix(mat, *domain);
00736         // Apply
00737         sdsolver.apply(lhs,rhs);
00738         }else
00739         {
00740           solver->apply(lhs, rhs);
00741           ++solver;
00742           
00743         }
00744       //Add relaxed correction to from SuperLU to v
00745       std::for_each(domain->begin(), domain->end(), Adder(v, x, lhs, relax));
00746       
00747     }
00748     
00749     Adder adder(v, x, lhs, relax);
00750     adder.axpy();
00751     delete[] lhs;
00752     delete[] rhs;
00753     
00754   }   
00755     
00756   template<class M, class X, class TM, bool onTheFly, class TA>
00757   template<typename T, typename A, int n>
00758   SeqOverlappingSchwarz<M,X,TM,onTheFly,TA>::Assigner<BlockVector<FieldVector<T,n>,A> >::Assigner(const M& mat_, T* rhs_, const BlockVector<FieldVector<T,n>,A>& b_,
00759                                                                                       const BlockVector<FieldVector<T,n>,A>& x_)
00760     : mat(&mat_), rhs(rhs_), b(&b_), x(&x_), i(0)
00761   {}
00762 
00763   template<class M, class X, class TM, bool onTheFly, class TA>
00764   template<typename T, typename A, int n>
00765   void SeqOverlappingSchwarz<M,X,TM,onTheFly,TA>::Assigner<BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
00766   {
00767     // The current index.
00768     size_type starti = i;
00769     
00770     //assign right hand side of current domainindex block
00771     // rhs is an array of doubles!
00772     // rhs[starti] = b[domainindex]
00773     for(size_type j=0; j<n; ++j, ++i)
00774       rhs[i]=(*b)[domainIndex][j];
00775     
00776     // loop over all Matrix row entries and calculate defect.
00777     typedef typename M::ConstColIterator col_iterator;
00778     typedef typename subdomain_type::const_iterator domain_iterator;
00779 
00780     // calculate defect for current row index block
00781     for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
00782         typename X::block_type tmp;
00783         // T = A(domainIndex,k)*x(k)
00784         (*col).mv((*x)[col.index()], tmp);
00785         i-=n;
00786         //rhs[starti] -= tmp
00787         for(size_type j=0; j<n; ++j, ++i)
00788           rhs[i]-=tmp[j];
00789     }
00790     assert(starti+static_cast<size_type>(n)==i);
00791   }
00792   namespace
00793   {
00794     
00795   template<typename T, typename A, int n>
00796   AdditiveAdder<BlockVector<FieldVector<T,n>,A> >::AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v_, 
00797                  BlockVector<FieldVector<T,n>,A>& x_, const T* t_, const T& relax_)
00798     : t(t_), v(&v_), x(&x_), relax(relax_), i(0)
00799   {}
00800 
00801   template<typename T, typename A, int n>
00802   void AdditiveAdder<BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
00803   {
00804     for(size_type j=0; j<n; ++j, ++i)
00805       (*v)[domainIndex][j]+=t[i];
00806   }
00807 
00808  
00809   template<typename T, typename A, int n>
00810   void AdditiveAdder<BlockVector<FieldVector<T,n>,A> >::axpy()
00811   {
00812     x->axpy(relax,*v);
00813   }
00814 
00815  
00816   template<typename T, typename A, int n>
00817   MultiplicativeAdder<BlockVector<FieldVector<T,n>,A> >
00818   ::MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v_, 
00819                         BlockVector<FieldVector<T,n>,A>& x_, const T* t_, const T& relax_)
00820     : t(t_), v(&v_), x(&x_), relax(relax_), i(0)
00821   {}
00822 
00823   
00824   template<typename T, typename A, int n>
00825   void MultiplicativeAdder<BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
00826   {
00827     for(size_type j=0; j<n; ++j, ++i)
00828       (*x)[domainIndex][j]+=relax*t[i];
00829   }
00830 
00831   
00832   template<typename T, typename A, int n>
00833   void MultiplicativeAdder<BlockVector<FieldVector<T,n>,A> >::axpy()
00834   {
00835   }
00836   }
00837   
00839 #endif  
00840 }
00841 
00842 #endif
Generated on Sat Apr 24 11:13:47 2010 for dune-istl by  doxygen 1.6.3