00001 #ifndef DUNE_MATRIX_HH
00002 #define DUNE_MATRIX_HH
00003
00008 #include <vector>
00009 #include <dune/istl/bvector.hh>
00010
00011 namespace Dune {
00012
00022 template<class T, class A=ISTLAllocator>
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 BlockVector<T> row_type;
00038
00040 typedef typename A::size_type size_type;
00041
00043 typedef typename BlockVector<row_type>::iterator RowIterator;
00044
00046 typedef typename row_type::iterator ColIterator;
00047
00049 typedef typename BlockVector<row_type>::const_iterator 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)
00066 {
00067 for (size_type i=0; i<rows; i++)
00068 data_[i].resize(cols);
00069 }
00070
00075 void setSize(size_type rows, size_type cols) {
00076 data_.resize(rows);
00077 for (size_type i=0; i<rows; i++)
00078 data_[i].resize(cols);
00079 cols_ = cols;
00080 }
00081
00083 RowIterator begin()
00084 {
00085 return data_.begin();
00086 }
00087
00089 RowIterator end()
00090 {
00091 return data_.end();
00092 }
00093
00095 RowIterator rbegin()
00096 {
00097 return data_.rbegin();
00098 }
00099
00101 RowIterator rend()
00102 {
00103 return data_.rend();
00104 }
00105
00107 ConstRowIterator begin() const
00108 {
00109 return data_.begin();
00110 }
00111
00113 ConstRowIterator end() const
00114 {
00115 return data_.end();
00116 }
00117
00119 ConstRowIterator rbegin() const
00120 {
00121 return data_.rbegin();
00122 }
00123
00125 ConstRowIterator rend() const
00126 {
00127 return data_.rend();
00128 }
00129
00131 Matrix& operator= (const field_type& t)
00132 {
00133 typedef typename BlockVector<row_type,allocator_type>::size_type
00134 size_type;
00135
00136 for (size_type i=0; i<data_.size(); i++)
00137 data_[i] = t;
00138 return *this;
00139 }
00140
00142 row_type& operator[](size_type row) {
00143 #ifdef DUNE_ISTL_WITH_CHECKING
00144 if (row<0)
00145 DUNE_THROW(ISTLError, "Can't access negative rows!");
00146 if (row>=N())
00147 DUNE_THROW(ISTLError, "Row index out of range!");
00148 #endif
00149 return data_[row];
00150 }
00151
00153 const row_type& operator[](size_type row) const {
00154 #ifdef DUNE_ISTL_WITH_CHECKING
00155 if (row<0)
00156 DUNE_THROW(ISTLError, "Can't access negative rows!");
00157 if (row>=N())
00158 DUNE_THROW(ISTLError, "Row index out of range!");
00159 #endif
00160 return data_[row];
00161 }
00162
00164 size_type N() const {
00165 return data_.size();
00166 }
00167
00169 size_type M() const {
00170 return cols_;
00171 }
00172
00174 size_type rowdim() const {
00175 #ifdef DUNE_ISTL_WITH_CHECKING
00176 if (M()==0)
00177 DUNE_THROW(ISTLError, "Can't compute rowdim() when there are no columns!");
00178 #endif
00179 size_type dim = 0;
00180 for (size_type i=0; i<data_.size(); i++)
00181 dim += data_[i][0].rowdim();
00182
00183 return dim;
00184 }
00185
00187 size_type coldim() const {
00188 #ifdef DUNE_ISTL_WITH_CHECKING
00189 if (N()==0)
00190 DUNE_THROW(ISTLError, "Can't compute coldim() when there are no rows!");
00191 #endif
00192 size_type dim = 0;
00193 for (size_type i=0; i<data_[0].size(); i++)
00194 dim += data_[0][i].coldim();
00195
00196 return dim;
00197 }
00198
00200 size_type rowdim(size_type r) const {
00201 #ifdef DUNE_ISTL_WITH_CHECKING
00202 if (r<0 || r>=N())
00203 DUNE_THROW(ISTLError, "Rowdim for nonexisting row " << r << " requested!");
00204 if (M()==0)
00205 DUNE_THROW(ISTLError, "Can't compute rowdim() when there are no columns!");
00206 #endif
00207 return data_[r][0].rowdim();
00208 }
00209
00211 size_type coldim(size_type c) const {
00212 #ifdef DUNE_ISTL_WITH_CHECKING
00213 if (c<0 || c>=M())
00214 DUNE_THROW(ISTLError, "Coldim for nonexisting column " << c << " requested!");
00215 if (N()==0)
00216 DUNE_THROW(ISTLError, "Can't compute coldim() when there are no rows!");
00217 #endif
00218 return data_[0][c].coldim();
00219 }
00220
00222 Matrix<T>& operator*=(const field_type& scalar) {
00223 for (size_type row=0; row<data_.size(); row++)
00224 for (size_type col=0; col<cols_; col++)
00225 (*this)[row][col] *= scalar;
00226
00227 return (*this);
00228 }
00229
00231 Matrix<T>& operator/=(const field_type& scalar) {
00232 for (size_type row=0; row<data_.size(); row++)
00233 for (size_type col=0; col<cols_; col++)
00234 (*this)[row][col] /= scalar;
00235
00236 return (*this);
00237 }
00238
00244 Matrix& operator+= (const Matrix& b) {
00245 #ifdef DUNE_ISTL_WITH_CHECKING
00246 if(N()!=b.N() || M() != b.M())
00247 DUNE_THROW(RangeError, "Matrix sizes do not match!");
00248 #endif
00249 for (size_type row=0; row<data_.size(); row++)
00250 (*this)[row] += b[row];
00251
00252 return (*this);
00253 }
00254
00260 Matrix& operator-= (const Matrix& b) {
00261 #ifdef DUNE_ISTL_WITH_CHECKING
00262 if(N()!=b.N() || M() != b.M())
00263 DUNE_THROW(RangeError, "Matrix sizes do not match!");
00264 #endif
00265 for (size_type row=0; row<data_.size(); row++)
00266 (*this)[row] -= b[row];
00267
00268 return (*this);
00269 }
00270
00272 Matrix transpose() const {
00273 Matrix out(N(), M());
00274 for (size_type i=0; i<M(); i++)
00275 for (size_type j=0; j<N(); j++)
00276 out[j][i] = (*this)[i][j];
00277
00278 return out;
00279 }
00280
00282 template <class X, class Y>
00283 Y transposedMult(const Y& vec) {
00284 #ifdef DUNE_ISTL_WITH_CHECKING
00285 if (N()!=vec.size())
00286 DUNE_THROW(ISTLError, "Vector size doesn't match the number of matrix rows!");
00287 #endif
00288 Y out(M());
00289 out = 0;
00290
00291 for (size_type i=0; i<out.size(); i++ ) {
00292 for ( size_type j=0; j<vec.size(); j++ )
00293 out[i] += (*this)[j][i]*vec[j];
00294 }
00295
00296 return out;
00297 }
00298
00300 friend Matrix<T> operator*(const Matrix<T>& m1, const Matrix<T>& m2) {
00301 Matrix<T> out(m1.N(), m2.M());
00302 out.clear();
00303
00304 for (size_type i=0; i<out.N(); i++ ) {
00305 for ( size_type j=0; j<out.M(); j++ )
00306 for (size_type k=0; k<m1.M(); k++)
00307 out[i][j] += m1[i][k]*m2[k][j];
00308 }
00309
00310 return out;
00311 }
00312
00314 template <class X, class Y>
00315 friend Y operator*(const Matrix<T>& m, const X& vec) {
00316 #ifdef DUNE_ISTL_WITH_CHECKING
00317 if (M()!=vec.size())
00318 DUNE_THROW(ISTLError, "Vector size doesn't match the number of matrix columns!");
00319 #endif
00320 Y out(m.N());
00321 out = 0;
00322
00323 for (size_type i=0; i<out.size(); i++ ) {
00324 for ( size_type j=0; j<vec.size(); j++ )
00325 out[i] += m[i][j]*vec[j];
00326 }
00327
00328 return out;
00329 }
00330
00332 template <class X, class Y>
00333 void mv(const X& x, Y& y) const
00334 {
00335 #ifdef DUNE_ISTL_WITH_CHECKING
00336 if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00337 if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00338 #endif
00339
00340 for (size_type i=0; i<data_.size(); i++) {
00341 y[i]=0;
00342 for (size_type j=0; j<cols_; j++)
00343 (*this)[i][j].umv(x[j], y[i]);
00344
00345 }
00346
00347 }
00349 template <class X, class Y>
00350 void umv(const X& x, Y& y) const
00351 {
00352 #ifdef DUNE_ISTL_WITH_CHECKING
00353 if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00354 if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00355 #endif
00356
00357 for (size_type i=0; i<data_.size(); i++) {
00358
00359 for (size_type j=0; j<cols_; j++)
00360 (*this)[i][j].umv(x[j], y[i]);
00361
00362 }
00363
00364 }
00365
00367 template<class X, class Y>
00368 void mmv (const X& x, Y& y) const
00369 {
00370 #ifdef DUNE_ISTL_WITH_CHECKING
00371 if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00372 if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00373 #endif
00374 ConstRowIterator endi=end();
00375 for (ConstRowIterator i=begin(); i!=endi; ++i)
00376 {
00377 ConstColIterator endj = (*i).end();
00378 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00379 (*j).mmv(x[j.index()],y[i.index()]);
00380 }
00381 }
00382
00384 template <class X, class Y>
00385 void usmv(const field_type& alpha, const X& x, Y& y) const
00386 {
00387 #ifdef DUNE_ISTL_WITH_CHECKING
00388 if (x.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00389 if (y.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00390 #endif
00391
00392 for (size_type i=0; i<data_.size(); i++) {
00393
00394 for (size_type j=0; j<cols_; j++)
00395 (*this)[i][j].usmv(alpha, x[j], y[i]);
00396
00397 }
00398
00399 }
00400
00402 template<class X, class Y>
00403 void umtv (const X& x, Y& y) const
00404 {
00405 #ifdef DUNE_ISTL_WITH_CHECKING
00406 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00407 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00408 #endif
00409 ConstRowIterator endi=end();
00410 for (ConstRowIterator i=begin(); i!=endi; ++i)
00411 {
00412 ConstColIterator endj = (*i).end();
00413 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00414 (*j).umtv(x[i.index()],y[j.index()]);
00415 }
00416 }
00417
00419 template<class X, class Y>
00420 void mmtv (const X& x, Y& y) const
00421 {
00422 #ifdef DUNE_ISTL_WITH_CHECKING
00423 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00424 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00425 #endif
00426 ConstRowIterator endi=end();
00427 for (ConstRowIterator i=begin(); i!=endi; ++i)
00428 {
00429 ConstColIterator endj = (*i).end();
00430 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00431 (*j).mmtv(x[i.index()],y[j.index()]);
00432 }
00433 }
00434
00436 template<class X, class Y>
00437 void usmtv (const field_type& alpha, const X& x, Y& y) const
00438 {
00439 #ifdef DUNE_ISTL_WITH_CHECKING
00440 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00441 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00442 #endif
00443 ConstRowIterator endi=end();
00444 for (ConstRowIterator i=begin(); i!=endi; ++i)
00445 {
00446 ConstColIterator endj = (*i).end();
00447 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00448 (*j).usmtv(alpha,x[i.index()],y[j.index()]);
00449 }
00450 }
00451
00453 template<class X, class Y>
00454 void umhv (const X& x, Y& y) const
00455 {
00456 #ifdef DUNE_ISTL_WITH_CHECKING
00457 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00458 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00459 #endif
00460 ConstRowIterator endi=end();
00461 for (ConstRowIterator i=begin(); i!=endi; ++i)
00462 {
00463 ConstColIterator endj = (*i).end();
00464 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00465 (*j).umhv(x[i.index()],y[j.index()]);
00466 }
00467 }
00468
00470 template<class X, class Y>
00471 void mmhv (const X& x, Y& y) const
00472 {
00473 #ifdef DUNE_ISTL_WITH_CHECKING
00474 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00475 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00476 #endif
00477 ConstRowIterator endi=end();
00478 for (ConstRowIterator i=begin(); i!=endi; ++i)
00479 {
00480 ConstColIterator endj = (*i).end();
00481 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00482 (*j).mmhv(x[i.index()],y[j.index()]);
00483 }
00484 }
00485
00487 template<class X, class Y>
00488 void usmhv (const field_type& alpha, const X& x, Y& y) const
00489 {
00490 #ifdef DUNE_ISTL_WITH_CHECKING
00491 if (x.N()!=N()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00492 if (y.N()!=M()) DUNE_THROW(ISTLError,"vector/matrix size mismatch!");
00493 #endif
00494 ConstRowIterator endi=end();
00495 for (ConstRowIterator i=begin(); i!=endi; ++i)
00496 {
00497 ConstColIterator endj = (*i).end();
00498 for (ConstColIterator j=(*i).begin(); j!=endj; ++j)
00499 (*j).usmhv(alpha,x[i.index()],y[j.index()]);
00500 }
00501 }
00502
00503
00504
00506 double frobenius_norm () const
00507 {
00508 return std::sqrt(frobenius_norm2());
00509 }
00510
00512 double frobenius_norm2 () const
00513 {
00514 double sum=0;
00515 for (size_type i=0; i<N(); ++i)
00516 for (size_type j=0; j<M(); ++j)
00517 sum += data_[i][j].frobenius_norm2();
00518 return sum;
00519 }
00520
00522 double infinity_norm () const
00523 {
00524 double max=0;
00525 for (size_type i=0; i<N(); ++i) {
00526 double sum=0;
00527 for (size_type j=0; j<M(); j++)
00528 sum += data_[i][j].infinity_norm();
00529 max = std::max(max,sum);
00530 }
00531 return max;
00532 }
00533
00535 double infinity_norm_real () const
00536 {
00537 double max=0;
00538 for (size_type i=0; i<N(); ++i) {
00539 double sum=0;
00540 for (size_type j=0; j<M(); j++)
00541 sum += data_[i][j].infinity_norm_real();
00542 max = std::max(max,sum);
00543 }
00544 return max;
00545 }
00546
00547
00548
00550 bool exists (size_type i, size_type j) const
00551 {
00552 #ifdef DUNE_ISTL_WITH_CHECKING
00553 if (i<0 || i>=N()) DUNE_THROW(ISTLError,"row index out of range");
00554 if (j<0 || i>=M()) DUNE_THROW(ISTLError,"column index out of range");
00555 #endif
00556 return true;
00557 }
00558
00559 protected:
00560
00561 BlockVector<row_type, allocator_type> data_;
00562
00563 size_type cols_;
00564 };
00566 }
00567
00568 #endif