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