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 // dsyev declaration (in liblapack)
00014 extern "C" {
00015 
00016 /*
00017     *
00018     **  purpose
00019     **  =======
00020     **
00021     **  xsyev computes all eigenvalues and, optionally, eigenvectors of a
00022     **  BASE DATA TYPE symmetric matrix a.
00023     **
00024     **  arguments
00025     **  =========
00026     **
00027     **  jobz    (input) char
00028     **          = 'n':  compute eigenvalues only;
00029     **          = 'v':  compute eigenvalues and eigenvectors.
00030     **
00031     **  uplo    (input) char
00032     **          = 'u':  upper triangle of a is stored;
00033     **          = 'l':  lower triangle of a is stored.
00034     **
00035     **  n       (input) long int
00036     **          the order of the matrix a.  n >= 0.
00037     **
00038     **  a       (input/output) BASE DATA TYPE array, dimension (lda, n)
00039     **          on entry, the symmetric matrix a.  if uplo = 'u', the
00040     **          leading n-by-n upper triangular part of a contains the
00041     **          upper triangular part of the matrix a.  if uplo = 'l',
00042     **          the leading n-by-n lower triangular part of a contains
00043     **          the lower triangular part of the matrix a.
00044     **          on exit, if jobz = 'v', then if info = 0, a contains the
00045     **          orthonormal eigenvectors of the matrix a.
00046     **          if jobz = 'n', then on exit the lower triangle (if uplo='l')
00047     **          or the upper triangle (if uplo='u') of a, including the
00048     **          diagonal, is destroyed.
00049     **
00050     **  lda     (input) long int
00051     **          the leading dimension of the array a.  lda >= max(1,n).
00052     **
00053     **  w       (output) BASE DATA TYPE array, dimension (n)
00054     **          if info = 0, the eigenvalues in ascending order.
00055     **
00056     **
00057     **
00058     **  info    (output) long int
00059     **          = 0:  successful exit
00060     **          < 0:  if info = -i, the i-th argument had an illegal value
00061     **          > 0:  if info = i, the algorithm failed to converge; i
00062     **                off-diagonal elements of an intermediate tridiagonal
00063     **                form did not converge to zero.
00064     **
00065 **/
00066 extern void dsyev_(const char* jobz, const char* uplo, const long
00067                    int* n, double* a, const long int* lda, double* w,
00068                    double* work, const long int* lwork, long int* info);
00069 } // end extern C 
00070 #endif
00071 
00072 namespace Dune {
00073 
00079 namespace FMatrixHelp {
00080 
00086 template <typename K> 
00087 static void eigenValues(const FieldMatrix<K, 1, 1>& matrix,
00088                         FieldVector<K, 1>& eigenvalues)  
00089 {
00090   eigenvalues[0] = matrix[0][0];
00091 }
00092 
00098 template <typename K> 
00099 static void eigenValues(const FieldMatrix<K, 2, 2>& matrix,
00100                         FieldVector<K, 2>& eigenvalues)  
00101 {
00102   const K detM = matrix[0][0] * matrix[1][1] - matrix[1][0] * matrix[0][1];
00103   const K p = 0.5 * (matrix[0][0] + matrix [1][1]);
00104   K q = p * p - detM;
00105   if( q < 0 && q > -1e-14 ) q = 0;
00106   if (p < 0 || q < 0)
00107   {
00108     std::cout << p << " p | q " << q << "\n";
00109     std::cout << matrix << std::endl;
00110     std::cout << "something went wrong in Eigenvalues for matrix!" << std::endl;
00111     assert(false);
00112     abort();
00113   }
00114 
00115   // get square root 
00116   q = std :: sqrt(q);
00117 
00118   // store eigenvalues in ascending order 
00119   eigenvalues[0] = p - q;
00120   eigenvalues[1] = p + q;
00121 }
00122 
00130 template <int dim, typename K> 
00131 static void eigenValues(const FieldMatrix<K, dim, dim>& matrix,
00132                         FieldVector<K, dim>& eigenvalues)  
00133 {
00134 #if HAVE_LAPACK 
00135   {
00136     const long int N = dim ;
00137     const char jobz = 'n'; // only calculate eigenvalues  
00138     const char uplo = 'u'; // use upper triangular matrix 
00139 
00140     // length of matrix vector 
00141     const long int w = N * N ;
00142 
00143     // matrix to put into dsyev 
00144     double matrixVector[dim * dim]; 
00145 
00146     // copy matrix  
00147     int row = 0;
00148     for(int i=0; i<dim; ++i) 
00149     {
00150       for(int j=0; j<dim; ++j, ++row) 
00151       {
00152         matrixVector[ row ] = matrix[ i ][ j ];
00153       }
00154     }
00155 
00156     // working memory 
00157     double workSpace[dim * dim]; 
00158 
00159     // return value information 
00160     long int info = 0;
00161 
00162     // call LAPACK dsyev 
00163     dsyev_(&jobz, &uplo, &N, &matrixVector[0], &N, 
00164            &eigenvalues[0], &workSpace[0], &w, &info);
00165 
00166     if( info != 0 ) 
00167     {
00168       std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
00169       DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
00170     }
00171   }
00172 #else 
00173   DUNE_THROW(NotImplemented,"LAPACK is not available, therefore no eigenvalue calculation");
00174 #endif
00175 }
00176 
00177 } // end namespace FMatrixHelp 
00178 
00181 } // end namespace Dune 
00182 #endif
Generated on Mon Apr 26 10:45:21 2010 for dune-common by  doxygen 1.6.3