Skip to main content

LinalgBackend

Trait LinalgBackend 

Source
pub trait LinalgBackend: TensorBackend {
    // Required methods
    fn cholesky(&mut self, input: &Tensor) -> Result<Tensor>;
    fn triangular_solve(
        &mut self,
        a: &Tensor,
        b: &Tensor,
        left_side: bool,
        lower: bool,
        transpose_a: bool,
        unit_diagonal: bool,
    ) -> Result<Tensor>;
    fn lu(&mut self, input: &Tensor) -> Result<Vec<Tensor>>;
    fn full_piv_lu(&mut self, input: &Tensor) -> Result<Vec<Tensor>>;
    fn full_piv_lu_solve(
        &mut self,
        a: &Tensor,
        b: &Tensor,
        transpose_a: bool,
    ) -> Result<Tensor>;
    fn svd(&mut self, input: &Tensor) -> Result<Vec<Tensor>>;
    fn qr(&mut self, input: &Tensor) -> Result<Vec<Tensor>>;
    fn eigh(&mut self, input: &Tensor) -> Result<Vec<Tensor>>;
    fn eig(&mut self, input: &Tensor) -> Result<Vec<Tensor>>;
    fn solve(&mut self, a: &Tensor, b: &Tensor) -> Result<Tensor>;

    // Provided method
    fn svd_read(&mut self, _input: TensorView<'_>) -> Result<Vec<Tensor>> { ... }
}
Expand description

Backend surface required by the linalg extension runtime.

§Examples

use tenferro_linalg::backend::LinalgBackend;
use tenferro_cpu::CpuBackend;

fn accepts_linalg_backend<B: LinalgBackend>(_backend: &mut B) {}

let mut backend = CpuBackend::new();
accepts_linalg_backend(&mut backend);

Required Methods§

Source

fn cholesky(&mut self, input: &Tensor) -> Result<Tensor>

Compute a Cholesky factorization.

Source

fn triangular_solve( &mut self, a: &Tensor, b: &Tensor, left_side: bool, lower: bool, transpose_a: bool, unit_diagonal: bool, ) -> Result<Tensor>

Solve a triangular linear system with explicit side, triangle, transpose, and unit-diagonal flags.

Source

fn lu(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Compute public LU outputs (P, L, U, parity).

Source

fn full_piv_lu(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Compute complete-pivot LU outputs (P, L, U, Q, parity).

The reconstruction convention is A = P^T * L * U * Q, equivalently P * A * Q^T = L * U. parity is a scalar real tensor containing +1 or -1: F32 for F32/C32 inputs and F64 for F64/C64 inputs.

Source

fn full_piv_lu_solve( &mut self, a: &Tensor, b: &Tensor, transpose_a: bool, ) -> Result<Tensor>

Solve a linear system through the complete-pivot LU path.

With transpose_a = false, this solves A * x = b. With transpose_a = true, this solves A^T * x = b.

Source

fn svd(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Compute public SVD outputs (U, S, Vt).

Source

fn qr(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Compute public QR outputs (Q, R).

QR is thin: for an m x n input, Q has shape m x min(m, n) and R has shape min(m, n) x n.

Source

fn eigh(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Compute public Hermitian eigendecomposition outputs (values, vectors).

The returned vector order is [values, vectors], where values has shape [n] and vectors has shape [n, n].

Source

fn eig(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Compute public general eigendecomposition outputs (values, vectors).

Source

fn solve(&mut self, a: &Tensor, b: &Tensor) -> Result<Tensor>

Solve a dense linear system.

Provided Methods§

Source

fn svd_read(&mut self, _input: TensorView<'_>) -> Result<Vec<Tensor>>

Compute a singular value decomposition from a borrowed tensor view.

Backends may canonicalize the view inside the same placement family, but must not silently transfer between CPU and GPU memory.

§Examples
use tenferro_linalg::LinalgBackend;
use tenferro_cpu::CpuBackend;
use tenferro_tensor::{TensorView, TypedTensor};

let input = TypedTensor::<f64>::from_vec_col_major(
    vec![2, 2],
    vec![1.0, 0.0, 0.0, 2.0],
)?;
let outputs = CpuBackend::new().svd_read(TensorView::F64(input.as_view()))?;
assert_eq!(outputs[1].shape(), &[2]);

Implementations on Foreign Types§

Source§

impl LinalgBackend for EagerBackend

Source§

fn cholesky(&mut self, input: &Tensor) -> Result<Tensor>

Source§

fn triangular_solve( &mut self, a: &Tensor, b: &Tensor, left_side: bool, lower: bool, transpose_a: bool, unit_diagonal: bool, ) -> Result<Tensor>

Source§

fn lu(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn full_piv_lu(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn full_piv_lu_solve( &mut self, a: &Tensor, b: &Tensor, transpose_a: bool, ) -> Result<Tensor>

Source§

fn svd(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn svd_read(&mut self, input: TensorView<'_>) -> Result<Vec<Tensor>>

Source§

fn qr(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn eigh(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn eig(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn solve(&mut self, a: &Tensor, b: &Tensor) -> Result<Tensor>

Source§

impl LinalgBackend for CpuBackend

Source§

fn cholesky(&mut self, input: &Tensor) -> Result<Tensor>

Source§

fn triangular_solve( &mut self, a: &Tensor, b: &Tensor, left_side: bool, lower: bool, transpose_a: bool, unit_diagonal: bool, ) -> Result<Tensor>

Source§

fn lu(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn full_piv_lu(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn full_piv_lu_solve( &mut self, a: &Tensor, b: &Tensor, transpose_a: bool, ) -> Result<Tensor>

Source§

fn svd(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn svd_read(&mut self, input: TensorView<'_>) -> Result<Vec<Tensor>>

Source§

fn qr(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn eigh(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn eig(&mut self, input: &Tensor) -> Result<Vec<Tensor>>

Source§

fn solve(&mut self, a: &Tensor, b: &Tensor) -> Result<Tensor>

Implementors§