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     
00044     OverlappingSchwarzInitializer(InitializerList& il,
00045                                   const IndexSet& indices);
00046     
00047     
00048     void addRowNnz(const Iter& row);
00049     
00050     void allocate();
00051         
00052     void countEntries(const Iter& row, const CIter& col) const;
00053 
00054     void calcColstart() const;
00055     
00056     void copyValue(const Iter& row, const CIter& col) const;
00057     
00058     void createMatrix() const;
00059   private:
00060     class IndexMap
00061     {
00062     public:
00063       typedef std::size_t size_type;
00064       typedef std::map<size_type,std::size_t> Map;      
00065       typedef typename Map::iterator iterator;
00066       typedef typename Map::const_iterator const_iterator;
00067       
00068       IndexMap();
00069       
00070       void insert(size_type grow);
00071       
00072       const_iterator find(size_type grow) const;
00073       
00074       iterator find(size_type grow);
00075       
00076       iterator end();
00077       
00078       const_iterator end() const;
00079       
00080     private:
00081       std::map<size_type,std::size_t> map_;
00082       std::size_t row;
00083     };
00084     
00085     
00086     typedef typename InitializerList::iterator InitIterator;
00087     typedef typename IndexSet::const_iterator IndexIteratur;
00088     InitializerList* initializers;
00089     const IndexSet *indices;
00090     std::vector<IndexMap> indexMaps;
00091   };
00092 
00096   struct AdditiveSchwarzMode
00097   {};
00098 
00102   struct MultiplicativeSchwarzMode
00103   {};  
00104 
00105   struct SymmetricMultiplicativeSchwarzMode
00106   {};  
00107 
00108   namespace
00109   {
00110     template<typename T>
00111     struct AdditiveAdder
00112     {
00113     };
00114     
00115     template<typename T, typename A, int n>
00116     struct AdditiveAdder<BlockVector<FieldVector<T,n>,A> >
00117     {
00118       AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x, const T* t, const T& relax_);
00119       void operator()(const std::size_t& domain);
00120       void axpy();
00121       
00122     private:
00123       const T* t;
00124       BlockVector<FieldVector<T,n>,A>* v;
00125       BlockVector<FieldVector<T,n>,A>* x;
00126       T relax;
00127       int i;
00128     };
00129 
00130     template<typename T>
00131     struct MultiplicativeAdder
00132     {
00133     };
00134     
00135     template<typename T, typename A, int n>
00136     struct MultiplicativeAdder<BlockVector<FieldVector<T,n>,A> >
00137     {
00138       MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x, const T* t, const T& relax_);
00139       void operator()(const std::size_t& domain);
00140       void axpy();
00141       
00142     private:
00143       const T* t;
00144       BlockVector<FieldVector<T,n>,A>* v;
00145       BlockVector<FieldVector<T,n>,A>* x;
00146       T relax;
00147       int i;
00148     };
00149 
00150     
00151     template<typename T, class X>
00152     struct AdderSelector
00153     {};
00154     
00155     template<class X>
00156     struct AdderSelector<AdditiveSchwarzMode,X>
00157     {
00158       typedef AdditiveAdder<X> Adder;
00159     };
00160     
00161     template<class X>
00162     struct AdderSelector<MultiplicativeSchwarzMode,X>
00163     {
00164       typedef MultiplicativeAdder<X> Adder;
00165     };    
00166 
00167     template<class X>
00168     struct AdderSelector<SymmetricMultiplicativeSchwarzMode,X>
00169     {
00170       typedef MultiplicativeAdder<X> Adder;
00171     };    
00172 
00173     template<typename T1, typename T2, bool forward>
00174     struct ApplyHelper
00175     {
00176       typedef T1 solver_vector;
00177       typedef typename solver_vector::iterator solver_iterator;
00178       typedef T2 subdomain_vector;
00179       typedef typename subdomain_vector::const_iterator domain_iterator;
00180     
00181       static solver_iterator begin(solver_vector& sv)
00182       {
00183         return sv.begin();
00184       }
00185       
00186       static solver_iterator end(solver_vector& sv)
00187       {
00188         return sv.end();
00189       }
00190       static domain_iterator begin(const subdomain_vector& sv)
00191       {
00192         return sv.begin();
00193       }
00194       
00195       static domain_iterator end(const subdomain_vector& sv)
00196       {
00197         return sv.end();
00198       }
00199     };
00200 
00201     template<typename T1, typename T2>
00202     struct ApplyHelper<T1,T2,false>
00203     {
00204       typedef T1 solver_vector;
00205       typedef typename solver_vector::reverse_iterator solver_iterator;
00206       typedef T2 subdomain_vector;
00207       typedef typename subdomain_vector::const_reverse_iterator domain_iterator;
00208     
00209       static solver_iterator begin(solver_vector& sv)
00210       {
00211         return sv.rbegin();
00212       }
00213       
00214       static solver_iterator end(solver_vector& sv)
00215       {
00216         return sv.rend();
00217       }
00218       static domain_iterator begin(const subdomain_vector& sv)
00219       {
00220         return sv.rbegin();
00221       }
00222       
00223       static domain_iterator end(const subdomain_vector& sv)
00224       {
00225         return sv.rend();
00226       }
00227     };
00228     
00229     template<class T>
00230     struct Applier
00231     {
00232       typedef T smoother;
00233       typedef typename smoother::range_type range_type;
00234       
00235       static void apply(smoother& sm, range_type& v, const range_type& b)
00236       {
00237         sm.template apply<true>(v, b);
00238       }
00239     };
00240     
00241     template<class M, class X, class TA>
00242     struct Applier<SeqOverlappingSchwarz<M,X,SymmetricMultiplicativeSchwarzMode, TA> >
00243     {
00244       typedef SeqOverlappingSchwarz<M,X,SymmetricMultiplicativeSchwarzMode, TA> smoother;
00245       typedef typename smoother::range_type range_type;
00246       
00247       static void apply(smoother& sm, range_type& v, const range_type& b)
00248       {
00249         sm.template apply<true>(v, b);
00250         sm.template apply<false>(v, b);
00251       }
00252     };
00253   }
00254   
00258   template<class M, class X, class TM=AdditiveSchwarzMode, class TA=std::allocator<X> >
00259   class SeqOverlappingSchwarz
00260     : public Preconditioner<X,X>
00261   {
00262   public:
00266     typedef M matrix_type;
00267 
00271     typedef X domain_type;
00272 
00276     typedef X range_type;
00277 
00284     typedef TM Mode;
00285 
00289     typedef typename X::field_type field_type;
00290 
00292     typedef typename matrix_type::size_type size_type;
00293 
00295     typedef TA allocator;
00296         
00298     typedef std::set<size_type, std::less<size_type>, typename TA::template rebind< std::less<size_type> >::other> subdomain_type;
00299 
00301     typedef std::vector<subdomain_type, typename TA::template rebind<subdomain_type>::other> subdomain_vector;
00302 
00304     typedef SLList<size_type, typename TA::template rebind<size_type>::other> subdomain_list;
00305       
00307     typedef std::vector<subdomain_list, typename TA::template rebind<subdomain_list>::other > rowtodomain_vector;
00308 
00310     typedef SuperLU<matrix_type> slu;
00311     
00313     typedef std::vector<slu, typename TA::template rebind<slu>::other> slu_vector;
00314 
00315     enum{
00317       category = SolverCategory::sequential
00318         };
00319     
00327     SeqOverlappingSchwarz(const matrix_type& mat, const subdomain_vector& subDomains,
00328                           field_type relaxationFactor=1);
00329 
00335     SeqOverlappingSchwarz(const matrix_type& mat, const rowtodomain_vector& rowToDomain,
00336                           field_type relaxationFactor=1);
00337                           
00343     virtual void pre (X& x, X& b) {}
00344 
00350     virtual void apply (X& v, const X& d);
00351 
00352     template<bool forward>
00353     void apply(X& v, const X& d);
00354     
00360     virtual void post (X& x) {}
00361   
00362   private:
00363     const M& mat;
00364     slu_vector solvers;
00365     subdomain_vector subDomains;
00366     field_type relax;
00367     
00368     int maxlength;
00369 
00370     void initialize(const rowtodomain_vector& rowToDomain, const matrix_type& mat);
00371 
00372     template<typename T>
00373     struct Assigner
00374     {
00375     };
00376     
00377     template<typename T, typename A, int n>
00378     struct Assigner<BlockVector<FieldVector<T,n>,A> >
00379     {
00380       Assigner(const M& mat, T* x, const BlockVector<FieldVector<T,n>,A>& b, const BlockVector<FieldVector<T,n>,A>& x);
00381       void operator()(const std::size_t& domain);
00382     private:
00383       const M* mat;
00384       T* rhs;
00385       const BlockVector<FieldVector<T,n>,A>* b;
00386       const BlockVector<FieldVector<T,n>,A>* x;
00387       int i;
00388     };
00389 
00390   };
00391   
00392   
00393   template<class I, class S>
00394   OverlappingSchwarzInitializer<I,S>::OverlappingSchwarzInitializer(InitializerList& il,
00395                                                                     const IndexSet& idx)
00396     : initializers(&il), indices(&idx), indexMaps(il.size())
00397   {
00398   }
00399   
00400   
00401   template<class I, class S>
00402   void OverlappingSchwarzInitializer<I,S>::addRowNnz(const Iter& row)
00403   {
00404     typedef typename IndexSet::value_type::const_iterator iterator;
00405     for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain){
00406       (*initializers)[*domain].addRowNnz(row);
00407       indexMaps[*domain].insert(row.index());
00408     }
00409   }
00410   
00411   template<class I, class S>
00412   void OverlappingSchwarzInitializer<I,S>::allocate()
00413   {
00414     std::for_each(initializers->begin(), initializers->end(),
00415                   std::mem_fun_ref(&AtomInitializer::allocateMatrixStorage));
00416     std::for_each(initializers->begin(), initializers->end(), 
00417                   std::mem_fun_ref(&AtomInitializer::allocateMarker));
00418   }
00419   
00420   template<class I, class S>
00421   void OverlappingSchwarzInitializer<I,S>::countEntries(const Iter& row, const CIter& col) const
00422   {
00423     typedef typename IndexSet::value_type::const_iterator iterator;
00424     for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain){
00425       typename std::map<std::size_t,std::size_t>::const_iterator v = indexMaps[*domain].find(col.index());
00426       if(v!= indexMaps[*domain].end()){
00427         (*initializers)[*domain].countEntries(indexMaps[*domain].find(col.index())->second);
00428       }
00429     }
00430   }
00431 
00432   template<class I, class S>
00433   void OverlappingSchwarzInitializer<I,S>::calcColstart() const
00434   {
00435     std::for_each(initializers->begin(), initializers->end(),
00436                   std::mem_fun_ref(&AtomInitializer::calcColstart));
00437   }
00438 
00439   template<class I, class S>
00440   void OverlappingSchwarzInitializer<I,S>::copyValue(const Iter& row, const CIter& col) const
00441   {
00442     typedef typename IndexSet::value_type::const_iterator iterator;
00443     for(iterator domain=(*indices)[row.index()].begin(); domain!= (*indices)[row.index()].end(); ++domain){
00444       typename std::map<std::size_t,std::size_t>::const_iterator v = indexMaps[*domain].find(col.index());
00445       if(v!= indexMaps[*domain].end()){
00446         (*initializers)[*domain].copyValue(col, indexMaps[*domain].find(row.index())->second, 
00447                                            indexMaps[*domain].find(col.index())->second);
00448       }
00449     }
00450   }
00451     
00452   template<class I, class S>
00453   void OverlappingSchwarzInitializer<I,S>::createMatrix() const
00454   {
00455     std::for_each(initializers->begin(), initializers->end(),
00456                   std::mem_fun_ref(&AtomInitializer::createMatrix));
00457   }
00458 
00459   template<class I, class S>
00460   OverlappingSchwarzInitializer<I,S>::IndexMap::IndexMap()
00461     : row(0)
00462   {}
00463 
00464   template<class I, class S>
00465   void OverlappingSchwarzInitializer<I,S>::IndexMap::insert(size_type grow)
00466   {
00467     assert(map_.find(grow)==map_.end());
00468     map_.insert(std::make_pair(grow, row++));
00469   }
00470 
00471   template<class I, class S>
00472   typename OverlappingSchwarzInitializer<I,S>::IndexMap::const_iterator 
00473   OverlappingSchwarzInitializer<I,S>::IndexMap::find(size_type grow) const
00474   {
00475     return map_.find(grow);
00476   }
00477   
00478   template<class I, class S>
00479   typename OverlappingSchwarzInitializer<I,S>::IndexMap::iterator 
00480   OverlappingSchwarzInitializer<I,S>::IndexMap::find(size_type grow)
00481   {
00482     return map_.find(grow);
00483   }
00484 
00485   template<class I, class S>
00486   typename OverlappingSchwarzInitializer<I,S>::IndexMap::const_iterator 
00487   OverlappingSchwarzInitializer<I,S>::IndexMap::end() const
00488   {
00489     return map_.end();
00490   }
00491   
00492   template<class I, class S>
00493   typename OverlappingSchwarzInitializer<I,S>::IndexMap::iterator 
00494   OverlappingSchwarzInitializer<I,S>::IndexMap::end()
00495   {
00496     return map_.end();
00497   }
00498 
00499   template<class M, class X, class TM, class TA>
00500   SeqOverlappingSchwarz<M,X,TM,TA>::SeqOverlappingSchwarz(const matrix_type& mat_, const rowtodomain_vector& rowToDomain,
00501                         field_type relaxationFactor)
00502     : mat(mat_), relax(relaxationFactor)
00503   {
00504     typedef typename rowtodomain_vector::const_iterator RowDomainIterator;
00505     typedef typename subdomain_list::const_iterator DomainIterator;
00506 #ifdef DUNE_ISTL_WITH_CHECKING
00507     assert(rowToDomain.size()==mat.N());
00508     assert(rowToDomain.size()==mat.M());
00509     
00510     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
00511       assert(iter->size()>0);
00512     
00513 #endif
00514     // calculate the number of domains
00515     int domains=0;
00516     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
00517       for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
00518         domains=std::max(domains, *d);
00519     ++domains;
00520     
00521     solvers.resize(domains);
00522     subDomains.resize(domains);
00523     
00524     // initialize subdomains to row mapping from row to subdomain mapping
00525     int row=0;
00526     for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter, ++row)
00527       for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
00528         subDomains[*d].insert(row);
00529 
00530 #ifdef DUNE_ISTL_WITH_CHECKING
00531     int i=0;
00532     typedef typename subdomain_vector::const_iterator iterator;
00533     for(iterator iter=subDomains.begin(); iter != subDomains.end(); ++iter){
00534       typedef typename subdomain_type::const_iterator entry_iterator;
00535       std::cout<<"domain "<<i++<<":";
00536       for(entry_iterator entry = iter->begin(); entry != iter->end(); ++entry){
00537         std::cout<<" "<<*entry;
00538       }
00539       std::cout<<std::endl;
00540     }
00541 #endif
00542     initialize(rowToDomain, mat);
00543   }
00544   
00545   template<class M, class X, class TM, class TA>
00546   SeqOverlappingSchwarz<M,X,TM,TA>::SeqOverlappingSchwarz(const matrix_type& mat_,
00547                                                        const subdomain_vector& sd,
00548                                                        field_type relaxationFactor)
00549     :  mat(mat_), solvers(sd.size()), subDomains(sd), relax(relaxationFactor)
00550   {
00551     typedef typename subdomain_vector::const_iterator DomainIterator;
00552 
00553 #ifdef DUNE_ISTL_WITH_CHECKING
00554     int i=0;
00555     
00556     for(DomainIterator d=sd.begin(); d != sd.end(); ++d,++i){
00557       //std::cout<<i<<": "<<d->size()<<std::endl;
00558       assert(d->size()>0);
00559       typedef typename DomainIterator::value_type::const_iterator entry_iterator;
00560       std::cout<<"domain "<<i<<":";
00561       for(entry_iterator entry = d->begin(); entry != d->end(); ++entry){
00562         std::cout<<" "<<*entry;
00563       }
00564       std::cout<<std::endl;
00565     }
00566     
00567 #endif
00568     
00569     // Create a row to subdomain mapping
00570     rowtodomain_vector rowToDomain(mat.N());
00571 
00572     int domainId=0;
00573     
00574     for(DomainIterator domain=sd.begin(); domain != sd.end(); ++domain, ++domainId){
00575       typedef typename subdomain_type::const_iterator iterator;
00576       for(iterator row=domain->begin(); row != domain->end(); ++row)
00577         rowToDomain[*row].push_back(domainId);
00578     }
00579 
00580     initialize(rowToDomain, mat);
00581   }
00582 
00589   template<class M>
00590   struct SeqOverlappingSchwarzDomainSize {};
00591   
00592   template<typename T, typename A, int n, int m>
00593   struct SeqOverlappingSchwarzDomainSize<BCRSMatrix<FieldMatrix<T,n,m>,A > >
00594   {
00595     template<class Domain>
00596     static int size(const Domain & d)
00597     {
00598       assert(n==m);
00599       return m*d.size();
00600     }
00601   };
00602   
00603   template<class M, class X, class TM, class TA>
00604   void SeqOverlappingSchwarz<M,X,TM,TA>::initialize(const rowtodomain_vector& rowToDomain, const matrix_type& mat)
00605   {
00606     typedef typename std::vector<SuperMatrixInitializer<matrix_type> >::iterator InitializerIterator;
00607     typedef typename subdomain_vector::const_iterator DomainIterator;
00608     typedef typename slu_vector::iterator SolverIterator;
00609     // initialize the initializers
00610     DomainIterator domain=subDomains.begin();
00611     
00612     // Create the initializers list.
00613     std::vector<SuperMatrixInitializer<matrix_type> > initializers(subDomains.size());    
00614 
00615     SolverIterator solver=solvers.begin();
00616     for(InitializerIterator initializer=initializers.begin(); initializer!=initializers.end();
00617         ++initializer, ++solver, ++domain){
00618       solver->mat.N_=SeqOverlappingSchwarzDomainSize<M>::size(*domain);
00619       solver->mat.M_=SeqOverlappingSchwarzDomainSize<M>::size(*domain);
00620       //solver->setVerbosity(true);
00621       *initializer=SuperMatrixInitializer<matrix_type>(solver->mat);
00622     }
00623 
00624         // Set up the supermatrices according to the subdomains
00625     typedef OverlappingSchwarzInitializer<std::vector<SuperMatrixInitializer<matrix_type> >,
00626       rowtodomain_vector > Initializer;
00627 
00628     Initializer initializer(initializers, rowToDomain);
00629     copyToSuperMatrix(initializer, mat);
00630     if(solvers.size()==1)
00631       assert(solvers[0].mat==mat);
00632     
00633     /*    for(SolverIterator solver=solvers.begin(); solver!=solvers.end(); ++solver)
00634           dPrint_CompCol_Matrix("superlu", &static_cast<SuperMatrix&>(solver->mat)); */
00635     
00636     // Calculate the LU decompositions
00637     std::for_each(solvers.begin(), solvers.end(), std::mem_fun_ref(&slu::decompose));
00638     maxlength=0;
00639     for(SolverIterator solver=solvers.begin(); solver!=solvers.end(); ++solver){
00640       assert(solver->mat.N()==solver->mat.M());
00641       maxlength=std::max(maxlength, solver->mat.N());
00642       //writeCompColMatrixToMatlab(static_cast<SuperLUMatrix<M>&>(solver->mat), std::cout);
00643     }
00644 
00645   }
00646   
00647   template<class M, class X, class TM, class TA>
00648   void SeqOverlappingSchwarz<M,X,TM,TA>::apply(X& x, const X& b)
00649   {
00650     Applier<SeqOverlappingSchwarz>::apply(*this, x, b);
00651   }
00652 
00653   template<class M, class X, class TM, class TA>
00654   template<bool forward>
00655   void SeqOverlappingSchwarz<M,X,TM,TA>::apply(X& x, const X& b)
00656   {
00657     typedef typename X::block_type block;
00658     typedef slu_vector solver_vector;
00659     typedef typename ApplyHelper<solver_vector,subdomain_vector,forward>::solver_iterator iterator;
00660     typedef typename ApplyHelper<solver_vector,subdomain_vector,forward>::domain_iterator
00661       domain_iterator;
00662     
00663     field_type *lhs=new field_type[maxlength];
00664     field_type *rhs=new field_type[maxlength];
00665     for(int i=0; i<maxlength; ++i)
00666       lhs[i]=0;
00667 
00668     
00669     domain_iterator domain=ApplyHelper<solver_vector,subdomain_vector,forward>::begin(subDomains);
00670     
00671     X v(x); // temporary for the update
00672     v=0;
00673     
00674     typedef typename AdderSelector<TM,X>::Adder Adder;
00675     for(iterator solver=ApplyHelper<solver_vector,subdomain_vector,forward>::begin(solvers); 
00676         solver != ApplyHelper<solver_vector,subdomain_vector,forward>::end(solvers); ++solver, ++domain){
00677       //Copy rhs to C-array for SuperLU
00678       std::for_each(domain->begin(), domain->end(), Assigner<X>(mat, rhs, b, x));
00679       
00680       solver->apply(lhs, rhs);
00681       //Add relaxed correction to from SuperLU to v
00682       std::for_each(domain->begin(), domain->end(), Adder(v, x, lhs, relax));
00683     }
00684     
00685     Adder adder(v, x, lhs, relax);
00686     adder.axpy();
00687     delete[] lhs;
00688     delete[] rhs;
00689     
00690   }   
00691     
00692   template<class M, class X, class TM, class TA>
00693   template<typename T, typename A, int n>
00694   SeqOverlappingSchwarz<M,X,TM,TA>::Assigner<BlockVector<FieldVector<T,n>,A> >::Assigner(const M& mat_, T* rhs_, const BlockVector<FieldVector<T,n>,A>& b_,
00695                                                                                       const BlockVector<FieldVector<T,n>,A>& x_)
00696     : mat(&mat_), rhs(rhs_), b(&b_), x(&x_), i(0)
00697   {}
00698 
00699   template<class M, class X, class TM, class TA>
00700   template<typename T, typename A, int n>
00701   void SeqOverlappingSchwarz<M,X,TM,TA>::Assigner<BlockVector<FieldVector<T,n>,A> >::operator()(const std::size_t& domainIndex)
00702   {
00703     int starti;
00704     starti = i;
00705     
00706     for(int j=0; j<n; ++j, ++i)
00707       rhs[i]=(*b)[domainIndex][j];
00708     
00709     // loop over all Matrix row entries and calculate defect.
00710     typedef typename M::ConstColIterator col_iterator;
00711     typedef typename subdomain_type::const_iterator domain_iterator;
00712 
00713     for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
00714         typename X::block_type tmp;
00715         (*col).mv((*x)[col.index()], tmp);
00716         i-=n;
00717         for(int j=0; j<n; ++j, ++i)
00718           rhs[i]-=tmp[j];
00719     }
00720     assert(starti+n==i);
00721   }
00722   namespace
00723   {
00724     
00725   template<typename T, typename A, int n>
00726   AdditiveAdder<BlockVector<FieldVector<T,n>,A> >::AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v_, 
00727                  BlockVector<FieldVector<T,n>,A>& x_, const T* t_, const T& relax_)
00728     : t(t_), v(&v_), x(&x_), relax(relax_), i(0)
00729   {}
00730 
00731   template<typename T, typename A, int n>
00732   void AdditiveAdder<BlockVector<FieldVector<T,n>,A> >::operator()(const std::size_t& domainIndex)
00733   {
00734     for(int j=0; j<n; ++j, ++i)
00735       (*v)[domainIndex][j]+=t[i];
00736   }
00737 
00738  
00739   template<typename T, typename A, int n>
00740   void AdditiveAdder<BlockVector<FieldVector<T,n>,A> >::axpy()
00741   {
00742     x->axpy(relax,*v);
00743   }
00744 
00745  
00746   template<typename T, typename A, int n>
00747   MultiplicativeAdder<BlockVector<FieldVector<T,n>,A> >
00748   ::MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v_, 
00749                         BlockVector<FieldVector<T,n>,A>& x_, const T* t_, const T& relax_)
00750     : t(t_), v(&v_), x(&x_), relax(relax_), i(0)
00751   {}
00752 
00753   
00754   template<typename T, typename A, int n>
00755   void MultiplicativeAdder<BlockVector<FieldVector<T,n>,A> >::operator()(const std::size_t& domainIndex)
00756   {
00757     for(int j=0; j<n; ++j, ++i)
00758       (*x)[domainIndex][j]+=relax*t[i];
00759   }
00760 
00761   
00762   template<typename T, typename A, int n>
00763   void MultiplicativeAdder<BlockVector<FieldVector<T,n>,A> >::axpy()
00764   {
00765   }
00766   }
00767   
00769 #endif  
00770 }
00771 
00772 #endif

Generated on 9 Apr 2008 with Doxygen (ver 1.5.2) [logfile].