nalgebra-lapack: add doc + fix warnings.
This commit is contained in:
parent
b22eb91a16
commit
b84c7e10df
|
@ -36,17 +36,37 @@ impl<N: CholeskyScalar + Zero, D: Dim> Cholesky<N, D>
|
||||||
Some(Cholesky { l: m })
|
Some(Cholesky { l: m })
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/// Retrieves the lower-triangular factor of the cholesky decomposition.
|
||||||
pub fn unpack(mut self) -> MatrixN<N, D> {
|
pub fn unpack(mut self) -> MatrixN<N, D> {
|
||||||
self.l.fill_upper_triangle(Zero::zero(), 1);
|
self.l.fill_upper_triangle(Zero::zero(), 1);
|
||||||
self.l
|
self.l
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/// Retrieves the lower-triangular factor of che cholesky decomposition, without zeroing-out
|
||||||
|
/// its strict upper-triangular part.
|
||||||
|
///
|
||||||
|
/// This is an allocation-less version of `self.l()`. The values of the strict upper-triangular
|
||||||
|
/// part are garbage and should be ignored by further computations.
|
||||||
|
pub fn unpack_dirty(self) -> MatrixN<N, D> {
|
||||||
|
self.l
|
||||||
|
}
|
||||||
|
|
||||||
|
/// Retrieves the lower-triangular factor of the cholesky decomposition.
|
||||||
pub fn l(&self) -> MatrixN<N, D> {
|
pub fn l(&self) -> MatrixN<N, D> {
|
||||||
let mut res = self.l.clone();
|
let mut res = self.l.clone();
|
||||||
res.fill_upper_triangle(Zero::zero(), 1);
|
res.fill_upper_triangle(Zero::zero(), 1);
|
||||||
res
|
res
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/// Retrieves the lower-triangular factor of the cholesky decomposition, without zeroing-out
|
||||||
|
/// its strict upper-triangular part.
|
||||||
|
///
|
||||||
|
/// This is an allocation-less version of `self.l()`. The values of the strict upper-triangular
|
||||||
|
/// part are garbage and should be ignored by further computations.
|
||||||
|
pub fn l_dirty(&self) -> &MatrixN<N, D> {
|
||||||
|
&self.l
|
||||||
|
}
|
||||||
|
|
||||||
/// Solves the symmetric-definite-positive linear system `self * x = b`, where `x` is the
|
/// Solves the symmetric-definite-positive linear system `self * x = b`, where `x` is the
|
||||||
/// unknown to be determined.
|
/// unknown to be determined.
|
||||||
pub fn solve<R2: Dim, C2: Dim, S2>(&self, b: &Matrix<N, R2, C2, S2>) -> Option<MatrixMN<N, R2, C2>>
|
pub fn solve<R2: Dim, C2: Dim, S2>(&self, b: &Matrix<N, R2, C2, S2>) -> Option<MatrixMN<N, R2, C2>>
|
||||||
|
@ -110,8 +130,11 @@ impl<N: CholeskyScalar + Zero, D: Dim> Cholesky<N, D>
|
||||||
/// Trait implemented by floats (`f32`, `f64`) and complex floats (`Complex<f32>`, `Complex<f64>`)
|
/// Trait implemented by floats (`f32`, `f64`) and complex floats (`Complex<f32>`, `Complex<f64>`)
|
||||||
/// supported by the cholesky decompotition.
|
/// supported by the cholesky decompotition.
|
||||||
pub trait CholeskyScalar: Scalar {
|
pub trait CholeskyScalar: Scalar {
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xpotrf(uplo: u8, n: i32, a: &mut [Self], lda: i32, info: &mut i32);
|
fn xpotrf(uplo: u8, n: i32, a: &mut [Self], lda: i32, info: &mut i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xpotrs(uplo: u8, n: i32, nrhs: i32, a: &[Self], lda: i32, b: &mut [Self], ldb: i32, info: &mut i32);
|
fn xpotrs(uplo: u8, n: i32, nrhs: i32, a: &[Self], lda: i32, b: &mut [Self], ldb: i32, info: &mut i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xpotri(uplo: u8, n: i32, a: &mut [Self], lda: i32, info: &mut i32);
|
fn xpotri(uplo: u8, n: i32, a: &mut [Self], lda: i32, info: &mut i32);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
|
@ -15,8 +15,11 @@ use lapack::fortran as interface;
|
||||||
pub struct Eigen<N: Scalar, D: Dim>
|
pub struct Eigen<N: Scalar, D: Dim>
|
||||||
where DefaultAllocator: Allocator<N, D> +
|
where DefaultAllocator: Allocator<N, D> +
|
||||||
Allocator<N, D, D> {
|
Allocator<N, D, D> {
|
||||||
|
/// The eigenvalues of the decomposed matrix.
|
||||||
pub eigenvalues: VectorN<N, D>,
|
pub eigenvalues: VectorN<N, D>,
|
||||||
|
/// The (right) eigenvectors of the decomposed matrix.
|
||||||
pub eigenvectors: Option<MatrixN<N, D>>,
|
pub eigenvectors: Option<MatrixN<N, D>>,
|
||||||
|
/// The left eigenvectors of the decomposed matrix.
|
||||||
pub left_eigenvectors: Option<MatrixN<N, D>>
|
pub left_eigenvectors: Option<MatrixN<N, D>>
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -125,7 +128,7 @@ impl<N: EigenScalar + Real, D: Dim> Eigen<N, D>
|
||||||
where DefaultAllocator: Allocator<Complex<N>, D> {
|
where DefaultAllocator: Allocator<Complex<N>, D> {
|
||||||
assert!(m.is_square(), "Unable to compute the eigenvalue decomposition of a non-square matrix.");
|
assert!(m.is_square(), "Unable to compute the eigenvalue decomposition of a non-square matrix.");
|
||||||
|
|
||||||
let (nrows, ncols) = m.data.shape();
|
let nrows = m.data.shape().0;
|
||||||
let n = nrows.value();
|
let n = nrows.value();
|
||||||
|
|
||||||
let lda = n as i32;
|
let lda = n as i32;
|
||||||
|
@ -181,11 +184,15 @@ impl<N: EigenScalar + Real, D: Dim> Eigen<N, D>
|
||||||
* Lapack functions dispatch.
|
* Lapack functions dispatch.
|
||||||
*
|
*
|
||||||
*/
|
*/
|
||||||
|
/// Trait implemented by scalar type for which Lapack funtion exist to compute the
|
||||||
|
/// eigendecomposition.
|
||||||
pub trait EigenScalar: Scalar {
|
pub trait EigenScalar: Scalar {
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xgeev(jobvl: u8, jobvr: u8, n: i32, a: &mut [Self], lda: i32,
|
fn xgeev(jobvl: u8, jobvr: u8, n: i32, a: &mut [Self], lda: i32,
|
||||||
wr: &mut [Self], wi: &mut [Self],
|
wr: &mut [Self], wi: &mut [Self],
|
||||||
vl: &mut [Self], ldvl: i32, vr: &mut [Self], ldvr: i32,
|
vl: &mut [Self], ldvl: i32, vr: &mut [Self], ldvr: i32,
|
||||||
work: &mut [Self], lwork: i32, info: &mut i32);
|
work: &mut [Self], lwork: i32, info: &mut i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xgeev_work_size(jobvl: u8, jobvr: u8, n: i32, a: &mut [Self], lda: i32,
|
fn xgeev_work_size(jobvl: u8, jobvr: u8, n: i32, a: &mut [Self], lda: i32,
|
||||||
wr: &mut [Self], wi: &mut [Self], vl: &mut [Self], ldvl: i32,
|
wr: &mut [Self], wi: &mut [Self], vl: &mut [Self], ldvl: i32,
|
||||||
vr: &mut [Self], ldvr: i32, info: &mut i32) -> i32;
|
vr: &mut [Self], ldvr: i32, info: &mut i32) -> i32;
|
||||||
|
|
|
@ -94,9 +94,12 @@ pub trait HessenbergScalar: Scalar {
|
||||||
tau: &mut [Self], info: &mut i32) -> i32;
|
tau: &mut [Self], info: &mut i32) -> i32;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/// Trait implemented by scalars for which Lapack implements the hessenberg decomposition.
|
||||||
pub trait HessenbergReal: HessenbergScalar {
|
pub trait HessenbergReal: HessenbergScalar {
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xorghr(n: i32, ilo: i32, ihi: i32, a: &mut [Self], lda: i32, tau: &[Self],
|
fn xorghr(n: i32, ilo: i32, ihi: i32, a: &mut [Self], lda: i32, tau: &[Self],
|
||||||
work: &mut [Self], lwork: i32, info: &mut i32);
|
work: &mut [Self], lwork: i32, info: &mut i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xorghr_work_size(n: i32, ilo: i32, ihi: i32, a: &mut [Self], lda: i32,
|
fn xorghr_work_size(n: i32, ilo: i32, ihi: i32, a: &mut [Self], lda: i32,
|
||||||
tau: &[Self], info: &mut i32) -> i32;
|
tau: &[Self], info: &mut i32) -> i32;
|
||||||
}
|
}
|
||||||
|
|
|
@ -1,10 +1,69 @@
|
||||||
|
//! # nalgebra-lapack
|
||||||
|
//!
|
||||||
|
//! Rust library for linear algebra using nalgebra and LAPACK.
|
||||||
|
//!
|
||||||
|
//! ## Documentation
|
||||||
|
//!
|
||||||
|
//! Documentation is available [here](https://docs.rs/nalgebra-lapack/).
|
||||||
|
//!
|
||||||
|
//! ## License
|
||||||
|
//!
|
||||||
|
//! MIT
|
||||||
|
//!
|
||||||
|
//! ## Cargo features to select lapack provider
|
||||||
|
//!
|
||||||
|
//! Like the [lapack crate](https://crates.io/crates/lapack) from which this
|
||||||
|
//! behavior is inherited, nalgebra-lapack uses [cargo
|
||||||
|
//! features](http://doc.crates.io/manifest.html#the-[features]-section) to select
|
||||||
|
//! which lapack provider (or implementation) is used. Command line arguments to
|
||||||
|
//! cargo are the easiest way to do this, and the best provider depends on your
|
||||||
|
//! particular system. In some cases, the providers can be further tuned with
|
||||||
|
//! environment variables.
|
||||||
|
//!
|
||||||
|
//! Below are given examples of how to invoke `cargo build` on two different systems
|
||||||
|
//! using two different providers. The `--no-default-features --features "provider"`
|
||||||
|
//! arguments will be consistent for other `cargo` commands.
|
||||||
|
//!
|
||||||
|
//! ### Ubuntu
|
||||||
|
//!
|
||||||
|
//! As tested on Ubuntu 12.04, do this to build the lapack package against
|
||||||
|
//! the system installation of netlib without LAPACKE (note the E) or
|
||||||
|
//! CBLAS:
|
||||||
|
//!
|
||||||
|
//! sudo apt-get install gfortran libblas3gf liblapack3gf
|
||||||
|
//! export CARGO_FEATURE_SYSTEM_NETLIB=1
|
||||||
|
//! export CARGO_FEATURE_EXCLUDE_LAPACKE=1
|
||||||
|
//! export CARGO_FEATURE_EXCLUDE_CBLAS=1
|
||||||
|
//!
|
||||||
|
//! export CARGO_FEATURES='--no-default-features --features netlib'
|
||||||
|
//! cargo build ${CARGO_FEATURES}
|
||||||
|
//!
|
||||||
|
//! ### Mac OS X
|
||||||
|
//!
|
||||||
|
//! On Mac OS X, do this to use Apple's Accelerate framework:
|
||||||
|
//!
|
||||||
|
//! export CARGO_FEATURES='--no-default-features --features accelerate'
|
||||||
|
//! cargo build ${CARGO_FEATURES}
|
||||||
|
//!
|
||||||
|
//! [version-img]: https://img.shields.io/crates/v/nalgebra-lapack.svg
|
||||||
|
//! [version-url]: https://crates.io/crates/nalgebra-lapack
|
||||||
|
//! [status-img]: https://travis-ci.org/strawlab/nalgebra-lapack.svg?branch=master
|
||||||
|
//! [status-url]: https://travis-ci.org/strawlab/nalgebra-lapack
|
||||||
|
//! [doc-img]: https://docs.rs/nalgebra-lapack/badge.svg
|
||||||
|
//! [doc-url]: https://docs.rs/nalgebra-lapack/
|
||||||
|
//!
|
||||||
|
//! ## Contributors
|
||||||
|
//! This integration of LAPACK on nalgebra was
|
||||||
|
//! [initiated](https://github.com/strawlab/nalgebra-lapack) by Andrew Straw. It
|
||||||
|
//! then became officially supported and integrated to the main nalgebra
|
||||||
|
//! repository.
|
||||||
|
|
||||||
#![deny(non_camel_case_types)]
|
#![deny(non_camel_case_types)]
|
||||||
#![deny(unused_parens)]
|
#![deny(unused_parens)]
|
||||||
#![deny(non_upper_case_globals)]
|
#![deny(non_upper_case_globals)]
|
||||||
#![deny(unused_qualifications)]
|
#![deny(unused_qualifications)]
|
||||||
#![deny(unused_results)]
|
#![deny(unused_results)]
|
||||||
#![warn(missing_docs)]
|
#![deny(missing_docs)]
|
||||||
#![doc(html_root_url = "http://nalgebra.org/rustdoc")]
|
#![doc(html_root_url = "http://nalgebra.org/rustdoc")]
|
||||||
|
|
||||||
extern crate num_traits as num;
|
extern crate num_traits as num;
|
||||||
|
|
|
@ -33,6 +33,7 @@ impl<N: LUScalar, R: Dim, C: Dim> LU<N, R, C>
|
||||||
Allocator<N, DimMinimum<R, C>, C> +
|
Allocator<N, DimMinimum<R, C>, C> +
|
||||||
Allocator<i32, DimMinimum<R, C>> {
|
Allocator<i32, DimMinimum<R, C>> {
|
||||||
|
|
||||||
|
/// Computes the LU decomposition with partial (row) pivoting of `matrix`.
|
||||||
pub fn new(mut m: MatrixMN<N, R, C>) -> Self {
|
pub fn new(mut m: MatrixMN<N, R, C>) -> Self {
|
||||||
let (nrows, ncols) = m.data.shape();
|
let (nrows, ncols) = m.data.shape();
|
||||||
let min_nrows_ncols = nrows.min(ncols);
|
let min_nrows_ncols = nrows.min(ncols);
|
||||||
|
@ -238,13 +239,19 @@ impl<N: LUScalar, D: Dim> LU<N, D, D>
|
||||||
* Lapack functions dispatch.
|
* Lapack functions dispatch.
|
||||||
*
|
*
|
||||||
*/
|
*/
|
||||||
|
/// Trait implemented by scalars for which Lapack implements the LU decomposition.
|
||||||
pub trait LUScalar: Scalar {
|
pub trait LUScalar: Scalar {
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xgetrf(m: i32, n: i32, a: &mut [Self], lda: i32, ipiv: &mut [i32], info: &mut i32);
|
fn xgetrf(m: i32, n: i32, a: &mut [Self], lda: i32, ipiv: &mut [i32], info: &mut i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xlaswp(n: i32, a: &mut [Self], lda: i32, k1: i32, k2: i32, ipiv: &[i32], incx: i32);
|
fn xlaswp(n: i32, a: &mut [Self], lda: i32, k1: i32, k2: i32, ipiv: &[i32], incx: i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xgetrs(trans: u8, n: i32, nrhs: i32, a: &[Self], lda: i32, ipiv: &[i32],
|
fn xgetrs(trans: u8, n: i32, nrhs: i32, a: &[Self], lda: i32, ipiv: &[i32],
|
||||||
b: &mut [Self], ldb: i32, info: &mut i32);
|
b: &mut [Self], ldb: i32, info: &mut i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xgetri(n: i32, a: &mut [Self], lda: i32, ipiv: &[i32],
|
fn xgetri(n: i32, a: &mut [Self], lda: i32, ipiv: &[i32],
|
||||||
work: &mut [Self], lwork: i32, info: &mut i32);
|
work: &mut [Self], lwork: i32, info: &mut i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xgetri_work_size(n: i32, a: &mut [Self], lda: i32, ipiv: &[i32], info: &mut i32) -> i32;
|
fn xgetri_work_size(n: i32, a: &mut [Self], lda: i32, ipiv: &[i32], info: &mut i32) -> i32;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
|
@ -102,6 +102,8 @@ impl<N: QRReal + Zero, R: DimMin<C>, C: Dim> QR<N, R, C>
|
||||||
* Lapack functions dispatch.
|
* Lapack functions dispatch.
|
||||||
*
|
*
|
||||||
*/
|
*/
|
||||||
|
/// Trait implemented by scalar types for which Lapack funtion exist to compute the
|
||||||
|
/// QR decomposition.
|
||||||
pub trait QRScalar: Scalar {
|
pub trait QRScalar: Scalar {
|
||||||
fn xgeqrf(m: i32, n: i32, a: &mut [Self], lda: i32, tau: &mut [Self],
|
fn xgeqrf(m: i32, n: i32, a: &mut [Self], lda: i32, tau: &mut [Self],
|
||||||
work: &mut [Self], lwork: i32, info: &mut i32);
|
work: &mut [Self], lwork: i32, info: &mut i32);
|
||||||
|
@ -110,10 +112,14 @@ pub trait QRScalar: Scalar {
|
||||||
tau: &mut [Self], info: &mut i32) -> i32;
|
tau: &mut [Self], info: &mut i32) -> i32;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/// Trait implemented by reals for which Lapack funtion exist to compute the
|
||||||
|
/// QR decomposition.
|
||||||
pub trait QRReal: QRScalar {
|
pub trait QRReal: QRScalar {
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xorgqr(m: i32, n: i32, k: i32, a: &mut [Self], lda: i32, tau: &[Self], work: &mut [Self],
|
fn xorgqr(m: i32, n: i32, k: i32, a: &mut [Self], lda: i32, tau: &[Self], work: &mut [Self],
|
||||||
lwork: i32, info: &mut i32);
|
lwork: i32, info: &mut i32);
|
||||||
|
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xorgqr_work_size(m: i32, n: i32, k: i32, a: &mut [Self], lda: i32,
|
fn xorgqr_work_size(m: i32, n: i32, k: i32, a: &mut [Self], lda: i32,
|
||||||
tau: &[Self], info: &mut i32) -> i32;
|
tau: &[Self], info: &mut i32) -> i32;
|
||||||
}
|
}
|
||||||
|
|
|
@ -23,7 +23,7 @@ pub struct RealSchur<N: Scalar, D: Dim>
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
impl<N: EigenScalar + Real, D: Dim> RealSchur<N, D>
|
impl<N: RealSchurScalar + Real, D: Dim> RealSchur<N, D>
|
||||||
where DefaultAllocator: Allocator<N, D, D> +
|
where DefaultAllocator: Allocator<N, D, D> +
|
||||||
Allocator<N, D> {
|
Allocator<N, D> {
|
||||||
/// Computes the eigenvalues and real Schur foorm of the matrix `m`.
|
/// Computes the eigenvalues and real Schur foorm of the matrix `m`.
|
||||||
|
@ -106,7 +106,9 @@ impl<N: EigenScalar + Real, D: Dim> RealSchur<N, D>
|
||||||
* Lapack functions dispatch.
|
* Lapack functions dispatch.
|
||||||
*
|
*
|
||||||
*/
|
*/
|
||||||
pub trait EigenScalar: Scalar {
|
/// Trait implemented by scalars for which Lapack implements the Real Schur decomposition.
|
||||||
|
pub trait RealSchurScalar: Scalar {
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xgees(jobvs: u8,
|
fn xgees(jobvs: u8,
|
||||||
sort: u8,
|
sort: u8,
|
||||||
// select: ???
|
// select: ???
|
||||||
|
@ -122,7 +124,8 @@ pub trait EigenScalar: Scalar {
|
||||||
lwork: i32,
|
lwork: i32,
|
||||||
bwork: &mut [i32],
|
bwork: &mut [i32],
|
||||||
info: &mut i32);
|
info: &mut i32);
|
||||||
|
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xgees_work_size(jobvs: u8,
|
fn xgees_work_size(jobvs: u8,
|
||||||
sort: u8,
|
sort: u8,
|
||||||
// select: ???
|
// select: ???
|
||||||
|
@ -141,7 +144,7 @@ pub trait EigenScalar: Scalar {
|
||||||
|
|
||||||
macro_rules! real_eigensystem_scalar_impl (
|
macro_rules! real_eigensystem_scalar_impl (
|
||||||
($N: ty, $xgees: path) => (
|
($N: ty, $xgees: path) => (
|
||||||
impl EigenScalar for $N {
|
impl RealSchurScalar for $N {
|
||||||
#[inline]
|
#[inline]
|
||||||
fn xgees(jobvs: u8,
|
fn xgees(jobvs: u8,
|
||||||
sort: u8,
|
sort: u8,
|
||||||
|
|
|
@ -15,9 +15,12 @@ pub struct SVD<N: Scalar, R: DimMin<C>, C: Dim>
|
||||||
where DefaultAllocator: Allocator<N, R, R> +
|
where DefaultAllocator: Allocator<N, R, R> +
|
||||||
Allocator<N, DimMinimum<R, C>> +
|
Allocator<N, DimMinimum<R, C>> +
|
||||||
Allocator<N, C, C> {
|
Allocator<N, C, C> {
|
||||||
|
/// The left-singular vectors `U` of this SVD.
|
||||||
pub u: MatrixN<N, R>,
|
pub u: MatrixN<N, R>,
|
||||||
pub s: VectorN<N, DimMinimum<R, C>>,
|
/// The right-singular vectors `V^t` of this SVD.
|
||||||
pub vt: MatrixN<N, C>
|
pub vt: MatrixN<N, C>,
|
||||||
|
/// The singular values of this SVD.
|
||||||
|
pub singular_values: VectorN<N, DimMinimum<R, C>>
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
@ -37,6 +40,7 @@ impl<N: SVDScalar<R, C>, R: DimMin<C>, C: Dim> SVD<N, R, C>
|
||||||
Allocator<N, R, C> +
|
Allocator<N, R, C> +
|
||||||
Allocator<N, DimMinimum<R, C>> +
|
Allocator<N, DimMinimum<R, C>> +
|
||||||
Allocator<N, C, C> {
|
Allocator<N, C, C> {
|
||||||
|
/// Computes the Singular Value Decomposition of `matrix`.
|
||||||
pub fn new(m: MatrixMN<N, R, C>) -> Option<Self> {
|
pub fn new(m: MatrixMN<N, R, C>) -> Option<Self> {
|
||||||
N::compute(m)
|
N::compute(m)
|
||||||
}
|
}
|
||||||
|
@ -87,7 +91,7 @@ macro_rules! svd_impl(
|
||||||
ldvt as i32, &mut work, lwork, &mut iwork, &mut info);
|
ldvt as i32, &mut work, lwork, &mut iwork, &mut info);
|
||||||
lapack_check!(info);
|
lapack_check!(info);
|
||||||
|
|
||||||
Some(SVD { u: u, s: s, vt: vt })
|
Some(SVD { u: u, singular_values: s, vt: vt })
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -108,7 +112,7 @@ macro_rules! svd_impl(
|
||||||
/// Useful if some components (e.g. some singular values) of this decomposition have
|
/// Useful if some components (e.g. some singular values) of this decomposition have
|
||||||
/// been manually changed by the user.
|
/// been manually changed by the user.
|
||||||
#[inline]
|
#[inline]
|
||||||
pub fn matrix(&self) -> MatrixMN<$t, R, C> {
|
pub fn recompose(self) -> MatrixMN<$t, R, C> {
|
||||||
let nrows = self.u.data.shape().0;
|
let nrows = self.u.data.shape().0;
|
||||||
let ncols = self.vt.data.shape().1;
|
let ncols = self.vt.data.shape().1;
|
||||||
let min_nrows_ncols = nrows.min(ncols);
|
let min_nrows_ncols = nrows.min(ncols);
|
||||||
|
@ -120,13 +124,13 @@ macro_rules! svd_impl(
|
||||||
sres.copy_from(&self.vt.rows_generic(0, min_nrows_ncols));
|
sres.copy_from(&self.vt.rows_generic(0, min_nrows_ncols));
|
||||||
|
|
||||||
for i in 0 .. min_nrows_ncols.value() {
|
for i in 0 .. min_nrows_ncols.value() {
|
||||||
let eigval = self.s[i];
|
let eigval = self.singular_values[i];
|
||||||
let mut row = sres.row_mut(i);
|
let mut row = sres.row_mut(i);
|
||||||
row *= eigval;
|
row *= eigval;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
&self.u * res
|
self.u * res
|
||||||
}
|
}
|
||||||
|
|
||||||
/// Computes the pseudo-inverse of the decomposed matrix.
|
/// Computes the pseudo-inverse of the decomposed matrix.
|
||||||
|
@ -145,7 +149,7 @@ macro_rules! svd_impl(
|
||||||
self.u.columns_generic(0, min_nrows_ncols).transpose_to(&mut sres);
|
self.u.columns_generic(0, min_nrows_ncols).transpose_to(&mut sres);
|
||||||
|
|
||||||
for i in 0 .. min_nrows_ncols.value() {
|
for i in 0 .. min_nrows_ncols.value() {
|
||||||
let eigval = self.s[i];
|
let eigval = self.singular_values[i];
|
||||||
let mut row = sres.row_mut(i);
|
let mut row = sres.row_mut(i);
|
||||||
|
|
||||||
if eigval.abs() > epsilon {
|
if eigval.abs() > epsilon {
|
||||||
|
@ -168,7 +172,7 @@ macro_rules! svd_impl(
|
||||||
pub fn rank(&self, epsilon: $t) -> usize {
|
pub fn rank(&self, epsilon: $t) -> usize {
|
||||||
let mut i = 0;
|
let mut i = 0;
|
||||||
|
|
||||||
for e in self.s.as_slice().iter() {
|
for e in self.singular_values.as_slice().iter() {
|
||||||
if e.abs() > epsilon {
|
if e.abs() > epsilon {
|
||||||
i += 1;
|
i += 1;
|
||||||
}
|
}
|
||||||
|
|
|
@ -15,8 +15,11 @@ use lapack::fortran as interface;
|
||||||
pub struct SymmetricEigen<N: Scalar, D: Dim>
|
pub struct SymmetricEigen<N: Scalar, D: Dim>
|
||||||
where DefaultAllocator: Allocator<N, D> +
|
where DefaultAllocator: Allocator<N, D> +
|
||||||
Allocator<N, D, D> {
|
Allocator<N, D, D> {
|
||||||
pub eigenvalues: VectorN<N, D>,
|
/// The eigenvectors of the decomposed matrix.
|
||||||
pub eigenvectors: MatrixN<N, D>,
|
pub eigenvectors: MatrixN<N, D>,
|
||||||
|
|
||||||
|
/// The unsorted eigenvalues of the decomposed matrix.
|
||||||
|
pub eigenvalues: VectorN<N, D>,
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
@ -48,7 +51,7 @@ impl<N: SymmetricEigenScalar + Real, D: Dim> SymmetricEigen<N, D>
|
||||||
|
|
||||||
let jobz = if eigenvectors { b'V' } else { b'N' };
|
let jobz = if eigenvectors { b'V' } else { b'N' };
|
||||||
|
|
||||||
let (nrows, ncols) = m.data.shape();
|
let nrows = m.data.shape().0;
|
||||||
let n = nrows.value();
|
let n = nrows.value();
|
||||||
|
|
||||||
let lda = n as i32;
|
let lda = n as i32;
|
||||||
|
@ -71,14 +74,14 @@ impl<N: SymmetricEigenScalar + Real, D: Dim> SymmetricEigen<N, D>
|
||||||
/// Computes only the eigenvalues of the input matrix.
|
/// Computes only the eigenvalues of the input matrix.
|
||||||
///
|
///
|
||||||
/// Panics if the method does not converge.
|
/// Panics if the method does not converge.
|
||||||
pub fn eigenvalues(mut m: MatrixN<N, D>) -> VectorN<N, D> {
|
pub fn eigenvalues(m: MatrixN<N, D>) -> VectorN<N, D> {
|
||||||
Self::do_decompose(m, false).expect("SymmetricEigen eigenvalues: convergence failure.").0
|
Self::do_decompose(m, false).expect("SymmetricEigen eigenvalues: convergence failure.").0
|
||||||
}
|
}
|
||||||
|
|
||||||
/// Computes only the eigenvalues of the input matrix.
|
/// Computes only the eigenvalues of the input matrix.
|
||||||
///
|
///
|
||||||
/// Returns `None` if the method does not converge.
|
/// Returns `None` if the method does not converge.
|
||||||
pub fn try_eigenvalues(mut m: MatrixN<N, D>) -> Option<VectorN<N, D>> {
|
pub fn try_eigenvalues(m: MatrixN<N, D>) -> Option<VectorN<N, D>> {
|
||||||
Self::do_decompose(m, false).map(|res| res.0)
|
Self::do_decompose(m, false).map(|res| res.0)
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -113,9 +116,13 @@ impl<N: SymmetricEigenScalar + Real, D: Dim> SymmetricEigen<N, D>
|
||||||
* Lapack functions dispatch.
|
* Lapack functions dispatch.
|
||||||
*
|
*
|
||||||
*/
|
*/
|
||||||
|
/// Trait implemented by scalars for which Lapack implements the eigendecomposition of symmetric
|
||||||
|
/// real matrices.
|
||||||
pub trait SymmetricEigenScalar: Scalar {
|
pub trait SymmetricEigenScalar: Scalar {
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xsyev(jobz: u8, uplo: u8, n: i32, a: &mut [Self], lda: i32, w: &mut [Self], work: &mut [Self],
|
fn xsyev(jobz: u8, uplo: u8, n: i32, a: &mut [Self], lda: i32, w: &mut [Self], work: &mut [Self],
|
||||||
lwork: i32, info: &mut i32);
|
lwork: i32, info: &mut i32);
|
||||||
|
#[allow(missing_docs)]
|
||||||
fn xsyev_work_size(jobz: u8, uplo: u8, n: i32, a: &mut [Self], lda: i32, info: &mut i32) -> i32;
|
fn xsyev_work_size(jobz: u8, uplo: u8, n: i32, a: &mut [Self], lda: i32, info: &mut i32) -> i32;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
Loading…
Reference in New Issue