fmatrixev.hh

Go to the documentation of this file.
00001 #ifndef DUNE_FMATRIXEIGENVALUES_HH 
00002 #define DUNE_FMATRIXEIGENVALUES_HH 
00003 
00004 #include <iostream>
00005 #include <cmath>
00006 #include <cassert>
00007 
00008 #include <dune/common/exceptions.hh>
00009 #include <dune/common/fvector.hh>
00010 #include <dune/common/fmatrix.hh>
00011 
00012 #if HAVE_LAPACK 
00013 
00014 // symetric matrices
00015 #define DSYEV_FORTRAN FC_FUNC (dsyev, DSYEV)
00016 
00017 // nonsymetric matrices
00018 #define DGEEV_FORTRAN FC_FUNC (dgeev, DGEEV)
00019 
00020 // dsyev declaration (in liblapack)
00021 extern "C" {
00022 
00023 /*
00024     *
00025     **  purpose
00026     **  =======
00027     **
00028     **  xsyev computes all eigenvalues and, optionally, eigenvectors of a
00029     **  BASE DATA TYPE symmetric matrix a.
00030     **
00031     **  arguments
00032     **  =========
00033     **
00034     **  jobz    (input) char
00035     **          = 'n':  compute eigenvalues only;
00036     **          = 'v':  compute eigenvalues and eigenvectors.
00037     **
00038     **  uplo    (input) char
00039     **          = 'u':  upper triangle of a is stored;
00040     **          = 'l':  lower triangle of a is stored.
00041     **
00042     **  n       (input) long int
00043     **          the order of the matrix a.  n >= 0.
00044     **
00045     **  a       (input/output) BASE DATA TYPE array, dimension (lda, n)
00046     **          on entry, the symmetric matrix a.  if uplo = 'u', the
00047     **          leading n-by-n upper triangular part of a contains the
00048     **          upper triangular part of the matrix a.  if uplo = 'l',
00049     **          the leading n-by-n lower triangular part of a contains
00050     **          the lower triangular part of the matrix a.
00051     **          on exit, if jobz = 'v', then if info = 0, a contains the
00052     **          orthonormal eigenvectors of the matrix a.
00053     **          if jobz = 'n', then on exit the lower triangle (if uplo='l')
00054     **          or the upper triangle (if uplo='u') of a, including the
00055     **          diagonal, is destroyed.
00056     **
00057     **  lda     (input) long int
00058     **          the leading dimension of the array a.  lda >= max(1,n).
00059     **
00060     **  w       (output) BASE DATA TYPE array, dimension (n)
00061     **          if info = 0, the eigenvalues in ascending order.
00062     **
00063     **
00064     **
00065     **  info    (output) long int
00066     **          = 0:  successful exit
00067     **          < 0:  if info = -i, the i-th argument had an illegal value
00068     **          > 0:  if info = i, the algorithm failed to converge; i
00069     **                off-diagonal elements of an intermediate tridiagonal
00070     **                form did not converge to zero.
00071     **
00072 **/
00073 extern void DSYEV_FORTRAN(const char* jobz, const char* uplo, const long
00074     int* n, double* a, const long int* lda, double* w,
00075     double* work, const long int* lwork, long int* info);
00076     
00077 } // end extern C 
00078 #endif
00079 
00080 namespace Dune {
00081 
00087 namespace FMatrixHelp {
00088 
00094 template <typename K> 
00095 static void eigenValues(const FieldMatrix<K, 1, 1>& matrix,
00096                         FieldVector<K, 1>& eigenvalues)  
00097 {
00098   eigenvalues[0] = matrix[0][0];
00099 }
00100 
00106 template <typename K>
00107 static void eigenValues(const FieldMatrix<K, 2, 2>& matrix,
00108                         FieldVector<K, 2>& eigenvalues)  
00109 {
00110   const K detM = matrix[0][0] * matrix[1][1] - matrix[1][0] * matrix[0][1];
00111   const K p = 0.5 * (matrix[0][0] + matrix [1][1]);
00112   K q = p * p - detM;
00113   if( q < 0 && q > -1e-14 ) q = 0;
00114   if (p < 0 || q < 0)
00115   {
00116     std::cout << p << " p | q " << q << "\n";
00117     std::cout << matrix << std::endl;
00118     std::cout << "something went wrong in Eigenvalues for matrix!" << std::endl;
00119     assert(false);
00120     abort();
00121   }
00122 
00123   // get square root 
00124   q = std :: sqrt(q);
00125 
00126   // store eigenvalues in ascending order 
00127   eigenvalues[0] = p - q;
00128   eigenvalues[1] = p + q;
00129 }
00130 
00131 #if HAVE_LAPACK || defined DOXYGEN
00132 
00139 template <int dim, typename K> 
00140 static void eigenValues(const FieldMatrix<K, dim, dim>& matrix,
00141                         FieldVector<K, dim>& eigenvalues)  
00142 {
00143   {
00144     const long int N = dim ;
00145     const char jobz = 'n'; // only calculate eigenvalues  
00146     const char uplo = 'u'; // use upper triangular matrix 
00147 
00148     // length of matrix vector 
00149     const long int w = N * N ;
00150 
00151     // matrix to put into dsyev 
00152     double matrixVector[dim * dim]; 
00153 
00154     // copy matrix  
00155     int row = 0;
00156     for(int i=0; i<dim; ++i) 
00157     {
00158       for(int j=0; j<dim; ++j, ++row) 
00159       {
00160         matrixVector[ row ] = matrix[ i ][ j ];
00161       }
00162     }
00163 
00164     // working memory 
00165     double workSpace[dim * dim]; 
00166 
00167     // return value information 
00168     long int info = 0;
00169 
00170     // call LAPACK dsyev 
00171     DSYEV_FORTRAN(&jobz, &uplo, &N, &matrixVector[0], &N, 
00172         &eigenvalues[0], &workSpace[0], &w, &info);
00173 
00174     if( info != 0 ) 
00175     {
00176       std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
00177       DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
00178     }
00179   }
00180 }
00181 #endif
00182 
00183 } // end namespace FMatrixHelp 
00184 
00187 } // end namespace Dune 
00188 #endif

Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].