[][src]Struct sprs::CsMatBase

pub struct CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr = I> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
{ /* fields omitted */ }

Compressed matrix in the CSR or CSC format, with sorted indices.

This sparse matrix format is the preferred format for performing arithmetic operations. Constructing a sparse matrix directly in this format requires a deep knowledge of its internals. For easier matrix construction, the triplet format is preferred.

The CsMatBase type is parameterized by the scalar type N, the indexing type I, the indexing storage backend types IptrStorage and IndStorage, and the value storage backend type DataStorage. Convenient aliases are available to specify frequent variants: CsMat refers to a sparse matrix that owns its data, similar to Vec<T>; CsMatView refers to a sparse matrix that borrows its data, similar to & [T]; and CsMatViewMut refers to a sparse matrix borrowing its data, with a mutable borrow for its values. No mutable borrow is allowed for the structure of the matrix, allowing the invariants to be preserved.

Additionaly, the type aliases CsMatI, CsMatViewI and CsMatViewMutI can be used to choose an index type different from the default usize.

Storage format

In the compressed storage format, the non-zero values of a sparse matrix are stored as the row and column location of the non-zero values, with a compression along the rows (CSR) or columns (CSC) indices. The dimension along which the storage is compressed is referred to as the outer dimension, the other dimension is called the inner dimension. For clarity, the remaining explanation will assume a CSR matrix, but the information stands for CSC matrices as well.

Indptr

An index pointer array indptr of size corresponding to the number of rows stores the cumulative sum of non-zero elements for each row. For instance, the number of non-zero elements of the i-th row can be obtained by computing indptr[i + 1] - indptr[i]. The total number of non-zero elements is thus nnz = indptr[nb_rows + 1]. This index pointer array can then be used to efficiently index the indices and data array, which respectively contain the column indices and the values of the non-zero elements.

Indices and data

The non-zero locations and values are stored in arrays of size nnz, indices and data. For row i, the non-zeros are located in the slices indices[indptr[i]..indptr[i+1]] and data[indptr[i]..indptr[i+1]]. We require and enforce sorted indices for each row.

Construction

A sparse matrix can be directly constructed by providing its index pointer, indices and data arrays. The coherence of the provided structure is then verified.

For situations where the compressed structure is hard to figure out up front, the triplet format can be used. A matrix in the triplet format can then be efficiently converted to a CsMat.

Alternately, a sparse matrix can be constructed from other sparse matrices using vstack, hstack or bmat.

Implementations

impl<N, I: SpIndex, Iptr: SpIndex> CsMatBase<N, I, Vec<Iptr, Global>, Vec<I, Global>, Vec<N, Global>, Iptr>[src]

pub fn eye(dim: usize) -> Self where
    N: Num + Clone
[src]

Identity matrix, stored as a CSR matrix.

use sprs::{CsMat, CsVec};
let eye = CsMat::eye(5);
assert!(eye.is_csr());
let x = CsVec::new(5, vec![0, 2, 4], vec![1., 2., 3.]);
let y = &eye * &x;
assert_eq!(x, y);

pub fn eye_csc(dim: usize) -> CsMatI<N, I, Iptr> where
    N: Num + Clone
[src]

Identity matrix, stored as a CSC matrix.

use sprs::{CsMat, CsVec};
let eye = CsMat::eye_csc(5);
assert!(eye.is_csc());
let x = CsVec::new(5, vec![0, 2, 4], vec![1., 2., 3.]);
let y = &eye * &x;
assert_eq!(x, y);

pub fn empty(
    storage: CompressedStorage,
    inner_size: usize
) -> CsMatI<N, I, Iptr>
[src]

Create an empty CsMat for building purposes

pub fn zero(shape: Shape) -> CsMatI<N, I, Iptr>[src]

Create a new CsMat representing the zero matrix. Hence it has no non-zero elements.

pub fn reserve_outer_dim(&mut self, outer_dim_additional: usize)[src]

Reserve the storage for the given additional number of nonzero data

pub fn reserve_nnz(&mut self, nnz_additional: usize)[src]

Reserve the storage for the given additional number of nonzero data

pub fn reserve_outer_dim_exact(&mut self, outer_dim_lim: usize)[src]

Reserve the storage for the given number of nonzero data

pub fn reserve_nnz_exact(&mut self, nnz_lim: usize)[src]

Reserve the storage for the given number of nonzero data

pub fn try_new(
    shape: Shape,
    indptr: Vec<Iptr>,
    indices: Vec<I>,
    data: Vec<N>
) -> Result<CsMatI<N, I, Iptr>, SprsError> where
    N: Copy
[src]

Try create an owned CSR matrix from moved data.

An owned CSC matrix can be created with try_new_csc().

If necessary, the indices will be sorted in place.

pub fn new(
    shape: Shape,
    indptr: Vec<Iptr>,
    indices: Vec<I>,
    data: Vec<N>
) -> CsMatI<N, I, Iptr> where
    N: Copy
[src]

Create an owned CSR matrix from moved data.

An owned CSC matrix can be created with new_csc().

If necessary, the indices will be sorted in place.

Panics

  • if indptr does not correspond to the number of rows.
  • if indices and data don't have exactly indptr[rows] elements.
  • if indices contains values greater or equal to the number of columns.

pub fn try_new_csc(
    shape: Shape,
    indptr: Vec<Iptr>,
    indices: Vec<I>,
    data: Vec<N>
) -> Result<Self, SprsError> where
    N: Copy
[src]

Try create an owned CSC matrix from moved data.

An owned CSC matrix can be created with new_csc().

If necessary, the indices will be sorted in place.

pub fn new_csc(
    shape: Shape,
    indptr: Vec<Iptr>,
    indices: Vec<I>,
    data: Vec<N>
) -> Self where
    N: Copy
[src]

Create an owned CSC matrix from moved data.

An owned CSC matrix can be created with new_csc().

If necessary, the indices will be sorted in place.

Panics

  • if indptr does not correspond to the number of rows.
  • if indices and data don't have exactly indptr[rows] elements.
  • if indices contains values greater or equal to the number of columns.

pub fn csr_from_dense(
    m: ArrayView<'_, N, Ix2>,
    epsilon: N
) -> CsMatI<N, I, Iptr> where
    N: Num + Clone + PartialOrd + Signed
[src]

Create a CSR matrix from a dense matrix, ignoring elements lower than epsilon.

If epsilon is negative, it will be clamped to zero.

pub fn csc_from_dense(
    m: ArrayView<'_, N, Ix2>,
    epsilon: N
) -> CsMatI<N, I, Iptr> where
    N: Num + Clone + PartialOrd + Signed
[src]

Create a CSC matrix from a dense matrix, ignoring elements lower than epsilon.

If epsilon is negative, it will be clamped to zero.

pub fn append_outer(self, data: &[N]) -> Self where
    N: Clone + Num
[src]

Append an outer dim to an existing matrix, compressing it in the process

pub fn append_outer_csvec(self, vec: CsVecBase<&[I], &[N]>) -> Self where
    N: Clone
[src]

Append an outer dim to an existing matrix, provided by a sparse vector

pub fn insert(&mut self, row: usize, col: usize, val: N)[src]

Insert an element in the matrix. If the element is already present, its value is overwritten.

Warning: this is not an efficient operation, as it requires a non-constant lookup followed by two Vec insertions.

The insertion will be efficient, however, if the elements are inserted according to the matrix's order, eg following the row order for a CSR matrix.

impl<'a, N: 'a, I: 'a + SpIndex, Iptr: 'a + SpIndex> CsMatBase<N, I, &'a [Iptr], &'a [I], &'a [N], Iptr>[src]

Constructor methods for sparse matrix views

These constructors can be used to create views over non-matrix data such as slices.

pub fn new_view(
    storage: CompressedStorage,
    shape: Shape,
    indptr: &'a [Iptr],
    indices: &'a [I],
    data: &'a [N]
) -> Result<CsMatViewI<'a, N, I, Iptr>, SprsError>
[src]

Create a borrowed CsMat matrix from sliced data, checking their validity

pub unsafe fn new_view_raw(
    storage: CompressedStorage,
    shape: Shape,
    indptr: *const Iptr,
    indices: *const I,
    data: *const N
) -> CsMatViewI<'a, N, I, Iptr>
[src]

Create a borrowed CsMat matrix from raw data, without checking their validity

Safety

This is unsafe because algorithms are free to assume that properties guaranteed by check_compressed_structure are enforced. For instance, non out-of-bounds indices can be relied upon to perform unchecked slice access.

pub fn middle_outer_views(
    &self,
    i: usize,
    count: usize
) -> CsMatViewI<'a, N, I, Iptr>
[src]

Get a view into count contiguous outer dimensions, starting from i.

eg this gets the rows from i to i + count in a CSR matrix

pub fn iter_rbr(&self) -> CsIter<'a, N, I, Iptr>

Notable traits for CsIter<'a, N, I, Iptr>

impl<'a, N, I, Iptr> Iterator for CsIter<'a, N, I, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    N: 'a, 
type Item = (&'a N, (I, I));
[src]

Get an iterator that yields the non-zero locations and values stored in this matrix, in the fastest iteration order.

This method will yield the correct lifetime for iterating over a sparse matrix view.

impl<N, I, Iptr, IptrStorage, IndStorage, DataStorage> CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

pub fn storage(&self) -> CompressedStorage[src]

The underlying storage of this matrix

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

The number of rows of this matrix

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

The number of cols of this matrix

pub fn shape(&self) -> Shape[src]

The shape of the matrix. Equivalent to let shape = (a.rows(), a.cols()).

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

The number of non-zero elements this matrix stores. This is often relevant for the complexity of most sparse matrix algorithms, which are often linear in the number of non-zeros.

pub fn density(&self) -> f64[src]

The density of the sparse matrix, defined as the number of non-zero elements divided by the maximum number of elements

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

Number of outer dimensions, that ie equal to self.rows() for a CSR matrix, and equal to self.cols() for a CSC matrix

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

Number of inner dimensions, that ie equal to self.cols() for a CSR matrix, and equal to self.rows() for a CSC matrix

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

Access the element located at row i and column j. Will return None if there is no non-zero element at this location.

This access is logarithmic in the number of non-zeros in the corresponding outer slice. It is therefore advisable not to rely on this for algorithms, and prefer outer_iterator() which accesses elements in storage order.

pub fn indptr(&self) -> &[Iptr][src]

The array of offsets in the indices() and data() slices. The elements of the slice at outer dimension i are available between the elements indptr[i] and indptr[i+1] in the indices() and data() slices.

Example

use sprs::{CsMat};
let eye : CsMat<f64> = CsMat::eye(5);
// get the element of row 3
// there is only one element in this row, with a column index of 3
// and a value of 1.
let loc = eye.indptr()[3];
assert_eq!(eye.indptr()[4], loc + 1);
assert_eq!(loc, 3);
assert_eq!(eye.indices()[loc], 3);
assert_eq!(eye.data()[loc], 1.);

pub fn indices(&self) -> &[I][src]

The inner dimension location for each non-zero value. See the documentation of indptr() for more explanations.

pub fn data(&self) -> &[N][src]

The non-zero values. See the documentation of indptr() for more explanations.

pub fn into_raw_storage(self) -> (IptrStorage, IndStorage, DataStorage)[src]

Destruct the matrix object and recycle its storage containers.

Example

use sprs::{CsMat};
let (indptr, indices, data) = CsMat::<i32>::eye(3).into_raw_storage();
assert_eq!(indptr, vec![0, 1, 2, 3]);
assert_eq!(indices, vec![0, 1, 2]);
assert_eq!(data, vec![1, 1, 1]);

pub fn is_csc(&self) -> bool[src]

Test whether the matrix is in CSC storage

pub fn is_csr(&self) -> bool[src]

Test whether the matrix is in CSR storage

pub fn transpose_mut(&mut self)[src]

Transpose a matrix in place No allocation required (this is simply a storage order change)

pub fn transpose_into(self) -> Self[src]

Transpose a matrix in place No allocation required (this is simply a storage order change)

pub fn transpose_view(&self) -> CsMatViewI<'_, N, I, Iptr>[src]

Transposed view of this matrix No allocation required (this is simply a storage order change)

pub fn to_owned(&self) -> CsMatI<N, I, Iptr> where
    N: Clone
[src]

Get an owned version of this matrix. If the matrix was already owned, this will make a deep copy.

pub fn to_inner_onehot(&self) -> CsMatI<N, I, Iptr> where
    N: Copy + Clone + Float + PartialOrd
[src]

Generate a one-hot matrix, compressing the inner dimension.

Returns a matrix with the same size, the same CSR/CSC type, and a single value of 1.0 within each populated inner vector.

See CsMatBase::into_csc and CsMatBase::into_csr if you need to prepare a matrix for one-hot compression.

pub fn to_other_types<I2, N2, Iptr2>(&self) -> CsMatI<N2, I2, Iptr2> where
    N: Clone + Into<N2>,
    I2: SpIndex,
    Iptr2: SpIndex
[src]

Clone the matrix with another integer type for indptr and indices

Panics

If the indices or indptr values cannot be represented by the requested integer type.

pub fn view(&self) -> CsMatViewI<'_, N, I, Iptr>[src]

Return a view into the current matrix

pub fn structure_view(&self) -> CsStructureViewI<'_, I, Iptr>[src]

pub fn to_dense(&self) -> Array<N, Ix2> where
    N: Clone + Zero
[src]

pub fn outer_iterator(&self) -> OuterIterator<'_, N, I, Iptr>

Notable traits for OuterIterator<'iter, N, I, Iptr>

impl<'iter, N: 'iter, I: 'iter + SpIndex, Iptr: 'iter + SpIndex> Iterator for OuterIterator<'iter, N, I, Iptr> type Item = CsVecBase<&'iter [I], &'iter [N]>;
[src]

Return an outer iterator for the matrix

This can be used for iterating over the rows (resp. cols) of a CSR (resp. CSC) matrix.

use sprs::{CsMat};
let eye = CsMat::eye(5);
for (row_ind, row_vec) in eye.outer_iterator().enumerate() {
    let (col_ind, &val): (_, &f64) = row_vec.iter().next().unwrap();
    assert_eq!(row_ind, col_ind);
    assert_eq!(val, 1.);
}

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

Get the max number of nnz for each outer dim

pub fn degrees(&self) -> Vec<usize>[src]

Get the degrees of each vertex on a symmetric matrix

The nonzero pattern of a symmetric matrix can be interpreted as an undirected graph. In such a graph, a vertex i is connected to another vertex j if there is a corresponding nonzero entry in the matrix at location (i, j).

This function returns a vector containing the degree of each vertex, that is to say the number of neighbor of each vertex. We do not count diagonal entries as a neighbor.

pub fn outer_view(&self, i: usize) -> Option<CsVecViewI<'_, N, I>>[src]

Get a view into the i-th outer dimension (eg i-th row for a CSR matrix)

pub fn diag(&self) -> CsVecI<N, I> where
    N: Clone
[src]

Get the diagonal of a sparse matrix

pub fn outer_block_iter(
    &self,
    block_size: usize
) -> ChunkOuterBlocks<'_, N, I, Iptr>
[src]

Iteration on outer blocks of size block_size

pub fn map<F, N2>(&self, f: F) -> CsMatI<N2, I, Iptr> where
    F: FnMut(&N) -> N2, 
[src]

Return a new sparse matrix with the same sparsity pattern, with all non-zero values mapped by the function f.

pub fn get_outer_inner(&self, outer_ind: usize, inner_ind: usize) -> Option<&N>[src]

Access an element given its outer_ind and inner_ind. Will return None if there is no non-zero element at this location.

This access is logarithmic in the number of non-zeros in the corresponding outer slice. It is therefore advisable not to rely on this for algorithms, and prefer outer_iterator() which accesses elements in storage order.

pub fn nnz_index(&self, row: usize, col: usize) -> Option<NnzIndex>[src]

Find the non-zero index of the element specified by row and col

Searching this index is logarithmic in the number of non-zeros in the corresponding outer slice. Once it is available, the NnzIndex enables retrieving the data with O(1) complexity.

pub fn nnz_index_outer_inner(
    &self,
    outer_ind: usize,
    inner_ind: usize
) -> Option<NnzIndex>
[src]

Find the non-zero index of the element specified by outer_ind and inner_ind.

Searching this index is logarithmic in the number of non-zeros in the corresponding outer slice.

pub fn check_compressed_structure(&self) -> Result<(), SprsError>[src]

Check the structure of CsMat components This will ensure that:

  • indptr is of length outer_dim() + 1
  • indices and data have the same length, nnz == indptr[outer_dims()]
  • indptr is sorted
  • indptr values do not exceed usize::MAX / 2, as that would mean indices and indptr would take more space than the addressable memory
  • indices is sorted for each outer slice
  • indices are lower than inner_dims()

pub fn iter(&self) -> CsIter<'_, N, I, Iptr>

Notable traits for CsIter<'a, N, I, Iptr>

impl<'a, N, I, Iptr> Iterator for CsIter<'a, N, I, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    N: 'a, 
type Item = (&'a N, (I, I));
[src]

Get an iterator that yields the non-zero locations and values stored in this matrix, in the fastest iteration order.

impl<N, I, Iptr, IptrStorage, IndStorage, DataStorage> CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    N: Default,
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

pub fn to_other_storage(&self) -> CsMatI<N, I, Iptr> where
    N: Clone
[src]

Create a matrix mathematically equal to this one, but with the opposed storage (a CSC matrix will be converted to CSR, and vice versa)

pub fn to_csc(&self) -> CsMatI<N, I, Iptr> where
    N: Clone
[src]

Create a new CSC matrix equivalent to this one. A new matrix will be created even if this matrix was already CSC.

pub fn to_csr(&self) -> CsMatI<N, I, Iptr> where
    N: Clone
[src]

Create a new CSR matrix equivalent to this one. A new matrix will be created even if this matrix was already CSR.

impl<N, I, Iptr> CsMatBase<N, I, Vec<Iptr, Global>, Vec<I, Global>, Vec<N, Global>, Iptr> where
    N: Default,
    I: SpIndex,
    Iptr: SpIndex
[src]

pub fn into_csc(self) -> CsMatI<N, I, Iptr> where
    N: Clone
[src]

Create a new CSC matrix equivalent to this one. If this matrix is CSR, it is converted to CSC If this matrix is CSC, it is returned by value

pub fn into_csr(self) -> CsMatI<N, I, Iptr> where
    N: Clone
[src]

Create a new CSR matrix equivalent to this one. If this matrix is CSC, it is converted to CSR If this matrix is CSR, it is returned by value

impl<N, I, Iptr, IptrStorage, IndStorage, DataStorage> CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: DerefMut<Target = [N]>, 
[src]

pub fn data_mut(&mut self) -> &mut [N][src]

Mutable access to the non zero values

This enables changing the values without changing the matrix's structure. To also change the matrix's structure, see modify

pub fn scale(&mut self, val: N) where
    N: Num + Copy
[src]

Sparse matrix self-multiplication by a scalar

pub fn outer_view_mut(&mut self, i: usize) -> Option<CsVecViewMutI<'_, N, I>>[src]

Get a mutable view into the i-th outer dimension (eg i-th row for a CSR matrix)

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

Get a mutable reference to the element located at row i and column j. Will return None if there is no non-zero element at this location.

This access is logarithmic in the number of non-zeros in the corresponding outer slice. It is therefore advisable not to rely on this for algorithms, and prefer outer_iterator_mut() which accesses elements in storage order. TODO: outer_iterator_mut is not yet implemented

pub fn get_outer_inner_mut(
    &mut self,
    outer_ind: usize,
    inner_ind: usize
) -> Option<&mut N>
[src]

Get a mutable reference to an element given its outer_ind and inner_ind. Will return None if there is no non-zero element at this location.

This access is logarithmic in the number of non-zeros in the corresponding outer slice. It is therefore advisable not to rely on this for algorithms, and prefer outer_iterator_mut() which accesses elements in storage order.

pub fn set(&mut self, row: usize, col: usize, val: N)[src]

Set the value of the non-zero element located at (row, col)

Panics

  • on out-of-bounds access
  • if no non-zero element exists at the given location

pub fn map_inplace<F>(&mut self, f: F) where
    F: FnMut(&N) -> N, 
[src]

Apply a function to every non-zero element

pub fn outer_iterator_mut(&mut self) -> OuterIteratorMut<'_, N, I, Iptr>

Notable traits for OuterIteratorMut<'iter, N, I, Iptr>

impl<'iter, N: 'iter, I: 'iter + SpIndex, Iptr: 'iter + SpIndex> Iterator for OuterIteratorMut<'iter, N, I, Iptr> type Item = CsVecViewMutI<'iter, N, I>;
[src]

Return a mutable outer iterator for the matrix

This iterator yields mutable sparse vector views for each outer dimension. Only the non-zero values can be modified, the structure is kept immutable.

impl<N, I, Iptr, IptrStorage, IndStorage, DataStorage> CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: DerefMut<Target = [Iptr]>,
    IndStorage: DerefMut<Target = [I]>,
    DataStorage: DerefMut<Target = [N]>, 
[src]

pub fn modify<F>(&mut self, f: F) where
    F: FnMut(&mut [Iptr], &mut [I], &mut [N]), 
[src]

Modify the matrix's structure without changing its nonzero count.

The coherence of the structure will be checked afterwards.

Panics

If the resulting matrix breaks the CsMat invariants (sorted indices, no out of bounds indices).

Example

use sprs::CsMat;
// |   1   |
// | 1     |
// |   1 1 |
let mut mat = CsMat::new_csc((3, 3),
                                  vec![0, 1, 3, 4],
                                  vec![1, 0, 2, 2],
                                  vec![1.; 4]);

// | 1 2   |
// | 1     |
// |   1   |
mat.modify(|indptr, indices, data| {
    indptr[1] = 2;
    indptr[2] = 4;
    indices[0] = 0;
    indices[1] = 1;
    indices[2] = 0;
    data[2] = 2.;
});

impl<'a, N: 'a, I: 'a + SpIndex, Iptr: 'a + SpIndex> CsMatBase<N, I, Vec<Iptr>, &'a [I], &'a [N], Iptr>[src]

pub unsafe fn new_vecview_raw(
    storage: CompressedStorage,
    nrows: usize,
    ncols: usize,
    indptr: *const Iptr,
    indices: *const I,
    data: *const N
) -> CsMatBase<N, I, Array2<Iptr>, &'a [I], &'a [N], Iptr>
[src]

Create a borrowed row or column CsMat matrix from raw data, without checking their validity

Safety

This is unsafe because algorithms are free to assume that properties guaranteed by check_compressed_structure are enforced. For instance, non out-of-bounds indices can be relied upon to perform unchecked slice access.

Trait Implementations

impl<N, I, Iptr, IS1, DS1, ISptr1, IS2, ISptr2, DS2> AbsDiffEq<CsMatBase<N, I, ISptr2, IS2, DS2, Iptr>> for CsMatBase<N, I, ISptr1, IS1, DS1, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    CsMatBase<N, I, ISptr1, IS1, DS1, Iptr>: PartialEq<CsMatBase<N, I, ISptr2, IS2, DS2, Iptr>>,
    IS1: Deref<Target = [I]>,
    IS2: Deref<Target = [I]>,
    ISptr1: Deref<Target = [Iptr]>,
    ISptr2: Deref<Target = [Iptr]>,
    DS1: Deref<Target = [N]>,
    DS2: Deref<Target = [N]>,
    N: AbsDiffEq,
    N::Epsilon: Clone,
    N: Zero
[src]

type Epsilon = N::Epsilon

Used for specifying relative comparisons.

impl<'a, 'b, N, I, Iptr, IpS, IS, DS, DS2> Add<&'b ArrayBase<DS2, Dim<[usize; 2]>>> for &'a CsMatBase<N, I, IpS, IS, DS, Iptr> where
    N: 'a + Copy + Num + Default,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpS: 'a + Deref<Target = [Iptr]>,
    IS: 'a + Deref<Target = [I]>,
    DS: 'a + Deref<Target = [N]>,
    DS2: 'b + Data<Elem = N>, 
[src]

type Output = Array<N, Ix2>

The resulting type after applying the + operator.

impl<'a, 'b, N, I, Iptr, IpStorage, IStorage, DStorage, IpS2, IS2, DS2> Add<&'b CsMatBase<N, I, IpS2, IS2, DS2, Iptr>> for &'a CsMatBase<N, I, IpStorage, IStorage, DStorage, Iptr> where
    N: 'a + Copy + Num + Default,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [N]>,
    IpS2: 'a + Deref<Target = [Iptr]>,
    IS2: 'a + Deref<Target = [I]>,
    DS2: 'a + Deref<Target = [N]>, 
[src]

type Output = CsMatI<N, I, Iptr>

The resulting type after applying the + operator.

impl<N: Clone, I: Clone, IptrStorage: Clone, IndStorage: Clone, DataStorage: Clone, Iptr: Clone> Clone for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

impl<N: Copy, I: Copy, IptrStorage: Copy, IndStorage: Copy, DataStorage: Copy, Iptr: Copy> Copy for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

impl<N: Debug, I: Debug, IptrStorage: Debug, IndStorage: Debug, DataStorage: Debug, Iptr: Debug> Debug for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

impl<'de, N, I, IptrStorage, IndStorage, DataStorage, Iptr> Deserialize<'de> for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>,
    IptrStorage: Deserialize<'de>,
    IndStorage: Deserialize<'de>,
    DataStorage: Deserialize<'de>, 
[src]

impl<'a, I, Iptr, IpStorage, IStorage, DStorage, T> DivAssign<T> for CsMatBase<T, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + DerefMut<Target = [T]>,
    T: DivAssign<T> + Clone
[src]

impl<'a, 'b, N, I, Iptr, IpS, IS, DS, DS2> Dot<ArrayBase<DS2, Dim<[usize; 1]>>> for CsMatBase<N, I, IpS, IS, DS, Iptr> where
    N: 'a + Copy + Num + Default,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpS: 'a + Deref<Target = [Iptr]>,
    IS: 'a + Deref<Target = [I]>,
    DS: 'a + Deref<Target = [N]>,
    DS2: 'b + Data<Elem = N>, 
[src]

type Output = Array<N, Ix1>

The result of the operation. Read more

impl<'a, 'b, N, I, Iptr, IpS, IS, DS, DS2> Dot<ArrayBase<DS2, Dim<[usize; 2]>>> for CsMatBase<N, I, IpS, IS, DS, Iptr> where
    N: 'a + Copy + Num + Default,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpS: 'a + Deref<Target = [Iptr]>,
    IS: 'a + Deref<Target = [I]>,
    DS: 'a + Deref<Target = [N]>,
    DS2: 'b + Data<Elem = N>, 
[src]

type Output = Array<N, Ix2>

The result of the operation. Read more

impl<'a, 'b, N, I, IpS, IS, DS, DS2> Dot<CsMatBase<N, I, IpS, IS, DS, I>> for ArrayBase<DS2, Ix2> where
    N: 'a + Copy + Num + Default + Debug,
    I: 'a + SpIndex,
    IpS: 'a + Deref<Target = [I]>,
    IS: 'a + Deref<Target = [I]>,
    DS: 'a + Deref<Target = [N]>,
    DS2: 'b + Data<Elem = N>, 
[src]

type Output = Array<N, Ix2>

The result of the operation. Read more

impl<N: Eq, I: Eq, IptrStorage: Eq, IndStorage: Eq, DataStorage: Eq, Iptr: Eq> Eq for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

impl<N: Hash, I: Hash, IptrStorage: Hash, IndStorage: Hash, DataStorage: Hash, Iptr: Hash> Hash for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

impl<N, I, Iptr, IpS, IS, DS> Index<[usize; 2]> for CsMatBase<N, I, IpS, IS, DS, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IpS: Deref<Target = [Iptr]>,
    IS: Deref<Target = [I]>,
    DS: Deref<Target = [N]>, 
[src]

type Output = N

The returned type after indexing.

impl<N, I, Iptr, IpS, IS, DS> IndexMut<[usize; 2]> for CsMatBase<N, I, IpS, IS, DS, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IpS: Deref<Target = [Iptr]>,
    IS: Deref<Target = [I]>,
    DS: DerefMut<Target = [N]>, 
[src]

impl<'a, N, I, IpS, IS, DS, Iptr> IntoIterator for &'a CsMatBase<N, I, IpS, IS, DS, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    N: 'a,
    IpS: Deref<Target = [Iptr]>,
    IS: Deref<Target = [I]>,
    DS: Deref<Target = [N]>, 
[src]

type Item = (&'a N, (I, I))

The type of the elements being iterated over.

type IntoIter = CsIter<'a, N, I, Iptr>

Which kind of iterator are we turning this into?

impl<'a, 'b, N, I, Iptr, IpS, IS, DS, DS2> Mul<&'b ArrayBase<DS2, Dim<[usize; 1]>>> for &'a CsMatBase<N, I, IpS, IS, DS, Iptr> where
    N: 'a + Copy + Num + Default,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpS: 'a + Deref<Target = [Iptr]>,
    IS: 'a + Deref<Target = [I]>,
    DS: 'a + Deref<Target = [N]>,
    DS2: 'b + Data<Elem = N>, 
[src]

type Output = Array<N, Ix1>

The resulting type after applying the * operator.

impl<'a, 'b, N, I, Iptr, IpS, IS, DS, DS2> Mul<&'b ArrayBase<DS2, Dim<[usize; 2]>>> for &'a CsMatBase<N, I, IpS, IS, DS, Iptr> where
    N: 'a + Copy + Num + Default,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpS: 'a + Deref<Target = [Iptr]>,
    IS: 'a + Deref<Target = [I]>,
    DS: 'a + Deref<Target = [N]>,
    DS2: 'b + Data<Elem = N>, 
[src]

type Output = Array<N, Ix2>

The resulting type after applying the * operator.

impl<'a, 'b, N, I, Iptr, IpS1, IS1, DS1, IpS2, IS2, DS2> Mul<&'b CsMatBase<N, I, IpS2, IS2, DS2, Iptr>> for &'a CsMatBase<N, I, IpS1, IS1, DS1, Iptr> where
    N: 'a + Copy + Num + Default + AddAssign + Send + Sync,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpS1: 'a + Deref<Target = [Iptr]>,
    IS1: 'a + Deref<Target = [I]>,
    DS1: 'a + Deref<Target = [N]>,
    IpS2: 'b + Deref<Target = [Iptr]>,
    IS2: 'b + Deref<Target = [I]>,
    DS2: 'b + Deref<Target = [N]>, 
[src]

type Output = CsMatI<N, I, Iptr>

The resulting type after applying the * operator.

impl<'a, 'b, N, I, Iptr, IS1, DS1, IpS2, IS2, DS2> Mul<&'b CsMatBase<N, I, IpS2, IS2, DS2, Iptr>> for &'a CsVecBase<IS1, DS1> where
    N: 'a + Copy + Num + Default + AddAssign + Send + Sync,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IS1: 'a + Deref<Target = [I]>,
    DS1: 'a + Deref<Target = [N]>,
    IpS2: 'b + Deref<Target = [Iptr]>,
    IS2: 'b + Deref<Target = [I]>,
    DS2: 'b + Deref<Target = [N]>, 
[src]

type Output = CsVecI<N, I>

The resulting type after applying the * operator.

impl<'a, 'b, N, I, Iptr, IpS1, IS1, DS1, IS2, DS2> Mul<&'b CsVecBase<IS2, DS2>> for &'a CsMatBase<N, I, IpS1, IS1, DS1, Iptr> where
    N: Copy + Num + Default + Sum + AddAssign + Send + Sync,
    I: SpIndex,
    Iptr: SpIndex,
    IpS1: Deref<Target = [Iptr]>,
    IS1: Deref<Target = [I]>,
    DS1: Deref<Target = [N]>,
    IS2: Deref<Target = [I]>,
    DS2: Deref<Target = [N]>, 
[src]

type Output = CsVecI<N, I>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<f32> for &'a CsMatBase<f32, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [f32]>, 
[src]

type Output = CsMatI<f32, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<f64> for &'a CsMatBase<f64, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [f64]>, 
[src]

type Output = CsMatI<f64, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<i16> for &'a CsMatBase<i16, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [i16]>, 
[src]

type Output = CsMatI<i16, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<i32> for &'a CsMatBase<i32, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [i32]>, 
[src]

type Output = CsMatI<i32, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<i64> for &'a CsMatBase<i64, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [i64]>, 
[src]

type Output = CsMatI<i64, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<i8> for &'a CsMatBase<i8, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [i8]>, 
[src]

type Output = CsMatI<i8, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<isize> for &'a CsMatBase<isize, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [isize]>, 
[src]

type Output = CsMatI<isize, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<u16> for &'a CsMatBase<u16, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [u16]>, 
[src]

type Output = CsMatI<u16, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<u32> for &'a CsMatBase<u32, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [u32]>, 
[src]

type Output = CsMatI<u32, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<u64> for &'a CsMatBase<u64, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [u64]>, 
[src]

type Output = CsMatI<u64, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<u8> for &'a CsMatBase<u8, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [u8]>, 
[src]

type Output = CsMatI<u8, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage> Mul<usize> for &'a CsMatBase<usize, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [usize]>, 
[src]

type Output = CsMatI<usize, I, Iptr>

The resulting type after applying the * operator.

impl<'a, I, Iptr, IpStorage, IStorage, DStorage, T> MulAssign<T> for CsMatBase<T, I, IpStorage, IStorage, DStorage, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + DerefMut<Target = [T]>,
    T: MulAssign<T> + Clone
[src]

impl<N: PartialEq, I: PartialEq, IptrStorage: PartialEq, IndStorage: PartialEq, DataStorage: PartialEq, Iptr: PartialEq> PartialEq<CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr>> for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

impl<N, I, Iptr, IS1, DS1, ISptr1, IS2, ISptr2, DS2> RelativeEq<CsMatBase<N, I, ISptr2, IS2, DS2, Iptr>> for CsMatBase<N, I, ISptr1, IS1, DS1, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    CsMatBase<N, I, ISptr1, IS1, DS1, Iptr>: PartialEq<CsMatBase<N, I, ISptr2, IS2, DS2, Iptr>>,
    IS1: Deref<Target = [I]>,
    IS2: Deref<Target = [I]>,
    ISptr1: Deref<Target = [Iptr]>,
    ISptr2: Deref<Target = [Iptr]>,
    DS1: Deref<Target = [N]>,
    DS2: Deref<Target = [N]>,
    N: RelativeEq,
    N::Epsilon: Clone,
    N: Zero
[src]

impl<N, I, IptrStorage, IndStorage, DataStorage, Iptr> Serialize for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>,
    IptrStorage: Serialize,
    IndStorage: Serialize,
    DataStorage: Serialize
[src]

impl<N, I, Iptr, IpS, IS, DS> SparseMat for CsMatBase<N, I, IpS, IS, DS, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IpS: Deref<Target = [Iptr]>,
    IS: Deref<Target = [I]>,
    DS: Deref<Target = [N]>, 
[src]

impl<'a, N, I, Iptr, IpS, IS, DS> SparseMat for &'a CsMatBase<N, I, IpS, IS, DS, Iptr> where
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    N: 'a,
    IpS: Deref<Target = [Iptr]>,
    IS: Deref<Target = [I]>,
    DS: Deref<Target = [N]>, 
[src]

impl<N, I, IptrStorage, IndStorage, DataStorage, Iptr> StructuralEq for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

impl<N, I, IptrStorage, IndStorage, DataStorage, Iptr> StructuralPartialEq for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    IptrStorage: Deref<Target = [Iptr]>,
    IndStorage: Deref<Target = [I]>,
    DataStorage: Deref<Target = [N]>, 
[src]

impl<'a, 'b, N, I, Iptr, IpStorage, IStorage, DStorage, Mat> Sub<&'b Mat> for &'a CsMatBase<N, I, IpStorage, IStorage, DStorage, Iptr> where
    N: 'a + Copy + Num + Default,
    I: 'a + SpIndex,
    Iptr: 'a + SpIndex,
    IpStorage: 'a + Deref<Target = [Iptr]>,
    IStorage: 'a + Deref<Target = [I]>,
    DStorage: 'a + Deref<Target = [N]>,
    Mat: SpMatView<N, I, Iptr>, 
[src]

type Output = CsMatI<N, I, Iptr>

The resulting type after applying the - operator.

impl<N, I, Iptr, IS1, DS1, ISptr1, IS2, ISptr2, DS2> UlpsEq<CsMatBase<N, I, ISptr2, IS2, DS2, Iptr>> for CsMatBase<N, I, ISptr1, IS1, DS1, Iptr> where
    I: SpIndex,
    Iptr: SpIndex,
    CsMatBase<N, I, ISptr1, IS1, DS1, Iptr>: PartialEq<CsMatBase<N, I, ISptr2, IS2, DS2, Iptr>>,
    IS1: Deref<Target = [I]>,
    IS2: Deref<Target = [I]>,
    ISptr1: Deref<Target = [Iptr]>,
    ISptr2: Deref<Target = [Iptr]>,
    DS1: Deref<Target = [N]>,
    DS2: Deref<Target = [N]>,
    N: UlpsEq,
    N::Epsilon: Clone,
    N: Zero
[src]

Auto Trait Implementations

impl<N, I, IptrStorage, IndStorage, DataStorage, Iptr> RefUnwindSafe for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    DataStorage: RefUnwindSafe,
    IndStorage: RefUnwindSafe,
    IptrStorage: RefUnwindSafe
[src]

impl<N, I, IptrStorage, IndStorage, DataStorage, Iptr> Send for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    DataStorage: Send,
    IndStorage: Send,
    IptrStorage: Send
[src]

impl<N, I, IptrStorage, IndStorage, DataStorage, Iptr> Sync for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    DataStorage: Sync,
    IndStorage: Sync,
    IptrStorage: Sync
[src]

impl<N, I, IptrStorage, IndStorage, DataStorage, Iptr> Unpin for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    DataStorage: Unpin,
    IndStorage: Unpin,
    IptrStorage: Unpin
[src]

impl<N, I, IptrStorage, IndStorage, DataStorage, Iptr> UnwindSafe for CsMatBase<N, I, IptrStorage, IndStorage, DataStorage, Iptr> where
    DataStorage: UnwindSafe,
    IndStorage: UnwindSafe,
    IptrStorage: UnwindSafe
[src]

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> Pointable for T

type Init = T

The type for initializers.

impl<SS, SP> SupersetOf<SS> for SP where
    SS: SubsetOf<SP>, 
[src]

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

type Owned = T

The resulting type after obtaining ownership.

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.