fmatrix.hh

Go to the documentation of this file.
00001 // $Id: fmatrix.hh 5077 2008-02-06 14:05:09Z robertk $
00002 #ifndef DUNE_FMATRIX_HH
00003 #define DUNE_FMATRIX_HH
00004 
00005 #include<cmath>
00006 #include<cstddef>
00007 #include<complex>
00008 #include<iostream>
00009 #include "exceptions.hh"
00010 #include "fvector.hh"
00011 #include "precision.hh"
00012 
00013 namespace Dune {
00014    
00026   template<class K, int n, int m> class FieldMatrix;
00027   
00029   class FMatrixError : public Exception {};
00030 
00031   // conjugate komplex does nothing for non-complex types
00032   template<class K>
00033   inline K fm_ck (const K& k)
00034   {
00035         return k;
00036   }
00037 
00038   // conjugate komplex
00039   template<class K>
00040   inline std::complex<K> fm_ck (const std::complex<K>& c)
00041   {
00042         return std::complex<K>(c.real(),-c.imag());
00043   }
00044 
00053 #ifdef DUNE_EXPRESSIONTEMPLATES
00054   template<class K, int n, int m>
00055   class FieldMatrix : ExprTmpl::Matrix< FieldMatrix<K,n,m> >
00056 #else
00057   template<class K, int n, int m>
00058   class FieldMatrix
00059 #endif
00060   {
00061   public:
00062         // standard constructor and everything is sufficient ...
00063 
00064         //===== type definitions and constants
00065 
00067         typedef K field_type;
00068 
00070         typedef K block_type;
00071     
00073     typedef std::size_t size_type;
00074     
00076         enum {
00078           blocklevel = 1
00079         };
00080 
00082         typedef FieldVector<K,m> row_type; 
00083 
00085         enum {
00087           rows = n, 
00089           cols = m
00090         };
00091 
00092         //===== constructors
00095         FieldMatrix () {}
00096 
00099         FieldMatrix (const K& k)
00100         {
00101           for (size_type i=0; i<n; i++) p[i] = k;
00102         }
00103 
00104         //===== random access interface to rows of the matrix
00105 
00107         row_type& operator[] (size_type i)
00108         {
00109 #ifdef DUNE_FMatrix_WITH_CHECKING
00110           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
00111 #endif
00112           return p[i];
00113         }
00114 
00116         const row_type& operator[] (size_type i) const
00117         {
00118 #ifdef DUNE_FMatrix_WITH_CHECKING
00119           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
00120 #endif
00121           return p[i];
00122         }
00123 
00124 
00125         //===== iterator interface to rows of the matrix
00127     typedef FieldIterator<FieldMatrix<K,n,m>,row_type> Iterator;
00129     typedef Iterator iterator;
00131         typedef Iterator RowIterator;
00133         typedef typename row_type::Iterator ColIterator;
00134 
00136         Iterator begin ()
00137         {
00138           return Iterator(*this,0);
00139         }
00140           
00142         Iterator end ()
00143         {
00144           return Iterator(*this,n);
00145         }
00146 
00148         Iterator rbegin ()
00149         {
00150           return Iterator(*this,n-1);
00151         }
00152           
00154         Iterator rend ()
00155         {
00156           return Iterator(*this,-1);
00157         }
00158 
00160     typedef FieldIterator<const FieldMatrix<K,n,m>,const row_type> ConstIterator;
00162     typedef ConstIterator const_iterator;
00164         typedef ConstIterator ConstRowIterator;
00166         typedef typename row_type::ConstIterator ConstColIterator;
00167 
00169         ConstIterator begin () const
00170         {
00171           return ConstIterator(*this,0);
00172         }
00173           
00175         ConstIterator end () const
00176         {
00177           return ConstIterator(*this,n);
00178         }
00179 
00181         ConstIterator rbegin () const
00182         {
00183           return ConstIterator(*this,n-1);
00184         }
00185           
00187         ConstIterator rend () const
00188         {
00189           return ConstIterator(*this,-1);
00190         }
00191 
00192         //===== assignment from scalar
00193         FieldMatrix& operator= (const K& k)
00194         {
00195             for (size_type i=0; i<n; i++)
00196                 p[i] = k;
00197           return *this;   
00198         }
00199 
00200         //===== vector space arithmetic
00201 
00203         FieldMatrix& operator+= (const FieldMatrix& y)
00204         {
00205             for (size_type i=0; i<n; i++)
00206                 p[i] += y.p[i];
00207           return *this;
00208         }
00209 
00211         FieldMatrix& operator-= (const FieldMatrix& y)
00212         {
00213             for (size_type i=0; i<n; i++)
00214                 p[i] -= y.p[i];
00215           return *this;
00216         }
00217 
00219         FieldMatrix& operator*= (const K& k)
00220         {
00221             for (size_type i=0; i<n; i++)
00222                 p[i] *= k;
00223           return *this;
00224         }
00225 
00227         FieldMatrix& operator/= (const K& k)
00228         {
00229             for (size_type i=0; i<n; i++)
00230                 p[i] /= k;
00231           return *this;
00232         }
00233 
00234         //===== linear maps
00235    
00237         template<class X, class Y>
00238         void mv (const X& x, Y& y) const
00239         {
00240 #ifdef DUNE_FMatrix_WITH_CHECKING
00241           if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00242           if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00243 #endif
00244     for (size_type i=0; i<n; ++i)
00245     {
00246       y[i] = 0;
00247       for (size_type j=0; j<m; j++)
00248         y[i] += (*this)[i][j] * x[j];
00249     }
00250         }
00251 
00253         template<class X, class Y>
00254         void umv (const X& x, Y& y) const
00255         {
00256 #ifdef DUNE_FMatrix_WITH_CHECKING
00257           if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00258           if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00259 #endif
00260           for (size_type i=0; i<n; i++)
00261               for (size_type j=0; j<m; j++)
00262                   y[i] += (*this)[i][j] * x[j];
00263         }
00264 
00266         template<class X, class Y>
00267         void umtv (const X& x, Y& y) const
00268         {
00269 #ifdef DUNE_FMatrix_WITH_CHECKING
00270           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00271           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00272 #endif
00273           
00274           for (size_type i=0; i<n; i++)
00275                   for (size_type j=0; j<m; j++)
00276                         y[j] += p[i][j]*x[i];
00277         }
00278 
00280         template<class X, class Y>
00281         void umhv (const X& x, Y& y) const
00282         {
00283 #ifdef DUNE_FMatrix_WITH_CHECKING
00284           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00285           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00286 #endif
00287           
00288           for (size_type i=0; i<n; i++)
00289                   for (size_type j=0; j<m; j++)
00290                         y[j] += fm_ck(p[i][j])*x[i];
00291         }
00292 
00294         template<class X, class Y>
00295         void mmv (const X& x, Y& y) const
00296         {
00297 #ifdef DUNE_FMatrix_WITH_CHECKING
00298           if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00299           if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00300 #endif
00301           for (size_type i=0; i<n; i++)
00302               for (size_type j=0; j<m; j++)
00303                   y[i] -= (*this)[i][j] * x[j];
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           
00315           for (size_type i=0; i<n; i++)
00316                   for (size_type j=0; j<m; j++)
00317                         y[j] -= p[i][j]*x[i];
00318         }
00319 
00321         template<class X, class Y>
00322         void mmhv (const X& x, Y& y) const
00323         {
00324 #ifdef DUNE_FMatrix_WITH_CHECKING
00325           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00326           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00327 #endif
00328           
00329           for (size_type i=0; i<n; i++)
00330                   for (size_type j=0; j<m; j++)
00331                         y[j] -= fm_ck(p[i][j])*x[i];
00332         }
00333 
00335         template<class X, class Y>
00336         void usmv (const K& alpha, const X& x, Y& y) const
00337         {
00338 #ifdef DUNE_FMatrix_WITH_CHECKING
00339           if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00340           if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00341 #endif
00342           for (size_type i=0; i<n; i++)
00343               for (size_type j=0; j<m; j++)
00344                   y[i] += alpha * (*this)[i][j] * x[j];
00345         }
00346 
00348         template<class X, class Y>
00349         void usmtv (const K& alpha, const X& x, Y& y) const
00350         {
00351 #ifdef DUNE_FMatrix_WITH_CHECKING
00352           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00353           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00354 #endif
00355           
00356           for (size_type i=0; i<n; i++)
00357                   for (size_type j=0; j<m; j++)
00358                         y[j] += alpha*p[i][j]*x[i];
00359         }
00360 
00362         template<class X, class Y>
00363         void usmhv (const K& alpha, const X& x, Y& y) const
00364         {
00365 #ifdef DUNE_FMatrix_WITH_CHECKING
00366           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00367           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00368 #endif
00369           
00370           for (size_type i=0; i<n; i++)
00371                   for (size_type j=0; j<m; j++)
00372                         y[j] += alpha*fm_ck(p[i][j])*x[i];
00373         }
00374 
00375         //===== norms
00376 
00378     double frobenius_norm () const
00379         {
00380           double sum=0;
00381           for (size_type i=0; i<n; ++i) sum += p[i].two_norm2();
00382           return sqrt(sum);
00383         }
00384 
00386     double frobenius_norm2 () const
00387         {
00388           double sum=0;
00389           for (size_type i=0; i<n; ++i) sum += p[i].two_norm2();
00390           return sum;
00391         }
00392 
00394     double infinity_norm () const
00395         {
00396           double max=0;
00397           for (size_type i=0; i<n; ++i) max = std::max(max,p[i].one_norm());
00398           return max;
00399         }
00400 
00402         double infinity_norm_real () const
00403         {
00404           double max=0;
00405           for (size_type i=0; i<n; ++i) max = std::max(max,p[i].one_norm_real());
00406           return max;
00407         }
00408 
00409         //===== solve
00410 
00415         template<class V>
00416         void solve (V& x, const V& b) const;
00417 
00422       void invert();
00423 
00425       K determinant () const;
00426 
00428         FieldMatrix& leftmultiply (const FieldMatrix<K,n,n>& M)
00429         {
00430             FieldMatrix<K,n,m> C(*this);
00431             
00432             for (size_type i=0; i<n; i++)
00433                 for (size_type j=0; j<m; j++) {
00434                     (*this)[i][j] = 0;
00435                     for (size_type k=0; k<n; k++)
00436                         (*this)[i][j] += M[i][k]*C[k][j];
00437                 }
00438 
00439           return *this;
00440         }
00441 
00443         FieldMatrix& rightmultiply (const FieldMatrix<K,m,m>& M)
00444         {
00445             FieldMatrix<K,n,m> C(*this);
00446 
00447             for (size_type i=0; i<n; i++)
00448                 for (size_type j=0; j<m; j++) {
00449                     (*this)[i][j] = 0;
00450                     for (size_type k=0; k<m; k++)
00451                         (*this)[i][j] += C[i][k]*M[k][j];
00452                 }
00453             return *this;
00454         }
00455 
00456 
00457         //===== sizes
00458 
00460         size_type N () const
00461         {
00462           return n;
00463         }
00464 
00466         size_type M () const
00467         {
00468           return m;
00469         }
00470 
00472         size_type rowdim (size_type r) const
00473         {
00474           return 1;
00475         }
00476 
00478         size_type coldim (size_type c) const
00479         {
00480           return 1;
00481         }
00482 
00484         size_type rowdim () const
00485         {
00486           return n;
00487         }
00488 
00490         size_type coldim () const
00491         {
00492           return m;
00493         }
00494 
00495         //===== query
00496         
00498         bool exists (size_type i, size_type j) const
00499         {
00500 #ifdef DUNE_FMatrix_WITH_CHECKING
00501           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"row index out of range");
00502           if (j<0 || j>=m) DUNE_THROW(FMatrixError,"column index out of range");
00503 #endif
00504           return true;
00505         }
00506 
00507         //===== conversion operator
00508 
00510       friend std::ostream& operator<< (std::ostream& s, const FieldMatrix<K,n,m>& a)
00511       {
00512           for (size_type i=0; i<n; i++)
00513               s << a.p[i] << std::endl;
00514           return s;
00515       }
00516 
00517   private:
00518         // the data, very simply a built in array with row-wise ordering
00519         row_type p[(n > 0) ? n : 1]; 
00520 
00521     struct ElimPivot
00522     {
00523       ElimPivot(size_type pivot[n]);
00524       
00525       void swap(int i, int j);
00526       
00527       template<typename T>
00528       void operator()(const T&, int k, int i)
00529       {}
00530       
00531       size_type* pivot_;
00532     };
00533 
00534     template<typename V>
00535     struct Elim
00536     {
00537       Elim(V& rhs);
00538       
00539       void swap(int i, int j);
00540       
00541       void operator()(const typename V::field_type& factor, int k, int i);
00542 
00543       V* rhs_;
00544     };
00545     
00546     template<class Func>
00547     void luDecomposition(FieldMatrix<K,n,n>& A, Func func) const;
00548   };
00549 
00550   template<typename K, int n, int m>
00551   FieldMatrix<K,n,m>::ElimPivot::ElimPivot(size_type pivot[n])
00552     : pivot_(pivot)
00553   {
00554     for(int i=0; i < n; ++i) pivot[i]=i;
00555   }
00556 
00557   template<typename K, int n, int m>
00558   void FieldMatrix<K,n,m>::ElimPivot::swap(int i, int j)
00559   {
00560     pivot_[i]=j;
00561   }
00562   
00563   template<typename K, int n, int m>
00564   template<typename V>
00565   FieldMatrix<K,n,m>::Elim<V>::Elim(V& rhs)
00566     : rhs_(&rhs)
00567   {}
00568   
00569    template<typename K, int n, int m>
00570    template<typename V>
00571    void FieldMatrix<K,n,m>::Elim<V>::swap(int i, int j)
00572    {
00573      std::swap((*rhs_)[i], (*rhs_)[j]);
00574    }
00575 
00576   template<typename K, int n, int m>
00577   template<typename V>
00578   void FieldMatrix<K,n,m>::
00579   Elim<V>::operator()(const typename V::field_type& factor, int k, int i)
00580   {
00581     (*rhs_)[k] -= factor*(*rhs_)[i];
00582   }
00583   template<typename K, int n, int m>
00584   template<typename Func>
00585   inline void FieldMatrix<K,n,m>::luDecomposition(FieldMatrix<K,n,n>& A, Func func) const
00586   {
00587     double norm=A.infinity_norm_real(); // for relative thresholds
00588     double pivthres = std::max(FMatrixPrecision<>::absolute_limit(),norm*FMatrixPrecision<>::pivoting_limit());
00589     double singthres = std::max(FMatrixPrecision<>::absolute_limit(),norm*FMatrixPrecision<>::singular_limit());
00590   
00591     // LU decomposition of A in A
00592     for (int i=0; i<n; i++)  // loop over all rows
00593       {
00594         double pivmax=fvmeta_absreal(A[i][i]);
00595       
00596         // pivoting ?
00597         if (pivmax<pivthres)
00598           {
00599             // compute maximum of column
00600             int imax=i; double abs;
00601             for (int k=i+1; k<n; k++)
00602               if ((abs=fvmeta_absreal(A[k][i]))>pivmax)
00603                 {
00604                   pivmax = abs; imax = k;
00605                 }
00606             // swap rows
00607             if (imax!=i){
00608               for (int j=0; j<n; j++)
00609                 std::swap(A[i][j],A[imax][j]);
00610               func.swap(i, imax); // swap the pivot or rhs
00611             }
00612           }
00613         
00614         // singular ?
00615         if (pivmax<singthres)
00616           DUNE_THROW(FMatrixError,"matrix is singular");                  
00617         
00618         // eliminate
00619         for (int k=i+1; k<n; k++)
00620           {
00621             K factor = A[k][i]/A[i][i];
00622             A[k][i] = factor;
00623             for (int j=i+1; j<n; j++)
00624               A[k][j] -= factor*A[i][j];
00625             func(factor, k, i);
00626           }
00627       }
00628   }
00629 
00630     template <class K, int n, int m>
00631     template <class V>
00632     inline void FieldMatrix<K,n,m>::solve(V& x, const V& b) const
00633     {
00634         // never mind those ifs, because they get optimized away
00635         if (n!=m)
00636             DUNE_THROW(FMatrixError, "Can't solve for a " << n << "x" << m << " matrix!");
00637 
00638         // no need to implement the case 1x1, because the whole matrix class is
00639         // specialized for this
00640         
00641         if (n==2) {
00642             
00643 #ifdef DUNE_FMatrix_WITH_CHECKING
00644             K detinv = p[0][0]*p[1][1]-p[0][1]*p[1][0];
00645             if (fvmeta_absreal(detinv)<FMatrixPrecision<>::absolute_limit())
00646                 DUNE_THROW(FMatrixError,"matrix is singular");
00647             detinv = 1/detinv;
00648 #else
00649             K detinv = 1.0/(p[0][0]*p[1][1]-p[0][1]*p[1][0]);
00650 #endif
00651             
00652             x[0] = detinv*(p[1][1]*b[0]-p[0][1]*b[1]);
00653             x[1] = detinv*(p[0][0]*b[1]-p[1][0]*b[0]);
00654 
00655         } else if (n==3) {
00656 
00657             K d = determinant();
00658 #ifdef DUNE_FMatrix_WITH_CHECKING
00659             if (fvmeta_absreal(d)<FMatrixPrecision<>::absolute_limit())
00660                 DUNE_THROW(FMatrixError,"matrix is singular");
00661 #endif
00662 
00663             x[0] = (b[0]*p[1][1]*p[2][2] - b[0]*p[2][1]*p[1][2]
00664                     - b[1] *p[0][1]*p[2][2] + b[1]*p[2][1]*p[0][2]
00665                     + b[2] *p[0][1]*p[1][2] - b[2]*p[1][1]*p[0][2]) / d;
00666 
00667             x[1] = (p[0][0]*b[1]*p[2][2] - p[0][0]*b[2]*p[1][2]
00668                     - p[1][0] *b[0]*p[2][2] + p[1][0]*b[2]*p[0][2]
00669                     + p[2][0] *b[0]*p[1][2] - p[2][0]*b[1]*p[0][2]) / d;
00670 
00671             x[2] = (p[0][0]*p[1][1]*b[2] - p[0][0]*p[2][1]*b[1]
00672                     - p[1][0] *p[0][1]*b[2] + p[1][0]*p[2][1]*b[0]
00673                     + p[2][0] *p[0][1]*b[1] - p[2][0]*p[1][1]*b[0]) / d;
00674 
00675         } else {
00676 
00677           V& rhs = x; // use x to store rhs
00678           rhs = b; // copy data
00679           Elim<V> elim(rhs);
00680           FieldMatrix<K,n,n> A(*this);
00681           
00682           luDecomposition(A, elim);
00683           
00684           // backsolve
00685           for(int i=n-1; i>=0; i--){
00686             for (int j=i+1; j<n; j++)
00687               rhs[i] -= A[i][j]*x[j];
00688             x[i] = rhs[i]/A[i][i];
00689           }
00690         }       
00691     }
00692 
00693     template <class K, int n, int m>
00694     inline void FieldMatrix<K,n,m>::invert()
00695     {
00696         // never mind those ifs, because they get optimized away
00697         if (n!=m)
00698             DUNE_THROW(FMatrixError, "Can't invert a " << n << "x" << m << " matrix!");
00699 
00700         // no need to implement the case 1x1, because the whole matrix class is
00701         // specialized for this
00702 
00703         if (n==2) {
00704 
00705             K detinv = p[0][0]*p[1][1]-p[0][1]*p[1][0];
00706 #ifdef DUNE_FMatrix_WITH_CHECKING
00707             if (fvmeta_absreal(detinv)<FMatrixPrecision<>::absolute_limit())
00708                 DUNE_THROW(FMatrixError,"matrix is singular");            
00709 #endif
00710             detinv = 1/detinv;
00711 
00712             K temp=p[0][0];
00713             p[0][0] =  p[1][1]*detinv;
00714             p[0][1] = -p[0][1]*detinv;
00715             p[1][0] = -p[1][0]*detinv;
00716             p[1][1] =  temp*detinv;
00717 
00718         } else {
00719 
00720             FieldMatrix<K,n,n> A(*this);
00721             size_type pivot[n];
00722             luDecomposition(A, ElimPivot(pivot));
00723             FieldMatrix<K,n,m>& L=A;
00724             FieldMatrix<K,n,m>& U=A;
00725             
00726             // initialize inverse
00727             *this=0;
00728             
00729             for(size_type i=0; i<n; ++i)
00730               p[i][i]=1;
00731             
00732             // L Y = I; multiple right hand sides
00733             for (size_type i=0; i<n; i++){
00734               for (size_type j=0; j<i; j++)
00735                 for (size_type k=0; k<n; k++)
00736                   p[i][k] -= L[i][j]*p[j][k];
00737             }
00738   
00739             // U A^{-1} = Y
00740             for (size_type i=n; i>0;){
00741               --i;
00742               for (size_type k=0; k<n; k++){
00743                 for (size_type j=i+1; j<n; j++)
00744                   p[i][k] -= U[i][j]*p[j][k];
00745                 p[i][k] /= U[i][i];
00746               }
00747             }
00748 
00749             for(size_type i=n; i>0; ){
00750               --i;
00751               if(i!=pivot[i])
00752                 for(size_type j=0; j<n; ++j)
00753                   std::swap(p[j][pivot[i]], p[j][i]);
00754             }
00755         }
00756     }
00757 
00758     // implementation of the determinant 
00759     template <class K, int n, int m>
00760     inline K FieldMatrix<K,n,m>::determinant() const
00761     {
00762         // never mind those ifs, because they get optimized away
00763         if (n!=m)
00764             DUNE_THROW(FMatrixError, "There is no determinant for a " << n << "x" << m << " matrix!");
00765 
00766         // no need to implement the case 1x1, because the whole matrix class is
00767         // specialized for this
00768 
00769         if (n==2)
00770             return p[0][0]*p[1][1] - p[0][1]*p[1][0]; 
00771 
00772         if (n==3) {
00773              // code generated by maple 
00774             K t4  = p[0][0] * p[1][1];
00775             K t6  = p[0][0] * p[1][2];
00776             K t8  = p[0][1] * p[1][0];
00777             K t10 = p[0][2] * p[1][0];
00778             K t12 = p[0][1] * p[2][0];
00779             K t14 = p[0][2] * p[2][0];
00780         
00781             return (t4*p[2][2]-t6*p[2][1]-t8*p[2][2]+
00782                     t10*p[2][1]+t12*p[1][2]-t14*p[1][1]);
00783 
00784         }
00785         
00786         DUNE_THROW(FMatrixError, "No implementation of determinantMatrix "
00787                    << "for FieldMatrix<" << n << "," << m << "> !");
00788 
00789     }
00790 
00791 
00794   template<class K>
00795   class FieldMatrix<K,1,1>
00796   {
00797   public:
00798         // standard constructor and everything is sufficient ...
00799 
00800         //===== type definitions and constants
00801 
00803         typedef K field_type;
00804 
00806         typedef K block_type;
00807 
00809     typedef std::size_t size_type;
00810     
00812         enum {
00815           blocklevel = 1
00816         };
00817 
00819         typedef FieldVector<K,1> row_type; 
00820 
00822         enum {
00825           rows = 1,
00826       n = 1,
00829           cols = 1,
00830       m = 1
00831         };
00832 
00833         //===== constructors
00836         FieldMatrix () {}
00837 
00840         FieldMatrix (const K& k)
00841         {
00842             a = k;
00843         }
00844 
00845         //===== random access interface to rows of the matrix
00846 
00848         row_type& operator[] (size_type i)
00849         {
00850 #ifdef DUNE_FMatrix_WITH_CHECKING
00851           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
00852 #endif
00853           return a;
00854         }
00855 
00857         const row_type& operator[] (size_type i) const
00858         {
00859 #ifdef DUNE_FMatrix_WITH_CHECKING
00860           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
00861 #endif
00862           return a;
00863         }
00864 
00865         //===== iterator interface to rows of the matrix
00867     typedef FieldIterator<FieldMatrix<K,n,m>,row_type> Iterator;
00869     typedef Iterator iterator;
00871         typedef Iterator RowIterator;
00873         typedef typename row_type::Iterator ColIterator;
00874 
00876         Iterator begin ()
00877         {
00878           return Iterator(*this,0);
00879         }
00880           
00882         Iterator end ()
00883         {
00884           return Iterator(*this,n);
00885         }
00886 
00888         Iterator rbegin ()
00889         {
00890           return Iterator(*this,n-1);
00891         }
00892           
00894         Iterator rend ()
00895         {
00896           return Iterator(*this,-1);
00897         }
00898 
00900     typedef FieldIterator<const FieldMatrix<K,n,m>,const row_type> ConstIterator;
00902     typedef ConstIterator const_iterator;
00904         typedef ConstIterator ConstRowIterator;
00906         typedef typename row_type::ConstIterator ConstColIterator;
00907 
00909         ConstIterator begin () const
00910         {
00911           return ConstIterator(*this,0);
00912         }
00913           
00915         ConstIterator end () const
00916         {
00917           return ConstIterator(*this,n);
00918         }
00919 
00921         ConstIterator rbegin () const
00922         {
00923           return ConstIterator(*this,n-1);
00924         }
00925           
00927         ConstIterator rend () const
00928         {
00929           return ConstIterator(*this,-1);
00930         }
00931 
00932         //===== assignment from scalar
00933 
00934         FieldMatrix& operator= (const K& k)
00935         {
00936           a[0] = k;
00937           return *this;   
00938         }
00939 
00940         //===== vector space arithmetic
00941 
00943         FieldMatrix& operator+= (const K& y)
00944         {
00945           a[0] += y;
00946           return *this;
00947         }
00948 
00950         FieldMatrix& operator-= (const K& y)
00951         {
00952           a[0] -= y;
00953           return *this;
00954         }
00955 
00957         FieldMatrix& operator*= (const K& k)
00958         {
00959           a[0] *= k;
00960           return *this;
00961         }
00962 
00964         FieldMatrix& operator/= (const K& k)
00965         {
00966           a[0] /= k;
00967           return *this;
00968         }
00969 
00970         //===== linear maps
00971    
00973         void mv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
00974         {
00975           y.p = a[0] * x.p;
00976         }
00977 
00979         void umv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
00980         {
00981           y.p += a[0] * x.p;
00982         }
00983 
00985         void umtv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
00986         {
00987           y.p += a[0] * x.p;
00988         }
00989 
00991         void umhv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
00992         {
00993           y.p += fm_ck(a[0]) * x.p;
00994         }
00995 
00997         void mmv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
00998         {
00999           y.p -= a[0] * x.p;
01000         }
01001 
01003         void mmtv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01004         {
01005           y.p -= a[0] * x.p;
01006         }
01007 
01009         void mmhv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01010         {
01011           y.p -= fm_ck(a[0]) * x.p;
01012         }
01013 
01015         void usmv (const K& alpha, const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01016         {
01017           y.p += alpha * a[0] * x.p;
01018         }
01019 
01021         void usmtv (const K& alpha, const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01022         {
01023           y.p += alpha * a[0] * x.p;
01024         }
01025 
01027         void usmhv (const K& alpha, const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01028         {
01029           y.p += alpha * fm_ck(a[0]) * x.p;
01030         }
01031 
01032         //===== norms
01033 
01035     double frobenius_norm () const
01036         {
01037           return sqrt(fvmeta_abs2(a[0]));
01038         }
01039 
01041     double frobenius_norm2 () const
01042         {
01043           return fvmeta_abs2(a[0]);
01044         }
01045 
01047     double infinity_norm () const
01048         {
01049             return std::abs(a[0]);
01050         }
01051 
01053         double infinity_norm_real () const
01054         {
01055           return fvmeta_abs_real(a[0]);
01056         }
01057 
01058         //===== solve
01059 
01061         void solve (FieldVector<K,1>& x, const FieldVector<K,1>& b) const
01062         {
01063 #ifdef DUNE_FMatrix_WITH_CHECKING
01064         if (fvmeta_absreal(a[0][0])<FMatrixPrecision<>::absolute_limit())
01065           DUNE_THROW(FMatrixError,"matrix is singular");                  
01066 #endif
01067           x.p = b.p/a[0];
01068         }
01069 
01071         void invert ()
01072         {
01073 #ifdef DUNE_FMatrix_WITH_CHECKING
01074             if (fvmeta_absreal(a[0][0])<FMatrixPrecision<>::absolute_limit())
01075                 DUNE_THROW(FMatrixError,"matrix is singular");            
01076 #endif
01077           a[0] = 1/a[0];
01078         }
01079 
01081     K determinant () const
01082     {
01083       return std::abs(a[0]);
01084     }
01085 
01087         FieldMatrix& leftmultiply (const FieldMatrix& M)
01088         {
01089           a[0] *= M.a[0];
01090           return *this;
01091         }
01092 
01094         FieldMatrix& rightmultiply (const FieldMatrix& M)
01095         {
01096           a[0] *= M.a[0];
01097           return *this;
01098         }
01099 
01100 
01101         //===== sizes
01102 
01104         size_type N () const
01105         {
01106           return 1;
01107         }
01108 
01110         size_type M () const
01111         {
01112           return 1;
01113         }
01114 
01116         size_type rowdim (size_type r) const
01117         {
01118           return 1;
01119         }
01120 
01122         size_type coldim (size_type c) const
01123         {
01124           return 1;
01125         }
01126 
01128         size_type rowdim () const
01129         {
01130           return 1;
01131         }
01132 
01134         size_type coldim () const
01135         {
01136           return 1;
01137         }
01138 
01139         //===== query
01140         
01142         bool exists (size_type i, size_type j) const 
01143         {
01144           return i==0 && j==0;
01145         }
01146 
01147         //===== conversion operator
01148 
01149         operator K () const {return a[0];}
01150 
01151         private:
01152         // the data, just a single row with a single scalar
01153     row_type a;
01154   };
01155 
01156 namespace FMatrixHelp {
01157 
01158 
01160 template <typename K>
01161 static inline K invertMatrix (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
01162 {
01163   inverse[0][0] = 1.0/matrix[0][0];
01164   return matrix[0][0];
01165 }
01166 
01168 template <typename K>
01169 static inline K invertMatrix_retTransposed (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
01170 {
01171   return invertMatrix(matrix,inverse); 
01172 }
01173 
01174 
01176 template <typename K>
01177 static inline K invertMatrix (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
01178 {
01179   // code generated by maple 
01180   K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
01181   K det_1 = 1.0/det;
01182   inverse[0][0] =   matrix[1][1] * det_1;
01183   inverse[0][1] = - matrix[0][1] * det_1;
01184   inverse[1][0] = - matrix[1][0] * det_1;
01185   inverse[1][1] =   matrix[0][0] * det_1;
01186   return det;
01187 }
01188 
01191 template <typename K>
01192 static inline K invertMatrix_retTransposed (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
01193 {
01194   // code generated by maple 
01195   K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
01196   K det_1 = 1.0/det;
01197   inverse[0][0] =   matrix[1][1] * det_1;
01198   inverse[1][0] = - matrix[0][1] * det_1;
01199   inverse[0][1] = - matrix[1][0] * det_1;
01200   inverse[1][1] =   matrix[0][0] * det_1;
01201   return det;
01202 }
01203 
01205 template <typename K>
01206 static inline K invertMatrix (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
01207 {
01208   // code generated by maple 
01209   K t4  = matrix[0][0] * matrix[1][1];
01210   K t6  = matrix[0][0] * matrix[1][2];
01211   K t8  = matrix[0][1] * matrix[1][0];
01212   K t10 = matrix[0][2] * matrix[1][0];
01213   K t12 = matrix[0][1] * matrix[2][0];
01214   K t14 = matrix[0][2] * matrix[2][0];
01215 
01216   K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
01217            t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
01218   K t17 = 1.0/det;
01219 
01220   inverse[0][0] =  (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
01221   inverse[0][1] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
01222   inverse[0][2] =  (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
01223   inverse[1][0] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
01224   inverse[1][1] =  (matrix[0][0] * matrix[2][2] - t14) * t17;
01225   inverse[1][2] = -(t6-t10) * t17;
01226   inverse[2][0] =  (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
01227   inverse[2][1] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
01228   inverse[2][2] =  (t4-t8) * t17;
01229 
01230   return det;
01231 }
01232 
01234 template <typename K>
01235 static inline K invertMatrix_retTransposed (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
01236 {
01237   // code generated by maple 
01238   K t4  = matrix[0][0] * matrix[1][1];
01239   K t6  = matrix[0][0] * matrix[1][2];
01240   K t8  = matrix[0][1] * matrix[1][0];
01241   K t10 = matrix[0][2] * matrix[1][0];
01242   K t12 = matrix[0][1] * matrix[2][0];
01243   K t14 = matrix[0][2] * matrix[2][0];
01244 
01245   K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
01246            t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
01247   K t17 = 1.0/det;
01248 
01249   inverse[0][0] =  (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
01250   inverse[1][0] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
01251   inverse[2][0] =  (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
01252   inverse[0][1] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
01253   inverse[1][1] =  (matrix[0][0] * matrix[2][2] - t14) * t17;
01254   inverse[2][1] = -(t6-t10) * t17;
01255   inverse[0][2] =  (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
01256   inverse[1][2] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
01257   inverse[2][2] =  (t4-t8) * t17;
01258 
01259   return det;
01260 }
01261 
01263 template <typename K, int rows,int cols>
01264 static inline void multTransposedMatrix(const FieldMatrix<K,rows,cols> &matrix, FieldMatrix<K,cols,cols>& ret)
01265 {
01266   typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
01267   
01268   for(size_type i=0; i<cols; i++)
01269     for(size_type j=0; j<cols; j++)
01270       { ret[i][j]=0.0;
01271         for(size_type k=0; k<rows; k++)
01272           ret[i][j]+=matrix[k][i]*matrix[k][j];
01273       }
01274 }
01275 
01277 template <typename K, int dim>
01278 static inline void multAssign(const FieldMatrix<K,dim,dim> &matrix, const FieldVector<K,dim> & x, FieldVector<K,dim> & ret) 
01279 {
01280   typedef typename FieldMatrix<K,dim,dim>::size_type size_type;
01281   
01282   for(size_type i=0; i<dim; i++)
01283   {
01284     ret[i] = 0.0;
01285     for(size_type j=0; j<dim; j++)
01286     {
01287       ret[i] += matrix[i][j]*x[j];
01288     }
01289   }
01290 }
01291 
01293 template <typename K, int rows,int cols>
01294 static inline void multAssign(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x, FieldVector<K,rows> & ret) 
01295 {
01296   typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
01297   
01298   for(size_type i=0; i<rows; i++)
01299   {
01300     ret[i] = 0.0;
01301     for(size_type j=0; j<cols; j++)
01302     {
01303       ret[i] += matrix[i][j]*x[j];
01304     }
01305   }
01306 }
01307 
01309 template <typename K, int dim>
01310 static inline void multAssignTransposed( const FieldMatrix<K,dim,dim> &matrix, const FieldVector<K,dim> & x, FieldVector<K,dim> & ret) 
01311 {
01312   typedef typename FieldMatrix<K,dim,dim>::size_type size_type;
01313   
01314   for(size_type i=0; i<dim; i++)
01315   {
01316     ret[i] = 0.0;
01317     for(size_type j=0; j<dim; j++)
01318     {
01319       ret[i] += matrix[j][i]*x[j];
01320     }
01321   }
01322 }
01323 
01325 template <typename K, int dim>
01326 static inline FieldVector<K,dim> mult(const FieldMatrix<K,dim,dim> &matrix, const FieldVector<K,dim> & x) 
01327 {
01328   FieldVector<K,dim> ret;
01329   multAssign(matrix,x,ret);
01330   return ret; 
01331 }
01332 
01333 
01334 
01336 template <typename K, int rows, int cols>
01337 static inline FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x) 
01338 {
01339   FieldVector<K,cols> ret;
01340   typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
01341   for(size_type i=0; i<cols; i++)
01342   {
01343     ret[i] = 0.0;
01344     for(size_type j=0; j<rows; j++)
01345     {
01346       ret[i] += matrix[j][i]*x[j];
01347     }
01348   }
01349   return ret; 
01350 }
01351 
01352 } // end namespace FMatrixHelp 
01353 
01354 #ifdef DUNE_EXPRESSIONTEMPLATES
01355 template <class K, int N, int M>
01356 struct BlockType< FieldMatrix<K,N,M> >
01357 {
01358   typedef K type;
01359 };
01360 
01361 template <class K, int N, int M>
01362 struct FieldType< FieldMatrix<K,N,M> >
01363 {
01364   typedef K type;
01365 };
01366 #endif // DUNE_EXPRESSIONTEMPLATES
01367 
01370 } // end namespace
01371 
01372 #endif

Generated on 6 Nov 2008 with Doxygen (ver 1.5.6) [logfile].