Dune Core Modules (2.7.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
17#include <dune/common/math.hh>
18
19namespace 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, 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 using real_type = typename FieldTraits<K>::real_type;
100 const K pi = MathematicalConstants<K>::pi();
101 K p1 = matrix[0][1]*matrix[0][1] + matrix[0][2]*matrix[0][2] + matrix[1][2]*matrix[1][2];
102
103 if (p1 <= 1e-8)
104 {
105 // A is diagonal.
106 eigenvalues[0] = matrix[0][0];
107 eigenvalues[1] = matrix[1][1];
108 eigenvalues[2] = matrix[2][2];
109 }
110 else
111 {
112 // q = trace(A)/3
113 K q = 0;
114 for (int i=0; i<3; i++)
115 q += matrix[i][i]/3.0;
116
117 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;
118 K p = sqrt(p2 / 6);
119 // B = (1 / p) * (A - q * I); // I is the identity matrix
121 for (int i=0; i<3; i++)
122 for (int j=0; j<3; j++)
123 B[i][j] = (real_type(1.0)/p) * (matrix[i][j] - q*(i==j));
124
125 K r = B.determinant() / 2.0;
126
127 // In exact arithmetic for a symmetric matrix -1 <= r <= 1
128 // but computation error can leave it slightly outside this range.
129 K phi;
130 if (r <= -1)
131 phi = pi / 3.0;
132 else if (r >= 1)
133 phi = 0;
134 else
135 phi = acos(r) / 3;
136
137 // the eigenvalues satisfy eig[2] <= eig[1] <= eig[0]
138 eigenvalues[2] = q + 2 * p * cos(phi);
139 eigenvalues[0] = q + 2 * p * cos(phi + (2*pi/3));
140 eigenvalues[1] = 3 * q - eigenvalues[0] - eigenvalues[2]; // since trace(matrix) = eig1 + eig2 + eig3
141 }
142 }
143
153 template <int dim, typename K>
155 FieldVector<K, dim>& eigenValues,
156 FieldMatrix<K, dim, dim>& eigenVectors,
157 const char& jobz)
158 {
159 {
160 const long int N = dim ;
161 const char uplo = 'u'; // use upper triangular matrix
162
163 // length of matrix vector, LWORK >= max(1,3*N-1)
164 const long int lwork = 3*N -1 ;
165
166 // matrix to put into dsyev
167 double matrixVector[dim * dim];
168
169 // copy matrix
170 int row = 0;
171 for(int i=0; i<dim; ++i)
172 {
173 for(int j=0; j<dim; ++j, ++row)
174 {
175 matrixVector[ row ] = matrix[ i ][ j ];
176 }
177 }
178
179 // working memory
180 double workSpace[lwork];
181
182 // return value information
183 long int info = 0;
184
185 // call LAPACK routine (see fmatrixev.cc)
186 eigenValuesLapackCall(&jobz, &uplo, &N, &matrixVector[0], &N,
187 &eigenValues[0], &workSpace[0], &lwork, &info);
188
189 // restore eigenvectors matrix
190 if (jobz=='v'){
191 row = 0;
192 for(int i=0; i<dim; ++i)
193 {
194 for(int j=0; j<dim; ++j, ++row)
195 {
196 eigenVectors[ i ][ j ] = matrixVector[ row ];
197 }
198 }
199 }
200
201 if( info != 0 )
202 {
203 std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
204 DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
205 }
206 }
207 }
208
216 template <int dim, typename K>
217 static void eigenValues(const FieldMatrix<K, dim, dim>& matrix,
218 FieldVector<K, dim>& eigenValues)
219 {
221 // calculate only eigenValues
222 eigenValuesVectorsLapack(matrix,eigenValues,dummy,'n');
223 }
224
233 template <int dim, typename K>
235 FieldVector<K, dim>& eigenValues,
236 FieldMatrix<K, dim, dim>& eigenVectors)
237 {
238 // calculate eigenValues and eigenVectors
239 eigenValuesVectorsLapack(matrix,eigenValues,eigenVectors,'v');
240 }
241
249 template <int dim, typename K, class C>
251 FieldVector<C, dim>& eigenValues)
252 {
253 {
254 const long int N = dim ;
255 const char jobvl = 'n';
256 const char jobvr = 'n';
257
258 // matrix to put into dgeev
259 double matrixVector[dim * dim];
260
261 // copy matrix
262 int row = 0;
263 for(int i=0; i<dim; ++i)
264 {
265 for(int j=0; j<dim; ++j, ++row)
266 {
267 matrixVector[ row ] = matrix[ i ][ j ];
268 }
269 }
270
271 // working memory
272 double eigenR[dim];
273 double eigenI[dim];
274 double work[3*dim];
275
276 // return value information
277 long int info = 0;
278 const long int lwork = 3*dim;
279
280 // call LAPACK routine (see fmatrixev_ext.cc)
281 eigenValuesNonsymLapackCall(&jobvl, &jobvr, &N, &matrixVector[0], &N,
282 &eigenR[0], &eigenI[0], nullptr, &N, nullptr, &N, &work[0],
283 &lwork, &info);
284
285 if( info != 0 )
286 {
287 std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
288 DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
289 }
290 for (int i=0; i<N; ++i) {
291 eigenValues[i].real = eigenR[i];
292 eigenValues[i].imag = eigenI[i];
293 }
294 }
295
296 }
297
298 } // end namespace FMatrixHelp
299
302} // end namespace Dune
303#endif
field_type determinant(bool doPivoting=true) const
calculates the determinant of this matrix
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:96
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 non-symmetric field matrix
Definition: fmatrixev.hh:250
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
static void eigenValuesVectors(const FieldMatrix< K, dim, dim > &matrix, FieldVector< K, dim > &eigenValues, FieldMatrix< K, dim, dim > &eigenVectors)
calculates the eigenvalues and eigenvectors of a symmetric field matrix
Definition: fmatrixev.hh:234
static void eigenValuesVectorsLapack(const FieldMatrix< K, dim, dim > &matrix, FieldVector< K, dim > &eigenValues, FieldMatrix< K, dim, dim > &eigenVectors, const char &jobz)
calculates the eigenvalues and/or eigenvectors of a symmetric field matrix
Definition: fmatrixev.hh:154
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:14
static const Field pi()
Archimedes' constant.
Definition: math.hh:46
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)