overlappingschwarz.hh

Go to the documentation of this file.
00001 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
00002 // vi: set et ts=8 sw=2 sts=2:
00003 #ifndef DUNE_OVERLAPPINGSCHWARZ_HH
00004 #define DUNE_OVERLAPPINGSCHWARZ_HH
00005 #include<cassert>
00006 #include<algorithm>
00007 #include<functional>
00008 #include<vector>
00009 #include<set>
00010 #include<dune/common/dynmatrix.hh>
00011 #include<dune/common/sllist.hh>
00012 #include"preconditioners.hh"
00013 #include"superlu.hh"
00014 #include"bvector.hh"
00015 #include"bcrsmatrix.hh"
00016 #include"ilusubdomainsolver.hh"
00017 
00018 namespace Dune
00019 {
00020 
00032   template<class M, class X, class TM, class TD, class TA>
00033   class SeqOverlappingSchwarz;
00034   
00038   template<class I, class S, class D>
00039   class OverlappingSchwarzInitializer
00040   {
00041   public:
00043     typedef D subdomain_vector;
00044 
00045     typedef I InitializerList;
00046     typedef typename InitializerList::value_type AtomInitializer;
00047     typedef typename AtomInitializer::Matrix Matrix;
00048     typedef typename Matrix::const_iterator Iter;
00049     typedef typename Matrix::row_type::const_iterator CIter;
00050     
00051     typedef S IndexSet;
00052     typedef typename IndexSet::size_type size_type;
00053     
00054     OverlappingSchwarzInitializer(InitializerList& il,
00055                                   const IndexSet& indices,
00056                                   const subdomain_vector& domains);
00057     
00058     
00059     void addRowNnz(const Iter& row);
00060     
00061     void allocate();
00062         
00063     void countEntries(const Iter& row, const CIter& col) const;
00064 
00065     void calcColstart() const;
00066     
00067     void copyValue(const Iter& row, const CIter& col) const;
00068     
00069     void createMatrix() const;
00070   private:
00071     class IndexMap
00072     {
00073     public:
00074       typedef typename S::size_type size_type;
00075       typedef std::map<size_type,size_type> Map;      
00076       typedef typename Map::iterator iterator;
00077       typedef typename Map::const_iterator const_iterator;
00078       
00079       IndexMap();
00080       
00081       void insert(size_type grow);
00082       
00083       const_iterator find(size_type grow) const;
00084       
00085       iterator find(size_type grow);
00086 
00087       iterator begin();
00088       
00089       const_iterator begin()const;
00090 
00091       iterator end();
00092       
00093       const_iterator end() const;
00094             
00095     private:
00096       std::map<size_type,size_type> map_;
00097       size_type row;
00098     };
00099     
00100     
00101     typedef typename InitializerList::iterator InitIterator;
00102     typedef typename IndexSet::const_iterator IndexIteratur;
00103     InitializerList* initializers;
00104     const IndexSet *indices;
00105     mutable std::vector<IndexMap> indexMaps;
00106     const subdomain_vector& domains;
00107   };
00108 
00112   struct AdditiveSchwarzMode
00113   {};
00114 
00118   struct MultiplicativeSchwarzMode
00119   {};
00120   
00125   struct SymmetricMultiplicativeSchwarzMode
00126   {};  
00127 
00132   template<class M, class X, class Y>
00133   class DynamicMatrixSubdomainSolver;
00134 
00135   // Specialization for BCRSMatrix
00136   template<class K, int n, class Al, class X, class Y>
00137   class DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >
00138   {
00139     typedef BCRSMatrix< FieldMatrix<K,n,n>, Al> M;
00140   public:
00142     typedef typename Dune::remove_const<M>::type matrix_type;
00143     typedef K field_type;
00144     typedef typename Dune::remove_const<M>::type rilu_type;
00146     typedef X domain_type;
00148     typedef Y range_type;
00149     
00154     void apply (DynamicVector<field_type>& v, DynamicVector<field_type>& d)
00155     {
00156       assert(v.size() > 0);
00157       assert(v.size() == d.size());
00158       assert(A.rows() <= v.size());
00159       assert(A.cols() <= v.size());
00160       size_t sz = A.rows();
00161       v.resize(sz);
00162       d.resize(sz);
00163       A.solve(v,d);
00164       v.resize(v.capacity());
00165       d.resize(d.capacity());
00166     }
00167 
00175     template<class S>
00176     void setSubMatrix(const M& BCRS, S& rowset)
00177     {
00178       size_t sz = rowset.size();
00179       A.resize(sz*n,sz*n);
00180       typename DynamicMatrix<K>::RowIterator rIt = A.begin();
00181       typedef typename S::const_iterator SIter;
00182       size_t r = 0, c = 0;
00183       for(SIter rowIdx = rowset.begin(), rowEnd=rowset.end();
00184           rowIdx!= rowEnd; ++rowIdx, r++)
00185       {
00186         c = 0;
00187         for(SIter colIdx = rowset.begin(), colEnd=rowset.end(); 
00188             colIdx!= colEnd; ++colIdx, c++)
00189         {
00190           if (BCRS[*rowIdx].find(*colIdx) == BCRS[*rowIdx].end())
00191             continue;
00192           for (size_t i=0; i<n; i++)
00193           {
00194             for (size_t j=0; j<n; j++)
00195             {
00196               A[r*n+i][c*n+j] = BCRS[*rowIdx][*colIdx][i][j];
00197             }
00198           }
00199         }
00200       }
00201     }
00202   private:
00203     DynamicMatrix<K> A;
00204   };
00205 
00206   template<typename T>
00207   struct OverlappingAssigner
00208   {
00209   };
00210     
00211   // specialization for DynamicMatrix
00212   template<class K, int n, class Al, class X, class Y>
00213   class OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
00214   {
00215   public:
00216     typedef BCRSMatrix< FieldMatrix<K,n,n>, Al> matrix_type;
00217     typedef K field_type;
00218     typedef Y range_type;
00219     typedef typename range_type::block_type block_type;
00220     typedef typename matrix_type::size_type size_type;
00221     
00229     OverlappingAssigner(std::size_t maxlength, const BCRSMatrix<FieldMatrix<K,n,n>, Al>& mat_, const X& b_, Y& x_);
00230 
00234     inline
00235     void deallocate();
00236     
00237     /*
00238      * @brief Resets the local index to zero.
00239      */
00240     inline
00241     void resetIndexForNextDomain();
00242     
00247     inline
00248     DynamicVector<K> & lhs();
00249     
00254     inline
00255     DynamicVector<K> & rhs();
00256     
00261     inline
00262     void relaxResult(field_type relax);
00263     
00268     void operator()(const size_type& domainIndex);
00269     
00277     inline
00278     void assignResult(block_type& res);
00279     
00280   private:
00284     const matrix_type* mat;
00286     // we need a pointer, because we have to avoid deep copies
00287     DynamicVector<field_type> * rhs_;
00289     // we need a pointer, because we have to avoid deep copies
00290     DynamicVector<field_type> * lhs_;
00292     const range_type* b;
00294     range_type* x;
00296     std::size_t i;
00298     std::size_t maxlength_;
00299   };
00300 
00301 #if HAVE_SUPERLU
00302   template<typename T, typename A, int n, int m>
00303   struct OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
00304     {
00305       typedef BCRSMatrix<FieldMatrix<T,n,m>,A> matrix_type;
00306       typedef typename SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::range_type range_type;
00307       typedef typename range_type::field_type field_type;
00308       typedef typename range_type::block_type block_type;
00309       
00310       typedef typename matrix_type::size_type size_type;
00311       
00319       OverlappingAssigner(std::size_t maxlength, const BCRSMatrix<FieldMatrix<T,n,m>,A>& mat, 
00320                           const range_type& b, range_type& x);
00326       void deallocate();
00327       
00328       /*
00329        * @brief Resets the local index to zero.
00330        */
00331       void resetIndexForNextDomain();
00332       
00337       field_type* lhs();
00338             
00343       field_type* rhs();
00344       
00349       void relaxResult(field_type relax);
00350 
00355       void operator()(const size_type& domain);
00356 
00364       void assignResult(block_type& res);
00365       
00366     private:
00370       const matrix_type* mat;
00372       field_type* rhs_;
00374       field_type* lhs_;
00376       const range_type* b;
00378       range_type* x;
00380       std::size_t i;
00382       std::size_t maxlength_;
00383     };
00384 
00385 #endif
00386 
00387   template<class M, class X, class Y>
00388   class OverlappingAssignerILUBase
00389   {
00390   public:
00391     typedef M matrix_type;
00392     
00393     typedef typename M::field_type field_type;
00394   
00395     typedef typename Y::block_type block_type;
00396       
00397     typedef typename matrix_type::size_type size_type;
00405     OverlappingAssignerILUBase(std::size_t maxlength, const M& mat, 
00406                         const Y& b, X& x);
00412     void deallocate();
00413     
00414     /*
00415      * @brief Resets the local index to zero.
00416      */
00417     void resetIndexForNextDomain();
00418     
00419     /*
00420      * @brief Resets the local index to zero.
00421      */
00422     X& lhs();
00423     
00428     Y& rhs();
00429     
00434     void relaxResult(field_type relax);
00435     
00440     void operator()(const size_type& domain);
00441     
00449     void assignResult(block_type& res);
00450 
00451   private:
00455     const M* mat;
00457     X* lhs_;
00459     Y* rhs_;
00461     const Y* b;
00463     X* x;
00465     size_type i;
00466   };
00467   
00468   // specialization for ILU0
00469   template<class M, class X, class Y>
00470   class OverlappingAssigner<ILU0SubdomainSolver<M,X,Y> >
00471     : public OverlappingAssignerILUBase<M,X,Y>
00472   {
00473   public:
00481     OverlappingAssigner(std::size_t maxlength, const M& mat, 
00482                         const Y& b, X& x)
00483       : OverlappingAssignerILUBase<M,X,Y>(maxlength, mat,b,x)
00484     {}
00485   };
00486 
00487   // specialization for ILUN
00488   template<class M, class X, class Y>
00489   class OverlappingAssigner<ILUNSubdomainSolver<M,X,Y> >
00490     : public OverlappingAssignerILUBase<M,X,Y>
00491   {
00492   public:
00500     OverlappingAssigner(std::size_t maxlength, const M& mat, 
00501                         const Y& b, X& x)
00502       : OverlappingAssignerILUBase<M,X,Y>(maxlength, mat,b,x)
00503     {}
00504   };
00505 
00506     template<typename S, typename T>
00507     struct AdditiveAdder
00508     {
00509     };
00510     
00511     template<typename S, typename T, typename A, int n>
00512     struct AdditiveAdder<S, BlockVector<FieldVector<T,n>,A> >
00513     {
00514       typedef typename A::size_type size_type;
00515       AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x,
00516                     OverlappingAssigner<S>& assigner, const T& relax_);
00517       void operator()(const size_type& domain);
00518       void axpy();
00519       
00520     private:
00521       BlockVector<FieldVector<T,n>,A>* v;
00522       BlockVector<FieldVector<T,n>,A>* x;
00523       OverlappingAssigner<S>* assigner;
00524       T relax;
00525     };
00526 
00527     template<typename S,typename T>
00528     struct MultiplicativeAdder
00529     {
00530     };
00531     
00532     template<typename S, typename T, typename A, int n>
00533     struct MultiplicativeAdder<S, BlockVector<FieldVector<T,n>,A> >
00534     {
00535       typedef typename A::size_type size_type;
00536       MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x, 
00537                           OverlappingAssigner<S>& assigner_, const T& relax_);
00538       void operator()(const size_type& domain);
00539       void axpy();
00540       
00541     private:
00542       BlockVector<FieldVector<T,n>,A>* x;
00543       OverlappingAssigner<S>* assigner;
00544       T relax;
00545     };
00546 
00556     template<typename T, class X, class S>
00557     struct AdderSelector
00558     {};
00559     
00560   template<class X, class S>
00561     struct AdderSelector<AdditiveSchwarzMode,X,S>
00562     {
00563       typedef AdditiveAdder<S,X> Adder;
00564     };
00565     
00566   template<class X, class S>
00567     struct AdderSelector<MultiplicativeSchwarzMode,X,S>
00568     {
00569       typedef MultiplicativeAdder<S,X> Adder;
00570     };    
00571 
00572   template<class X, class S>
00573     struct AdderSelector<SymmetricMultiplicativeSchwarzMode,X,S>
00574     {
00575       typedef MultiplicativeAdder<S,X> Adder;
00576     };    
00577 
00589     template<typename T1, typename T2, bool forward>
00590     struct IteratorDirectionSelector
00591     {
00592       typedef T1 solver_vector;
00593       typedef typename solver_vector::iterator solver_iterator;
00594       typedef T2 subdomain_vector;
00595       typedef typename subdomain_vector::const_iterator domain_iterator;
00596     
00597       static solver_iterator begin(solver_vector& sv)
00598       {
00599         return sv.begin();
00600       }
00601       
00602       static solver_iterator end(solver_vector& sv)
00603       {
00604         return sv.end();
00605       }
00606       static domain_iterator begin(const subdomain_vector& sv)
00607       {
00608         return sv.begin();
00609       }
00610       
00611       static domain_iterator end(const subdomain_vector& sv)
00612       {
00613         return sv.end();
00614       }
00615     };
00616 
00617     template<typename T1, typename T2>
00618     struct IteratorDirectionSelector<T1,T2,false>
00619     {
00620       typedef T1 solver_vector;
00621       typedef typename solver_vector::reverse_iterator solver_iterator;
00622       typedef T2 subdomain_vector;
00623       typedef typename subdomain_vector::const_reverse_iterator domain_iterator;
00624     
00625       static solver_iterator begin(solver_vector& sv)
00626       {
00627         return sv.rbegin();
00628       }
00629       
00630       static solver_iterator end(solver_vector& sv)
00631       {
00632         return sv.rend();
00633       }
00634       static domain_iterator begin(const subdomain_vector& sv)
00635       {
00636         return sv.rbegin();
00637       }
00638       
00639       static domain_iterator end(const subdomain_vector& sv)
00640       {
00641         return sv.rend();
00642       }
00643     };
00644     
00653     template<class T>
00654     struct SeqOverlappingSchwarzApplier
00655     {
00656       typedef T smoother;
00657       typedef typename smoother::range_type range_type;
00658       
00659       static void apply(smoother& sm, range_type& v, const range_type& b)
00660       {
00661         sm.template apply<true>(v, b);
00662       }
00663     };
00664     
00665     template<class M, class X, class TD, class TA>
00666     struct SeqOverlappingSchwarzApplier<SeqOverlappingSchwarz<M,X,SymmetricMultiplicativeSchwarzMode,TD,TA> >
00667     {
00668       typedef SeqOverlappingSchwarz<M,X,SymmetricMultiplicativeSchwarzMode,TD,TA> smoother;
00669       typedef typename smoother::range_type range_type;
00670       
00671       static void apply(smoother& sm, range_type& v, const range_type& b)
00672       {
00673         sm.template apply<true>(v, b);
00674         sm.template apply<false>(v, b);
00675       }
00676     };
00677 
00678   template<class T>
00679   struct SeqOverlappingSchwarzAssembler
00680   {};
00681   
00682 
00683   template<class K, int n, class Al, class X, class Y>
00684   struct SeqOverlappingSchwarzAssembler< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
00685   {
00686     typedef BCRSMatrix< FieldMatrix<K,n,n>, Al> matrix_type;
00687     template<class RowToDomain, class Solvers, class SubDomains>
00688     static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
00689                                              Solvers& solvers, const SubDomains& domains,
00690                                              bool onTheFly);
00691   };
00692 
00693 #if HAVE_SUPERLU
00694   template<class T>
00695   struct SeqOverlappingSchwarzAssembler<SuperLU<T> >
00696   {
00697     typedef T matrix_type;
00698     template<class RowToDomain, class Solvers, class SubDomains>
00699     static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
00700                                              Solvers& solvers, const SubDomains& domains,
00701                                              bool onTheFly);
00702   };
00703 #endif
00704 
00705   template<class M,class X, class Y>
00706   struct SeqOverlappingSchwarzAssemblerILUBase
00707   {
00708     typedef M matrix_type;
00709     template<class RowToDomain, class Solvers, class SubDomains>
00710     static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
00711                                              Solvers& solvers, const SubDomains& domains,
00712                                              bool onTheFly);
00713   };
00714 
00715   template<class M,class X, class Y>
00716   struct SeqOverlappingSchwarzAssembler<ILU0SubdomainSolver<M,X,Y> >
00717     : public SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>
00718   {};
00719 
00720   template<class M,class X, class Y>
00721   struct SeqOverlappingSchwarzAssembler<ILUNSubdomainSolver<M,X,Y> >
00722     : public SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>
00723   {};
00724   
00735   template<class M, class X, class TM=AdditiveSchwarzMode, 
00736            class TD=ILU0SubdomainSolver<M,X,X>, class TA=std::allocator<X> >
00737   class SeqOverlappingSchwarz
00738     : public Preconditioner<X,X>
00739   {
00740   public:
00744     typedef M matrix_type;
00745 
00749     typedef X domain_type;
00750 
00754     typedef X range_type;
00755 
00762     typedef TM Mode;
00763 
00767     typedef typename X::field_type field_type;
00768 
00770     typedef typename matrix_type::size_type size_type;
00771 
00773     typedef TA allocator;
00774         
00776     typedef std::set<size_type, std::less<size_type>, 
00777                      typename TA::template rebind<size_type>::other> 
00778       subdomain_type;
00779 
00781     typedef std::vector<subdomain_type, typename TA::template rebind<subdomain_type>::other> subdomain_vector;
00782 
00784     typedef SLList<size_type, typename TA::template rebind<size_type>::other> subdomain_list;
00785       
00787     typedef std::vector<subdomain_list, typename TA::template rebind<subdomain_list>::other > rowtodomain_vector;
00788 
00790     typedef TD slu;
00791     
00793     typedef std::vector<slu, typename TA::template rebind<slu>::other> slu_vector;
00794 
00795     enum{
00797       category = SolverCategory::sequential
00798         };
00799     
00813     SeqOverlappingSchwarz(const matrix_type& mat, const subdomain_vector& subDomains,
00814                           field_type relaxationFactor=1, bool onTheFly_=true);
00815 
00827     SeqOverlappingSchwarz(const matrix_type& mat, const rowtodomain_vector& rowToDomain,
00828                           field_type relaxationFactor=1, bool onTheFly_=true);
00829                           
00835     virtual void pre (X& x, X& b) {}
00836 
00842     virtual void apply (X& v, const X& d);
00843 
00844     template<bool forward>
00845     void apply(X& v, const X& d);
00846     
00852     virtual void post (X& x) {
00853       Dune::dverb<<" avg nnz over subdomain is "<<nnz<<std::endl;
00854     }
00855   
00856   private:
00857     const M& mat;
00858     slu_vector solvers;
00859     subdomain_vector subDomains;
00860     field_type relax;
00861     
00862     typename M::size_type maxlength;
00863     std::size_t nnz;
00864 
00865     bool onTheFly;
00866   };
00867   
00868 
00869     
00870   template<class I, class S, class D>
00871   OverlappingSchwarzInitializer<I,S,D>::OverlappingSchwarzInitializer(InitializerList& il,
00872                                                                     const IndexSet& idx,
00873                                                                     const subdomain_vector& domains_)
00874     : initializers(&il), indices(&idx), indexMaps(il.size()), domains(domains_)
00875   {
00876   }
00877   
00878   
00879   template<class I, class S, class D>
00880   void OverlappingSchwarzInitializer<I,S,D>::addRowNnz(const Iter& row)
00881   {
00882     typedef typename IndexSet::value_type::const_iterator iterator;
00883     for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain){
00884       (*initializers)[*domain].addRowNnz(row, domains[*domain]);
00885       indexMaps[*domain].insert(row.index());
00886     }
00887   }
00888   
00889   template<class I, class S, class D>
00890   void OverlappingSchwarzInitializer<I,S,D>::allocate()
00891   {
00892     std::for_each(initializers->begin(), initializers->end(),
00893                   std::mem_fun_ref(&AtomInitializer::allocateMatrixStorage));
00894     std::for_each(initializers->begin(), initializers->end(), 
00895                   std::mem_fun_ref(&AtomInitializer::allocateMarker));
00896   }
00897   
00898   template<class I, class S, class D>
00899   void OverlappingSchwarzInitializer<I,S,D>::countEntries(const Iter& row, const CIter& col) const
00900   {
00901     typedef typename IndexSet::value_type::const_iterator iterator;
00902     for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain){
00903       typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
00904       if(v!= indexMaps[*domain].end()){
00905         (*initializers)[*domain].countEntries(indexMaps[*domain].find(col.index())->second);
00906       }
00907     }
00908   }
00909 
00910   template<class I, class S, class D>
00911   void OverlappingSchwarzInitializer<I,S,D>::calcColstart() const
00912   {
00913     std::for_each(initializers->begin(), initializers->end(),
00914                   std::mem_fun_ref(&AtomInitializer::calcColstart));
00915   }
00916 
00917   template<class I, class S, class D>
00918   void OverlappingSchwarzInitializer<I,S,D>::copyValue(const Iter& row, const CIter& col) const
00919   {
00920     typedef typename IndexSet::value_type::const_iterator iterator;
00921     for(iterator domain=(*indices)[row.index()].begin(); domain!= (*indices)[row.index()].end(); ++domain){
00922       typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
00923       if(v!= indexMaps[*domain].end()){
00924         assert(indexMaps[*domain].end()!=indexMaps[*domain].find(row.index()));
00925         (*initializers)[*domain].copyValue(col, indexMaps[*domain].find(row.index())->second, 
00926                                            v->second);
00927       }
00928     }
00929   }
00930     
00931   template<class I, class S, class D>
00932   void OverlappingSchwarzInitializer<I,S,D>::createMatrix() const
00933   {
00934     indexMaps.clear();
00935     indexMaps.swap(std::vector<IndexMap>(indexMaps));
00936     std::for_each(initializers->begin(), initializers->end(),
00937                   std::mem_fun_ref(&AtomInitializer::createMatrix));
00938   }
00939 
00940   template<class I, class S, class D>
00941   OverlappingSchwarzInitializer<I,S,D>::IndexMap::IndexMap()
00942     : row(0)
00943   {}
00944 
00945   template<class I, class S, class D>
00946   void OverlappingSchwarzInitializer<I,S,D>::IndexMap::insert(size_type grow)
00947   {
00948     assert(map_.find(grow)==map_.end());
00949     map_.insert(std::make_pair(grow, row++));
00950   }
00951 
00952   template<class I, class S, class D>
00953   typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator 
00954   OverlappingSchwarzInitializer<I,S,D>::IndexMap::find(size_type grow) const
00955   {
00956     return map_.find(grow);
00957   }
00958   
00959   template<class I, class S, class D>
00960   typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator 
00961   OverlappingSchwarzInitializer<I,S,D>::IndexMap::find(size_type grow)
00962   {
00963     return map_.find(grow);
00964   }
00965 
00966   template<class I, class S, class D>
00967   typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator 
00968   OverlappingSchwarzInitializer<I,S,D>::IndexMap::end() const
00969   {
00970     return map_.end();
00971   }
00972   
00973   template<class I, class S, class D>
00974   typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator 
00975   OverlappingSchwarzInitializer<I,S,D>::IndexMap::end()
00976   {
00977     return map_.end();
00978   }
00979 
00980   template<class I, class S, class D>
00981   typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator 
00982   OverlappingSchwarzInitializer<I,S,D>::IndexMap::begin() const
00983   {
00984     return map_.begin();
00985   }
00986   
00987   template<class I, class S, class D>
00988   typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator 
00989   OverlappingSchwarzInitializer<I,S,D>::IndexMap::begin()
00990   {
00991     return map_.begin();
00992   }
00993 
00994   template<class M, class X, class TM, class TD, class TA>
00995   SeqOverlappingSchwarz<M,X,TM,TD,TA>::SeqOverlappingSchwarz(const matrix_type& mat_, const rowtodomain_vector& rowToDomain,
00996                                                           field_type relaxationFactor, bool fly)
00997     : mat(mat_), relax(relaxationFactor), onTheFly(fly)
00998   {
00999     typedef typename rowtodomain_vector::const_iterator RowDomainIterator;
01000     typedef typename subdomain_list::const_iterator DomainIterator;
01001 #ifdef DUNE_ISTL_WITH_CHECKING
01002     assert(rowToDomain.size()==mat.N());
01003     assert(rowToDomain.size()==mat.M());
01004     
01005     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
01006       assert(iter->size()>0);
01007     
01008 #endif
01009     // calculate the number of domains
01010     size_type domains=0;
01011     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
01012       for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
01013         domains=std::max(domains, *d);
01014     ++domains;
01015     
01016     solvers.resize(domains);
01017     subDomains.resize(domains);
01018     
01019     // initialize subdomains to row mapping from row to subdomain mapping
01020     size_type row=0;
01021     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter, ++row)
01022       for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
01023         subDomains[*d].insert(row);
01024 
01025 #ifdef DUNE_ISTL_WITH_CHECKING
01026     size_type i=0;
01027     typedef typename subdomain_vector::const_iterator iterator;
01028     for(iterator iter=subDomains.begin(); iter != subDomains.end(); ++iter){
01029       typedef typename subdomain_type::const_iterator entry_iterator;
01030       Dune::dvverb<<"domain "<<i++<<":";
01031       for(entry_iterator entry = iter->begin(); entry != iter->end(); ++entry){
01032         Dune::dvverb<<" "<<*entry;
01033       }
01034       Dune::dvverb<<std::endl;
01035     }
01036 #endif
01037     maxlength = SeqOverlappingSchwarzAssembler<slu>
01038       ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
01039   }
01040   
01041   template<class M, class X, class TM, class TD, class TA>
01042   SeqOverlappingSchwarz<M,X,TM,TD,TA>::SeqOverlappingSchwarz(const matrix_type& mat_,
01043                                                           const subdomain_vector& sd,
01044                                                           field_type relaxationFactor,
01045                                                           bool fly)
01046     :  mat(mat_), solvers(sd.size()), subDomains(sd), relax(relaxationFactor),
01047        onTheFly(fly)
01048   {
01049     typedef typename subdomain_vector::const_iterator DomainIterator;
01050 
01051 #ifdef DUNE_ISTL_WITH_CHECKING
01052     size_type i=0;
01053     
01054     for(DomainIterator d=sd.begin(); d != sd.end(); ++d,++i){
01055       //std::cout<<i<<": "<<d->size()<<std::endl;
01056       assert(d->size()>0);
01057       typedef typename DomainIterator::value_type::const_iterator entry_iterator;
01058       Dune::dvverb<<"domain "<<i<<":";
01059       for(entry_iterator entry = d->begin(); entry != d->end(); ++entry){
01060         Dune::dvverb<<" "<<*entry;
01061       }
01062       Dune::dvverb<<std::endl;
01063     }
01064     
01065 #endif
01066     
01067     // Create a row to subdomain mapping
01068     rowtodomain_vector rowToDomain(mat.N());
01069 
01070     size_type domainId=0;
01071     
01072     for(DomainIterator domain=sd.begin(); domain != sd.end(); ++domain, ++domainId){
01073       typedef typename subdomain_type::const_iterator iterator;
01074       for(iterator row=domain->begin(); row != domain->end(); ++row)
01075         rowToDomain[*row].push_back(domainId);
01076     }
01077 
01078     maxlength = SeqOverlappingSchwarzAssembler<slu>
01079       ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
01080   }
01081 
01088   template<class M>
01089   struct SeqOverlappingSchwarzDomainSize {};
01090   
01091   template<typename T, typename A, int n, int m>
01092   struct SeqOverlappingSchwarzDomainSize<BCRSMatrix<FieldMatrix<T,n,m>,A > >
01093   {
01094     template<class Domain>
01095     static int size(const Domain & d)
01096     {
01097       assert(n==m);
01098       return m*d.size();
01099     }
01100   };
01101   
01102   template<class K, int n, class Al, class X, class Y>
01103   template<class RowToDomain, class Solvers, class SubDomains>
01104   std::size_t
01105   SeqOverlappingSchwarzAssembler< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >::
01106   assembleLocalProblems(const RowToDomain& rowToDomain, 
01107     const matrix_type& mat,
01108     Solvers& solvers,
01109     const SubDomains& subDomains,
01110     bool onTheFly)
01111   {
01112     typedef typename SubDomains::const_iterator DomainIterator;
01113     std::size_t maxlength = 0;
01114     
01115     assert(onTheFly);
01116     
01117     for(DomainIterator domain=subDomains.begin();domain!=subDomains.end();++domain)
01118       maxlength=std::max(maxlength, domain->size());
01119     maxlength*=n;
01120     
01121     return maxlength;
01122   }
01123 
01124 #if HAVE_SUPERLU
01125   template<class T>
01126   template<class RowToDomain, class Solvers, class SubDomains>
01127   std::size_t SeqOverlappingSchwarzAssembler<SuperLU<T> >::assembleLocalProblems(const RowToDomain& rowToDomain, 
01128                                                                                 const matrix_type& mat,
01129                                                                                 Solvers& solvers,
01130                                                                                 const SubDomains& subDomains,
01131                                                                                 bool onTheFly)
01132   {
01133     typedef typename std::vector<SuperMatrixInitializer<matrix_type> >::iterator InitializerIterator;
01134     typedef typename SubDomains::const_iterator DomainIterator;
01135     typedef typename Solvers::iterator SolverIterator;
01136     std::size_t maxlength = 0;
01137     
01138     if(onTheFly){
01139       for(DomainIterator domain=subDomains.begin();domain!=subDomains.end();++domain)
01140         maxlength=std::max(maxlength, domain->size());
01141       maxlength*=mat[0].begin()->N();
01142     }else{
01143       // initialize the initializers
01144       DomainIterator domain=subDomains.begin();
01145     
01146       // Create the initializers list.
01147       std::vector<SuperMatrixInitializer<matrix_type> > initializers(subDomains.size());    
01148       
01149       SolverIterator solver=solvers.begin();
01150       for(InitializerIterator initializer=initializers.begin(); initializer!=initializers.end();
01151           ++initializer, ++solver, ++domain){
01152         solver->mat.N_=SeqOverlappingSchwarzDomainSize<matrix_type>::size(*domain);
01153         solver->mat.M_=SeqOverlappingSchwarzDomainSize<matrix_type>::size(*domain);
01154         //solver->setVerbosity(true);
01155         *initializer=SuperMatrixInitializer<matrix_type>(solver->mat);
01156       }
01157 
01158       // Set up the supermatrices according to the subdomains
01159       typedef OverlappingSchwarzInitializer<std::vector<SuperMatrixInitializer<matrix_type> >,
01160         RowToDomain, SubDomains> Initializer;
01161 
01162       Initializer initializer(initializers, rowToDomain, subDomains);
01163       copyToSuperMatrix(initializer, mat);
01164       if(solvers.size()==1)
01165         assert(solvers[0].mat==mat);
01166     
01167       /*    for(SolverIterator solver=solvers.begin(); solver!=solvers.end(); ++solver)
01168             dPrint_CompCol_Matrix("superlu", &static_cast<SuperMatrix&>(solver->mat)); */
01169       
01170       // Calculate the LU decompositions
01171       std::for_each(solvers.begin(), solvers.end(), std::mem_fun_ref(&SuperLU<matrix_type>::decompose));
01172       for(SolverIterator solver=solvers.begin(); solver!=solvers.end(); ++solver){
01173         assert(solver->mat.N()==solver->mat.M());
01174         maxlength=std::max(maxlength, solver->mat.N());
01175         //writeCompColMatrixToMatlab(static_cast<SuperLUMatrix<M>&>(solver->mat), std::cout);
01176       }
01177     }
01178     return maxlength;
01179   }
01180 
01181 #endif
01182 
01183   template<class M,class X,class Y>
01184   template<class RowToDomain, class Solvers, class SubDomains>
01185   std::size_t SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>::assembleLocalProblems(const RowToDomain& rowToDomain, 
01186                                                                                                  const matrix_type& mat,
01187                                                                                                  Solvers& solvers,
01188                                                                                                  const SubDomains& subDomains,
01189                                                                                                  bool onTheFly)
01190   {
01191     typedef typename SubDomains::const_iterator DomainIterator;
01192     typedef typename Solvers::iterator SolverIterator;
01193     std::size_t maxlength = 0;
01194     
01195     if(onTheFly){
01196       for(DomainIterator domain=subDomains.begin();domain!=subDomains.end();++domain)
01197         maxlength=std::max(maxlength, domain->size());
01198     }else{
01199       // initialize the solvers of the local prolems.
01200       SolverIterator solver=solvers.begin();
01201       for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end();
01202           ++domain, ++solver){
01203         solver->setSubMatrix(mat, *domain);
01204         maxlength=std::max(maxlength, domain->size());
01205       }
01206     }
01207     
01208     return maxlength;
01209       
01210   }
01211 
01212   
01213   template<class M, class X, class TM, class TD, class TA>
01214   void SeqOverlappingSchwarz<M,X,TM,TD,TA>::apply(X& x, const X& b)
01215   {
01216     SeqOverlappingSchwarzApplier<SeqOverlappingSchwarz>::apply(*this, x, b);
01217   }
01218   
01219   template<class M, class X, class TM, class TD, class TA>
01220   template<bool forward>
01221   void SeqOverlappingSchwarz<M,X,TM,TD,TA>::apply(X& x, const X& b)
01222   {
01223     typedef typename X::block_type block;
01224     typedef slu_vector solver_vector;
01225     typedef typename IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::solver_iterator iterator;
01226     typedef typename IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::domain_iterator
01227       domain_iterator;
01228     
01229     OverlappingAssigner<TD> assigner(maxlength, mat, b, x);
01230         
01231     domain_iterator domain=IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::begin(subDomains);
01232     iterator solver = IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::begin(solvers);
01233     X v(x); // temporary for the update
01234     v=0;
01235     
01236     typedef typename AdderSelector<TM,X,TD >::Adder Adder;    
01237     Adder adder(v, x, assigner, relax);
01238     
01239     nnz=0;
01240     std::size_t no=0;
01241     for(;domain != IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::end(subDomains); ++domain){
01242       //Copy rhs to C-array for SuperLU
01243       std::for_each(domain->begin(), domain->end(), assigner);
01244       assigner.resetIndexForNextDomain();
01245       if(onTheFly){
01246         // Create the subdomain solver
01247         slu sdsolver;
01248         sdsolver.setSubMatrix(mat, *domain);
01249         // Apply
01250         sdsolver.apply(assigner.lhs(), assigner.rhs());
01251         //nnz+=sdsolver.nnz();
01252       }else{
01253         solver->apply(assigner.lhs(), assigner.rhs());
01254         //nnz+=solver->nnz();
01255         ++solver;
01256       }
01257       ++no;
01258       //Add relaxed correction to from SuperLU to v
01259       std::for_each(domain->begin(), domain->end(), adder);
01260       assigner.resetIndexForNextDomain();
01261       
01262     }
01263     nnz/=no;
01264     
01265     adder.axpy();
01266     assigner.deallocate();
01267   }
01268   
01269   template<class K, int n, class Al, class X, class Y>
01270   OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
01271   ::OverlappingAssigner(std::size_t maxlength, const BCRSMatrix<FieldMatrix<K,n,n>, Al>& mat_, 
01272     const X& b_, Y& x_) :
01273     mat(&mat_), 
01274     rhs_( new DynamicVector<field_type>(maxlength, 42) ),
01275     lhs_( new DynamicVector<field_type>(maxlength, -42) ),
01276     b(&b_),
01277     x(&x_),
01278     i(0),
01279     maxlength_(maxlength)
01280   {}
01281 
01282   template<class K, int n, class Al, class X, class Y>
01283   void
01284   OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
01285   ::deallocate()
01286   {
01287     delete rhs_;
01288     delete lhs_;
01289   }
01290     
01291   template<class K, int n, class Al, class X, class Y>
01292   void
01293   OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
01294   ::resetIndexForNextDomain()
01295   {
01296     i=0;
01297   }
01298   
01299   template<class K, int n, class Al, class X, class Y>
01300   DynamicVector<K> &
01301   OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
01302   ::lhs()
01303   {
01304     return *lhs_;
01305   }
01306     
01307   template<class K, int n, class Al, class X, class Y>
01308   DynamicVector<K> &   
01309   OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
01310   ::rhs()
01311   {
01312     return *rhs_;
01313   }
01314     
01315   template<class K, int n, class Al, class X, class Y>
01316   void
01317   OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
01318   ::relaxResult(field_type relax)
01319   {
01320     lhs() *= relax;
01321   }
01322   
01323   template<class K, int n, class Al, class X, class Y>
01324   void
01325   OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
01326   ::operator()(const size_type& domainIndex)
01327   {
01328     lhs() = 0.0;
01329 #if 0
01330     //assign right hand side of current domainindex block
01331     for(size_type j=0; j<n; ++j, ++i){
01332       assert(i<maxlength_);
01333       rhs()[i]=(*b)[domainIndex][j];
01334     }
01335     
01336     // loop over all Matrix row entries and calculate defect.
01337     typedef typename matrix_type::ConstColIterator col_iterator;
01338     
01339     // calculate defect for current row index block
01340     for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
01341       block_type tmp(0.0);
01342       (*col).mv((*x)[col.index()], tmp);
01343       i-=n;
01344       for(size_type j=0; j<n; ++j, ++i){
01345         assert(i<maxlength_);
01346         rhs()[i]-=tmp[j];
01347       }
01348     }
01349 #else
01350     //assign right hand side of current domainindex block
01351     for(size_type j=0; j<n; ++j, ++i){
01352       assert(i<maxlength_);
01353       rhs()[i]=(*b)[domainIndex][j];
01354       
01355       // loop over all Matrix row entries and calculate defect.
01356       typedef typename matrix_type::ConstColIterator col_iterator;
01357       
01358       // calculate defect for current row index block
01359       for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
01360         for(size_type k=0; k<n; ++k){
01361           rhs()[i]-=(*col)[j][k] * (*x)[col.index()][k];
01362         }
01363       }
01364     }
01365 #endif 
01366   }
01367 
01368   template<class K, int n, class Al, class X, class Y>
01369   void
01370   OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
01371   ::assignResult(block_type& res)
01372   {
01373     // assign the result of the local solve to the global vector
01374     for(size_type j=0; j<n; ++j, ++i){
01375       assert(i<maxlength_);
01376       res[j]+=lhs()[i];
01377     }
01378   }
01379 
01380 #if HAVE_SUPERLU
01381 
01382   template<typename T, typename A, int n, int m>
01383   OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
01384   ::OverlappingAssigner(std::size_t maxlength, 
01385                         const BCRSMatrix<FieldMatrix<T,n,m>,A>& mat_, 
01386                         const range_type& b_,
01387                         range_type& x_)
01388     : mat(&mat_), 
01389       b(&b_), 
01390       x(&x_), i(0), maxlength_(maxlength)
01391   {
01392     rhs_ = new field_type[maxlength];
01393     lhs_ = new field_type[maxlength];
01394     
01395   }
01396   
01397   template<typename T, typename A, int n, int m>
01398   void OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::deallocate()
01399   {
01400     delete[] rhs_;
01401     delete[] lhs_;
01402   }
01403   
01404   template<typename T, typename A, int n, int m>
01405   void OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::operator()(const size_type& domainIndex)
01406   {
01407     //assign right hand side of current domainindex block
01408     // rhs is an array of doubles!
01409     // rhs[starti] = b[domainindex]
01410     for(size_type j=0; j<n; ++j, ++i){
01411       assert(i<maxlength_);
01412       rhs_[i]=(*b)[domainIndex][j];
01413     }
01414     
01415     
01416     // loop over all Matrix row entries and calculate defect.
01417     typedef typename matrix_type::ConstColIterator col_iterator;
01418 
01419     // calculate defect for current row index block
01420     for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
01421         block_type tmp;
01422         (*col).mv((*x)[col.index()], tmp);
01423         i-=n;
01424         for(size_type j=0; j<n; ++j, ++i){
01425           assert(i<maxlength_);
01426           rhs_[i]-=tmp[j];
01427         }
01428         
01429     }
01430     
01431   }
01432   template<typename T, typename A, int n, int m>
01433   void OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::relaxResult(field_type relax)
01434   {
01435     for(size_type j=i+n; i<j; ++i){
01436       assert(i<maxlength_);
01437       lhs_[i]*=relax;
01438     }
01439     i-=n;
01440   }
01441   
01442   template<typename T, typename A, int n, int m>
01443   void OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::assignResult(block_type& res)
01444   {
01445     // assign the result of the local solve to the global vector
01446     for(size_type j=0; j<n; ++j, ++i){
01447       assert(i<maxlength_);
01448       res[j]+=lhs_[i];
01449     }
01450   }
01451   
01452   template<typename T, typename A, int n, int m>
01453   void OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::resetIndexForNextDomain()
01454   {
01455     i=0;
01456   }
01457   
01458   template<typename T, typename A, int n, int m>
01459   typename OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::field_type* 
01460   OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::lhs()
01461   {
01462     return lhs_;
01463   }
01464   
01465   template<typename T, typename A, int n, int m>
01466   typename OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::field_type* 
01467   OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >::rhs()
01468   {
01469     return rhs_;
01470   }
01471 
01472 #endif
01473 
01474   template<class M, class X, class Y>
01475   OverlappingAssignerILUBase<M,X,Y>::OverlappingAssignerILUBase(std::size_t maxlength, 
01476                                                                         const M& mat_,
01477                                                                         const Y& b_,
01478                                                                         X& x_)
01479     : mat(&mat_), 
01480       b(&b_), 
01481       x(&x_), i(0)
01482   {
01483     rhs_= new Y(maxlength);
01484     lhs_ = new X(maxlength);
01485   }
01486   
01487   template<class M, class X, class Y>
01488   void OverlappingAssignerILUBase<M,X,Y>::deallocate()
01489   {
01490     delete rhs_;
01491     delete lhs_;
01492   }
01493   
01494   template<class M, class X, class Y>
01495   void OverlappingAssignerILUBase<M,X,Y>::operator()(const size_type& domainIndex)
01496   {
01497     (*rhs_)[i]=(*b)[domainIndex];
01498     
01499      // loop over all Matrix row entries and calculate defect.
01500     typedef typename matrix_type::ConstColIterator col_iterator;
01501 
01502     // calculate defect for current row index block
01503     for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
01504       (*col).mmv((*x)[col.index()], (*rhs_)[i]);
01505     }
01506     // Goto next local index
01507     ++i;
01508   }
01509   
01510   template<class M, class X, class Y>
01511   void OverlappingAssignerILUBase<M,X,Y>::relaxResult(field_type relax)
01512   {
01513     (*lhs_)[i]*=relax;
01514   }
01515   
01516   template<class M, class X, class Y>
01517   void OverlappingAssignerILUBase<M,X,Y>::assignResult(block_type& res)
01518   {
01519     res+=(*lhs_)[i++];
01520   }
01521 
01522   template<class M, class X, class Y>
01523   X& OverlappingAssignerILUBase<M,X,Y>::lhs()
01524   {
01525     return *lhs_;
01526   }
01527   
01528   template<class M, class X, class Y>
01529   Y& OverlappingAssignerILUBase<M,X,Y>::rhs()
01530   {
01531     return *rhs_;
01532   }
01533 
01534   template<class M, class X, class Y>
01535   void OverlappingAssignerILUBase<M,X,Y>::resetIndexForNextDomain()
01536   {
01537     i=0;
01538   }
01539   
01540   template<typename S, typename T, typename A, int n>
01541   AdditiveAdder<S,BlockVector<FieldVector<T,n>,A> >::AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v_, 
01542                                                                    BlockVector<FieldVector<T,n>,A>& x_, 
01543                                                                    OverlappingAssigner<S>& assigner_, 
01544                                                                    const T& relax_)
01545     : v(&v_), x(&x_), assigner(&assigner_), relax(relax_)
01546   {}
01547 
01548   template<typename S, typename T, typename A, int n>
01549   void AdditiveAdder<S,BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
01550   {
01551     // add the result of the local solve to the current update
01552     assigner->assignResult((*v)[domainIndex]);
01553   }
01554 
01555  
01556   template<typename S, typename T, typename A, int n>
01557   void AdditiveAdder<S,BlockVector<FieldVector<T,n>,A> >::axpy()
01558   {
01559     // relax the update and add it to the current guess.
01560     x->axpy(relax,*v);
01561   }
01562 
01563  
01564   template<typename S, typename T, typename A, int n>
01565   MultiplicativeAdder<S,BlockVector<FieldVector<T,n>,A> >
01566   ::MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v_, 
01567                         BlockVector<FieldVector<T,n>,A>& x_, 
01568                         OverlappingAssigner<S>& assigner_, const T& relax_)
01569     : x(&x_), assigner(&assigner_), relax(relax_)
01570   {}
01571 
01572   
01573   template<typename S,typename T, typename A, int n>
01574   void MultiplicativeAdder<S,BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
01575   {
01576     // add the result of the local solve to the current guess
01577     assigner->relaxResult(relax);
01578     assigner->assignResult((*x)[domainIndex]);
01579   }
01580 
01581   
01582   template<typename S,typename T, typename A, int n>
01583   void MultiplicativeAdder<S,BlockVector<FieldVector<T,n>,A> >::axpy()
01584   {
01585     // nothing to do, as the corrections already relaxed and added in operator()
01586   }
01587   
01588   
01590 }
01591 
01592 #endif

Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].