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 extern "C" {
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
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 }
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
00116 q = std :: sqrt(q);
00117
00118
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';
00138 const char uplo = 'u';
00139
00140
00141 const long int w = N * N ;
00142
00143
00144 double matrixVector[dim * dim];
00145
00146
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
00157 double workSpace[dim * dim];
00158
00159
00160 long int info = 0;
00161
00162
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 }
00178
00181 }
00182 #endif