fmatrix.hh

Go to the documentation of this file.
00001 // $Id: fmatrix.hh 5684 2009-11-04 08:52:28Z mnolte $
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 #include "static_assert.hh"
00013 
00014 namespace Dune {
00015    
00027   template<class K, int n, int m> class FieldMatrix;
00028 
00029   template<class K, int n, int m, typename T>
00030   void istl_assign_to_fmatrix(FieldMatrix<K,n,m>& f, const T& t)
00031   {
00032     DUNE_THROW(NotImplemented, "You need to specialise this function for type T!");
00033   }
00034 
00035   namespace
00036   {
00037     template<bool b>
00038     struct Assigner
00039     {
00040       template<class K, int n, int m, class T>
00041       static void assign(FieldMatrix<K,n,m>& fm, const T& t)
00042       {
00043         istl_assign_to_fmatrix(fm, t);
00044       }
00045       
00046     };
00047     
00048 
00049     template<>
00050     struct Assigner<true>
00051     {
00052       template<class K, int n, int m, class T>
00053       static void assign(FieldMatrix<K,n,m>& fm, const T& t)
00054       {
00055         fm = static_cast<const K>(t);
00056       }
00057     };  
00058   }
00059   
00061   class FMatrixError : public Exception {};
00062 
00063   // conjugate komplex does nothing for non-complex types
00064   template<class K>
00065   inline K fm_ck (const K& k)
00066   {
00067         return k;
00068   }
00069 
00070   // conjugate komplex
00071   template<class K>
00072   inline std::complex<K> fm_ck (const std::complex<K>& c)
00073   {
00074         return std::complex<K>(c.real(),-c.imag());
00075   }
00076 
00085 #ifdef DUNE_EXPRESSIONTEMPLATES
00086   template<class K, int n, int m>
00087   class FieldMatrix : ExprTmpl::Matrix< FieldMatrix<K,n,m> >
00088 #else
00089   template<class K, int n, int m>
00090   class FieldMatrix
00091 #endif
00092   {
00093   public:
00094         // standard constructor and everything is sufficient ...
00095 
00096         //===== type definitions and constants
00097 
00099         typedef K field_type;
00100 
00102         typedef K block_type;
00103     
00105     typedef std::size_t size_type;
00106     
00108         enum {
00110           blocklevel = 1
00111         };
00112 
00114         typedef FieldVector<K,m> row_type; 
00115 
00117         enum {
00119           rows = n, 
00121           cols = m
00122         };
00123 
00124         //===== constructors
00127         FieldMatrix () {}
00128 
00131         FieldMatrix (const K& k)
00132         {
00133           for (size_type i=0; i<n; i++) p[i] = k;
00134         }
00135 
00136     template<typename T>
00137     explicit FieldMatrix( const T& t)
00138     {
00139       Assigner<Conversion<T,K>::exists>::assign(*this, t);
00140     }
00141     
00142         //===== random access interface to rows of the matrix
00143 
00145         row_type& operator[] (size_type i)
00146         {
00147 #ifdef DUNE_FMatrix_WITH_CHECKING
00148           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
00149 #endif
00150           return p[i];
00151         }
00152 
00154         const row_type& operator[] (size_type i) const
00155         {
00156 #ifdef DUNE_FMatrix_WITH_CHECKING
00157           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
00158 #endif
00159           return p[i];
00160         }
00161 
00162 
00163         //===== iterator interface to rows of the matrix
00165     typedef FieldIterator<FieldMatrix<K,n,m>,row_type> Iterator;
00167     typedef Iterator iterator;
00169         typedef Iterator RowIterator;
00171         typedef typename row_type::Iterator ColIterator;
00172 
00174         Iterator begin ()
00175         {
00176           return Iterator(*this,0);
00177         }
00178           
00180         Iterator end ()
00181         {
00182           return Iterator(*this,n);
00183         }
00184 
00186         Iterator rbegin ()
00187         {
00188           return Iterator(*this,n-1);
00189         }
00190           
00192         Iterator rend ()
00193         {
00194           return Iterator(*this,-1);
00195         }
00196 
00198     typedef FieldIterator<const FieldMatrix<K,n,m>,const row_type> ConstIterator;
00200     typedef ConstIterator const_iterator;
00202         typedef ConstIterator ConstRowIterator;
00204         typedef typename row_type::ConstIterator ConstColIterator;
00205 
00207         ConstIterator begin () const
00208         {
00209           return ConstIterator(*this,0);
00210         }
00211           
00213         ConstIterator end () const
00214         {
00215           return ConstIterator(*this,n);
00216         }
00217 
00219         ConstIterator rbegin () const
00220         {
00221           return ConstIterator(*this,n-1);
00222         }
00223           
00225         ConstIterator rend () const
00226         {
00227           return ConstIterator(*this,-1);
00228         }
00229 
00230         //===== assignment from scalar
00231         FieldMatrix& operator= (const K& k)
00232         {
00233             for (size_type i=0; i<n; i++)
00234                 p[i] = k;
00235           return *this;   
00236         }
00237 
00238     template<typename T>
00239     FieldMatrix& operator= ( const T& t)
00240     {
00241       Assigner<Conversion<T,K>::exists>::assign(*this, t);
00242       return *this;
00243     }
00244         //===== vector space arithmetic
00245 
00247         FieldMatrix& operator+= (const FieldMatrix& y)
00248         {
00249             for (size_type i=0; i<n; i++)
00250                 p[i] += y.p[i];
00251           return *this;
00252         }
00253 
00255         FieldMatrix& operator-= (const FieldMatrix& y)
00256         {
00257             for (size_type i=0; i<n; i++)
00258                 p[i] -= y.p[i];
00259           return *this;
00260         }
00261 
00263         FieldMatrix& operator*= (const K& k)
00264         {
00265             for (size_type i=0; i<n; i++)
00266                 p[i] *= k;
00267           return *this;
00268         }
00269 
00271         FieldMatrix& operator/= (const K& k)
00272         {
00273             for (size_type i=0; i<n; i++)
00274                 p[i] /= k;
00275           return *this;
00276         }
00277 
00279         FieldMatrix &axpy ( const K &k, const FieldMatrix &y )
00280         {
00281           for( size_type i = 0; i < n; ++i )
00282             p[ i ].axpy( k, y[ i ] );
00283           return *this;
00284         }
00285 
00286         //===== linear maps
00287    
00289         template<class X, class Y>
00290         void mv (const X& x, Y& y) const
00291         {
00292 #ifdef DUNE_FMatrix_WITH_CHECKING
00293           if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00294           if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00295 #endif
00296     for (size_type i=0; i<n; ++i)
00297     {
00298       y[i] = 0;
00299       for (size_type j=0; j<m; j++)
00300         y[i] += (*this)[i][j] * x[j];
00301     }
00302         }
00303 
00304 
00306     template< class X, class Y >
00307     void mtv ( const X &x, Y &y ) const
00308     {
00309 #ifdef DUNE_FMatrix_WITH_CHECKING
00310       assert( &x != &y );
00311       if( x.N() != N() )
00312         DUNE_THROW( FMatrixError, "Index out of range." );
00313       if( y.N() != M() )
00314         DUNE_THROW( FMatrixError, "Index out of range." );
00315 #endif
00316       for( size_type i = 0; i < cols; ++i )
00317       {
00318         y[ i ] = 0;
00319         for( size_type j = 0; j < rows; ++j )
00320           y[ i ] += (*this)[ j ][ i ] * x[ j ];
00321       }
00322     }
00323 
00324 
00326         template<class X, class Y>
00327         void umv (const X& x, Y& y) const
00328         {
00329 #ifdef DUNE_FMatrix_WITH_CHECKING
00330           if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00331           if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00332 #endif
00333           for (size_type i=0; i<n; i++)
00334               for (size_type j=0; j<m; j++)
00335                   y[i] += (*this)[i][j] * x[j];
00336         }
00337 
00339         template<class X, class Y>
00340         void umtv (const X& x, Y& y) const
00341         {
00342 #ifdef DUNE_FMatrix_WITH_CHECKING
00343           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00344           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00345 #endif
00346           
00347           for (size_type i=0; i<n; i++)
00348                   for (size_type j=0; j<m; j++)
00349                         y[j] += p[i][j]*x[i];
00350         }
00351 
00353         template<class X, class Y>
00354         void umhv (const X& x, Y& y) const
00355         {
00356 #ifdef DUNE_FMatrix_WITH_CHECKING
00357           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00358           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00359 #endif
00360           
00361           for (size_type i=0; i<n; i++)
00362                   for (size_type j=0; j<m; j++)
00363                         y[j] += fm_ck(p[i][j])*x[i];
00364         }
00365 
00367         template<class X, class Y>
00368         void mmv (const X& x, Y& y) const
00369         {
00370 #ifdef DUNE_FMatrix_WITH_CHECKING
00371           if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00372           if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00373 #endif
00374           for (size_type i=0; i<n; i++)
00375               for (size_type j=0; j<m; j++)
00376                   y[i] -= (*this)[i][j] * x[j];
00377         }
00378 
00380         template<class X, class Y>
00381         void mmtv (const X& x, Y& y) const
00382         {
00383 #ifdef DUNE_FMatrix_WITH_CHECKING
00384           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00385           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00386 #endif
00387           
00388           for (size_type i=0; i<n; i++)
00389                   for (size_type j=0; j<m; j++)
00390                         y[j] -= p[i][j]*x[i];
00391         }
00392 
00394         template<class X, class Y>
00395         void mmhv (const X& x, Y& y) const
00396         {
00397 #ifdef DUNE_FMatrix_WITH_CHECKING
00398           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00399           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00400 #endif
00401           
00402           for (size_type i=0; i<n; i++)
00403                   for (size_type j=0; j<m; j++)
00404                         y[j] -= fm_ck(p[i][j])*x[i];
00405         }
00406 
00408         template<class X, class Y>
00409         void usmv (const K& alpha, const X& x, Y& y) const
00410         {
00411 #ifdef DUNE_FMatrix_WITH_CHECKING
00412           if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00413           if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00414 #endif
00415           for (size_type i=0; i<n; i++)
00416               for (size_type j=0; j<m; j++)
00417                   y[i] += alpha * (*this)[i][j] * x[j];
00418         }
00419 
00421         template<class X, class Y>
00422         void usmtv (const K& alpha, const X& x, Y& y) const
00423         {
00424 #ifdef DUNE_FMatrix_WITH_CHECKING
00425           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00426           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00427 #endif
00428           
00429           for (size_type i=0; i<n; i++)
00430                   for (size_type j=0; j<m; j++)
00431                         y[j] += alpha*p[i][j]*x[i];
00432         }
00433 
00435         template<class X, class Y>
00436         void usmhv (const K& alpha, const X& x, Y& y) const
00437         {
00438 #ifdef DUNE_FMatrix_WITH_CHECKING
00439           if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
00440           if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
00441 #endif
00442           
00443           for (size_type i=0; i<n; i++)
00444                   for (size_type j=0; j<m; j++)
00445                         y[j] += alpha*fm_ck(p[i][j])*x[i];
00446         }
00447 
00448         //===== norms
00449 
00451     double frobenius_norm () const
00452         {
00453           double sum=0;
00454           for (size_type i=0; i<n; ++i) sum += p[i].two_norm2();
00455           return sqrt(sum);
00456         }
00457 
00459     double frobenius_norm2 () const
00460         {
00461           double sum=0;
00462           for (size_type i=0; i<n; ++i) sum += p[i].two_norm2();
00463           return sum;
00464         }
00465 
00467     double infinity_norm () const
00468         {
00469           double max=0;
00470           for (size_type i=0; i<n; ++i) max = std::max(max,p[i].one_norm());
00471           return max;
00472         }
00473 
00475         double infinity_norm_real () const
00476         {
00477           double max=0;
00478           for (size_type i=0; i<n; ++i) max = std::max(max,p[i].one_norm_real());
00479           return max;
00480         }
00481 
00482         //===== solve
00483 
00488         template<class V>
00489         void solve (V& x, const V& b) const;
00490 
00495       void invert();
00496 
00498       K determinant () const;
00499 
00501         FieldMatrix& leftmultiply (const FieldMatrix<K,n,n>& M)
00502         {
00503             FieldMatrix<K,n,m> C(*this);
00504             
00505             for (size_type i=0; i<n; i++)
00506                 for (size_type j=0; j<m; j++) {
00507                     (*this)[i][j] = 0;
00508                     for (size_type k=0; k<n; k++)
00509                         (*this)[i][j] += M[i][k]*C[k][j];
00510                 }
00511 
00512           return *this;
00513         }
00514 
00516         FieldMatrix& rightmultiply (const FieldMatrix<K,m,m>& M)
00517         {
00518             FieldMatrix<K,n,m> C(*this);
00519 
00520             for (size_type i=0; i<n; i++)
00521                 for (size_type j=0; j<m; j++) {
00522                     (*this)[i][j] = 0;
00523                     for (size_type k=0; k<m; k++)
00524                         (*this)[i][j] += C[i][k]*M[k][j];
00525                 }
00526             return *this;
00527         }
00528 
00529 
00530         //===== sizes
00531 
00533         size_type N () const
00534         {
00535           return n;
00536         }
00537 
00539         size_type M () const
00540         {
00541           return m;
00542         }
00543 
00544         //===== query
00545         
00547         bool exists (size_type i, size_type j) const
00548         {
00549 #ifdef DUNE_FMatrix_WITH_CHECKING
00550           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"row index out of range");
00551           if (j<0 || j>=m) DUNE_THROW(FMatrixError,"column index out of range");
00552 #endif
00553           return true;
00554         }
00555 
00556         //===== conversion operator
00557 
00559       friend std::ostream& operator<< (std::ostream& s, const FieldMatrix<K,n,m>& a)
00560       {
00561           for (size_type i=0; i<n; i++)
00562               s << a.p[i] << std::endl;
00563           return s;
00564       }
00565 
00566   private:
00567         // the data, very simply a built in array with row-wise ordering
00568         row_type p[(n > 0) ? n : 1]; 
00569 
00570     struct ElimPivot
00571     {
00572       ElimPivot(size_type pivot[n]);
00573       
00574       void swap(int i, int j);
00575       
00576       template<typename T>
00577       void operator()(const T&, int k, int i)
00578       {}
00579       
00580       size_type* pivot_;
00581     };
00582 
00583     template<typename V>
00584     struct Elim
00585     {
00586       Elim(V& rhs);
00587       
00588       void swap(int i, int j);
00589       
00590       void operator()(const typename V::field_type& factor, int k, int i);
00591 
00592       V* rhs_;
00593     };
00594     
00595     template<class Func>
00596     void luDecomposition(FieldMatrix<K,n,n>& A, Func func) const;
00597   };
00598 
00599   template<typename K, int n, int m>
00600   FieldMatrix<K,n,m>::ElimPivot::ElimPivot(size_type pivot[n])
00601     : pivot_(pivot)
00602   {
00603     for(int i=0; i < n; ++i) pivot[i]=i;
00604   }
00605 
00606   template<typename K, int n, int m>
00607   void FieldMatrix<K,n,m>::ElimPivot::swap(int i, int j)
00608   {
00609     pivot_[i]=j;
00610   }
00611   
00612   template<typename K, int n, int m>
00613   template<typename V>
00614   FieldMatrix<K,n,m>::Elim<V>::Elim(V& rhs)
00615     : rhs_(&rhs)
00616   {}
00617   
00618    template<typename K, int n, int m>
00619    template<typename V>
00620    void FieldMatrix<K,n,m>::Elim<V>::swap(int i, int j)
00621    {
00622      std::swap((*rhs_)[i], (*rhs_)[j]);
00623    }
00624 
00625   template<typename K, int n, int m>
00626   template<typename V>
00627   void FieldMatrix<K,n,m>::
00628   Elim<V>::operator()(const typename V::field_type& factor, int k, int i)
00629   {
00630     (*rhs_)[k] -= factor*(*rhs_)[i];
00631   }
00632   template<typename K, int n, int m>
00633   template<typename Func>
00634   inline void FieldMatrix<K,n,m>::luDecomposition(FieldMatrix<K,n,n>& A, Func func) const
00635   {
00636     double norm=A.infinity_norm_real(); // for relative thresholds
00637     double pivthres = std::max(FMatrixPrecision<>::absolute_limit(),norm*FMatrixPrecision<>::pivoting_limit());
00638     double singthres = std::max(FMatrixPrecision<>::absolute_limit(),norm*FMatrixPrecision<>::singular_limit());
00639   
00640     // LU decomposition of A in A
00641     for (int i=0; i<n; i++)  // loop over all rows
00642       {
00643         double pivmax=fvmeta_absreal(A[i][i]);
00644       
00645         // pivoting ?
00646         if (pivmax<pivthres)
00647           {
00648             // compute maximum of column
00649             int imax=i; double abs;
00650             for (int k=i+1; k<n; k++)
00651               if ((abs=fvmeta_absreal(A[k][i]))>pivmax)
00652                 {
00653                   pivmax = abs; imax = k;
00654                 }
00655             // swap rows
00656             if (imax!=i){
00657               for (int j=0; j<n; j++)
00658                 std::swap(A[i][j],A[imax][j]);
00659               func.swap(i, imax); // swap the pivot or rhs
00660             }
00661           }
00662         
00663         // singular ?
00664         if (pivmax<singthres)
00665           DUNE_THROW(FMatrixError,"matrix is singular");                  
00666         
00667         // eliminate
00668         for (int k=i+1; k<n; k++)
00669           {
00670             K factor = A[k][i]/A[i][i];
00671             A[k][i] = factor;
00672             for (int j=i+1; j<n; j++)
00673               A[k][j] -= factor*A[i][j];
00674             func(factor, k, i);
00675           }
00676       }
00677   }
00678 
00679     template <class K, int n, int m>
00680     template <class V>
00681     inline void FieldMatrix<K,n,m>::solve(V& x, const V& b) const
00682     {
00683         // never mind those ifs, because they get optimized away
00684         if (n!=m)
00685             DUNE_THROW(FMatrixError, "Can't solve for a " << n << "x" << m << " matrix!");
00686 
00687         // no need to implement the case 1x1, because the whole matrix class is
00688         // specialized for this
00689         
00690         if (n==2) {
00691             
00692 #ifdef DUNE_FMatrix_WITH_CHECKING
00693             K detinv = p[0][0]*p[1][1]-p[0][1]*p[1][0];
00694             if (fvmeta_absreal(detinv)<FMatrixPrecision<>::absolute_limit())
00695                 DUNE_THROW(FMatrixError,"matrix is singular");
00696             detinv = 1/detinv;
00697 #else
00698             K detinv = 1.0/(p[0][0]*p[1][1]-p[0][1]*p[1][0]);
00699 #endif
00700             
00701             x[0] = detinv*(p[1][1]*b[0]-p[0][1]*b[1]);
00702             x[1] = detinv*(p[0][0]*b[1]-p[1][0]*b[0]);
00703 
00704         } else if (n==3) {
00705 
00706             K d = determinant();
00707 #ifdef DUNE_FMatrix_WITH_CHECKING
00708             if (fvmeta_absreal(d)<FMatrixPrecision<>::absolute_limit())
00709                 DUNE_THROW(FMatrixError,"matrix is singular");
00710 #endif
00711 
00712             x[0] = (b[0]*p[1][1]*p[2][2] - b[0]*p[2][1]*p[1][2]
00713                     - b[1] *p[0][1]*p[2][2] + b[1]*p[2][1]*p[0][2]
00714                     + b[2] *p[0][1]*p[1][2] - b[2]*p[1][1]*p[0][2]) / d;
00715 
00716             x[1] = (p[0][0]*b[1]*p[2][2] - p[0][0]*b[2]*p[1][2]
00717                     - p[1][0] *b[0]*p[2][2] + p[1][0]*b[2]*p[0][2]
00718                     + p[2][0] *b[0]*p[1][2] - p[2][0]*b[1]*p[0][2]) / d;
00719 
00720             x[2] = (p[0][0]*p[1][1]*b[2] - p[0][0]*p[2][1]*b[1]
00721                     - p[1][0] *p[0][1]*b[2] + p[1][0]*p[2][1]*b[0]
00722                     + p[2][0] *p[0][1]*b[1] - p[2][0]*p[1][1]*b[0]) / d;
00723 
00724         } else {
00725 
00726           V& rhs = x; // use x to store rhs
00727           rhs = b; // copy data
00728           Elim<V> elim(rhs);
00729           FieldMatrix<K,n,n> A(*this);
00730           
00731           luDecomposition(A, elim);
00732           
00733           // backsolve
00734           for(int i=n-1; i>=0; i--){
00735             for (int j=i+1; j<n; j++)
00736               rhs[i] -= A[i][j]*x[j];
00737             x[i] = rhs[i]/A[i][i];
00738           }
00739         }       
00740     }
00741 
00742     template <class K, int n, int m>
00743     inline void FieldMatrix<K,n,m>::invert()
00744     {
00745         // never mind those ifs, because they get optimized away
00746         if (n!=m)
00747             DUNE_THROW(FMatrixError, "Can't invert a " << n << "x" << m << " matrix!");
00748 
00749         // no need to implement the case 1x1, because the whole matrix class is
00750         // specialized for this
00751 
00752         if (n==2) {
00753 
00754             K detinv = p[0][0]*p[1][1]-p[0][1]*p[1][0];
00755 #ifdef DUNE_FMatrix_WITH_CHECKING
00756             if (fvmeta_absreal(detinv)<FMatrixPrecision<>::absolute_limit())
00757                 DUNE_THROW(FMatrixError,"matrix is singular");            
00758 #endif
00759             detinv = 1/detinv;
00760 
00761             K temp=p[0][0];
00762             p[0][0] =  p[1][1]*detinv;
00763             p[0][1] = -p[0][1]*detinv;
00764             p[1][0] = -p[1][0]*detinv;
00765             p[1][1] =  temp*detinv;
00766 
00767         } else {
00768 
00769             FieldMatrix<K,n,n> A(*this);
00770             size_type pivot[n];
00771             luDecomposition(A, ElimPivot(pivot));
00772             FieldMatrix<K,n,m>& L=A;
00773             FieldMatrix<K,n,m>& U=A;
00774             
00775             // initialize inverse
00776             *this=K();
00777             
00778             for(size_type i=0; i<n; ++i)
00779               p[i][i]=1;
00780             
00781             // L Y = I; multiple right hand sides
00782             for (size_type i=0; i<n; i++){
00783               for (size_type j=0; j<i; j++)
00784                 for (size_type k=0; k<n; k++)
00785                   p[i][k] -= L[i][j]*p[j][k];
00786             }
00787   
00788             // U A^{-1} = Y
00789             for (size_type i=n; i>0;){
00790               --i;
00791               for (size_type k=0; k<n; k++){
00792                 for (size_type j=i+1; j<n; j++)
00793                   p[i][k] -= U[i][j]*p[j][k];
00794                 p[i][k] /= U[i][i];
00795               }
00796             }
00797 
00798             for(size_type i=n; i>0; ){
00799               --i;
00800               if(i!=pivot[i])
00801                 for(size_type j=0; j<n; ++j)
00802                   std::swap(p[j][pivot[i]], p[j][i]);
00803             }
00804         }
00805     }
00806 
00807     // implementation of the determinant 
00808     template <class K, int n, int m>
00809     inline K FieldMatrix<K,n,m>::determinant() const
00810     {
00811         // never mind those ifs, because they get optimized away
00812         if (n!=m)
00813             DUNE_THROW(FMatrixError, "There is no determinant for a " << n << "x" << m << " matrix!");
00814 
00815         // no need to implement the case 1x1, because the whole matrix class is
00816         // specialized for this
00817 
00818         if (n==2)
00819             return p[0][0]*p[1][1] - p[0][1]*p[1][0]; 
00820 
00821         if (n==3) {
00822              // code generated by maple 
00823             K t4  = p[0][0] * p[1][1];
00824             K t6  = p[0][0] * p[1][2];
00825             K t8  = p[0][1] * p[1][0];
00826             K t10 = p[0][2] * p[1][0];
00827             K t12 = p[0][1] * p[2][0];
00828             K t14 = p[0][2] * p[2][0];
00829         
00830             return (t4*p[2][2]-t6*p[2][1]-t8*p[2][2]+
00831                     t10*p[2][1]+t12*p[1][2]-t14*p[1][1]);
00832 
00833         }
00834         
00835         DUNE_THROW(FMatrixError, "No implementation of determinantMatrix "
00836                    << "for FieldMatrix<" << n << "," << m << "> !");
00837 
00838     }
00839 
00840 
00843   template<class K>
00844   class FieldMatrix<K,1,1>
00845   {
00846   public:
00847         // standard constructor and everything is sufficient ...
00848 
00849         //===== type definitions and constants
00850 
00852         typedef K field_type;
00853 
00855         typedef K block_type;
00856 
00858     typedef std::size_t size_type;
00859     
00861         enum {
00864           blocklevel = 1
00865         };
00866 
00868         typedef FieldVector<K,1> row_type; 
00869 
00871         enum {
00874           rows = 1,
00875       n = 1,
00878           cols = 1,
00879       m = 1
00880         };
00881 
00882         //===== constructors
00885         FieldMatrix () {}
00886 
00889         FieldMatrix (const K& k)
00890         {
00891             a = k;
00892         }
00893     template<typename T>
00894     explicit FieldMatrix( const T& t)
00895     {
00896       Assigner<Conversion<T,K>::exists>::assign(*this, t);
00897     }
00898         //===== random access interface to rows of the matrix
00899 
00901         row_type& operator[] (size_type i)
00902         {
00903 #ifdef DUNE_FMatrix_WITH_CHECKING
00904           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
00905 #endif
00906           return a;
00907         }
00908 
00910         const row_type& operator[] (size_type i) const
00911         {
00912 #ifdef DUNE_FMatrix_WITH_CHECKING
00913           if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
00914 #endif
00915           return a;
00916         }
00917 
00918         //===== iterator interface to rows of the matrix
00920     typedef FieldIterator<FieldMatrix<K,n,m>,row_type> Iterator;
00922     typedef Iterator iterator;
00924         typedef Iterator RowIterator;
00926         typedef typename row_type::Iterator ColIterator;
00927 
00929         Iterator begin ()
00930         {
00931           return Iterator(*this,0);
00932         }
00933           
00935         Iterator end ()
00936         {
00937           return Iterator(*this,n);
00938         }
00939 
00941         Iterator rbegin ()
00942         {
00943           return Iterator(*this,n-1);
00944         }
00945           
00947         Iterator rend ()
00948         {
00949           return Iterator(*this,-1);
00950         }
00951 
00953     typedef FieldIterator<const FieldMatrix<K,n,m>,const row_type> ConstIterator;
00955     typedef ConstIterator const_iterator;
00957         typedef ConstIterator ConstRowIterator;
00959         typedef typename row_type::ConstIterator ConstColIterator;
00960 
00962         ConstIterator begin () const
00963         {
00964           return ConstIterator(*this,0);
00965         }
00966           
00968         ConstIterator end () const
00969         {
00970           return ConstIterator(*this,n);
00971         }
00972 
00974         ConstIterator rbegin () const
00975         {
00976           return ConstIterator(*this,n-1);
00977         }
00978           
00980         ConstIterator rend () const
00981         {
00982           return ConstIterator(*this,-1);
00983         }
00984 
00985         //===== assignment from scalar
00986 
00987         FieldMatrix& operator= (const K& k)
00988         {
00989           a[0] = k;
00990           return *this;   
00991         }
00992 
00993     template<typename T>
00994     FieldMatrix& operator= ( const T& t)
00995     {
00996       Assigner<Conversion<T,K>::exists>::assign(*this, t);
00997       return *this;
00998     }
00999 
01000         //===== vector space arithmetic
01001 
01003         FieldMatrix& operator+= (const K& y)
01004         {
01005           a[0] += y;
01006           return *this;
01007         }
01008 
01010         FieldMatrix& operator-= (const K& y)
01011         {
01012           a[0] -= y;
01013           return *this;
01014         }
01015 
01017         FieldMatrix& operator*= (const K& k)
01018         {
01019           a[0] *= k;
01020           return *this;
01021         }
01022 
01024         FieldMatrix& operator/= (const K& k)
01025         {
01026           a[0] /= k;
01027           return *this;
01028         }
01029 
01031         FieldMatrix &axpy ( const K &k, const FieldMatrix &y )
01032         {
01033           a[ 0 ] += k * y.a[ 0 ];
01034           return *this;
01035         }
01036 
01037         //===== linear maps
01038    
01040         void mv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01041         {
01042           y.p = a[0] * x.p;
01043         }
01044 
01046         void mtv ( const FieldVector< K, 1 > &x, FieldVector< K, 1 > &y ) const
01047         {
01048           y.p = a[ 0 ] * x.p;
01049         }
01050 
01052         void umv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01053         {
01054           y.p += a[0] * x.p;
01055         }
01056 
01058         void umtv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01059         {
01060           y.p += a[0] * x.p;
01061         }
01062 
01064         void umhv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01065         {
01066           y.p += fm_ck(a[0]) * x.p;
01067         }
01068 
01070         void mmv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01071         {
01072           y.p -= a[0] * x.p;
01073         }
01074 
01076         void mmtv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01077         {
01078           y.p -= a[0] * x.p;
01079         }
01080 
01082         void mmhv (const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01083         {
01084           y.p -= fm_ck(a[0]) * x.p;
01085         }
01086 
01088         void usmv (const K& alpha, const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01089         {
01090           y.p += alpha * a[0] * x.p;
01091         }
01092 
01094         void usmtv (const K& alpha, const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01095         {
01096           y.p += alpha * a[0] * x.p;
01097         }
01098 
01100         void usmhv (const K& alpha, const FieldVector<K,1>& x, FieldVector<K,1>& y) const
01101         {
01102           y.p += alpha * fm_ck(a[0]) * x.p;
01103         }
01104 
01105         //===== norms
01106 
01108     double frobenius_norm () const
01109         {
01110           return sqrt(fvmeta_abs2(a[0]));
01111         }
01112 
01114     double frobenius_norm2 () const
01115         {
01116           return fvmeta_abs2(a[0]);
01117         }
01118 
01120     double infinity_norm () const
01121         {
01122             return std::abs(a[0]);
01123         }
01124 
01126         double infinity_norm_real () const
01127         {
01128           return fvmeta_abs_real(a[0]);
01129         }
01130 
01131         //===== solve
01132 
01134         void solve (FieldVector<K,1>& x, const FieldVector<K,1>& b) const
01135         {
01136 #ifdef DUNE_FMatrix_WITH_CHECKING
01137         if (fvmeta_absreal(a[0][0])<FMatrixPrecision<>::absolute_limit())
01138           DUNE_THROW(FMatrixError,"matrix is singular");                  
01139 #endif
01140           x.p = b.p/a[0];
01141         }
01142 
01144         void invert ()
01145         {
01146 #ifdef DUNE_FMatrix_WITH_CHECKING
01147             if (fvmeta_absreal(a[0][0])<FMatrixPrecision<>::absolute_limit())
01148                 DUNE_THROW(FMatrixError,"matrix is singular");            
01149 #endif
01150           a[0] = 1/a[0];
01151         }
01152 
01154     K determinant () const
01155     {
01156       return a[ 0 ];
01157     }
01158 
01160         FieldMatrix& leftmultiply (const FieldMatrix& M)
01161         {
01162           a[0] *= M.a[0];
01163           return *this;
01164         }
01165 
01167         FieldMatrix& rightmultiply (const FieldMatrix& M)
01168         {
01169           a[0] *= M.a[0];
01170           return *this;
01171         }
01172 
01173 
01174         //===== sizes
01175 
01177         size_type N () const
01178         {
01179           return 1;
01180         }
01181 
01183         size_type M () const
01184         {
01185           return 1;
01186         }
01187 
01189         size_type rowdim (size_type r) const
01190         {
01191           return 1;
01192         }
01193 
01195         size_type coldim (size_type c) const
01196         {
01197           return 1;
01198         }
01199 
01201         size_type rowdim () const
01202         {
01203           return 1;
01204         }
01205 
01207         size_type coldim () const
01208         {
01209           return 1;
01210         }
01211 
01212         //===== query
01213         
01215         bool exists (size_type i, size_type j) const 
01216         {
01217           return i==0 && j==0;
01218         }
01219 
01220         //===== conversion operator
01221 
01222         operator K () const {return a[0];}
01223 
01224         private:
01225         // the data, just a single row with a single scalar
01226     row_type a;
01227     
01228   };
01229 
01230 namespace FMatrixHelp {
01231 
01232 
01234 template <typename K>
01235 static inline K invertMatrix (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
01236 {
01237   inverse[0][0] = 1.0/matrix[0][0];
01238   return matrix[0][0];
01239 }
01240 
01242 template <typename K>
01243 static inline K invertMatrix_retTransposed (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
01244 {
01245   return invertMatrix(matrix,inverse); 
01246 }
01247 
01248 
01250 template <typename K>
01251 static inline K invertMatrix (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
01252 {
01253   // code generated by maple 
01254   K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
01255   K det_1 = 1.0/det;
01256   inverse[0][0] =   matrix[1][1] * det_1;
01257   inverse[0][1] = - matrix[0][1] * det_1;
01258   inverse[1][0] = - matrix[1][0] * det_1;
01259   inverse[1][1] =   matrix[0][0] * det_1;
01260   return det;
01261 }
01262 
01265 template <typename K>
01266 static inline K invertMatrix_retTransposed (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
01267 {
01268   // code generated by maple 
01269   K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
01270   K det_1 = 1.0/det;
01271   inverse[0][0] =   matrix[1][1] * det_1;
01272   inverse[1][0] = - matrix[0][1] * det_1;
01273   inverse[0][1] = - matrix[1][0] * det_1;
01274   inverse[1][1] =   matrix[0][0] * det_1;
01275   return det;
01276 }
01277 
01279 template <typename K>
01280 static inline K invertMatrix (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
01281 {
01282   // code generated by maple 
01283   K t4  = matrix[0][0] * matrix[1][1];
01284   K t6  = matrix[0][0] * matrix[1][2];
01285   K t8  = matrix[0][1] * matrix[1][0];
01286   K t10 = matrix[0][2] * matrix[1][0];
01287   K t12 = matrix[0][1] * matrix[2][0];
01288   K t14 = matrix[0][2] * matrix[2][0];
01289 
01290   K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
01291            t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
01292   K t17 = 1.0/det;
01293 
01294   inverse[0][0] =  (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
01295   inverse[0][1] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
01296   inverse[0][2] =  (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
01297   inverse[1][0] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
01298   inverse[1][1] =  (matrix[0][0] * matrix[2][2] - t14) * t17;
01299   inverse[1][2] = -(t6-t10) * t17;
01300   inverse[2][0] =  (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
01301   inverse[2][1] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
01302   inverse[2][2] =  (t4-t8) * t17;
01303 
01304   return det;
01305 }
01306 
01308 template <typename K>
01309 static inline K invertMatrix_retTransposed (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
01310 {
01311   // code generated by maple 
01312   K t4  = matrix[0][0] * matrix[1][1];
01313   K t6  = matrix[0][0] * matrix[1][2];
01314   K t8  = matrix[0][1] * matrix[1][0];
01315   K t10 = matrix[0][2] * matrix[1][0];
01316   K t12 = matrix[0][1] * matrix[2][0];
01317   K t14 = matrix[0][2] * matrix[2][0];
01318 
01319   K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
01320            t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
01321   K t17 = 1.0/det;
01322 
01323   inverse[0][0] =  (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
01324   inverse[1][0] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
01325   inverse[2][0] =  (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
01326   inverse[0][1] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
01327   inverse[1][1] =  (matrix[0][0] * matrix[2][2] - t14) * t17;
01328   inverse[2][1] = -(t6-t10) * t17;
01329   inverse[0][2] =  (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
01330   inverse[1][2] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
01331   inverse[2][2] =  (t4-t8) * t17;
01332 
01333   return det;
01334 }
01335 
01337 template< class K, int m, int n, int p >
01338 static inline void multMatrix ( const FieldMatrix< K, m, n > &A,
01339                                 const FieldMatrix< K, n, p > &B,
01340                                 FieldMatrix< K, m, p > &ret )
01341 {
01342   typedef typename FieldMatrix< K, m, p > :: size_type size_type;
01343 
01344   for( size_type i = 0; i < m; ++i )
01345   {
01346     for( size_type j = 0; j < p; ++j )
01347     {
01348       ret[ i ][ j ] = K( 0 );
01349       for( size_type k = 0; k < n; ++k )
01350         ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ];
01351     }
01352   }
01353 }
01354 
01356 template <typename K, int rows,int cols>
01357 static inline void multTransposedMatrix(const FieldMatrix<K,rows,cols> &matrix, FieldMatrix<K,cols,cols>& ret)
01358 {
01359   typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
01360   
01361   for(size_type i=0; i<cols; i++)
01362     for(size_type j=0; j<cols; j++)
01363       { ret[i][j]=0.0;
01364         for(size_type k=0; k<rows; k++)
01365           ret[i][j]+=matrix[k][i]*matrix[k][j];
01366       }
01367 }
01368 
01370 template <typename K, int rows,int cols>
01371 static inline void multAssign(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x, FieldVector<K,rows> & ret) 
01372 {
01373   typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
01374   
01375   for(size_type i=0; i<rows; ++i)
01376   {
01377     ret[i] = 0.0;
01378     for(size_type j=0; j<cols; ++j)
01379     {
01380       ret[i] += matrix[i][j]*x[j];
01381     }
01382   }
01383 }
01384 
01386 template <typename K, int rows, int cols>
01387 static inline void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret) 
01388 {
01389   typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
01390   
01391   for(size_type i=0; i<cols; ++i)
01392   {
01393     ret[i] = 0.0;
01394     for(size_type j=0; j<rows; ++j)
01395       ret[i] += matrix[j][i]*x[j];
01396   }
01397 }
01398 
01400 template <typename K, int rows,int cols>
01401 static inline FieldVector<K,rows> mult(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x) 
01402 {
01403   FieldVector<K,rows> ret;
01404   multAssign(matrix,x,ret);
01405   return ret;
01406 }
01407 
01409 template <typename K, int rows, int cols>
01410 static inline FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x) 
01411 {
01412   FieldVector<K,cols> ret;
01413   multAssignTransposed( matrix, x, ret );
01414   return ret; 
01415 }
01416 
01417 } // end namespace FMatrixHelp 
01418 
01419 #ifdef DUNE_EXPRESSIONTEMPLATES
01420 template <class K, int N, int M>
01421 struct BlockType< FieldMatrix<K,N,M> >
01422 {
01423   typedef K type;
01424 };
01425 
01426 template <class K, int N, int M>
01427 struct FieldType< FieldMatrix<K,N,M> >
01428 {
01429   typedef K type;
01430 };
01431 #endif // DUNE_EXPRESSIONTEMPLATES
01432 
01435 } // end namespace
01436 
01437 #endif

Generated on Sun Nov 15 22:28:13 2009 for dune-common by  doxygen 1.5.6