Dune Core Modules (2.6.0)

fmatrixev.hh File Reference

Eigenvalue computations for the FieldMatrix class. More...

#include <iostream>
#include <cmath>
#include <cassert>
#include <dune/common/exceptions.hh>
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <dune/common/math.hh>

Go to the source code of this file.

Namespaces

namespace  Dune
 Dune namespace.
 

Functions

template<typename K >
static void Dune::FMatrixHelp::eigenValues (const FieldMatrix< K, 1, 1 > &matrix, FieldVector< K, 1 > &eigenvalues)
 calculates the eigenvalues of a symmetric field matrix More...
 
template<typename K >
static void Dune::FMatrixHelp::eigenValues (const FieldMatrix< K, 2, 2 > &matrix, FieldVector< K, 2 > &eigenvalues)
 calculates the eigenvalues of a symmetric field matrix More...
 
template<typename K >
static void Dune::FMatrixHelp::eigenValues (const FieldMatrix< K, 3, 3 > &matrix, FieldVector< K, 3 > &eigenvalues)
 Calculates the eigenvalues of a symmetric 3x3 field matrix. More...
 
template<int dim, typename K >
static void Dune::FMatrixHelp::eigenValues (const FieldMatrix< K, dim, dim > &matrix, FieldVector< K, dim > &eigenvalues)
 calculates the eigenvalues of a symmetric field matrix More...
 
template<int dim, typename K , class C >
static void Dune::FMatrixHelp::eigenValuesNonSym (const FieldMatrix< K, dim, dim > &matrix, FieldVector< C, dim > &eigenValues)
 calculates the eigenvalues of a symmetric field matrix More...
 

Detailed Description

Eigenvalue computations for the FieldMatrix class.

Function Documentation

◆ eigenValues() [1/4]

template<typename K >
static void Dune::FMatrixHelp::eigenValues ( const FieldMatrix< K, 1, 1 > &  matrix,
FieldVector< K, 1 > &  eigenvalues 
)
static

calculates the eigenvalues of a symmetric field matrix

Parameters
[in]matrixmatrix eigenvalues are calculated for
[out]eigenvaluesFieldVector that contains eigenvalues in ascending order

References Dune::FMatrixHelp::eigenValues().

Referenced by Dune::FMatrixHelp::eigenValues().

◆ eigenValues() [2/4]

template<typename K >
static void Dune::FMatrixHelp::eigenValues ( const FieldMatrix< K, 2, 2 > &  matrix,
FieldVector< K, 2 > &  eigenvalues 
)
static

calculates the eigenvalues of a symmetric field matrix

Parameters
[in]matrixmatrix eigenvalues are calculated for
[out]eigenvaluesFieldVector that contains eigenvalues in ascending order

References DUNE_THROW, and Dune::FMatrixHelp::eigenValues().

◆ eigenValues() [3/4]

template<typename K >
static void Dune::FMatrixHelp::eigenValues ( const FieldMatrix< K, 3, 3 > &  matrix,
FieldVector< K, 3 > &  eigenvalues 
)
static

Calculates the eigenvalues of a symmetric 3x3 field matrix.

Parameters
[in]matrixmatrix eigenvalues are calculated for
[out]eigenvaluesEigenvalues in ascending order
Note
If the input matrix is not symmetric the behavior of this method is undefined.

This implementation was adapted from the pseudo-code (Python?) implementation found on http://en.wikipedia.org/wiki/Eigenvalue_algorithm (retrieved late August 2014). Wikipedia claims to have taken it from Smith, Oliver K. (April 1961), Eigenvalues of a symmetric 3 × 3 matrix., Communications of the ACM 4 (4): 168, doi:10.1145/355578.366316

References Dune::DenseMatrix< MAT >::determinant(), and Dune::FMatrixHelp::eigenValues().

◆ eigenValues() [4/4]

template<int dim, typename K >
static void Dune::FMatrixHelp::eigenValues ( const FieldMatrix< K, dim, dim > &  matrix,
FieldVector< K, dim > &  eigenvalues 
)
static

calculates the eigenvalues of a symmetric field matrix

Parameters
[in]matrixmatrix eigenvalues are calculated for
[out]eigenvaluesFieldVector that contains eigenvalues in ascending order
Note
LAPACK::dsyev is used to calculate the eigenvalues

References Dune::FMatrixHelp::eigenValues().

◆ eigenValuesNonSym()

template<int dim, typename K , class C >
static void Dune::FMatrixHelp::eigenValuesNonSym ( const FieldMatrix< K, dim, dim > &  matrix,
FieldVector< C, dim > &  eigenValues 
)
static

calculates the eigenvalues of a symmetric field matrix

Parameters
[in]matrixmatrix eigenvalues are calculated for
[out]eigenValuesFieldVector that contains eigenvalues in ascending order
Note
LAPACK::dgeev is used to calculate the eigen values

References Dune::FMatrixHelp::eigenValuesNonSym().

Referenced by Dune::FMatrixHelp::eigenValuesNonSym().

Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)