Thanks to visit codestin.com
Credit goes to docs.rs

arrayfire/lapack/
mod.rs

1use super::core::{
2    af_array, AfError, Array, FloatingPoint, HasAfEnum, MatProp, NormType, HANDLE_ERROR,
3};
4
5use libc::{c_double, c_int, c_uint};
6
7extern "C" {
8    fn af_svd(u: *mut af_array, s: *mut af_array, vt: *mut af_array, input: af_array) -> c_int;
9    fn af_svd_inplace(
10        u: *mut af_array,
11        s: *mut af_array,
12        vt: *mut af_array,
13        input: af_array,
14    ) -> c_int;
15    fn af_lu(
16        lower: *mut af_array,
17        upper: *mut af_array,
18        pivot: *mut af_array,
19        input: af_array,
20    ) -> c_int;
21    fn af_lu_inplace(pivot: *mut af_array, input: af_array, is_lapack_piv: bool) -> c_int;
22    fn af_qr(q: *mut af_array, r: *mut af_array, tau: *mut af_array, input: af_array) -> c_int;
23    fn af_qr_inplace(tau: *mut af_array, input: af_array) -> c_int;
24    fn af_cholesky(out: *mut af_array, info: *mut c_int, input: af_array, is_upper: bool) -> c_int;
25    fn af_cholesky_inplace(info: *mut c_int, input: af_array, is_upper: bool) -> c_int;
26    fn af_solve(x: *mut af_array, a: af_array, b: af_array, options: c_uint) -> c_int;
27    fn af_solve_lu(
28        x: *mut af_array,
29        a: af_array,
30        piv: af_array,
31        b: af_array,
32        options: c_uint,
33    ) -> c_int;
34    fn af_inverse(out: *mut af_array, input: af_array, options: c_uint) -> c_int;
35    fn af_rank(rank: *mut c_uint, input: af_array, tol: c_double) -> c_int;
36    fn af_det(det_real: *mut c_double, det_imag: *mut c_double, input: af_array) -> c_int;
37    fn af_norm(
38        out: *mut c_double,
39        input: af_array,
40        ntype: c_uint,
41        p: c_double,
42        q: c_double,
43    ) -> c_int;
44    fn af_is_lapack_available(out: *mut bool) -> c_int;
45    fn af_pinverse(out: *mut af_array, input: af_array, tol: c_double, options: c_uint) -> c_int;
46}
47
48/// Perform Singular Value Decomposition
49///
50/// This function factorizes a matrix A into two unitary matrices U and Vt, and a diagonal matrix S
51/// such that
52///
53/// A = U∗S∗Vt
54///
55/// 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
56/// x N
57///
58/// # Parameters
59///
60/// - `in` is the input matrix
61///
62/// # Return Values
63///
64/// A triplet of Arrays.
65///
66/// The first Array is the output array containing U
67///
68/// The second Array is the output array containing the diagonal values of sigma, (singular values of the input matrix))
69///
70/// The third Array is the output array containing V ^ H
71pub fn svd<T>(input: &Array<T>) -> (Array<T>, Array<T::BaseType>, Array<T>)
72where
73    T: HasAfEnum + FloatingPoint,
74    T::BaseType: HasAfEnum,
75{
76    unsafe {
77        let mut u: af_array = std::ptr::null_mut();
78        let mut s: af_array = std::ptr::null_mut();
79        let mut vt: af_array = std::ptr::null_mut();
80        let err_val = af_svd(
81            &mut u as *mut af_array,
82            &mut s as *mut af_array,
83            &mut vt as *mut af_array,
84            input.get(),
85        );
86        HANDLE_ERROR(AfError::from(err_val));
87        (u.into(), s.into(), vt.into())
88    }
89}
90
91/// Perform Singular Value Decomposition inplace
92///
93/// This function factorizes a matrix A into two unitary matrices U and Vt, and a diagonal matrix S
94/// such that
95///
96/// A = U∗S∗Vt
97///
98/// 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
99/// x N
100///
101/// # Parameters
102///
103/// - `in` is the input/output matrix. This will contain random data after the function call is
104/// complete.
105///
106/// # Return Values
107///
108/// A triplet of Arrays.
109///
110/// The first Array is the output array containing U
111///
112/// The second Array is the output array containing the diagonal values of sigma, (singular values of the input matrix))
113///
114/// The third Array is the output array containing V ^ H
115pub fn svd_inplace<T>(input: &mut Array<T>) -> (Array<T>, Array<T::BaseType>, Array<T>)
116where
117    T: HasAfEnum + FloatingPoint,
118    T::BaseType: HasAfEnum,
119{
120    unsafe {
121        let mut u: af_array = std::ptr::null_mut();
122        let mut s: af_array = std::ptr::null_mut();
123        let mut vt: af_array = std::ptr::null_mut();
124        let err_val = af_svd_inplace(
125            &mut u as *mut af_array,
126            &mut s as *mut af_array,
127            &mut vt as *mut af_array,
128            input.get(),
129        );
130        HANDLE_ERROR(AfError::from(err_val));
131        (u.into(), s.into(), vt.into())
132    }
133}
134
135/// Perform LU decomposition
136///
137/// # Parameters
138///
139/// - `input` is the input matrix
140///
141/// # Return Values
142///
143/// A triplet of Arrays.
144///
145/// The first Array will contain the lower triangular matrix of the LU decomposition.
146///
147/// The second Array will contain the lower triangular matrix of the LU decomposition.
148///
149/// The third Array will contain the permutation indices to map the input to the decomposition.
150pub fn lu<T>(input: &Array<T>) -> (Array<T>, Array<T>, Array<i32>)
151where
152    T: HasAfEnum + FloatingPoint,
153{
154    unsafe {
155        let mut lower: af_array = std::ptr::null_mut();
156        let mut upper: af_array = std::ptr::null_mut();
157        let mut pivot: af_array = std::ptr::null_mut();
158        let err_val = af_lu(
159            &mut lower as *mut af_array,
160            &mut upper as *mut af_array,
161            &mut pivot as *mut af_array,
162            input.get(),
163        );
164        HANDLE_ERROR(AfError::from(err_val));
165        (lower.into(), upper.into(), pivot.into())
166    }
167}
168
169/// Perform inplace LU decomposition
170///
171/// # Parameters
172///
173/// - `input` contains the input matrix on entry and packed LU decomposition on exit
174/// - `is_lapack_pic` specified if the pivot is returned in original LAPACK compliant format
175///
176/// # Return Values
177///
178/// An Array with permutation indices to map the input to the decomposition. Since, the input
179/// matrix is modified in place, only pivot values are returned.
180pub fn lu_inplace<T>(input: &mut Array<T>, is_lapack_piv: bool) -> Array<i32>
181where
182    T: HasAfEnum + FloatingPoint,
183{
184    unsafe {
185        let mut pivot: af_array = std::ptr::null_mut();
186        let err_val = af_lu_inplace(&mut pivot as *mut af_array, input.get(), is_lapack_piv);
187        HANDLE_ERROR(AfError::from(err_val));
188        pivot.into()
189    }
190}
191
192/// Perform QR decomposition
193///
194/// # Parameters
195///
196/// - `input` is the input matrix
197///
198/// # Return Values
199///
200/// A triplet of Arrays.
201///
202/// The first Array is the orthogonal matrix from QR decomposition
203///
204/// The second Array is the upper triangular matrix from QR decomposition
205///
206/// The third Array will contain additional information needed for solving a least squares problem
207/// using q and r
208pub fn qr<T>(input: &Array<T>) -> (Array<T>, Array<T>, Array<T>)
209where
210    T: HasAfEnum + FloatingPoint,
211{
212    unsafe {
213        let mut q: af_array = std::ptr::null_mut();
214        let mut r: af_array = std::ptr::null_mut();
215        let mut tau: af_array = std::ptr::null_mut();
216        let err_val = af_qr(
217            &mut q as *mut af_array,
218            &mut r as *mut af_array,
219            &mut tau as *mut af_array,
220            input.get(),
221        );
222        HANDLE_ERROR(AfError::from(err_val));
223        (q.into(), r.into(), tau.into())
224    }
225}
226
227/// Perform inplace QR decomposition
228///
229/// # Parameters
230///
231/// - `input` contains the input matrix on entry, and packed QR decomposition on exit
232///
233/// # Return Values
234///
235/// An Array with additional information needed for unpacking the data.
236pub fn qr_inplace<T>(input: &mut Array<T>) -> Array<T>
237where
238    T: HasAfEnum + FloatingPoint,
239{
240    unsafe {
241        let mut tau: af_array = std::ptr::null_mut();
242        let err_val = af_qr_inplace(&mut tau as *mut af_array, input.get());
243        HANDLE_ERROR(AfError::from(err_val));
244        tau.into()
245    }
246}
247
248/// Perform Cholesky decomposition
249///
250/// # Parameters
251///
252/// - `input` is the input matrix
253/// - `is_upper` is a boolean to indicate if the output has to be upper or lower triangular matrix
254///
255/// # Return Values
256///
257/// A tuple of an Array and signed 32-bit integer.
258///
259/// The Array contains the triangular matrix (multiply it with conjugate transpose to reproduce the input).
260///
261/// If the integer is 0, it means the cholesky decomposition passed. Otherwise, it will contain the rank at
262/// which the decomposition failed.
263pub fn cholesky<T>(input: &Array<T>, is_upper: bool) -> (Array<T>, i32)
264where
265    T: HasAfEnum + FloatingPoint,
266{
267    unsafe {
268        let mut temp: af_array = std::ptr::null_mut();
269        let mut info: i32 = 0;
270        let err_val = af_cholesky(
271            &mut temp as *mut af_array,
272            &mut info as *mut c_int,
273            input.get(),
274            is_upper,
275        );
276        HANDLE_ERROR(AfError::from(err_val));
277        (temp.into(), info)
278    }
279}
280
281/// Perform inplace Cholesky decomposition
282///
283/// # Parameters
284///
285/// - `input` contains the input matrix on entry, and triangular matrix on exit.
286/// - `is_upper` is a boolean to indicate if the output has to be upper or lower triangular matrix
287///
288/// # Return Values
289///
290/// A signed 32-bit integer. If the integer is 0, it means the cholesky decomposition passed. Otherwise,
291/// it will contain the rank at which the decomposition failed.
292pub fn cholesky_inplace<T>(input: &mut Array<T>, is_upper: bool) -> i32
293where
294    T: HasAfEnum + FloatingPoint,
295{
296    let mut info: i32 = 0;
297    unsafe {
298        let err_val = af_cholesky_inplace(&mut info as *mut c_int, input.get(), is_upper);
299        HANDLE_ERROR(AfError::from(err_val));
300    }
301    info
302}
303
304/// Solve a system of equations
305///
306/// # Parameters
307///
308/// - `a` is the coefficient matrix
309/// - `b` has the measured values
310/// - `options` determine the various properties of matrix a
311///
312/// The `options` parameter currently needs to be either `NONE`, `LOWER` or `UPPER`, other values are not supported yet.
313///
314/// # Return Values
315///
316/// An Array which is the matrix of unknown variables
317pub fn solve<T>(a: &Array<T>, b: &Array<T>, options: MatProp) -> Array<T>
318where
319    T: HasAfEnum + FloatingPoint,
320{
321    unsafe {
322        let mut temp: af_array = std::ptr::null_mut();
323        let err_val = af_solve(
324            &mut temp as *mut af_array,
325            a.get(),
326            b.get(),
327            options as c_uint,
328        );
329        HANDLE_ERROR(AfError::from(err_val));
330        temp.into()
331    }
332}
333
334/// Solve a system of equations
335///
336/// # Parameters
337///
338/// - `a` is the output matrix from packed LU decomposition of the coefficient matrix
339/// - `piv` is the pivot array from packed LU decomposition of the coefficient matrix
340/// - `b` has the measured values
341/// - `options` determine the various properties of matrix a
342///
343/// The `options` parameter currently needs to be `NONE`, other values are not supported yet.
344///
345/// # Return Values
346///
347/// An Array which is the matrix of unknown variables
348pub fn solve_lu<T>(a: &Array<T>, piv: &Array<i32>, b: &Array<T>, options: MatProp) -> Array<T>
349where
350    T: HasAfEnum + FloatingPoint,
351{
352    unsafe {
353        let mut temp: af_array = std::ptr::null_mut();
354        let err_val = af_solve_lu(
355            &mut temp as *mut af_array,
356            a.get(),
357            piv.get(),
358            b.get(),
359            options as c_uint,
360        );
361        HANDLE_ERROR(AfError::from(err_val));
362        temp.into()
363    }
364}
365
366/// Compute inverse of a matrix
367///
368/// # Parameters
369///
370/// - `input` is the input matrix
371/// - `options` determine various properties of input matrix
372///
373/// The parameter `options` currently take only the value `NONE`.
374///
375/// # Return Values
376///
377/// An Array with values of the inverse of input matrix.
378pub fn inverse<T>(input: &Array<T>, options: MatProp) -> Array<T>
379where
380    T: HasAfEnum + FloatingPoint,
381{
382    unsafe {
383        let mut temp: af_array = std::ptr::null_mut();
384        let err_val = af_inverse(&mut temp as *mut af_array, input.get(), options as c_uint);
385        HANDLE_ERROR(AfError::from(err_val));
386        temp.into()
387    }
388}
389
390/// Find rank of a matrix
391///
392/// # Parameters
393///
394/// - `input` is the input matrix
395/// - `tol` is the tolerance value
396///
397/// # Return Values
398///
399/// An unsigned 32-bit integer which is the rank of the input matrix.
400pub fn rank<T>(input: &Array<T>, tol: f64) -> u32
401where
402    T: HasAfEnum + FloatingPoint,
403{
404    let mut temp: u32 = 0;
405    unsafe {
406        let err_val = af_rank(&mut temp as *mut c_uint, input.get(), tol);
407        HANDLE_ERROR(AfError::from(err_val));
408    }
409    temp
410}
411
412/// Find the determinant of the matrix
413///
414/// # Parameters
415///
416/// - `input` is the input matrix
417///
418/// # Return Values
419///
420/// A tuple of 32-bit floating point values.
421///
422/// If the input matrix is non-complex type, only first values of tuple contains the result.
423pub fn det<T>(input: &Array<T>) -> (f64, f64)
424where
425    T: HasAfEnum + FloatingPoint,
426{
427    let mut real: f64 = 0.0;
428    let mut imag: f64 = 0.0;
429    unsafe {
430        let err_val = af_det(
431            &mut real as *mut c_double,
432            &mut imag as *mut c_double,
433            input.get(),
434        );
435        HANDLE_ERROR(AfError::from(err_val));
436    }
437    (real, imag)
438}
439
440/// Find the norm of a matrix
441///
442/// # Parameters
443///
444/// - `input` is the input matrix
445/// - `ntype` is specifies the required norm type using enum [NormType](./enum.NormType.html)
446/// - `p` specifies the value of *P* when `ntype` is one of VECTOR_P, MATRIX_L_PQ. It is ignored
447/// for other values of `ntype`
448/// - `q` specifies the value of *Q* when `ntype` is MATRIX_L_PQ. This parameter is ignored if
449/// `ntype` is anything else.
450///
451/// # Return Values
452///
453/// A 64-bit floating point value that contains the norm of input matrix.
454pub fn norm<T>(input: &Array<T>, ntype: NormType, p: f64, q: f64) -> f64
455where
456    T: HasAfEnum + FloatingPoint,
457{
458    let mut out: f64 = 0.0;
459    unsafe {
460        let err_val = af_norm(
461            &mut out as *mut c_double,
462            input.get(),
463            ntype as c_uint,
464            p,
465            q,
466        );
467        HANDLE_ERROR(AfError::from(err_val));
468    }
469    out
470}
471
472/// Function to check if lapack support is available
473///
474/// # Parameters
475///
476/// None
477///
478/// # Return Values
479///
480/// Return a boolean indicating if ArrayFire was compiled with lapack support
481pub fn is_lapack_available() -> bool {
482    let mut temp: bool = false;
483    unsafe {
484        af_is_lapack_available(&mut temp as *mut bool);
485    }
486    temp
487}
488
489/// Psuedo Inverse of Matrix
490///
491/// # Parameters
492///
493/// - `input` is input matrix
494/// - `tolerance` defines the lower threshold for singular values from SVD
495/// - `option` must be [MatProp::NONE](./enum.MatProp.html) (more options might be supported in the future)
496///
497/// Notes:
498///
499/// - Tolerance is not the actual lower threshold, but it is passed in as a
500///   parameter to the calculation of the actual threshold relative to the shape and contents of input.
501/// - First, try setting tolerance to 1e-6 for single precision and 1e-12 for double.
502///
503/// # Return
504///
505/// Pseudo Inverse matrix for the input matrix
506pub fn pinverse<T>(input: &Array<T>, tolerance: f64, option: MatProp) -> Array<T>
507where
508    T: HasAfEnum + FloatingPoint,
509{
510    unsafe {
511        let mut out: af_array = std::ptr::null_mut();
512        let err_val = af_pinverse(
513            &mut out as *mut af_array,
514            input.get(),
515            tolerance,
516            option as c_uint,
517        );
518        HANDLE_ERROR(AfError::from(err_val));
519        out.into()
520    }
521}