Dune Core Modules (2.5.1)

fmatrixev.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_FMATRIXEIGENVALUES_HH
4 #define DUNE_FMATRIXEIGENVALUES_HH
5 
10 #include <iostream>
11 #include <cmath>
12 #include <cassert>
13 
15 #include <dune/common/fvector.hh>
16 #include <dune/common/fmatrix.hh>
17 
18 namespace Dune {
19 
25  namespace FMatrixHelp {
26 
27  // defined in fmatrixev.cc
28  extern void eigenValuesLapackCall(
29  const char* jobz, const char* uplo, const long
30  int* n, double* a, const long int* lda, double* w,
31  double* work, const long int* lwork, long int* info);
32 
33  extern void eigenValuesNonsymLapackCall(
34  const char* jobvl, const char* jobvr, const long
35  int* n, double* a, const long int* lda, double* wr, double* wi, double* vl,
36  const long int* ldvl, double* vr, const long int* ldvr, double* work,
37  const long int* lwork, const long int* info);
38 
44  template <typename K>
45  static void eigenValues(const FieldMatrix<K, 1, 1>& matrix,
46  FieldVector<K, 1>& eigenvalues)
47  {
48  eigenvalues[0] = matrix[0][0];
49  }
50 
56  template <typename K>
57  static void eigenValues(const FieldMatrix<K, 2, 2>& matrix,
58  FieldVector<K, 2>& eigenvalues)
59  {
60  const K detM = matrix[0][0] * matrix[1][1] - matrix[1][0] * matrix[0][1];
61  const K p = 0.5 * (matrix[0][0] + matrix [1][1]);
62  K q = p * p - detM;
63  if( q < 0 && q > -1e-14 ) q = 0;
64  if (q < 0)
65  {
66  std::cout << matrix << std::endl;
67  // Complex eigenvalues are either caused by non-symmetric matrices or by round-off errors
68  DUNE_THROW(MathError, "Complex eigenvalue detected (which this implementation cannot handle).");
69  }
70 
71  // get square root
72  q = std :: sqrt(q);
73 
74  // store eigenvalues in ascending order
75  eigenvalues[0] = p - q;
76  eigenvalues[1] = p + q;
77  }
78 
91  template <typename K>
92  static void eigenValues(const FieldMatrix<K, 3, 3>& matrix,
93  FieldVector<K, 3>& eigenvalues)
94  {
95  K p1 = matrix[0][1]*matrix[0][1] + matrix[0][2]*matrix[0][2] + matrix[1][2]*matrix[1][2];
96 
97  if (p1 <= 1e-8)
98  {
99  // A is diagonal.
100  eigenvalues[0] = matrix[0][0];
101  eigenvalues[1] = matrix[1][1];
102  eigenvalues[2] = matrix[2][2];
103  }
104  else
105  {
106  // q = trace(A)/3
107  K q = 0;
108  for (int i=0; i<3; i++)
109  q += matrix[i][i]/3.0;
110 
111  K p2 = (matrix[0][0] - q)*(matrix[0][0] - q) + (matrix[1][1] - q)*(matrix[1][1] - q) + (matrix[2][2] - q)*(matrix[2][2] - q) + 2 * p1;
112  K p = std::sqrt(p2 / 6);
113  // B = (1 / p) * (A - q * I); // I is the identity matrix
115  for (int i=0; i<3; i++)
116  for (int j=0; j<3; j++)
117  B[i][j] = (1/p) * (matrix[i][j] - q*(i==j));
118 
119  K r = B.determinant() / 2.0;
120 
121  // In exact arithmetic for a symmetric matrix -1 <= r <= 1
122  // but computation error can leave it slightly outside this range.
123  K phi;
124  if (r <= -1)
125  phi = M_PI / 3.0;
126  else if (r >= 1)
127  phi = 0;
128  else
129  phi = std::acos(r) / 3;
130 
131  // the eigenvalues satisfy eig[2] <= eig[1] <= eig[0]
132  eigenvalues[2] = q + 2 * p * cos(phi);
133  eigenvalues[0] = q + 2 * p * cos(phi + (2*M_PI/3));
134  eigenvalues[1] = 3 * q - eigenvalues[0] - eigenvalues[2]; // since trace(matrix) = eig1 + eig2 + eig3
135  }
136  }
137 
145  template <int dim, typename K>
146  static void eigenValues(const FieldMatrix<K, dim, dim>& matrix,
147  FieldVector<K, dim>& eigenvalues)
148  {
149  {
150  const long int N = dim ;
151  const char jobz = 'n'; // only calculate eigenvalues
152  const char uplo = 'u'; // use upper triangular matrix
153 
154  // length of matrix vector
155  const long int w = N * N ;
156 
157  // matrix to put into dsyev
158  double matrixVector[dim * dim];
159 
160  // copy matrix
161  int row = 0;
162  for(int i=0; i<dim; ++i)
163  {
164  for(int j=0; j<dim; ++j, ++row)
165  {
166  matrixVector[ row ] = matrix[ i ][ j ];
167  }
168  }
169 
170  // working memory
171  double workSpace[dim * dim];
172 
173  // return value information
174  long int info = 0;
175 
176  // call LAPACK routine (see fmatrixev.cc)
177  eigenValuesLapackCall(&jobz, &uplo, &N, &matrixVector[0], &N,
178  &eigenvalues[0], &workSpace[0], &w, &info);
179 
180  if( info != 0 )
181  {
182  std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
183  DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
184  }
185  }
186  }
194  template <int dim, typename K, class C>
195  static void eigenValuesNonSym(const FieldMatrix<K, dim, dim>& matrix,
196  FieldVector<C, dim>& eigenValues)
197  {
198  {
199  const long int N = dim ;
200  const char jobvl = 'n';
201  const char jobvr = 'n';
202 
203  // matrix to put into dgeev
204  double matrixVector[dim * dim];
205 
206  // copy matrix
207  int row = 0;
208  for(int i=0; i<dim; ++i)
209  {
210  for(int j=0; j<dim; ++j, ++row)
211  {
212  matrixVector[ row ] = matrix[ i ][ j ];
213  }
214  }
215 
216  // working memory
217  double eigenR[dim];
218  double eigenI[dim];
219  double work[3*dim];
220 
221  // return value information
222  long int info = 0;
223  long int lwork = 3*dim;
224 
225  // call LAPACK routine (see fmatrixev_ext.cc)
226  eigenValuesNonsymLapackCall(&jobvl, &jobvr, &N, &matrixVector[0], &N,
227  &eigenR[0], &eigenI[0], 0, &N, 0, &N, &work[0],
228  &lwork, &info);
229 
230  if( info != 0 )
231  {
232  std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
233  DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
234  }
235  for (int i=0; i<N; ++i) {
236  eigenValues[i].real = eigenR[i];
237  eigenValues[i].imag = eigenI[i];
238  }
239  }
240 
241  }
242 
243  } // end namespace FMatrixHelp
244 
247 } // end namespace Dune
248 #endif
field_type determinant() const
calculates the determinant of this matrix
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:279
Default exception class for mathematical errors.
Definition: exceptions.hh:239
A few common exception classes.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
static void eigenValuesNonSym(const FieldMatrix< K, dim, dim > &matrix, FieldVector< C, dim > &eigenValues)
calculates the eigenvalues of a symmetric field matrix
Definition: fmatrixev.hh:195
static void eigenValues(const FieldMatrix< K, 1, 1 > &matrix, FieldVector< K, 1 > &eigenvalues)
calculates the eigenvalues of a symmetric field matrix
Definition: fmatrixev.hh:45
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Dune namespace.
Definition: alignment.hh:11
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 15, 22:30, 2024)