5#ifndef DUNE_FMATRIXEIGENVALUES_HH
6#define DUNE_FMATRIXEIGENVALUES_HH
29 namespace FMatrixHelp {
33 extern void eigenValuesLapackCall(
34 const char* jobz,
const char* uplo,
const long
35 int* n,
double* a,
const long int* lda,
double* w,
36 double* work,
const long int* lwork,
long int* info);
38 extern void eigenValuesNonsymLapackCall(
39 const char* jobvl,
const char* jobvr,
const long
40 int* n,
double* a,
const long int* lda,
double* wr,
double* wi,
double* vl,
41 const long int* ldvl,
double* vr,
const long int* ldvr,
double* work,
42 const long int* lwork,
long int* info);
44 extern void eigenValuesLapackCall(
45 const char* jobz,
const char* uplo,
const long
46 int* n,
float* a,
const long int* lda,
float* w,
47 float* work,
const long int* lwork,
long int* info);
49 extern void eigenValuesNonsymLapackCall(
50 const char* jobvl,
const char* jobvr,
const long
51 int* n,
float* a,
const long int* lda,
float* wr,
float* wi,
float* vl,
52 const long int* ldvl,
float* vr,
const long int* ldvr,
float* work,
53 const long int* lwork,
long int* info);
59 enum Jobs { OnlyEigenvalues=0, EigenvaluesEigenvectors=1 };
62 template<
typename K,
int dim>
63 using EVDummy = FieldMatrix<K, dim, dim>;
67 inline FieldVector<K,3> crossProduct(
const FieldVector<K,3>& vec0,
const FieldVector<K,3>& vec1) {
68 return {vec0[1]*vec1[2] - vec0[2]*vec1[1], vec0[2]*vec1[0] - vec0[0]*vec1[2], vec0[0]*vec1[1] - vec0[1]*vec1[0]};
72 static void eigenValues2dImpl(
const FieldMatrix<K, 2, 2>& matrix,
73 FieldVector<K, 2>& eigenvalues)
76 const K p = 0.5 * (matrix[0][0] + matrix [1][1]);
77 const K p2 = p - matrix[1][1];
78 K q = p2 * p2 + matrix[1][0] * matrix[0][1];
79 if( q < 0 && q > -1e-14 ) q = 0;
82 std::cout << matrix << std::endl;
84 DUNE_THROW(MathError,
"Complex eigenvalue detected (which this implementation cannot handle).");
91 eigenvalues[0] = p - q;
92 eigenvalues[1] = p + q;
102 template <
typename K>
103 static K eigenValues3dImpl(
const FieldMatrix<K, 3, 3>& matrix,
104 FieldVector<K, 3>& eigenvalues)
108 using real_type =
typename FieldTraits<K>::real_type;
110 K p1 = matrix[0][1]*matrix[0][1] + matrix[0][2]*matrix[0][2] + matrix[1][2]*matrix[1][2];
112 if (p1 <= std::numeric_limits<K>::epsilon()) {
114 eigenvalues[0] = matrix[0][0];
115 eigenvalues[1] = matrix[1][1];
116 eigenvalues[2] = matrix[2][2];
117 std::sort(eigenvalues.begin(), eigenvalues.end());
125 for (
int i=0; i<3; i++)
126 q += matrix[i][i] / 3.0;
128 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.0 * p1;
131 FieldMatrix<K,3,3> B;
132 for (
int i=0; i<3; i++)
133 for (
int j=0; j<3; j++)
134 B[i][j] = (real_type(1.0)/p) * (matrix[i][j] - q*(i==j));
136 K r = B.determinant() / 2.0;
144 r = clamp<K>(r, -1.0, 1.0);
145 K phi = acos(r) / 3.0;
148 eigenvalues[2] = q + 2 * p * cos(phi);
149 eigenvalues[0] = q + 2 * p * cos(phi + (2*pi/3));
150 eigenvalues[1] = 3 * q - eigenvalues[0] - eigenvalues[2];
159 void orthoComp(
const FieldVector<K,3>& evec0, FieldVector<K,3>& u, FieldVector<K,3>& v) {
161 if(abs(evec0[0]) > abs(evec0[1])) {
163 FieldVector<K,2> temp = {evec0[0], evec0[2]};
164 auto L = 1.0 / temp.two_norm();
165 u = L * FieldVector<K,3>({-evec0[2], 0.0, evec0[0]});
169 FieldVector<K,2> temp = {evec0[1], evec0[2]};
170 auto L = 1.0 / temp.two_norm();
171 u = L * FieldVector<K,3>({0.0, evec0[2], -evec0[1]});
173 v = crossProduct(evec0, u);
178 void eig0(
const FieldMatrix<K,3,3>& matrix, K eval0, FieldVector<K,3>& evec0) {
184 using Vector = FieldVector<K,3>;
185 Vector row0 = {matrix[0][0]-eval0, matrix[0][1], matrix[0][2]};
186 Vector row1 = {matrix[1][0], matrix[1][1]-eval0, matrix[1][2]};
187 Vector row2 = {matrix[2][0], matrix[2][1], matrix[2][2]-eval0};
189 Vector r0xr1 = crossProduct(row0, row1);
190 Vector r0xr2 = crossProduct(row0, row2);
191 Vector r1xr2 = crossProduct(row1, row2);
192 auto d0 = r0xr1.two_norm();
193 auto d1 = r0xr2.two_norm();
194 auto d2 = r1xr2.two_norm();
215 void eig1(
const FieldMatrix<K,3,3>& matrix,
const FieldVector<K,3>& evec0, FieldVector<K,3>& evec1, K eval1) {
216 using Vector = FieldVector<K,3>;
220 orthoComp(evec0, u, v);
244 auto m00 = u.dot(Au) - eval1;
245 auto m01 = u.dot(Av);
246 auto m11 = v.dot(Av) - eval1;
253 using std::abs, std::sqrt,
std::max;
254 auto absM00 = abs(m00);
255 auto absM01 = abs(m01);
256 auto absM11 = abs(m11);
257 if(absM00 >= absM11) {
258 auto maxAbsComp =
max(absM00, absM01);
259 if(maxAbsComp > 0.0) {
260 if(absM00 >= absM01) {
262 m00 = 1.0 / sqrt(1.0 + m01*m01);
267 m01 = 1.0 / sqrt(1.0 + m00*m00);
270 evec1 = m01*u - m00*v;
276 auto maxAbsComp =
max(absM11, absM01);
277 if(maxAbsComp > 0.0) {
278 if(absM11 >= absM01) {
280 m11 = 1.0 / sqrt(1.0 + m01*m01);
285 m01 = 1.0 / sqrt(1.0 + m11*m11);
288 evec1 = m11*u - m01*v;
296 template<Jobs Tag,
typename K>
297 static void eigenValuesVectorsImpl(
const FieldMatrix<K, 1, 1>& matrix,
298 FieldVector<K, 1>& eigenValues,
299 FieldMatrix<K, 1, 1>& eigenVectors)
302 if constexpr(Tag==EigenvaluesEigenvectors)
303 eigenVectors[0] = {1.0};
308 template <Jobs Tag,
typename K>
309 static void eigenValuesVectorsImpl(
const FieldMatrix<K, 2, 2>& matrix,
310 FieldVector<K, 2>& eigenValues,
311 FieldMatrix<K, 2, 2>& eigenVectors)
314 Impl::eigenValues2dImpl(matrix, eigenValues);
323 if constexpr(Tag==EigenvaluesEigenvectors) {
326 FieldMatrix<K,2,2> temp = matrix;
329 if(temp.infinity_norm() <= 1e-14) {
330 eigenVectors[0] = {1.0, 0.0};
331 eigenVectors[1] = {0.0, 1.0};
336 FieldVector<K,2> ev0 = {matrix[0][0]-
eigenValues[1], matrix[1][0]};
337 FieldVector<K,2> ev1 = {matrix[0][1], matrix[1][1]-
eigenValues[1]};
338 eigenVectors[0] = (ev0.two_norm2() >= ev1.two_norm2()) ? ev0/ev0.two_norm() : ev1/ev1.two_norm();
344 eigenVectors[1] = (ev0.two_norm2() >= ev1.two_norm2()) ? ev0/ev0.two_norm() : ev1/ev1.two_norm();
350 template <Jobs Tag,
typename K>
351 static void eigenValuesVectorsImpl(
const FieldMatrix<K, 3, 3>& matrix,
352 FieldVector<K, 3>& eigenValues,
353 FieldMatrix<K, 3, 3>& eigenVectors)
355 using Vector = FieldVector<K,3>;
356 using Matrix = FieldMatrix<K,3,3>;
363 K maxAbsElement = (isnormal(matrix.infinity_norm())) ? matrix.infinity_norm() : K(1.0);
364 Matrix scaledMatrix = matrix / maxAbsElement;
365 K r = Impl::eigenValues3dImpl(scaledMatrix, eigenValues);
367 if constexpr(Tag==EigenvaluesEigenvectors) {
368 K offDiagNorm = Vector{scaledMatrix[0][1],scaledMatrix[0][2],scaledMatrix[1][2]}.two_norm2();
369 if (offDiagNorm <= std::numeric_limits<K>::epsilon())
371 eigenValues = {scaledMatrix[0][0], scaledMatrix[1][1], scaledMatrix[2][2]};
372 eigenVectors = {{1.0, 0.0, 0.0}, {0.0, 1.0, 0.0}, {0.0, 0.0, 1.0}};
376 if (eigenValues[0] > eigenValues[1])
378 std::swap(eigenValues[0], eigenValues[1]);
379 std::swap(eigenVectors[0], eigenVectors[1]);
381 if (eigenValues[1] > eigenValues[2])
383 std::swap(eigenValues[1], eigenValues[2]);
384 std::swap(eigenVectors[1], eigenVectors[2]);
386 if (eigenValues[0] > eigenValues[1])
388 std::swap(eigenValues[0], eigenValues[1]);
389 std::swap(eigenVectors[0], eigenVectors[1]);
398 Vector eval(eigenValues);
400 Impl::eig0(scaledMatrix, eval[2], evec[2]);
401 Impl::eig1(scaledMatrix, evec[2], evec[1], eval[1]);
402 evec[0] = Impl::crossProduct(evec[1], evec[2]);
405 Impl::eig0(scaledMatrix, eval[0], evec[0]);
406 Impl::eig1(scaledMatrix, evec[0], evec[1], eval[1]);
407 evec[2] = Impl::crossProduct(evec[0], evec[1]);
410 using EVPair = std::pair<K, Vector>;
411 std::vector<EVPair> pairs;
412 for(std::size_t i=0; i<=2; ++i)
413 pairs.push_back(EVPair(eval[i], evec[i]));
414 auto comp = [](EVPair x, EVPair y){
return x.first < y.first; };
415 std::sort(pairs.begin(), pairs.end(), comp);
416 for(std::size_t i=0; i<=2; ++i){
418 eigenVectors[i] = pairs[i].second;
427 template <Jobs Tag,
int dim,
typename K>
428 static void eigenValuesVectorsLapackImpl(
const FieldMatrix<K, dim, dim>& matrix,
429 FieldVector<K, dim>& eigenValues,
430 FieldMatrix<K, dim, dim>& eigenVectors)
436 const char jobz =
"nv"[Tag];
438 const long int N = dim ;
439 const char uplo =
'u';
442 const long int lwork = 3*N -1 ;
444 constexpr bool isKLapackType = std::is_same_v<K,double> || std::is_same_v<K,float>;
445 using LapackNumType = std::conditional_t<isKLapackType, K, double>;
448 LapackNumType matrixVector[dim * dim];
452 for(
int i=0; i<dim; ++i)
454 for(
int j=0; j<dim; ++j, ++row)
456 matrixVector[ row ] = matrix[ i ][ j ];
461 LapackNumType workSpace[lwork];
466 if constexpr (isKLapackType){
469 ev =
new LapackNumType[dim];
473 eigenValuesLapackCall(&jobz, &uplo, &N, &matrixVector[0], &N,
474 ev, &workSpace[0], &lwork, &info);
476 if constexpr (!isKLapackType){
477 for(
size_t i=0;i<dim;++i)
478 eigenValues[i] = ev[i];
483 if (Tag==Jobs::EigenvaluesEigenvectors){
485 for(
int i=0; i<dim; ++i)
487 for(
int j=0; j<dim; ++j, ++row)
489 eigenVectors[ i ][ j ] = matrixVector[ row ];
496 std::cerr <<
"For matrix " << matrix <<
" eigenvalue calculation failed! " << std::endl;
497 DUNE_THROW(InvalidStateException,
"eigenValues: Eigenvalue calculation failed!");
500 DUNE_THROW(NotImplemented,
"LAPACK not found!");
506 template <Jobs Tag,
int dim,
typename K>
507 static void eigenValuesVectorsImpl(
const FieldMatrix<K, dim, dim>& matrix,
508 FieldVector<K, dim>& eigenValues,
509 FieldMatrix<K, dim, dim>& eigenVectors)
511 eigenValuesVectorsLapackImpl<Tag>(matrix,eigenValues,eigenVectors);
522 template <
int dim,
typename K>
527 Impl::eigenValuesVectorsImpl<Impl::Jobs::OnlyEigenvalues>(matrix,
eigenValues, dummy);
538 template <
int dim,
typename K>
543 Impl::eigenValuesVectorsImpl<Impl::Jobs::EigenvaluesEigenvectors>(matrix,
eigenValues, eigenVectors);
553 template <
int dim,
typename K>
558 Impl::eigenValuesVectorsLapackImpl<Impl::Jobs::EigenvaluesEigenvectors>(matrix,
eigenValues, dummy);
569 template <
int dim,
typename K>
574 Impl::eigenValuesVectorsLapackImpl<Impl::Jobs::EigenvaluesEigenvectors>(matrix,
eigenValues, eigenVectors);
584 template <
int dim,
typename K,
class C>
590 const long int N = dim ;
591 const char jobvl =
'n';
592 const char jobvr =
'n';
594 constexpr bool isKLapackType = std::is_same_v<K,double> || std::is_same_v<K,float>;
595 using LapackNumType = std::conditional_t<isKLapackType, K, double>;
598 LapackNumType matrixVector[dim * dim];
602 for(
int i=0; i<dim; ++i)
604 for(
int j=0; j<dim; ++j, ++row)
606 matrixVector[ row ] = matrix[ i ][ j ];
611 LapackNumType eigenR[dim];
612 LapackNumType eigenI[dim];
613 LapackNumType work[3*dim];
617 const long int lwork = 3*dim;
620 eigenValuesNonsymLapackCall(&jobvl, &jobvr, &N, &matrixVector[0], &N,
621 &eigenR[0], &eigenI[0],
nullptr, &N,
nullptr, &N, &work[0],
626 std::cerr <<
"For matrix " << matrix <<
" eigenvalue calculation failed! " << std::endl;
629 for (
int i=0; i<N; ++i) {
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:281
Default exception for dummy implementations.
Definition: exceptions.hh:263
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:585
static void eigenValues(const FieldMatrix< K, dim, dim > &matrix, FieldVector< K, dim > &eigenValues)
calculates the eigenvalues of a symmetric field matrix
Definition: fmatrixev.hh:523
static void eigenValuesLapack(const FieldMatrix< K, dim, dim > &matrix, FieldVector< K, dim > &eigenValues)
calculates the eigenvalues of a symmetric field matrix
Definition: fmatrixev.hh:554
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:539
static void eigenValuesVectorsLapack(const FieldMatrix< K, dim, dim > &matrix, FieldVector< K, dim > &eigenValues, FieldMatrix< K, dim, dim > &eigenVectors)
calculates the eigenvalues and -vectors of a symmetric field matrix
Definition: fmatrixev.hh:570
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:218
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
static const Field pi()
Archimedes' constant.
Definition: math.hh:48