pub enum Tensor {
F32(TypedTensor<f32>),
F64(TypedTensor<f64>),
C32(TypedTensor<Complex<f32>>),
C64(TypedTensor<Complex<f64>>),
}Expand description
Variants§
F32(TypedTensor<f32>)
F64(TypedTensor<f64>)
C32(TypedTensor<Complex<f32>>)
C64(TypedTensor<Complex<f64>>)
Implementations§
Source§impl Tensor
impl Tensor
Sourcepub fn from_vec<T>(shape: Vec<usize>, data: Vec<T>) -> Tensorwhere
T: TensorScalar,
pub fn from_vec<T>(shape: Vec<usize>, data: Vec<T>) -> Tensorwhere
T: TensorScalar,
Create a tensor from a shape and flat data.
This is the Tensor-level equivalent of TypedTensor::<T>::from_vec.
§Examples
use tenferro_tensor::Tensor;
let t = Tensor::from_vec(vec![2, 3], vec![1.0_f64, 2.0, 3.0, 4.0, 5.0, 6.0]);
assert_eq!(t.shape(), &[2, 3]);
assert_eq!(t.as_slice::<f64>().unwrap(), &[1.0, 2.0, 3.0, 4.0, 5.0, 6.0]);Sourcepub fn as_slice<T>(&self) -> Option<&[T]>where
T: TensorScalar,
pub fn as_slice<T>(&self) -> Option<&[T]>where
T: TensorScalar,
Try to borrow the host data as a typed slice.
Returns None if the tensor dtype does not match T.
§Examples
use tenferro_tensor::{Tensor, TypedTensor};
let t = Tensor::F64(TypedTensor::from_vec(vec![3], vec![1.0, 2.0, 3.0]));
assert_eq!(t.as_slice::<f64>(), Some([1.0, 2.0, 3.0].as_slice()));
assert_eq!(t.as_slice::<f32>(), None);Sourcepub fn svd(
&self,
ctx: &mut impl TensorBackend,
) -> Result<(Tensor, Tensor, Tensor), Error>
pub fn svd( &self, ctx: &mut impl TensorBackend, ) -> Result<(Tensor, Tensor, Tensor), Error>
Singular value decomposition: A = U diag(S) Vt.
Returns (U, S, Vt) using the thin/economy SVD.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![3, 2], vec![1.0_f64, 2.0, 3.0, 4.0, 5.0, 6.0]);
let (u, s, vt) = a.svd(&mut ctx).unwrap();
assert_eq!(u.shape(), &[3, 2]);
assert_eq!(s.shape(), &[2]);
assert_eq!(vt.shape(), &[2, 2]);Sourcepub fn qr(
&self,
ctx: &mut impl TensorBackend,
) -> Result<(Tensor, Tensor), Error>
pub fn qr( &self, ctx: &mut impl TensorBackend, ) -> Result<(Tensor, Tensor), Error>
QR decomposition: A = Q R.
Returns (Q, R) using the thin/economy QR decomposition.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![3, 2], vec![1.0_f64, 2.0, 3.0, 4.0, 5.0, 6.0]);
let (q, r) = a.qr(&mut ctx).unwrap();
assert_eq!(q.shape(), &[3, 2]);
assert_eq!(r.shape(), &[2, 2]);Sourcepub fn lu(
&self,
ctx: &mut impl TensorBackend,
) -> Result<(Tensor, Tensor, Tensor, Tensor), Error>
pub fn lu( &self, ctx: &mut impl TensorBackend, ) -> Result<(Tensor, Tensor, Tensor, Tensor), Error>
LU decomposition with partial pivoting: P A = L U.
Returns (P, L, U, parity).
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 2], vec![0.0_f64, 1.0, 1.0, 0.0]);
let (p, l, u, parity) = a.lu(&mut ctx).unwrap();
assert_eq!(p.shape(), &[2, 2]);
assert_eq!(l.shape(), &[2, 2]);
assert_eq!(u.shape(), &[2, 2]);
assert_eq!(parity.shape(), &[] as &[usize]);Sourcepub fn cholesky(&self, ctx: &mut impl TensorBackend) -> Result<Tensor, Error>
pub fn cholesky(&self, ctx: &mut impl TensorBackend) -> Result<Tensor, Error>
Cholesky decomposition: A = L L^T or A = L L^H for complex inputs.
Returns the lower-triangular factor L.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 2], vec![4.0_f64, 1.0, 1.0, 3.0]);
let l = a.cholesky(&mut ctx).unwrap();
assert_eq!(l.shape(), &[2, 2]);Sourcepub fn eigh(
&self,
ctx: &mut impl TensorBackend,
) -> Result<(Tensor, Tensor), Error>
pub fn eigh( &self, ctx: &mut impl TensorBackend, ) -> Result<(Tensor, Tensor), Error>
Symmetric or Hermitian eigendecomposition: A = V diag(W) V^T.
Returns (eigenvalues, eigenvectors).
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 2], vec![4.0_f64, 1.0, 1.0, 3.0]);
let (w, v) = a.eigh(&mut ctx).unwrap();
assert_eq!(w.shape(), &[2]);
assert_eq!(v.shape(), &[2, 2]);Sourcepub fn eig(
&self,
ctx: &mut impl TensorBackend,
) -> Result<(Tensor, Tensor), Error>
pub fn eig( &self, ctx: &mut impl TensorBackend, ) -> Result<(Tensor, Tensor), Error>
General eigendecomposition.
Returns (eigenvalues, eigenvectors).
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 2], vec![1.0_f64, 0.0, 0.0, 3.0]);
let (w, v) = a.eig(&mut ctx).unwrap();
assert_eq!(w.shape(), &[2]);
assert_eq!(v.shape(), &[2, 2]);Sourcepub fn solve(
&self,
b: &Tensor,
ctx: &mut impl TensorBackend,
) -> Result<Tensor, Error>
pub fn solve( &self, b: &Tensor, ctx: &mut impl TensorBackend, ) -> Result<Tensor, Error>
Solve A x = b for x.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 2], vec![2.0_f64, 1.0, 1.0, 2.0]);
let b = Tensor::from_vec(vec![2, 1], vec![1.0_f64, 0.0]);
let x = a.solve(&b, &mut ctx).unwrap();
assert_eq!(x.shape(), &[2, 1]);Sourcepub fn triangular_solve(
&self,
b: &Tensor,
left_side: bool,
lower: bool,
transpose_a: bool,
unit_diagonal: bool,
ctx: &mut impl TensorBackend,
) -> Result<Tensor, Error>
pub fn triangular_solve( &self, b: &Tensor, left_side: bool, lower: bool, transpose_a: bool, unit_diagonal: bool, ctx: &mut impl TensorBackend, ) -> Result<Tensor, Error>
Solve a triangular system.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 2], vec![2.0_f64, 1.0, 0.0, 3.0]);
let b = Tensor::from_vec(vec![2, 1], vec![2.0_f64, 7.0]);
let x = a
.triangular_solve(&b, true, true, false, false, &mut ctx)
.unwrap();
assert_eq!(x.shape(), &[2, 1]);Sourcepub fn add(
&self,
other: &Tensor,
ctx: &mut impl TensorBackend,
) -> Result<Tensor, Error>
pub fn add( &self, other: &Tensor, ctx: &mut impl TensorBackend, ) -> Result<Tensor, Error>
Elementwise addition.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![3], vec![1.0_f64, 2.0, 3.0]);
let b = Tensor::from_vec(vec![3], vec![4.0_f64, 5.0, 6.0]);
let c = a.add(&b, &mut ctx).unwrap();
assert_eq!(c.as_slice::<f64>().unwrap(), &[5.0, 7.0, 9.0]);Sourcepub fn mul(
&self,
other: &Tensor,
ctx: &mut impl TensorBackend,
) -> Result<Tensor, Error>
pub fn mul( &self, other: &Tensor, ctx: &mut impl TensorBackend, ) -> Result<Tensor, Error>
Elementwise multiplication.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![3], vec![1.0_f64, 2.0, 3.0]);
let b = Tensor::from_vec(vec![3], vec![4.0_f64, 5.0, 6.0]);
let c = a.mul(&b, &mut ctx).unwrap();
assert_eq!(c.as_slice::<f64>().unwrap(), &[4.0, 10.0, 18.0]);Sourcepub fn neg(&self, ctx: &mut impl TensorBackend) -> Result<Tensor, Error>
pub fn neg(&self, ctx: &mut impl TensorBackend) -> Result<Tensor, Error>
Negation.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![3], vec![1.0_f64, -2.0, 3.0]);
let b = a.neg(&mut ctx).unwrap();
assert_eq!(b.as_slice::<f64>().unwrap(), &[-1.0, 2.0, -3.0]);Sourcepub fn transpose(
&self,
perm: &[usize],
ctx: &mut impl TensorBackend,
) -> Result<Tensor, Error>
pub fn transpose( &self, perm: &[usize], ctx: &mut impl TensorBackend, ) -> Result<Tensor, Error>
Transpose with an explicit permutation.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 2], vec![1.0_f64, 2.0, 3.0, 4.0]);
let b = a.transpose(&[1, 0], &mut ctx).unwrap();
assert_eq!(b.shape(), &[2, 2]);
assert_eq!(b.as_slice::<f64>().unwrap(), &[1.0, 3.0, 2.0, 4.0]);Sourcepub fn reshape(
&self,
shape: &[usize],
ctx: &mut impl TensorBackend,
) -> Result<Tensor, Error>
pub fn reshape( &self, shape: &[usize], ctx: &mut impl TensorBackend, ) -> Result<Tensor, Error>
Reshape to a new shape with the same number of elements.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 3], vec![1.0_f64, 2.0, 3.0, 4.0, 5.0, 6.0]);
let b = a.reshape(&[3, 2], &mut ctx).unwrap();
assert_eq!(b.shape(), &[3, 2]);
assert_eq!(b.as_slice::<f64>().unwrap(), &[1.0, 2.0, 3.0, 4.0, 5.0, 6.0]);Sourcepub fn reduce_sum(
&self,
axes: &[usize],
ctx: &mut impl TensorBackend,
) -> Result<Tensor, Error>
pub fn reduce_sum( &self, axes: &[usize], ctx: &mut impl TensorBackend, ) -> Result<Tensor, Error>
Reduce sum over the specified axes.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 3], vec![1.0_f64, 2.0, 3.0, 4.0, 5.0, 6.0]);
let b = a.reduce_sum(&[1], &mut ctx).unwrap();
assert_eq!(b.shape(), &[2]);
assert_eq!(b.as_slice::<f64>().unwrap(), &[9.0, 12.0]);Sourcepub fn matmul(
&self,
other: &Tensor,
ctx: &mut impl TensorBackend,
) -> Result<Tensor, Error>
pub fn matmul( &self, other: &Tensor, ctx: &mut impl TensorBackend, ) -> Result<Tensor, Error>
Matrix multiplication for rank-2 tensors.
This is a convenience wrapper around dot_general.
§Examples
use tenferro_tensor::{Tensor, TensorBackend, cpu::CpuBackend};
let mut ctx = CpuBackend::new();
let a = Tensor::from_vec(vec![2, 3], vec![1.0_f64, 2.0, 3.0, 4.0, 5.0, 6.0]);
let b = Tensor::from_vec(vec![3, 2], vec![1.0_f64, 2.0, 3.0, 4.0, 5.0, 6.0]);
let c = a.matmul(&b, &mut ctx).unwrap();
assert_eq!(c.shape(), &[2, 2]);
assert_eq!(c.as_slice::<f64>().unwrap(), &[22.0, 28.0, 49.0, 64.0]);Trait Implementations§
Source§impl From<TypedTensor<Complex<f32>>> for Tensor
Wrap a Complex32 TypedTensor into the corresponding Tensor
variant.
impl From<TypedTensor<Complex<f32>>> for Tensor
Wrap a Complex32 TypedTensor into the corresponding Tensor
variant.
§Examples
use num_complex::Complex32;
use tenferro_tensor::{Tensor, TypedTensor};
let typed = TypedTensor::from_vec(
vec![1],
vec![Complex32::new(1.0, 2.0)],
);
let tensor: Tensor = typed.into();
assert_eq!(tensor.shape(), &[1]);Source§impl From<TypedTensor<Complex<f64>>> for Tensor
Wrap a Complex64 TypedTensor into the corresponding Tensor
variant.
impl From<TypedTensor<Complex<f64>>> for Tensor
Wrap a Complex64 TypedTensor into the corresponding Tensor
variant.
§Examples
use num_complex::Complex64;
use tenferro_tensor::{Tensor, TypedTensor};
let typed = TypedTensor::from_vec(
vec![1],
vec![Complex64::new(1.0, 2.0)],
);
let tensor: Tensor = typed.into();
assert_eq!(tensor.shape(), &[1]);Source§impl From<TypedTensor<f32>> for Tensor
Wrap an f32 TypedTensor into the corresponding Tensor variant.
impl From<TypedTensor<f32>> for Tensor
Wrap an f32 TypedTensor into the corresponding Tensor variant.
§Examples
use tenferro_tensor::{Tensor, TypedTensor};
let typed = TypedTensor::from_vec(vec![2], vec![1.0_f32, 2.0]);
let tensor: Tensor = typed.into();
assert_eq!(tensor.shape(), &[2]);Source§impl From<TypedTensor<f64>> for Tensor
Wrap an f64 TypedTensor into the corresponding Tensor variant.
impl From<TypedTensor<f64>> for Tensor
Wrap an f64 TypedTensor into the corresponding Tensor variant.
§Examples
use tenferro_tensor::{Tensor, TypedTensor};
let typed = TypedTensor::from_vec(vec![2], vec![1.0_f64, 2.0]);
let tensor: Tensor = typed.into();
assert_eq!(tensor.shape(), &[2]);Auto Trait Implementations§
impl Freeze for Tensor
impl RefUnwindSafe for Tensor
impl Send for Tensor
impl Sync for Tensor
impl Unpin for Tensor
impl UnsafeUnpin for Tensor
impl UnwindSafe for Tensor
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
§impl<T> DistributionExt for Twhere
T: ?Sized,
impl<T> DistributionExt for Twhere
T: ?Sized,
fn rand<T>(&self, rng: &mut (impl Rng + ?Sized)) -> Twhere
Self: Distribution<T>,
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
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 moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
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