ndarray_linalg/lib.rs
1//! The `ndarray-linalg` crate provides linear algebra functionalities for `ArrayBase`, the n-dimensional array data structure provided by [`ndarray`](https://github.com/rust-ndarray/ndarray).
2//!
3//! `ndarray-linalg` leverages [LAPACK](http://www.netlib.org/lapack/)'s routines using the bindings provided by [blas-lapack-rs/lapack](https://github.com/blas-lapack-rs/lapack).
4//!
5//! Linear algebra methods
6//! -----------------------
7//! - Decomposition methods:
8//! - [QR decomposition](qr/index.html)
9//! - [Cholesky/LU decomposition](cholesky/index.html)
10//! - [Eigenvalue decomposition](eig/index.html)
11//! - [Eigenvalue decomposition for Hermite matrices](eigh/index.html)
12//! - [**S**ingular **V**alue **D**ecomposition](svd/index.html)
13//! - Solution of linear systems:
14//! - [General matrices](solve/index.html)
15//! - [Triangular matrices](triangular/index.html)
16//! - [Hermitian/real symmetric matrices](solveh/index.html)
17//! - [Tridiagonal matrices](tridiagonal/index.html)
18//! - [Inverse matrix computation](solve/trait.Inverse.html)
19//!
20//! Naming Convention
21//! -----------------------
22//! Each routine is usually exposed as a trait, implemented by the relevant types.
23//!
24//! For each routine there might be multiple "variants": different traits corresponding to the different ownership possibilities of the array you intend to work on.
25//!
26//! For example, if you are interested in the QR decomposition of a square matrix, you can use:
27//! - [QRSquare](qr/trait.QRSquare.html), if you hold an immutable reference (i.e. `&self`) to the matrix you want to decompose;
28//! - [QRSquareInplace](qr/trait.QRSquareInplace.html), if you hold a mutable reference (i.e. `&mut self`) to the matrix you want to decompose;
29//! - [QRSquareInto](qr/trait.QRSquareInto.html), if you can pass the matrix you want to decompose by value (e.g. `self`).
30//!
31//! Depending on the algorithm, each variant might require more or less copy operations of the underlying data.
32//!
33//! Details are provided in the description of each routine.
34//!
35//! Utilities
36//! -----------
37//! - [Assertions for array](index.html#macros)
38//! - [Random matrix generators](generate/index.html)
39//! - [Scalar trait](types/trait.Scalar.html)
40
41#![allow(
42 clippy::module_inception,
43 clippy::many_single_char_names,
44 clippy::type_complexity,
45 clippy::ptr_arg
46)]
47
48#[macro_use]
49extern crate ndarray;
50
51pub mod assert;
52pub mod cholesky;
53pub mod convert;
54pub mod diagonal;
55pub mod eig;
56pub mod eigh;
57pub mod error;
58pub mod generate;
59pub mod inner;
60pub mod krylov;
61pub mod layout;
62pub mod least_squares;
63pub mod lobpcg;
64pub mod norm;
65pub mod operator;
66pub mod opnorm;
67pub mod qr;
68pub mod solve;
69pub mod solveh;
70pub mod svd;
71pub mod svddc;
72pub mod trace;
73pub mod triangular;
74pub mod tridiagonal;
75pub mod types;
76
77pub use crate::assert::*;
78pub use crate::cholesky::*;
79pub use crate::convert::*;
80pub use crate::diagonal::*;
81pub use crate::eig::*;
82pub use crate::eigh::*;
83pub use crate::generate::*;
84pub use crate::inner::*;
85pub use crate::layout::*;
86pub use crate::least_squares::*;
87pub use crate::lobpcg::{TruncatedEig, TruncatedOrder, TruncatedSvd};
88pub use crate::norm::*;
89pub use crate::operator::*;
90pub use crate::opnorm::*;
91pub use crate::qr::*;
92pub use crate::solve::*;
93pub use crate::solveh::*;
94pub use crate::svd::*;
95pub use crate::svddc::*;
96pub use crate::trace::*;
97pub use crate::triangular::*;
98pub use crate::tridiagonal::*;
99pub use crate::types::*;