Fix all tests and the ColPivQR::solve.

This commit is contained in:
Crozet Sébastien 2021-02-25 12:06:04 +01:00
parent 63a34528e0
commit 308d95386e
6 changed files with 198 additions and 174 deletions

View File

@ -1,98 +1,94 @@
use num::Zero;
#[cfg(feature = "serde-serialize")] #[cfg(feature = "serde-serialize")]
use serde::{Deserialize, Serialize}; use serde::{Deserialize, Serialize};
use num::Zero;
use alga::general::ComplexField;
use crate::allocator::{Allocator, Reallocator}; use crate::allocator::{Allocator, Reallocator};
use crate::base::{DefaultAllocator, Matrix, MatrixMN, MatrixN, Unit, VectorN}; use crate::base::{DefaultAllocator, Matrix, MatrixMN, MatrixN, Unit, VectorN};
use crate::constraint::{SameNumberOfRows, ShapeConstraint}; use crate::constraint::{SameNumberOfRows, ShapeConstraint};
use crate::dimension::{Dim, DimMin, DimMinimum, U1}; use crate::dimension::{Dim, DimMin, DimMinimum, U1};
use crate::storage::{Storage, StorageMut}; use crate::storage::{Storage, StorageMut};
use crate::ComplexField;
use crate::geometry::Reflection; use crate::geometry::Reflection;
use crate::linalg::householder; use crate::linalg::{householder, PermutationSequence};
use crate::linalg::PermutationSequence;
/// The QRP decomposition of a general matrix. /// The QR decomposition (with column pivoting) of a general matrix.
#[cfg_attr(feature = "serde-serialize", derive(Serialize, Deserialize))] #[cfg_attr(feature = "serde-serialize", derive(Serialize, Deserialize))]
#[cfg_attr( #[cfg_attr(
feature = "serde-serialize", feature = "serde-serialize",
serde(bound( serde(bound(serialize = "DefaultAllocator: Allocator<N, R, C> +
serialize = "DefaultAllocator: Allocator<N, R, C> +
Allocator<N, DimMinimum<R, C>>, Allocator<N, DimMinimum<R, C>>,
MatrixMN<N, R, C>: Serialize, MatrixMN<N, R, C>: Serialize,
PermutationSequence<DimMinimum<R, C>>: Serialize, PermutationSequence<DimMinimum<R, C>>: Serialize,
VectorN<N, DimMinimum<R, C>>: Serialize" VectorN<N, DimMinimum<R, C>>: Serialize"))
))
)] )]
#[cfg_attr( #[cfg_attr(
feature = "serde-serialize", feature = "serde-serialize",
serde(bound( serde(bound(deserialize = "DefaultAllocator: Allocator<N, R, C> +
deserialize = "DefaultAllocator: Allocator<N, R, C> +
Allocator<N, DimMinimum<R, C>>, Allocator<N, DimMinimum<R, C>>,
MatrixMN<N, R, C>: Deserialize<'de>, MatrixMN<N, R, C>: Deserialize<'de>,
PermutationSequence<DimMinimum<R, C>>: Deserialize<'de>, PermutationSequence<DimMinimum<R, C>>: Deserialize<'de>,
VectorN<N, DimMinimum<R, C>>: Deserialize<'de>" VectorN<N, DimMinimum<R, C>>: Deserialize<'de>"))
))
)] )]
#[derive(Clone, Debug)] #[derive(Clone, Debug)]
pub struct QRP<N: ComplexField, R: DimMin<C>, C: Dim> pub struct ColPivQR<N: ComplexField, R: DimMin<C>, C: Dim>
where DefaultAllocator: Allocator<N, R, C> + Allocator<N, DimMinimum<R, C>> + where
Allocator<(usize, usize), DimMinimum<R, C>>, DefaultAllocator: Allocator<N, R, C>
+ Allocator<N, DimMinimum<R, C>>
+ Allocator<(usize, usize), DimMinimum<R, C>>,
{ {
qrp: MatrixMN<N, R, C>, col_piv_qr: MatrixMN<N, R, C>,
p: PermutationSequence<DimMinimum<R, C>>, p: PermutationSequence<DimMinimum<R, C>>,
diag: VectorN<N, DimMinimum<R, C>>, diag: VectorN<N, DimMinimum<R, C>>,
} }
impl<N: ComplexField, R: DimMin<C>, C: Dim> Copy for QRP<N, R, C> impl<N: ComplexField, R: DimMin<C>, C: Dim> Copy for ColPivQR<N, R, C>
where where
DefaultAllocator: Allocator<N, R, C> + Allocator<N, DimMinimum<R, C>> + DefaultAllocator: Allocator<N, R, C>
Allocator<(usize, usize), DimMinimum<R, C>>, + Allocator<N, DimMinimum<R, C>>
+ Allocator<(usize, usize), DimMinimum<R, C>>,
MatrixMN<N, R, C>: Copy, MatrixMN<N, R, C>: Copy,
PermutationSequence<DimMinimum<R, C>>: Copy, PermutationSequence<DimMinimum<R, C>>: Copy,
VectorN<N, DimMinimum<R, C>>: Copy, VectorN<N, DimMinimum<R, C>>: Copy,
{}
impl<N: ComplexField, R: DimMin<C>, C: Dim> QRP<N, R, C>
where DefaultAllocator: Allocator<N, R, C> + Allocator<N, R> + Allocator<N, DimMinimum<R, C>> + Allocator<(usize, usize), DimMinimum<R, C>>
{ {
/// Computes the QRP decomposition using householder reflections. }
impl<N: ComplexField, R: DimMin<C>, C: Dim> ColPivQR<N, R, C>
where
DefaultAllocator: Allocator<N, R, C>
+ Allocator<N, R>
+ Allocator<N, DimMinimum<R, C>>
+ Allocator<(usize, usize), DimMinimum<R, C>>,
{
/// Computes the ColPivQR decomposition using householder reflections.
pub fn new(mut matrix: MatrixMN<N, R, C>) -> Self { pub fn new(mut matrix: MatrixMN<N, R, C>) -> Self {
let (nrows, ncols) = matrix.data.shape(); let (nrows, ncols) = matrix.data.shape();
let min_nrows_ncols = nrows.min(ncols); let min_nrows_ncols = nrows.min(ncols);
let mut p = PermutationSequence::identity_generic(min_nrows_ncols); let mut p = PermutationSequence::identity_generic(min_nrows_ncols);
let mut diag = unsafe { MatrixMN::new_uninitialized_generic(min_nrows_ncols, U1) }; let mut diag = unsafe { MatrixMN::new_uninitialized_generic(min_nrows_ncols, U1) };
println!("diag: {:?}", &diag);
if min_nrows_ncols.value() == 0 { if min_nrows_ncols.value() == 0 {
return QRP { return ColPivQR {
qrp: matrix, col_piv_qr: matrix,
p: p, p,
diag: diag, diag,
}; };
} }
for ite in 0..min_nrows_ncols.value() { for i in 0..min_nrows_ncols.value() {
let mut col_norm = Vec::new(); let piv = matrix.slice_range(i.., i..).icamax_full();
for column in matrix.column_iter() { let col_piv = piv.1 + i;
col_norm.push(column.norm_squared()); matrix.swap_columns(i, col_piv);
} p.append_permutation(i, col_piv);
let piv = matrix.slice_range(ite.., ite..).icamax_full(); householder::clear_column_unchecked(&mut matrix, &mut diag[i], i, 0, None);
let col_piv = piv.1 + ite;
matrix.swap_columns(ite, col_piv);
p.append_permutation(ite, col_piv);
householder::clear_column_unchecked(&mut matrix, &mut diag[ite], ite, 0, None);
println!("matrix: {:?}", &matrix.data);
} }
QRP { ColPivQR {
qrp: matrix, col_piv_qr: matrix,
p: p, p,
diag: diag, diag,
} }
} }
@ -102,8 +98,11 @@ where DefaultAllocator: Allocator<N, R, C> + Allocator<N, R> + Allocator<N, DimM
where where
DefaultAllocator: Allocator<N, DimMinimum<R, C>, C>, DefaultAllocator: Allocator<N, DimMinimum<R, C>, C>,
{ {
let (nrows, ncols) = self.qrp.data.shape(); let (nrows, ncols) = self.col_piv_qr.data.shape();
let mut res = self.qrp.rows_generic(0, nrows.min(ncols)).upper_triangle(); let mut res = self
.col_piv_qr
.rows_generic(0, nrows.min(ncols))
.upper_triangle();
res.set_partial_diagonal(self.diag.iter().map(|e| N::from_real(e.modulus()))); res.set_partial_diagonal(self.diag.iter().map(|e| N::from_real(e.modulus())));
res res
} }
@ -116,8 +115,10 @@ where DefaultAllocator: Allocator<N, R, C> + Allocator<N, R> + Allocator<N, DimM
where where
DefaultAllocator: Reallocator<N, R, C, DimMinimum<R, C>, C>, DefaultAllocator: Reallocator<N, R, C, DimMinimum<R, C>, C>,
{ {
let (nrows, ncols) = self.qrp.data.shape(); let (nrows, ncols) = self.col_piv_qr.data.shape();
let mut res = self.qrp.resize_generic(nrows.min(ncols), ncols, N::zero()); let mut res = self
.col_piv_qr
.resize_generic(nrows.min(ncols), ncols, N::zero());
res.fill_lower_triangle(N::zero(), 1); res.fill_lower_triangle(N::zero(), 1);
res.set_partial_diagonal(self.diag.iter().map(|e| N::from_real(e.modulus()))); res.set_partial_diagonal(self.diag.iter().map(|e| N::from_real(e.modulus())));
res res
@ -125,8 +126,10 @@ where DefaultAllocator: Allocator<N, R, C> + Allocator<N, R> + Allocator<N, DimM
/// Computes the orthogonal matrix `Q` of this decomposition. /// Computes the orthogonal matrix `Q` of this decomposition.
pub fn q(&self) -> MatrixMN<N, R, DimMinimum<R, C>> pub fn q(&self) -> MatrixMN<N, R, DimMinimum<R, C>>
where DefaultAllocator: Allocator<N, R, DimMinimum<R, C>> { where
let (nrows, ncols) = self.qrp.data.shape(); DefaultAllocator: Allocator<N, R, DimMinimum<R, C>>,
{
let (nrows, ncols) = self.col_piv_qr.data.shape();
// NOTE: we could build the identity matrix and call q_mul on it. // NOTE: we could build the identity matrix and call q_mul on it.
// Instead we don't so that we take in account the matrix sparseness. // Instead we don't so that we take in account the matrix sparseness.
@ -134,8 +137,8 @@ where DefaultAllocator: Allocator<N, R, C> + Allocator<N, R> + Allocator<N, DimM
let dim = self.diag.len(); let dim = self.diag.len();
for i in (0..dim).rev() { for i in (0..dim).rev() {
let axis = self.qrp.slice_range(i.., i); let axis = self.col_piv_qr.slice_range(i.., i);
// FIXME: sometimes, the axis might have a zero magnitude. // TODO: sometimes, the axis might have a zero magnitude.
let refl = Reflection::new(Unit::new_unchecked(axis), N::zero()); let refl = Reflection::new(Unit::new_unchecked(axis), N::zero());
let mut res_rows = res.slice_range_mut(i.., i..); let mut res_rows = res.slice_range_mut(i.., i..);
@ -160,25 +163,27 @@ where DefaultAllocator: Allocator<N, R, C> + Allocator<N, R> + Allocator<N, DimM
) )
where where
DimMinimum<R, C>: DimMin<C, Output = DimMinimum<R, C>>, DimMinimum<R, C>: DimMin<C, Output = DimMinimum<R, C>>,
DefaultAllocator: DefaultAllocator: Allocator<N, R, DimMinimum<R, C>>
Allocator<N, R, DimMinimum<R, C>> + Reallocator<N, R, C, DimMinimum<R, C>, C> + Allocator<(usize, usize), DimMinimum<R, C>>, + Reallocator<N, R, C, DimMinimum<R, C>, C>
+ Allocator<(usize, usize), DimMinimum<R, C>>,
{ {
(self.q(), self.r(), self.p) (self.q(), self.r(), self.p)
} }
#[doc(hidden)] #[doc(hidden)]
pub fn qrp_internal(&self) -> &MatrixMN<N, R, C> { pub fn col_piv_qr_internal(&self) -> &MatrixMN<N, R, C> {
&self.qrp &self.col_piv_qr
} }
/// Multiplies the provided matrix by the transpose of the `Q` matrix of this decomposition. /// Multiplies the provided matrix by the transpose of the `Q` matrix of this decomposition.
pub fn q_tr_mul<R2: Dim, C2: Dim, S2>(&self, rhs: &mut Matrix<N, R2, C2, S2>) pub fn q_tr_mul<R2: Dim, C2: Dim, S2>(&self, rhs: &mut Matrix<N, R2, C2, S2>)
// FIXME: do we need a static constraint on the number of rows of rhs? where
where S2: StorageMut<N, R2, C2> { S2: StorageMut<N, R2, C2>,
{
let dim = self.diag.len(); let dim = self.diag.len();
for i in 0..dim { for i in 0..dim {
let axis = self.qrp.slice_range(i.., i); let axis = self.col_piv_qr.slice_range(i.., i);
let refl = Reflection::new(Unit::new_unchecked(axis), N::zero()); let refl = Reflection::new(Unit::new_unchecked(axis), N::zero());
let mut rhs_rows = rhs.rows_range_mut(i..); let mut rhs_rows = rhs.rows_range_mut(i..);
@ -187,9 +192,10 @@ where DefaultAllocator: Allocator<N, R, C> + Allocator<N, R> + Allocator<N, DimM
} }
} }
impl<N: ComplexField, D: DimMin<D, Output = D>> QRP<N, D, D> impl<N: ComplexField, D: DimMin<D, Output = D>> ColPivQR<N, D, D>
where DefaultAllocator: Allocator<N, D, D> + Allocator<N, D> + where
Allocator<(usize, usize), DimMinimum<D, D>> DefaultAllocator:
Allocator<N, D, D> + Allocator<N, D> + Allocator<(usize, usize), DimMinimum<D, D>>,
{ {
/// Solves the linear system `self * x = b`, where `x` is the unknown to be determined. /// Solves the linear system `self * x = b`, where `x` is the unknown to be determined.
/// ///
@ -222,20 +228,23 @@ where DefaultAllocator: Allocator<N, D, D> + Allocator<N, D> +
ShapeConstraint: SameNumberOfRows<R2, D>, ShapeConstraint: SameNumberOfRows<R2, D>,
{ {
assert_eq!( assert_eq!(
self.qrp.nrows(), self.col_piv_qr.nrows(),
b.nrows(), b.nrows(),
"QRP solve matrix dimension mismatch." "ColPivQR solve matrix dimension mismatch."
); );
assert!( assert!(
self.qrp.is_square(), self.col_piv_qr.is_square(),
"QRP solve: unable to solve a non-square system." "ColPivQR solve: unable to solve a non-square system."
); );
self.q_tr_mul(b); self.q_tr_mul(b);
self.solve_upper_triangular_mut(b) let solved = self.solve_upper_triangular_mut(b);
self.p.inv_permute_rows(b);
solved
} }
// FIXME: duplicate code from the `solve` module. // TODO: duplicate code from the `solve` module.
fn solve_upper_triangular_mut<R2: Dim, C2: Dim, S2>( fn solve_upper_triangular_mut<R2: Dim, C2: Dim, S2>(
&self, &self,
b: &mut Matrix<N, R2, C2, S2>, b: &mut Matrix<N, R2, C2, S2>,
@ -244,7 +253,7 @@ where DefaultAllocator: Allocator<N, D, D> + Allocator<N, D> +
S2: StorageMut<N, R2, C2>, S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>, ShapeConstraint: SameNumberOfRows<R2, D>,
{ {
let dim = self.qrp.nrows(); let dim = self.col_piv_qr.nrows();
for k in 0..b.ncols() { for k in 0..b.ncols() {
let mut b = b.column_mut(k); let mut b = b.column_mut(k);
@ -263,7 +272,7 @@ where DefaultAllocator: Allocator<N, D, D> + Allocator<N, D> +
} }
b.rows_range_mut(..i) b.rows_range_mut(..i)
.axpy(-coeff, &self.qrp.slice_range(..i, i), N::one()); .axpy(-coeff, &self.col_piv_qr.slice_range(..i, i), N::one());
} }
} }
@ -275,12 +284,12 @@ where DefaultAllocator: Allocator<N, D, D> + Allocator<N, D> +
/// Returns `None` if the decomposed matrix is not invertible. /// Returns `None` if the decomposed matrix is not invertible.
pub fn try_inverse(&self) -> Option<MatrixN<N, D>> { pub fn try_inverse(&self) -> Option<MatrixN<N, D>> {
assert!( assert!(
self.qrp.is_square(), self.col_piv_qr.is_square(),
"QRP inverse: unable to compute the inverse of a non-square matrix." "ColPivQR inverse: unable to compute the inverse of a non-square matrix."
); );
// FIXME: is there a less naive method ? // TODO: is there a less naive method ?
let (nrows, ncols) = self.qrp.data.shape(); let (nrows, ncols) = self.col_piv_qr.data.shape();
let mut res = MatrixN::identity_generic(nrows, ncols); let mut res = MatrixN::identity_generic(nrows, ncols);
if self.solve_mut(&mut res) { if self.solve_mut(&mut res) {
@ -293,8 +302,8 @@ where DefaultAllocator: Allocator<N, D, D> + Allocator<N, D> +
/// Indicates if the decomposed matrix is invertible. /// Indicates if the decomposed matrix is invertible.
pub fn is_invertible(&self) -> bool { pub fn is_invertible(&self) -> bool {
assert!( assert!(
self.qrp.is_square(), self.col_piv_qr.is_square(),
"QRP: unable to test the invertibility of a non-square matrix." "ColPivQR: unable to test the invertibility of a non-square matrix."
); );
for i in 0..self.diag.len() { for i in 0..self.diag.len() {
@ -306,25 +315,32 @@ where DefaultAllocator: Allocator<N, D, D> + Allocator<N, D> +
true true
} }
/// Computes the determinant of the decomposed matrix. /// Computes the determinant of the decomposed matrix.
pub fn determinant(&self) -> N { pub fn determinant(&self) -> N {
let dim = self.qrp.nrows(); let dim = self.col_piv_qr.nrows();
assert!(self.qrp.is_square(), "QRP determinant: unable to compute the determinant of a non-square matrix."); assert!(
self.col_piv_qr.is_square(),
"ColPivQR determinant: unable to compute the determinant of a non-square matrix."
);
let mut res = N::one(); let mut res = N::one();
for i in 0 .. dim { for i in 0..dim {
res *= unsafe { *self.diag.vget_unchecked(i) }; res *= unsafe { *self.diag.vget_unchecked(i) };
} }
res * self.p.determinant() res * self.p.determinant()
} }
} }
impl<N: ComplexField, R: DimMin<C>, C: Dim, S: Storage<N, R, C>> Matrix<N, R, C, S> impl<N: ComplexField, R: DimMin<C>, C: Dim, S: Storage<N, R, C>> Matrix<N, R, C, S>
where DefaultAllocator: Allocator<N, R, C> + Allocator<N, R> + Allocator<N, DimMinimum<R, C>> + Allocator<(usize, usize), DimMinimum<R, C>> where
DefaultAllocator: Allocator<N, R, C>
+ Allocator<N, R>
+ Allocator<N, DimMinimum<R, C>>
+ Allocator<(usize, usize), DimMinimum<R, C>>,
{ {
/// Computes the QRP decomposition of this matrix. /// Computes the QR decomposition (with column pivoting) of this matrix.
pub fn qrp(self) -> QRP<N, R, C> { pub fn col_piv_qr(self) -> ColPivQR<N, R, C> {
QRP::new(self.into_owned()) ColPivQR::new(self.into_owned())
} }
} }

View File

@ -34,7 +34,7 @@ pub fn reflection_axis_mut<N: ComplexField, D: Dim, S: StorageMut<N, D>>(
if !factor.is_zero() { if !factor.is_zero() {
column.unscale_mut(factor.sqrt()); column.unscale_mut(factor.sqrt());
(signed_norm, true) (-signed_norm, true)
} else { } else {
// TODO: not sure why we don't have a - sign here. // TODO: not sure why we don't have a - sign here.
(signed_norm, false) (signed_norm, false)

View File

@ -19,6 +19,7 @@ mod inverse;
mod lu; mod lu;
mod permutation_sequence; mod permutation_sequence;
mod qr; mod qr;
mod col_piv_qr;
mod schur; mod schur;
mod solve; mod solve;
mod svd; mod svd;
@ -39,6 +40,7 @@ pub use self::hessenberg::*;
pub use self::lu::*; pub use self::lu::*;
pub use self::permutation_sequence::*; pub use self::permutation_sequence::*;
pub use self::qr::*; pub use self::qr::*;
pub use self::col_piv_qr::*;
pub use self::schur::*; pub use self::schur::*;
pub use self::svd::*; pub use self::svd::*;
pub use self::symmetric_eigen::*; pub use self::symmetric_eigen::*;

View File

@ -60,8 +60,8 @@ where
return QR { qr: matrix, diag }; return QR { qr: matrix, diag };
} }
for ite in 0..min_nrows_ncols.value() { for i in 0..min_nrows_ncols.value() {
householder::clear_column_unchecked(&mut matrix, &mut diag[ite], ite, 0, None); householder::clear_column_unchecked(&mut matrix, &mut diag[i], i, 0, None);
} }
QR { qr: matrix, diag } QR { qr: matrix, diag }

View File

@ -3,19 +3,21 @@
use na::Matrix4; use na::Matrix4;
#[test] #[test]
fn qrp() { fn col_piv_qr() {
let m = Matrix4::new ( let m = Matrix4::new(
1.0, -1.0, 2.0, 1.0, 1.0, -1.0, 2.0, 1.0, -1.0, 3.0, -1.0, -1.0, 3.0, -5.0, 5.0, 3.0, 1.0, 2.0, 1.0, -2.0,
-1.0, 3.0, -1.0, -1.0, );
3.0, -5.0, 5.0, 3.0, let col_piv_qr = m.col_piv_qr();
1.0, 2.0, 1.0, -2.0); assert!(relative_eq!(
let qrp = m.qrp(); col_piv_qr.determinant(),
assert!(relative_eq!(qrp.determinant(), 0.0, epsilon = 1.0e-7)); 0.0,
epsilon = 1.0e-7
));
let (q, r, p) = qrp.unpack(); let (q, r, p) = col_piv_qr.unpack();
let mut qr = q * r; let mut qr = q * r;
p.inv_permute_columns(& mut qr); p.inv_permute_columns(&mut qr);
assert!(relative_eq!(m, qr, epsilon = 1.0e-7)); assert!(relative_eq!(m, qr, epsilon = 1.0e-7));
} }
@ -31,64 +33,68 @@ mod quickcheck_tests {
use crate::core::helper::{RandScalar, RandComplex}; use crate::core::helper::{RandScalar, RandComplex};
quickcheck! { quickcheck! {
fn qrp(m: DMatrix<$scalar>) -> bool { fn col_piv_qr(m: DMatrix<$scalar>) -> bool {
let m = m.map(|e| e.0); let m = m.map(|e| e.0);
let qrp = m.clone().qrp(); let col_piv_qr = m.clone().col_piv_qr();
let q = qrp.q(); let (q, r, p) = col_piv_qr.unpack();
let r = qrp.r(); let mut qr = &q * &r;
p.inv_permute_columns(&mut qr);
println!("m: {}", m); println!("m: {}", m);
println!("qrp: {}", &q * &r); println!("col_piv_qr: {}", &q * &r);
relative_eq!(m, &q * r, epsilon = 1.0e-7) && relative_eq!(m, &qr, epsilon = 1.0e-7) &&
q.is_orthogonal(1.0e-7) q.is_orthogonal(1.0e-7)
} }
fn qrp_static_5_3(m: Matrix5x3<$scalar>) -> bool { fn col_piv_qr_static_5_3(m: Matrix5x3<$scalar>) -> bool {
let m = m.map(|e| e.0); let m = m.map(|e| e.0);
let qrp = m.qrp(); let col_piv_qr = m.col_piv_qr();
let q = qrp.q(); let (q, r, p) = col_piv_qr.unpack();
let r = qrp.r(); let mut qr = q * r;
p.inv_permute_columns(&mut qr);
relative_eq!(m, q * r, epsilon = 1.0e-7) && relative_eq!(m, qr, epsilon = 1.0e-7) &&
q.is_orthogonal(1.0e-7) q.is_orthogonal(1.0e-7)
} }
fn qrp_static_3_5(m: Matrix3x5<$scalar>) -> bool { fn col_piv_qr_static_3_5(m: Matrix3x5<$scalar>) -> bool {
let m = m.map(|e| e.0); let m = m.map(|e| e.0);
let qrp = m.qrp(); let col_piv_qr = m.col_piv_qr();
let q = qrp.q(); let (q, r, p) = col_piv_qr.unpack();
let r = qrp.r(); let mut qr = q * r;
p.inv_permute_columns(&mut qr);
relative_eq!(m, q * r, epsilon = 1.0e-7) && relative_eq!(m, qr, epsilon = 1.0e-7) &&
q.is_orthogonal(1.0e-7)
}
fn qrp_static_square(m: Matrix4<$scalar>) -> bool {
let m = m.map(|e| e.0);
let qrp = m.qrp();
let q = qrp.q();
let r = qrp.r();
println!("{}{}{}{}", q, r, q * r, m);
relative_eq!(m, q * r, epsilon = 1.0e-7) &&
q.is_orthogonal(1.0e-7) q.is_orthogonal(1.0e-7)
} }
fn qrp_solve(n: usize, nb: usize) -> bool { fn col_piv_qr_static_square(m: Matrix4<$scalar>) -> bool {
let m = m.map(|e| e.0);
let col_piv_qr = m.col_piv_qr();
let (q, r, p) = col_piv_qr.unpack();
let mut qr = q * r;
p.inv_permute_columns(&mut qr);
println!("{}{}{}{}", q, r, qr, m);
relative_eq!(m, qr, epsilon = 1.0e-7) &&
q.is_orthogonal(1.0e-7)
}
fn col_piv_qr_solve(n: usize, nb: usize) -> bool {
if n != 0 && nb != 0 { if n != 0 && nb != 0 {
let n = cmp::min(n, 50); // To avoid slowing down the test too much. let n = cmp::min(n, 50); // To avoid slowing down the test too much.
let nb = cmp::min(nb, 50); // To avoid slowing down the test too much. let nb = cmp::min(nb, 50); // To avoid slowing down the test too much.
let m = DMatrix::<$scalar>::new_random(n, n).map(|e| e.0); let m = DMatrix::<$scalar>::new_random(n, n).map(|e| e.0);
let qrp = m.clone().qrp(); let col_piv_qr = m.clone().col_piv_qr();
let b1 = DVector::<$scalar>::new_random(n).map(|e| e.0); let b1 = DVector::<$scalar>::new_random(n).map(|e| e.0);
let b2 = DMatrix::<$scalar>::new_random(n, nb).map(|e| e.0); let b2 = DMatrix::<$scalar>::new_random(n, nb).map(|e| e.0);
if qrp.is_invertible() { if col_piv_qr.is_invertible() {
let sol1 = qrp.solve(&b1).unwrap(); let sol1 = col_piv_qr.solve(&b1).unwrap();
let sol2 = qrp.solve(&b2).unwrap(); let sol2 = col_piv_qr.solve(&b2).unwrap();
return relative_eq!(&m * sol1, b1, epsilon = 1.0e-6) && return relative_eq!(&m * sol1, b1, epsilon = 1.0e-6) &&
relative_eq!(&m * sol2, b2, epsilon = 1.0e-6) relative_eq!(&m * sol2, b2, epsilon = 1.0e-6)
@ -98,15 +104,15 @@ mod quickcheck_tests {
return true; return true;
} }
fn qrp_solve_static(m: Matrix4<$scalar>) -> bool { fn col_piv_qr_solve_static(m: Matrix4<$scalar>) -> bool {
let m = m.map(|e| e.0); let m = m.map(|e| e.0);
let qrp = m.qrp(); let col_piv_qr = m.col_piv_qr();
let b1 = Vector4::<$scalar>::new_random().map(|e| e.0); let b1 = Vector4::<$scalar>::new_random().map(|e| e.0);
let b2 = Matrix4x3::<$scalar>::new_random().map(|e| e.0); let b2 = Matrix4x3::<$scalar>::new_random().map(|e| e.0);
if qrp.is_invertible() { if col_piv_qr.is_invertible() {
let sol1 = qrp.solve(&b1).unwrap(); let sol1 = col_piv_qr.solve(&b1).unwrap();
let sol2 = qrp.solve(&b2).unwrap(); let sol2 = col_piv_qr.solve(&b2).unwrap();
relative_eq!(m * sol1, b1, epsilon = 1.0e-6) && relative_eq!(m * sol1, b1, epsilon = 1.0e-6) &&
relative_eq!(m * sol2, b2, epsilon = 1.0e-6) relative_eq!(m * sol2, b2, epsilon = 1.0e-6)
@ -116,11 +122,11 @@ mod quickcheck_tests {
} }
} }
fn qrp_inverse(n: usize) -> bool { fn col_piv_qr_inverse(n: usize) -> bool {
let n = cmp::max(1, cmp::min(n, 15)); // To avoid slowing down the test too much. let n = cmp::max(1, cmp::min(n, 15)); // To avoid slowing down the test too much.
let m = DMatrix::<$scalar>::new_random(n, n).map(|e| e.0); let m = DMatrix::<$scalar>::new_random(n, n).map(|e| e.0);
if let Some(m1) = m.clone().qrp().try_inverse() { if let Some(m1) = m.clone().col_piv_qr().try_inverse() {
let id1 = &m * &m1; let id1 = &m * &m1;
let id2 = &m1 * &m; let id2 = &m1 * &m;
@ -131,11 +137,11 @@ mod quickcheck_tests {
} }
} }
fn qrp_inverse_static(m: Matrix4<$scalar>) -> bool { fn col_piv_qr_inverse_static(m: Matrix4<$scalar>) -> bool {
let m = m.map(|e| e.0); let m = m.map(|e| e.0);
let qrp = m.qrp(); let col_piv_qr = m.col_piv_qr();
if let Some(m1) = qrp.try_inverse() { if let Some(m1) = col_piv_qr.try_inverse() {
let id1 = &m * &m1; let id1 = &m * &m1;
let id2 = &m1 * &m; let id2 = &m1 * &m;
@ -153,4 +159,3 @@ mod quickcheck_tests {
gen_tests!(complex, RandComplex<f64>); gen_tests!(complex, RandComplex<f64>);
gen_tests!(f64, RandScalar<f64>); gen_tests!(f64, RandScalar<f64>);
} }

View File

@ -1,6 +1,7 @@
mod balancing; mod balancing;
mod bidiagonal; mod bidiagonal;
mod cholesky; mod cholesky;
mod col_piv_qr;
mod convolution; mod convolution;
mod eigen; mod eigen;
mod exp; mod exp;