Skip to content

Commit 2167299

Browse files
committed
Added svd, svd_inplace functions
1 parent d369255 commit 2167299

File tree

2 files changed

+84
-4
lines changed

2 files changed

+84
-4
lines changed

src/lapack/mod.rs

Lines changed: 82 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -12,6 +12,8 @@ type AfArray = self::libc::c_longlong;
1212

1313
#[allow(dead_code)]
1414
extern {
15+
fn af_svd(u: MutAfArray, s: MutAfArray, vt: MutAfArray, input: AfArray) -> c_int;
16+
fn af_svd_inplace(u: MutAfArray, s: MutAfArray, vt: MutAfArray, input: AfArray) -> c_int;
1517
fn af_lu(lower: MutAfArray, upper: MutAfArray, pivot: MutAfArray, input: AfArray) -> c_int;
1618
fn af_lu_inplace(pivot: MutAfArray, input: AfArray, is_lapack_piv: c_int) -> c_int;
1719
fn af_qr(q: MutAfArray, r: MutAfArray, tau: MutAfArray, input: AfArray) -> c_int;
@@ -26,6 +28,83 @@ extern {
2628
fn af_norm(out: MutDouble, input: AfArray, ntype: uint8_t, p: c_double, q: c_double) -> c_int;
2729
}
2830

31+
/// Perform Singular Value Decomposition
32+
///
33+
/// This function factorizes a matrix A into two unitary matrices U and Vt, and a diagonal matrix S
34+
/// such that
35+
///
36+
/// A = U∗S∗Vt
37+
///
38+
/// If A has M rows and N columns, U is of the size M x M , V is of size N x N, and S is of size M
39+
/// x N
40+
///
41+
/// # Parameters
42+
///
43+
/// - `in` is the input matrix
44+
///
45+
/// # Return Values
46+
///
47+
/// A triplet of Arrays.
48+
///
49+
/// The first Array is the output array containing U
50+
///
51+
/// The second Array is the output array containing the diagonal values of sigma, (singular values of the input matrix))
52+
///
53+
/// The third Array is the output array containing V ^ H
54+
#[allow(unused_mut)]
55+
pub fn svd(input: &Array) -> Result<(Array, Array, Array), AfError> {
56+
unsafe {
57+
let mut u: i64 = 0;
58+
let mut s: i64 = 0;
59+
let mut vt: i64 = 0;
60+
let err_val = af_svd(&mut u as MutAfArray, &mut s as MutAfArray, &mut vt as MutAfArray,
61+
input.get() as AfArray);
62+
match err_val {
63+
0 => Ok((Array::from(u), Array::from(s), Array::from(vt))),
64+
_ => Err(AfError::from(err_val)),
65+
}
66+
}
67+
}
68+
69+
/// Perform Singular Value Decomposition inplace
70+
///
71+
/// This function factorizes a matrix A into two unitary matrices U and Vt, and a diagonal matrix S
72+
/// such that
73+
///
74+
/// A = U∗S∗Vt
75+
///
76+
/// If A has M rows and N columns, U is of the size M x M , V is of size N x N, and S is of size M
77+
/// x N
78+
///
79+
/// # Parameters
80+
///
81+
/// - `in` is the input/output matrix. This will contain random data after the function call is
82+
/// complete.
83+
///
84+
/// # Return Values
85+
///
86+
/// A triplet of Arrays.
87+
///
88+
/// The first Array is the output array containing U
89+
///
90+
/// The second Array is the output array containing the diagonal values of sigma, (singular values of the input matrix))
91+
///
92+
/// The third Array is the output array containing V ^ H
93+
#[allow(unused_mut)]
94+
pub fn svd_inplace(input: &mut Array) -> Result<(Array, Array, Array), AfError> {
95+
unsafe {
96+
let mut u: i64 = 0;
97+
let mut s: i64 = 0;
98+
let mut vt: i64 = 0;
99+
let err_val = af_svd_inplace(&mut u as MutAfArray, &mut s as MutAfArray,
100+
&mut vt as MutAfArray, input.get() as AfArray);
101+
match err_val {
102+
0 => Ok((Array::from(u), Array::from(s), Array::from(vt))),
103+
_ => Err(AfError::from(err_val)),
104+
}
105+
}
106+
}
107+
29108
/// Perform LU decomposition
30109
///
31110
/// # Parameters
@@ -60,7 +139,7 @@ pub fn lu(input: &Array) -> Result<(Array, Array, Array), AfError> {
60139
///
61140
/// # Parameters
62141
///
63-
/// - `input` is the input matrix
142+
/// - `input` contains the input matrix on entry and packed LU decomposition on exit
64143
/// - `is_lapack_pic` specified if the pivot is returned in original LAPACK compliant format
65144
///
66145
/// # Return Values
@@ -115,7 +194,7 @@ pub fn qr(input: &Array) -> Result<(Array, Array, Array), AfError> {
115194
///
116195
/// # Parameters
117196
///
118-
/// - `input` is the input matrix
197+
/// - `input` contains the input matrix on entry, and packed QR decomposition on exit
119198
///
120199
/// # Return Values
121200
///
@@ -165,7 +244,7 @@ pub fn cholesky(input: &Array, is_upper: bool) -> Result<(Array, i32), AfError>
165244
///
166245
/// # Parameters
167246
///
168-
/// - `input` is the input matrix
247+
/// - `input` contains the input matrix on entry, and triangular matrix on exit.
169248
/// - `is_upper` is a boolean to indicate if the output has to be upper or lower triangular matrix
170249
///
171250
/// # Return Values

src/lib.rs

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -65,7 +65,8 @@ pub use image::{bilateral, mean_shift, medfilt, sobel};
6565
pub use image::{unwrap, wrap, sat, rgb2ycbcr, ycbcr2rgb};
6666
mod image;
6767

68-
pub use lapack::{lu, lu_inplace, qr, qr_inplace, cholesky, cholesky_inplace, solve, solve_lu, inverse, det, rank, norm};
68+
pub use lapack::{svd, lu, qr, cholesky, solve, solve_lu, inverse, det, rank, norm};
69+
pub use lapack::{svd_inplace, lu_inplace, qr_inplace, cholesky_inplace};
6970
mod lapack;
7071

7172
pub use signal::{approx1, approx2};

0 commit comments

Comments
 (0)