1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202
/*!
sprs is a sparse linear algebra library for Rust.
It features a sparse matrix type, [**`CsMat`**](struct.CsMatBase.html), and a sparse vector type,
[**`CsVec`**](struct.CsVecBase.html), both based on the
[compressed storage scheme](https://en.wikipedia.org/wiki/Sparse_matrix#Compressed_sparse_row_.28CSR.2C_CRS_or_Yale_format.29).
## Features
- sparse matrix/sparse matrix addition, multiplication.
- sparse vector/sparse vector addition, dot product.
- sparse matrix/dense matrix addition, multiplication.
- sparse triangular solves.
- powerful iteration over the sparse structure, enabling easy extension of the library.
- matrix construction using the [triplet format](struct.TriMatBase.html),
vertical and horizontal stacking, block construction.
- sparse cholesky solver in the separate crate `sprs-ldl`.
- fully generic integer type for the storage of indices, enabling compact
representations.
- planned interoperability with existing sparse solvers such as `SuiteSparse`.
## Quick Examples
Matrix construction:
```rust
use sprs::{CsMat, TriMat};
let mut a = TriMat::new((4, 4));
a.add_triplet(0, 0, 3.0_f64);
a.add_triplet(1, 2, 2.0);
a.add_triplet(3, 0, -2.0);
// This matrix type does not allow computations, and must to
// converted to a compatible sparse type, using for example
let b: CsMat<_> = a.to_csr();
```
Constructing matrix using the more efficient direct sparse constructor
```rust
use sprs::{CsMat, CsVec};
let eye : CsMat<f64> = CsMat::eye(3);
let a = CsMat::new_csc((3, 3),
vec![0, 2, 4, 5],
vec![0, 1, 0, 2, 2],
vec![1., 2., 3., 4., 5.]);
```
Matrix vector multiplication:
```rust
use sprs::{CsMat, CsVec};
let eye = CsMat::eye(5);
let x = CsVec::new(5, vec![0, 2, 4], vec![1., 2., 3.]);
let y = &eye * &x;
assert_eq!(x, y);
```
Matrix matrix multiplication, addition:
```rust
use sprs::{CsMat, CsVec};
let eye = CsMat::eye(3);
let a = CsMat::new_csc((3, 3),
vec![0, 2, 4, 5],
vec![0, 1, 0, 2, 2],
vec![1., 2., 3., 4., 5.]);
let b = &eye * &a;
assert_eq!(a, b.to_csc());
```
*/
#![allow(clippy::redundant_slicing)]
pub mod array_backend;
mod dense_vector;
pub mod errors;
pub mod indexing;
pub mod io;
mod mul_acc;
pub mod num_kinds;
pub mod num_matrixmarket;
mod range;
mod sparse;
pub mod stack;
pub type Ix1 = ndarray::Ix1;
pub type Ix2 = ndarray::Ix2;
pub use crate::indexing::SpIndex;
pub use crate::sparse::{
csmat::CsIter,
indptr::{IndPtr, IndPtrBase, IndPtrView},
kronecker::kronecker_product,
CsMat, CsMatBase, CsMatI, CsMatVecView, CsMatView, CsMatViewI,
CsMatViewMut, CsMatViewMutI, CsStructure, CsStructureI, CsStructureView,
CsStructureViewI, CsVec, CsVecBase, CsVecI, CsVecView, CsVecViewI,
CsVecViewMut, CsVecViewMutI, SparseMat, TriMat, TriMatBase, TriMatI,
TriMatIter, TriMatView, TriMatViewI, TriMatViewMut, TriMatViewMutI,
};
pub use crate::dense_vector::{DenseVector, DenseVectorMut};
pub use crate::mul_acc::MulAcc;
pub use crate::sparse::symmetric::is_symmetric;
pub use crate::sparse::permutation::{
perm_is_valid, transform_mat_papt, PermOwned, PermOwnedI, PermView,
PermViewI, Permutation,
};
pub use crate::sparse::CompressedStorage::{self, CSC, CSR};
pub use crate::sparse::binop;
pub use crate::sparse::linalg;
pub use crate::sparse::prod;
pub use crate::sparse::smmp;
pub use crate::sparse::special_mats;
pub use crate::sparse::visu;
pub mod vec {
pub use crate::sparse::{CsVec, CsVecBase, CsVecView, CsVecViewMut};
pub use crate::sparse::vec::{
IntoSparseVecIter, NnzEither, NnzIndex, NnzOrZip, SparseIterTools,
VecDim, VectorIterator, VectorIteratorMut,
};
}
pub use crate::sparse::construct::{bmat, hstack, vstack};
pub use crate::sparse::to_dense::assign_to_dense;
/// The shape of a matrix. This a 2-tuple with the first element indicating
/// the number of rows, and the second element indicating the number of
/// columns.
pub type Shape = (usize, usize); // FIXME: maybe we could use Ix2 here?
/// Configuration enum to ask for symmetry checks in algorithms
#[derive(Copy, Clone, Eq, PartialEq, Debug)]
pub enum SymmetryCheck {
CheckSymmetry,
DontCheckSymmetry,
}
pub use SymmetryCheck::*;
/// Configuration enum to ask for permutation checks in algorithms
#[derive(Copy, Clone, Eq, PartialEq, Debug)]
pub enum PermutationCheck {
CheckPerm,
DontCheckPerm,
}
pub use PermutationCheck::*;
/// The different kinds of fill-in-reduction algorithms supported by sprs
#[derive(Copy, Clone, Eq, PartialEq, Debug)]
#[non_exhaustive]
pub enum FillInReduction {
NoReduction,
ReverseCuthillMcKee,
#[allow(clippy::upper_case_acronyms)]
CAMDSuiteSparse,
}
#[cfg(feature = "approx")]
/// Traits for comparing vectors and matrices using the approx traits
///
/// Comparisons of sparse matrices with different storages might be slow.
/// It is advised to compare using the same storage order for efficiency
///
/// These traits requires the `approx` feature to be activated
pub mod approx {
pub use approx::{AbsDiffEq, RelativeEq, UlpsEq};
}
#[cfg(test)]
mod test_data;
#[cfg(test)]
mod test {
use super::CsMat;
#[test]
fn iter_rbr() {
let mat = CsMat::new(
(3, 3),
vec![0, 2, 3, 3],
vec![1, 2, 0],
vec![0.1, 0.2, 0.3],
);
let view = mat.view();
let mut iter = view.iter();
assert_eq!(iter.next(), Some((&0.1, (0, 1))));
assert_eq!(iter.next(), Some((&0.2, (0, 2))));
assert_eq!(iter.next(), Some((&0.3, (1, 0))));
assert_eq!(iter.next(), None);
}
}