peroxide/structure/matrix.rs
1//! Matrix for Scientific computation
2//!
3//! ## Declare matrix
4//!
5//! * You can declare matrix by various ways.
6//! * R's way - Default
7//! * MATLAB's way
8//! * Python's way
9//! * Other macro
10//!
11//! ### R's way
12//!
13//! * Description: Same as R - `matrix(Vec<f64>, Row, Col, Shape)`
14//! * Type: `matrix(Vec<T>, usize, usize, Shape) where T: std::convert::Into<f64> + Copy`
15//! * `Shape`: `Enum` for matrix shape - `Row` & `Col`
16//!
17//! ```rust
18//! #[macro_use]
19//! extern crate peroxide;
20//! use peroxide::fuga::*;
21//!
22//! fn main() {
23//! let a = matrix(c!(1,2,3,4), 2, 2, Row);
24//! a.print();
25//! // c[0] c[1]
26//! // r[0] 1 2
27//! // r[1] 3 4
28//!
29//! let b = matrix(c!(1,2,3,4), 2, 2, Col);
30//! b.print();
31//! // c[0] c[1]
32//! // r[0] 1 3
33//! // r[1] 2 4
34//! }
35//! ```
36//!
37//! ### MATLAB's way
38//!
39//! * Description: Similar to MATLAB (But should use `&str`)
40//! * Type: `ml_matrix(&str)`
41//!
42//! ```rust
43//! use peroxide::fuga::*;
44//!
45//! fn main() {
46//! let a = ml_matrix("1 2; 3 4");
47//! a.print();
48//! // c[0] c[1]
49//! // r[0] 1 2
50//! // r[1] 3 4
51//! }
52//! ```
53//!
54//! ### Python's way
55//!
56//! * Description: Declare matrix as vector of vectors.
57//! * Type: `py_matrix(Vec<Vec<T>>) where T: std::convert::Into<f64> + Copy`
58//!
59//! ```rust
60//! use peroxide::fuga::*;
61//!
62//! fn main() {
63//! let a = py_matrix(vec![vec![1, 2], vec![3, 4]]);
64//! a.print();
65//! // c[0] c[1]
66//! // r[0] 1 2
67//! // r[1] 3 4
68//! }
69//! ```
70//!
71//! ### Other macro
72//!
73//! * Description: R-like macro to declare matrix
74//! * For `R`,
75//!
76//! ```R
77//! # R
78//! a = matrix(1:4, nrow = 2, ncol = 2, byrow = T)
79//! print(a)
80//! # [,1] [,2]
81//! # [1,] 1 2
82//! # [2,] 3 4
83//! ```
84//!
85//! * For `Peroxide`,
86//!
87//! ```rust
88//! #[macro_use]
89//! extern crate peroxide;
90//! use peroxide::fuga::*;
91//!
92//! fn main() {
93//! let a = matrix!(1;4;1, 2, 2, Row);
94//! a.print();
95//! // c[0] c[1]
96//! // r[0] 1 2
97//! // r[1] 3 4
98//! }
99//! ```
100//!
101//! ## Basic Method for Matrix
102//!
103//! There are some useful methods for `Matrix`
104//!
105//! * `row(&self, index: usize) -> Vec<f64>` : Extract specific row as `Vec<f64>`
106//! * `col(&self, index: usize) -> Vec<f64>` : Extract specific column as `Vec<f64>`
107//! * `diag(&self) -> Vec<f64>`: Extract diagonal components as `Vec<f64>`
108//! * `swap(&self, usize, usize, Shape)`: Swap two rows or columns (unsafe function)
109//! * `subs_col(&mut self, usize, Vec<f64>)`: Substitute column with `Vec<f64>`
110//! * `subs_row(&mut self, usize, Vec<f64>)`: Substitute row with `Vec<f64>`
111//!
112//! ```rust
113//! #[macro_use]
114//! extern crate peroxide;
115//! use peroxide::fuga::*;
116//!
117//! fn main() {
118//! let mut a = ml_matrix("1 2; 3 4");
119//!
120//! a.row(0).print(); // [1, 2]
121//! a.col(0).print(); // [1, 3]
122//! a.diag().print(); // [1, 4]
123//! unsafe {
124//! a.swap(0, 1, Row);
125//! }
126//! a.print();
127//! // c[0] c[1]
128//! // r[0] 3 4
129//! // r[1] 1 2
130//!
131//! let mut b = ml_matrix("1 2;3 4");
132//! b.subs_col(0, &c!(5, 6));
133//! b.subs_row(1, &c!(7, 8));
134//! b.print();
135//! // c[0] c[1]
136//! // r[0] 5 2
137//! // r[1] 7 8
138//! }
139//! ```
140//!
141//! ## Read & Write
142//!
143//! In peroxide, we can write matrix to `csv`
144//!
145//! ### CSV (Not recommended)
146//!
147//! * `csv` feature should be required
148//! * `write(&self, file_path: &str)`: Write matrix to csv
149//! * `write_with_header(&self, file_path, header: Vec<&str>)`: Write with header
150//!
151//! ```rust
152//! use peroxide::fuga::*;
153//!
154//! fn main() {
155//! # #[cfg(feature="csv")]
156//! # {
157//! let a = ml_matrix("1 2;3 4");
158//! a.write("example_data/matrix.csv").expect("Can't write file");
159//!
160//! let b = ml_matrix("1 2; 3 4; 5 6");
161//! b.write_with_header("example_data/header.csv", vec!["odd", "even"])
162//! .expect("Can't write header file");
163//! # }
164//! println!("Complete!")
165//! }
166//! ```
167//!
168//! Also, you can read matrix from csv.
169//!
170//! * Type: `read(&str, bool, char) -> Result<Matrix, Box<Error>>`
171//! * Description: `read(file_path, is_header, delimiter)`
172//!
173//! ```no_run
174//! use peroxide::fuga::*;
175//!
176//! fn main() {
177//! # #[cfg(feature="csv")]
178//! # {
179//! // `csv` feature should be required
180//! let a = Matrix::read("example_data/matrix.csv", false, ',')
181//! .expect("Can't read matrix.csv file");
182//! a.print();
183//! # }
184//! // c[0] c[1]
185//! // r[0] 1 2
186//! // r[1] 3 4
187//! }
188//! ```
189//!
190//! ### Convert to DataFrame (Recommended)
191//!
192//! To write columns or rows, `DataFrame` and `nc` feature could be the best choice.
193//!
194//! * `nc` feature should be required - `netcdf` or `libnetcdf` are prerequisites.
195//!
196//! ```no_run
197//! #[macro_use]
198//! extern crate peroxide;
199//! use peroxide::fuga::*;
200//!
201//! fn main() {
202//! let a = matrix(c!(1,2,3,4,5,6), 3, 2, Col);
203//!
204//! // Construct DataFrame
205//! let mut df = DataFrame::new(vec![]);
206//! df.push("x", Series::new(a.col(0)));
207//! df.push("y", Series::new(a.col(1)));
208//!
209//! // Write nc file (`nc` feature should be required)
210//! # #[cfg(feature="nc")]
211//! # {
212//! df.write_nc("data.nc").expect("Can't write data.nc");
213//!
214//! // Read nc file (`nc` feature should be required)
215//! let dg = DataFrame::read_nc("data.nc").expect("Can't read data.nc");
216//! let x: Vec<f64> = dg["x"].to_vec();
217//! let y: Vec<f64> = dg["y"].to_vec();
218//!
219//! assert_eq!(a.col(0), x);
220//! assert_eq!(a.col(1), y);
221//! # }
222//! }
223//! ```
224//!
225//! ## Concatenation
226//!
227//! There are two options to concatenate matrices.
228//!
229//! * `cbind`: Concatenate two matrices by column direction.
230//! * `rbind`: Concatenate two matrices by row direction.
231//!
232//! ```rust
233//! use peroxide::fuga::*;
234//!
235//! fn main() -> Result<(), Box<dyn Error>> {
236//! let a = ml_matrix("1 2;3 4");
237//! let b = ml_matrix("5 6;7 8");
238//!
239//! cbind(a.clone(), b.clone())?.print();
240//! // c[0] c[1] c[2] c[3]
241//! // r[0] 1 2 5 7
242//! // r[1] 3 4 6 8
243//!
244//! rbind(a, b)?.print();
245//! // c[0] c[1]
246//! // r[0] 1 2
247//! // r[1] 3 4
248//! // r[2] 5 6
249//! // r[3] 7 8
250//!
251//! Ok(())
252//! }
253//! ```
254//!
255//! ## Matrix operations
256//!
257//! * In peroxide, can use basic operations between matrices. I'll show you by examples.
258//!
259//! ```rust
260//! #[macro_use]
261//! extern crate peroxide;
262//! use peroxide::fuga::*;
263//!
264//! fn main() {
265//! let a = matrix!(1;4;1, 2, 2, Row);
266//! (a.clone() + 1).print(); // -, *, / are also available
267//! // c[0] c[1]
268//! // r[0] 2 3
269//! // r[1] 4 5
270//!
271//! let b = matrix!(5;8;1, 2, 2, Row);
272//! (a.clone() + b.clone()).print(); // - is also available
273//! // c[0] c[1]
274//! // r[0] 6 8
275//! // r[1] 10 12
276//!
277//! (a.clone() * b.clone()).print(); // Matrix multiplication
278//! // c[0] c[1]
279//! // r[0] 19 22
280//! // r[1] 43 50
281//! }
282//! ```
283//!
284//! * `clone` is too annoying - We can use reference operations!
285//!
286//! ```rust
287//! use peroxide::fuga::*;
288//!
289//! fn main() {
290//! let a = ml_matrix("1 2;3 4");
291//! let b = ml_matrix("5 6;7 8");
292//!
293//! (&a + 1).print();
294//! (&a + &b).print();
295//! (&a - &b).print();
296//! (&a * &b).print();
297//! }
298//! ```
299//!
300//! ## Extract & modify components
301//!
302//! * In peroxide, matrix data is saved as linear structure.
303//! * But you can use two-dimensional index to extract or modify components.
304//!
305//! ```rust
306//! #[macro_use]
307//! extern crate peroxide;
308//! use peroxide::fuga::*;
309//!
310//! fn main() {
311//! let mut a = matrix!(1;4;1, 2, 2, Row);
312//! a[(0,0)].print(); // 1
313//! a[(0,0)] = 2f64; // Modify component
314//! a.print();
315//! // c[0] c[1]
316//! // r[0] 2 2
317//! // r[1] 3 4
318//! }
319//! ```
320//!
321//! ## Conversion to `Vec<f64>`
322//!
323//! * Just use `row` or `col` method (I already showed at Basic method section).
324//!
325//! ```rust
326//! #[macro_use]
327//! extern crate peroxide;
328//! use peroxide::fuga::*;
329//!
330//! fn main() {
331//! let a = matrix!(1;4;1, 2, 2, Row);
332//! a.row(0).print(); // [1, 2]
333//! assert_eq!(a.row(0), vec![1f64, 2f64]);
334//! }
335//! ```
336//!
337//! ## Useful constructor
338//!
339//! * `zeros(usize, usize)`: Construct matrix which elements are all zero
340//! * `eye(usize)`: Identity matrix
341//! * `rand(usize, usize)`: Construct random uniform matrix (from 0 to 1)
342//! * `rand_with_rng(usize, usize, &mut Rng)`: Construct random matrix with user-defined RNG
343//!
344//! ```rust
345//! use peroxide::fuga::*;
346//!
347//! fn main() {
348//! let a = zeros(2, 2);
349//! assert_eq!(a, ml_matrix("0 0;0 0"));
350//!
351//! let b = eye(2);
352//! assert_eq!(b, ml_matrix("1 0;0 1"));
353//!
354//! let c = rand(2, 2);
355//! c.print(); // Random 2x2 matrix
356//! }
357//! ```
358//! # Linear Algebra
359//!
360//! ## Transpose
361//!
362//! * Caution: Transpose does not consume the original value.
363//!
364//! ```rust
365//! #[macro_use]
366//! extern crate peroxide;
367//! use peroxide::fuga::*;
368//!
369//! fn main() {
370//! let a = matrix!(1;4;1, 2, 2, Row);
371//! a.transpose().print();
372//! // Or you can use shorter one
373//! a.t().print();
374//! // c[0] c[1]
375//! // r[0] 1 3
376//! // r[1] 2 4
377//! }
378//! ```
379//!
380//! ## LU Decomposition
381//!
382//! * Peroxide uses **complete pivoting** for LU decomposition - Very stable
383//! * Since there are lots of causes to generate error, you should use `Option`
384//! * `lu` returns `PQLU`
385//! * `PQLU` has four field - `p`, `q`, `l` , `u`
386//! * `p` means row permutations
387//! * `q` means column permutations
388//! * `l` means lower triangular matrix
389//! * `u` menas upper triangular matrix
390//! * The structure of `PQLU` is as follows:
391//!
392//! ```rust
393//! use peroxide::fuga::*;
394//!
395//! #[derive(Debug, Clone)]
396//! pub struct PQLU {
397//! pub p: Vec<usize>,
398//! pub q: Vec<usize>,
399//! pub l: Matrix,
400//! pub u: Matrix,
401//! }
402//! ```
403//!
404//! * Example of LU decomposition:
405//!
406//! ```rust
407//! #[macro_use]
408//! extern crate peroxide;
409//! use peroxide::fuga::*;
410//!
411//! fn main() {
412//! let a = matrix(c!(1,2,3,4), 2, 2, Row);
413//! let pqlu = a.lu();
414//! let (p,q,l,u) = (pqlu.p, pqlu.q, pqlu.l, pqlu.u);
415//! assert_eq!(p, vec![1]); // swap 0 & 1 (Row)
416//! assert_eq!(q, vec![1]); // swap 0 & 1 (Col)
417//! assert_eq!(l, matrix(c!(1,0,0.5,1),2,2,Row));
418//! // c[0] c[1]
419//! // r[0] 1 0
420//! // r[1] 0.5 1
421//! assert_eq!(u, matrix(c!(4,3,0,-0.5),2,2,Row));
422//! // c[0] c[1]
423//! // r[0] 4 3
424//! // r[1] 0 -0.5
425//! }
426//! ```
427//!
428//! ## Determinant
429//!
430//! * Peroxide uses LU decomposition to obtain determinant ($ \mathcal{O}(n^3) $)
431//!
432//! ```rust
433//! #[macro_use]
434//! extern crate peroxide;
435//! use peroxide::fuga::*;
436//!
437//! fn main() {
438//! let a = matrix!(1;4;1, 2, 2, Row);
439//! assert_eq!(a.det(), -2f64);
440//! }
441//! ```
442//!
443//! ## Inverse matrix
444//!
445//! * Peroxide uses LU decomposition (via GECP) to obtain inverse matrix.
446//! * It needs two sub functions - `inv_l`, `inv_u`
447//! * For inverse of `L, U`, I use block partitioning. For example, for lower triangular matrix :
448//! $$ \begin{aligned} L &= \begin{pmatrix} L_1 & \mathbf{0} \\\ L_2 & L_3 \end{pmatrix} \\\ L^{-1} &= \begin{pmatrix} L_1^{-1} & \mathbf{0} \\\ -L_3^{-1}L_2 L_1^{-1} & L_3^{-1} \end{pmatrix} \end{aligned} $$
449//! ```rust
450//! #[macro_use]
451//! extern crate peroxide;
452//! use peroxide::fuga::*;
453//!
454//! fn main() {
455//! let a = matrix!(1;4;1, 2, 2, Row);
456//! a.inv().print();
457//! // c[0] c[1]
458//! // r[0] -2 1
459//! // r[1] 1.5 -0.5
460//! }
461//! ```
462//!
463//! ## Tips for LU, Det, Inverse
464//!
465//! * If you save `self.lu()` rather than the direct use of `self.det()` or `self.lu()` then you can get better performance (via memoization)
466//!
467//! ```rust
468//! use peroxide::fuga::*;
469//!
470//! fn main() {
471//! let a = ml_matrix("1 2;3 4");
472//! let pqlu = a.lu(); // Memoization of LU
473//! pqlu.det().print(); // Same as a.det() but do not need an additional LU
474//! pqlu.inv().print(); // Same as a.inv() but do not need an additional LU
475//! }
476//! ```
477//!
478//! ## QR Decomposition (`O3` feature only)
479//!
480//! * Use `dgeqrf` of LAPACK
481//! * Return `QR` structure.
482//!
483//! ```ignore
484//! pub struct QR {
485//! pub q: Matrix,
486//! pub r: Matrix,
487//! }
488//! ```
489//!
490//! * Example
491//!
492//! ```ignore
493//! use peroxide::fuga::*;
494//!
495//! let a = ml_matrix("
496//! 1 -1 4;
497//! 1 4 -2;
498//! 1 4 2;
499//! 1 -1 0
500//! ");
501//!
502//! // QR decomposition
503//! let qr = a.qr();
504//!
505//! qr.q().print();
506//! // c[0] c[1] c[2]
507//! // r[0] -0.5 0.5 -0.5
508//! // r[1] -0.5 -0.5000 0.5000
509//! // r[2] -0.5 -0.5000 -0.5
510//! // r[3] -0.5 0.5 0.5000
511//!
512//! qr.r().print();
513//! // c[0] c[1] c[2]
514//! // r[0] -2 -3 -2
515//! // r[1] 0 -5 2
516//! // r[2] 0 0 -4
517//! ```
518//!
519//! ## Singular Value Decomposition (`O3` feature only)
520//!
521//! * Use `dgesvd` of LAPACK
522//! * Return `SVD` structure
523//!
524//! ```no_run
525//! use peroxide::fuga::*;
526//!
527//! pub struct SVD {
528//! pub s: Vec<f64>,
529//! pub u: Matrix,
530//! pub vt: Matrix,
531//! }
532//! ```
533//!
534//! * Example
535//!
536//! ```
537//! use peroxide::fuga::*;
538//!
539//! fn main() {
540//! let a = ml_matrix("3 2 2;2 3 -2");
541//! #[cfg(feature="O3")]
542//! {
543//! // Full SVD
544//! let svd = a.svd();
545//! assert!(eq_vec(&vec![5f64, 3f64], &svd.s, 1e-7));
546//!
547//! // Or Truncated SVD
548//! let svd2 = svd.truncated();
549//! assert!(eq_vec(&vec![5f64, 3f64], &svd2.s, 1e-7));
550//! }
551//! a.print();
552//! }
553//! ```
554//!
555//! ## Cholesky Decomposition
556//!
557//! * Use `dpotrf` of LAPACK
558//! * Return Matrix (But there can be panic! - Not symmetric or Not positive definite)
559//! * Example
560//!
561//! ```rust
562//! extern crate peroxide;
563//! use peroxide::fuga::*;
564//!
565//! fn main() {
566//! let a = ml_matrix("1 2;2 5");
567//! #[cfg(feature = "O3")]
568//! {
569//! let u = a.cholesky(Upper);
570//! assert_eq!(u, ml_matrix("1 2;0 1"));
571//!
572//! let l = a.cholesky(Lower);
573//! assert_eq!(l, ml_matrix("1 0;2 1"));
574//! }
575//! a.print();
576//! }
577//! ```
578//!
579//! ## Moore-Penrose Pseudo Inverse
580//!
581//! * $ X^\dagger = \left(X^T X\right)^{-1} X^T $
582//! * For `O3` feature, peroxide use SVD to obtain pseudo inverse
583//!
584//! ```rust
585//! #[macro_use]
586//! extern crate peroxide;
587//! use peroxide::fuga::*;
588//!
589//! fn main() {
590//! let a = matrix!(1;4;1, 2, 2, Row);
591//! let pinv_a = a.pseudo_inv();
592//! let inv_a = a.inv();
593//!
594//! assert_eq!(inv_a, pinv_a); // Nearly equal (not actually equal)
595//! }
596//! ```
597
598#[cfg(feature = "csv")]
599extern crate csv;
600
601#[cfg(feature = "csv")]
602use self::csv::{ReaderBuilder, StringRecord, WriterBuilder};
603
604#[cfg(feature = "O3")]
605extern crate blas;
606#[cfg(feature = "O3")]
607extern crate lapack;
608#[cfg(feature = "O3")]
609use blas::{daxpy, dgemm, dgemv};
610#[cfg(feature = "O3")]
611use lapack::{dgecon, dgeqrf, dgesvd, dgetrf, dgetri, dgetrs, dorgqr, dpotrf};
612#[cfg(feature = "parallel")]
613use rayon::iter::{IntoParallelIterator, IntoParallelRefIterator, ParallelIterator};
614#[cfg(feature = "serde")]
615use serde::{Deserialize, Serialize};
616
617pub use self::Shape::{Col, Row};
618use crate::numerical::eigen::{eigen, EigenMethod};
619use crate::structure::dataframe::{Series, TypedVector};
620#[cfg(feature = "parallel")]
621use crate::traits::math::{ParallelInnerProduct, ParallelNormed};
622use crate::traits::sugar::ScalableMut;
623use crate::traits::{
624 fp::{FPMatrix, FPVector},
625 general::Algorithm,
626 math::{InnerProduct, LinearOp, MatrixProduct, Norm, Normed, Vector},
627 matrix::{Form, LinearAlgebra, MatrixTrait, SolveKind, PQLU, QR, SVD, UPLO, WAZD},
628 mutable::MutMatrix,
629};
630#[allow(unused_imports)]
631use crate::util::{
632 low_level::{copy_vec_ptr, swap_vec_ptr},
633 non_macro::{cbind, eye, rbind, zeros},
634 useful::{nearly_eq, tab},
635};
636use peroxide_num::{ExpLogOps, Numeric, PowOps, TrigOps};
637use std::cmp::{max, min};
638pub use std::error::Error;
639use std::fmt;
640use std::ops::{Add, Div, Index, IndexMut, Mul, Neg, Sub};
641
642pub type Perms = Vec<(usize, usize)>;
643
644/// To select matrices' binding.
645///
646/// Row - Row binding
647/// Col - Column binding
648///
649/// # Examples
650/// ```
651/// use peroxide::fuga::*;
652///
653/// let a = matrix(vec![1,2,3,4], 2, 2, Row);
654/// let b = matrix(vec![1,2,3,4], 2, 2, Col);
655/// println!("{}", a); // Similar to [[1,2],[3,4]]
656/// println!("{}", b); // Similar to [[1,3],[2,4]]
657/// ```
658#[derive(Default, Debug, PartialEq, Clone, Copy)]
659#[cfg_attr(feature = "serde", derive(Serialize, Deserialize))]
660#[cfg_attr(
661 feature = "rkyv",
662 derive(rkyv::Archive, rkyv::Serialize, rkyv::Deserialize)
663)]
664pub enum Shape {
665 #[default]
666 Col,
667 Row,
668}
669
670/// Print for Shape
671impl fmt::Display for Shape {
672 fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
673 let to_display = match self {
674 Row => "Row",
675 Col => "Col",
676 };
677 write!(f, "{}", to_display)
678 }
679}
680
681/// R-like matrix structure
682///
683/// # Examples
684///
685/// ```
686/// use peroxide::fuga::*;
687///
688/// let a = Matrix {
689/// data: vec![1f64,2f64,3f64,4f64],
690/// row: 2,
691/// col: 2,
692/// shape: Row,
693/// }; // [[1,2],[3,4]]
694/// ```
695#[derive(Debug, Clone, Default)]
696#[cfg_attr(feature = "serde", derive(Serialize, Deserialize))]
697#[cfg_attr(
698 feature = "rkyv",
699 derive(rkyv::Archive, rkyv::Serialize, rkyv::Deserialize)
700)]
701pub struct Matrix {
702 pub data: Vec<f64>,
703 pub row: usize,
704 pub col: usize,
705 pub shape: Shape,
706}
707
708// =============================================================================
709// Various matrix constructor
710// =============================================================================
711
712/// R-like matrix constructor
713///
714/// # Examples
715/// ```
716/// #[macro_use]
717/// extern crate peroxide;
718/// use peroxide::fuga::*;
719///
720/// fn main() {
721/// let a = matrix(c!(1,2,3,4), 2, 2, Row);
722/// a.print(); // Print matrix
723/// }
724/// ```
725pub fn matrix<T>(v: Vec<T>, r: usize, c: usize, shape: Shape) -> Matrix
726where
727 T: Into<f64>,
728{
729 Matrix {
730 data: v.into_iter().map(|t| t.into()).collect::<Vec<f64>>(),
731 row: r,
732 col: c,
733 shape,
734 }
735}
736
737/// R-like matrix constructor (Explicit ver.)
738pub fn r_matrix<T>(v: Vec<T>, r: usize, c: usize, shape: Shape) -> Matrix
739where
740 T: Into<f64>,
741{
742 matrix(v, r, c, shape)
743}
744
745/// Python-like matrix constructor
746///
747/// # Examples
748/// ```
749/// #[macro_use]
750/// extern crate peroxide;
751/// use peroxide::fuga::*;
752///
753/// fn main() {
754/// let a = py_matrix(vec![c!(1,2), c!(3,4)]);
755/// let b = matrix(c!(1,2,3,4), 2, 2, Row);
756/// assert_eq!(a, b);
757/// }
758/// ```
759pub fn py_matrix<T>(v: Vec<Vec<T>>) -> Matrix
760where
761 T: Into<f64> + Copy,
762{
763 let r = v.len();
764 let c = v[0].len();
765 let mut data = vec![0f64; r * c];
766 for i in 0..r {
767 for j in 0..c {
768 let idx = i * c + j;
769 data[idx] = v[i][j].into();
770 }
771 }
772 matrix(data, r, c, Row)
773}
774
775/// Matlab-like matrix constructor
776///
777/// # Examples
778/// ```
779/// #[macro_use]
780/// extern crate peroxide;
781/// use peroxide::fuga::*;
782///
783/// fn main() {
784/// let a = ml_matrix("1 2; 3 4");
785/// let b = matrix(c!(1,2,3,4), 2, 2, Row);
786/// assert_eq!(a, b);
787/// }
788/// ```
789pub fn ml_matrix(s: &str) -> Matrix where {
790 let str_rows: Vec<&str> = s.split(';').collect();
791 let r = str_rows.len();
792 let str_data = str_rows
793 .into_iter()
794 .map(|x| x.trim().split(' ').collect::<Vec<&str>>())
795 .collect::<Vec<Vec<&str>>>();
796 let c = str_data[0].len();
797 let data = str_data
798 .into_iter()
799 .flat_map(|t| {
800 t.into_iter()
801 .map(|x| x.parse::<f64>().unwrap())
802 .collect::<Vec<f64>>()
803 })
804 .collect::<Vec<f64>>();
805 matrix(data, r, c, Row)
806}
807
808/// Pretty Print
809impl fmt::Display for Matrix {
810 fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
811 write!(f, "{}", self.spread())
812 }
813}
814
815/// PartialEq implements
816impl PartialEq for Matrix {
817 fn eq(&self, other: &Matrix) -> bool {
818 if self.shape == other.shape {
819 self.data
820 .clone()
821 .into_iter()
822 .zip(other.data.clone())
823 .all(|(x, y)| nearly_eq(x, y))
824 && self.row == other.row
825 } else {
826 self.eq(&other.change_shape())
827 }
828 }
829}
830
831/// Main matrix structure
832#[allow(dead_code)]
833impl MatrixTrait for Matrix {
834 type Scalar = f64;
835
836 /// Raw pointer for `self.data`
837 fn ptr(&self) -> *const f64 {
838 &self.data[0] as *const f64
839 }
840
841 /// Raw mutable pointer for `self.data`
842 fn mut_ptr(&mut self) -> *mut f64 {
843 &mut self.data[0] as *mut f64
844 }
845
846 /// Slice of `self.data`
847 ///
848 /// # Examples
849 /// ```
850 /// use peroxide::fuga::*;
851 ///
852 /// let a = matrix(vec![1,2,3,4], 2, 2, Col);
853 /// let b = a.as_slice();
854 /// assert_eq!(b, &[1f64,2f64,3f64,4f64]);
855 /// ```
856 fn as_slice(&self) -> &[f64] {
857 &self.data[..]
858 }
859
860 /// Mutable slice of `self.data`
861 ///
862 /// # Examples
863 /// ```
864 /// use peroxide::fuga::*;
865 ///
866 /// let mut a = matrix(vec![1,2,3,4], 2, 2, Col);
867 /// let mut b = a.as_mut_slice();
868 /// b[0] = 5f64;
869 /// assert_eq!(b, &[5f64, 2f64, 3f64, 4f64]);
870 /// assert_eq!(a, matrix(vec![5,2,3,4], 2, 2, Col));
871 /// ```
872 fn as_mut_slice(&mut self) -> &mut [f64] {
873 &mut self.data[..]
874 }
875
876 /// Change Bindings
877 ///
878 /// `Row` -> `Col` or `Col` -> `Row`
879 ///
880 /// # Examples
881 /// ```
882 /// use peroxide::fuga::*;
883 ///
884 /// let a = matrix(vec![1,2,3,4],2,2,Row);
885 /// assert_eq!(a.shape, Row);
886 /// let b = a.change_shape();
887 /// assert_eq!(b.shape, Col);
888 /// ```
889 fn change_shape(&self) -> Self {
890 let r = self.row;
891 let c = self.col;
892 assert_eq!(r * c, self.data.len());
893 let l = r * c - 1;
894 let mut data: Vec<f64> = self.data.clone();
895 let ref_data = &self.data;
896
897 match self.shape {
898 Row => {
899 for i in 0..l {
900 let s = (i * c) % l;
901 data[i] = ref_data[s];
902 }
903 data[l] = ref_data[l];
904 matrix(data, r, c, Col)
905 }
906 Col => {
907 for i in 0..l {
908 let s = (i * r) % l;
909 data[i] = ref_data[s];
910 }
911 data[l] = ref_data[l];
912 matrix(data, r, c, Row)
913 }
914 }
915 }
916
917 /// Change Bindings Mutably
918 ///
919 /// `Row` -> `Col` or `Col` -> `Row`
920 ///
921 /// # Examples
922 /// ```
923 /// use peroxide::fuga::*;
924 ///
925 /// let a = matrix(vec![1,2,3,4],2,2,Row);
926 /// assert_eq!(a.shape, Row);
927 /// let b = a.change_shape();
928 /// assert_eq!(b.shape, Col);
929 /// ```
930 fn change_shape_mut(&mut self) {
931 let r = self.row;
932 let c = self.col;
933 assert_eq!(r * c, self.data.len());
934 let l = r * c - 1;
935 let ref_data = self.data.clone();
936
937 match self.shape {
938 Row => {
939 for i in 0..l {
940 let s = (i * c) % l;
941 self.data[i] = ref_data[s];
942 }
943 self.data[l] = ref_data[l];
944 self.shape = Col;
945 }
946 Col => {
947 for i in 0..l {
948 let s = (i * r) % l;
949 self.data[i] = ref_data[s];
950 }
951 self.data[l] = ref_data[l];
952 self.shape = Row;
953 }
954 }
955 }
956
957 /// Spread data(1D vector) to 2D formatted String
958 ///
959 /// # Examples
960 /// ```
961 /// use peroxide::fuga::*;
962 ///
963 /// let a = matrix(vec![1,2,3,4],2,2,Row);
964 /// println!("{}", a.spread()); // same as println!("{}", a);
965 /// // Result:
966 /// // c[0] c[1]
967 /// // r[0] 1 3
968 /// // r[1] 2 4
969 /// ```
970 fn spread(&self) -> String {
971 assert_eq!(self.row * self.col, self.data.len());
972 let r = self.row;
973 let c = self.col;
974 let mut key_row = 20usize;
975 let mut key_col = 20usize;
976
977 if r > 100 || c > 100 || (r > 20 && c > 20) {
978 let part = if r <= 10 {
979 key_row = r;
980 key_col = 100;
981 self.take_col(100)
982 } else if c <= 10 {
983 key_row = 100;
984 key_col = c;
985 self.take_row(100)
986 } else {
987 self.take_row(20).take_col(20)
988 };
989 return format!(
990 "Result is too large to print - {}x{}\nonly print {}x{} parts:\n{}",
991 self.row,
992 self.col,
993 key_row,
994 key_col,
995 part.spread()
996 );
997 }
998
999 // Find maximum length of data
1000 let sample = self.data.clone();
1001 let mut space: usize = sample
1002 .into_iter()
1003 .map(
1004 |x| min(format!("{:.4}", x).len(), x.to_string().len()), // Choose minimum of approx vs normal
1005 )
1006 .fold(0, max)
1007 + 1;
1008
1009 if space < 5 {
1010 space = 5;
1011 }
1012
1013 let mut result = String::new();
1014
1015 result.push_str(&tab("", 5));
1016 for i in 0..c {
1017 result.push_str(&tab(&format!("c[{}]", i), space)); // Header
1018 }
1019 result.push('\n');
1020
1021 for i in 0..r {
1022 result.push_str(&tab(&format!("r[{}]", i), 5));
1023 for j in 0..c {
1024 let st1 = format!("{:.4}", self[(i, j)]); // Round at fourth position
1025 let st2 = self[(i, j)].to_string(); // Normal string
1026 let mut st = st2.clone();
1027
1028 // Select more small thing
1029 if st1.len() < st2.len() {
1030 st = st1;
1031 }
1032
1033 result.push_str(&tab(&st, space));
1034 }
1035 if i == (r - 1) {
1036 break;
1037 }
1038 result.push('\n');
1039 }
1040
1041 result
1042 }
1043
1044 /// Extract Column
1045 ///
1046 /// # Examples
1047 /// ```
1048 /// #[macro_use]
1049 /// extern crate peroxide;
1050 /// use peroxide::fuga::*;
1051 ///
1052 /// fn main() {
1053 /// let a = matrix(c!(1,2,3,4), 2, 2, Row);
1054 /// assert_eq!(a.col(0), c!(1,3));
1055 /// }
1056 /// ```
1057 fn col(&self, index: usize) -> Vec<f64> {
1058 assert!(index < self.col);
1059 let mut container: Vec<f64> = vec![0f64; self.row];
1060 for i in 0..self.row {
1061 container[i] = self[(i, index)];
1062 }
1063 container
1064 }
1065
1066 /// Extract Row
1067 ///
1068 /// # Examples
1069 /// ```
1070 /// #[macro_use]
1071 /// extern crate peroxide;
1072 /// use peroxide::fuga::*;
1073 ///
1074 /// fn main() {
1075 /// let a = matrix(c!(1,2,3,4), 2, 2, Row);
1076 /// assert_eq!(a.row(0), c!(1,2));
1077 /// }
1078 /// ```
1079 fn row(&self, index: usize) -> Vec<f64> {
1080 assert!(index < self.row);
1081 let mut container: Vec<f64> = vec![0f64; self.col];
1082 for i in 0..self.col {
1083 container[i] = self[(index, i)];
1084 }
1085 container
1086 }
1087
1088 /// Extract diagonal components
1089 ///
1090 /// # Examples
1091 /// ```
1092 /// #[macro_use]
1093 /// extern crate peroxide;
1094 /// use peroxide::fuga::*;
1095 ///
1096 /// fn main() {
1097 /// let a = matrix!(1;4;1, 2, 2, Row);
1098 /// assert_eq!(a.diag(), c!(1,4));
1099 /// }
1100 /// ```
1101 fn diag(&self) -> Vec<f64> {
1102 let mut container = vec![0f64; self.row];
1103 let r = self.row;
1104 let c = self.col;
1105 assert_eq!(r, c);
1106 let c2 = c + 1;
1107 for i in 0..r {
1108 container[i] = self.data[i * c2];
1109 }
1110 container
1111 }
1112
1113 /// Transpose
1114 ///
1115 /// # Examples
1116 /// ```
1117 /// use peroxide::fuga::*;
1118 ///
1119 /// let a = matrix(vec![1,2,3,4], 2, 2, Row);
1120 /// println!("{}", a); // [[1,3],[2,4]]
1121 /// ```
1122 fn transpose(&self) -> Self {
1123 match self.shape {
1124 Row => matrix(self.data.clone(), self.col, self.row, Col),
1125 Col => matrix(self.data.clone(), self.col, self.row, Row),
1126 }
1127 }
1128
1129 /// R-like transpose function
1130 ///
1131 /// # Examples
1132 /// ```
1133 /// #[macro_use]
1134 /// extern crate peroxide;
1135 /// use peroxide::fuga::*;
1136 ///
1137 /// fn main() {
1138 /// let a = matrix!(1;4;1, 2, 2, Row);
1139 /// assert_eq!(a.transpose(), a.t());
1140 /// }
1141 /// ```
1142 fn t(&self) -> Self {
1143 self.transpose()
1144 }
1145
1146 /// Substitute Col
1147 #[inline]
1148 fn subs_col(&mut self, idx: usize, v: &[f64]) {
1149 for i in 0..self.row {
1150 self[(i, idx)] = v[i];
1151 }
1152 }
1153
1154 /// Substitute Row
1155 #[inline]
1156 fn subs_row(&mut self, idx: usize, v: &[f64]) {
1157 for j in 0..self.col {
1158 self[(idx, j)] = v[j];
1159 }
1160 }
1161
1162 /// From index operations
1163 fn from_index<F, G>(f: F, size: (usize, usize)) -> Matrix
1164 where
1165 F: Fn(usize, usize) -> G + Copy,
1166 G: Into<f64>,
1167 {
1168 let row = size.0;
1169 let col = size.1;
1170
1171 let mut mat = matrix(vec![0f64; row * col], row, col, Row);
1172
1173 for i in 0..row {
1174 for j in 0..col {
1175 mat[(i, j)] = f(i, j).into();
1176 }
1177 }
1178 mat
1179 }
1180
1181 /// Matrix to `Vec<Vec<f64>>`
1182 ///
1183 /// To send `Matrix` to `inline-python`
1184 fn to_vec(&self) -> Vec<Vec<f64>> {
1185 let mut result = vec![vec![0f64; self.col]; self.row];
1186 for i in 0..self.row {
1187 result[i] = self.row(i);
1188 }
1189 result
1190 }
1191
1192 fn to_diag(&self) -> Matrix {
1193 assert_eq!(self.row, self.col, "Should be square matrix");
1194 let mut result = matrix(vec![0f64; self.row * self.col], self.row, self.col, Row);
1195 let diag = self.diag();
1196 for i in 0..self.row {
1197 result[(i, i)] = diag[i];
1198 }
1199 result
1200 }
1201
1202 /// Submatrix
1203 ///
1204 /// # Description
1205 /// Return below elements of matrix to new matrix
1206 ///
1207 /// $$
1208 /// \begin{pmatrix}
1209 /// \\ddots & & & & \\\\
1210 /// & start & \\cdots & end.1 & \\\\
1211 /// & \\vdots & \\ddots & \\vdots & \\\\
1212 /// & end.0 & \\cdots & end & \\\\
1213 /// & & & & \\ddots
1214 /// \end{pmatrix}
1215 /// $$
1216 ///
1217 /// # Examples
1218 /// ```
1219 /// extern crate peroxide;
1220 /// use peroxide::fuga::*;
1221 ///
1222 /// fn main() {
1223 /// let a = ml_matrix("1 2 3;4 5 6;7 8 9");
1224 /// let b = ml_matrix("5 6;8 9");
1225 /// let c = a.submat((1, 1), (2, 2));
1226 /// assert_eq!(b, c);
1227 /// }
1228 /// ```
1229 fn submat(&self, start: (usize, usize), end: (usize, usize)) -> Matrix {
1230 let row = end.0 + 1 - start.0;
1231 let col = end.1 + 1 - start.1;
1232 let mut result = matrix(vec![0f64; row * col], row, col, self.shape);
1233 for i in 0..row {
1234 for j in 0..col {
1235 result[(i, j)] = self[(start.0 + i, start.1 + j)];
1236 }
1237 }
1238 result
1239 }
1240
1241 /// Substitute matrix to specific position
1242 ///
1243 /// # Description
1244 /// Substitute below elements of matrix
1245 ///
1246 /// $$
1247 /// \begin{pmatrix}
1248 /// \\ddots & & & & \\\\
1249 /// & start & \\cdots & end.1 & \\\\
1250 /// & \\vdots & \\ddots & \\vdots & \\\\
1251 /// & end.0 & \\cdots & end & \\\\
1252 /// & & & & \\ddots
1253 /// \end{pmatrix}
1254 /// $$
1255 ///
1256 /// # Examples
1257 /// ```
1258 /// extern crate peroxide;
1259 /// use peroxide::fuga::*;
1260 ///
1261 /// fn main() {
1262 /// let mut a = ml_matrix("1 2 3;4 5 6;7 8 9");
1263 /// let b = ml_matrix("1 2;3 4");
1264 /// let c = ml_matrix("1 2 3;4 1 2;7 3 4");
1265 /// a.subs_mat((1,1), (2,2), &b);
1266 /// assert_eq!(a, c);
1267 /// }
1268 /// ```
1269 fn subs_mat(&mut self, start: (usize, usize), end: (usize, usize), m: &Matrix) {
1270 let row = end.0 - start.0 + 1;
1271 let col = end.1 - start.1 + 1;
1272 for i in 0..row {
1273 for j in 0..col {
1274 self[(start.0 + i, start.1 + j)] = m[(i, j)];
1275 }
1276 }
1277 }
1278}
1279
1280impl Matrix {
1281 /// Write to CSV
1282 ///
1283 /// # Examples
1284 /// ```
1285 /// #[macro_use]
1286 /// extern crate peroxide;
1287 /// use peroxide::fuga::*;
1288 ///
1289 /// fn main() -> Result<(), Box<dyn Error>> {
1290 /// let a = matrix(c!(1,2,3,3,2,1), 3, 2, Col);
1291 /// a.write("example_data/test.csv")?;
1292 ///
1293 /// Ok(())
1294 /// }
1295 /// ```
1296 #[cfg(feature = "csv")]
1297 pub fn write(&self, file_path: &str) -> Result<(), Box<dyn Error>> {
1298 let mut wtr = WriterBuilder::new().from_path(file_path)?;
1299 let r = self.row;
1300 let c = self.col;
1301 for i in 0..r {
1302 let mut record: Vec<String> = vec!["".to_string(); c];
1303 for j in 0..c {
1304 record[j] = self[(i, j)].to_string();
1305 }
1306 wtr.write_record(record)?;
1307 }
1308 wtr.flush()?;
1309 Ok(())
1310 }
1311
1312 /// Write to CSV (with round option)
1313 ///
1314 /// # Examples
1315 /// ```
1316 /// #[macro_use]
1317 /// extern crate peroxide;
1318 /// use peroxide::fuga::*;
1319 ///
1320 /// fn main() -> Result<(), Box<dyn Error>> {
1321 /// let a = matrix(c!(1,2,3,3,2,1), 3, 2, Col);
1322 /// a.write_round("example_data/test.csv", 0)?;
1323 ///
1324 /// Ok(())
1325 /// }
1326 /// ```
1327 #[cfg(feature = "csv")]
1328 pub fn write_round(&self, file_path: &str, round: usize) -> Result<(), Box<dyn Error>> {
1329 let mut wtr = WriterBuilder::new().from_path(file_path)?;
1330 let r = self.row;
1331 let c = self.col;
1332 for i in 0..r {
1333 let mut record: Vec<String> = vec!["".to_string(); c];
1334 for j in 0..c {
1335 record[j] = format!("{:.*}", round, self[(i, j)]);
1336 }
1337 wtr.write_record(record)?;
1338 }
1339 wtr.flush()?;
1340 Ok(())
1341 }
1342
1343 #[cfg(feature = "csv")]
1344 pub fn write_with_header(
1345 &self,
1346 file_path: &str,
1347 header: Vec<&str>,
1348 ) -> Result<(), Box<dyn Error>> {
1349 let mut wtr = WriterBuilder::new().from_path(file_path)?;
1350 let r = self.row;
1351 let c = self.col;
1352 assert_eq!(c, header.len());
1353 wtr.write_record(header)?;
1354 for i in 0..r {
1355 let mut record: Vec<String> = vec!["".to_string(); c];
1356 for j in 0..c {
1357 record[j] = self[(i, j)].to_string();
1358 }
1359 wtr.write_record(record)?;
1360 }
1361 wtr.flush()?;
1362 Ok(())
1363 }
1364
1365 #[cfg(feature = "csv")]
1366 pub fn write_with_header_round(
1367 &self,
1368 file_path: &str,
1369 header: Vec<&str>,
1370 round: usize,
1371 ) -> Result<(), Box<dyn Error>> {
1372 let mut wtr = WriterBuilder::new().from_path(file_path)?;
1373 let r = self.row;
1374 let c = self.col;
1375 assert_eq!(c, header.len());
1376 wtr.write_record(header)?;
1377 for i in 0..r {
1378 let mut record: Vec<String> = vec!["".to_string(); c];
1379 for j in 0..c {
1380 record[j] = format!("{:.*}", round, self[(i, j)]);
1381 }
1382 wtr.write_record(record)?;
1383 }
1384 wtr.flush()?;
1385 Ok(())
1386 }
1387
1388 /// Read from CSV
1389 ///
1390 /// # Examples
1391 /// ```
1392 /// #[macro_use]
1393 /// extern crate peroxide;
1394 /// use peroxide::fuga::*;
1395 /// use std::process;
1396 ///
1397 /// fn main() -> Result<(), Box<dyn Error>> {
1398 /// let a = matrix(c!(1,2,3,3,2,1), 3, 2, Col);
1399 /// a.write_round("example_data/test.csv", 0)?;
1400 ///
1401 /// let b = Matrix::read("example_data/test.csv", false, ','); // header = false, delimiter = ','
1402 /// match b {
1403 /// Ok(mat) => println!("{}", mat),
1404 /// Err(err) => {
1405 /// println!("{}", err);
1406 /// process::exit(1);
1407 /// }
1408 /// }
1409 ///
1410 /// Ok(())
1411 /// }
1412 /// ```
1413 #[cfg(feature = "csv")]
1414 pub fn read(file_path: &str, header: bool, delimiter: char) -> Result<Matrix, Box<dyn Error>> {
1415 let mut rdr = ReaderBuilder::new()
1416 .has_headers(header)
1417 .delimiter(delimiter as u8)
1418 .from_path(file_path)?;
1419
1420 let records = rdr
1421 .records()
1422 .collect::<Result<Vec<StringRecord>, csv::Error>>()?;
1423 let result = records;
1424
1425 let l_row = result.len();
1426 let l_col = result[0].len();
1427
1428 let mut m = matrix(vec![0f64; l_row * l_col], l_row, l_col, Row);
1429
1430 for i in 0..l_row {
1431 for j in 0..l_col {
1432 m[(i, j)] = result[i][j].parse::<f64>().unwrap();
1433 }
1434 }
1435
1436 Ok(m)
1437 }
1438 /// Matrix from series
1439 ///
1440 /// # Example
1441 /// ```rust
1442 /// #[macro_use]
1443 /// extern crate peroxide;
1444 /// use peroxide::fuga::*;
1445 ///
1446 /// fn main() {
1447 /// let a = Series::new(c!(1,2,3,4));
1448 /// let m = Matrix::from_series(&a, 2, 2, Row);
1449 ///
1450 /// assert_eq!(m, matrix(c!(1,2,3,4), 2, 2, Row));
1451 /// }
1452 /// ```
1453 pub fn from_series(series: &Series, row: usize, col: usize, shape: Shape) -> Self {
1454 let v: Vec<f64> = series.to_vec();
1455 matrix(v, row, col, shape)
1456 }
1457
1458 /// From index operations in parallel
1459 #[cfg(feature = "parallel")]
1460 fn par_from_index<F, G>(f: F, size: (usize, usize)) -> Matrix
1461 where
1462 F: Fn(usize, usize) -> G + Copy + Send + Sync,
1463 G: Into<f64>,
1464 {
1465 let row = size.0;
1466 let col = size.1;
1467
1468 let data = (0..row)
1469 .into_par_iter()
1470 .flat_map(|i| {
1471 (0..col)
1472 .into_par_iter()
1473 .map(|j| f(i, j).into())
1474 .collect::<Vec<f64>>()
1475 })
1476 .collect::<Vec<f64>>();
1477 matrix(data, row, col, Row)
1478 }
1479}
1480
1481// =============================================================================
1482// Mathematics for Matrix
1483// =============================================================================
1484impl Vector for Matrix {
1485 type Scalar = f64;
1486
1487 fn add_vec(&self, other: &Self) -> Self {
1488 assert_eq!(self.row, other.row);
1489 assert_eq!(self.col, other.col);
1490
1491 match () {
1492 #[cfg(feature = "O3")]
1493 () => {
1494 if self.shape != other.shape {
1495 return self.add(&other.change_shape());
1496 }
1497 let x = &self.data;
1498 let mut y = other.data.clone();
1499 let n_i32 = x.len() as i32;
1500 let a_f64 = 1f64;
1501 unsafe {
1502 daxpy(n_i32, a_f64, x, 1, &mut y, 1);
1503 }
1504 matrix(y, self.row, self.col, self.shape)
1505 }
1506 _ => {
1507 let mut result = matrix(self.data.clone(), self.row, self.col, self.shape);
1508 for i in 0..self.row {
1509 for j in 0..self.col {
1510 result[(i, j)] += other[(i, j)];
1511 }
1512 }
1513 result
1514 }
1515 }
1516 }
1517
1518 fn sub_vec(&self, other: &Self) -> Self {
1519 assert_eq!(self.row, other.row);
1520 assert_eq!(self.col, other.col);
1521 match () {
1522 #[cfg(feature = "O3")]
1523 () => {
1524 if self.shape != other.shape {
1525 return self.sub(&other.change_shape());
1526 }
1527 let x = &other.data;
1528 let mut y = self.data.clone();
1529 let n_i32 = x.len() as i32;
1530 let a_f64 = -1f64;
1531 unsafe {
1532 daxpy(n_i32, a_f64, x, 1, &mut y, 1);
1533 }
1534 matrix(y, self.row, self.col, self.shape)
1535 }
1536 _ => {
1537 let mut result = matrix(self.data.clone(), self.row, self.col, self.shape);
1538 for i in 0..self.row {
1539 for j in 0..self.col {
1540 result[(i, j)] -= other[(i, j)];
1541 }
1542 }
1543 result
1544 }
1545 }
1546 }
1547
1548 fn mul_scalar(&self, other: Self::Scalar) -> Self {
1549 match () {
1550 #[cfg(feature = "O3")]
1551 () => {
1552 let x = &self.data;
1553 let mut y = vec![0f64; x.len()];
1554 let a_f64 = other;
1555 let n_i32 = x.len() as i32;
1556
1557 unsafe {
1558 daxpy(n_i32, a_f64, x, 1, &mut y, 1);
1559 }
1560 matrix(y, self.row, self.col, self.shape)
1561 }
1562 _ => {
1563 let scalar = other;
1564 self.fmap(|x| x * scalar)
1565 }
1566 }
1567 }
1568}
1569
1570impl Normed for Matrix {
1571 type UnsignedScalar = f64;
1572
1573 /// Norm of Matrix
1574 ///
1575 /// # Example
1576 /// ```
1577 /// use peroxide::fuga::*;
1578 ///
1579 /// let a = ml_matrix("1 2;3 4");
1580 /// assert_eq!(a.norm(Norm::F), (1f64 + 4f64 + 9f64 + 16f64).sqrt());
1581 /// ```
1582 fn norm(&self, kind: Norm) -> f64 {
1583 match kind {
1584 Norm::F => {
1585 let mut s = 0f64;
1586 for i in 0..self.data.len() {
1587 s += self.data[i].powi(2);
1588 }
1589 s.sqrt()
1590 }
1591 Norm::Lpq(p, q) => {
1592 let mut s = 0f64;
1593 for j in 0..self.col {
1594 let mut s_row = 0f64;
1595 for i in 0..self.row {
1596 s_row += self[(i, j)].powi(p as i32);
1597 }
1598 s += s_row.powf(q / p);
1599 }
1600 s.powf(1f64 / q)
1601 }
1602 Norm::L1 => {
1603 let mut m = f64::MIN;
1604 match self.shape {
1605 Row => self.change_shape().norm(Norm::L1),
1606 Col => {
1607 for c in 0..self.col {
1608 let s = self.col(c).iter().sum();
1609 if s > m {
1610 m = s;
1611 }
1612 }
1613 m
1614 }
1615 }
1616 }
1617 Norm::LInf => {
1618 let mut m = f64::MIN;
1619 match self.shape {
1620 Col => self.change_shape().norm(Norm::LInf),
1621 Row => {
1622 for r in 0..self.row {
1623 let s = self.row(r).iter().sum();
1624 if s > m {
1625 m = s;
1626 }
1627 }
1628 m
1629 }
1630 }
1631 }
1632 Norm::L2 => {
1633 let a = &self.t() * self;
1634 let eig = eigen(&a, EigenMethod::Jacobi);
1635 eig.eigenvalue.norm(Norm::LInf)
1636 }
1637 Norm::Lp(_) => unimplemented!(),
1638 }
1639 }
1640 fn normalize(&self, _kind: Norm) -> Self
1641 where
1642 Self: Sized,
1643 {
1644 unimplemented!()
1645 }
1646}
1647
1648#[cfg(feature = "parallel")]
1649impl ParallelNormed for Matrix {
1650 type UnsignedScalar = f64;
1651
1652 /// Parallel version of norm
1653 ///
1654 /// # Example
1655 /// ```
1656 /// use peroxide::fuga::*;
1657 ///
1658 /// let a = ml_matrix("1 2;3 4");
1659 /// assert_eq!(a.par_norm(Norm::F), (1f64 + 4f64 + 9f64 + 16f64).sqrt());
1660 /// ```
1661 ///
1662 fn par_norm(&self, kind: Norm) -> f64 {
1663 match kind {
1664 Norm::F => {
1665 let mut s = 0f64;
1666 for i in 0..self.data.len() {
1667 s += self.data[i].powi(2);
1668 }
1669 s.sqrt()
1670 }
1671 Norm::Lpq(p, q) => {
1672 let s = (0..self.col)
1673 .into_par_iter()
1674 .map(|j| {
1675 let s_row = (0..self.row)
1676 .into_iter()
1677 .map(|i| self[(i, j)].powi(p as i32))
1678 .sum::<f64>();
1679 s_row.powf(q / p)
1680 })
1681 .sum::<f64>();
1682 s.powf(1f64 / q)
1683 }
1684 Norm::L1 => {
1685 let mut m = f64::MIN;
1686 match self.shape {
1687 Row => self.change_shape().norm(Norm::L1),
1688 Col => {
1689 for c in 0..self.col {
1690 let s = self.col(c).par_iter().sum();
1691 if s > m {
1692 m = s;
1693 }
1694 }
1695 m
1696 }
1697 }
1698 }
1699 Norm::LInf => {
1700 let mut m = f64::MIN;
1701 match self.shape {
1702 Col => self.change_shape().norm(Norm::LInf),
1703 Row => {
1704 for r in 0..self.row {
1705 let s = self.row(r).iter().sum();
1706 if s > m {
1707 m = s;
1708 }
1709 }
1710 m
1711 }
1712 }
1713 }
1714 Norm::L2 => {
1715 let a = &self.t() * self;
1716 let eig = eigen(&a, EigenMethod::Jacobi);
1717 eig.eigenvalue.norm(Norm::LInf)
1718 }
1719 Norm::Lp(_) => unimplemented!(),
1720 }
1721 }
1722}
1723
1724/// Frobenius inner product
1725impl InnerProduct for Matrix {
1726 fn dot(&self, rhs: &Self) -> f64 {
1727 if self.shape == rhs.shape {
1728 self.data.dot(&rhs.data)
1729 } else {
1730 self.data.dot(&rhs.change_shape().data)
1731 }
1732 }
1733}
1734
1735/// Frobenius inner product in parallel
1736#[cfg(feature = "parallel")]
1737impl ParallelInnerProduct for Matrix {
1738 fn par_dot(&self, rhs: &Self) -> f64 {
1739 if self.shape == rhs.shape {
1740 self.data.par_dot(&rhs.data)
1741 } else {
1742 self.data.par_dot(&rhs.change_shape().data)
1743 }
1744 }
1745}
1746
1747/// Matrix as Linear operator for Vector
1748#[allow(non_snake_case)]
1749impl LinearOp<Vec<f64>, Vec<f64>> for Matrix {
1750 fn apply(&self, other: &Vec<f64>) -> Vec<f64> {
1751 match () {
1752 #[cfg(feature = "O3")]
1753 () => {
1754 let x = other;
1755 let mut y = vec![0f64; self.row];
1756 let A = &self.data;
1757 let m_i32 = self.row as i32;
1758 let n_i32 = self.col as i32;
1759 match self.shape {
1760 Row => unsafe {
1761 dgemv(b'T', n_i32, m_i32, 1f64, A, n_i32, x, 1, 0f64, &mut y, 1);
1762 },
1763 Col => unsafe {
1764 dgemv(b'N', m_i32, n_i32, 1f64, A, m_i32, x, 1, 0f64, &mut y, 1);
1765 },
1766 }
1767 y
1768 }
1769 _ => {
1770 assert_eq!(self.col, other.len());
1771 let mut c = vec![0f64; self.row];
1772 gemv(1f64, self, other, 0f64, &mut c);
1773 c
1774 }
1775 }
1776 }
1777}
1778
1779impl MatrixProduct for Matrix {
1780 fn kronecker(&self, other: &Self) -> Self {
1781 let r1 = self.row;
1782 let c1 = self.col;
1783
1784 let mut result = self[(0, 0)] * other;
1785
1786 for j in 1..c1 {
1787 let n = self[(0, j)] * other;
1788 result = cbind(result, n).unwrap();
1789 }
1790
1791 for i in 1..r1 {
1792 let mut m = self[(i, 0)] * other;
1793 for j in 1..c1 {
1794 let n = self[(i, j)] * other;
1795 m = cbind(m, n).unwrap();
1796 }
1797 result = rbind(result, m).unwrap();
1798 }
1799 result
1800 }
1801
1802 fn hadamard(&self, other: &Self) -> Matrix {
1803 assert_eq!(self.row, other.row);
1804 assert_eq!(self.col, other.col);
1805
1806 let r = self.row;
1807 let c = self.col;
1808
1809 let mut m = matrix(vec![0f64; r * c], r, c, self.shape);
1810 for i in 0..r {
1811 for j in 0..c {
1812 m[(i, j)] = self[(i, j)] * other[(i, j)]
1813 }
1814 }
1815 m
1816 }
1817}
1818
1819// =============================================================================
1820// Common Properties of Matrix & Vec<f64>
1821// =============================================================================
1822/// `Matrix` to `Vec<f64>`
1823impl Into<Vec<f64>> for Matrix {
1824 fn into(self) -> Vec<f64> {
1825 self.data
1826 }
1827}
1828
1829/// `&Matrix` to `&Vec<f64>`
1830impl<'a> Into<&'a Vec<f64>> for &'a Matrix {
1831 fn into(self) -> &'a Vec<f64> {
1832 &self.data
1833 }
1834}
1835
1836/// `Vec<f64>` to `Matrix`
1837impl Into<Matrix> for Vec<f64> {
1838 fn into(self) -> Matrix {
1839 let l = self.len();
1840 matrix(self, l, 1, Col)
1841 }
1842}
1843
1844impl Into<Matrix> for &Vec<f64> {
1845 fn into(self) -> Matrix {
1846 let l = self.len();
1847 matrix(self.clone(), l, 1, Col)
1848 }
1849}
1850
1851// =============================================================================
1852// Standard Operation for Matrix (ADD)
1853// =============================================================================
1854
1855/// Element-wise addition of Matrix
1856///
1857/// # Caution
1858/// > You should remember ownership.
1859/// > If you use Matrix `a,b` then you can't use them after.
1860impl Add<Matrix> for Matrix {
1861 type Output = Self;
1862
1863 fn add(self, other: Self) -> Self {
1864 assert_eq!(&self.row, &other.row);
1865 assert_eq!(&self.col, &other.col);
1866
1867 match () {
1868 #[cfg(feature = "O3")]
1869 () => {
1870 if self.shape != other.shape {
1871 return self.add(other.change_shape());
1872 }
1873 let x = &self.data;
1874 let mut y = other.data.clone();
1875 let n_i32 = x.len() as i32;
1876 let a_f64 = 1f64;
1877 unsafe {
1878 daxpy(n_i32, a_f64, x, 1, &mut y, 1);
1879 }
1880 matrix(y, self.row, self.col, self.shape)
1881 }
1882 _ => {
1883 let mut result = matrix(self.data.clone(), self.row, self.col, self.shape);
1884 for i in 0..self.row {
1885 for j in 0..self.col {
1886 result[(i, j)] += other[(i, j)];
1887 }
1888 }
1889 result
1890 }
1891 }
1892 }
1893}
1894
1895impl<'a, 'b> Add<&'b Matrix> for &'a Matrix {
1896 type Output = Matrix;
1897
1898 fn add(self, rhs: &'b Matrix) -> Self::Output {
1899 self.add_vec(rhs)
1900 }
1901}
1902
1903/// Element-wise addition between Matrix & f64
1904///
1905/// # Examples
1906/// ```
1907/// #[macro_use]
1908/// extern crate peroxide;
1909/// use peroxide::fuga::*;
1910///
1911/// fn main() {
1912/// let a = matrix!(1;4;1, 2, 2, Row);
1913/// assert_eq!(a + 1, matrix!(2;5;1, 2, 2, Row));
1914/// }
1915/// ```
1916impl<T> Add<T> for Matrix
1917where
1918 T: Into<f64> + Copy,
1919{
1920 type Output = Self;
1921 fn add(self, other: T) -> Self {
1922 match () {
1923 #[cfg(feature = "O3")]
1924 () => {
1925 let x = &self.data;
1926 let mut y = vec![other.into(); x.len()];
1927 let n_i32 = x.len() as i32;
1928 unsafe {
1929 daxpy(n_i32, 1f64, x, 1, &mut y, 1);
1930 }
1931 matrix(y, self.row, self.col, self.shape)
1932 }
1933 _ => self.fmap(|x| x + other.into()),
1934 }
1935 }
1936}
1937
1938/// Element-wise addition between &Matrix & f64
1939impl<'a, T> Add<T> for &'a Matrix
1940where
1941 T: Into<f64> + Copy,
1942{
1943 type Output = Matrix;
1944 fn add(self, other: T) -> Self::Output {
1945 match () {
1946 #[cfg(feature = "O3")]
1947 () => {
1948 let x = &self.data;
1949 let mut y = vec![other.into(); x.len()];
1950 let n_i32 = x.len() as i32;
1951 unsafe {
1952 daxpy(n_i32, 1f64, x, 1, &mut y, 1);
1953 }
1954 matrix(y, self.row, self.col, self.shape)
1955 }
1956 _ => self.fmap(|x| x + other.into()),
1957 }
1958 }
1959}
1960
1961/// Element-wise addition between f64 & matrix
1962///
1963/// # Examples
1964/// ```
1965/// #[macro_use]
1966/// extern crate peroxide;
1967/// use peroxide::fuga::*;
1968///
1969/// fn main() {
1970/// let a = matrix!(1;4;1, 2, 2, Row);
1971/// assert_eq!(1f64 + a, matrix!(2;5;1, 2, 2, Row));
1972/// }
1973/// ```
1974impl Add<Matrix> for f64 {
1975 type Output = Matrix;
1976
1977 fn add(self, other: Matrix) -> Matrix {
1978 other.add(self)
1979 }
1980}
1981
1982/// Element-wise addition between f64 & &Matrix
1983impl<'a> Add<&'a Matrix> for f64 {
1984 type Output = Matrix;
1985
1986 fn add(self, other: &'a Matrix) -> Self::Output {
1987 other.add(self)
1988 }
1989}
1990
1991/// Element-wise addition between i32 & Matrix
1992///
1993/// # Examples
1994/// ```
1995/// #[macro_use]
1996/// extern crate peroxide;
1997/// use peroxide::fuga::*;
1998///
1999/// fn main() {
2000/// let a = matrix!(1;4;1, 2, 2, Row);
2001/// assert_eq!(1 + a, matrix!(2;5;1, 2, 2, Row));
2002/// }
2003/// ```
2004impl Add<Matrix> for i32 {
2005 type Output = Matrix;
2006
2007 fn add(self, other: Matrix) -> Matrix {
2008 other.add(self)
2009 }
2010}
2011
2012/// Element-wise addition between i32 & &Matrix
2013impl<'a> Add<&'a Matrix> for i32 {
2014 type Output = Matrix;
2015
2016 fn add(self, other: &'a Matrix) -> Self::Output {
2017 other.add(self)
2018 }
2019}
2020
2021/// Element-wise addition between usize & matrix
2022///
2023/// # Examples
2024/// ```
2025/// #[macro_use]
2026/// extern crate peroxide;
2027/// use peroxide::fuga::*;
2028///
2029/// fn main() {
2030/// let a = matrix!(1;4;1, 2, 2, Row);
2031/// assert_eq!(1usize + a, matrix!(2;5;1, 2, 2, Row));
2032/// }
2033/// ```
2034impl Add<Matrix> for usize {
2035 type Output = Matrix;
2036
2037 fn add(self, other: Matrix) -> Matrix {
2038 other.add(self as f64)
2039 }
2040}
2041
2042/// Element-wise addition between usize & &Matrix
2043impl<'a> Add<&'a Matrix> for usize {
2044 type Output = Matrix;
2045
2046 fn add(self, other: &'a Matrix) -> Self::Output {
2047 other.add(self as f64)
2048 }
2049}
2050
2051// =============================================================================
2052// Standard Operation for Matrix (Neg)
2053// =============================================================================
2054/// Negation of Matrix
2055///
2056/// # Examples
2057/// ```
2058/// extern crate peroxide;
2059/// use peroxide::fuga::*;
2060///
2061/// let a = matrix(vec![1,2,3,4],2,2,Row);
2062/// println!("{}", -a); // [[-1,-2],[-3,-4]]
2063/// ```
2064impl Neg for Matrix {
2065 type Output = Self;
2066
2067 fn neg(self) -> Self {
2068 match () {
2069 #[cfg(feature = "O3")]
2070 () => {
2071 let x = &self.data;
2072 let mut y = vec![0f64; x.len()];
2073 let n_i32 = x.len() as i32;
2074 unsafe {
2075 daxpy(n_i32, -1f64, x, 1, &mut y, 1);
2076 }
2077 matrix(y, self.row, self.col, self.shape)
2078 }
2079 _ => matrix(
2080 self.data.into_iter().map(|x: f64| -x).collect::<Vec<f64>>(),
2081 self.row,
2082 self.col,
2083 self.shape,
2084 ),
2085 }
2086 }
2087}
2088
2089/// Negation of &'a Matrix
2090impl<'a> Neg for &'a Matrix {
2091 type Output = Matrix;
2092
2093 fn neg(self) -> Self::Output {
2094 match () {
2095 #[cfg(feature = "O3")]
2096 () => {
2097 let x = &self.data;
2098 let mut y = vec![0f64; x.len()];
2099 let n_i32 = x.len() as i32;
2100 unsafe {
2101 daxpy(n_i32, -1f64, x, 1, &mut y, 1);
2102 }
2103 matrix(y, self.row, self.col, self.shape)
2104 }
2105 _ => matrix(
2106 self.data
2107 .clone()
2108 .into_iter()
2109 .map(|x: f64| -x)
2110 .collect::<Vec<f64>>(),
2111 self.row,
2112 self.col,
2113 self.shape,
2114 ),
2115 }
2116 }
2117}
2118
2119// =============================================================================
2120// Standard Operation for Matrix (Sub)
2121// =============================================================================
2122/// Subtraction between Matrix
2123///
2124/// # Examples
2125/// ```
2126/// extern crate peroxide;
2127/// use peroxide::fuga::*;
2128///
2129/// let a = matrix(vec![1,2,3,4],2,2,Row);
2130/// let b = matrix(vec![1,2,3,4],2,2,Col);
2131/// println!("{}", a - b); // [[0, -1], [1, 0]]
2132/// ```
2133impl Sub<Matrix> for Matrix {
2134 type Output = Self;
2135
2136 fn sub(self, other: Self) -> Self {
2137 assert_eq!(&self.row, &other.row);
2138 assert_eq!(&self.col, &other.col);
2139 match () {
2140 #[cfg(feature = "O3")]
2141 () => {
2142 if self.shape != other.shape {
2143 return self.sub(other.change_shape());
2144 }
2145 let x = &other.data;
2146 let mut y = self.data.clone();
2147 let n_i32 = x.len() as i32;
2148 let a_f64 = -1f64;
2149 unsafe {
2150 daxpy(n_i32, a_f64, x, 1, &mut y, 1);
2151 }
2152 matrix(y, self.row, self.col, self.shape)
2153 }
2154 _ => {
2155 let mut result = matrix(self.data.clone(), self.row, self.col, self.shape);
2156 for i in 0..self.row {
2157 for j in 0..self.col {
2158 result[(i, j)] -= other[(i, j)];
2159 }
2160 }
2161 result
2162 }
2163 }
2164 }
2165}
2166
2167impl<'a, 'b> Sub<&'b Matrix> for &'a Matrix {
2168 type Output = Matrix;
2169
2170 fn sub(self, rhs: &'b Matrix) -> Matrix {
2171 self.sub_vec(rhs)
2172 }
2173}
2174
2175/// Subtraction between Matrix & f64
2176impl<T> Sub<T> for Matrix
2177where
2178 T: Into<f64> + Copy,
2179{
2180 type Output = Self;
2181
2182 fn sub(self, other: T) -> Self {
2183 match () {
2184 #[cfg(feature = "O3")]
2185 () => {
2186 let mut y = self.data;
2187 let x = vec![other.into(); y.len()];
2188 let n_i32 = y.len() as i32;
2189 unsafe {
2190 daxpy(n_i32, -1f64, &x, 1, &mut y, 1);
2191 }
2192 matrix(y, self.row, self.col, self.shape)
2193 }
2194 _ => self.fmap(|x| x - other.into()),
2195 }
2196 }
2197}
2198
2199/// Subtraction between &Matrix & f64
2200impl<'a, T> Sub<T> for &'a Matrix
2201where
2202 T: Into<f64> + Copy,
2203{
2204 type Output = Matrix;
2205
2206 fn sub(self, other: T) -> Self::Output {
2207 match () {
2208 #[cfg(feature = "O3")]
2209 () => {
2210 let mut y = self.data.clone();
2211 let x = vec![other.into(); y.len()];
2212 let n_i32 = y.len() as i32;
2213 unsafe {
2214 daxpy(n_i32, -1f64, &x, 1, &mut y, 1);
2215 }
2216 matrix(y, self.row, self.col, self.shape)
2217 }
2218 _ => self.fmap(|x| x - other.into()),
2219 }
2220 }
2221}
2222
2223/// Subtraction Matrix with f64
2224///
2225/// # Examples
2226/// ```
2227/// #[macro_use]
2228/// extern crate peroxide;
2229/// use peroxide::fuga::*;
2230///
2231/// fn main() {
2232/// let a = matrix(vec![1,2,3,4],2,2,Row);
2233/// assert_eq!(a - 1f64, matrix!(0;3;1, 2, 2, Row));
2234/// }
2235/// ```
2236impl Sub<Matrix> for f64 {
2237 type Output = Matrix;
2238
2239 fn sub(self, other: Matrix) -> Matrix {
2240 -other.sub(self)
2241 }
2242}
2243
2244impl<'a> Sub<&'a Matrix> for f64 {
2245 type Output = Matrix;
2246
2247 fn sub(self, other: &'a Matrix) -> Self::Output {
2248 -other.sub(self)
2249 }
2250}
2251
2252impl Sub<Matrix> for i32 {
2253 type Output = Matrix;
2254
2255 fn sub(self, other: Matrix) -> Matrix {
2256 -other.sub(self)
2257 }
2258}
2259
2260impl<'a> Sub<&'a Matrix> for i32 {
2261 type Output = Matrix;
2262
2263 fn sub(self, other: &'a Matrix) -> Self::Output {
2264 -other.sub(self)
2265 }
2266}
2267
2268impl Sub<Matrix> for usize {
2269 type Output = Matrix;
2270
2271 fn sub(self, other: Matrix) -> Matrix {
2272 -other.sub(self as f64)
2273 }
2274}
2275
2276impl<'a> Sub<&'a Matrix> for usize {
2277 type Output = Matrix;
2278
2279 fn sub(self, other: &'a Matrix) -> Self::Output {
2280 -other.sub(self as f64)
2281 }
2282}
2283
2284// =============================================================================
2285// Multiplication for Matrix
2286// =============================================================================
2287/// Element-wise multiplication between Matrix vs f64
2288impl Mul<f64> for Matrix {
2289 type Output = Self;
2290
2291 fn mul(self, other: f64) -> Self {
2292 match () {
2293 #[cfg(feature = "O3")]
2294 () => {
2295 let x = &self.data;
2296 let mut y = vec![0f64; x.len()];
2297 let a_f64 = other;
2298 let n_i32 = x.len() as i32;
2299
2300 unsafe {
2301 daxpy(n_i32, a_f64, x, 1, &mut y, 1);
2302 }
2303 matrix(y, self.row, self.col, self.shape)
2304 }
2305 _ => self.fmap(|x| x * other),
2306 }
2307 }
2308}
2309
2310impl Mul<i64> for Matrix {
2311 type Output = Self;
2312
2313 fn mul(self, other: i64) -> Self {
2314 self.mul(other as f64)
2315 }
2316}
2317
2318impl Mul<i32> for Matrix {
2319 type Output = Self;
2320
2321 fn mul(self, other: i32) -> Self {
2322 self.mul(other as f64)
2323 }
2324}
2325
2326impl Mul<usize> for Matrix {
2327 type Output = Self;
2328
2329 fn mul(self, other: usize) -> Self {
2330 self.mul(other as f64)
2331 }
2332}
2333
2334// impl<'a> Mul<i64> for &'a Matrix {
2335// type Output = Matrix;
2336//
2337// fn mul(self, other: i64) -> Self::Output {
2338// self.mul(other as f64)
2339// }
2340// }
2341//
2342// impl<'a> Mul<i32> for &'a Matrix {
2343// type Output = Matrix;
2344//
2345// fn mul(self, other: i32) -> Self::Output {
2346// self.mul(other as f64)
2347// }
2348// }
2349//
2350// impl<'a> Mul<usize> for &'a Matrix {
2351// type Output = Matrix;
2352//
2353// fn mul(self, other: usize) -> Self::Output {
2354// self.mul(other as f64)
2355// }
2356// }
2357
2358impl Mul<Matrix> for f64 {
2359 type Output = Matrix;
2360
2361 fn mul(self, other: Matrix) -> Matrix {
2362 other.mul(self)
2363 }
2364}
2365
2366impl Mul<Matrix> for i64 {
2367 type Output = Matrix;
2368
2369 fn mul(self, other: Matrix) -> Matrix {
2370 other.mul(self as f64)
2371 }
2372}
2373
2374impl Mul<Matrix> for i32 {
2375 type Output = Matrix;
2376
2377 fn mul(self, other: Matrix) -> Matrix {
2378 other.mul(self)
2379 }
2380}
2381
2382impl Mul<Matrix> for usize {
2383 type Output = Matrix;
2384
2385 fn mul(self, other: Matrix) -> Matrix {
2386 other.mul(self as f64)
2387 }
2388}
2389
2390impl<'a> Mul<&'a Matrix> for f64 {
2391 type Output = Matrix;
2392
2393 fn mul(self, other: &'a Matrix) -> Matrix {
2394 other.mul_scalar(self)
2395 }
2396}
2397
2398impl<'a> Mul<&'a Matrix> for i64 {
2399 type Output = Matrix;
2400
2401 fn mul(self, other: &'a Matrix) -> Matrix {
2402 other.mul_scalar(self as f64)
2403 }
2404}
2405
2406impl<'a> Mul<&'a Matrix> for i32 {
2407 type Output = Matrix;
2408
2409 fn mul(self, other: &'a Matrix) -> Matrix {
2410 other.mul_scalar(self as f64)
2411 }
2412}
2413
2414impl<'a> Mul<&'a Matrix> for usize {
2415 type Output = Matrix;
2416
2417 fn mul(self, other: &'a Matrix) -> Matrix {
2418 other.mul_scalar(self as f64)
2419 }
2420}
2421
2422/// Matrix Multiplication
2423///
2424/// # Examples
2425/// ```
2426/// #[macro_use]
2427/// extern crate peroxide;
2428/// use peroxide::fuga::*;
2429///
2430/// fn main() {
2431/// let a = matrix!(1;4;1, 2, 2, Row);
2432/// let b = matrix!(1;4;1, 2, 2, Col);
2433/// assert_eq!(a * b, matrix(c!(5, 11, 11, 25), 2, 2, Row));
2434/// }
2435/// ```
2436impl Mul<Matrix> for Matrix {
2437 type Output = Self;
2438
2439 fn mul(self, other: Self) -> Self {
2440 match () {
2441 #[cfg(feature = "O3")]
2442 () => blas_mul(&self, &other),
2443 _ => matmul(&self, &other),
2444 }
2445 }
2446}
2447
2448impl<'a, 'b> Mul<&'b Matrix> for &'a Matrix {
2449 type Output = Matrix;
2450
2451 fn mul(self, other: &'b Matrix) -> Self::Output {
2452 match () {
2453 #[cfg(feature = "O3")]
2454 () => blas_mul(self, other),
2455 _ => matmul(self, other),
2456 }
2457 }
2458}
2459
2460#[allow(non_snake_case)]
2461impl Mul<Vec<f64>> for Matrix {
2462 type Output = Vec<f64>;
2463
2464 fn mul(self, other: Vec<f64>) -> Self::Output {
2465 self.apply(&other)
2466 }
2467}
2468
2469#[allow(non_snake_case)]
2470impl<'a, 'b> Mul<&'b Vec<f64>> for &'a Matrix {
2471 type Output = Vec<f64>;
2472
2473 fn mul(self, other: &'b Vec<f64>) -> Self::Output {
2474 self.apply(other)
2475 }
2476}
2477
2478/// Matrix multiplication for `Vec<f64>` vs `Matrix`
2479///
2480/// # Examples
2481/// ```
2482/// #[macro_use]
2483/// extern crate peroxide;
2484/// use peroxide::fuga::*;
2485///
2486/// fn main() {
2487/// let a = matrix!(1;4;1, 2, 2, Row);
2488/// let v = c!(1,2);
2489/// assert_eq!(v * a, c!(7, 10));
2490/// }
2491/// ```
2492impl Mul<Matrix> for Vec<f64> {
2493 type Output = Vec<f64>;
2494
2495 fn mul(self, other: Matrix) -> Self::Output {
2496 assert_eq!(self.len(), other.row);
2497 let mut c = vec![0f64; other.col];
2498 gevm(1f64, &self, &other, 0f64, &mut c);
2499 c
2500 }
2501}
2502
2503impl<'a, 'b> Mul<&'b Matrix> for &'a Vec<f64> {
2504 type Output = Vec<f64>;
2505
2506 fn mul(self, other: &'b Matrix) -> Self::Output {
2507 assert_eq!(self.len(), other.row);
2508 let mut c = vec![0f64; other.col];
2509 gevm(1f64, self, other, 0f64, &mut c);
2510 c
2511 }
2512}
2513
2514// =============================================================================
2515// Standard Operation for Matrix (DIV)
2516// =============================================================================
2517/// Element-wise division between Matrix vs f64
2518impl Div<f64> for Matrix {
2519 type Output = Self;
2520
2521 fn div(self, other: f64) -> Self {
2522 match () {
2523 #[cfg(feature = "O3")]
2524 () => {
2525 let x = &self.data;
2526 let mut y = vec![0f64; x.len()];
2527 let a_f64 = other;
2528 let n_i32 = x.len() as i32;
2529
2530 unsafe {
2531 daxpy(n_i32, 1f64 / a_f64, x, 1, &mut y, 1);
2532 }
2533 matrix(y, self.row, self.col, self.shape)
2534 }
2535 _ => self.fmap(|x| x / other),
2536 }
2537 }
2538}
2539
2540impl Div<i64> for Matrix {
2541 type Output = Self;
2542
2543 fn div(self, other: i64) -> Self {
2544 self.div(other as f64)
2545 }
2546}
2547
2548impl Div<i32> for Matrix {
2549 type Output = Self;
2550
2551 fn div(self, other: i32) -> Self {
2552 self.div(other as f64)
2553 }
2554}
2555
2556impl Div<usize> for Matrix {
2557 type Output = Self;
2558
2559 fn div(self, other: usize) -> Self {
2560 self.div(other as f64)
2561 }
2562}
2563
2564impl<'a> Div<f64> for &'a Matrix {
2565 type Output = Matrix;
2566
2567 fn div(self, other: f64) -> Self::Output {
2568 match () {
2569 #[cfg(feature = "O3")]
2570 () => {
2571 let x = &self.data;
2572 let mut y = vec![0f64; x.len()];
2573 let a_f64 = other;
2574 let n_i32 = x.len() as i32;
2575
2576 unsafe {
2577 daxpy(n_i32, 1f64 / a_f64, x, 1, &mut y, 1);
2578 }
2579 matrix(y, self.row, self.col, self.shape)
2580 }
2581 _ => self.fmap(|x| x / other),
2582 }
2583 }
2584}
2585
2586impl<'a> Div<i64> for &'a Matrix {
2587 type Output = Matrix;
2588
2589 fn div(self, other: i64) -> Self::Output {
2590 self.div(other as f64)
2591 }
2592}
2593
2594impl<'a> Div<i32> for &'a Matrix {
2595 type Output = Matrix;
2596
2597 fn div(self, other: i32) -> Self::Output {
2598 self.div(other as f64)
2599 }
2600}
2601
2602impl<'a> Div<usize> for &'a Matrix {
2603 type Output = Matrix;
2604
2605 fn div(self, other: usize) -> Self::Output {
2606 self.div(other as f64)
2607 }
2608}
2609
2610/// Index for Matrix
2611///
2612/// `(usize, usize) -> f64`
2613///
2614/// # Examples
2615/// ```
2616/// extern crate peroxide;
2617/// use peroxide::fuga::*;
2618///
2619/// let a = matrix(vec![1,2,3,4],2,2,Row);
2620/// assert_eq!(a[(0,1)], 2f64);
2621/// ```
2622impl Index<(usize, usize)> for Matrix {
2623 type Output = f64;
2624
2625 fn index(&self, pair: (usize, usize)) -> &f64 {
2626 let p = self.ptr();
2627 let i = pair.0;
2628 let j = pair.1;
2629 assert!(i < self.row && j < self.col, "Index out of range");
2630 match self.shape {
2631 Row => unsafe { &*p.add(i * self.col + j) },
2632 Col => unsafe { &*p.add(i + j * self.row) },
2633 }
2634 }
2635}
2636
2637/// IndexMut for Matrix (Assign)
2638///
2639/// `(usize, usize) -> f64`
2640///
2641/// # Examples
2642/// ```
2643/// #[macro_use]
2644/// extern crate peroxide;
2645/// use peroxide::fuga::*;
2646///
2647/// fn main() {
2648/// let mut a = matrix!(1;4;1, 2, 2, Row);
2649/// a[(1,1)] = 10.0;
2650/// assert_eq!(a, matrix(c!(1,2,3,10), 2, 2, Row));
2651/// }
2652/// ```
2653impl IndexMut<(usize, usize)> for Matrix {
2654 fn index_mut(&mut self, pair: (usize, usize)) -> &mut f64 {
2655 let i = pair.0;
2656 let j = pair.1;
2657 let r = self.row;
2658 let c = self.col;
2659 assert!(i < self.row && j < self.col, "Index out of range");
2660 let p = self.mut_ptr();
2661 match self.shape {
2662 Row => {
2663 let idx = i * c + j;
2664 unsafe { &mut *p.add(idx) }
2665 }
2666 Col => {
2667 let idx = i + j * r;
2668 unsafe { &mut *p.add(idx) }
2669 }
2670 }
2671 }
2672}
2673
2674// =============================================================================
2675// Functional Programming Tools (Hand-written)
2676// =============================================================================
2677
2678impl FPMatrix for Matrix {
2679 type Scalar = f64;
2680
2681 fn take_row(&self, n: usize) -> Self {
2682 if n >= self.row {
2683 return self.clone();
2684 }
2685 match self.shape {
2686 Row => {
2687 let new_data = self
2688 .data
2689 .clone()
2690 .into_iter()
2691 .take(n * self.col)
2692 .collect::<Vec<f64>>();
2693 matrix(new_data, n, self.col, Row)
2694 }
2695 Col => {
2696 let mut temp_data: Vec<f64> = Vec::new();
2697 for i in 0..n {
2698 temp_data.extend(self.row(i));
2699 }
2700 matrix(temp_data, n, self.col, Row)
2701 }
2702 }
2703 }
2704
2705 fn take_col(&self, n: usize) -> Self {
2706 if n >= self.col {
2707 return self.clone();
2708 }
2709 match self.shape {
2710 Col => {
2711 let new_data = self
2712 .data
2713 .clone()
2714 .into_iter()
2715 .take(n * self.row)
2716 .collect::<Vec<f64>>();
2717 matrix(new_data, self.row, n, Col)
2718 }
2719 Row => {
2720 let mut temp_data: Vec<f64> = Vec::new();
2721 for i in 0..n {
2722 temp_data.extend(self.col(i));
2723 }
2724 matrix(temp_data, self.row, n, Col)
2725 }
2726 }
2727 }
2728
2729 fn skip_row(&self, n: usize) -> Self {
2730 assert!(n < self.row, "Skip range is larger than row of matrix");
2731
2732 let mut temp_data: Vec<f64> = Vec::new();
2733 for i in n..self.row {
2734 temp_data.extend(self.row(i));
2735 }
2736 matrix(temp_data, self.row - n, self.col, Row)
2737 }
2738
2739 fn skip_col(&self, n: usize) -> Self {
2740 assert!(n < self.col, "Skip range is larger than col of matrix");
2741
2742 let mut temp_data: Vec<f64> = Vec::new();
2743 for i in n..self.col {
2744 temp_data.extend(self.col(i));
2745 }
2746 matrix(temp_data, self.row, self.col - n, Col)
2747 }
2748
2749 fn fmap<F>(&self, f: F) -> Matrix
2750 where
2751 F: Fn(f64) -> f64,
2752 {
2753 let result = self.data.iter().map(|x| f(*x)).collect::<Vec<f64>>();
2754 matrix(result, self.row, self.col, self.shape)
2755 }
2756
2757 /// Column map
2758 ///
2759 /// # Example
2760 /// ```rust
2761 /// use peroxide::fuga::*;
2762 ///
2763 /// let x = ml_matrix("1 2;3 4;5 6");
2764 /// let y = x.col_map(|c| c.fmap(|t| t - c.mean()));
2765 ///
2766 /// assert_eq!(y, ml_matrix("-2 -2;0 0;2 2"));
2767 /// ```
2768 fn col_map<F>(&self, f: F) -> Matrix
2769 where
2770 F: Fn(Vec<f64>) -> Vec<f64>,
2771 {
2772 let mut result = matrix(vec![0f64; self.row * self.col], self.row, self.col, Col);
2773
2774 for i in 0..self.col {
2775 result.subs_col(i, &f(self.col(i)));
2776 }
2777
2778 result
2779 }
2780
2781 /// Row map
2782 ///
2783 /// # Example
2784 /// ```rust
2785 /// use peroxide::fuga::*;
2786 ///
2787 /// let x = ml_matrix("1 2 3;4 5 6");
2788 /// let y = x.row_map(|r| r.fmap(|t| t - r.mean()));
2789 ///
2790 /// assert_eq!(y, ml_matrix("-1 0 1;-1 0 1"));
2791 /// ```
2792 fn row_map<F>(&self, f: F) -> Matrix
2793 where
2794 F: Fn(Vec<f64>) -> Vec<f64>,
2795 {
2796 let mut result = matrix(vec![0f64; self.row * self.col], self.row, self.col, Row);
2797
2798 for i in 0..self.row {
2799 result.subs_row(i, &f(self.row(i)));
2800 }
2801
2802 result
2803 }
2804
2805 fn col_mut_map<F>(&mut self, f: F)
2806 where
2807 F: Fn(Vec<f64>) -> Vec<f64>,
2808 {
2809 for i in 0..self.col {
2810 unsafe {
2811 let mut p = self.col_mut(i);
2812 let fv = f(self.col(i));
2813 for j in 0..p.len() {
2814 *p[j] = fv[j];
2815 }
2816 }
2817 }
2818 }
2819
2820 fn row_mut_map<F>(&mut self, f: F)
2821 where
2822 F: Fn(Vec<f64>) -> Vec<f64>,
2823 {
2824 for i in 0..self.col {
2825 unsafe {
2826 let mut p = self.row_mut(i);
2827 let fv = f(self.row(i));
2828 for j in 0..p.len() {
2829 *p[j] = fv[j];
2830 }
2831 }
2832 }
2833 }
2834
2835 fn reduce<F, T>(&self, init: T, f: F) -> f64
2836 where
2837 F: Fn(f64, f64) -> f64,
2838 T: Into<f64>,
2839 {
2840 self.data.iter().fold(init.into(), |x, y| f(x, *y))
2841 }
2842
2843 fn zip_with<F>(&self, f: F, other: &Matrix) -> Self
2844 where
2845 F: Fn(f64, f64) -> f64,
2846 {
2847 assert_eq!(self.data.len(), other.data.len());
2848 let mut a = other.clone();
2849 if self.shape != other.shape {
2850 a = a.change_shape();
2851 }
2852 let result = self
2853 .data
2854 .iter()
2855 .zip(a.data.iter())
2856 .map(|(x, y)| f(*x, *y))
2857 .collect::<Vec<f64>>();
2858 matrix(result, self.row, self.col, self.shape)
2859 }
2860
2861 fn col_reduce<F>(&self, f: F) -> Vec<f64>
2862 where
2863 F: Fn(Vec<f64>) -> f64,
2864 {
2865 let mut v = vec![0f64; self.col];
2866 for i in 0..self.col {
2867 v[i] = f(self.col(i));
2868 }
2869 v
2870 }
2871
2872 fn row_reduce<F>(&self, f: F) -> Vec<f64>
2873 where
2874 F: Fn(Vec<f64>) -> f64,
2875 {
2876 let mut v = vec![0f64; self.row];
2877 for i in 0..self.row {
2878 v[i] = f(self.row(i));
2879 }
2880 v
2881 }
2882}
2883
2884// =============================================================================
2885// Linear Algebra
2886// =============================================================================
2887pub fn diag(n: usize) -> Matrix {
2888 let mut v: Vec<f64> = vec![0f64; n * n];
2889 for i in 0..n {
2890 let idx = i * (n + 1);
2891 v[idx] = 1f64;
2892 }
2893 matrix(v, n, n, Row)
2894}
2895
2896impl PQLU<Matrix> {
2897 pub fn extract(&self) -> (Vec<usize>, Vec<usize>, Matrix, Matrix) {
2898 (
2899 self.p.clone(),
2900 self.q.clone(),
2901 self.l.clone(),
2902 self.u.clone(),
2903 )
2904 }
2905
2906 pub fn det(&self) -> f64 {
2907 // sgn of perms
2908 let mut sgn_p = 1f64;
2909 let mut sgn_q = 1f64;
2910 for (i, &j) in self.p.iter().enumerate() {
2911 if i != j {
2912 sgn_p *= -1f64;
2913 }
2914 }
2915 for (i, &j) in self.q.iter().enumerate() {
2916 if i != j {
2917 sgn_q *= -1f64;
2918 }
2919 }
2920
2921 self.u.diag().reduce(1f64, |x, y| x * y) * sgn_p * sgn_q
2922 }
2923
2924 pub fn inv(&self) -> Matrix {
2925 let (p, q, l, u) = self.extract();
2926 let mut m = inv_u(u) * inv_l(l);
2927 // Q = Q1 Q2 Q3 ..
2928 for (idx1, idx2) in q.into_iter().enumerate().rev() {
2929 unsafe {
2930 m.swap(idx1, idx2, Row);
2931 }
2932 }
2933 // P = Pn-1 .. P3 P2 P1
2934 for (idx1, idx2) in p.into_iter().enumerate().rev() {
2935 unsafe {
2936 m.swap(idx1, idx2, Col);
2937 }
2938 }
2939 m
2940 }
2941}
2942
2943impl SVD<Matrix> {
2944 pub fn u(&self) -> &Matrix {
2945 &self.u
2946 }
2947
2948 pub fn vt(&self) -> &Matrix {
2949 &self.vt
2950 }
2951
2952 pub fn s_mat(&self) -> Matrix {
2953 let mut mat = zeros(self.u.col, self.vt.row);
2954 for i in 0..mat.row.min(mat.col) {
2955 mat[(i, i)] = self.s[i];
2956 }
2957 mat
2958 }
2959
2960 /// Generated Truncated SVD
2961 ///
2962 /// ```rust
2963 /// extern crate peroxide;
2964 /// use peroxide::fuga::*;
2965 ///
2966 /// fn main() {
2967 /// let x = ml_matrix("1 2 3;4 5 6");
2968 /// # #[cfg(feature = "O3")] {
2969 /// // Full SVD
2970 /// let svd = x.svd();
2971 /// svd.u().print(); // m x m matrix
2972 /// svd.s_mat().print(); // m x n matrix
2973 /// svd.vt().print(); // n x n matrix
2974 ///
2975 /// // Truncated SVD
2976 /// let svd2 = svd.truncated();
2977 /// svd2.u().print(); // m x p matrix
2978 /// svd2.s_mat().print(); // p x p matrix
2979 /// svd2.vt().print(); // p x n matrix
2980 /// # }
2981 /// }
2982 /// ```
2983 pub fn truncated(&self) -> Self {
2984 let mut s: Vec<f64> = vec![];
2985 let mut u = matrix::<f64>(vec![], self.u.row, 0, Col);
2986 let mut vt = matrix::<f64>(vec![], 0, self.vt.col, Row);
2987 for (i, sig) in self.s.iter().enumerate() {
2988 if *sig == 0f64 {
2989 continue;
2990 } else {
2991 s.push(*sig);
2992 u.add_col_mut(&self.u.col(i));
2993 vt.add_row_mut(&self.vt.row(i));
2994 }
2995 }
2996
2997 SVD { s, u, vt }
2998 }
2999}
3000
3001impl LinearAlgebra<Matrix> for Matrix {
3002 /// Backward Substitution for Upper Triangular
3003 fn back_subs(&self, b: &[f64]) -> Vec<f64> {
3004 let n = self.col;
3005 let mut y = vec![0f64; n];
3006 y[n - 1] = b[n - 1] / self[(n - 1, n - 1)];
3007 for i in (0..n - 1).rev() {
3008 let mut s = 0f64;
3009 for j in i + 1..n {
3010 s += self[(i, j)] * y[j];
3011 }
3012 y[i] = 1f64 / self[(i, i)] * (b[i] - s);
3013 }
3014 y
3015 }
3016
3017 /// Forward substitution for Lower Triangular
3018 fn forward_subs(&self, b: &[f64]) -> Vec<f64> {
3019 let n = self.col;
3020 let mut y = vec![0f64; n];
3021 y[0] = b[0] / self[(0, 0)];
3022 for i in 1..n {
3023 let mut s = 0f64;
3024 for j in 0..i {
3025 s += self[(i, j)] * y[j];
3026 }
3027 y[i] = 1f64 / self[(i, i)] * (b[i] - s);
3028 }
3029 y
3030 }
3031
3032 /// LU Decomposition Implements (Complete Pivot)
3033 ///
3034 /// # Description
3035 /// It use complete pivoting LU decomposition.
3036 /// You can get two permutations, and LU matrices.
3037 ///
3038 /// # Caution
3039 /// It returns `Option<PQLU>` - You should unwrap to obtain real value.
3040 /// - `PQLU` has four field - `p`, `q`, `l`, `u`.
3041 /// - `p`, `q` are permutations.
3042 /// - `l`, `u` are matrices.
3043 ///
3044 /// # Examples
3045 /// ```
3046 /// use peroxide::fuga::*;
3047 ///
3048 /// let a = matrix(vec![1,2,3,4], 2, 2, Row);
3049 /// let pqlu = a.lu();
3050 /// let (p,q,l,u) = (pqlu.p, pqlu.q, pqlu.l, pqlu.u);
3051 /// assert_eq!(p, vec![1]); // swap 0 & 1 (Row)
3052 /// assert_eq!(q, vec![1]); // swap 0 & 1 (Col)
3053 /// assert_eq!(l, matrix(vec![1.0,0.0,0.5,1.0],2,2,Row));
3054 /// assert_eq!(u, matrix(vec![4.0,3.0,0.0,-0.5],2,2,Row));
3055 /// ```
3056 fn lu(&self) -> PQLU<Matrix> {
3057 assert_eq!(self.col, self.row);
3058 let n = self.row;
3059 let len: usize = n * n;
3060
3061 let mut l = eye(n);
3062 let mut u = matrix(vec![0f64; len], n, n, self.shape);
3063
3064 let mut temp = self.clone();
3065 let (p, q) = gecp(&mut temp);
3066 for i in 0..n {
3067 for j in 0..i {
3068 // Inverse multiplier
3069 l[(i, j)] = -temp[(i, j)];
3070 }
3071 for j in i..n {
3072 u[(i, j)] = temp[(i, j)];
3073 }
3074 }
3075 // Pivoting L
3076 for i in 0..n - 1 {
3077 unsafe {
3078 let l_i = l.col_mut(i);
3079 for j in i + 1..l.col - 1 {
3080 let dst = p[j];
3081 std::ptr::swap(l_i[j], l_i[dst]);
3082 }
3083 }
3084 }
3085 PQLU { p, q, l, u }
3086 }
3087
3088 fn waz(&self, d_form: Form) -> Option<WAZD<Matrix>> {
3089 match d_form {
3090 Form::Diagonal => {
3091 let n = self.row;
3092 let mut w = eye(n);
3093 let mut z = eye(n);
3094 let mut d = eye(n);
3095 let mut q = vec![0f64; n];
3096 let mut p = vec![0f64; n];
3097
3098 for i in 0..n {
3099 let m_i = self.col(i);
3100 let pq = w.col(i).dot(&m_i);
3101 d[(i, i)] = pq;
3102 if pq == 0f64 {
3103 return None;
3104 }
3105 for j in i + 1..n {
3106 q[j] = w.col(j).dot(&m_i) / pq;
3107 p[j] = z.col(j).dot(&self.row(i)) / pq;
3108 }
3109 for j in i + 1..n {
3110 for k in 0..i + 1 {
3111 w[(k, j)] -= q[j] * w[(k, i)];
3112 z[(k, j)] -= p[j] * z[(k, i)];
3113 }
3114 }
3115 }
3116 Some(WAZD { w, z, d })
3117 }
3118 Form::Identity => {
3119 let n = self.row;
3120 let mut w = eye(n);
3121 let mut z = eye(n);
3122 let mut p = zeros(n, n);
3123 let mut q = zeros(n, n);
3124
3125 for i in 0..n {
3126 let m_i = self.col(i);
3127 let p_ii = w.col(i).dot(&m_i);
3128 p[(i, i)] = p_ii;
3129 if p_ii == 0f64 {
3130 return None;
3131 }
3132 for j in i + 1..n {
3133 q[(i, j)] = w.col(j).dot(&m_i) / p_ii;
3134 p[(i, j)] = z.col(j).dot(&self.row(i)) / p_ii;
3135 for k in 0..j {
3136 w[(k, j)] -= q[(i, j)] * w[(k, i)];
3137 z[(k, j)] -= p[(i, j)] * z[(k, i)];
3138 }
3139 }
3140 unsafe {
3141 let col_ptr = z.col_mut(i);
3142 col_ptr.into_iter().for_each(|x| *x /= p_ii);
3143 }
3144 }
3145 Some(WAZD { w, z, d: eye(n) })
3146 }
3147 }
3148 }
3149
3150 /// QR Decomposition
3151 ///
3152 /// Translation of [RosettaCode#Python](https://rosettacode.org/wiki/QR_decomposition#Python)
3153 ///
3154 /// # Example
3155 /// ```
3156 /// use peroxide::fuga::*;
3157 ///
3158 /// let a = ml_matrix("12 -51 4;6 167 -68; -4 24 -41");
3159 /// let qr = a.qr();
3160 /// let r = ml_matrix("-14 -21 14; 0 -175 70; 0 0 -35");
3161 /// #[cfg(feature="O3")]
3162 /// {
3163 /// assert_eq!(r, qr.r);
3164 /// }
3165 /// qr.r.print();
3166 /// ```
3167 #[allow(non_snake_case)]
3168 fn qr(&self) -> QR<Matrix> {
3169 match () {
3170 #[cfg(feature = "O3")]
3171 () => {
3172 let opt_dgeqrf = lapack_dgeqrf(self);
3173 match opt_dgeqrf {
3174 None => panic!("Serious problem in QR decomposition"),
3175 Some(dgeqrf) => {
3176 let q = dgeqrf.get_Q();
3177 let r = dgeqrf.get_R();
3178 QR { q, r }
3179 }
3180 }
3181 }
3182 _ => {
3183 let m = self.row;
3184 let n = self.col;
3185
3186 let mut r = self.clone();
3187 let mut q = eye(m);
3188 let sub = if m == n { 1 } else { 0 };
3189 for i in 0..n - sub {
3190 let mut H = eye(m);
3191 let hh = gen_householder(&self.col(i).skip(i));
3192 for j in i..m {
3193 for k in i..m {
3194 H[(j, k)] = hh[(j - i, k - i)];
3195 }
3196 }
3197 q = &q * &H;
3198 r = &H * &r;
3199 }
3200
3201 QR { q, r }
3202 }
3203 }
3204 }
3205
3206 /// Singular Value Decomposition
3207 ///
3208 /// # Examples
3209 /// ```
3210 /// use peroxide::fuga::*;
3211 ///
3212 /// let a = ml_matrix("3 2 2;2 3 -2");
3213 /// #[cfg(feature="O3")]
3214 /// {
3215 /// let svd = a.svd();
3216 /// assert!(eq_vec(&vec![5f64, 3f64], &svd.s, 1e-7));
3217 /// }
3218 /// a.print();
3219 /// ```
3220 fn svd(&self) -> SVD<Matrix> {
3221 match () {
3222 #[cfg(feature = "O3")]
3223 () => {
3224 let opt_dgesvd = lapack_dgesvd(self);
3225 match opt_dgesvd {
3226 None => panic!("Illegal value in LAPACK SVD"),
3227 Some(dgesvd) => match dgesvd.status {
3228 SVD_STATUS::Diverge(i) => {
3229 panic!("Divergent occurs in SVD - {} iterations", i)
3230 }
3231 SVD_STATUS::Success => SVD {
3232 s: dgesvd.s,
3233 u: dgesvd.u,
3234 vt: dgesvd.vt,
3235 },
3236 },
3237 }
3238 }
3239 _ => {
3240 unimplemented!()
3241 }
3242 }
3243 }
3244
3245 /// Cholesky Decomposition
3246 ///
3247 /// # Examples
3248 /// ```
3249 /// use peroxide::fuga::*;
3250 ///
3251 /// let a = ml_matrix("1 2;2 5");
3252 /// #[cfg(feature = "O3")]
3253 /// {
3254 /// let u = a.cholesky(Upper);
3255 /// let l = a.cholesky(Lower);
3256 ///
3257 /// assert_eq!(u, ml_matrix("1 2;0 1"));
3258 /// assert_eq!(l, ml_matrix("1 0;2 1"));
3259 /// }
3260 /// a.print();
3261 /// ```
3262 #[cfg(feature = "O3")]
3263 fn cholesky(&self, uplo: UPLO) -> Matrix {
3264 match () {
3265 #[cfg(feature = "O3")]
3266 () => {
3267 if !self.is_symmetric() {
3268 panic!("Cholesky Error: Matrix is not symmetric!");
3269 }
3270 let dpotrf = lapack_dpotrf(self, uplo);
3271 match dpotrf {
3272 None => panic!("Cholesky Error: Not symmetric or not positive definite."),
3273 Some(x) => {
3274 match x.status {
3275 POSITIVE_STATUS::Failed(i) => panic!("Cholesky Error: the leading minor of order {} is not positive definite", i),
3276 POSITIVE_STATUS::Success => {
3277 match uplo {
3278 UPLO::Upper => x.get_U().unwrap(),
3279 UPLO::Lower => x.get_L().unwrap()
3280 }
3281 }
3282 }
3283 }
3284 }
3285 }
3286 _ => {
3287 unimplemented!()
3288 }
3289 }
3290 }
3291
3292 /// Reduced Row Echelon Form
3293 ///
3294 /// Implementation of [RosettaCode](https://rosettacode.org/wiki/Reduced_row_echelon_form)
3295 fn rref(&self) -> Matrix {
3296 let mut lead = 0usize;
3297 let mut result = self.clone();
3298 'outer: for r in 0..self.row {
3299 if self.col <= lead {
3300 break;
3301 }
3302 let mut i = r;
3303 while result[(i, lead)] == 0f64 {
3304 i += 1;
3305 if self.row == i {
3306 i = r;
3307 lead += 1;
3308 if self.col == lead {
3309 break 'outer;
3310 }
3311 }
3312 }
3313 unsafe {
3314 result.swap(i, r, Row);
3315 }
3316 let tmp = result[(r, lead)];
3317 if tmp != 0f64 {
3318 unsafe {
3319 result.row_mut(r).iter_mut().for_each(|t| *(*t) /= tmp);
3320 }
3321 }
3322 for j in 0..result.row {
3323 if j != r {
3324 let tmp1 = result.row(r).mul_scalar(result[(j, lead)]);
3325 let tmp2 = result.row(j).sub_vec(&tmp1);
3326 result.subs_row(j, &tmp2);
3327 }
3328 }
3329 lead += 1;
3330 }
3331 result
3332 }
3333
3334 /// Determinant
3335 ///
3336 /// # Examples
3337 /// ```
3338 /// #[macro_use]
3339 /// extern crate peroxide;
3340 /// use peroxide::fuga::*;
3341 ///
3342 /// fn main() {
3343 /// let a = matrix!(1;4;1, 2, 2, Row);
3344 /// assert_eq!(a.det(), -2f64);
3345 /// }
3346 /// ```
3347 fn det(&self) -> f64 {
3348 assert_eq!(self.row, self.col);
3349 match () {
3350 #[cfg(feature = "O3")]
3351 () => {
3352 let opt_dgrf = lapack_dgetrf(self);
3353 match opt_dgrf {
3354 None => f64::NAN,
3355 Some(dgrf) => match dgrf.status {
3356 LAPACK_STATUS::Singular => 0f64,
3357 LAPACK_STATUS::NonSingular => {
3358 let mat = &dgrf.fact_mat;
3359 let ipiv = &dgrf.ipiv;
3360 let n = mat.col;
3361 let mut sgn = 1i32;
3362 let mut d = 1f64;
3363 for i in 0..n {
3364 d *= mat[(i, i)];
3365 }
3366 for i in 0..ipiv.len() {
3367 if ipiv[i] - 1 != i as i32 {
3368 sgn *= -1;
3369 }
3370 }
3371 (sgn as f64) * d
3372 }
3373 },
3374 }
3375 }
3376 _ => self.lu().det(),
3377 }
3378 }
3379
3380 /// Block Partition
3381 ///
3382 /// # Examples
3383 /// ```
3384 /// #[macro_use]
3385 /// extern crate peroxide;
3386 /// use peroxide::fuga::*;
3387 ///
3388 /// fn main() {
3389 /// let a = matrix!(1;16;1, 4, 4, Row);
3390 /// let (m1, m2, m3, m4) = a.block();
3391 /// assert_eq!(m1, matrix(c!(1,2,5,6), 2, 2, Row));
3392 /// assert_eq!(m2, matrix(c!(3,4,7,8), 2, 2, Row));
3393 /// assert_eq!(m3, matrix(c!(9,10,13,14), 2, 2, Row));
3394 /// assert_eq!(m4, matrix(c!(11,12,15,16), 2, 2, Row));
3395 ///
3396 /// let b = matrix!(1;16;1, 4, 4, Col);
3397 /// let (m1, m2, m3, m4) = b.block();
3398 /// assert_eq!(m1, matrix(c!(1,2,5,6), 2, 2, Col));
3399 /// assert_eq!(m3, matrix(c!(3,4,7,8), 2, 2, Col));
3400 /// assert_eq!(m2, matrix(c!(9,10,13,14), 2, 2, Col));
3401 /// assert_eq!(m4, matrix(c!(11,12,15,16), 2, 2, Col));
3402 /// }
3403 /// ```
3404 fn block(&self) -> (Self, Self, Self, Self) {
3405 let r = self.row;
3406 let c = self.col;
3407 let l_r = self.row / 2;
3408 let l_c = self.col / 2;
3409 let r_l = r - l_r;
3410 let c_l = c - l_c;
3411
3412 let mut m1 = matrix(vec![0f64; l_r * l_c], l_r, l_c, self.shape);
3413 let mut m2 = matrix(vec![0f64; l_r * c_l], l_r, c_l, self.shape);
3414 let mut m3 = matrix(vec![0f64; r_l * l_c], r_l, l_c, self.shape);
3415 let mut m4 = matrix(vec![0f64; r_l * c_l], r_l, c_l, self.shape);
3416
3417 for idx_row in 0..r {
3418 for idx_col in 0..c {
3419 match (idx_row, idx_col) {
3420 (i, j) if (i < l_r) && (j < l_c) => {
3421 m1[(i, j)] = self[(i, j)];
3422 }
3423 (i, j) if (i < l_r) && (j >= l_c) => {
3424 m2[(i, j - l_c)] = self[(i, j)];
3425 }
3426 (i, j) if (i >= l_r) && (j < l_c) => {
3427 m3[(i - l_r, j)] = self[(i, j)];
3428 }
3429 (i, j) if (i >= l_r) && (j >= l_c) => {
3430 m4[(i - l_r, j - l_c)] = self[(i, j)];
3431 }
3432 _ => (),
3433 }
3434 }
3435 }
3436 (m1, m2, m3, m4)
3437 }
3438
3439 /// Inverse of Matrix
3440 ///
3441 /// # Caution
3442 ///
3443 /// `inv` function returns `Option<Matrix>`
3444 /// Thus, you should use pattern matching or `unwrap` to obtain inverse.
3445 ///
3446 /// # Examples
3447 /// ```
3448 /// #[macro_use]
3449 /// extern crate peroxide;
3450 /// use peroxide::fuga::*;
3451 ///
3452 /// fn main() {
3453 /// // Non-singular
3454 /// let a = matrix!(1;4;1, 2, 2, Row);
3455 /// assert_eq!(a.inv(), matrix(c!(-2,1,1.5,-0.5),2,2,Row));
3456 ///
3457 /// // Singular
3458 /// //let b = matrix!(1;9;1, 3, 3, Row);
3459 /// //let c = b.inv(); // Can compile but..
3460 /// //let d = b.det();
3461 /// //assert_eq!(d, 0f64);
3462 /// }
3463 /// ```
3464 fn inv(&self) -> Self {
3465 match () {
3466 #[cfg(feature = "O3")]
3467 () => {
3468 let opt_dgrf = lapack_dgetrf(self);
3469 match opt_dgrf {
3470 None => panic!("Singular matrix (opt_dgrf)"),
3471 Some(dgrf) => match dgrf.status {
3472 LAPACK_STATUS::NonSingular => lapack_dgetri(&dgrf).unwrap(),
3473 LAPACK_STATUS::Singular => {
3474 panic!("Singular matrix (LAPACK_STATUS Singular)")
3475 }
3476 },
3477 }
3478 }
3479 _ => self.lu().inv(),
3480 // _ => {
3481 // match self.lu() {
3482 // None => None,
3483 // Some(lu) => Some(lu.inv())
3484 // }
3485 // }
3486 }
3487 }
3488
3489 /// Moore-Penrose Pseudo inverse
3490 ///
3491 /// # Description
3492 /// `$X^\dagger = (X^T X)^{-1} X^T$`
3493 ///
3494 /// # Examples
3495 /// ```
3496 /// #[macro_use]
3497 /// extern crate peroxide;
3498 /// use peroxide::fuga::*;
3499 ///
3500 /// fn main() {
3501 /// let a = matrix!(1;4;1, 2, 2, Row);
3502 /// let inv_a = a.inv();
3503 /// let pse_a = a.pseudo_inv();
3504 ///
3505 /// assert_eq!(inv_a, pse_a); // Nearly equal
3506 /// }
3507 /// ```
3508 fn pseudo_inv(&self) -> Self {
3509 match () {
3510 #[cfg(feature = "O3")]
3511 () => {
3512 let svd = self.svd();
3513 let row = svd.vt.row;
3514 let col = svd.u.row;
3515 let mut sp = zeros(row, col);
3516 for i in 0..row.min(col) {
3517 sp[(i, i)] = 1f64 / svd.s[i];
3518 }
3519 svd.vt.t() * sp * svd.u.t()
3520 }
3521 _ => {
3522 let xt = self.t();
3523 let xtx = &xt * self;
3524 xtx.inv() * xt
3525 }
3526 }
3527 }
3528
3529 /// Solve with Vector
3530 ///
3531 /// # Solve options
3532 ///
3533 /// * LU: Gaussian elimination with Complete pivoting LU (GECP)
3534 /// * WAZ: Solve with WAZ decomposition
3535 ///
3536 /// # Reference
3537 ///
3538 /// * Biswa Nath Datta, *Numerical Linear Algebra and Applications, Second Edition*
3539 /// * Ke Chen, *Matrix Preconditioning Techniques and Applications*, Cambridge Monographs on Applied and Computational Mathematics
3540 fn solve(&self, b: &[f64], sk: SolveKind) -> Vec<f64> {
3541 match sk {
3542 #[cfg(feature = "O3")]
3543 SolveKind::LU => {
3544 let opt_dgrf = lapack_dgetrf(self);
3545 match opt_dgrf {
3546 None => panic!("Try solve for Singluar matrix"),
3547 Some(dgrf) => match dgrf.status {
3548 LAPACK_STATUS::Singular => panic!("Try solve for Singluar matrix"),
3549 LAPACK_STATUS::NonSingular => {
3550 let b = b.to_vec();
3551 lapack_dgetrs(&dgrf, &(b.into())).unwrap().into()
3552 }
3553 },
3554 }
3555 }
3556 #[cfg(not(feature = "O3"))]
3557 SolveKind::LU => {
3558 let lu = self.lu();
3559 let (p, q, l, u) = lu.extract();
3560 let mut v = b.to_vec();
3561 v.swap_with_perm(&p.into_iter().enumerate().collect());
3562 let z = l.forward_subs(&v);
3563 let mut y = u.back_subs(&z);
3564 y.swap_with_perm(&q.into_iter().enumerate().rev().collect());
3565 y
3566 }
3567 SolveKind::WAZ => {
3568 let wazd = match self.waz(Form::Identity) {
3569 None => panic!("Can't solve by WAZ with Singular matrix!"),
3570 Some(obj) => obj,
3571 };
3572 let x = &wazd.w.t() * &b.to_vec();
3573 let x = &wazd.z * &x;
3574 x
3575 }
3576 }
3577 }
3578
3579 fn solve_mat(&self, m: &Matrix, sk: SolveKind) -> Matrix {
3580 match sk {
3581 #[cfg(feature = "O3")]
3582 SolveKind::LU => {
3583 let opt_dgrf = lapack_dgetrf(self);
3584 match opt_dgrf {
3585 None => panic!("Try solve for Singluar matrix"),
3586 Some(dgrf) => match dgrf.status {
3587 LAPACK_STATUS::Singular => panic!("Try solve for Singluar matrix"),
3588 LAPACK_STATUS::NonSingular => lapack_dgetrs(&dgrf, m).unwrap(),
3589 },
3590 }
3591 }
3592 #[cfg(not(feature = "O3"))]
3593 SolveKind::LU => {
3594 let lu = self.lu();
3595 let (p, q, l, u) = lu.extract();
3596 let mut x = matrix(vec![0f64; self.col * m.col], self.col, m.col, Col);
3597 for i in 0..m.col {
3598 let mut v = m.col(i).clone();
3599 for (r, &s) in p.iter().enumerate() {
3600 v.swap(r, s);
3601 }
3602 let z = l.forward_subs(&v);
3603 let mut y = u.back_subs(&z);
3604 for (r, &s) in q.iter().enumerate() {
3605 y.swap(r, s);
3606 }
3607 unsafe {
3608 let mut c = x.col_mut(i);
3609 copy_vec_ptr(&mut c, &y);
3610 }
3611 }
3612 x
3613 }
3614 SolveKind::WAZ => {
3615 let wazd = match self.waz(Form::Identity) {
3616 None => panic!("Try solve for Singular matrix"),
3617 Some(obj) => obj,
3618 };
3619 let x = &wazd.w.t() * m;
3620 let x = &wazd.z * &x;
3621 x
3622 }
3623 }
3624 }
3625
3626 fn is_symmetric(&self) -> bool {
3627 if self.row != self.col {
3628 return false;
3629 }
3630
3631 for i in 0..self.row {
3632 for j in i..self.col {
3633 if !nearly_eq(self[(i, j)], self[(j, i)]) {
3634 return false;
3635 }
3636 }
3637 }
3638 true
3639 }
3640}
3641
3642impl MutMatrix for Matrix {
3643 type Scalar = f64;
3644
3645 unsafe fn col_mut(&mut self, idx: usize) -> Vec<*mut f64> {
3646 assert!(idx < self.col, "Index out of range");
3647 match self.shape {
3648 Col => {
3649 let mut v: Vec<*mut f64> = vec![&mut 0f64; self.row];
3650 let start_idx = idx * self.row;
3651 let p = self.mut_ptr();
3652 for (i, j) in (start_idx..start_idx + v.len()).enumerate() {
3653 v[i] = p.add(j);
3654 }
3655 v
3656 }
3657 Row => {
3658 let mut v: Vec<*mut f64> = vec![&mut 0f64; self.row];
3659 let p = self.mut_ptr();
3660 for i in 0..v.len() {
3661 v[i] = p.add(idx + i * self.col);
3662 }
3663 v
3664 }
3665 }
3666 }
3667
3668 unsafe fn row_mut(&mut self, idx: usize) -> Vec<*mut f64> {
3669 assert!(idx < self.row, "Index out of range");
3670 match self.shape {
3671 Row => {
3672 let mut v: Vec<*mut f64> = vec![&mut 0f64; self.col];
3673 let start_idx = idx * self.col;
3674 let p = self.mut_ptr();
3675 for (i, j) in (start_idx..start_idx + v.len()).enumerate() {
3676 v[i] = p.add(j);
3677 }
3678 v
3679 }
3680 Col => {
3681 let mut v: Vec<*mut f64> = vec![&mut 0f64; self.col];
3682 let p = self.mut_ptr();
3683 for i in 0..v.len() {
3684 v[i] = p.add(idx + i * self.row);
3685 }
3686 v
3687 }
3688 }
3689 }
3690
3691 unsafe fn swap(&mut self, idx1: usize, idx2: usize, shape: Shape) {
3692 match shape {
3693 Col => swap_vec_ptr(&mut self.col_mut(idx1), &mut self.col_mut(idx2)),
3694 Row => swap_vec_ptr(&mut self.row_mut(idx1), &mut self.row_mut(idx2)),
3695 }
3696 }
3697
3698 unsafe fn swap_with_perm(&mut self, p: &Vec<(usize, usize)>, shape: Shape) {
3699 for (i, j) in p.iter() {
3700 self.swap(*i, *j, shape)
3701 }
3702 }
3703}
3704
3705// =============================================================================
3706// Matrix as Numeric<f64>
3707// =============================================================================
3708impl Div<Matrix> for f64 {
3709 type Output = Matrix;
3710 fn div(self, m: Matrix) -> Matrix {
3711 m.fmap(|x| self / x)
3712 }
3713}
3714
3715impl Div<Matrix> for f32 {
3716 type Output = Matrix;
3717 fn div(self, m: Matrix) -> Matrix {
3718 m.fmap(|x| self as f64 / x)
3719 }
3720}
3721
3722impl Div<Matrix> for Matrix {
3723 type Output = Matrix;
3724
3725 fn div(self, m: Matrix) -> Matrix {
3726 self.mul(m.inv())
3727 }
3728}
3729
3730impl ExpLogOps for Matrix {
3731 type Float = f64;
3732 fn exp(&self) -> Self {
3733 self.fmap(|x| x.exp())
3734 }
3735 fn ln(&self) -> Self {
3736 self.fmap(|x| x.ln())
3737 }
3738 fn log(&self, base: Self::Float) -> Self {
3739 self.fmap(|x| x.log(base))
3740 }
3741 fn log2(&self) -> Self {
3742 self.fmap(|x| x.log2())
3743 }
3744 fn log10(&self) -> Self {
3745 self.fmap(|x| x.log10())
3746 }
3747}
3748
3749impl PowOps for Matrix {
3750 type Float = f64;
3751
3752 fn powi(&self, n: i32) -> Self {
3753 self.fmap(|x| x.powi(n))
3754 }
3755
3756 fn powf(&self, f: Self::Float) -> Self {
3757 self.fmap(|x| x.powf(f))
3758 }
3759
3760 fn pow(&self, _f: Self) -> Self {
3761 unimplemented!()
3762 }
3763
3764 fn sqrt(&self) -> Self {
3765 self.fmap(|x| x.sqrt())
3766 }
3767}
3768
3769impl TrigOps for Matrix {
3770 fn sin_cos(&self) -> (Self, Self) {
3771 let (sin, cos) = self.data.iter().map(|x| x.sin_cos()).unzip();
3772 (
3773 matrix(sin, self.row, self.col, self.shape),
3774 matrix(cos, self.row, self.col, self.shape),
3775 )
3776 }
3777
3778 fn sin(&self) -> Self {
3779 self.fmap(|x| x.sin())
3780 }
3781
3782 fn cos(&self) -> Self {
3783 self.fmap(|x| x.cos())
3784 }
3785
3786 fn tan(&self) -> Self {
3787 self.fmap(|x| x.tan())
3788 }
3789
3790 fn sinh(&self) -> Self {
3791 self.fmap(|x| x.sinh())
3792 }
3793
3794 fn cosh(&self) -> Self {
3795 self.fmap(|x| x.cosh())
3796 }
3797
3798 fn tanh(&self) -> Self {
3799 self.fmap(|x| x.tanh())
3800 }
3801
3802 fn asin(&self) -> Self {
3803 self.fmap(|x| x.asin())
3804 }
3805
3806 fn acos(&self) -> Self {
3807 self.fmap(|x| x.acos())
3808 }
3809
3810 fn atan(&self) -> Self {
3811 self.fmap(|x| x.atan())
3812 }
3813
3814 fn asinh(&self) -> Self {
3815 self.fmap(|x| x.asinh())
3816 }
3817
3818 fn acosh(&self) -> Self {
3819 self.fmap(|x| x.acosh())
3820 }
3821
3822 fn atanh(&self) -> Self {
3823 self.fmap(|x| x.atanh())
3824 }
3825}
3826
3827impl Numeric<f64> for Matrix {}
3828
3829// =============================================================================
3830// Back-end Utils
3831// =============================================================================
3832
3833/// Combine separated matrix to one matrix
3834///
3835/// # Examples
3836/// ```
3837/// #[macro_use]
3838/// extern crate peroxide;
3839/// use peroxide::fuga::*;
3840///
3841/// fn main() {
3842/// let a = matrix!(1;16;1, 4, 4, Row);
3843/// let (m1, m2, m3, m4) = a.block();
3844/// let m = combine(m1,m2,m3,m4);
3845/// assert_eq!(m, a);
3846///
3847/// let b = matrix!(1;16;1, 4, 4, Col);
3848/// let (n1, n2, n3, n4) = b.block();
3849/// let n = combine(n1,n2,n3,n4);
3850/// assert_eq!(n, b);
3851/// }
3852/// ```
3853pub fn combine(m1: Matrix, m2: Matrix, m3: Matrix, m4: Matrix) -> Matrix {
3854 let l_r = m1.row;
3855 let l_c = m1.col;
3856 let c_l = m2.col;
3857 let r_l = m3.row;
3858
3859 let r = l_r + r_l;
3860 let c = l_c + c_l;
3861
3862 let mut m = matrix(vec![0f64; r * c], r, c, m1.shape);
3863
3864 for idx_row in 0..r {
3865 for idx_col in 0..c {
3866 match (idx_row, idx_col) {
3867 (i, j) if (i < l_r) && (j < l_c) => {
3868 m[(i, j)] = m1[(i, j)];
3869 }
3870 (i, j) if (i < l_r) && (j >= l_c) => {
3871 m[(i, j)] = m2[(i, j - l_c)];
3872 }
3873 (i, j) if (i >= l_r) && (j < l_c) => {
3874 m[(i, j)] = m3[(i - l_r, j)];
3875 }
3876 (i, j) if (i >= l_r) && (j >= l_c) => {
3877 m[(i, j)] = m4[(i - l_r, j - l_c)];
3878 }
3879 _ => (),
3880 }
3881 }
3882 }
3883 m
3884}
3885
3886/// Inverse of Lower matrix
3887///
3888/// # Examples
3889/// ```
3890/// #[macro_use]
3891/// extern crate peroxide;
3892/// use peroxide::fuga::*;
3893///
3894/// fn main() {
3895/// let a = matrix(c!(1,0,2,1), 2, 2, Row);
3896/// assert_eq!(inv_l(a), matrix(c!(1,0,-2,1), 2, 2, Row));
3897///
3898/// let b = matrix(c!(1,0,0,2,1,0,4,3,1), 3, 3, Row);
3899/// assert_eq!(inv_l(b), matrix(c!(1,0,0,-2,1,0,2,-3,1), 3, 3, Row));
3900/// }
3901/// ```
3902pub fn inv_l(l: Matrix) -> Matrix {
3903 let mut m = l.clone();
3904
3905 match l.row {
3906 1 => l,
3907 2 => {
3908 m[(1, 0)] = -m[(1, 0)];
3909 m
3910 }
3911 _ => {
3912 let (l1, l2, l3, l4) = l.block();
3913
3914 let m1 = inv_l(l1);
3915 let m2 = l2;
3916 let m4 = inv_l(l4);
3917 let m3 = -(&(&m4 * &l3) * &m1);
3918
3919 combine(m1, m2, m3, m4)
3920 }
3921 }
3922}
3923
3924/// Inverse of upper triangular matrix
3925///
3926/// # Examples
3927/// ```
3928/// #[macro_use]
3929/// extern crate peroxide;
3930/// use peroxide::fuga::*;
3931///
3932/// fn main() {
3933/// let u = matrix(c!(2,2,0,1), 2, 2, Row);
3934/// assert_eq!(inv_u(u), matrix(c!(0.5,-1,0,1), 2, 2, Row));
3935/// }
3936/// ```
3937pub fn inv_u(u: Matrix) -> Matrix {
3938 let mut w = u.clone();
3939
3940 match u.row {
3941 1 => {
3942 w[(0, 0)] = 1f64 / w[(0, 0)];
3943 w
3944 }
3945 2 => {
3946 let a = w[(0, 0)];
3947 let b = w[(0, 1)];
3948 let c = w[(1, 1)];
3949 let d = a * c;
3950
3951 w[(0, 0)] = 1f64 / a;
3952 w[(0, 1)] = -b / d;
3953 w[(1, 1)] = 1f64 / c;
3954 w
3955 }
3956 _ => {
3957 let (u1, u2, u3, u4) = u.block();
3958 let m1 = inv_u(u1);
3959 let m3 = u3;
3960 let m4 = inv_u(u4);
3961 let m2 = -(m1.clone() * u2 * m4.clone());
3962
3963 combine(m1, m2, m3, m4)
3964 }
3965 }
3966}
3967
3968/// Matrix multiply back-ends
3969fn matmul(a: &Matrix, b: &Matrix) -> Matrix {
3970 assert_eq!(a.col, b.row);
3971 let mut c = matrix(vec![0f64; a.row * b.col], a.row, b.col, a.shape);
3972 gemm(1f64, a, b, 0f64, &mut c);
3973 c
3974}
3975
3976/// GEMM wrapper for Matrixmultiply
3977///
3978/// # Examples
3979/// ```
3980/// #[macro_use]
3981/// extern crate peroxide;
3982/// use peroxide::prelude::*;
3983///
3984/// fn main() {
3985/// let a = ml_matrix("1 2 3;4 5 6");
3986/// let b = ml_matrix("1 2;3 4;5 6");
3987/// let mut c1 = zeros(2, 2);
3988/// let mut c2 = matrix(vec![1f64; 9], 3, 3, Col);
3989///
3990/// gemm(1f64, &a, &b, 0f64, &mut c1);
3991/// gemm(1f64, &b, &a, 2f64, &mut c2);
3992///
3993/// assert_eq!(c1, ml_matrix("22 28; 49 64"));
3994/// assert_eq!(c2, ml_matrix("11 14 17;21 28 35;31 42 53"));
3995/// }
3996/// ```
3997pub fn gemm(alpha: f64, a: &Matrix, b: &Matrix, beta: f64, c: &mut Matrix) {
3998 let m = a.row;
3999 let k = a.col;
4000 let n = b.col;
4001 let (rsa, csa) = match a.shape {
4002 Row => (a.col as isize, 1isize),
4003 Col => (1isize, a.row as isize),
4004 };
4005 let (rsb, csb) = match b.shape {
4006 Row => (b.col as isize, 1isize),
4007 Col => (1isize, b.row as isize),
4008 };
4009 let (rsc, csc) = match c.shape {
4010 Row => (c.col as isize, 1isize),
4011 Col => (1isize, c.row as isize),
4012 };
4013
4014 unsafe {
4015 matrixmultiply::dgemm(
4016 m,
4017 k,
4018 n,
4019 alpha,
4020 a.ptr(),
4021 rsa,
4022 csa,
4023 b.ptr(),
4024 rsb,
4025 csb,
4026 beta,
4027 c.mut_ptr(),
4028 rsc,
4029 csc,
4030 )
4031 }
4032}
4033
4034/// General Matrix-Vector multiplication
4035///
4036/// # Example
4037/// ```
4038/// #[macro_use]
4039/// extern crate peroxide;
4040/// use peroxide::fuga::*;
4041///
4042/// fn main() {
4043/// let a = ml_matrix("1 2 3; 4 5 6");
4044/// let b = c!(1, 2, 3);
4045/// let mut c = vec![0f64; 2];
4046/// gemv(1f64, &a, &b, 0f64, &mut c);
4047/// assert_eq!(c, c!(14, 32));
4048/// }
4049/// ```
4050pub fn gemv(alpha: f64, a: &Matrix, b: &Vec<f64>, beta: f64, c: &mut Vec<f64>) {
4051 let m = a.row;
4052 let k = a.col;
4053 let n = 1usize;
4054 let (rsa, csa) = match a.shape {
4055 Row => (a.col as isize, 1isize),
4056 Col => (1isize, a.row as isize),
4057 };
4058 let (rsb, csb) = (1isize, 1isize);
4059 let (rsc, csc) = (1isize, 1isize);
4060
4061 unsafe {
4062 matrixmultiply::dgemm(
4063 m,
4064 k,
4065 n,
4066 alpha,
4067 a.ptr(),
4068 rsa,
4069 csa,
4070 b.as_ptr(),
4071 rsb,
4072 csb,
4073 beta,
4074 c.as_mut_ptr(),
4075 rsc,
4076 csc,
4077 )
4078 }
4079}
4080
4081/// General Vector-Matrix multiplication
4082///
4083/// # Example
4084/// ```
4085/// #[macro_use]
4086/// extern crate peroxide;
4087/// use peroxide::fuga::*;
4088///
4089/// fn main() {
4090/// let a = c!(1, 2);
4091/// let b = ml_matrix("1 2 3; 4 5 6");
4092/// let mut c = vec![0f64; 3];
4093/// gevm(1f64, &a, &b, 0f64, &mut c);
4094/// assert_eq!(c, c!(9, 12, 15));
4095/// }
4096/// ```
4097pub fn gevm(alpha: f64, a: &Vec<f64>, b: &Matrix, beta: f64, c: &mut Vec<f64>) {
4098 let m = 1usize;
4099 let k = a.len();
4100 let n = b.col;
4101 let (rsa, csa) = (1isize, 1isize);
4102 let (rsb, csb) = match b.shape {
4103 Row => (b.col as isize, 1isize),
4104 Col => (1isize, b.row as isize),
4105 };
4106 let (rsc, csc) = (1isize, 1isize);
4107
4108 unsafe {
4109 matrixmultiply::dgemm(
4110 m,
4111 k,
4112 n,
4113 alpha,
4114 a.as_ptr(),
4115 rsa,
4116 csa,
4117 b.ptr(),
4118 rsb,
4119 csb,
4120 beta,
4121 c.as_mut_ptr(),
4122 rsc,
4123 csc,
4124 )
4125 }
4126}
4127
4128//fn matmul(a: &Matrix, b: &Matrix) -> Matrix {
4129// match (a.row, a.col) {
4130// (p, q) if p <= 100 && q <= 100 => {
4131// let r_self = a.row;
4132// let c_self = a.col;
4133// let new_other = b;
4134// let r_other = new_other.row;
4135// let c_other = new_other.col;
4136//
4137// assert_eq!(c_self, r_other);
4138//
4139// let r_new = r_self;
4140// let c_new = c_other;
4141//
4142// let mut result = matrix(vec![0f64; r_new * c_new], r_new, c_new, a.shape);
4143//
4144// for i in 0..r_new {
4145// for j in 0..c_new {
4146// let mut s = 0f64;
4147// for k in 0..c_self {
4148// s += a[(i, k)] * new_other[(k, j)];
4149// }
4150// result[(i, j)] = s;
4151// }
4152// }
4153// result
4154// }
4155// _ => {
4156// let (a1, a2, a3, a4) = a.block();
4157// let (b1, b2, b3, b4) = b.block();
4158//
4159// let m1 = matmul(&a1, &b1) + matmul(&a2, &b3);
4160// let m2 = matmul(&a1, &b2) + matmul(&a2, &b4);
4161// let m3 = matmul(&a3, &b1) + matmul(&a4, &b3);
4162// let m4 = matmul(&a3, &b2) + matmul(&a4, &b4);
4163//
4164// combine(m1, m2, m3, m4)
4165// }
4166// }
4167//}
4168
4169// =============================================================================
4170// BLAS & LAPACK Area
4171// =============================================================================
4172
4173/// Matrix multiplication with BLAS
4174///
4175/// * m1: m x k matrix
4176/// * m2: k x n matrix
4177/// * result: m x n matrix
4178#[cfg(feature = "O3")]
4179pub fn blas_mul(m1: &Matrix, m2: &Matrix) -> Matrix {
4180 let m = m1.row;
4181 let k = m1.col;
4182 assert_eq!(k, m2.row);
4183 let n = m2.col;
4184
4185 let m_i32 = m as i32;
4186 let n_i32 = n as i32;
4187 let k_i32 = k as i32;
4188
4189 let a = &m1.data;
4190 let b = &m2.data;
4191 let mut c = vec![0f64; m * n];
4192
4193 match (m1.shape, m2.shape) {
4194 (Row, Row) => {
4195 unsafe {
4196 dgemm(
4197 b'N', b'N', n_i32, m_i32, k_i32, 1f64, b, n_i32, a, k_i32, 0f64, &mut c, n_i32,
4198 );
4199 }
4200 matrix(c, m, n, Row)
4201 }
4202 (Row, Col) => {
4203 unsafe {
4204 dgemm(
4205 b'T', b'N', n_i32, m_i32, k_i32, 1f64, b, k_i32, a, k_i32, 0f64, &mut c, n_i32,
4206 );
4207 }
4208 matrix(c, m, n, Row)
4209 }
4210 (Col, Col) => {
4211 unsafe {
4212 // (nxk) x (kxm) = n x m
4213 dgemm(
4214 b'N', b'N', m_i32, n_i32, k_i32, 1f64, a, m_i32, b, k_i32, 0f64, &mut c, m_i32,
4215 );
4216 }
4217 matrix(c, m, n, Col)
4218 }
4219 (Col, Row) => {
4220 unsafe {
4221 dgemm(
4222 b'N', b'T', m_i32, n_i32, k_i32, 1f64, a, m_i32, b, n_i32, 0f64, &mut c, m_i32,
4223 );
4224 }
4225 matrix(c, m, n, Col)
4226 }
4227 }
4228}
4229
4230#[allow(non_camel_case_types)]
4231#[derive(Debug, Copy, Clone, Eq, PartialEq)]
4232pub enum LAPACK_STATUS {
4233 Singular,
4234 NonSingular,
4235}
4236
4237#[allow(non_camel_case_types)]
4238#[derive(Debug, Copy, Clone, Eq, PartialEq)]
4239pub enum SVD_STATUS {
4240 Success,
4241 Diverge(i32),
4242}
4243
4244#[allow(non_camel_case_types)]
4245#[derive(Debug, Copy, Clone, Eq, PartialEq)]
4246pub enum POSITIVE_STATUS {
4247 Success,
4248 Failed(i32),
4249}
4250
4251/// Temporary data structure from `dgetrf`
4252#[derive(Debug, Clone)]
4253pub struct DGETRF {
4254 pub fact_mat: Matrix,
4255 pub ipiv: Vec<i32>,
4256 pub status: LAPACK_STATUS,
4257}
4258
4259/// Temporary data structure from `dgeqrf`
4260#[derive(Debug, Clone)]
4261pub struct DGEQRF {
4262 pub fact_mat: Matrix,
4263 pub tau: Vec<f64>,
4264 pub status: LAPACK_STATUS,
4265}
4266
4267#[derive(Debug, Clone)]
4268pub struct DGESVD {
4269 pub s: Vec<f64>,
4270 pub u: Matrix,
4271 pub vt: Matrix,
4272 pub status: SVD_STATUS,
4273}
4274
4275#[derive(Debug, Clone)]
4276pub struct DPOTRF {
4277 pub fact_mat: Matrix,
4278 pub uplo: UPLO,
4279 pub status: POSITIVE_STATUS,
4280}
4281
4282///// Temporary data structure from `dgeev`
4283//#[derive(Debug, Clone)]
4284//pub struct DGEEV {
4285// pub fact_mat: Matrix,
4286// pub tau: Vec<f64>,
4287// pub status: LAPACK_STATUS,
4288//}
4289
4290/// Peroxide version of `dgetrf`
4291#[cfg(feature = "O3")]
4292pub fn lapack_dgetrf(mat: &Matrix) -> Option<DGETRF> {
4293 let m = mat.row;
4294 let n = mat.col;
4295 let m_i32 = m as i32;
4296 let n_i32 = n as i32;
4297 // Should column major
4298 let mut a = match mat.shape {
4299 Row => mat.change_shape().data.clone(),
4300 Col => mat.data.clone(),
4301 };
4302
4303 let mut info = 0i32;
4304 let mut ipiv: Vec<i32> = vec![0i32; min(m, n)];
4305
4306 unsafe {
4307 dgetrf(m_i32, n_i32, &mut a, m_i32, &mut ipiv, &mut info);
4308 }
4309
4310 if info < 0 {
4311 None
4312 } else if info == 0 {
4313 Some(DGETRF {
4314 fact_mat: matrix(a, m, n, Col),
4315 ipiv,
4316 status: LAPACK_STATUS::NonSingular,
4317 })
4318 } else {
4319 Some(DGETRF {
4320 fact_mat: matrix(a, m, n, Col),
4321 ipiv,
4322 status: LAPACK_STATUS::Singular,
4323 })
4324 }
4325}
4326
4327/// Peroxide version of `dgetri`
4328#[cfg(feature = "O3")]
4329pub fn lapack_dgetri(dgrf: &DGETRF) -> Option<Matrix> {
4330 let mut result = dgrf.fact_mat.clone();
4331 let ipiv = &dgrf.ipiv;
4332 let (n, lda) = (result.col as i32, result.row as i32);
4333 let mut info = 0i32;
4334 let mut work = vec![0f64; result.col];
4335
4336 // Workspace Query
4337 unsafe {
4338 dgetri(n, &mut result.data, lda, ipiv, &mut work, -1, &mut info);
4339 }
4340
4341 let optimal_lwork = work[0] as usize;
4342 let mut optimal_work = vec![0f64; optimal_lwork];
4343
4344 // Real dgetri
4345 unsafe {
4346 dgetri(
4347 n,
4348 &mut result.data,
4349 lda,
4350 ipiv,
4351 &mut optimal_work,
4352 optimal_lwork as i32,
4353 &mut info,
4354 );
4355 }
4356
4357 if info == 0 {
4358 Some(result)
4359 } else {
4360 None
4361 }
4362}
4363
4364/// Peroxide version of `dgetrs`
4365#[allow(non_snake_case)]
4366#[cfg(feature = "O3")]
4367pub fn lapack_dgetrs(dgrf: &DGETRF, b: &Matrix) -> Option<Matrix> {
4368 match b.shape {
4369 Row => lapack_dgetrs(dgrf, &b.change_shape()),
4370 Col => {
4371 let A = &dgrf.fact_mat;
4372 let mut b_vec = b.data.clone();
4373 let ipiv = &dgrf.ipiv;
4374 let n = A.col as i32;
4375 let nrhs = b.col as i32;
4376 let lda = A.row as i32;
4377 let ldb = b.row as i32;
4378 let mut info = 0i32;
4379
4380 unsafe {
4381 dgetrs(
4382 b'N', n, nrhs, &A.data, lda, ipiv, &mut b_vec, ldb, &mut info,
4383 );
4384 }
4385
4386 if info == 0 {
4387 Some(matrix(b_vec, A.col, b.col, Col))
4388 } else {
4389 None
4390 }
4391 }
4392 }
4393}
4394
4395/// Peroxide version of `dgeqrf`
4396#[allow(non_snake_case)]
4397#[cfg(feature = "O3")]
4398pub fn lapack_dgeqrf(mat: &Matrix) -> Option<DGEQRF> {
4399 match mat.shape {
4400 Row => lapack_dgeqrf(&mat.change_shape()),
4401 Col => {
4402 let m = mat.row as i32;
4403 let n = mat.col as i32;
4404 let mut A = mat.clone();
4405 let mut tau = vec![0f64; min(mat.row, mat.col)];
4406 let mut work = vec![0f64; mat.col];
4407 let mut info = 0i32;
4408
4409 // Workspace query
4410 unsafe {
4411 dgeqrf(m, n, &mut A.data, m, &mut tau, &mut work, -1, &mut info);
4412 }
4413
4414 let optimal_lwork = work[0] as usize;
4415 let mut optimal_work = vec![0f64; optimal_lwork];
4416
4417 // Real dgeqrf
4418 unsafe {
4419 dgeqrf(
4420 m,
4421 n,
4422 &mut A.data,
4423 m,
4424 &mut tau,
4425 &mut optimal_work,
4426 optimal_lwork as i32,
4427 &mut info,
4428 );
4429 }
4430
4431 if info == 0 {
4432 Some(DGEQRF {
4433 fact_mat: A,
4434 tau,
4435 status: LAPACK_STATUS::NonSingular,
4436 })
4437 } else if info > 0 {
4438 Some(DGEQRF {
4439 fact_mat: A,
4440 tau,
4441 status: LAPACK_STATUS::Singular,
4442 })
4443 } else {
4444 None
4445 }
4446 }
4447 }
4448}
4449
4450#[allow(non_snake_case)]
4451#[cfg(feature = "O3")]
4452pub fn lapack_dgesvd(mat: &Matrix) -> Option<DGESVD> {
4453 match mat.shape {
4454 Row => lapack_dgesvd(&mat.change_shape()),
4455 Col => {
4456 let jobu = b'A';
4457 let jobvt = b'A';
4458 let m = mat.row as i32;
4459 let n = mat.col as i32;
4460 let mut A = mat.clone();
4461 let lda = m;
4462 let mut s = vec![0f64; m.min(n) as usize];
4463 let ldu = m;
4464 let mut u = vec![0f64; (ldu * m) as usize];
4465 let ldvt = n;
4466 let mut vt = vec![0f64; (ldvt * n) as usize];
4467 let mut work = vec![0f64; mat.col];
4468 let lwork = -1i32;
4469 let mut info = 0i32;
4470
4471 // Workspace query
4472 unsafe {
4473 dgesvd(
4474 jobu,
4475 jobvt,
4476 m,
4477 n,
4478 &mut A.data,
4479 lda,
4480 &mut s,
4481 &mut u,
4482 ldu,
4483 &mut vt,
4484 ldvt,
4485 &mut work,
4486 lwork,
4487 &mut info,
4488 );
4489 }
4490
4491 let optimal_lwork = work[0] as usize;
4492 let mut optimal_work = vec![0f64; optimal_lwork];
4493
4494 // Real dgesvd
4495 unsafe {
4496 dgesvd(
4497 jobu,
4498 jobvt,
4499 m,
4500 n,
4501 &mut A.data,
4502 lda,
4503 &mut s,
4504 &mut u,
4505 ldu,
4506 &mut vt,
4507 ldvt,
4508 &mut optimal_work,
4509 optimal_lwork as i32,
4510 &mut info,
4511 )
4512 }
4513
4514 if info == 0 {
4515 Some(DGESVD {
4516 s,
4517 u: matrix(u, m as usize, m as usize, Col),
4518 vt: matrix(vt, n as usize, n as usize, Col),
4519 status: SVD_STATUS::Success,
4520 })
4521 } else if info < 0 {
4522 None
4523 } else {
4524 Some(DGESVD {
4525 s,
4526 u: matrix(u, m as usize, m as usize, Col),
4527 vt: matrix(vt, n as usize, n as usize, Col),
4528 status: SVD_STATUS::Diverge(info),
4529 })
4530 }
4531 }
4532 }
4533}
4534
4535#[allow(non_snake_case)]
4536#[cfg(feature = "O3")]
4537pub fn lapack_dpotrf(mat: &Matrix, UPLO: UPLO) -> Option<DPOTRF> {
4538 match mat.shape {
4539 Row => lapack_dpotrf(&mat.change_shape(), UPLO),
4540 Col => {
4541 let lda = mat.row as i32;
4542 let N = mat.col as i32;
4543 let mut A = mat.clone();
4544 let mut info = 0i32;
4545 let uplo = match UPLO {
4546 UPLO::Upper => b'U',
4547 UPLO::Lower => b'L',
4548 };
4549
4550 unsafe { dpotrf(uplo, N, &mut A.data, lda, &mut info) }
4551
4552 if info == 0 {
4553 Some(DPOTRF {
4554 fact_mat: matrix(A.data, mat.row, mat.col, Col),
4555 uplo: UPLO,
4556 status: POSITIVE_STATUS::Success,
4557 })
4558 } else if info > 0 {
4559 Some(DPOTRF {
4560 fact_mat: matrix(A.data, mat.row, mat.col, Col),
4561 uplo: UPLO,
4562 status: POSITIVE_STATUS::Failed(info),
4563 })
4564 } else {
4565 None
4566 }
4567 }
4568 }
4569}
4570
4571#[allow(non_snake_case)]
4572#[cfg(feature = "O3")]
4573impl DGETRF {
4574 pub fn get_P(&self) -> Vec<i32> {
4575 self.ipiv.clone()
4576 }
4577
4578 pub fn get_L(&self) -> Matrix {
4579 let mut l = self.fact_mat.clone();
4580 for i in 0..l.row {
4581 l[(i, i)] = 1f64;
4582 for j in i + 1..l.col {
4583 l[(i, j)] = 0f64;
4584 }
4585 }
4586 l
4587 }
4588
4589 pub fn get_U(&self) -> Matrix {
4590 let mut u = self.fact_mat.clone();
4591 for i in 0..u.row {
4592 for j in 0..i {
4593 u[(i, j)] = 0f64;
4594 }
4595 }
4596 u
4597 }
4598
4599 pub fn get_cond(&self) -> Option<f64> {
4600 let A = &self.fact_mat;
4601 let lda = A.row as i32;
4602 let n = A.col as i32;
4603 let anorm = A.norm(Norm::L1);
4604 let mut work = vec![0f64; 4 * A.col];
4605 let mut iwork = vec![0i32; A.col];
4606 let mut info = 0i32;
4607 let mut rcond = 0f64;
4608
4609 unsafe {
4610 dgecon(
4611 b'1', n, &A.data, lda, anorm, &mut rcond, &mut work, &mut iwork, &mut info,
4612 );
4613 }
4614
4615 if info == 0 {
4616 Some(rcond)
4617 } else {
4618 None
4619 }
4620 }
4621}
4622
4623#[allow(non_snake_case)]
4624#[cfg(feature = "O3")]
4625impl DGEQRF {
4626 pub fn get_Q(&self) -> Matrix {
4627 let mut A = self.fact_mat.clone();
4628 let m = A.row as i32;
4629 let n = A.col as i32;
4630 let k = min(m, n);
4631 let lda = m;
4632 let tau = &self.tau;
4633 let mut lwork = -1i32;
4634 let mut work = vec![0f64; 1];
4635 let mut info = 0i32;
4636
4637 // Optimize
4638 unsafe {
4639 dorgqr(m, n, k, &mut A.data, lda, tau, &mut work, lwork, &mut info);
4640 }
4641
4642 lwork = work[0] as i32;
4643 work = vec![0f64; lwork as usize];
4644
4645 // Real dorgqr
4646 unsafe {
4647 dorgqr(m, n, k, &mut A.data, lda, tau, &mut work, lwork, &mut info);
4648 }
4649
4650 A
4651 }
4652
4653 pub fn get_R(&self) -> Matrix {
4654 let qr = &self.fact_mat;
4655 let row = min(qr.row, qr.col);
4656 let col = qr.col;
4657 let mut result = matrix(vec![0f64; row * col], row, col, qr.shape);
4658 for i in 0..row {
4659 for j in i..col {
4660 result[(i, j)] = qr[(i, j)];
4661 }
4662 }
4663 result
4664 }
4665}
4666
4667#[allow(non_snake_case)]
4668impl DPOTRF {
4669 pub fn get_U(&self) -> Option<Matrix> {
4670 if self.uplo == UPLO::Lower {
4671 return None;
4672 }
4673
4674 let mat = &self.fact_mat;
4675 let n = mat.col;
4676 let mut result = matrix(vec![0f64; n.pow(2)], n, n, mat.shape);
4677 for i in 0..n {
4678 for j in i..n {
4679 result[(i, j)] = mat[(i, j)];
4680 }
4681 }
4682 Some(result)
4683 }
4684
4685 pub fn get_L(&self) -> Option<Matrix> {
4686 if self.uplo == UPLO::Upper {
4687 return None;
4688 }
4689
4690 let mat = &self.fact_mat;
4691 let n = mat.col;
4692 let mut result = matrix(vec![0f64; n.pow(2)], n, n, mat.shape);
4693
4694 for i in 0..n {
4695 for j in 0..i + 1 {
4696 result[(i, j)] = mat[(i, j)];
4697 }
4698 }
4699 Some(result)
4700 }
4701}
4702
4703#[allow(non_snake_case)]
4704pub fn gen_householder(a: &Vec<f64>) -> Matrix {
4705 let mut v = a.fmap(|t| t / (a[0] + a.norm(Norm::L2) * a[0].signum()));
4706 v[0] = 1f64;
4707 let mut H = eye(a.len());
4708 let vt: Matrix = v.clone().into();
4709 H = H - 2f64 / v.dot(&v) * (&vt * &vt.t());
4710 H
4711}
4712
4713/// LU via Gaussian Elimination with Partial Pivoting
4714#[allow(dead_code)]
4715fn gepp(m: &mut Matrix) -> Vec<usize> {
4716 let mut r = vec![0usize; m.col - 1];
4717 for k in 0..(m.col - 1) {
4718 // Find the pivot row
4719 let r_k = m
4720 .col(k)
4721 .into_iter()
4722 .skip(k)
4723 .enumerate()
4724 .max_by(|x1, x2| x1.1.abs().partial_cmp(&x2.1.abs()).unwrap())
4725 .unwrap()
4726 .0
4727 + k;
4728 r[k] = r_k;
4729
4730 // Interchange the rows r_k and k
4731 for j in k..m.col {
4732 unsafe {
4733 std::ptr::swap(&mut m[(k, j)], &mut m[(r_k, j)]);
4734 println!("Swap! k:{}, r_k:{}", k, r_k);
4735 }
4736 }
4737 // Form the multipliers
4738 for i in k + 1..m.col {
4739 m[(i, k)] = -m[(i, k)] / m[(k, k)];
4740 }
4741 // Update the entries
4742 for i in k + 1..m.col {
4743 for j in k + 1..m.col {
4744 m[(i, j)] += m[(i, k)] * m[(k, j)];
4745 }
4746 }
4747 }
4748 r
4749}
4750
4751/// LU via Gauss Elimination with Complete Pivoting
4752fn gecp(m: &mut Matrix) -> (Vec<usize>, Vec<usize>) {
4753 let n = m.col;
4754 let mut r = vec![0usize; n - 1];
4755 let mut s = vec![0usize; n - 1];
4756 for k in 0..n - 1 {
4757 // Find pivot
4758 let (r_k, s_k) = match m.shape {
4759 Col => {
4760 let mut row_ics = 0usize;
4761 let mut col_ics = 0usize;
4762 let mut max_val = 0f64;
4763 for i in k..n {
4764 let c = m
4765 .col(i)
4766 .into_iter()
4767 .skip(k)
4768 .enumerate()
4769 .max_by(|x1, x2| x1.1.abs().partial_cmp(&x2.1.abs()).unwrap())
4770 .unwrap();
4771 let c_ics = c.0 + k;
4772 let c_val = c.1.abs();
4773 if c_val > max_val {
4774 row_ics = c_ics;
4775 col_ics = i;
4776 max_val = c_val;
4777 }
4778 }
4779 (row_ics, col_ics)
4780 }
4781 Row => {
4782 let mut row_ics = 0usize;
4783 let mut col_ics = 0usize;
4784 let mut max_val = 0f64;
4785 for i in k..n {
4786 let c = m
4787 .row(i)
4788 .into_iter()
4789 .skip(k)
4790 .enumerate()
4791 .max_by(|x1, x2| x1.1.abs().partial_cmp(&x2.1.abs()).unwrap())
4792 .unwrap();
4793 let c_ics = c.0 + k;
4794 let c_val = c.1.abs();
4795 if c_val > max_val {
4796 col_ics = c_ics;
4797 row_ics = i;
4798 max_val = c_val;
4799 }
4800 }
4801 (row_ics, col_ics)
4802 }
4803 };
4804 r[k] = r_k;
4805 s[k] = s_k;
4806
4807 // Interchange rows
4808 for j in k..n {
4809 unsafe {
4810 std::ptr::swap(&mut m[(k, j)], &mut m[(r_k, j)]);
4811 }
4812 }
4813
4814 // Interchange cols
4815 for i in 0..n {
4816 unsafe {
4817 std::ptr::swap(&mut m[(i, k)], &mut m[(i, s_k)]);
4818 }
4819 }
4820
4821 // Form the multipliers
4822 for i in k + 1..n {
4823 m[(i, k)] = -m[(i, k)] / m[(k, k)];
4824 for j in k + 1..n {
4825 m[(i, j)] += m[(i, k)] * m[(k, j)];
4826 }
4827 }
4828 }
4829 (r, s)
4830}