[][src]Struct mathru::algebra::linear::matrix::matrix::Matrix

pub struct Matrix<T> { /* fields omitted */ }

Implementations

impl<T> Matrix<T>[src]

pub fn iter(&self) -> MatrixIterator<'_, T>

Notable traits for MatrixIterator<'a, T>

impl<'a, T> Iterator for MatrixIterator<'a, T> where
    T: Field + Scalar
type Item = &'a T;
[src]

pub fn iter_mut(&mut self) -> MatrixIteratorMut<'_, T>

Notable traits for MatrixIteratorMut<'a, T>

impl<'a, T> Iterator for MatrixIteratorMut<'a, T> where
    T: Field + Scalar
type Item = &'a mut T;
[src]

pub fn row_into_iter(&self) -> MatrixRowIntoIterator<'_, T>

Notable traits for MatrixRowIntoIterator<'a, T>

impl<'a, T> Iterator for MatrixRowIntoIterator<'a, T> where
    T: Field + Scalar
type Item = Vector<T>;
[src]

pub fn column_into_iter(&self) -> MatrixColumnIntoIterator<'_, T>

Notable traits for MatrixColumnIntoIterator<'a, T>

impl<'a, T> Iterator for MatrixColumnIntoIterator<'a, T> where
    T: Field + Scalar
type Item = Vector<T>;
[src]

impl<T> Matrix<T> where
    T: Clone
[src]

pub fn apply_mut(self: Matrix<T>, f: &dyn Fn(&T) -> T) -> Matrix<T>[src]

Applies the function f on every element in the matrix

pub fn apply(self: &Matrix<T>, f: &dyn Fn(&T) -> T) -> Matrix<T>[src]

impl<T> Matrix<T> where
    T: Field + Scalar + Power
[src]

pub fn det(&self) -> T[src]

Calculates the determinant

Example

use mathru::algebra::linear::Matrix;

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

pub fn trace(&self) -> T where
    T: Field + Scalar
[src]

Calculates the trace of a matrix

Arguments

self: square matrix

Panics

if it is not a square matrix

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, -2.0, 3.0, -7.0]);
let tr: f64 = a.trace();

assert_eq!(-6.0, tr);

impl<T> Matrix<T> where
    T: Field + Scalar
[src]

pub fn get_column(&self, i: usize) -> Vector<T>[src]

pub fn get_row(&self, i: usize) -> Vector<T>[src]

return row vector

i: row

pub fn set_column(&mut self, column: &Vector<T>, i: usize)[src]

set column

pub fn set_row(&mut self, row: &Vector<T>, i: usize)[src]

set row

Arguments

  • 'row'
  • 'i'

Panics

impl<T> Matrix<T> where
    T: Field + Scalar + Power
[src]

pub fn givens<'d, 'e>(m: usize, i: usize, j: usize, c: T, s: T) -> Self[src]

pub fn givens_cosine_sine_pair(a: T, b: T) -> (T, T)[src]

function [c,s] = Givens(a,b) Givens rotation computation Determines cosine-sine pair (c,s) so that [c s;-s c]'*[a;b] = [r;0] GVL4: Algorithm 5.1.3

impl<T> Matrix<T> where
    T: Field + Scalar + Power
[src]

pub fn householder(v: &Vector<T>, k: usize) -> Self[src]

Returns the householder matrix

Arguments

v: Column vector k: index 0 <= k < m

Panics

if index out of bounds

impl<T> Matrix<T> where
    T: Field + Scalar + Power
[src]

pub fn dec_sv(&self) -> (Self, Self, Self)[src]

Computes the singular value decomposition

M = U * S * V*

Return

(u, s, v)

Example

use mathru::algebra::linear::Matrix;

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

let (u, s, v): (Matrix<f64>, Matrix<f64>, Matrix<f64>) = a.dec_sv();

pub fn rot(f: T, g: T) -> (T, T, T)[src]

pub fn householder_bidiag(&self) -> (Self, Self, Self)[src]

self is an m times n matrix with m >= n A = UBV^{T} U \in T^{m \times n} B \in T^{n \times n} V \in T^{n \times n}

impl<T> Matrix<T> where
    T: Field + Scalar
[src]

pub fn get_slice(
    &self,
    row_s: usize,
    row_e: usize,
    column_s: usize,
    column_e: usize
) -> Matrix<T>
[src]

Returns a slice of the matrix

Arugments

0 <= row_s < m
0 <= row_e < m
0 <= column_s < n
0 <= column_e <= n

row_s: start row
row_e: end row
column_s: start column
column_e: end column

Example

use mathru::algebra::linear::{Matrix};

let mut a: Matrix<f64> = matrix![1.0, -2.0; 3.0, -7.0];
a = a.get_slice(0, 0, 1, 1);

let a_ref: Matrix<f64> = Matrix::new(1, 1, vec![-2.0]);

assert_eq!(a_ref, a);

pub fn set_slice(self, slice: &Self, row: usize, column: usize) -> Matrix<T>[src]

Replaces parts of the matrix with the given values

Arugments

0 <= row < m
0 <= column < n

Example

use mathru::algebra::linear::{Matrix};

let mut a: Matrix<f64> = matrix![   1.0, 0.0;
                                    3.0, -7.0];

let b: Matrix<f64> = matrix![2.0, -1.0];
a = a.set_slice(&b, 0, 0);

let a_updated: Matrix<f64> = matrix![   2.0, -1.0;
                                        3.0, -7.0];

assert_eq!(a_updated, a);

impl<T> Matrix<T>[src]

pub fn get_mut(&mut self, i: usize, j: usize) -> &mut T[src]

Returns the mutual element a_ij from the matrix

Example

use mathru::algebra::linear::{Matrix};

let mut a: Matrix<f64> = matrix![1.0, 0.0; 3.0, -7.0];
*a.get_mut(1, 0) = -8.0;

let a_updated: Matrix<f64> = matrix![1.0, 0.0; -8.0, -7.0];
assert_eq!(a_updated, a);

pub fn get(&self, i: usize, j: usize) -> &T[src]

Returns the element a_ij from the matrix

Example

use mathru::algebra::linear::{Matrix};

let a: Matrix<f64> = matrix![   1.0, 0.0;
                                3.0, -7.0];

let a_ref: f64 = 3.0;
let element: f64 = *a.get(1, 0);

assert_eq!(a_ref, element);

impl<T> Matrix<T> where
    T: Clone + Copy
[src]

pub fn new(m: usize, n: usize, data: Vec<T>) -> Self[src]

Creates a new Matrix object

Fortran like, column wise

[ 0, 1, 2] 3, 4, 5, 6, 7, 8 ] => vec![ 0, 3, 6, 1, 4, 7, 2, 5, 8]

impl<T> Matrix<T>[src]

pub fn convert_to_vec(self) -> Vec<T>[src]

impl<T> Matrix<T> where
    T: Scalar + Clone + Copy
[src]

pub fn new_random(m: usize, n: usize) -> Matrix<T>[src]

impl<T> Matrix<T> where
    T: Field + Scalar
[src]

pub fn zero(m: usize, n: usize) -> Self[src]

Returns the zero matrix(additive neutral element)

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = &a + &Matrix::zero(2, 2);

impl<T> Matrix<T> where
    T: Field + Scalar
[src]

pub fn one(size: usize) -> Self[src]

Returns the eye matrix(multiplicative neutral element)

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = &a * &Matrix::one(2);

pub fn ones(m: usize, n: usize) -> Self[src]

impl<T> Matrix<T> where
    T: Field + Scalar + Power
[src]

pub fn pinv(&self) -> Matrix<T>[src]

Calculates the pseudo inverse matrix

A^+ = (A^TA)^-1A^T

impl<T> Matrix<T>[src]

pub fn dim(&self) -> (usize, usize)[src]

Returns the matrix dimension

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(4, 2, vec![1.0, 0.0, 3.0, 0.0, 1.0, -7.0, 0.5, 0.25]);
let (m, n) = a.dim();

assert_eq!(4, m);
assert_eq!(2, n);

pub fn nrows(&self) -> usize[src]

Returns the number of rows

Example

use mathru::algebra::linear::{Matrix};

let a: Matrix<f64> = Matrix::new(4, 2, vec![1.0, 0.0, 3.0, 0.0, 1.0, -7.0, 0.5, 0.25]);
let m: usize = a.nrows();

assert_eq!(4, m);

pub fn ncols(&self) -> usize[src]

Returns the number of columns

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(4, 2, vec![1.0, 0.0, 3.0, 0.0, 1.0, -7.0, 0.5, 0.25]);
let n: usize = a.ncols();

assert_eq!(2, n);

impl<T> Matrix<T> where
    T: Real
[src]

pub fn dec_eigen(self) -> EigenDec<T>[src]

Computes the eigenvalues of a real matrix

Arguments

Example

use mathru::algebra::linear::{matrix::EigenDec, Matrix, Vector};

let a: Matrix<f64> = Matrix::new(3, 3, vec![1.0, -3.0, 3.0, 3.0, -5.0, 3.0, 6.0, -6.0, 4.0]);
let eigen: EigenDec<f64> = a.dec_eigen();

pub fn eigenvalue_r(&self) -> Vector<T>[src]

pub fn eigenvector_r(&self, value: &Vector<T>) -> Matrix<T>[src]

impl<T> Matrix<T> where
    T: Field + Scalar + Power
[src]

pub fn dec_hessenberg(&self) -> HessenbergDec<T>[src]

Decomposes self in to the M

q * h * q^T = self

Arguments

Return

(q, h)

Panics

if M is not a square matrix

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(3, 3, vec![1.0, 5.0, 3.0, 1.0, 0.0, -7.0, 3.0, 8.0, 9.0]);
let (q, h): (Matrix<f64>, Matrix<f64>) = a.dec_hessenberg().qh();

impl<T> Matrix<T> where
    T: Field + Scalar
[src]

pub fn dec_lu<'a>(&'a self) -> Result<LUDec<T>, ()>[src]

Decomposes the matrix into a upper and a lower matrix

PA = LU

Arguments

Example

use mathru::algebra::linear::Matrix;

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

let (l, u, p): (Matrix<f64>, Matrix<f64>, Matrix<f64>) = a.dec_lu().unwrap().lup();

impl<T> Matrix<T> where
    T: Field + Scalar + Power
[src]

pub fn dec_qr<'a>(&'a self) -> QRDec<T>[src]

QR Decomposition with Givens rotations

A = QR
Q is an orthogonal matrix
R is an upper triangular matrix

Panics

if A is not a square matrix

Example

use mathru::algebra::linear::Matrix;

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

let (q, r): (Matrix<f64>, Matrix<f64>) = a.dec_qr().qr();

impl<T> Matrix<T> where
    T: Field + Scalar
[src]

pub fn inv_r(&self) -> Result<Matrix<T>, ()>[src]

impl<T> Matrix<T> where
    T: Field + Scalar + Power
[src]

pub fn dec_cholesky(&self) -> Result<CholeskyDec<T>, ()>[src]

Decomposes the symetric, positive definite quadractic matrix A into a lower triangular matrix L A = L L^T

Arguments

A has to be symetric and postive definite

Panics

If A is not a quadratic matrix

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = matrix![   2.0, -1.0, 0.0;
                               -1.0, 2.0, -1.0;
                                0.0, -1.0,  2.0];

let l: (Matrix<f64>) = a.dec_cholesky().unwrap().l();

Trait Implementations

impl AbsDiffEq<Matrix<f32>> for Matrix<f32>[src]

type Epsilon = f32

Used for specifying relative comparisons.

impl AbsDiffEq<Matrix<f64>> for Matrix<f64>[src]

type Epsilon = f64

Used for specifying relative comparisons.

impl<'a, 'b, T> Add<&'b Matrix<T>> for &'a Matrix<T> where
    T: Field + Scalar
[src]

Adds two matrices

type Output = Matrix<T>

The resulting type after applying the + operator.

pub fn add(self, rhs: &'b Matrix<T>) -> Self::Output[src]

Adds two matrices

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);

let c: Matrix<f64> = &b + &a;

impl<'a, 'b, T> Add<&'b T> for &'a Matrix<T> where
    T: Field + Scalar
[src]

Add scalar to matrix

type Output = Matrix<T>

The resulting type after applying the + operator.

pub fn add(self, rhs: &T) -> Self::Output[src]

Add a scalar to the matrix

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = &a + &-4.0;

impl<T> Add<Matrix<T>> for Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the + operator.

pub fn add(self, rhs: Self) -> Self::Output[src]

Adds two matrices

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);

let c: Matrix<f64> = a + b;

impl<T> Add<T> for Matrix<T> where
    T: Field + Scalar
[src]

Add scalar to matrix

type Output = Matrix<T>

The resulting type after applying the + operator.

pub fn add(self, rhs: T) -> Self::Output[src]

Add a scalar to the matrix

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = a + -4.0;

impl<T: Clone> Clone for Matrix<T>[src]

impl<T: Debug> Debug for Matrix<T>[src]

impl<'de, T> Deserialize<'de> for Matrix<T> where
    T: Deserialize<'de>, 
[src]

impl<T> Display for Matrix<T> where
    T: Display
[src]

impl<'a, 'b, T> Div<&'b T> for &'a Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the / operator.

pub fn div(self, m: &'b T) -> Matrix<T>[src]

Divide all matrix with a scalar

Example

use mathru::algebra::linear::Matrix;

let res_ref: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let a: Matrix<f64> = Matrix::new(2, 2, vec![4.0, 0.0, 12.0, -28.0]);

impl<T> Div<T> for Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the / operator.

pub fn div(self, m: T) -> Matrix<T>[src]

Divides all matrix element with a scalar

Example

use mathru::algebra::linear::Matrix;

let res_ref: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let f: f64 = 7.0;
let a: Matrix<f64> = Matrix::new(2, 2, vec![7.0, 0.0, 21.0, -49.0]);

impl<T> From<Vector<T>> for Matrix<T> where
    T: Field + Scalar
[src]

impl<T> Identity<Addition> for Matrix<T> where
    T: Identity<Addition>, 
[src]

pub fn id() -> Self[src]

Returns the additive neutral element)

Example

use mathru::algebra::linear::Matrix;

impl<T> IntoIterator for Matrix<T> where
    T: Field + Scalar
[src]

type Item = T

The type of the elements being iterated over.

type IntoIter = MatrixIntoIterator<T>

Which kind of iterator are we turning this into?

impl<T> Inverse<T> for Matrix<T> where
    T: Field + Scalar
[src]

pub fn inv(&self) -> Result<Matrix<T>, ()>[src]

Inverse Matrix

Example

use mathru::algebra::linear::{matrix::*, Matrix};

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b_inv: Matrix<f64> = a.inv().unwrap();

impl<'a, 'b, T> Mul<&'b Matrix<T>> for &'a Vector<T> where
    T: Field + Scalar
[src]

type Output = Vector<T>

The resulting type after applying the * operator.

impl<'a, 'b, T> Mul<&'b Matrix<T>> for &'a Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the * operator.

pub fn mul(self, rhs: &'b Matrix<T>) -> Self::Output[src]

Multiplies two matrices

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let res_ref: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, -18.0, 49.0]);
assert_eq!(res_ref, &a * &b);

impl<'a, 'b, T> Mul<&'b T> for &'a Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the * operator.

pub fn mul(self, m: &'b T) -> Matrix<T>[src]

Multiplies a matrix with a scalar

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let res_ref: Matrix<f64> = Matrix::new(2, 2, vec![4.0, 0.0, 12.0, -28.0]);

assert_eq!(res_ref, &a * &4.0);

impl<'a, 'b, T> Mul<&'b Vector<T>> for &'a Matrix<T> where
    T: Field + Scalar
[src]

Multiplies matrix by vector.

type Output = Vector<T>

The resulting type after applying the * operator.

impl<T> Mul<Matrix<T>> for Vector<T> where
    T: Field + Scalar
[src]

type Output = Vector<T>

The resulting type after applying the * operator.

impl<T> Mul<Matrix<T>> for Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the * operator.

pub fn mul(self, rhs: Self) -> Self::Output[src]

Multiplies two matrices

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let res_ref: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, -18.0, 49.0]);
assert_eq!(res_ref, a * b);

impl<T> Mul<T> for Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the * operator.

pub fn mul(self, m: T) -> Matrix<T>[src]

Multiplies a matrix with a scalar

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let f: f64 = 7.0;
let res_ref: Matrix<f64> = Matrix::new(2, 2, vec![7.0, 0.0, 21.0, -49.0]);

assert_eq!(res_ref, a * f);

impl<T> Mul<Vector<T>> for Matrix<T> where
    T: Field + Scalar
[src]

type Output = Vector<T>

The resulting type after applying the * operator.

impl<T> PartialEq<Matrix<T>> for Matrix<T> where
    T: PartialEq
[src]

pub fn eq(&self, other: &Self) -> bool[src]

Checks if two matrices are equal

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);

assert_eq!(true, a == b);

impl RelativeEq<Matrix<f32>> for Matrix<f32>[src]

pub fn relative_eq(
    &self,
    other: &Matrix<f32>,
    epsilon: Self::Epsilon,
    max_relative: Self::Epsilon
) -> bool
[src]

A test for equality that uses a relative comparison if the values are far apart.

impl RelativeEq<Matrix<f64>> for Matrix<f64>[src]

pub fn relative_eq(
    &self,
    other: &Matrix<f64>,
    epsilon: Self::Epsilon,
    max_relative: Self::Epsilon
) -> bool
[src]

A test for equality that uses a relative comparison if the values are far apart.

impl<T> Serialize for Matrix<T> where
    T: Serialize
[src]

impl<T> Solve<Matrix<T>> for LUDec<T> where
    T: Field + Scalar
[src]

impl<T> Solve<Matrix<T>> for Matrix<T> where
    T: Field + Scalar
[src]

impl<T> Solve<Vector<T>> for Matrix<T> where
    T: Field + Scalar
[src]

pub fn solve(&self, rhs: &Vector<T>) -> Result<Vector<T>, ()>[src]

Solves Ax = y where A \in R^{m * n}, x \in R^n, y \in R^m

impl<'a, 'b, T> Sub<&'b Matrix<T>> for &'a Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the - operator.

impl<'a, 'b, T> Sub<&'b T> for &'a Matrix<T> where
    T: Field + Scalar
[src]

Subtracts scalar from all matrix elements

type Output = Matrix<T>

The resulting type after applying the - operator.

pub fn sub(self, rhs: &T) -> Self::Output[src]

Subtracts a scalar value from all matrix elements

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = &a - &-4.0;

impl<T> Sub<Matrix<T>> for Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the - operator.

pub fn sub(self, rhs: Self) -> Self::Output[src]

Subtracts two matrices

A = (a_{ij}) \in T^{m \times n} B = (b_{ij}) \in T^{m \times n} A - B = ( a_{ij} - b_{ij} )

Arguments

rhs:

Example

use mathru::algebra::linear::Matrix;

let a: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);
let b: Matrix<f64> = Matrix::new(2, 2, vec![1.0, 0.0, 3.0, -7.0]);

let c: Matrix<f64> = a - b;

impl<T> Sub<T> for Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

The resulting type after applying the - operator.

pub fn sub(self, rhs: T) -> Self::Output[src]

Subtracts a scalar from all matrix elements

Example

use mathru::algebra::linear::Matrix;

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

impl<T> Substitute<Matrix<T>> for Matrix<T> where
    T: Field + Scalar
[src]

impl<T> Substitute<Vector<T>> for Matrix<T> where
    T: Field + Scalar
[src]

impl<T> Transpose for Matrix<T> where
    T: Field + Scalar
[src]

type Output = Matrix<T>

pub fn transpose(self) -> Matrix<T>[src]

Function to transpose a matrix without allocating memory for the transposed matrix

catanzaro.name/papers/PPoPP-2014.pdf TODO

Example

use mathru::algebra::linear::Matrix;
use mathru::algebra::linear::matrix::Transpose;

let mut uut: Matrix<f64> = Matrix::new(4, 2, vec![1.0, 0.0, 3.0, 0.0, 1.0, -7.0, 0.5, 0.25]);
uut = uut.transpose();

let refer: Matrix<f64> = Matrix::new(2, 4, vec![1.0, 1.0, 0.0, -7.0, 3.0, 0.5, 0.0, 0.25]);
assert_eq!(refer, uut);

Auto Trait Implementations

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

impl<T> Any for T where
    T: 'static + ?Sized
[src]

impl<T> Borrow<T> for T where
    T: ?Sized
[src]

impl<T> BorrowMut<T> for T where
    T: ?Sized
[src]

impl<T> DeserializeOwned for T where
    T: for<'de> Deserialize<'de>, 
[src]

impl<T> From<T> for T[src]

impl<T, U> Into<U> for T where
    U: From<T>, 
[src]

impl<T> ToOwned for T where
    T: Clone
[src]

type Owned = T

The resulting type after obtaining ownership.

impl<T> ToString for T where
    T: Display + ?Sized
[src]

impl<T, U> TryFrom<U> for T where
    U: Into<T>, 
[src]

type Error = Infallible

The type returned in the event of a conversion error.

impl<T, U> TryInto<U> for T where
    U: TryFrom<T>, 
[src]

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

The type returned in the event of a conversion error.

impl<V, T> VZip<V> for T where
    V: MultiLane<T>,