Struct sprs::TriMatBase [−][src]
pub struct TriMatBase<IStorage, DStorage> { /* fields omitted */ }
Sparse matrix in the triplet format.
Sparse matrices in the triplet format use three arrays of equal sizes (accessible through the
methods row_inds
, col_inds
, data
), the first one
storing the row indices of non-zero values, the second storing the
corresponding column indices and the last array storing the corresponding
scalar value. If a non-zero location is repeated in the arrays, the
non-zero value is taken as the sum of the corresponding scalar entries.
This format is useful for iteratively building a sparse matrix, since the various non-zero entries can be specified in any order, or even partially as is common in physics with partial derivatives equations.
This format cannot be used for arithmetic operations. Arithmetic operations
are more efficient in the compressed format.
A matrix in the triplet format can be converted to the compressed format
using the methods to_csc
and to_csr
.
The TriMatBase
type is parameterized by the storage type for the row and
column indices, IStorage
, and by the storage type for the non-zero values
DStorage
. Convenient aliases are availaible to specify frequent variant:
TriMat
refers to a triplet matrix owning the storage of its indices and
and values, TriMatView
refers to a triplet matrix with slices to store
its indices and values, while TriMatViewMut
refers to a a triplet matrix
using mutable slices.
Additionaly, the type aliases TriMatI
, TriMatViewI
and
TriMatViewMutI
can be used to choose an index type different from the
default usize
.
Methods
impl<N, I: SpIndex> TriMatBase<Vec<I>, Vec<N>>
[src]
impl<N, I: SpIndex> TriMatBase<Vec<I>, Vec<N>>
pub fn new(shape: (usize, usize)) -> TriMatI<N, I>
[src]
pub fn new(shape: (usize, usize)) -> TriMatI<N, I>
Create a new triplet matrix of shape (nb_rows, nb_cols)
pub fn with_capacity(shape: (usize, usize), cap: usize) -> TriMatI<N, I>
[src]
pub fn with_capacity(shape: (usize, usize), cap: usize) -> TriMatI<N, I>
Create a new triplet matrix of shape (nb_rows, nb_cols)
, and
pre-allocate cap
elements on the backing storage
pub fn from_triplets(
shape: (usize, usize),
row_inds: Vec<I>,
col_inds: Vec<I>,
data: Vec<N>
) -> TriMatI<N, I>
[src]
pub fn from_triplets(
shape: (usize, usize),
row_inds: Vec<I>,
col_inds: Vec<I>,
data: Vec<N>
) -> TriMatI<N, I>
Create a triplet matrix from its raw components. All arrays should have the same length.
Panics
- if the arrays don't have the same length
- if either the row or column indices are out of bounds.
pub fn add_triplet(&mut self, row: usize, col: usize, val: N)
[src]
pub fn add_triplet(&mut self, row: usize, col: usize, val: N)
Append a non-zero triplet to this matrix.
pub fn reserve(&mut self, cap: usize)
[src]
pub fn reserve(&mut self, cap: usize)
Reserve cap
additional non-zeros
pub fn reserve_exact(&mut self, cap: usize)
[src]
pub fn reserve_exact(&mut self, cap: usize)
Reserve exactly cap
non-zeros
impl<N, I: SpIndex, IStorage, DStorage> TriMatBase<IStorage, DStorage> where
IStorage: Deref<Target = [I]>,
DStorage: Deref<Target = [N]>,
[src]
impl<N, I: SpIndex, IStorage, DStorage> TriMatBase<IStorage, DStorage> where
IStorage: Deref<Target = [I]>,
DStorage: Deref<Target = [N]>,
pub fn rows(&self) -> usize
[src]
pub fn rows(&self) -> usize
The number of rows of the matrix
pub fn cols(&self) -> usize
[src]
pub fn cols(&self) -> usize
The number of cols of the matrix
pub fn shape(&self) -> (usize, usize)
[src]
pub fn shape(&self) -> (usize, usize)
The shape of the matrix, as a (rows, cols)
tuple
pub fn nnz(&self) -> usize
[src]
pub fn nnz(&self) -> usize
The number of non-zero entries
pub fn row_inds(&self) -> &[I]
[src]
pub fn row_inds(&self) -> &[I]
The non-zero row indices
pub fn col_inds(&self) -> &[I]
[src]
pub fn col_inds(&self) -> &[I]
The non-zero column indices
pub fn data(&self) -> &[N]
[src]
pub fn data(&self) -> &[N]
The non-zero values
pub fn find_locations(&self, row: usize, col: usize) -> Vec<TripletIndex>
[src]
pub fn find_locations(&self, row: usize, col: usize) -> Vec<TripletIndex>
Find all non-zero entries at the location given by row
and col
pub fn transpose_view(&self) -> TriMatViewI<N, I>
[src]
pub fn transpose_view(&self) -> TriMatViewI<N, I>
Get a transposed view of this matrix
ⓘImportant traits for TriMatIter<RI, CI, DI>pub fn triplet_iter(&self) -> TriMatIter<Iter<I>, Iter<I>, Iter<N>>
[src]
pub fn triplet_iter(&self) -> TriMatIter<Iter<I>, Iter<I>, Iter<N>>
Get an iterator over non-zero elements stored by this matrix
pub fn to_csc(&self) -> CsMatI<N, I> where
N: Clone + Num,
[src]
pub fn to_csc(&self) -> CsMatI<N, I> where
N: Clone + Num,
Create a CSC matrix from this triplet matrix
pub fn to_csr(&self) -> CsMatI<N, I> where
N: Clone + Num,
[src]
pub fn to_csr(&self) -> CsMatI<N, I> where
N: Clone + Num,
Create a CSR matrix from this triplet matrix
pub fn view(&self) -> TriMatViewI<N, I>
[src]
pub fn view(&self) -> TriMatViewI<N, I>
impl<'a, N, I: SpIndex> TriMatBase<&'a [I], &'a [N]>
[src]
impl<'a, N, I: SpIndex> TriMatBase<&'a [I], &'a [N]>
ⓘImportant traits for TriMatIter<RI, CI, DI>pub fn triplet_iter_rbr(
&self
) -> TriMatIter<Iter<'a, I>, Iter<'a, I>, Iter<'a, N>>
[src]
pub fn triplet_iter_rbr(
&self
) -> TriMatIter<Iter<'a, I>, Iter<'a, I>, Iter<'a, N>>
Get an iterator over non-zero elements stored by this matrix
Reborrowing version of triplet_iter()
.
impl<N, I: SpIndex, IStorage, DStorage> TriMatBase<IStorage, DStorage> where
IStorage: DerefMut<Target = [I]>,
DStorage: DerefMut<Target = [N]>,
[src]
impl<N, I: SpIndex, IStorage, DStorage> TriMatBase<IStorage, DStorage> where
IStorage: DerefMut<Target = [I]>,
DStorage: DerefMut<Target = [N]>,
pub fn set_triplet(
&mut self,
TripletIndex: TripletIndex,
row: usize,
col: usize,
val: N
)
[src]
pub fn set_triplet(
&mut self,
TripletIndex: TripletIndex,
row: usize,
col: usize,
val: N
)
Replace a non-zero value at the given index. Indices can be obtained using find_locations.
pub fn view_mut(&mut self) -> TriMatViewMutI<N, I>
[src]
pub fn view_mut(&mut self) -> TriMatViewMutI<N, I>
Trait Implementations
impl<'a, N, I, IS, DS> IntoIterator for &'a TriMatBase<IS, DS> where
I: 'a + SpIndex,
N: 'a,
IS: Deref<Target = [I]>,
DS: Deref<Target = [N]>,
[src]
impl<'a, N, I, IS, DS> IntoIterator for &'a TriMatBase<IS, DS> where
I: 'a + SpIndex,
N: 'a,
IS: Deref<Target = [I]>,
DS: Deref<Target = [N]>,
type Item = (&'a N, (I, I))
The type of the elements being iterated over.
type IntoIter = TriMatIter<Iter<'a, I>, Iter<'a, I>, Iter<'a, N>>
Which kind of iterator are we turning this into?
fn into_iter(self) -> Self::IntoIter
[src]
fn into_iter(self) -> Self::IntoIter
Creates an iterator from a value. Read more
impl<N, I, IS, DS> SparseMat for TriMatBase<IS, DS> where
I: SpIndex,
IS: Deref<Target = [I]>,
DS: Deref<Target = [N]>,
[src]
impl<N, I, IS, DS> SparseMat for TriMatBase<IS, DS> where
I: SpIndex,
IS: Deref<Target = [I]>,
DS: Deref<Target = [N]>,
fn rows(&self) -> usize
[src]
fn rows(&self) -> usize
The number of rows of this matrix
fn cols(&self) -> usize
[src]
fn cols(&self) -> usize
The number of columns of this matrix
fn nnz(&self) -> usize
[src]
fn nnz(&self) -> usize
The number of nonzeros of this matrix
impl<'a, N, I, IS, DS> SparseMat for &'a TriMatBase<IS, DS> where
I: 'a + SpIndex,
N: 'a,
IS: Deref<Target = [I]>,
DS: Deref<Target = [N]>,
[src]
impl<'a, N, I, IS, DS> SparseMat for &'a TriMatBase<IS, DS> where
I: 'a + SpIndex,
N: 'a,
IS: Deref<Target = [I]>,
DS: Deref<Target = [N]>,
fn rows(&self) -> usize
[src]
fn rows(&self) -> usize
The number of rows of this matrix
fn cols(&self) -> usize
[src]
fn cols(&self) -> usize
The number of columns of this matrix
fn nnz(&self) -> usize
[src]
fn nnz(&self) -> usize
The number of nonzeros of this matrix
impl<IStorage: PartialEq, DStorage: PartialEq> PartialEq for TriMatBase<IStorage, DStorage>
[src]
impl<IStorage: PartialEq, DStorage: PartialEq> PartialEq for TriMatBase<IStorage, DStorage>
fn eq(&self, other: &TriMatBase<IStorage, DStorage>) -> bool
[src]
fn eq(&self, other: &TriMatBase<IStorage, DStorage>) -> bool
This method tests for self
and other
values to be equal, and is used by ==
. Read more
fn ne(&self, other: &TriMatBase<IStorage, DStorage>) -> bool
[src]
fn ne(&self, other: &TriMatBase<IStorage, DStorage>) -> bool
This method tests for !=
.
impl<IStorage: Debug, DStorage: Debug> Debug for TriMatBase<IStorage, DStorage>
[src]
impl<IStorage: Debug, DStorage: Debug> Debug for TriMatBase<IStorage, DStorage>
Auto Trait Implementations
impl<IStorage, DStorage> Send for TriMatBase<IStorage, DStorage> where
DStorage: Send,
IStorage: Send,
impl<IStorage, DStorage> Send for TriMatBase<IStorage, DStorage> where
DStorage: Send,
IStorage: Send,
impl<IStorage, DStorage> Sync for TriMatBase<IStorage, DStorage> where
DStorage: Sync,
IStorage: Sync,
impl<IStorage, DStorage> Sync for TriMatBase<IStorage, DStorage> where
DStorage: Sync,
IStorage: Sync,