dune-istl  2.1.1
matrix.hh
Go to the documentation of this file.
00001 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 4 -*-
00002 // vi: set et ts=8 sw=4 sts=4:
00003 #ifndef DUNE_MATRIX_HH
00004 #define DUNE_MATRIX_HH
00005 
00010 #include <vector>
00011 #include <memory>
00012 
00013 #include <dune/istl/vbvector.hh>
00014 
00015 namespace Dune {
00016 
00022     template<class T, class A=std::allocator<T> >
00023 class Matrix
00024 {
00025 public:
00026 
00028     typedef typename T::field_type field_type;
00029 
00031     typedef T block_type;
00032 
00034     typedef A allocator_type;
00035 
00037     typedef typename VariableBlockVector<T,A>::window_type row_type;
00038 
00040     typedef typename A::size_type size_type;
00041 
00043     typedef typename VariableBlockVector<T,A>::Iterator RowIterator;
00044 
00046     typedef typename row_type::iterator ColIterator;
00047 
00049     typedef typename VariableBlockVector<T,A>::ConstIterator ConstRowIterator;
00050 
00052     typedef typename row_type::const_iterator ConstColIterator;
00053 
00054     enum {
00056         blocklevel = T::blocklevel+1
00057     };
00058 
00060     Matrix() : data_(0), cols_(0) 
00061     {}
00062 
00065     Matrix(size_type rows, size_type cols) : data_(rows,cols), cols_(cols)
00066     {}
00067 
00072     void setSize(size_type rows, size_type cols) {
00073         data_.resize(rows,cols);
00074         cols_ = cols;
00075     }
00076 
00078     RowIterator begin()
00079     {
00080         return data_.begin();
00081     }
00082 
00084     RowIterator end()
00085     {
00086         return data_.end();
00087     }
00088 
00093     RowIterator rbegin() DUNE_DEPRECATED
00094     {
00095         return beforeEnd();
00096     }
00097 
00100     RowIterator beforeEnd ()
00101     {
00102         return data_.beforeEnd();
00103     }
00104 
00109     RowIterator rend ()  DUNE_DEPRECATED
00110     {
00111         return beforeBegin();
00112     }
00113 
00116     RowIterator beforeBegin ()
00117     {
00118         return data_.beforeBegin();
00119     }
00120 
00122     ConstRowIterator begin() const
00123     {
00124         return data_.begin();
00125     }
00126 
00128     ConstRowIterator end() const
00129     {
00130         return data_.end();
00131     }
00132 
00137     ConstRowIterator rbegin() const DUNE_DEPRECATED
00138     {
00139         return beforeEnd();
00140     }
00141 
00144     ConstRowIterator beforeEnd() const
00145     {
00146         return data_.beforeEnd();
00147     }
00148 
00153     ConstRowIterator rend () const DUNE_DEPRECATED
00154     {
00155         return beforeBegin();
00156     }
00157     
00160     ConstRowIterator beforeBegin () const
00161     {
00162         return data_.beforeBegin();
00163     }
00164 
00166     Matrix& operator= (const field_type& t)
00167     {
00168         data_ = t;
00169       return *this;
00170     }
00171 
00173     row_type& operator[](size_type row) {
00174 #ifdef DUNE_ISTL_WITH_CHECKING
00175         if (row<0)
00176             DUNE_THROW(ISTLError, "Can't access negative rows!");
00177         if (row>=N())
00178             DUNE_THROW(ISTLError, "Row index out of range!");
00179 #endif
00180         return data_[row];
00181     }
00182 
00184     const row_type& operator[](size_type row) const {
00185 #ifdef DUNE_ISTL_WITH_CHECKING
00186         if (row<0)
00187             DUNE_THROW(ISTLError, "Can't access negative rows!");
00188         if (row>=N())
00189             DUNE_THROW(ISTLError, "Row index out of range!");
00190 #endif
00191         return data_[row];
00192     }
00193 
00195     size_type N() const {
00196         return data_.N();
00197     }
00198 
00200     size_type M() const {
00201         return cols_;
00202     }
00203 
00205     size_type rowdim() const {
00206 #ifdef DUNE_ISTL_WITH_CHECKING
00207         if (M()==0)
00208             DUNE_THROW(ISTLError, "Can't compute rowdim() when there are no columns!");
00209 #endif
00210         size_type dim = 0;
00211         for (size_type i=0; i<data_.N(); i++)
00212             dim += data_[i][0].rowdim();
00213 
00214         return dim;
00215     }
00216 
00218     size_type coldim() const {
00219 #ifdef DUNE_ISTL_WITH_CHECKING
00220         if (N()==0)
00221             DUNE_THROW(ISTLError, "Can't compute coldim() when there are no rows!");
00222 #endif
00223         size_type dim = 0;
00224         for (size_type i=0; i<data_[0].size(); i++)
00225             dim += data_[0][i].coldim();
00226 
00227         return dim;
00228     }
00229 
00231     size_type rowdim(size_type r) const {
00232 #ifdef DUNE_ISTL_WITH_CHECKING
00233         if (r<0 || r>=N())
00234             DUNE_THROW(ISTLError, "Rowdim for nonexisting row " << r << " requested!");
00235         if (M()==0)
00236             DUNE_THROW(ISTLError, "Can't compute rowdim() when there are no columns!");
00237 #endif
00238         return data_[r][0].rowdim();
00239     }
00240 
00242     size_type coldim(size_type c) const {
00243 #ifdef DUNE_ISTL_WITH_CHECKING
00244         if (c<0 || c>=M())
00245             DUNE_THROW(ISTLError, "Coldim for nonexisting column " << c << " requested!");
00246         if (N()==0)
00247             DUNE_THROW(ISTLError, "Can't compute coldim() when there are no rows!");
00248 #endif
00249         return data_[0][c].coldim();
00250     }
00251 
00253     Matrix<T>& operator*=(const field_type& scalar) {
00254         data_ *= scalar;
00255         return (*this);
00256     }
00257 
00259     Matrix<T>& operator/=(const field_type& scalar) {
00260         data_ /= scalar;
00261         return (*this);
00262     }
00263 
00269     Matrix& operator+= (const Matrix& b) {
00270 #ifdef DUNE_ISTL_WITH_CHECKING
00271         if(N()!=b.N() || M() != b.M())
00272             DUNE_THROW(RangeError, "Matrix sizes do not match!");
00273 #endif
00274         data_ += b.data_;
00275         return (*this);
00276     }
00277 
00283     Matrix& operator-= (const Matrix& b) {
00284 #ifdef DUNE_ISTL_WITH_CHECKING
00285         if(N()!=b.N() || M() != b.M())
00286             DUNE_THROW(RangeError, "Matrix sizes do not match!");
00287 #endif
00288         data_ -= b.data_;
00289         return (*this);
00290     }
00291 
00293     Matrix transpose() const {
00294         Matrix out(N(), M());
00295         for (size_type i=0; i<M(); i++)
00296             for (size_type j=0; j<N(); j++)
00297                 out[j][i] = (*this)[i][j];
00298 
00299         return out;
00300     }
00301 
00303     template <class X, class Y>
00304     Y transposedMult(const X& vec) {
00305 #ifdef DUNE_ISTL_WITH_CHECKING
00306         if (N()!=vec.size())
00307             DUNE_THROW(ISTLError, "Vector size doesn't match the number of matrix rows!");
00308 #endif
00309         Y out(M());
00310         out = 0;
00311 
00312         for (size_type i=0; i<out.size(); i++ ) {
00313             for ( size_type j=0; j<vec.size(); j++ ) 
00314                 out[i] += (*this)[j][i]*vec[j];
00315         }
00316 
00317         return out;
00318     }
00319 
00321     friend Matrix<T> operator*(const Matrix<T>& m1, const Matrix<T>& m2) {
00322         Matrix<T> out(m1.N(), m2.M());
00323         out = 0;
00324         
00325         for (size_type i=0; i<out.N(); i++ ) {
00326             for ( size_type j=0; j<out.M(); j++ ) 
00327                 for (size_type k=0; k<m1.M(); k++)
00328                     out[i][j] += m1[i][k]*m2[k][j];
00329         }
00330 
00331         return out;
00332     }
00333 
00335     template <class X, class Y>
00336     friend Y operator*(const Matrix<T>& m, const X& vec) {
00337 #ifdef DUNE_ISTL_WITH_CHECKING
00338         if (m.M()!=vec.size())
00339             DUNE_THROW(ISTLError, "Vector size doesn't match the number of matrix columns!");
00340 #endif
00341         Y out(m.N());
00342         out = 0;
00343 
00344         for (size_type i=0; i<out.size(); i++ ) {
00345             for ( size_type j=0; j<vec.size(); j++ ) 
00346                 out[i] += m[i][j]*vec[j];
00347         }
00348 
00349         return out;
00350     }
00351 
00353     template <class X, class Y>
00354     void mv(const X& x, Y& y) const
00355     {
00356 #ifdef DUNE_ISTL_WITH_CHECKING
00357         if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00358         if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00359 #endif
00360 
00361         for (size_type i=0; i<data_.N(); i++) {
00362           y[i]=0;
00363             for (size_type j=0; j<cols_; j++)
00364                 (*this)[i][j].umv(x[j], y[i]);
00365 
00366         }
00367 
00368     }
00369 
00371     template<class X, class Y>
00372     void mtv (const X& x, Y& y) const
00373     {
00374 #ifdef DUNE_ISTL_WITH_CHECKING
00375         if (x.N()!=N()) DUNE_THROW(ISTLError,"index out of range");
00376         if (y.N()!=M()) DUNE_THROW(ISTLError,"index out of range");
00377 #endif
00378         
00379         for(size_type i=0; i<y.N(); ++i)
00380             y[i]=0;
00381         umtv(x,y);
00382     }
00383 
00385     template <class X, class Y>
00386     void umv(const X& x, Y& y) const
00387     {
00388 #ifdef DUNE_ISTL_WITH_CHECKING
00389         if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00390         if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00391 #endif
00392 
00393         for (size_type i=0; i<data_.N(); i++) {
00394 
00395             for (size_type j=0; j<cols_; j++)
00396                 (*this)[i][j].umv(x[j], y[i]);
00397 
00398         }
00399 
00400     }
00401 
00403     template<class X, class Y>
00404     void mmv (const X& x, Y& y) const
00405     {
00406 #ifdef DUNE_ISTL_WITH_CHECKING
00407         if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00408         if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00409 #endif
00410         ConstRowIterator endi=end();
00411         for (ConstRowIterator i=begin(); i!=endi; ++i)
00412             {
00413                 ConstColIterator endj = (*i).end();
00414                 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00415                     (*j).mmv(x[j.index()],y[i.index()]);
00416             }
00417     }
00418 
00420     template <class X, class Y>
00421     void usmv(const field_type& alpha, const X& x, Y& y) const
00422     {
00423 #ifdef DUNE_ISTL_WITH_CHECKING
00424         if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00425         if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00426 #endif
00427 
00428         for (size_type i=0; i<data_.N(); i++) {
00429 
00430             for (size_type j=0; j<cols_; j++)
00431                 (*this)[i][j].usmv(alpha, x[j], y[i]);
00432             
00433         }
00434         
00435     }
00436     
00438     template<class X, class Y>
00439     void umtv (const X& x, Y& y) const
00440     {
00441 #ifdef DUNE_ISTL_WITH_CHECKING
00442         if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00443         if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00444 #endif
00445         ConstRowIterator endi=end();
00446         for (ConstRowIterator i=begin(); i!=endi; ++i)
00447             {
00448                 ConstColIterator endj = (*i).end();
00449                 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00450                     (*j).umtv(x[i.index()],y[j.index()]);
00451             }
00452     }
00453     
00455     template<class X, class Y>
00456     void mmtv (const X& x, Y& y) const
00457     {
00458 #ifdef DUNE_ISTL_WITH_CHECKING
00459         if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00460         if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00461 #endif
00462         ConstRowIterator endi=end();
00463         for (ConstRowIterator i=begin(); i!=endi; ++i)
00464             {
00465                 ConstColIterator endj = (*i).end();
00466                 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00467                     (*j).mmtv(x[i.index()],y[j.index()]);
00468             }
00469     }
00470     
00472     template<class X, class Y>
00473     void usmtv (const field_type& alpha, const X& x, Y& y) const
00474     {
00475 #ifdef DUNE_ISTL_WITH_CHECKING
00476         if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00477         if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00478 #endif
00479         ConstRowIterator endi=end();
00480         for (ConstRowIterator i=begin(); i!=endi; ++i)
00481             {
00482                 ConstColIterator endj = (*i).end();
00483                 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00484                     (*j).usmtv(alpha,x[i.index()],y[j.index()]);
00485             }
00486     }
00487     
00489     template<class X, class Y>
00490     void umhv (const X& x, Y& y) const
00491     {
00492 #ifdef DUNE_ISTL_WITH_CHECKING
00493         if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00494         if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00495 #endif
00496         ConstRowIterator endi=end();
00497         for (ConstRowIterator i=begin(); i!=endi; ++i)
00498             {
00499                 ConstColIterator endj = (*i).end();
00500                 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00501                     (*j).umhv(x[i.index()],y[j.index()]);
00502             }
00503     }
00504     
00506     template<class X, class Y>
00507     void mmhv (const X& x, Y& y) const
00508     {
00509 #ifdef DUNE_ISTL_WITH_CHECKING
00510         if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00511         if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00512 #endif
00513         ConstRowIterator endi=end();
00514         for (ConstRowIterator i=begin(); i!=endi; ++i)
00515             {
00516                 ConstColIterator endj = (*i).end();
00517                 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00518                     (*j).mmhv(x[i.index()],y[j.index()]);
00519             }
00520     }
00521     
00523     template<class X, class Y>
00524     void usmhv (const field_type& alpha, const X& x, Y& y) const
00525     {
00526 #ifdef DUNE_ISTL_WITH_CHECKING
00527         if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00528         if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00529 #endif
00530         ConstRowIterator endi=end();
00531         for (ConstRowIterator i=begin(); i!=endi; ++i)
00532             {
00533                 ConstColIterator endj = (*i).end();
00534                 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00535                     (*j).usmhv(alpha,x[i.index()],y[j.index()]);
00536             }
00537     }
00538 
00539     //===== norms
00540     
00542     double frobenius_norm () const
00543     {
00544         return std::sqrt(frobenius_norm2());
00545     }
00546     
00548     double frobenius_norm2 () const
00549     {
00550         double sum=0;
00551         for (size_type i=0; i<N(); ++i) 
00552             for (size_type j=0; j<M(); ++j)
00553                 sum += data_[i][j].frobenius_norm2();
00554         return sum;
00555     }
00556     
00558     double infinity_norm () const
00559     {
00560         double max=0;
00561         for (size_type i=0; i<N(); ++i) {
00562             double sum=0;
00563             for (size_type j=0; j<M(); j++)
00564                 sum += data_[i][j].infinity_norm();
00565             max = std::max(max,sum);
00566         }
00567         return max;
00568     }
00569     
00571     double infinity_norm_real () const
00572     {
00573         double max=0;
00574         for (size_type i=0; i<N(); ++i) {
00575             double sum=0;
00576             for (size_type j=0; j<M(); j++)
00577                 sum += data_[i][j].infinity_norm_real();
00578             max = std::max(max,sum);
00579         }
00580         return max;
00581     }
00582     
00583     //===== query
00584     
00586     bool exists (size_type i, size_type j) const
00587     {
00588 #ifdef DUNE_ISTL_WITH_CHECKING
00589         if (i<0 || i>=N()) DUNE_THROW(ISTLError,"row index out of range");
00590         if (j<0 || i>=M()) DUNE_THROW(ISTLError,"column index out of range");
00591 #endif
00592         return true;
00593     }
00594 
00595 protected:
00596 
00602     VariableBlockVector<T,A> data_;
00603 
00609     size_type cols_;
00610 };
00612 } // end namespace Dune
00613 
00614 #endif