New wrapper for generalized eigenvalues and associated eigenvectors via LAPACK routines sggev/dggev
This commit is contained in:
parent
7e9345d91e
commit
7d8fb3d384
339
nalgebra-lapack/src/generalized_eigenvalues.rs
Normal file
339
nalgebra-lapack/src/generalized_eigenvalues.rs
Normal file
@ -0,0 +1,339 @@
|
||||
#[cfg(feature = "serde-serialize")]
|
||||
use serde::{Deserialize, Serialize};
|
||||
|
||||
use num::Zero;
|
||||
use num_complex::Complex;
|
||||
|
||||
use simba::scalar:: RealField;
|
||||
|
||||
use crate::ComplexHelper;
|
||||
use na::allocator::Allocator;
|
||||
use na::dimension::{Const, Dim};
|
||||
use na::{DefaultAllocator, Matrix, OMatrix, OVector, Scalar};
|
||||
|
||||
use lapack;
|
||||
|
||||
/// Generalized eigenvalues and generalized eigenvectors(left and right) of a pair of N*N square matrices.
|
||||
#[cfg_attr(feature = "serde-serialize", derive(Serialize, Deserialize))]
|
||||
#[cfg_attr(
|
||||
feature = "serde-serialize",
|
||||
serde(
|
||||
bound(serialize = "DefaultAllocator: Allocator<T, D, D> + Allocator<T, D>,
|
||||
OVector<T, D>: Serialize,
|
||||
OMatrix<T, D, D>: Serialize")
|
||||
)
|
||||
)]
|
||||
#[cfg_attr(
|
||||
feature = "serde-serialize",
|
||||
serde(
|
||||
bound(deserialize = "DefaultAllocator: Allocator<T, D, D> + Allocator<T, D>,
|
||||
OVector<T, D>: Serialize,
|
||||
OMatrix<T, D, D>: Deserialize<'de>")
|
||||
)
|
||||
)]
|
||||
#[derive(Clone, Debug)]
|
||||
pub struct GE<T: Scalar, D: Dim>
|
||||
where
|
||||
DefaultAllocator: Allocator<T, D> + Allocator<T, D, D>,
|
||||
{
|
||||
alphar: OVector<T, D>,
|
||||
alphai: OVector<T, D>,
|
||||
beta: OVector<T, D>,
|
||||
vsl: OMatrix<T, D, D>,
|
||||
vsr: OMatrix<T, D, D>,
|
||||
}
|
||||
|
||||
impl<T: Scalar + Copy, D: Dim> Copy for GE<T, D>
|
||||
where
|
||||
DefaultAllocator: Allocator<T, D, D> + Allocator<T, D>,
|
||||
OMatrix<T, D, D>: Copy,
|
||||
OVector<T, D>: Copy,
|
||||
{
|
||||
}
|
||||
|
||||
impl<T: GEScalar + RealField + Copy, D: Dim> GE<T, D>
|
||||
where
|
||||
DefaultAllocator: Allocator<T, D, D> + Allocator<T, D>,
|
||||
{
|
||||
/// Attempts to compute the generalized eigenvalues (and eigenvectors) via the raw returns from LAPACK's
|
||||
/// dggev and sggev routines
|
||||
///
|
||||
/// For each e in generalized eigenvalues and the associated eigenvectors e_l and e_r (left andf right)
|
||||
/// it satisfies e_l*a = e*e_l*b and a*e_r = e*b*e_r
|
||||
///
|
||||
/// Panics if the method did not converge.
|
||||
pub fn new(a: OMatrix<T, D, D>, b: OMatrix<T, D, D>) -> Self {
|
||||
Self::try_new(a, b).expect("Calculation of generalized eigenvalues failed.")
|
||||
}
|
||||
|
||||
/// Attempts to compute the generalized eigenvalues (and eigenvectors) via the raw returns from LAPACK's
|
||||
/// dggev and sggev routines
|
||||
///
|
||||
/// For each e in generalized eigenvalues and the associated eigenvectors e_l and e_r (left andf right)
|
||||
/// it satisfies e_l*a = e*e_l*b and a*e_r = e*b*e_r
|
||||
///
|
||||
/// Returns `None` if the method did not converge.
|
||||
pub fn try_new(mut a: OMatrix<T, D, D>, mut b: OMatrix<T, D, D>) -> Option<Self> {
|
||||
assert!(
|
||||
a.is_square() && b.is_square(),
|
||||
"Unable to compute the generalized eigenvalues of non-square matrices."
|
||||
);
|
||||
|
||||
assert!(
|
||||
a.shape_generic() == b.shape_generic(),
|
||||
"Unable to compute the generalized eigenvalues of two square matrices of different dimensions."
|
||||
);
|
||||
|
||||
let (nrows, ncols) = a.shape_generic();
|
||||
let n = nrows.value();
|
||||
|
||||
let mut info = 0;
|
||||
|
||||
let mut alphar = Matrix::zeros_generic(nrows, Const::<1>);
|
||||
let mut alphai = Matrix::zeros_generic(nrows, Const::<1>);
|
||||
let mut beta = Matrix::zeros_generic(nrows, Const::<1>);
|
||||
let mut vsl = Matrix::zeros_generic(nrows, ncols);
|
||||
let mut vsr = Matrix::zeros_generic(nrows, ncols);
|
||||
|
||||
let lwork = T::xggev_work_size(
|
||||
b'V',
|
||||
b'V',
|
||||
n as i32,
|
||||
a.as_mut_slice(),
|
||||
n as i32,
|
||||
b.as_mut_slice(),
|
||||
n as i32,
|
||||
alphar.as_mut_slice(),
|
||||
alphai.as_mut_slice(),
|
||||
beta.as_mut_slice(),
|
||||
vsl.as_mut_slice(),
|
||||
n as i32,
|
||||
vsr.as_mut_slice(),
|
||||
n as i32,
|
||||
&mut info,
|
||||
);
|
||||
lapack_check!(info);
|
||||
|
||||
let mut work = vec![T::zero(); lwork as usize];
|
||||
|
||||
T::xggev(
|
||||
b'V',
|
||||
b'V',
|
||||
n as i32,
|
||||
a.as_mut_slice(),
|
||||
n as i32,
|
||||
b.as_mut_slice(),
|
||||
n as i32,
|
||||
alphar.as_mut_slice(),
|
||||
alphai.as_mut_slice(),
|
||||
beta.as_mut_slice(),
|
||||
vsl.as_mut_slice(),
|
||||
n as i32,
|
||||
vsr.as_mut_slice(),
|
||||
n as i32,
|
||||
&mut work,
|
||||
lwork,
|
||||
&mut info,
|
||||
);
|
||||
lapack_check!(info);
|
||||
|
||||
Some(GE {
|
||||
alphar,
|
||||
alphai,
|
||||
beta,
|
||||
vsl,
|
||||
vsr,
|
||||
})
|
||||
}
|
||||
|
||||
/// Calculates the generalized eigenvectors (left and right) associated with the generalized eigenvalues
|
||||
pub fn eigenvectors(self) -> (OMatrix<Complex<T>, D, D>, OMatrix<Complex<T>, D, D>)
|
||||
where
|
||||
DefaultAllocator: Allocator<Complex<T>, D, D> + Allocator<Complex<T>, D>,
|
||||
{
|
||||
let n = self.vsl.shape().0;
|
||||
let mut l = self
|
||||
.vsl
|
||||
.clone()
|
||||
.map(|x| Complex::new(x, T::RealField::zero()));
|
||||
let mut r = self
|
||||
.vsr
|
||||
.clone()
|
||||
.map(|x| Complex::new(x, T::RealField::zero()));
|
||||
|
||||
let eigenvalues = &self.eigenvalues();
|
||||
|
||||
let mut ll;
|
||||
let mut c = 0;
|
||||
while c < n {
|
||||
if eigenvalues[c].im.abs() > T::RealField::default_epsilon() && c + 1 < n && {
|
||||
let e_conj = eigenvalues[c].conj();
|
||||
let e = eigenvalues[c + 1];
|
||||
((e_conj.re - e.re).abs() < T::RealField::default_epsilon())
|
||||
&& ((e_conj.im - e.im).abs() < T::RealField::default_epsilon())
|
||||
} {
|
||||
ll = l.column(c + 1).into_owned();
|
||||
l.column_mut(c).zip_apply(&ll, |r, i| {
|
||||
*r = Complex::new(r.re.clone(), i.re);
|
||||
});
|
||||
ll.copy_from(&l.column(c));
|
||||
l.column_mut(c + 1).zip_apply(&ll, |r, i| {
|
||||
*r = i.conj();
|
||||
});
|
||||
|
||||
ll.copy_from(&r.column(c + 1));
|
||||
r.column_mut(c).zip_apply(&ll, |r, i| {
|
||||
*r = Complex::new(r.re, i.re);
|
||||
});
|
||||
ll.copy_from(&r.column(c));
|
||||
r.column_mut(c + 1).zip_apply(&ll, |r, i| {
|
||||
*r = i.conj();
|
||||
});
|
||||
|
||||
c += 2;
|
||||
} else {
|
||||
c += 1;
|
||||
}
|
||||
}
|
||||
|
||||
(l, r)
|
||||
}
|
||||
|
||||
/// computes the generalized eigenvalues
|
||||
#[must_use]
|
||||
pub fn eigenvalues(&self) -> OVector<Complex<T>, D>
|
||||
where
|
||||
DefaultAllocator: Allocator<Complex<T>, D>,
|
||||
{
|
||||
let mut out = Matrix::zeros_generic(self.vsl.shape_generic().0, Const::<1>);
|
||||
|
||||
for i in 0..out.len() {
|
||||
out[i] = if self.beta[i].clone().abs() < T::RealField::default_epsilon() {
|
||||
Complex::zero()
|
||||
} else {
|
||||
let mut cr = self.alphar[i].clone();
|
||||
let mut ci = self.alphai[i].clone();
|
||||
let b = self.beta[i].clone();
|
||||
|
||||
if cr.clone().abs() < T::RealField::default_epsilon() {
|
||||
cr = T::RealField::zero()
|
||||
} else {
|
||||
cr = cr / b.clone()
|
||||
};
|
||||
|
||||
if ci.clone().abs() < T::RealField::default_epsilon() {
|
||||
ci = T::RealField::zero()
|
||||
} else {
|
||||
ci = ci / b
|
||||
};
|
||||
|
||||
Complex::new(cr, ci)
|
||||
}
|
||||
}
|
||||
|
||||
out
|
||||
}
|
||||
}
|
||||
|
||||
/*
|
||||
*
|
||||
* Lapack functions dispatch.
|
||||
*
|
||||
*/
|
||||
/// Trait implemented by scalars for which Lapack implements the RealField GE decomposition.
|
||||
pub trait GEScalar: Scalar {
|
||||
#[allow(missing_docs)]
|
||||
fn xggev(
|
||||
jobvsl: u8,
|
||||
jobvsr: u8,
|
||||
n: i32,
|
||||
a: &mut [Self],
|
||||
lda: i32,
|
||||
b: &mut [Self],
|
||||
ldb: i32,
|
||||
alphar: &mut [Self],
|
||||
alphai: &mut [Self],
|
||||
beta: &mut [Self],
|
||||
vsl: &mut [Self],
|
||||
ldvsl: i32,
|
||||
vsr: &mut [Self],
|
||||
ldvsr: i32,
|
||||
work: &mut [Self],
|
||||
lwork: i32,
|
||||
info: &mut i32,
|
||||
);
|
||||
|
||||
#[allow(missing_docs)]
|
||||
fn xggev_work_size(
|
||||
jobvsl: u8,
|
||||
jobvsr: u8,
|
||||
n: i32,
|
||||
a: &mut [Self],
|
||||
lda: i32,
|
||||
b: &mut [Self],
|
||||
ldb: i32,
|
||||
alphar: &mut [Self],
|
||||
alphai: &mut [Self],
|
||||
beta: &mut [Self],
|
||||
vsl: &mut [Self],
|
||||
ldvsl: i32,
|
||||
vsr: &mut [Self],
|
||||
ldvsr: i32,
|
||||
info: &mut i32,
|
||||
) -> i32;
|
||||
}
|
||||
|
||||
macro_rules! real_eigensystem_scalar_impl (
|
||||
($N: ty, $xggev: path) => (
|
||||
impl GEScalar for $N {
|
||||
#[inline]
|
||||
fn xggev(jobvsl: u8,
|
||||
jobvsr: u8,
|
||||
n: i32,
|
||||
a: &mut [$N],
|
||||
lda: i32,
|
||||
b: &mut [$N],
|
||||
ldb: i32,
|
||||
alphar: &mut [$N],
|
||||
alphai: &mut [$N],
|
||||
beta : &mut [$N],
|
||||
vsl: &mut [$N],
|
||||
ldvsl: i32,
|
||||
vsr: &mut [$N],
|
||||
ldvsr: i32,
|
||||
work: &mut [$N],
|
||||
lwork: i32,
|
||||
info: &mut i32) {
|
||||
unsafe { $xggev(jobvsl, jobvsr, n, a, lda, b, ldb, alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr, work, lwork, info); }
|
||||
}
|
||||
|
||||
|
||||
#[inline]
|
||||
fn xggev_work_size(jobvsl: u8,
|
||||
jobvsr: u8,
|
||||
n: i32,
|
||||
a: &mut [$N],
|
||||
lda: i32,
|
||||
b: &mut [$N],
|
||||
ldb: i32,
|
||||
alphar: &mut [$N],
|
||||
alphai: &mut [$N],
|
||||
beta : &mut [$N],
|
||||
vsl: &mut [$N],
|
||||
ldvsl: i32,
|
||||
vsr: &mut [$N],
|
||||
ldvsr: i32,
|
||||
info: &mut i32)
|
||||
-> i32 {
|
||||
let mut work = [ Zero::zero() ];
|
||||
let lwork = -1 as i32;
|
||||
|
||||
unsafe { $xggev(jobvsl, jobvsr, n, a, lda, b, ldb, alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr, &mut work, lwork, info); }
|
||||
ComplexHelper::real_part(work[0]) as i32
|
||||
}
|
||||
}
|
||||
)
|
||||
);
|
||||
|
||||
real_eigensystem_scalar_impl!(f32, lapack::sggev);
|
||||
real_eigensystem_scalar_impl!(f64, lapack::dggev);
|
@ -87,6 +87,7 @@ mod hessenberg;
|
||||
mod lu;
|
||||
mod qr;
|
||||
mod qz;
|
||||
mod generalized_eigenvalues;
|
||||
mod schur;
|
||||
mod svd;
|
||||
mod symmetric_eigen;
|
||||
@ -99,6 +100,7 @@ pub use self::hessenberg::Hessenberg;
|
||||
pub use self::lu::{LUScalar, LU};
|
||||
pub use self::qr::QR;
|
||||
pub use self::qz::QZ;
|
||||
pub use self::generalized_eigenvalues::GE;
|
||||
pub use self::schur::Schur;
|
||||
pub use self::svd::SVD;
|
||||
pub use self::symmetric_eigen::SymmetricEigen;
|
||||
|
59
nalgebra-lapack/tests/linalg/generalized_eigenvalues.rs
Normal file
59
nalgebra-lapack/tests/linalg/generalized_eigenvalues.rs
Normal file
@ -0,0 +1,59 @@
|
||||
use na::dimension::{Const, Dynamic};
|
||||
use na::{DMatrix, EuclideanNorm, Norm, OMatrix};
|
||||
use nl::GE;
|
||||
use num_complex::Complex;
|
||||
use simba::scalar::ComplexField;
|
||||
use std::cmp;
|
||||
|
||||
use crate::proptest::*;
|
||||
use proptest::{prop_assert, proptest};
|
||||
|
||||
proptest! {
|
||||
#[test]
|
||||
fn ge(n in PROPTEST_MATRIX_DIM) {
|
||||
let n = cmp::max(1, cmp::min(n, 10));
|
||||
let a = DMatrix::<f64>::new_random(n, n);
|
||||
let b = DMatrix::<f64>::new_random(n, n);
|
||||
let a_condition_no = a.clone().try_inverse().and_then(|x| Some(EuclideanNorm.norm(&x)* EuclideanNorm.norm(&a)));
|
||||
let b_condition_no = b.clone().try_inverse().and_then(|x| Some(EuclideanNorm.norm(&x)* EuclideanNorm.norm(&b)));
|
||||
|
||||
if a_condition_no.unwrap_or(200000.0) < 10.0 && b_condition_no.unwrap_or(200000.0) < 10.0 {
|
||||
let a_c =a.clone().map(|x| Complex::new(x, 0.0));
|
||||
let b_c = b.clone().map(|x| Complex::new(x, 0.0));
|
||||
|
||||
let ge = GE::new(a.clone(), b.clone());
|
||||
let (vsl,vsr) = ge.clone().eigenvectors();
|
||||
let eigenvalues = ge.clone().eigenvalues();
|
||||
|
||||
for i in 0..n {
|
||||
let left_eigenvector = &vsl.column(i);
|
||||
prop_assert!(relative_eq!((left_eigenvector.transpose()*&a_c - left_eigenvector.transpose()*&b_c*eigenvalues[i]).map(|x| x.modulus()), OMatrix::zeros_generic(Const::<1>,Dynamic::new(n)) ,epsilon = 1.0e-7));
|
||||
|
||||
let right_eigenvector = &vsr.column(i);
|
||||
prop_assert!(relative_eq!((&a_c*right_eigenvector - &b_c*right_eigenvector*eigenvalues[i]).map(|x| x.modulus()), OMatrix::zeros_generic(Dynamic::new(n), Const::<1>) ,epsilon = 1.0e-7));
|
||||
};
|
||||
};
|
||||
}
|
||||
|
||||
#[test]
|
||||
fn ge_static(a in matrix4(), b in matrix4()) {
|
||||
let a_condition_no = a.clone().try_inverse().and_then(|x| Some(EuclideanNorm.norm(&x)* EuclideanNorm.norm(&a)));
|
||||
let b_condition_no = b.clone().try_inverse().and_then(|x| Some(EuclideanNorm.norm(&x)* EuclideanNorm.norm(&b)));
|
||||
|
||||
if a_condition_no.unwrap_or(200000.0) < 10.0 && b_condition_no.unwrap_or(200000.0) < 10.0{
|
||||
let ge = GE::new(a.clone(), b.clone());
|
||||
let a_c =a.clone().map(|x| Complex::new(x, 0.0));
|
||||
let b_c = b.clone().map(|x| Complex::new(x, 0.0));
|
||||
let (vsl,vsr) = ge.eigenvectors();
|
||||
let eigenvalues = ge.eigenvalues();
|
||||
|
||||
for i in 0..4 {
|
||||
let left_eigenvector = &vsl.column(i);
|
||||
prop_assert!(relative_eq!((left_eigenvector.transpose()*&a_c - left_eigenvector.transpose()*&b_c*eigenvalues[i]).map(|x| x.modulus()), OMatrix::zeros_generic(Const::<1>,Const::<4>) ,epsilon = 1.0e-7));
|
||||
|
||||
let right_eigenvector = &vsr.column(i);
|
||||
prop_assert!(relative_eq!((&a_c*right_eigenvector - &b_c*right_eigenvector*eigenvalues[i]).map(|x| x.modulus()), OMatrix::zeros_generic(Const::<4>, Const::<1>) ,epsilon = 1.0e-7));
|
||||
};
|
||||
};
|
||||
}
|
||||
}
|
@ -2,6 +2,7 @@ mod cholesky;
|
||||
mod lu;
|
||||
mod qr;
|
||||
mod qz;
|
||||
mod generalized_eigenvalues;
|
||||
mod real_eigensystem;
|
||||
mod schur;
|
||||
mod svd;
|
||||
|
Loading…
Reference in New Issue
Block a user