dune-istl  2.1.1
scaledidmatrix.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_SCALED_IDENTITY_MATRIX_HH
00004 #define DUNE_SCALED_IDENTITY_MATRIX_HH
00005 
00012 #include<cmath>
00013 #include<cstddef>
00014 #include<complex>
00015 #include<iostream>
00016 #include <dune/common/exceptions.hh>
00017 #include <dune/common/fmatrix.hh>
00018 #include <dune/istl/diagonalmatrix.hh>
00019 
00020 namespace Dune {
00021 
00025 template<class K, int n>
00026 class ScaledIdentityMatrix
00027 {
00028     typedef DiagonalMatrixWrapper< ScaledIdentityMatrix<K,n> > WrapperType;
00029 
00030 public:
00031     //===== type definitions and constants
00032 
00034     typedef K field_type;
00035 
00037     typedef K block_type;
00038 
00040         typedef std::size_t size_type;
00041 
00043     enum {
00045             blocklevel = 1
00046     };
00047 
00049     typedef DiagonalRowVector<K,n> row_type;
00050     typedef row_type reference;
00051     typedef DiagonalRowVectorConst<K,n> const_row_type;
00052     typedef const_row_type const_reference;
00053 
00055     enum {
00057         rows = n,
00059         cols = n
00060     };
00061 
00062     //===== constructors
00065     ScaledIdentityMatrix () {}
00066 
00069     ScaledIdentityMatrix (const K& k)
00070         : p_(k)
00071     {}
00072 
00073     //===== assignment from scalar
00074     ScaledIdentityMatrix& operator= (const K& k)
00075     {
00076         p_ = k;
00077         return *this;
00078     }
00079 
00080     // check if matrix is identical to other matrix (not only identical values)
00081     bool identical(const ScaledIdentityMatrix<K,n>& other) const
00082     {
00083         return (this==&other);
00084     }
00085 
00086     //===== iterator interface to rows of the matrix
00088     typedef ContainerWrapperIterator<const WrapperType, reference, reference> Iterator;
00090     typedef Iterator iterator;
00092     typedef Iterator RowIterator;
00094     typedef typename row_type::Iterator ColIterator;
00095 
00097     Iterator begin ()
00098     {
00099         return Iterator(WrapperType(this),0);
00100     }
00101 
00103     Iterator end ()
00104     {
00105         return Iterator(WrapperType(this),n);
00106     }
00107 
00112     Iterator rbegin() DUNE_DEPRECATED
00113     {
00114         beforeEnd();
00115     }
00116 
00119     Iterator beforeEnd ()
00120     {
00121         return Iterator(WrapperType(this),n-1);
00122     }
00123 
00128     Iterator rend ()  DUNE_DEPRECATED
00129     {
00130         return beforeBegin();
00131     }
00132 
00135     Iterator beforeBegin ()
00136     {
00137         return Iterator(WrapperType(this),-1);
00138     }
00139 
00140 
00142     typedef ContainerWrapperIterator<const WrapperType, const_reference, const_reference> ConstIterator;
00144     typedef ConstIterator const_iterator;
00146     typedef ConstIterator ConstRowIterator;
00148     typedef typename const_row_type::ConstIterator ConstColIterator;
00149 
00151     ConstIterator begin () const
00152     {
00153         return ConstIterator(WrapperType(this),0);
00154     }
00155 
00157     ConstIterator end () const
00158     {
00159         return ConstIterator(WrapperType(this),n);
00160     }
00161 
00166     ConstIterator rbegin() const DUNE_DEPRECATED
00167     {
00168         beforeEnd();
00169     }
00170 
00173     ConstIterator beforeEnd() const
00174     {
00175         return ConstIterator(WrapperType(this),n-1);
00176     }
00177 
00182     ConstIterator rend () const DUNE_DEPRECATED
00183     {
00184         return beforeBegin();
00185     }
00186     
00189     ConstIterator beforeBegin () const
00190     {
00191         return ConstIterator(WrapperType(this),-1);
00192     }
00193 
00194     //===== vector space arithmetic
00195 
00197     ScaledIdentityMatrix& operator+= (const ScaledIdentityMatrix& y)
00198     {
00199         p_ += y.p_;
00200         return *this;
00201     }
00202 
00204     ScaledIdentityMatrix& operator-= (const ScaledIdentityMatrix& y)
00205     {
00206         p_ -= y.p_;
00207         return *this;
00208     }
00209 
00211     ScaledIdentityMatrix& operator+= (const K& k)
00212     {
00213         p_ += k;
00214         return *this;
00215     }
00216 
00218     ScaledIdentityMatrix& operator-= (const K& k)
00219     {
00220         p_ -= k;
00221         return *this;
00222     }
00224     ScaledIdentityMatrix& operator*= (const K& k)
00225     {
00226         p_ *= k;
00227         return *this;
00228     }
00229 
00231     ScaledIdentityMatrix& operator/= (const K& k)
00232     {
00233         p_ /= k;
00234         return *this;
00235     }
00236 
00237     //===== linear maps
00238 
00240     template<class X, class Y>
00241     void mv (const X& x, Y& y) const
00242     {
00243 #ifdef DUNE_FMatrix_WITH_CHECKING
00244         if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00245         if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00246 #endif
00247         for (size_type i=0; i<n; ++i)
00248             y[i] = p_ * x[i];
00249     }
00250 
00252     template<class X, class Y>
00253     void mtv (const X& x, Y& y) const
00254     {
00255         mv(x, y);
00256     }
00257 
00259     template<class X, class Y>
00260     void umv (const X& x, Y& y) const
00261     {
00262 #ifdef DUNE_FMatrix_WITH_CHECKING
00263         if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00264         if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00265 #endif
00266         for (size_type i=0; i<n; ++i)
00267             y[i] += p_ * x[i];
00268     }
00269 
00271     template<class X, class Y>
00272     void umtv (const X& x, Y& y) const
00273     {
00274 #ifdef DUNE_FMatrix_WITH_CHECKING
00275         if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00276         if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00277 #endif
00278         for (size_type i=0; i<n; ++i)
00279             y[i] += p_ * x[i];
00280     }
00281 
00283     template<class X, class Y>
00284     void umhv (const X& x, Y& y) const
00285     {
00286 #ifdef DUNE_FMatrix_WITH_CHECKING
00287         if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00288         if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00289 #endif
00290         for (size_type i=0; i<n; i++)
00291             y[i] += conjugateComplex(p_)*x[i];
00292     }
00293 
00295     template<class X, class Y>
00296     void mmv (const X& x, Y& y) const
00297     {
00298 #ifdef DUNE_FMatrix_WITH_CHECKING
00299         if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00300         if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00301 #endif
00302         for (size_type i=0; i<n; ++i)
00303             y[i] -= p_ * x[i];
00304     }
00305 
00307     template<class X, class Y>
00308     void mmtv (const X& x, Y& y) const
00309     {
00310 #ifdef DUNE_FMatrix_WITH_CHECKING
00311         if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00312         if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00313 #endif
00314         for (size_type i=0; i<n; ++i)
00315             y[i] -= p_ * x[i];
00316     }
00317 
00319     template<class X, class Y>
00320     void mmhv (const X& x, Y& y) const
00321     {
00322 #ifdef DUNE_FMatrix_WITH_CHECKING
00323         if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00324         if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00325 #endif
00326         for (size_type i=0; i<n; i++)
00327             y[i] -= conjugateComplex(p_)*x[i];
00328     }
00329 
00331     template<class X, class Y>
00332     void usmv (const K& alpha, const X& x, Y& y) const
00333     {
00334 #ifdef DUNE_FMatrix_WITH_CHECKING
00335         if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00336         if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00337 #endif
00338         for (size_type i=0; i<n; i++)
00339             y[i] += alpha * p_ * x[i];
00340     }
00341 
00343     template<class X, class Y>
00344     void usmtv (const K& alpha, const X& x, Y& y) const
00345     {
00346 #ifdef DUNE_FMatrix_WITH_CHECKING
00347         if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00348         if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00349 #endif
00350         for (size_type i=0; i<n; i++)
00351             y[i] += alpha * p_ * x[i];
00352     }
00353 
00355     template<class X, class Y>
00356     void usmhv (const K& alpha, const X& x, Y& y) const
00357     {
00358 #ifdef DUNE_FMatrix_WITH_CHECKING
00359         if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00360         if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00361 #endif
00362         for (size_type i=0; i<n; i++)
00363             y[i] += alpha * conjugateComplex(p_) * x[i];
00364     }
00365 
00366     //===== norms
00367 
00369     double frobenius_norm () const
00370     {
00371         return fvmeta::sqrt(n*p_*p_);
00372     }
00373 
00375     double frobenius_norm2 () const
00376     {
00377         return n*p_*p_;
00378     }
00379 
00381     double infinity_norm () const
00382     {
00383         return std::abs(p_);
00384     }
00385 
00387     double infinity_norm_real () const
00388     {
00389         return fvmeta::absreal(p_);
00390     }
00391 
00392     //===== solve
00393 
00396     template<class V>
00397     void solve (V& x, const V& b) const
00398     {
00399         for (int i=0; i<n; i++)
00400             x[i] = b[i]/p_;
00401     }
00402 
00405     void invert()
00406     {
00407         p_ = 1/p_;
00408     }
00409 
00411     K determinant () const {
00412         return std::pow(p_,n);
00413     }
00414 
00415     //===== sizes
00416 
00418     size_type N () const
00419     {
00420         return n;
00421     }
00422 
00424     size_type M () const
00425     {
00426         return n;
00427     }
00428 
00429     //===== query
00430 
00432     bool exists (size_type i, size_type j) const
00433     {
00434 #ifdef DUNE_FMatrix_WITH_CHECKING
00435         if (i<0 || i>=n) DUNE_THROW(FMatrixError,"row index out of range");
00436         if (j<0 || j>=n) DUNE_THROW(FMatrixError,"column index out of range");
00437 #endif
00438         return i==j;
00439     }
00440 
00441     //===== conversion operator
00442 
00444     friend std::ostream& operator<< (std::ostream& s, const ScaledIdentityMatrix<K,n>& a)
00445     {
00446         for (size_type i=0; i<n; i++) {
00447             for (size_type j=0; j<n; j++)
00448                 s << ((i==j) ? a.p_ : 0) << " ";
00449             s << std::endl;
00450         }
00451         return s;
00452     }
00453 
00455     reference operator[](size_type i)
00456     {
00457         return reference(const_cast<K*>(&p_), i);
00458     }
00459 
00461     const_reference operator[](size_type i) const
00462     {
00463         return const_reference(const_cast<K*>(&p_), i);
00464     }
00465 
00467     const K& diagonal(size_type i) const
00468     {
00469         return p_;
00470     }
00471 
00473     K& diagonal(size_type i)
00474     {
00475         return p_;
00476     }
00477 
00480     const K& scalar() const
00481     {
00482         return p_;
00483     }
00484 
00487     K& scalar()
00488     {
00489         return p_;
00490     }
00491 
00492 private:
00493     // the data, very simply a single number
00494     K p_;
00495 
00496 };
00497 
00498 template<class M, class K, int n>
00499 void istl_assign_to_fmatrix(DenseMatrix<M>& fm, const ScaledIdentityMatrix<K,n>& s)
00500 {
00501     fm = K();
00502     for(int i=0; i<n; ++i)
00503         fm[i][i] = s.scalar();
00504 }
00505 
00506 } // end namespace
00507 
00508 #endif