Dune Core Modules (2.6.0)

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 #include <dune/common/math.hh>
18 
19 namespace Dune {
20 
26  namespace FMatrixHelp {
27 
28  // defined in fmatrixev.cc
29  extern void eigenValuesLapackCall(
30  const char* jobz, const char* uplo, const long
31  int* n, double* a, const long int* lda, double* w,
32  double* work, const long int* lwork, long int* info);
33 
34  extern void eigenValuesNonsymLapackCall(
35  const char* jobvl, const char* jobvr, const long
36  int* n, double* a, const long int* lda, double* wr, double* wi, double* vl,
37  const long int* ldvl, double* vr, const long int* ldvr, double* work,
38  const long int* lwork, const long int* info);
39 
45  template <typename K>
46  static void eigenValues(const FieldMatrix<K, 1, 1>& matrix,
47  FieldVector<K, 1>& eigenvalues)
48  {
49  eigenvalues[0] = matrix[0][0];
50  }
51 
57  template <typename K>
58  static void eigenValues(const FieldMatrix<K, 2, 2>& matrix,
59  FieldVector<K, 2>& eigenvalues)
60  {
61  using std::sqrt;
62  const K detM = matrix[0][0] * matrix[1][1] - matrix[1][0] * matrix[0][1];
63  const K p = 0.5 * (matrix[0][0] + matrix [1][1]);
64  K q = p * p - detM;
65  if( q < 0 && q > -1e-14 ) q = 0;
66  if (q < 0)
67  {
68  std::cout << matrix << std::endl;
69  // Complex eigenvalues are either caused by non-symmetric matrices or by round-off errors
70  DUNE_THROW(MathError, "Complex eigenvalue detected (which this implementation cannot handle).");
71  }
72 
73  // get square root
74  q = sqrt(q);
75 
76  // store eigenvalues in ascending order
77  eigenvalues[0] = p - q;
78  eigenvalues[1] = p + q;
79  }
80 
93  template <typename K>
94  static void eigenValues(const FieldMatrix<K, 3, 3>& matrix,
95  FieldVector<K, 3>& eigenvalues)
96  {
97  using std::sqrt;
98  using std::acos;
99  const K pi = MathematicalConstants<K>::pi();
100  K p1 = matrix[0][1]*matrix[0][1] + matrix[0][2]*matrix[0][2] + matrix[1][2]*matrix[1][2];
101 
102  if (p1 <= 1e-8)
103  {
104  // A is diagonal.
105  eigenvalues[0] = matrix[0][0];
106  eigenvalues[1] = matrix[1][1];
107  eigenvalues[2] = matrix[2][2];
108  }
109  else
110  {
111  // q = trace(A)/3
112  K q = 0;
113  for (int i=0; i<3; i++)
114  q += matrix[i][i]/3.0;
115 
116  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;
117  K p = sqrt(p2 / 6);
118  // B = (1 / p) * (A - q * I); // I is the identity matrix
120  for (int i=0; i<3; i++)
121  for (int j=0; j<3; j++)
122  B[i][j] = (1/p) * (matrix[i][j] - q*(i==j));
123 
124  K r = B.determinant() / 2.0;
125 
126  // In exact arithmetic for a symmetric matrix -1 <= r <= 1
127  // but computation error can leave it slightly outside this range.
128  K phi;
129  if (r <= -1)
130  phi = pi / 3.0;
131  else if (r >= 1)
132  phi = 0;
133  else
134  phi = acos(r) / 3;
135 
136  // the eigenvalues satisfy eig[2] <= eig[1] <= eig[0]
137  eigenvalues[2] = q + 2 * p * cos(phi);
138  eigenvalues[0] = q + 2 * p * cos(phi + (2*pi/3));
139  eigenvalues[1] = 3 * q - eigenvalues[0] - eigenvalues[2]; // since trace(matrix) = eig1 + eig2 + eig3
140  }
141  }
142 
150  template <int dim, typename K>
151  static void eigenValues(const FieldMatrix<K, dim, dim>& matrix,
152  FieldVector<K, dim>& eigenvalues)
153  {
154  {
155  const long int N = dim ;
156  const char jobz = 'n'; // only calculate eigenvalues
157  const char uplo = 'u'; // use upper triangular matrix
158 
159  // length of matrix vector
160  const long int w = N * N ;
161 
162  // matrix to put into dsyev
163  double matrixVector[dim * dim];
164 
165  // copy matrix
166  int row = 0;
167  for(int i=0; i<dim; ++i)
168  {
169  for(int j=0; j<dim; ++j, ++row)
170  {
171  matrixVector[ row ] = matrix[ i ][ j ];
172  }
173  }
174 
175  // working memory
176  double workSpace[dim * dim];
177 
178  // return value information
179  long int info = 0;
180 
181  // call LAPACK routine (see fmatrixev.cc)
182  eigenValuesLapackCall(&jobz, &uplo, &N, &matrixVector[0], &N,
183  &eigenvalues[0], &workSpace[0], &w, &info);
184 
185  if( info != 0 )
186  {
187  std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
188  DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
189  }
190  }
191  }
199  template <int dim, typename K, class C>
200  static void eigenValuesNonSym(const FieldMatrix<K, dim, dim>& matrix,
201  FieldVector<C, dim>& eigenValues)
202  {
203  {
204  const long int N = dim ;
205  const char jobvl = 'n';
206  const char jobvr = 'n';
207 
208  // matrix to put into dgeev
209  double matrixVector[dim * dim];
210 
211  // copy matrix
212  int row = 0;
213  for(int i=0; i<dim; ++i)
214  {
215  for(int j=0; j<dim; ++j, ++row)
216  {
217  matrixVector[ row ] = matrix[ i ][ j ];
218  }
219  }
220 
221  // working memory
222  double eigenR[dim];
223  double eigenI[dim];
224  double work[3*dim];
225 
226  // return value information
227  long int info = 0;
228  long int lwork = 3*dim;
229 
230  // call LAPACK routine (see fmatrixev_ext.cc)
231  eigenValuesNonsymLapackCall(&jobvl, &jobvr, &N, &matrixVector[0], &N,
232  &eigenR[0], &eigenI[0], 0, &N, 0, &N, &work[0],
233  &lwork, &info);
234 
235  if( info != 0 )
236  {
237  std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
238  DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
239  }
240  for (int i=0; i<N; ++i) {
241  eigenValues[i].real = eigenR[i];
242  eigenValues[i].imag = eigenI[i];
243  }
244  }
245 
246  }
247 
248  } // end namespace FMatrixHelp
249 
252 } // end namespace Dune
253 #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:200
static void eigenValues(const FieldMatrix< K, 1, 1 > &matrix, FieldVector< K, 1 > &eigenvalues)
calculates the eigenvalues of a symmetric field matrix
Definition: fmatrixev.hh:46
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
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:10
Provides commonly used mathematical constants.
Definition: math.hh:24
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 27, 22:29, 2024)