Reference

Defined in xtensor-blas/xlinalg.hpp

The functions here are closely modeled after NumPy’s linalg package.

Matrix, vector and tensor products

template<class T, class O>
auto xt::linalg::dot(const xexpression<T> &xt, const xexpression<O> &xo)

Non-broadcasting dot function.

In the case of two 1D vectors, computes the vector dot product. In the case of complex vectors, computes the dot product without conjugating the first argument. If t or o is a 2D matrix, computes the matrix-times-vector product. If both t and o ar 2D matrices, computes the matrix-product.

Parameters:
  • t – input array

  • o – input array

Returns:

resulting array

template<class T, class O>
auto xt::linalg::vdot(const xexpression<T> &a, const xexpression<O> &b)

Computes the dot product for two vectors.

Behaves different from dot in the case of complex vectors. If vectors are complex, vdot conjugates the first argument t. Note: Unlike NumPy, xtensor-blas currently doesn’t flatten the input arguments.

Parameters:
  • t – input vector (1D)

  • o – input vector (1D)

Returns:

resulting array

template<class T, class O>
auto xt::linalg::outer(const xexpression<T> &a, const xexpression<O> &b)

Compute the outer product of two vectors.

Parameters:
  • t – input vector (1D)

  • o – input vector (1D)

Returns:

resulting array

template<class E>
auto xt::linalg::matrix_power(const xexpression<E> &A, long n)

Calculate matrix power A**n.

Parameters:
  • mat – The matrix

  • n – The exponent

Returns:

resulting array

template<class T, class E>
auto xt::linalg::kron(const xexpression<T> &a, const xexpression<E> &b)

Calculate the Kronecker product between two 2D xexpressions.

template<class T, class O>
auto xt::linalg::tensordot(const xexpression<T> &xa, const xexpression<O> &xb, std::size_t naxes = 2)

Compute tensor dot product along specified axes for arrays.

Compute the sum of products along the last naxes axes of a and first naxes axes of b.

Parameters:
  • xa – input array

  • xb – input array

  • naxes – the number of axes to sum over

Returns:

resulting array

template<class T, class O>
auto xt::linalg::tensordot(const xexpression<T> &xa, const xexpression<O> &xb, const std::vector<std::size_t> &ax_a, const std::vector<std::size_t> &ax_b)

Compute tensor dot product along specified axes for arrays.

Compute the sum of products along the axes ax_a for a and ax_b for b

Parameters:
  • xa – input array

  • xb – input array

  • ax_a – axes to sum over for a

  • ax_b – axes to sum over for b

Returns:

resulting array

Decompositions

template<class T>
auto xt::linalg::cholesky(const xexpression<T> &A)

Compute the Cholesky decomposition of A.

Returns:

the decomposed matrix

enum class xt::linalg::qrmode

Select the mode for the qr decomposition K = min(M, K)

Values:

enumerator reduced

return Q, R with dimensions (M, K), (K, N) (default)

enumerator complete

return Q, R with dimensions (M, M), (M, N)

enumerator r

return empty Q and R with dimensions (0, 0), (K, N)

enumerator raw

return H, Tau with dimensions (N, M), (K, 1)

template<class T>
auto xt::linalg::qr(const xexpression<T> &A, qrmode mode = qrmode::reduced)

Compute the QR decomposition of A.

Parameters:

t – The matrix to calculate Q and R for

Returns:

std::tuple with Q and R

template<class T>
auto xt::linalg::svd(const xexpression<T> &A, bool full_matrices = true, bool compute_uv = true)

Compute the SVD decomposition of A.

Returns:

tuple containing S, V, and D

Matrix eigenvalues

template<class E, std::enable_if_t<!xtl::is_complex<typename E::value_type>::value>* = nullptr>
auto xt::linalg::eig(const xexpression<E> &A)

Compute the eigenvalues and right eigenvectors of a square array.

Parameters:

Matrix – for which the eigenvalues and right eigenvectors will be computed

Returns:

(eigenvalues, eigenvectors) tuple. The first element corresponds to the eigenvalues, each repeated according to its multiplicity. The eigenvalues are not necessarily ordered. The second (1) element are the normalized (unit “length”) eigenvectors, such that the column v[:, i] corresponds to the eigenvalue w[i].

template<class E, std::enable_if_t<!xtl::is_complex<typename E::value_type>::value>* = nullptr>
auto xt::linalg::eigvals(const xexpression<E> &A)

Compute the eigenvalues of a square xexpression.

Parameters:

A – Matrix for which the eigenvalues are computed

Returns:

xtensor containing the eigenvalues.

template<class E, std::enable_if_t<!xtl::is_complex<typename E::value_type>::value>* = nullptr>
auto xt::linalg::eigh(const xexpression<E> &A, char UPLO = 'L')

Compute the eigenvalues and eigenvectors of a square Hermitian or real symmetric xexpression.

Parameters:

A – Matrix for which the eigenvalues and right eigenvectors are computed

Returns:

xtensor containing the eigenvalues.

template<class E, std::enable_if_t<!xtl::is_complex<typename E::value_type>::value>* = nullptr>
auto xt::linalg::eigh(const xexpression<E> &A, const xexpression<E> &B, const char UPLO = 'L')

Compute the generalized eigenvalues and eigenvectors of a square Hermitian or real symmetric xexpression.

Parameters:

A, B – Matrices for which the generalized eigenvalues and right eigenvectors are computed

Returns:

xtensor containing the eigenvalues.

template<class E, std::enable_if_t<!xtl::is_complex<typename E::value_type>::value>* = nullptr>
auto xt::linalg::eigvalsh(const xexpression<E> &A, char UPLO = 'L')

Compute the eigenvalues of a Hermitian or real symmetric matrix xexpression.

Parameters:

Matrix – for which the eigenvalues are computed

Returns:

xtensor containing the eigenvalues.

Norms and other numbers

enum class xt::linalg::normorder

Selects special norm orders.

Values:

enumerator frob

Frobenius norm.

enumerator nuc

Nuclear norm.

enumerator inf

Positive infinity norm.

enumerator neg_inf

Negative infinity norm.

template<class E>
auto xt::linalg::norm(const xexpression<E> &vec, int ord)

Calculate norm of vector, or matrix.

Parameters:
  • vec – input vector

  • ord – order of norm. This can be any integer for a vector or [-2,-1,1,2] for a matrix.

Template Parameters:

type – of xexpression

Returns:

scalar result

template<class E>
auto xt::linalg::norm(const xexpression<E> &vec, normorder ord)

Calculate matrix or vector norm using normorder.

Parameters:
  • vec – The matrix or vector to take the norm of

  • ord – normorder (frob, nuc, inf, neg_inf)

Returns:

norm value

template<class E>
E::value_type xt::linalg::norm(const xexpression<E> &vec)

Calculate default norm (2-norm for vector, Frobenius norm for matrix)

Parameters:

vec – Input vector or matrix

Returns:

norm

Warning

doxygenfunction: Unable to resolve function “xt::linalg::det” with arguments (const xexpression<E>&) in doxygen xml output for project “xtensor-blas” from directory: ../xml. Potential matches:

- template<class T> auto det(const xexpression<T> &A)

Warning

doxygenfunction: Unable to resolve function “xt::linalg::slogdet” with arguments (const xexpression<E>&) in doxygen xml output for project “xtensor-blas” from directory: ../xml. Potential matches:

- template<class T, std::enable_if_t<xtl::is_complex<typename T::value_type>::value, int> = 0> auto slogdet(const xexpression<T> &A)
template<class T>
int xt::linalg::matrix_rank(const xexpression<T> &m, double tol = -1.0)

Calculate the matrix rank of m.

If tol == -1, the tolerance is automatically computed.

Parameters:
  • m – matrix for which rank is calculated

  • tol – tolerance for finding rank

template<class T>
auto xt::linalg::trace(const xexpression<T> &M, int offset = 0, int axis1 = 0, int axis2 = 1)

Compute the trace of a xexpression.

Solving equations and inverting matrices

template<class E1, class E2>
auto xt::linalg::solve(const xexpression<E1> &A, const xexpression<E2> &b)

Solve a linear matrix equation, or system of linear scalar equations.

Computes the “exact” solution, x, of the well-determined, i.e., full rank, linear matrix equation ax = b.

Parameters:
  • a – Coefficient matrix

  • b – Ordinate or “dependent variable” values.

Returns:

Solution to the system a x = b. Returned shape is identical to b.

template<class T, class E>
auto xt::linalg::lstsq(const xexpression<T> &A, const xexpression<E> &b, double rcond = -1.0)

Calculate the least-squares solution to a linear matrix equation.

Parameters:
  • A – coefficient matrix

  • b – Ordinate, or dependent variable values. If b is two-dimensional, the least-squares solution is calculated for each of the K columns of b.

  • rcond – Cut-off ratio for small singular values of A. For the purposes of rank determination, singular values are treated as zero if they are smaller than rcond times the largest singular value of a.

Returns:

tuple containing (x, residuals, rank, s) where: x is the least squares solution. Note that the solution is always returned as a 2D matrix where the columns are the solutions (even for a 1D b). s Sums of residuals; squared Euclidean 2-norm for each column in b - a*x. If the rank of A is < N or M <= N, this is an empty xtensor. rank the rank of A s singular values of A

template<class E1>
auto xt::linalg::inv(const xexpression<E1> &A)

Compute the (multiplicative) inverse of a matrix.

Parameters:

A – xexpression to be inverted

Returns:

(Multiplicative) inverse of the matrix a.

template<class T>
auto xt::linalg::pinv(const xexpression<T> &A, double rcond = 1e-15)

Calculate Moore-Rose pseudo inverse using LAPACK SVD.

Other

template<class E1, class E2>
auto xt::linalg::cross(const xexpression<E1> &a, const xexpression<E2> &b)

Non-broadcasting cross product between two vectors.

Calculate cross product between two 1D vectors with 2- or 3 entries. If only two entries are available, the third entry is assumed to be 0.

Parameters:
  • a – input vector

  • b – input vector

Returns:

resulting array