Matrix

Struct Matrix 

Source
pub struct Matrix<T> { /* private fields */ }
Expand description

The Matrix struct.

Can be instantiated with any type.

Implementations§

Source§

impl<T> Matrix<T>
where T: Any + Float,

Source

pub fn qr_decomp(self) -> Result<(Matrix<T>, Matrix<T>), Error>

Compute the QR decomposition of the matrix.

Returns the tuple (Q,R).

§Examples
use rulinalg::matrix::Matrix;

let m = Matrix::new(3,3, vec![1.0,0.5,0.5,0.5,1.0,0.5,0.5,0.5,1.0]);

let (q, r) = m.qr_decomp().unwrap();
§Failures
  • Cannot compute the QR decomposition.
Source§

impl<T> Matrix<T>
where T: Any + Float,

Source

pub fn cholesky(&self) -> Result<Matrix<T>, Error>

Cholesky decomposition

Returns the cholesky decomposition of a positive definite matrix.

§Examples
use rulinalg::matrix::Matrix;

let m = Matrix::new(3,3, vec![1.0,0.5,0.5,0.5,1.0,0.5,0.5,0.5,1.0]);

let l = m.cholesky();
§Panics
  • The matrix is not square.
§Failures
  • Matrix is not positive definite.
Source§

impl<T> Matrix<T>
where T: Any + Float,

Source

pub fn bidiagonal_decomp( self, ) -> Result<(Matrix<T>, Matrix<T>, Matrix<T>), Error>

Converts matrix to bidiagonal form

Returns (B, U, V), where B is bidiagonal and self = U B V_T.

Note that if self has self.rows() > self.cols() the matrix will be transposed and then reduced - this will lead to a sub-diagonal instead of super-diagonal.

§Failures
  • The matrix cannot be reduced to bidiagonal form.
Source§

impl<T> Matrix<T>
where T: Any + Float + Signed,

Source

pub fn svd(self) -> Result<(Matrix<T>, Matrix<T>, Matrix<T>), Error>

Singular Value Decomposition

Computes the SVD using the Golub-Reinsch algorithm.

Returns Σ, U, V, such that self = U Σ VT. Σ is a diagonal matrix whose elements correspond to the non-negative singular values of the matrix. The singular values are ordered in non-increasing order. U and V have orthonormal columns, and each column represents the left and right singular vectors for the corresponding singular value in Σ, respectively.

If self has M rows and N columns, the dimensions of the returned matrices are as follows.

If M >= N:

  • Σ: N x N
  • U: M x N
  • V: N x N

If M < N:

  • Σ: M x M
  • U: M x M
  • V: N x M

Note: This version of the SVD is sometimes referred to as the ‘economy SVD’.

§Failures

This function may fail in some cases. The current decomposition whilst being efficient is fairly basic. Hopefully the algorithm can be made not to fail in the near future.

Source§

impl<T> Matrix<T>
where T: Any + Float,

Source

pub fn upper_hessenberg(self) -> Result<Matrix<T>, Error>

Returns H, where H is the upper hessenberg form.

If the transformation matrix is also required, you should use upper_hess_decomp.

§Examples
use rulinalg::matrix::Matrix;

let a = Matrix::new(4,4,vec![2.,0.,1.,1.,2.,0.,1.,2.,1.,2.,0.,0.,2.,0.,1.,1.]);
let h = a.upper_hessenberg();

println!("{:?}", h.expect("Could not get upper Hessenberg form.").data());
§Panics
  • The matrix is not square.
§Failures
  • The matrix cannot be reduced to upper hessenberg form.
Source

pub fn upper_hess_decomp(self) -> Result<(Matrix<T>, Matrix<T>), Error>

Returns (U,H), where H is the upper hessenberg form and U is the unitary transform matrix.

Note: The current transform matrix seems broken…

§Examples
use rulinalg::matrix::{Matrix, BaseMatrix};

let a = Matrix::new(3,3,vec![1.,2.,3.,4.,5.,6.,7.,8.,9.]);

// u is the transform, h is the upper hessenberg form.
let (u,h) = a.clone().upper_hess_decomp().expect("This matrix should decompose!");

println!("The hess : {:?}", h.data());
println!("Manual hess : {:?}", (u.transpose() * a * u).data());
§Panics
  • The matrix is not square.
§Failures
  • The matrix cannot be reduced to upper hessenberg form.
Source§

impl<T> Matrix<T>
where T: Any + Float,

Source

pub fn lup_decomp(&self) -> Result<(Matrix<T>, Matrix<T>, Matrix<T>), Error>

Computes L, U, and P for LUP decomposition.

Returns L,U, and P respectively.

§Examples
use rulinalg::matrix::Matrix;

let a = Matrix::new(3,3, vec![1.0,2.0,0.0,
                              0.0,3.0,4.0,
                              5.0, 1.0, 2.0]);

let (l,u,p) = a.lup_decomp().expect("This matrix should decompose!");
§Panics
  • Matrix is not square.
§Failures
  • Matrix cannot be LUP decomposed.
Source§

impl<T> Matrix<T>
where T: Any + Float + Signed,

Source

pub fn eigenvalues(&self) -> Result<Vec<T>, Error>

Eigenvalues of a square matrix.

Returns a Vec of eigenvalues.

§Examples
use rulinalg::matrix::Matrix;

let a = Matrix::new(4,4, (1..17).map(|v| v as f64).collect::<Vec<f64>>());
let e = a.eigenvalues().expect("We should be able to compute these eigenvalues!");
println!("{:?}", e);
§Panics
  • The matrix is not square.
§Failures
  • Eigenvalues cannot be computed.
Source

pub fn eigendecomp(&self) -> Result<(Vec<T>, Matrix<T>), Error>

Eigendecomposition of a square matrix.

Returns a Vec of eigenvalues, and a matrix with eigenvectors as the columns.

The eigenvectors are only gauranteed to be correct if the matrix is real-symmetric.

§Examples
use rulinalg::matrix::Matrix;

let a = Matrix::new(3,3,vec![3.,2.,4.,2.,0.,2.,4.,2.,3.]);

let (e, m) = a.eigendecomp().expect("We should be able to compute this eigendecomp!");
println!("{:?}", e);
println!("{:?}", m.data());
§Panics
  • The matrix is not square.
§Failures
  • The eigen decomposition can not be computed.
Source§

impl<T> Matrix<T>

Source

pub fn new<U>(rows: usize, cols: usize, data: U) -> Matrix<T>
where U: Into<Vec<T>>,

Constructor for Matrix struct.

Requires both the row and column dimensions.

§Examples
use rulinalg::matrix::{Matrix, BaseMatrix};

let mat = Matrix::new(2,2, vec![1.0,2.0,3.0,4.0]);

assert_eq!(mat.rows(), 2);
assert_eq!(mat.cols(), 2);
§Panics
  • The input data does not match the given dimensions.
Examples found in repository?
examples/k-means_generating_cluster.rs (lines 39-41)
12fn generate_data(centroids: &Matrix<f64>,
13                 points_per_centroid: usize,
14                 noise: f64)
15                 -> Matrix<f64> {
16    assert!(centroids.cols() > 0, "Centroids cannot be empty.");
17    assert!(centroids.rows() > 0, "Centroids cannot be empty.");
18    assert!(noise >= 0f64, "Noise must be non-negative.");
19    let mut raw_cluster_data = Vec::with_capacity(centroids.rows() * points_per_centroid *
20                                                  centroids.cols());
21
22    let mut rng = thread_rng();
23    let normal_rv = Normal::new(0f64, noise);
24
25    for _ in 0..points_per_centroid {
26        // Generate points from each centroid
27        for centroid in centroids.iter_rows() {
28            // Generate a point randomly around the centroid
29            let mut point = Vec::with_capacity(centroids.cols());
30            for feature in centroid {
31                point.push(feature + normal_rv.ind_sample(&mut rng));
32            }
33
34            // Push point to raw_cluster_data
35            raw_cluster_data.extend(point);
36        }
37    }
38
39    Matrix::new(centroids.rows() * points_per_centroid,
40                centroids.cols(),
41                raw_cluster_data)
42}
43
44fn main() {
45    println!("K-Means clustering example:");
46
47    const SAMPLES_PER_CENTROID: usize = 2000;
48
49    println!("Generating {0} samples from each centroids:",
50             SAMPLES_PER_CENTROID);
51    // Choose two cluster centers, at (-0.5, -0.5) and (0, 0.5).
52    let centroids = Matrix::new(2, 2, vec![-0.5, -0.5, 0.0, 0.5]);
53    println!("{}", centroids);
54
55    // Generate some data randomly around the centroids
56    let samples = generate_data(&centroids, SAMPLES_PER_CENTROID, 0.4);
57
58    // Create a new model with 2 clusters
59    let mut model = KMeansClassifier::new(2);
60
61    // Train the model
62    println!("Training the model...");
63    // Our train function returns a Result<(), E>
64    model.train(&samples).unwrap();
65
66    let centroids = model.centroids().as_ref().unwrap();
67    println!("Model Centroids:\n{:.3}", centroids);
68
69    // Predict the classes and partition into
70    println!("Classifying the samples...");
71    let classes = model.predict(&samples).unwrap();
72    let (first, second): (Vec<usize>, Vec<usize>) = classes.data().iter().partition(|&x| *x == 0);
73
74    println!("Samples closest to first centroid: {}", first.len());
75    println!("Samples closest to second centroid: {}", second.len());
76}
More examples
Hide additional examples
examples/svm-sign_learner.rs (lines 22-25)
17fn main() {
18    println!("Sign learner sample:");
19
20    println!("Training...");
21    // Training data
22    let inputs = Matrix::new(11, 1, vec![
23                             -0.1, -2., -9., -101., -666.7,
24                             0., 0.1, 1., 11., 99., 456.7
25                             ]);
26    let targets = Vector::new(vec![
27                              -1., -1., -1., -1., -1.,
28                              1., 1., 1., 1., 1., 1.
29                              ]);
30
31    // Trainee
32    let mut svm_mod = SVM::new(HyperTan::new(100., 0.), 0.3);
33    // Our train function returns a Result<(), E>
34    svm_mod.train(&inputs, &targets).unwrap();
35
36    println!("Evaluation...");
37    let mut hits = 0;
38    let mut misses = 0;
39    // Evaluation
40    //   Note: We could pass all input values at once to the `predict` method!
41    //         Here, we use a loop just to count and print logs.
42    for n in (-1000..1000).filter(|&x| x % 100 == 0) {
43        let nf = n as f64;
44        let input = Matrix::new(1, 1, vec![nf]);
45        let out = svm_mod.predict(&input).unwrap();
46        let res = if out[0] * nf > 0. {
47            hits += 1;
48            true
49        } else if nf == 0. {
50            hits += 1;
51            true
52        } else {
53            misses += 1;
54            false
55        };
56
57        println!("{} -> {}: {}", Matrix::data(&input)[0], out[0], res);
58    }
59
60    println!("Performance report:");
61    println!("Hits: {}, Misses: {}", hits, misses);
62    let hits_f = hits as f64;
63    let total = (hits + misses) as f64;
64    println!("Accuracy: {}", (hits_f / total) * 100.);
65}
examples/nnet-and_gate.rs (line 39)
15fn main() {
16    println!("AND gate learner sample:");
17
18    const THRESHOLD: f64 = 0.7;
19
20    const SAMPLES: usize = 10000;
21    println!("Generating {} training data and labels...", SAMPLES as u32);
22
23    let mut input_data = Vec::with_capacity(SAMPLES * 2);
24    let mut label_data = Vec::with_capacity(SAMPLES);
25
26    for _ in 0..SAMPLES {
27        // The two inputs are "signals" between 0 and 1
28        let Closed01(left) = random::<Closed01<f64>>();
29        let Closed01(right) = random::<Closed01<f64>>();
30        input_data.push(left);
31        input_data.push(right);
32        if left > THRESHOLD && right > THRESHOLD {
33            label_data.push(1.0);
34        } else {
35            label_data.push(0.0)
36        }
37    }
38
39    let inputs = Matrix::new(SAMPLES, 2, input_data);
40    let targets = Matrix::new(SAMPLES, 1, label_data);
41
42    let layers = &[2, 1];
43    let criterion = BCECriterion::new(Regularization::L2(0.));
44    let mut model = NeuralNet::new(layers, criterion, StochasticGD::default());
45
46    println!("Training...");
47    // Our train function returns a Result<(), E>
48    model.train(&inputs, &targets).unwrap();
49
50    let test_cases = vec![
51        0.0, 0.0,
52        0.0, 1.0,
53        1.0, 1.0,
54        1.0, 0.0,
55        ];
56    let expected = vec![
57        0.0,
58        0.0,
59        1.0,
60        0.0,
61        ];
62    let test_inputs = Matrix::new(test_cases.len() / 2, 2, test_cases);
63    let res = model.predict(&test_inputs).unwrap();
64
65    println!("Evaluation...");
66    let mut hits = 0;
67    let mut misses = 0;
68    // Evaluation
69    println!("Got\tExpected");
70    for (idx, prediction) in res.into_vec().iter().enumerate() {
71        println!("{:.2}\t{}", prediction, expected[idx]);
72        if (prediction - 0.5) * (expected[idx] - 0.5) > 0. {
73            hits += 1;
74        } else {
75            misses += 1;
76        }
77    }
78
79    println!("Hits: {}, Misses: {}", hits, misses);
80    let hits_f = hits as f64;
81    let total = (hits + misses) as f64;
82    println!("Accuracy: {}%", (hits_f / total) * 100.);
83}
Source

pub fn from_fn<F>(rows: usize, cols: usize, f: F) -> Matrix<T>
where F: FnMut(usize, usize) -> T,

Constructor for Matrix struct that takes a function f and constructs a new matrix such that A_ij = f(i, j), where i is the row index and j the column index.

Requires both the row and column dimensions as well as a generating function.

§Examples
use rulinalg::matrix::{Matrix, BaseMatrix};

// Let's assume you have an array of "things" for
// which you want to generate a distance matrix:
let things: [i32; 3] = [1, 2, 3];
let distances: Matrix<f64> = Matrix::from_fn(things.len(), things.len(), |col, row| {
    (things[col] - things[row]).abs().into()
});

assert_eq!(distances.rows(), 3);
assert_eq!(distances.cols(), 3);
assert_eq!(distances.data(), &vec![
    0.0, 1.0, 2.0,
    1.0, 0.0, 1.0,
    2.0, 1.0, 0.0,
]);
Source

pub fn data(&self) -> &Vec<T>

Returns a non-mutable reference to the underlying data.

Examples found in repository?
examples/svm-sign_learner.rs (line 57)
17fn main() {
18    println!("Sign learner sample:");
19
20    println!("Training...");
21    // Training data
22    let inputs = Matrix::new(11, 1, vec![
23                             -0.1, -2., -9., -101., -666.7,
24                             0., 0.1, 1., 11., 99., 456.7
25                             ]);
26    let targets = Vector::new(vec![
27                              -1., -1., -1., -1., -1.,
28                              1., 1., 1., 1., 1., 1.
29                              ]);
30
31    // Trainee
32    let mut svm_mod = SVM::new(HyperTan::new(100., 0.), 0.3);
33    // Our train function returns a Result<(), E>
34    svm_mod.train(&inputs, &targets).unwrap();
35
36    println!("Evaluation...");
37    let mut hits = 0;
38    let mut misses = 0;
39    // Evaluation
40    //   Note: We could pass all input values at once to the `predict` method!
41    //         Here, we use a loop just to count and print logs.
42    for n in (-1000..1000).filter(|&x| x % 100 == 0) {
43        let nf = n as f64;
44        let input = Matrix::new(1, 1, vec![nf]);
45        let out = svm_mod.predict(&input).unwrap();
46        let res = if out[0] * nf > 0. {
47            hits += 1;
48            true
49        } else if nf == 0. {
50            hits += 1;
51            true
52        } else {
53            misses += 1;
54            false
55        };
56
57        println!("{} -> {}: {}", Matrix::data(&input)[0], out[0], res);
58    }
59
60    println!("Performance report:");
61    println!("Hits: {}, Misses: {}", hits, misses);
62    let hits_f = hits as f64;
63    let total = (hits + misses) as f64;
64    println!("Accuracy: {}", (hits_f / total) * 100.);
65}
Source

pub fn mut_data(&mut self) -> &mut [T]

Returns a mutable slice of the underlying data.

Source

pub fn into_vec(self) -> Vec<T>

Consumes the Matrix and returns the Vec of data.

Examples found in repository?
examples/nnet-and_gate.rs (line 70)
15fn main() {
16    println!("AND gate learner sample:");
17
18    const THRESHOLD: f64 = 0.7;
19
20    const SAMPLES: usize = 10000;
21    println!("Generating {} training data and labels...", SAMPLES as u32);
22
23    let mut input_data = Vec::with_capacity(SAMPLES * 2);
24    let mut label_data = Vec::with_capacity(SAMPLES);
25
26    for _ in 0..SAMPLES {
27        // The two inputs are "signals" between 0 and 1
28        let Closed01(left) = random::<Closed01<f64>>();
29        let Closed01(right) = random::<Closed01<f64>>();
30        input_data.push(left);
31        input_data.push(right);
32        if left > THRESHOLD && right > THRESHOLD {
33            label_data.push(1.0);
34        } else {
35            label_data.push(0.0)
36        }
37    }
38
39    let inputs = Matrix::new(SAMPLES, 2, input_data);
40    let targets = Matrix::new(SAMPLES, 1, label_data);
41
42    let layers = &[2, 1];
43    let criterion = BCECriterion::new(Regularization::L2(0.));
44    let mut model = NeuralNet::new(layers, criterion, StochasticGD::default());
45
46    println!("Training...");
47    // Our train function returns a Result<(), E>
48    model.train(&inputs, &targets).unwrap();
49
50    let test_cases = vec![
51        0.0, 0.0,
52        0.0, 1.0,
53        1.0, 1.0,
54        1.0, 0.0,
55        ];
56    let expected = vec![
57        0.0,
58        0.0,
59        1.0,
60        0.0,
61        ];
62    let test_inputs = Matrix::new(test_cases.len() / 2, 2, test_cases);
63    let res = model.predict(&test_inputs).unwrap();
64
65    println!("Evaluation...");
66    let mut hits = 0;
67    let mut misses = 0;
68    // Evaluation
69    println!("Got\tExpected");
70    for (idx, prediction) in res.into_vec().iter().enumerate() {
71        println!("{:.2}\t{}", prediction, expected[idx]);
72        if (prediction - 0.5) * (expected[idx] - 0.5) > 0. {
73            hits += 1;
74        } else {
75            misses += 1;
76        }
77    }
78
79    println!("Hits: {}, Misses: {}", hits, misses);
80    let hits_f = hits as f64;
81    let total = (hits + misses) as f64;
82    println!("Accuracy: {}%", (hits_f / total) * 100.);
83}
Source§

impl<T> Matrix<T>
where T: Clone + Zero,

Source

pub fn zeros(rows: usize, cols: usize) -> Matrix<T>

Constructs matrix of all zeros.

Requires both the row and the column dimensions.

§Examples
use rulinalg::matrix::Matrix;

let mat = Matrix::<f64>::zeros(2,3);
Source

pub fn from_diag(diag: &[T]) -> Matrix<T>

Constructs matrix with given diagonal.

Requires slice of diagonal elements.

§Examples
use rulinalg::matrix::Matrix;

let mat = Matrix::from_diag(&vec![1.0,2.0,3.0,4.0]);
Source§

impl<T> Matrix<T>
where T: Clone + One,

Source

pub fn ones(rows: usize, cols: usize) -> Matrix<T>

Constructs matrix of all ones.

Requires both the row and the column dimensions.

§Examples
use rulinalg::matrix::Matrix;

let mat = Matrix::<f64>::ones(2,3);
Source§

impl<T> Matrix<T>
where T: Clone + Zero + One,

Source

pub fn identity(size: usize) -> Matrix<T>

Constructs the identity matrix.

Requires the size of the matrix.

§Examples
use rulinalg::matrix::Matrix;

let I = Matrix::<f64>::identity(4);
Source§

impl<T> Matrix<T>
where T: Float + FromPrimitive,

Source

pub fn mean(&self, axis: Axes) -> Vector<T>

The mean of the matrix along the specified axis.

  • Axis Row - Arithmetic mean of rows.
  • Axis Col - Arithmetic mean of columns.

Calling mean() on an empty matrix will return an empty matrix.

§Examples
use rulinalg::matrix::{Matrix, Axes};

let a = Matrix::<f64>::new(2,2, vec![1.0,2.0,3.0,4.0]);

let c = a.mean(Axes::Row);
assert_eq!(*c.data(), vec![2.0, 3.0]);

let d = a.mean(Axes::Col);
assert_eq!(*d.data(), vec![1.5, 3.5]);
Source

pub fn variance(&self, axis: Axes) -> Result<Vector<T>, Error>

The variance of the matrix along the specified axis.

  • Axis Row - Sample variance of rows.
  • Axis Col - Sample variance of columns.
§Examples
use rulinalg::matrix::{Matrix, Axes};

let a = Matrix::<f32>::new(2,2,vec![1.0,2.0,3.0,4.0]);

let c = a.variance(Axes::Row).unwrap();
assert_eq!(*c.data(), vec![2.0, 2.0]);

let d = a.variance(Axes::Col).unwrap();
assert_eq!(*d.data(), vec![0.5, 0.5]);
§Failures
  • There are one or fewer row/columns in the working axis.
Source§

impl<T> Matrix<T>
where T: Any + Float,

Source

pub fn solve(&self, y: Vector<T>) -> Result<Vector<T>, Error>

Solves the equation Ax = y.

Requires a Vector y as input.

§Examples
use rulinalg::matrix::Matrix;
use rulinalg::vector::Vector;

let a = Matrix::new(2,2, vec![2.0,3.0,1.0,2.0]);
let y = Vector::new(vec![13.0,8.0]);

let x = a.solve(y).unwrap();

assert_eq!(*x.data(), vec![2.0, 3.0]);
§Panics
  • The matrix column count and vector size are different.
  • The matrix is not square.
§Failures
  • The matrix cannot be decomposed into an LUP form to solve.
  • There is no valid solution as the matrix is singular.
Source

pub fn inverse(&self) -> Result<Matrix<T>, Error>

Computes the inverse of the matrix.

§Examples
use rulinalg::matrix::Matrix;

let a = Matrix::new(2,2, vec![2.,3.,1.,2.]);
let inv = a.inverse().expect("This matrix should have an inverse!");

let I = a * inv;

assert_eq!(*I.data(), vec![1.0,0.0,0.0,1.0]);
§Panics
  • The matrix is not square.
§Failures
  • The matrix could not be LUP decomposed.
  • The matrix has zero determinant.
Source

pub fn det(&self) -> T

Computes the determinant of the matrix.

§Examples
use rulinalg::matrix::Matrix;

let a = Matrix::new(3,3, vec![1.0,2.0,0.0,
                              0.0,3.0,4.0,
                              5.0, 1.0, 2.0]);

let det = a.det();
§Panics
  • The matrix is not square.

Trait Implementations§

Source§

impl<'a, 'b, T> Add<&'b Matrix<T>> for &'a Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between two matrices.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: &Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, 'c, T> Add<&'c Matrix<T>> for &'b MatrixSlice<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: &Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, 'c, T> Add<&'c Matrix<T>> for &'b MatrixSliceMut<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: &Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> Add<&'a Matrix<T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between two matrices.

This will reuse allocated memory from self.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: &Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<&'b Matrix<T>> for MatrixSlice<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: &Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<&'b Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: &Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, 'c, T> Add<&'c MatrixSlice<'a, T>> for &'b Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, s: &MatrixSlice<'_, T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<&'b MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, s: &MatrixSlice<'_, T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, 'c, T> Add<&'c MatrixSliceMut<'a, T>> for &'b Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, s: &MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<&'b MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, s: &MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<&'b T> for &'a Matrix<T>
where T: Copy + Add<Output = T>,

Scalar addition with matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, f: &T) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> Add<&'a T> for Matrix<T>
where T: Copy + Add<Output = T>,

Scalar addition with matrix.

Will reuse the memory allocated for the existing matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, f: &T) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> Add<Matrix<T>> for &'a Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between two matrices.

This will reuse allocated memory from m.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<Matrix<T>> for &'b MatrixSlice<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<Matrix<T>> for &'b MatrixSliceMut<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> Add<Matrix<T>> for MatrixSlice<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> Add<Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<MatrixSlice<'a, T>> for &'b Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, s: MatrixSlice<'_, T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> Add<MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, s: MatrixSlice<'_, T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, 'b, T> Add<MatrixSliceMut<'a, T>> for &'b Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, s: MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> Add<MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, s: MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> Add<T> for &'a Matrix<T>
where T: Copy + Add<Output = T>,

Scalar addition with matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, f: T) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<T> Add<T> for Matrix<T>
where T: Copy + Add<Output = T>,

Scalar addition with matrix.

Will reuse the memory allocated for the existing matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, f: T) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<T> Add for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition between two matrices.

This will reuse allocated memory from self.

Source§

type Output = Matrix<T>

The resulting type after applying the + operator.
Source§

fn add(self, m: Matrix<T>) -> Matrix<T>

Performs the + operation. Read more
Source§

impl<'a, T> AddAssign<&'a Matrix<T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition assignment between two matrices.

Source§

fn add_assign(&mut self, _rhs: &Matrix<T>)

Performs the += operation. Read more
Source§

impl<'a, 'b, T> AddAssign<&'b Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition assignment between two matrices.

Source§

fn add_assign(&mut self, _rhs: &Matrix<T>)

Performs the += operation. Read more
Source§

impl<'a, 'b, T> AddAssign<&'b MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition assignment between two matrices.

Source§

fn add_assign(&mut self, _rhs: &MatrixSlice<'_, T>)

Performs the += operation. Read more
Source§

impl<'a, 'b, T> AddAssign<&'b MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition assignment between two matrices.

Source§

fn add_assign(&mut self, _rhs: &MatrixSliceMut<'_, T>)

Performs the += operation. Read more
Source§

impl<'a, T> AddAssign<&'a T> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs addition assignment between a matrix and a scalar.

Source§

fn add_assign(&mut self, _rhs: &T)

Performs the += operation. Read more
Source§

impl<'a, T> AddAssign<Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Add<Output = T>,

Performs elementwise addition assignment between two matrices.

Source§

fn add_assign(&mut self, _rhs: Matrix<T>)

Performs the += operation. Read more
Source§

impl<'a, T> AddAssign<MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition assignment between two matrices.

Source§

fn add_assign(&mut self, _rhs: MatrixSlice<'_, T>)

Performs the += operation. Read more
Source§

impl<'a, T> AddAssign<MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition assignment between two matrices.

Source§

fn add_assign(&mut self, _rhs: MatrixSliceMut<'_, T>)

Performs the += operation. Read more
Source§

impl<T> AddAssign<T> for Matrix<T>
where T: Copy + Add<Output = T>,

Performs addition assignment between a matrix and a scalar.

Source§

fn add_assign(&mut self, _rhs: T)

Performs the += operation. Read more
Source§

impl<T> AddAssign for Matrix<T>
where T: Copy + Add<Output = T>,

Performs elementwise addition assignment between two matrices.

Source§

fn add_assign(&mut self, _rhs: Matrix<T>)

Performs the += operation. Read more
Source§

impl<T> BaseMatrix<T> for Matrix<T>

Source§

fn rows(&self) -> usize

Rows in the matrix.
Source§

fn cols(&self) -> usize

Columns in the matrix.
Source§

fn row_stride(&self) -> usize

Row stride in the matrix.
Source§

fn is_empty(&self) -> bool

Returns true if the matrix contais no elements
Source§

fn as_ptr(&self) -> *const T

Top left index of the matrix.
Source§

fn into_matrix(self) -> Matrix<T>
where T: Copy,

Convert the matrix struct into a owned Matrix.
Source§

fn sum(&self) -> T
where T: Copy + Zero<Output = T> + Add,

The sum of all elements in the matrix Read more
Source§

fn elemul(&self, m: &Matrix<T>) -> Matrix<T>
where T: Copy + Mul<Output = T>,

The elementwise product of two matrices. Read more
Source§

fn elediv(&self, m: &Matrix<T>) -> Matrix<T>
where T: Copy + Div<Output = T>,

The elementwise division of two matrices. Read more
Source§

fn vcat<S>(&self, m: &S) -> Matrix<T>
where T: Copy, S: BaseMatrix<T>,

Vertically concatenates two matrices. With self on top. Read more
Source§

fn as_slice(&self) -> MatrixSlice<'_, T>

Returns a MatrixSlice over the whole matrix. Read more
Source§

unsafe fn get_unchecked(&self, index: [usize; 2]) -> &T

Get a reference to a point in the matrix without bounds checking.
Source§

fn get_row(&self, index: usize) -> Option<&[T]>

Returns the row of a matrix at the given index. None if the index is out of bounds. Read more
Source§

unsafe fn get_row_unchecked(&self, index: usize) -> &[T]

Returns the row of a matrix at the given index without doing unbounds checking Read more
Source§

fn iter<'a>(&self) -> SliceIter<'a, T>
where T: 'a,

Returns an iterator over the matrix data. Read more
Source§

fn iter_rows(&self) -> Rows<'_, T>

Iterate over the rows of the matrix. Read more
Source§

fn iter_diag(&self, k: DiagOffset) -> Diagonal<'_, T, Self>

Iterate over diagonal entries Read more
Source§

fn sum_rows(&self) -> Vector<T>
where T: Copy + Zero<Output = T> + Add,

The sum of the rows of the matrix. Read more
Source§

fn sum_cols(&self) -> Vector<T>
where T: Copy + Zero<Output = T> + Add,

The sum of the columns of the matrix. Read more
Source§

fn select_rows<'a, I>(&self, rows: I) -> Matrix<T>
where T: Copy, I: IntoIterator<Item = &'a usize>, <I as IntoIterator>::IntoIter: ExactSizeIterator + Clone,

Select rows from matrix Read more
Source§

fn select_cols<'a, I>(&self, cols: I) -> Matrix<T>
where T: Copy, I: IntoIterator<Item = &'a usize>, <I as IntoIterator>::IntoIter: ExactSizeIterator + Clone,

Select columns from matrix Read more
Source§

fn select(&self, rows: &[usize], cols: &[usize]) -> Matrix<T>
where T: Copy,

Select block matrix from matrix Read more
Source§

fn hcat<S>(&self, m: &S) -> Matrix<T>
where T: Copy, S: BaseMatrix<T>,

Horizontally concatenates two matrices. With self on the left. Read more
Source§

fn diag(&self) -> Vector<T>
where T: Copy,

Extract the diagonal of the matrix Read more
Source§

fn transpose(&self) -> Matrix<T>
where T: Copy,

Tranposes the given matrix Read more
Source§

fn is_diag(&self) -> bool
where T: Zero + PartialEq,

Checks if matrix is diagonal. Read more
Source§

fn solve_u_triangular(&self, y: Vector<T>) -> Result<Vector<T>, Error>
where T: Any + Float,

Solves an upper triangular linear system. Read more
Source§

fn solve_l_triangular(&self, y: Vector<T>) -> Result<Vector<T>, Error>
where T: Any + Float,

Solves a lower triangular linear system. Read more
Source§

fn split_at( &self, mid: usize, axis: Axes, ) -> (MatrixSlice<'_, T>, MatrixSlice<'_, T>)

Split the matrix at the specified axis returning two MatrixSlices. Read more
Source§

fn sub_slice<'a>( &self, start: [usize; 2], rows: usize, cols: usize, ) -> MatrixSlice<'a, T>
where T: 'a,

Produce a MatrixSlice from an existing matrix. Read more
Source§

impl<T> BaseMatrixMut<T> for Matrix<T>

Source§

fn as_mut_ptr(&mut self) -> *mut T

Top left index of the slice.

Source§

fn as_mut_slice(&mut self) -> MatrixSliceMut<'_, T>

Returns a MatrixSliceMut over the whole matrix. Read more
Source§

unsafe fn get_unchecked_mut(&mut self, index: [usize; 2]) -> &mut T

Get a mutable reference to a point in the matrix without bounds checks.
Source§

fn iter_mut<'a>(&mut self) -> SliceIterMut<'a, T>
where T: 'a,

Returns a mutable iterator over the matrix. Read more
Source§

fn get_row_mut(&mut self, index: usize) -> Option<&mut [T]>

Returns a mutable reference to the row of a matrix at the given index. None if the index is out of bounds. Read more
Source§

unsafe fn get_row_unchecked_mut(&mut self, index: usize) -> &mut [T]

Returns a mutable reference to the row of a matrix at the given index without doing unbounds checking Read more
Source§

fn swap_rows(&mut self, a: usize, b: usize)

Swaps two rows in a matrix. Read more
Source§

fn swap_cols(&mut self, a: usize, b: usize)

Swaps two columns in a matrix. Read more
Source§

fn iter_rows_mut(&mut self) -> RowsMut<'_, T>

Iterate over the mutable rows of the matrix. Read more
Source§

fn iter_diag_mut(&mut self, k: DiagOffset) -> DiagonalMut<'_, T, Self>

Iterate over diagonal entries mutably Read more
Source§

fn set_to<M>(self, target: M)
where M: BaseMatrix<T>, T: Copy,

Sets the underlying matrix data to the target data. Read more
Source§

fn apply(self, f: &dyn Fn(T) -> T) -> Self
where T: Copy,

Applies a function to each element in the matrix. Read more
Source§

fn split_at_mut( &mut self, mid: usize, axis: Axes, ) -> (MatrixSliceMut<'_, T>, MatrixSliceMut<'_, T>)

Split the matrix at the specified axis returning two MatrixSliceMuts. Read more
Source§

fn sub_slice_mut<'a>( &mut self, start: [usize; 2], rows: usize, cols: usize, ) -> MatrixSliceMut<'a, T>
where T: 'a,

Produce a MatrixSliceMut from an existing matrix. Read more
Source§

impl<T> Clone for Matrix<T>
where T: Clone,

Source§

fn clone(&self) -> Matrix<T>

Clones the Matrix.

1.0.0 · Source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
Source§

impl CostFunc<Matrix<f64>> for CrossEntropyError

Source§

fn cost(outputs: &Matrix<f64>, targets: &Matrix<f64>) -> f64

The cost function.
Source§

fn grad_cost(outputs: &Matrix<f64>, targets: &Matrix<f64>) -> Matrix<f64>

The gradient of the cost function.
Source§

impl CostFunc<Matrix<f64>> for MeanSqError

Source§

fn cost(outputs: &Matrix<f64>, targets: &Matrix<f64>) -> f64

The cost function.
Source§

fn grad_cost(outputs: &Matrix<f64>, targets: &Matrix<f64>) -> Matrix<f64>

The gradient of the cost function.
Source§

impl<T> Debug for Matrix<T>
where T: Debug,

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result<(), Error>

Formats the value using the given formatter. Read more
Source§

impl<T> Display for Matrix<T>
where T: Display,

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result<(), Error>

Formats the Matrix for display.

Source§

impl<'a, 'b, T> Div<&'b T> for &'a Matrix<T>
where T: Copy + Div<Output = T>,

Scalar division with matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the / operator.
Source§

fn div(self, f: &T) -> Matrix<T>

Performs the / operation. Read more
Source§

impl<'a, T> Div<&'a T> for Matrix<T>
where T: Copy + Div<Output = T>,

Scalar division with matrix.

Will reuse the memory allocated for the existing matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the / operator.
Source§

fn div(self, f: &T) -> Matrix<T>

Performs the / operation. Read more
Source§

impl<'a, T> Div<T> for &'a Matrix<T>
where T: Copy + Div<Output = T>,

Scalar division with matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the / operator.
Source§

fn div(self, f: T) -> Matrix<T>

Performs the / operation. Read more
Source§

impl<T> Div<T> for Matrix<T>
where T: Copy + Div<Output = T>,

Scalar division with matrix.

Will reuse the memory allocated for the existing matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the / operator.
Source§

fn div(self, f: T) -> Matrix<T>

Performs the / operation. Read more
Source§

impl<'a, T> DivAssign<&'a T> for Matrix<T>
where T: Copy + Div<Output = T>,

Performs division assignment between a matrix and a scalar.

Source§

fn div_assign(&mut self, _rhs: &T)

Performs the /= operation. Read more
Source§

impl<T> DivAssign<T> for Matrix<T>
where T: Copy + Div<Output = T>,

Performs division assignment between a matrix and a scalar.

Source§

fn div_assign(&mut self, _rhs: T)

Performs the /= operation. Read more
Source§

impl<'a, T> From<MatrixSlice<'a, T>> for Matrix<T>
where T: Copy,

Source§

fn from(slice: MatrixSlice<'a, T>) -> Matrix<T>

Converts to this type from the input type.
Source§

impl<'a, T> From<MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy,

Source§

fn from(slice: MatrixSliceMut<'a, T>) -> Matrix<T>

Converts to this type from the input type.
Source§

impl<T> From<Vector<T>> for Matrix<T>

Source§

fn from(vector: Vector<T>) -> Matrix<T>

Converts to this type from the input type.
Source§

impl<'a, T> FromIterator<&'a [T]> for Matrix<T>
where T: 'a + Copy,

Creates a Matrix from an iterator over slices.

Each of the slices produced by the iterator will become a row in the matrix.

§Panics

Will panic if the iterators items do not have constant length.

§Examples

We can create a new matrix from some data.

use rulinalg::matrix::{Matrix, BaseMatrix};

let a : Matrix<f64> = vec![4f64; 16].chunks(4).collect();

assert_eq!(a.rows(), 4);
assert_eq!(a.cols(), 4);

We can also do more interesting things.

use rulinalg::matrix::{Matrix, BaseMatrix};

let a = Matrix::new(4,2, (0..8).collect::<Vec<usize>>());

// Here we skip the first row and take only those
// where the first entry is less than 6.
let b = a.iter_rows()
         .skip(1)
         .filter(|x| x[0] < 6)
         .collect::<Matrix<usize>>();

// We take the middle rows
assert_eq!(b.into_vec(), vec![2,3,4,5]);
Source§

fn from_iter<I>(iterable: I) -> Matrix<T>
where I: IntoIterator<Item = &'a [T]>,

Creates a value from an iterator. Read more
Source§

impl<T> Hash for Matrix<T>
where T: Hash,

Source§

fn hash<__H>(&self, state: &mut __H)
where __H: Hasher,

Feeds this value into the given Hasher. Read more
1.3.0 · Source§

fn hash_slice<H>(data: &[Self], state: &mut H)
where H: Hasher, Self: Sized,

Feeds a slice of this type into the given Hasher. Read more
Source§

impl<T> Index<[usize; 2]> for Matrix<T>

Indexes matrix.

Takes row index first then column.

Source§

type Output = T

The returned type after indexing.
Source§

fn index(&self, idx: [usize; 2]) -> &T

Performs the indexing (container[index]) operation. Read more
Source§

impl<T> IndexMut<[usize; 2]> for Matrix<T>

Indexes mutable matrix.

Takes row index first then column.

Source§

fn index_mut(&mut self, idx: [usize; 2]) -> &mut T

Performs the mutable indexing (container[index]) operation. Read more
Source§

impl<T: Float> Invertible<Matrix<T>> for MinMaxScaler<T>

Source§

fn inv_transform(&self, inputs: Matrix<T>) -> Result<Matrix<T>, Error>

Maps the inputs using the inverse of the fitted transform.
Source§

impl<T: Float + FromPrimitive> Invertible<Matrix<T>> for Standardizer<T>

Source§

fn inv_transform(&self, inputs: Matrix<T>) -> Result<Matrix<T>, Error>

Maps the inputs using the inverse of the fitted transform.
Source§

impl<T> Metric<T> for Matrix<T>
where T: Float,

Source§

fn norm(&self) -> T

Compute euclidean norm for matrix.

§Examples
use rulinalg::matrix::Matrix;
use rulinalg::Metric;

let a = Matrix::new(2,1, vec![3.0,4.0]);
let c = a.norm();

assert_eq!(c, 5.0);
Source§

impl<'a, 'b, T> Mul<&'b Matrix<T>> for &'a Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, 'c, T> Mul<&'c Matrix<T>> for &'b MatrixSlice<'a, T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, 'c, T> Mul<&'c Matrix<T>> for &'b MatrixSliceMut<'a, T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<&'a Matrix<T>> for Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<&'b Matrix<T>> for MatrixSlice<'a, T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<&'b Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, 'c, T> Mul<&'c MatrixSlice<'a, T>> for &'b Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &MatrixSlice<'_, T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<&'b MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &MatrixSlice<'_, T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, 'c, T> Mul<&'c MatrixSliceMut<'a, T>> for &'b Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<&'b MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<&'b T> for &'a Matrix<T>
where T: Copy + Mul<Output = T>,

Scalar multiplication with matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, f: &T) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<&'a T> for Matrix<T>
where T: Copy + Mul<Output = T>,

Scalar multiplication with matrix.

Will reuse the memory allocated for the existing matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, f: &T) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<&'b Vector<T>> for &'a Matrix<T>
where T: Zero<Output = T> + Mul<Output = T> + Add + Copy,

Multiplies matrix by vector.

Source§

type Output = Vector<T>

The resulting type after applying the * operator.
Source§

fn mul(self, v: &Vector<T>) -> Vector<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<&'a Vector<T>> for Matrix<T>
where T: Zero<Output = T> + Mul<Output = T> + Add + Copy,

Multiplies matrix by vector.

Source§

type Output = Vector<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: &Vector<T>) -> Vector<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<Matrix<T>> for &'a Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<Matrix<T>> for &'b MatrixSlice<'a, T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<Matrix<T>> for &'b MatrixSliceMut<'a, T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<Matrix<T>> for MatrixSlice<'a, T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<MatrixSlice<'a, T>> for &'b Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: MatrixSlice<'_, T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: MatrixSlice<'_, T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, 'b, T> Mul<MatrixSliceMut<'a, T>> for &'b Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<T> for &'a Matrix<T>
where T: Copy + Mul<Output = T>,

Scalar multiplication with matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, f: T) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<T> Mul<T> for Matrix<T>
where T: Copy + Mul<Output = T>,

Scalar multiplication with matrix.

Will reuse the memory allocated for the existing matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, f: T) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul<Vector<T>> for &'a Matrix<T>
where T: Zero<Output = T> + Mul<Output = T> + Add + Copy,

Multiplies matrix by vector.

Source§

type Output = Vector<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: Vector<T>) -> Vector<T>

Performs the * operation. Read more
Source§

impl<T> Mul<Vector<T>> for Matrix<T>
where T: Zero<Output = T> + Mul<Output = T> + Add + Copy,

Multiplies matrix by vector.

Source§

type Output = Vector<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: Vector<T>) -> Vector<T>

Performs the * operation. Read more
Source§

impl<'a, T> Mul for Matrix<T>
where T: Copy + Zero<Output = T> + Add + Mul<Output = T> + Any,

Multiplies two matrices together.

Source§

type Output = Matrix<T>

The resulting type after applying the * operator.
Source§

fn mul(self, m: Matrix<T>) -> Matrix<T>

Performs the * operation. Read more
Source§

impl<'a, T> MulAssign<&'a T> for Matrix<T>
where T: Copy + Mul<Output = T>,

Performs multiplication assignment between a matrix and a scalar.

Source§

fn mul_assign(&mut self, _rhs: &T)

Performs the *= operation. Read more
Source§

impl<T> MulAssign<T> for Matrix<T>
where T: Copy + Mul<Output = T>,

Performs multiplication assignment between a matrix and a scalar.

Source§

fn mul_assign(&mut self, _rhs: T)

Performs the *= operation. Read more
Source§

impl<'a, T> Neg for &'a Matrix<T>
where T: Neg<Output = T> + Copy,

Gets negative of matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn neg(self) -> Matrix<T>

Performs the unary - operation. Read more
Source§

impl<T> Neg for Matrix<T>
where T: Neg<Output = T> + Copy,

Gets negative of matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn neg(self) -> Matrix<T>

Performs the unary - operation. Read more
Source§

impl<T> PartialEq for Matrix<T>
where T: PartialEq,

Source§

fn eq(&self, other: &Matrix<T>) -> bool

Tests for self and other values to be equal, and is used by ==.
1.0.0 · Source§

fn ne(&self, other: &Rhs) -> bool

Tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
Source§

impl<'a, 'b, T> Sub<&'b Matrix<T>> for &'a Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between two matrices.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: &Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, 'c, T> Sub<&'c Matrix<T>> for &'b MatrixSlice<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: &Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, 'c, T> Sub<&'c Matrix<T>> for &'b MatrixSliceMut<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: &Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> Sub<&'a Matrix<T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between two matrices.

This will reuse allocated memory from self.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: &Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<&'b Matrix<T>> for MatrixSlice<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: &Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<&'b Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: &Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, 'c, T> Sub<&'c MatrixSlice<'a, T>> for &'b Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, s: &MatrixSlice<'_, T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<&'b MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, s: &MatrixSlice<'_, T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, 'c, T> Sub<&'c MatrixSliceMut<'a, T>> for &'b Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, s: &MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<&'b MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, s: &MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<&'b T> for &'a Matrix<T>
where T: Copy + Sub<Output = T>,

Scalar subtraction with matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, f: &T) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> Sub<&'a T> for Matrix<T>
where T: Copy + Sub<Output = T>,

Scalar subtraction with matrix.

Will reuse the memory allocated for the existing matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, f: &T) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> Sub<Matrix<T>> for &'a Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between two matrices.

This will reuse allocated memory from m.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<Matrix<T>> for &'b MatrixSlice<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<Matrix<T>> for &'b MatrixSliceMut<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> Sub<Matrix<T>> for MatrixSlice<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> Sub<Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<MatrixSlice<'a, T>> for &'b Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, s: MatrixSlice<'_, T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> Sub<MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, s: MatrixSlice<'_, T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, 'b, T> Sub<MatrixSliceMut<'a, T>> for &'b Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, s: MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> Sub<MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between Matrix and MatrixSlice.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, s: MatrixSliceMut<'_, T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> Sub<T> for &'a Matrix<T>
where T: Copy + Sub<Output = T>,

Scalar subtraction with matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, f: T) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<T> Sub<T> for Matrix<T>
where T: Copy + Sub<Output = T>,

Scalar subtraction with matrix.

Will reuse the memory allocated for the existing matrix.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, f: T) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<T> Sub for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction between two matrices.

This will reuse allocated memory from self.

Source§

type Output = Matrix<T>

The resulting type after applying the - operator.
Source§

fn sub(self, m: Matrix<T>) -> Matrix<T>

Performs the - operation. Read more
Source§

impl<'a, T> SubAssign<&'a Matrix<T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction assignment between two matrices.

Source§

fn sub_assign(&mut self, _rhs: &Matrix<T>)

Performs the -= operation. Read more
Source§

impl<'a, 'b, T> SubAssign<&'b Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction assignment between two matrices.

Source§

fn sub_assign(&mut self, _rhs: &Matrix<T>)

Performs the -= operation. Read more
Source§

impl<'a, 'b, T> SubAssign<&'b MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction assignment between two matrices.

Source§

fn sub_assign(&mut self, _rhs: &MatrixSlice<'_, T>)

Performs the -= operation. Read more
Source§

impl<'a, 'b, T> SubAssign<&'b MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction assignment between two matrices.

Source§

fn sub_assign(&mut self, _rhs: &MatrixSliceMut<'_, T>)

Performs the -= operation. Read more
Source§

impl<'a, T> SubAssign<&'a T> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs subtraction assignment between a matrix and a scalar.

Source§

fn sub_assign(&mut self, _rhs: &T)

Performs the -= operation. Read more
Source§

impl<'a, T> SubAssign<Matrix<T>> for MatrixSliceMut<'a, T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction assignment between two matrices.

Source§

fn sub_assign(&mut self, _rhs: Matrix<T>)

Performs the -= operation. Read more
Source§

impl<'a, T> SubAssign<MatrixSlice<'a, T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction assignment between two matrices.

Source§

fn sub_assign(&mut self, _rhs: MatrixSlice<'_, T>)

Performs the -= operation. Read more
Source§

impl<'a, T> SubAssign<MatrixSliceMut<'a, T>> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction assignment between two matrices.

Source§

fn sub_assign(&mut self, _rhs: MatrixSliceMut<'_, T>)

Performs the -= operation. Read more
Source§

impl<T> SubAssign<T> for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs subtraction assignment between a matrix and a scalar.

Source§

fn sub_assign(&mut self, _rhs: T)

Performs the -= operation. Read more
Source§

impl<T> SubAssign for Matrix<T>
where T: Copy + Sub<Output = T>,

Performs elementwise subtraction assignment between two matrices.

Source§

fn sub_assign(&mut self, _rhs: Matrix<T>)

Performs the -= operation. Read more
Source§

impl<T: Distribution> SupModel<Matrix<f64>, Matrix<f64>> for NaiveBayes<T>

Train and predict from the Naive Bayes model.

The input matrix must be rows made up of features. The target matrix should have indicator vectors in each row specifying the input class. e.g. [[1,0,0],[0,0,1]] shows class 1 first, then class 3.

Source§

fn train( &mut self, inputs: &Matrix<f64>, targets: &Matrix<f64>, ) -> LearningResult<()>

Train the model using inputs and targets.

Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Matrix<f64>>

Predict output from inputs.

Source§

impl<'a, T, A> SupModel<Matrix<f64>, Matrix<f64>> for NeuralNet<'a, T, A>
where T: Criterion, A: OptimAlgorithm<BaseNeuralNet<'a, T>>,

Supervised learning for the Neural Network.

The model is trained using back propagation.

Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Matrix<f64>>

Predict neural network output using forward propagation.

Source§

fn train( &mut self, inputs: &Matrix<f64>, targets: &Matrix<f64>, ) -> LearningResult<()>

Train the model using gradient optimization and back propagation.

Source§

impl<T: Kernel, U: MeanFunc> SupModel<Matrix<f64>, Vector<f64>> for GaussianProcess<T, U>

Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Vector<f64>>

Predict output from inputs.

Source§

fn train( &mut self, inputs: &Matrix<f64>, targets: &Vector<f64>, ) -> LearningResult<()>

Train the model using data and outputs.

Source§

impl<C: Criterion> SupModel<Matrix<f64>, Vector<f64>> for GenLinearModel<C>

Supervised model trait for the GLM.

Predictions are made from the model by computing g^-1(Xb).

The model is trained using Iteratively Re-weighted Least Squares.

Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Vector<f64>>

Predict output from inputs.

Source§

fn train( &mut self, inputs: &Matrix<f64>, targets: &Vector<f64>, ) -> LearningResult<()>

Train the model using inputs and targets.

Source§

impl SupModel<Matrix<f64>, Vector<f64>> for LinRegressor

Source§

fn train( &mut self, inputs: &Matrix<f64>, targets: &Vector<f64>, ) -> LearningResult<()>

Train the linear regression model.

Takes training data and output values as input.

§Examples
use rusty_machine::learning::lin_reg::LinRegressor;
use rusty_machine::linalg::Matrix;
use rusty_machine::linalg::Vector;
use rusty_machine::learning::SupModel;

let mut lin_mod = LinRegressor::default();
let inputs = Matrix::new(3,1, vec![2.0, 3.0, 4.0]);
let targets = Vector::new(vec![5.0, 6.0, 7.0]);

lin_mod.train(&inputs, &targets).unwrap();
Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Vector<f64>>

Predict output value from input data.

Model must be trained before prediction can be made.

Source§

impl<A> SupModel<Matrix<f64>, Vector<f64>> for LogisticRegressor<A>

Source§

fn train( &mut self, inputs: &Matrix<f64>, targets: &Vector<f64>, ) -> LearningResult<()>

Train the logistic regression model.

Takes training data and output values as input.

§Examples
use rusty_machine::learning::logistic_reg::LogisticRegressor;
use rusty_machine::linalg::Matrix;
use rusty_machine::linalg::Vector;
use rusty_machine::learning::SupModel;

let mut logistic_mod = LogisticRegressor::default();
let inputs = Matrix::new(3,2, vec![1.0, 2.0, 1.0, 3.0, 1.0, 4.0]);
let targets = Vector::new(vec![5.0, 6.0, 7.0]);

logistic_mod.train(&inputs, &targets).unwrap();
Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Vector<f64>>

Predict output value from input data.

Model must be trained before prediction can be made.

Source§

impl<K: Kernel> SupModel<Matrix<f64>, Vector<f64>> for SVM<K>

Train the model using the Pegasos algorithm and predict the model output from new data.

Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Vector<f64>>

Predict output from inputs.
Source§

fn train( &mut self, inputs: &Matrix<f64>, targets: &Vector<f64>, ) -> LearningResult<()>

Train the model using inputs and targets.
Source§

impl<T: Float> Transformer<Matrix<T>> for MinMaxScaler<T>

Source§

fn fit(&mut self, inputs: &Matrix<T>) -> Result<(), Error>

Fit Transformer to input data, and stores the transformation in the Transformer
Source§

fn transform(&mut self, inputs: Matrix<T>) -> Result<Matrix<T>, Error>

Transforms the inputs and stores the transformation in the Transformer
Source§

impl<R: Rng, T> Transformer<Matrix<T>> for Shuffler<R>

The Shuffler will transform the input Matrix by shuffling its rows in place.

Under the hood this uses a Fisher-Yates shuffle.

Source§

fn fit(&mut self, inputs: &Matrix<T>) -> Result<(), Error>

Fit Transformer to input data, and stores the transformation in the Transformer
Source§

fn transform(&mut self, inputs: Matrix<T>) -> LearningResult<Matrix<T>>

Transforms the inputs and stores the transformation in the Transformer
Source§

impl<T: Float + FromPrimitive> Transformer<Matrix<T>> for Standardizer<T>

Source§

fn fit(&mut self, inputs: &Matrix<T>) -> Result<(), Error>

Fit Transformer to input data, and stores the transformation in the Transformer
Source§

fn transform(&mut self, inputs: Matrix<T>) -> Result<Matrix<T>, Error>

Transforms the inputs and stores the transformation in the Transformer
Source§

impl UnSupModel<Matrix<f64>, Matrix<f64>> for GaussianMixtureModel

Source§

fn train(&mut self, inputs: &Matrix<f64>) -> LearningResult<()>

Train the model using inputs.

Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Matrix<f64>>

Predict output from inputs.

Source§

impl UnSupModel<Matrix<f64>, Vector<Option<usize>>> for DBSCAN

Source§

fn train(&mut self, inputs: &Matrix<f64>) -> LearningResult<()>

Train the classifier using input data.

Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Vector<Option<usize>>>

Predict output from inputs.
Source§

impl<InitAlg: Initializer> UnSupModel<Matrix<f64>, Vector<usize>> for KMeansClassifier<InitAlg>

Source§

fn predict(&self, inputs: &Matrix<f64>) -> LearningResult<Vector<usize>>

Predict classes from data.

Model must be trained.

Source§

fn train(&mut self, inputs: &Matrix<f64>) -> LearningResult<()>

Train the classifier using input data.

Source§

impl<T> Eq for Matrix<T>
where T: Eq,

Source§

impl<T> StructuralPartialEq for Matrix<T>

Auto Trait Implementations§

§

impl<T> Freeze for Matrix<T>

§

impl<T> RefUnwindSafe for Matrix<T>
where T: RefUnwindSafe,

§

impl<T> Send for Matrix<T>
where T: Send,

§

impl<T> Sync for Matrix<T>
where T: Sync,

§

impl<T> Unpin for Matrix<T>
where T: Unpin,

§

impl<T> UnwindSafe for Matrix<T>
where T: UnwindSafe,

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dest: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dest. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> ToOwned for T
where T: Clone,

Source§

type Owned = T

The resulting type after obtaining ownership.
Source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
Source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
Source§

impl<T> ToString for T
where T: Display + ?Sized,

Source§

fn to_string(&self) -> String

Converts the given value to a String. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.