Dune Core Modules (2.3.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 (p < 0 || q < 0)
65  {
66  std::cout << p << " p | q " << q << "\n";
67  std::cout << matrix << std::endl;
68  std::cout << "something went wrong in Eigenvalues for matrix!" << std::endl;
69  assert(false);
70  abort();
71  }
72 
73  // get square root
74  q = std :: sqrt(q);
75 
76  // store eigenvalues in ascending order
77  eigenvalues[0] = p - q;
78  eigenvalues[1] = p + q;
79  }
80 
88  template <int dim, typename K>
89  static void eigenValues(const FieldMatrix<K, dim, dim>& matrix,
90  FieldVector<K, dim>& eigenvalues)
91  {
92  {
93  const long int N = dim ;
94  const char jobz = 'n'; // only calculate eigenvalues
95  const char uplo = 'u'; // use upper triangular matrix
96 
97  // length of matrix vector
98  const long int w = N * N ;
99 
100  // matrix to put into dsyev
101  double matrixVector[dim * dim];
102 
103  // copy matrix
104  int row = 0;
105  for(int i=0; i<dim; ++i)
106  {
107  for(int j=0; j<dim; ++j, ++row)
108  {
109  matrixVector[ row ] = matrix[ i ][ j ];
110  }
111  }
112 
113  // working memory
114  double workSpace[dim * dim];
115 
116  // return value information
117  long int info = 0;
118 
119  // call LAPACK routine (see fmatrixev.cc)
120  eigenValuesLapackCall(&jobz, &uplo, &N, &matrixVector[0], &N,
121  &eigenvalues[0], &workSpace[0], &w, &info);
122 
123  if( info != 0 )
124  {
125  std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
126  DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
127  }
128  }
129  }
137  template <int dim, typename K, class C>
138  static void eigenValuesNonSym(const FieldMatrix<K, dim, dim>& matrix,
139  FieldVector<C, dim>& eigenValues)
140  {
141  {
142  const long int N = dim ;
143  const char jobvl = 'n';
144  const char jobvr = 'n';
145 
146  // matrix to put into dgeev
147  double matrixVector[dim * dim];
148 
149  // copy matrix
150  int row = 0;
151  for(int i=0; i<dim; ++i)
152  {
153  for(int j=0; j<dim; ++j, ++row)
154  {
155  matrixVector[ row ] = matrix[ i ][ j ];
156  }
157  }
158 
159  // working memory
160  double eigenR[dim];
161  double eigenI[dim];
162  double work[3*dim];
163 
164  // return value information
165  long int info = 0;
166  long int lwork = 3*dim;
167 
168  // call LAPACK routine (see fmatrixev_ext.cc)
169  eigenValuesNonsymLapackCall(&jobvl, &jobvr, &N, &matrixVector[0], &N,
170  &eigenR[0], &eigenI[0], 0, &N, 0, &N, &work[0],
171  &lwork, &info);
172 
173  if( info != 0 )
174  {
175  std::cerr << "For matrix " << matrix << " eigenvalue calculation failed! " << std::endl;
176  DUNE_THROW(InvalidStateException,"eigenValues: Eigenvalue calculation failed!");
177  }
178  for (int i=0; i<N; ++i) {
179  eigenValues[i].real = eigenR[i];
180  eigenValues[i].imag = eigenI[i];
181  }
182  }
183 
184  }
185 
186  } // end namespace FMatrixHelp
187 
190 } // end namespace Dune
191 #endif
A dense n x m matrix.
Definition: fmatrix.hh:67
vector space out of a tensor product of fields.
Definition: fvector.hh:92
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:307
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:138
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:244
Dune namespace.
Definition: alignment.hh:14
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)