Struct Array

Source
pub struct Array<Item, ArrayImpl, const NDIM: usize>(/* private fields */)
where
    Item: RlstBase,
    ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>;
Expand description

The basic tuple type defining an array.

Implementations§

Source§

impl<Item: RlstScalar + MatrixInverse, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + Shape<2> + RawAccessMut<Item = Item>> Array<Item, ArrayImpl, 2>

Source

pub fn into_inverse_alloc(self) -> RlstResult<()>

Return the matrix inverse.

Source§

impl<Item: RlstScalar + MatrixLu, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + RawAccessMut<Item = Item> + Shape<2>> Array<Item, ArrayImpl, 2>

Source

pub fn into_lu_alloc(self) -> RlstResult<LuDecomposition<Item, ArrayImpl>>

Compute the LU decomposition of a matrix.

The LU Decomposition of an (m,n) matrix A is defined by A = PLU, where P is an (m, m) permutation matrix, L is a (m, k) unit lower triangular matrix, and U is an (k, n) upper triangular matrix.

Source§

impl<Item: RlstScalar + MatrixPseudoInverse, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + Shape<2> + UnsafeRandomAccessMut<2, Item = Item> + RawAccessMut<Item = Item>> Array<Item, ArrayImpl, 2>

Source

pub fn into_pseudo_inverse_resize_alloc<ArrayImplPInv: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + Shape<2> + RawAccessMut<Item = Item> + UnsafeRandomAccessMut<2, Item = Item> + ResizeInPlace<2>>( self, pinv: Array<Item, ArrayImplPInv, 2>, tol: <Item as RlstScalar>::Real, ) -> RlstResult<()>

Compute the pseudoinverse into the array pinv.

This method dynamically allocates memory for the computation.

§Parameters
  • pinv: Array to store the pseudo-inverse in. Array is resized if necessary.
  • tol: The relative tolerance. Singular values smaller or equal tol will be discarded.
Source

pub fn into_pseudo_inverse_alloc<ArrayImplPInv: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + Shape<2> + RawAccessMut<Item = Item> + UnsafeRandomAccessMut<2, Item = Item>>( self, pinv: Array<Item, ArrayImplPInv, 2>, tol: <Item as RlstScalar>::Real, ) -> RlstResult<()>

Compute the pseudo inverse into the array pinv.

This method dynamically allocates memory for the computation.

§Parameters
  • pinv: Array to store the pseudo-inverse in. If self has shape [m, n] then pinv must have shape [n, m].
  • tol: The relative tolerance. Singular values smaller or equal tol * s\[0\] will be discarded, where s[0] is the largest singular value.
Source§

impl<Item: RlstScalar + MatrixQr, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + RawAccessMut<Item = Item> + Shape<2>> Array<Item, ArrayImpl, 2>

Source

pub fn into_qr_alloc(self) -> RlstResult<QrDecomposition<Item, ArrayImpl>>

Compute the QR decomposition of a given 2-dimensional array.

Source§

impl<Item: RlstScalar + MatrixSvd, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + Shape<2> + RawAccessMut<Item = Item>> Array<Item, ArrayImpl, 2>

Source

pub fn into_singular_values_alloc( self, singular_values: &mut [<Item as RlstScalar>::Real], ) -> RlstResult<()>

Compute the singular values of the matrix.

For a (m, n) matrix A the slice singular_values has length k=min(m, n).

This method allocates temporary memory during execution.

Source

pub fn into_svd_alloc<ArrayImplU: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + Shape<2> + RawAccessMut<Item = Item>, ArrayImplVt: UnsafeRandomAccessByValue<2, Item = Item> + Stride<2> + Shape<2> + RawAccessMut<Item = Item>>( self, u: Array<Item, ArrayImplU, 2>, vt: Array<Item, ArrayImplVt, 2>, singular_values: &mut [<Item as RlstScalar>::Real], mode: SvdMode, ) -> RlstResult<()>

Compute the singular value decomposition.

We assume that A is a (m, n) matrix and assume k=min(m, n). This method computes the singular value decomposition A = USVt.

§Parameters
  • u - Stores the matrix U. For the full SVD the shape needs to be (m, m). For the reduced SVD it needs to be (m, k).
  • vt - Stores the matrix Vt. For the full SVD the shape needs to be (n, n). For the reduced SVD it needs to be (k, n). Note that vt stores the complex conjugate transpose of the matrix of right singular vectors. Hence, the columns of vt.transpose().conj() will be the right singular vectors.
  • singular_values - Stores the k singular values of A.
  • mode - Choose between full SVD SvdMode::Full or reduced SVD SvdMode::Reduced.

This method allocates temporary memory during execution.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<ADIM, Item = Item> + Shape<ADIM>, const ADIM: usize> Array<Item, ArrayImpl, ADIM>

Source

pub fn insert_empty_axis<const NDIM: usize>( self, axis_position: AxisPosition, ) -> Array<Item, ArrayAppendAxis<Item, ArrayImpl, ADIM, NDIM>, NDIM>
where NumberType<ADIM>: IsSmallerByOne<NDIM>,

Insert empty axis

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2>> Array<Item, ArrayImpl, 2>

Source

pub fn row_iter(&self) -> RowIterator<'_, Item, ArrayImpl, 2>

Return a row iterator for a two-dimensional array.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + UnsafeRandomAccessMut<2, Item = Item>> Array<Item, ArrayImpl, 2>

Source

pub fn row_iter_mut(&mut self) -> RowIteratorMut<'_, Item, ArrayImpl, 2>

Return a mutable row iterator for a two-dimensional array.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2>> Array<Item, ArrayImpl, 2>

Source

pub fn col_iter(&self) -> ColIterator<'_, Item, ArrayImpl, 2>

Return a column iterator for a two-dimensional array.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + UnsafeRandomAccessMut<2, Item = Item>> Array<Item, ArrayImpl, 2>

Source

pub fn col_iter_mut(&mut self) -> ColIteratorMut<'_, Item, ArrayImpl, 2>

Return a mutable column iterator for a two-dimensional array.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1>> Array<Item, ArrayImpl, 1>

Source

pub fn len(&self) -> usize

Length of 1-d vector

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn is_empty(&self) -> bool

Return true of array is empty (that is one dimension is zero), otherwise false.

Source

pub fn get_diag<ArrayImplOther: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1> + UnsafeRandomAccessMut<1, Item = Item>>( &self, other: Array<Item, ArrayImplOther, 1>, )

Get the diagonal of an array.

Argument must be a 1d array of length self.shape().iter().min().

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + UnsafeRandomAccessMut<NDIM, Item = Item>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn set_diag<ArrayImplOther: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1> + UnsafeRandomAccessMut<1, Item = Item>>( &mut self, other: Array<Item, ArrayImplOther, 1>, )

Set the diagonal of an array.

Argument must be a 1d array of length self.shape().iter().min().

Source

pub fn fill_from<ArrayImplOther: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>>( &mut self, other: Array<Item, ArrayImplOther, NDIM>, )

Fill an array with values from another array.

Source

pub fn fill_from_resize<ArrayImplOther: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>>( &mut self, other: Array<Item, ArrayImplOther, NDIM>, )
where Self: ResizeInPlace<NDIM>,

Fill an array from another array and resize if necessary.

Source

pub fn fill_from_chunked<Other: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + ChunkedAccess<N, Item = Item>, const N: usize>( &mut self, other: Other, )

Fill an array with values from an other arrays using chunks of size N.

Source

pub fn fill_from_chunked_resize<Other: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + ChunkedAccess<N, Item = Item>, const N: usize>( &mut self, other: Other, )
where Self: ResizeInPlace<NDIM>,

Fill an array with values from an other arrays using chunks of size N. Resize if necessary.

Source§

impl<Item: RlstNum, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + UnsafeRandomAccessMut<NDIM, Item = Item>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn set_zero(&mut self)

Set all elements of an array to zero.

Source

pub fn set_one(&mut self)

Set all elements of an array to one.

Source

pub fn set_identity(&mut self)

Fill the diagonal of an array with the value 1 and all other elements zero.

Source

pub fn scale_inplace(&mut self, alpha: Item)

Multiply all array elements with the scalar alpha.

Source

pub fn sum_into<ArrayImplOther: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>>( &mut self, other: Array<Item, ArrayImplOther, NDIM>, )

Sum other array into array.

Source

pub fn cmp_mult_into<ArrayImplOther: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>>( &mut self, other: Array<Item, ArrayImplOther, NDIM>, )

Componentwise multiply other array into array.

Source

pub fn sum_into_chunked<ArrayImplOther: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + ChunkedAccess<N, Item = Item>, const N: usize>( &mut self, other: Array<Item, ArrayImplOther, NDIM>, )
where Self: ChunkedAccess<N, Item = Item>,

Chunked summation into array.

Source

pub fn cmp_mult_into_chunked<ArrayImplOther: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + ChunkedAccess<N, Item = Item>, const N: usize>( &mut self, other: Array<Item, ArrayImplOther, NDIM>, )
where Self: ChunkedAccess<N, Item = Item>,

Chunked componentwise multiplication into array.

Source§

impl<Item: RlstNum, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn trace(self) -> Item

Return the trace of an array.

Source

pub fn sum(self) -> Item

Return the sum of the elements.

Source§

impl<Item: RlstScalar, ArrayImpl: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1>> Array<Item, ArrayImpl, 1>
where Item::Real: Float,

Source

pub fn inner<ArrayImplOther: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1>>( &self, other: Array<Item, ArrayImplOther, 1>, ) -> Item

Compute the inner product between two vectors.

The inner product takes the complex conjugate of the other argument.

Source

pub fn norm_inf(self) -> <Item as RlstScalar>::Real

Compute the maximum (or inf) norm of a vector.

Source

pub fn norm_1(self) -> <Item as RlstScalar>::Real

Compute the 1-norm of a vector.

Source

pub fn norm_2(self) -> <Item as RlstScalar>::Real

Compute the 2-norm of a vector.

Source

pub fn cross<ArrayImplOther: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1>, ArrayImplRes: UnsafeRandomAccessByValue<1, Item = Item> + UnsafeRandomAccessMut<1, Item = Item> + UnsafeRandomAccessByRef<1, Item = Item> + Shape<1>>( &self, other: Array<Item, ArrayImplOther, 1>, res: Array<Item, ArrayImplRes, 1>, )

Compute the cross product with vector other and store into res.

Source§

impl<Item, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccessMut<Item = Item>> Array<Item, ArrayImpl, 2>
where Item: MatrixSvd + RlstScalar,

Source

pub fn norm_2_alloc(self) -> RlstResult<<Item as RlstScalar>::Real>

Compute the 2-norm of a matrix.

This method allocates temporary memory during execution.

Source§

impl<Item: RlstScalar, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2>> Array<Item, ArrayImpl, 2>

Source

pub fn norm_fro(self) -> Item::Real

Compute the Frobenius-norm of a matrix.

Source

pub fn norm_inf(self) -> Item::Real

Compute the inf-norm of a matrix.

Source

pub fn norm_1(self) -> Item::Real

Compute the 1-norm of a matrix.

Source§

impl<Item: RlstScalar, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccessMut<Item = Item>> Array<Item, ArrayImpl, 2>
where LuDecomposition<Item, ArrayImpl>: MatrixLuDecomposition<Item = Item, ArrayImpl = ArrayImpl>,

Source

pub fn into_solve_vec<ArrayImplMut: RawAccessMut<Item = Item> + UnsafeRandomAccessByValue<1, Item = Item> + UnsafeRandomAccessMut<1, Item = Item> + Shape<1> + Stride<1>>( self, trans: TransMode, rhs: Array<Item, ArrayImplMut, 1>, ) -> RlstResult<Array<Item, ArrayImplMut, 1>>

Solve a linear system with a single right-hand side.

The array is overwritten with the LU Decomposition and the right-hand side is overwritten with the solution. The solution is also returned.

Source

pub fn into_det<ArrayImplMut: RawAccessMut<Item = Item> + UnsafeRandomAccessByValue<2, Item = Item> + UnsafeRandomAccessMut<2, Item = Item> + Shape<2> + Stride<2>>( self, ) -> RlstResult<Item>

Compute the determinant of a matrix.

The array is overwritten by the determinant computation.

Source

pub fn into_solve_mat<ArrayImplMut: RawAccessMut<Item = Item> + UnsafeRandomAccessByValue<2, Item = Item> + UnsafeRandomAccessMut<2, Item = Item> + Shape<2> + Stride<2>>( self, trans: TransMode, rhs: Array<Item, ArrayImplMut, 2>, ) -> RlstResult<Array<Item, ArrayImplMut, 2>>

Solve a linear system with multiple right-hand sides.

The array is overwritten with the LU Decomposition and the right-hand side is overwritten with the solution. The solution is also returned.

Source§

impl<Item: RlstNum, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn cast<T: RlstNum>( self, ) -> Array<T, ArrayCast<Item, T, ArrayImpl, NDIM>, NDIM>

Cast array to type T.

Source§

impl<Item: RlstScalar, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn conj(self) -> Array<Item, ArrayConjugate<Item, ArrayImpl, NDIM>, NDIM>

Conjugate

Source§

impl<Item: RlstNum, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn scalar_mul( self, other: Item, ) -> Array<Item, ArrayScalarMult<Item, ArrayImpl, NDIM>, NDIM>

Multiplication by a scalar

Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = f32> + Shape<NDIM>, const NDIM: usize> Array<f32, ArrayImpl, NDIM>

Source

pub fn to_complex( self, ) -> Array<c32, ArrayToComplex<c32, ArrayImpl, NDIM>, NDIM>

Convert to complex

Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = f64> + Shape<NDIM>, const NDIM: usize> Array<f64, ArrayImpl, NDIM>

Source

pub fn to_complex( self, ) -> Array<c64, ArrayToComplex<c64, ArrayImpl, NDIM>, NDIM>

Convert to complex

Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = c64> + Shape<NDIM>, const NDIM: usize> Array<c64, ArrayImpl, NDIM>

Source

pub fn to_complex(self) -> Self

Convert to complex

Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = c32> + Shape<NDIM>, const NDIM: usize> Array<c32, ArrayImpl, NDIM>

Source

pub fn to_complex(self) -> Self

Convert to complex

Source§

impl<Item: RlstNum, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn permute_axes( self, permutation: [usize; NDIM], ) -> Array<Item, ArrayTranspose<Item, ArrayImpl, NDIM>, NDIM>

Permute axes of an array

Source

pub fn transpose( self, ) -> Array<Item, ArrayTranspose<Item, ArrayImpl, NDIM>, NDIM>

Transpose an array

Source§

impl<Item: RlstScalar + RandScalar, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + UnsafeRandomAccessMut<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn fill_from_standard_normal<R: Rng>(&mut self, rng: &mut R)

Fill an array with normally distributed random numbers.

Source

pub fn fill_from_equally_distributed<R: Rng>(&mut self, rng: &mut R)

Fill an array with equally distributed random numbers.

Source

pub fn fill_from_seed_equally_distributed(&mut self, seed: usize)

Fill with equally distributed numbers using a given seed.

Source

pub fn fill_from_seed_normally_distributed(&mut self, seed: usize)

Fill with normally distributed numbers using a given seed.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<ADIM, Item = Item> + Shape<ADIM>, const ADIM: usize> Array<Item, ArrayImpl, ADIM>

Source

pub fn slice<const NDIM: usize>( self, axis: usize, index: usize, ) -> Array<Item, ArraySlice<Item, ArrayImpl, ADIM, NDIM>, NDIM>

Create a slice from a given array.

Consider an array arr with shape [a0, a1, a2, a3, ...]. The function call arr.slice(2, 3) returns a one dimension smaller array indexed by [a0, a1, 3, a3, ...]. Hence, the dimension 2 has been fixed to always have the value 3.

§Examples

If arr is a matrix then the first column of the matrix is obtained from arr.slice(1, 0), while the third row of the matrix is obtained from arr.slice(0, 2).

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn into_subview( self, offset: [usize; NDIM], shape: [usize; NDIM], ) -> Array<Item, ArraySubView<Item, ArrayImpl, NDIM>, NDIM>

Move the array into a subview specified by an offset and shape of the subview.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn view(&self) -> Array<Item, ArrayView<'_, Item, ArrayImpl, NDIM>, NDIM>

Return a view onto the array.

Source

pub fn view_flat( &self, ) -> Array<Item, ArrayFlatView<'_, Item, ArrayImpl, NDIM>, 1>

Return a flattened 1d view onto the array. The view is flattened in column-major order.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + UnsafeRandomAccessMut<NDIM, Item = Item>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn view_mut( &mut self, ) -> Array<Item, ArrayViewMut<'_, Item, ArrayImpl, NDIM>, NDIM>

Return a mutable view onto the array.

Source

pub fn view_flat_mut( &mut self, ) -> Array<Item, ArrayFlatViewMut<'_, Item, ArrayImpl, NDIM>, 1>

Return a flattened 1d view onto the array. The view is flattened in column-major order.

Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Array<Item, ArrayImpl, NDIM>

Source

pub fn new(arr: ArrayImpl) -> Self

Instantiate a new array from an ArrayImpl structure.

Source

pub fn number_of_elements(&self) -> usize

Return the number of elements in the array.

Source§

impl<Item: RlstBase, const NDIM: usize> Array<Item, BaseArray<Item, VectorContainer<Item>, NDIM>, NDIM>

Source

pub fn from_shape(shape: [usize; NDIM]) -> Self

Create a new heap allocated array from a given shape.

Source

pub fn from_shape_with_stride( shape: [usize; NDIM], stride: [usize; NDIM], ) -> Self

Create a new heap allocated array from a given shape.

Source§

impl<'a, Item: RlstBase, const NDIM: usize> Array<Item, BaseArray<Item, SliceContainer<'a, Item>, NDIM>, NDIM>

Source

pub fn from_shape(slice: &'a [Item], shape: [usize; NDIM]) -> Self

Create a new array from a slice with a given shape.

The stride is automatically assumed to be column major.

Source

pub fn from_shape_with_stride( slice: &'a [Item], shape: [usize; NDIM], stride: [usize; NDIM], ) -> Self

Create a new array from a slice with a given shape and stride.

Source§

impl<'a, Item: RlstBase, const NDIM: usize> Array<Item, BaseArray<Item, SliceContainerMut<'a, Item>, NDIM>, NDIM>

Source

pub fn from_shape(slice: &'a mut [Item], shape: [usize; NDIM]) -> Self

Create a new array from a slice with a given shape.

The stride is automatically assumed to be column major.

Source

pub fn from_shape_with_stride( slice: &'a mut [Item], shape: [usize; NDIM], stride: [usize; NDIM], ) -> Self

Create a new array from a slice with a given shape and stride.

Trait Implementations§

Source§

impl<Item: RlstNum, ArrayImpl1: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, ArrayImpl2: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Add<Array<Item, ArrayImpl2, NDIM>> for Array<Item, ArrayImpl1, NDIM>

Source§

type Output = Array<Item, ArrayAddition<Item, ArrayImpl1, ArrayImpl2, NDIM>, NDIM>

The resulting type after applying the + operator.
Source§

fn add(self, rhs: Array<Item, ArrayImpl2, NDIM>) -> Self::Output

Performs the + operation. Read more
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + ChunkedAccess<N, Item = Item>, const NDIM: usize, const N: usize> ChunkedAccess<N> for Array<Item, ArrayImpl, NDIM>

Source§

type Item = Item

Item type
Source§

fn get_chunk(&self, chunk_index: usize) -> Option<DataChunk<Self::Item, N>>

Get chunk
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Debug for Array<Item, ArrayImpl, NDIM>

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> DefaultIterator for Array<Item, ArrayImpl, NDIM>

Source§

type Item = Item

Item type
Source§

type Iter<'a> = ArrayDefaultIterator<'a, Item, ArrayImpl, NDIM> where Self: 'a

Iterator type
Source§

fn iter(&self) -> Self::Iter<'_>

Get iterator
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + UnsafeRandomAccessMut<NDIM, Item = Item>, const NDIM: usize> DefaultIteratorMut for Array<Item, ArrayImpl, NDIM>

Source§

type Item = Item

Item type
Source§

type IterMut<'a> = ArrayDefaultIteratorMut<'a, Item, ArrayImpl, NDIM> where Self: 'a

Iterator
Source§

fn iter_mut(&mut self) -> Self::IterMut<'_>

Get mutable iterator
Source§

impl<Item: RlstNum, ArrayImpl1: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, ArrayImpl2: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Div<Array<Item, ArrayImpl2, NDIM>> for Array<Item, ArrayImpl1, NDIM>

Source§

type Output = Array<Item, CmpWiseDivision<Item, ArrayImpl1, ArrayImpl2, NDIM>, NDIM>

The resulting type after applying the / operator.
Source§

fn div(self, rhs: Array<Item, ArrayImpl2, NDIM>) -> Self::Output

Performs the / operation. Read more
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + UnsafeRandomAccessByRef<NDIM, Item = Item>, const NDIM: usize> Index<[usize; NDIM]> for Array<Item, ArrayImpl, NDIM>

Source§

type Output = Item

The returned type after indexing.
Source§

fn index(&self, index: [usize; NDIM]) -> &Self::Output

Performs the indexing (container[index]) operation. Read more
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + UnsafeRandomAccessByRef<NDIM, Item = Item> + UnsafeRandomAccessMut<NDIM, Item = Item>, const NDIM: usize> IndexMut<[usize; NDIM]> for Array<Item, ArrayImpl, NDIM>

Source§

fn index_mut(&mut self, index: [usize; NDIM]) -> &mut Self::Output

Performs the mutable indexing (container[index]) operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = c32> + Shape<NDIM>, const NDIM: usize> Mul<Array<Complex<f32>, ArrayImpl, NDIM>> for c32

Source§

type Output = Array<Complex<f32>, ArrayScalarMult<Complex<f32>, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<c32, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = c64> + Shape<NDIM>, const NDIM: usize> Mul<Array<Complex<f64>, ArrayImpl, NDIM>> for c64

Source§

type Output = Array<Complex<f64>, ArrayScalarMult<Complex<f64>, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<c64, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<Item: RlstNum, ArrayImpl1: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, ArrayImpl2: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Mul<Array<Item, ArrayImpl2, NDIM>> for Array<Item, ArrayImpl1, NDIM>

Source§

type Output = Array<Item, CmpWiseProduct<Item, ArrayImpl1, ArrayImpl2, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<Item, ArrayImpl2, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = f32> + Shape<NDIM>, const NDIM: usize> Mul<Array<f32, ArrayImpl, NDIM>> for f32

Source§

type Output = Array<f32, ArrayScalarMult<f32, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<f32, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = f64> + Shape<NDIM>, const NDIM: usize> Mul<Array<f64, ArrayImpl, NDIM>> for f64

Source§

type Output = Array<f64, ArrayScalarMult<f64, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<f64, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = i16> + Shape<NDIM>, const NDIM: usize> Mul<Array<i16, ArrayImpl, NDIM>> for i16

Source§

type Output = Array<i16, ArrayScalarMult<i16, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<i16, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = i32> + Shape<NDIM>, const NDIM: usize> Mul<Array<i32, ArrayImpl, NDIM>> for i32

Source§

type Output = Array<i32, ArrayScalarMult<i32, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<i32, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = i64> + Shape<NDIM>, const NDIM: usize> Mul<Array<i64, ArrayImpl, NDIM>> for i64

Source§

type Output = Array<i64, ArrayScalarMult<i64, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<i64, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = i8> + Shape<NDIM>, const NDIM: usize> Mul<Array<i8, ArrayImpl, NDIM>> for i8

Source§

type Output = Array<i8, ArrayScalarMult<i8, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<i8, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = u16> + Shape<NDIM>, const NDIM: usize> Mul<Array<u16, ArrayImpl, NDIM>> for u16

Source§

type Output = Array<u16, ArrayScalarMult<u16, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<u16, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = u32> + Shape<NDIM>, const NDIM: usize> Mul<Array<u32, ArrayImpl, NDIM>> for u32

Source§

type Output = Array<u32, ArrayScalarMult<u32, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<u32, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = u64> + Shape<NDIM>, const NDIM: usize> Mul<Array<u64, ArrayImpl, NDIM>> for u64

Source§

type Output = Array<u64, ArrayScalarMult<u64, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<u64, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = u8> + Shape<NDIM>, const NDIM: usize> Mul<Array<u8, ArrayImpl, NDIM>> for u8

Source§

type Output = Array<u8, ArrayScalarMult<u8, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<u8, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = usize> + Shape<NDIM>, const NDIM: usize> Mul<Array<usize, ArrayImpl, NDIM>> for usize

Source§

type Output = Array<usize, ArrayScalarMult<usize, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: Array<usize, ArrayImpl, NDIM>) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = c32> + Shape<NDIM>, const NDIM: usize> Mul<Complex<f32>> for Array<c32, ArrayImpl, NDIM>

Source§

type Output = Array<Complex<f32>, ArrayScalarMult<Complex<f32>, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: c32) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = c64> + Shape<NDIM>, const NDIM: usize> Mul<Complex<f64>> for Array<c64, ArrayImpl, NDIM>

Source§

type Output = Array<Complex<f64>, ArrayScalarMult<Complex<f64>, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: c64) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = f32> + Shape<NDIM>, const NDIM: usize> Mul<f32> for Array<f32, ArrayImpl, NDIM>

Source§

type Output = Array<f32, ArrayScalarMult<f32, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: f32) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = f64> + Shape<NDIM>, const NDIM: usize> Mul<f64> for Array<f64, ArrayImpl, NDIM>

Source§

type Output = Array<f64, ArrayScalarMult<f64, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: f64) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = i16> + Shape<NDIM>, const NDIM: usize> Mul<i16> for Array<i16, ArrayImpl, NDIM>

Source§

type Output = Array<i16, ArrayScalarMult<i16, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: i16) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = i32> + Shape<NDIM>, const NDIM: usize> Mul<i32> for Array<i32, ArrayImpl, NDIM>

Source§

type Output = Array<i32, ArrayScalarMult<i32, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: i32) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = i64> + Shape<NDIM>, const NDIM: usize> Mul<i64> for Array<i64, ArrayImpl, NDIM>

Source§

type Output = Array<i64, ArrayScalarMult<i64, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: i64) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = i8> + Shape<NDIM>, const NDIM: usize> Mul<i8> for Array<i8, ArrayImpl, NDIM>

Source§

type Output = Array<i8, ArrayScalarMult<i8, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: i8) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = u16> + Shape<NDIM>, const NDIM: usize> Mul<u16> for Array<u16, ArrayImpl, NDIM>

Source§

type Output = Array<u16, ArrayScalarMult<u16, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: u16) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = u32> + Shape<NDIM>, const NDIM: usize> Mul<u32> for Array<u32, ArrayImpl, NDIM>

Source§

type Output = Array<u32, ArrayScalarMult<u32, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: u32) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = u64> + Shape<NDIM>, const NDIM: usize> Mul<u64> for Array<u64, ArrayImpl, NDIM>

Source§

type Output = Array<u64, ArrayScalarMult<u64, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: u64) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = u8> + Shape<NDIM>, const NDIM: usize> Mul<u8> for Array<u8, ArrayImpl, NDIM>

Source§

type Output = Array<u8, ArrayScalarMult<u8, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: u8) -> Self::Output

Performs the * operation. Read more
Source§

impl<ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = usize> + Shape<NDIM>, const NDIM: usize> Mul<usize> for Array<usize, ArrayImpl, NDIM>

Source§

type Output = Array<usize, ArrayScalarMult<usize, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the * operator.
Source§

fn mul(self, rhs: usize) -> Self::Output

Performs the * operation. Read more
Source§

impl<Item: RlstScalar + Gemm, ArrayImpl: UnsafeRandomAccessByValue<1, Item = Item> + UnsafeRandomAccessMut<1, Item = Item> + Shape<1> + Stride<1> + RawAccessMut<Item = Item>, ArrayImplFirst: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1> + Stride<1> + RawAccess<Item = Item>, ArrayImplSecond: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccess<Item = Item>> MultInto<Array<Item, ArrayImplFirst, 1>, Array<Item, ArrayImplSecond, 2>> for Array<Item, ArrayImpl, 1>

Source§

type Item = Item

Item type
Source§

fn mult_into( self, transa: TransMode, transb: TransMode, alpha: Item, arr_a: Array<Item, ArrayImplFirst, 1>, arr_b: Array<Item, ArrayImplSecond, 2>, beta: Item, ) -> Self

Multiply into
Source§

fn simple_mult_into(self, arr_a: First, arr_b: Second) -> Self
where Self: Sized,

Multiply First * Second and sum into Self
Source§

impl<Item: RlstScalar + Gemm, ArrayImpl: UnsafeRandomAccessByValue<1, Item = Item> + UnsafeRandomAccessMut<1, Item = Item> + Shape<1> + Stride<1> + RawAccessMut<Item = Item>, ArrayImplFirst: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccess<Item = Item>, ArrayImplSecond: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1> + Stride<1> + RawAccess<Item = Item>> MultInto<Array<Item, ArrayImplFirst, 2>, Array<Item, ArrayImplSecond, 1>> for Array<Item, ArrayImpl, 1>

Source§

type Item = Item

Item type
Source§

fn mult_into( self, transa: TransMode, transb: TransMode, alpha: Item, arr_a: Array<Item, ArrayImplFirst, 2>, arr_b: Array<Item, ArrayImplSecond, 1>, beta: Item, ) -> Self

Multiply into
Source§

fn simple_mult_into(self, arr_a: First, arr_b: Second) -> Self
where Self: Sized,

Multiply First * Second and sum into Self
Source§

impl<Item: RlstScalar, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + UnsafeRandomAccessMut<2, Item = Item> + Shape<2> + Stride<2> + RawAccessMut<Item = Item>, ArrayImplFirst: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccess<Item = Item>, ArrayImplSecond: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccess<Item = Item>> MultInto<Array<Item, ArrayImplFirst, 2>, Array<Item, ArrayImplSecond, 2>> for Array<Item, ArrayImpl, 2>

Source§

type Item = Item

Item type
Source§

fn mult_into( self, transa: TransMode, transb: TransMode, alpha: Item, arr_a: Array<Item, ArrayImplFirst, 2>, arr_b: Array<Item, ArrayImplSecond, 2>, beta: Item, ) -> Self

Multiply into
Source§

fn simple_mult_into(self, arr_a: First, arr_b: Second) -> Self
where Self: Sized,

Multiply First * Second and sum into Self
Source§

impl<Item: RlstScalar + Gemm, ArrayImpl: UnsafeRandomAccessByValue<1, Item = Item> + UnsafeRandomAccessMut<1, Item = Item> + Shape<1> + Stride<1> + RawAccessMut<Item = Item> + ResizeInPlace<1>, ArrayImplFirst: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1> + Stride<1> + RawAccess<Item = Item>, ArrayImplSecond: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccess<Item = Item>> MultIntoResize<Array<Item, ArrayImplFirst, 1>, Array<Item, ArrayImplSecond, 2>> for Array<Item, ArrayImpl, 1>

Source§

type Item = Item

Item type
Source§

fn mult_into_resize( self, transa: TransMode, transb: TransMode, alpha: Item, arr_a: Array<Item, ArrayImplFirst, 1>, arr_b: Array<Item, ArrayImplSecond, 2>, beta: Item, ) -> Self

Multiply into with resize
Source§

fn simple_mult_into_resize(self, arr_a: First, arr_b: Second) -> Self
where Self: Sized,

Multiply First * Second and sum into Self. Allow to resize Self if necessary
Source§

impl<Item: RlstScalar + Gemm, ArrayImpl: UnsafeRandomAccessByValue<1, Item = Item> + UnsafeRandomAccessMut<1, Item = Item> + Shape<1> + Stride<1> + RawAccessMut<Item = Item> + ResizeInPlace<1>, ArrayImplFirst: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccess<Item = Item>, ArrayImplSecond: UnsafeRandomAccessByValue<1, Item = Item> + Shape<1> + Stride<1> + RawAccess<Item = Item>> MultIntoResize<Array<Item, ArrayImplFirst, 2>, Array<Item, ArrayImplSecond, 1>> for Array<Item, ArrayImpl, 1>

Source§

type Item = Item

Item type
Source§

fn mult_into_resize( self, transa: TransMode, transb: TransMode, alpha: Item, arr_a: Array<Item, ArrayImplFirst, 2>, arr_b: Array<Item, ArrayImplSecond, 1>, beta: Item, ) -> Self

Multiply into with resize
Source§

fn simple_mult_into_resize(self, arr_a: First, arr_b: Second) -> Self
where Self: Sized,

Multiply First * Second and sum into Self. Allow to resize Self if necessary
Source§

impl<Item: RlstScalar + Gemm, ArrayImpl: UnsafeRandomAccessByValue<2, Item = Item> + UnsafeRandomAccessMut<2, Item = Item> + Shape<2> + Stride<2> + RawAccessMut<Item = Item> + ResizeInPlace<2>, ArrayImplFirst: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccess<Item = Item>, ArrayImplSecond: UnsafeRandomAccessByValue<2, Item = Item> + Shape<2> + Stride<2> + RawAccess<Item = Item>> MultIntoResize<Array<Item, ArrayImplFirst, 2>, Array<Item, ArrayImplSecond, 2>> for Array<Item, ArrayImpl, 2>

MultIntoResize

Source§

type Item = Item

Item type
Source§

fn mult_into_resize( self, transa: TransMode, transb: TransMode, alpha: Item, arr_a: Array<Item, ArrayImplFirst, 2>, arr_b: Array<Item, ArrayImplSecond, 2>, beta: Item, ) -> Self

Multiply into with resize
Source§

fn simple_mult_into_resize(self, arr_a: First, arr_b: Second) -> Self
where Self: Sized,

Multiply First * Second and sum into Self. Allow to resize Self if necessary
Source§

impl<Item: RlstScalar, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Neg for Array<Item, ArrayImpl, NDIM>

Source§

type Output = Array<Item, ArrayScalarMult<Item, ArrayImpl, NDIM>, NDIM>

The resulting type after applying the - operator.
Source§

fn neg(self) -> Self::Output

Performs the unary - operation. Read more
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> NumberOfElements for Array<Item, ArrayImpl, NDIM>

Source§

fn number_of_elements(&self) -> usize

Return the number of elements.
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + RawAccess<Item = Item>, const NDIM: usize> RawAccess for Array<Item, ArrayImpl, NDIM>

Source§

type Item = Item

Item type
Source§

fn data(&self) -> &[Self::Item]

Get a slice of the whole data.
Source§

fn buff_ptr(&self) -> *const Self::Item

Return a pointer to the start of the underlying buffer.
Source§

fn offset(&self) -> usize

Offset of beginning of data from start ptr.
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + RawAccessMut<Item = Item>, const NDIM: usize> RawAccessMut for Array<Item, ArrayImpl, NDIM>

Source§

fn data_mut(&mut self) -> &mut [Self::Item]

Get a mutable slice of the whole data.
Source§

fn buff_ptr_mut(&mut self) -> *mut Self::Item

Return a pointer to the start of the underlying buffer.
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + ResizeInPlace<NDIM>, const NDIM: usize> ResizeInPlace<NDIM> for Array<Item, ArrayImpl, NDIM>

Source§

fn resize_in_place(&mut self, shape: [usize; NDIM])

Resize an operator in place
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Shape<NDIM> for Array<Item, ArrayImpl, NDIM>

Source§

fn shape(&self) -> [usize; NDIM]

Return the shape of the object.
Source§

fn is_empty(&self) -> bool

Return true if a dimension is 0.
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + Stride<NDIM>, const NDIM: usize> Stride<NDIM> for Array<Item, ArrayImpl, NDIM>

Source§

fn stride(&self) -> [usize; NDIM]

Return the stride of the object.
Source§

impl<Item: RlstNum, ArrayImpl1: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, ArrayImpl2: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> Sub<Array<Item, ArrayImpl2, NDIM>> for Array<Item, ArrayImpl1, NDIM>

Source§

type Output = Array<Item, ArraySubtraction<Item, ArrayImpl1, ArrayImpl2, NDIM>, NDIM>

The resulting type after applying the - operator.
Source§

fn sub(self, rhs: Array<Item, ArrayImpl2, NDIM>) -> Self::Output

Performs the - operation. Read more
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + UnsafeRandomAccessByRef<NDIM, Item = Item>, const NDIM: usize> UnsafeRandomAccessByRef<NDIM> for Array<Item, ArrayImpl, NDIM>

Source§

type Item = Item

Item type
Source§

unsafe fn get_unchecked(&self, multi_index: [usize; NDIM]) -> &Self::Item

Return a mutable reference to the element at position determined by multi_index. Read more
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>, const NDIM: usize> UnsafeRandomAccessByValue<NDIM> for Array<Item, ArrayImpl, NDIM>

Source§

type Item = Item

Item type
Source§

unsafe fn get_value_unchecked(&self, multi_index: [usize; NDIM]) -> Self::Item

Return the element at position determined by multi_index. Read more
Source§

impl<Item: RlstBase, ArrayImpl: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM> + UnsafeRandomAccessMut<NDIM, Item = Item>, const NDIM: usize> UnsafeRandomAccessMut<NDIM> for Array<Item, ArrayImpl, NDIM>

Source§

type Item = Item

Item type
Source§

unsafe fn get_unchecked_mut( &mut self, multi_index: [usize; NDIM], ) -> &mut Self::Item

Return a mutable reference to the element at position determined by multi_index. Read more

Auto Trait Implementations§

§

impl<Item, ArrayImpl, const NDIM: usize> Freeze for Array<Item, ArrayImpl, NDIM>
where ArrayImpl: Freeze,

§

impl<Item, ArrayImpl, const NDIM: usize> RefUnwindSafe for Array<Item, ArrayImpl, NDIM>
where ArrayImpl: RefUnwindSafe,

§

impl<Item, ArrayImpl, const NDIM: usize> Send for Array<Item, ArrayImpl, NDIM>
where ArrayImpl: Send,

§

impl<Item, ArrayImpl, const NDIM: usize> Sync for Array<Item, ArrayImpl, NDIM>
where ArrayImpl: Sync,

§

impl<Item, ArrayImpl, const NDIM: usize> Unpin for Array<Item, ArrayImpl, NDIM>
where ArrayImpl: Unpin,

§

impl<Item, ArrayImpl, const NDIM: usize> UnwindSafe for Array<Item, ArrayImpl, NDIM>
where ArrayImpl: UnwindSafe,

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> IntoEither for T

Source§

fn into_either(self, into_left: bool) -> Either<Self, Self>

Converts self into a Left variant of Either<Self, Self> if into_left is true. Converts self into a Right variant of Either<Self, Self> otherwise. Read more
Source§

fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
where F: FnOnce(&Self) -> bool,

Converts self into a Left variant of Either<Self, Self> if into_left(&self) returns true. Converts self into a Right variant of Either<Self, Self> otherwise. Read more
Source§

impl<Mat> PrettyPrint<Complex<f32>> for Mat
where Mat: RandomAccessByValue<2, Item = Complex<f32>> + Shape<2>,

Source§

fn pretty_print(&self)

Pretty printing
Source§

fn pretty_print_with_dimension(&self, rows: usize, cols: usize)

Pretty printing with dimension
Source§

fn pretty_print_advanced( &self, rows: usize, cols: usize, print_width: usize, mantissa: usize, exponent: usize, )

Pretty printing with advanced options
Source§

impl<Mat> PrettyPrint<Complex<f64>> for Mat
where Mat: RandomAccessByValue<2, Item = Complex<f64>> + Shape<2>,

Source§

fn pretty_print(&self)

Pretty printing
Source§

fn pretty_print_with_dimension(&self, rows: usize, cols: usize)

Pretty printing with dimension
Source§

fn pretty_print_advanced( &self, rows: usize, cols: usize, print_width: usize, mantissa: usize, exponent: usize, )

Pretty printing with advanced options
Source§

impl<Mat> PrettyPrint<f32> for Mat
where Mat: RandomAccessByValue<2, Item = f32> + Shape<2>,

Source§

fn pretty_print(&self)

Pretty printing
Source§

fn pretty_print_with_dimension(&self, rows: usize, cols: usize)

Pretty printing with dimension
Source§

fn pretty_print_advanced( &self, rows: usize, cols: usize, print_width: usize, mantissa: usize, exponent: usize, )

Pretty printing with advanced options
Source§

impl<Mat> PrettyPrint<f64> for Mat
where Mat: RandomAccessByValue<2, Item = f64> + Shape<2>,

Source§

fn pretty_print(&self)

Pretty printing
Source§

fn pretty_print_with_dimension(&self, rows: usize, cols: usize)

Pretty printing with dimension
Source§

fn pretty_print_advanced( &self, rows: usize, cols: usize, print_width: usize, mantissa: usize, exponent: usize, )

Pretty printing with advanced options
Source§

impl<Item, Mat, const NDIM: usize> RandomAccessByRef<NDIM> for Mat
where Item: RlstBase, Mat: UnsafeRandomAccessByRef<NDIM, Item = Item> + Shape<NDIM>,

Source§

fn get( &self, multi_index: [usize; NDIM], ) -> Option<&<Mat as UnsafeRandomAccessByRef<NDIM>>::Item>

Return a reference to the element at position determined by multi_index.
Source§

impl<Item, Mat, const NDIM: usize> RandomAccessByValue<NDIM> for Mat
where Item: RlstBase, Mat: UnsafeRandomAccessByValue<NDIM, Item = Item> + Shape<NDIM>,

Source§

fn get_value( &self, multi_index: [usize; NDIM], ) -> Option<<Mat as UnsafeRandomAccessByValue<NDIM>>::Item>

Return the element at position determined by multi_index.
Source§

impl<Item, Mat, const NDIM: usize> RandomAccessMut<NDIM> for Mat
where Item: RlstBase, Mat: UnsafeRandomAccessMut<NDIM, Item = Item> + Shape<NDIM>,

Source§

fn get_mut( &mut self, multi_index: [usize; NDIM], ) -> Option<&mut <Mat as UnsafeRandomAccessMut<NDIM>>::Item>

Return a mutable reference to the element at position determined by multi_index.
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

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

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
Source§

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

Source§

fn vzip(self) -> V