mptensor  v0.3.0
Parallel Library for Tensor Network Methods
template<template< typename > class Matrix, typename C >
int mptensor::svd (const Tensor< Matrix, C > &a, std::vector< double > &s)
 Singular value decomposition for rank-2 tensor (matrix) More...
 
template<template< typename > class Matrix, typename C >
int mptensor::svd (const Tensor< Matrix, C > &a, Tensor< Matrix, C > &u, std::vector< double > &s, Tensor< Matrix, C > &vt)
 Singular value decomposition for rank-2 tensor (matrix) More...
 
template<template< typename > class Matrix, typename C >
int mptensor::svd (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, std::vector< double > &s)
 Singular value decomposition for tensor. More...
 
template<template< typename > class Matrix, typename C >
int mptensor::svd (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, Tensor< Matrix, C > &u, std::vector< double > &s, Tensor< Matrix, C > &vt)
 Singular value decomposition for tensor. More...
 
template<template< typename > class Matrix, typename C >
int mptensor::psvd (const Tensor< Matrix, C > &a, std::vector< double > &s, const size_t target_rank)
 Partial SVD for rank-2 tensor (matrix) More...
 
template<template< typename > class Matrix, typename C >
int mptensor::psvd (const Tensor< Matrix, C > &a, Tensor< Matrix, C > &u, std::vector< double > &s, Tensor< Matrix, C > &vt, const size_t target_rank)
 Partial SVD for rank-2 tensor (matrix) More...
 
template<template< typename > class Matrix, typename C >
int mptensor::psvd (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, std::vector< double > &s, const size_t target_rank)
 Partial SVD for tensor. More...
 
template<template< typename > class Matrix, typename C >
int mptensor::psvd (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, Tensor< Matrix, C > &u, std::vector< double > &s, Tensor< Matrix, C > &vt, const size_t target_rank)
 Partial SVD for tensor. More...
 
template<template< typename > class Matrix, typename C >
int mptensor::qr (const Tensor< Matrix, C > &a, Tensor< Matrix, C > &q, Tensor< Matrix, C > &r)
 QR decomposition of matrix (rank-2 tensor) More...
 
template<template< typename > class Matrix, typename C >
int mptensor::qr (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, Tensor< Matrix, C > &q, Tensor< Matrix, C > &r)
 QR decomposition of tensor. More...
 
template<template< typename > class Matrix, typename C >
int mptensor::eigh (const Tensor< Matrix, C > &a, std::vector< double > &eigval, Tensor< Matrix, C > &eigvec)
 
template<template< typename > class Matrix, typename C >
int mptensor::eigh (const Tensor< Matrix, C > &a, std::vector< double > &eigval)
 
template<template< typename > class Matrix, typename C >
int mptensor::eigh (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, std::vector< double > &eigval, Tensor< Matrix, C > &eigvec)
 
template<template< typename > class Matrix, typename C >
int mptensor::eigh (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, std::vector< double > &w)
 Compute the eigenvalues of a complex Hermitian or real symmetric tensor. More...
 
template<template< typename > class Matrix, typename C >
int mptensor::eigh (const Tensor< Matrix, C > &a, const Axes &axes_row_a, const Axes &axes_col_a, const Tensor< Matrix, C > &b, const Axes &axes_row_b, const Axes &axes_col_b, std::vector< double > &eigval, Tensor< Matrix, C > &eigvec)
 
template<template< typename > class Matrix, typename C >
int mptensor::eig (const Tensor< Matrix, C > &a, std::vector< complex > &w, Tensor< Matrix, complex > &z)
 Compute the eigenvalues and eigenvectors of a general matrix (rank-2 tensor) More...
 
template<template< typename > class Matrix, typename C >
int mptensor::eig (const Tensor< Matrix, C > &a, std::vector< complex > &w)
 Compute only the eigenvalues of a general matrix (rank-2 tensor) More...
 
template<template< typename > class Matrix, typename C >
int mptensor::eig (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, std::vector< complex > &w, Tensor< Matrix, complex > &z)
 Compute the eigenvalues and eigenvectors of a general tensor. More...
 
template<template< typename > class Matrix, typename C >
int mptensor::eig (const Tensor< Matrix, C > &a, const Axes &axes_row, const Axes &axes_col, std::vector< complex > &w)
 Compute the eigenvalues of a general tensor. More...
 

Detailed Description

Functions to decompose a tensor into some tensors.

Function Documentation

◆ eig() [1/4]

template<template< typename > class Matrix, typename C >
int mptensor::eig ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
std::vector< complex > &  w 
)

Compute the eigenvalues of a general tensor.

Parameters
[in]aA tensor
[in]axes_rowAxes for left singular vectors.
[in]axes_colAxes for right singular vectors.
[out]wEigenvalues
Returns
Information from linear-algebra library.
Warning
This function may require huge memory because it uses LAPACK routines.

◆ eig() [2/4]

template<template< typename > class Matrix, typename C >
int mptensor::eig ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
std::vector< complex > &  w,
Tensor< Matrix, complex > &  z 
)

Compute the eigenvalues and eigenvectors of a general tensor.

For example,

eig(A, Axes(0,3), Axes(1,2), w, U)
int eig(const Tensor< Matrix, C > &a, std::vector< complex > &eigval, Tensor< Matrix, complex > &eigvec)
Compute the eigenvalues and eigenvectors of a general matrix (rank-2 tensor)
Definition: tensor_impl.hpp:2295
Index Axes
Definition: tensor.hpp:45

calculates the following decomposition.

\[ A_{abcd} = A_{(ad)(bc)} = \sum_i U_{adi} w_i (U^\dagger)_{ibc}. \]

Parameters
[in]aA tensor
[in]axes_rowAxes for left singular vectors.
[in]axes_colAxes for right singular vectors.
[out]wEigenvalues
[out]zTensor corresponds to Eigenvectors
Returns
Information from linear-algebra library.
Warning
This function may require huge memory because it uses LAPACK routines.

◆ eig() [3/4]

template<template< typename > class Matrix, typename C >
int mptensor::eig ( const Tensor< Matrix, C > &  a,
std::vector< complex > &  w 
)

Compute only the eigenvalues of a general matrix (rank-2 tensor)

Parameters
[in]aRank-2 tensor
[out]wEigenvalues
Returns
Information from linear-algebra library.
Warning
This function may require huge memory because it uses LAPACK routines.

◆ eig() [4/4]

template<template< typename > class Matrix, typename C >
int mptensor::eig ( const Tensor< Matrix, C > &  a,
std::vector< complex > &  w,
Tensor< Matrix, complex > &  z 
)

Compute the eigenvalues and eigenvectors of a general matrix (rank-2 tensor)

Parameters
[in]aRank-2 tensor
[out]wEigenvalues
[out]zEigenvectors
Returns
Information from linear-algebra library.
Warning
This function may require huge memory because it uses LAPACK routines.

◆ eigh() [1/5]

template<template< typename > class Matrix, typename C >
int mptensor::eigh ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
std::vector< double > &  w,
Tensor< Matrix, C > &  z 
)

Compute the eigenvalues and eigenvectors of a complex Hermitian or real symmetric tensor

For example,

eigh(A, Axes(0,3), Axes(1,2), w, U)
int eigh(const Tensor< Matrix, C > &a, std::vector< double > &eigval, Tensor< Matrix, C > &eigvec)
Definition: tensor_impl.hpp:2071

calculates the following decomposition.

\[ A_{abcd} = A_{(ad)(bc)} = \sum_i U_{adi} w_i (U^\dagger)_{ibc}. \]

Parameters
[in]aA tensor
[in]axes_rowAxes for left singular vectors.
[in]axes_colAxes for right singular vectors.
[out]wEigenvalues
[out]zTensor corresponds to Eigenvectors
Returns
Information from linear-algebra library.
Warning
Input tensor should be Hermite.

◆ eigh() [2/5]

template<template< typename > class Matrix, typename C >
int mptensor::eigh ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
std::vector< double > &  w 
)

Compute the eigenvalues of a complex Hermitian or real symmetric tensor.

Parameters
[in]aA tensor
[in]axes_rowAxes for left singular vectors.
[in]axes_colAxes for right singular vectors.
[out]wEigenvalues
Returns
Information from linear-algebra library.
Warning
Input tensor should be Hermite.

◆ eigh() [3/5]

template<template< typename > class Matrix, typename C >
int mptensor::eigh ( const Tensor< Matrix, C > &  a,
const Axes axes_row_a,
const Axes axes_col_a,
const Tensor< Matrix, C > &  b,
const Axes axes_row_b,
const Axes axes_col_b,
std::vector< double > &  w,
Tensor< Matrix, C > &  z 
)

Solve a generalized eigenvalue problem \( A \vec{v} = \lambda B \vec{v} \) for a complex Hermitian or real symmetric tensor

Parameters
[in]aA complex Hermitian or real symmetric tensor.
[in]axes_row_aAxes of tensor a to be row after matricization.
[in]axes_col_aAxes of tensor a to be column after matricization.
[in]bA complex Hermitian or real symmetric definitive positive tensor.
[in]axes_row_bAxes of tensor b to be row after matricization.
[in]axes_col_bAxes of tensor b to be column after matricization.
[out]wEigenvalues
[out]zTensor corresponds to Eigenvectors
Returns
Information from linear-algebra library.

For example,

eigh(A, Axes(0,3), Axes(1,2), B, Axes(1, 3), Axes(2, 0), w, U)

solve a generalized eigenvalue problem

\[ \sum_{bc} A_{abcd} U_{bck} = \sum_{bc} B_{cabd} U_{bck} w_k. \]

In the matricized representation, it is written as

\[ \sum_{bc} A'_{(ad)(bc)} U'_{(bc)(k)} = \sum_{bc} B_{(ad)(bc)} U'_{(bc)(k)} w_k, \]

where \( A \) and \( B \) are matricized as \( A_{abcd} = A'_{(ad)(bc)} \) and \( B_{cabd} = B'_{(ad)(bc)} \).

Warning
The matricies obtained by a and b should be square and have the same size.
Note
The shape of z is determined by axes_col_a.

◆ eigh() [4/5]

template<template< typename > class Matrix, typename C >
int mptensor::eigh ( const Tensor< Matrix, C > &  a,
std::vector< double > &  w 
)

Compute only the eigenvalues of a complex Hermitian or real symmetric matrix (rank-2 tensor)

Parameters
[in]aRank-2 tensor
[out]wEigenvalues
Returns
Information from linear-algebra library.

◆ eigh() [5/5]

template<template< typename > class Matrix, typename C >
int mptensor::eigh ( const Tensor< Matrix, C > &  a,
std::vector< double > &  w,
Tensor< Matrix, C > &  z 
)

Compute the eigenvalues and eigenvectors of a complex Hermitian or real symmetric matrix (rank-2 tensor)

Parameters
[in]aRank-2 tensor
[out]wEigenvalues
[out]zEigenvectors
Returns
Information from linear-algebra library.

◆ psvd() [1/4]

template<template< typename > class Matrix, typename C >
int mptensor::psvd ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
std::vector< double > &  s,
const size_t  target_rank 
)

Partial SVD for tensor.

This function only calculates singular values.

Parameters
[in]aTensor to be decomposed.
[in]axes_rowAxes for left singular vectors.
[in]axes_colAxes for right singular vectors.
[out]sSingular values.
[in]target_rankThe number of singular values to be calculated.
Returns
Information from linear-algebra library.
Attention
Computational cost is the same as full SVD.

◆ psvd() [2/4]

template<template< typename > class Matrix, typename C >
int mptensor::psvd ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
Tensor< Matrix, C > &  u,
std::vector< double > &  s,
Tensor< Matrix, C > &  vt,
const size_t  target_rank 
)

Partial SVD for tensor.

This may be useful for creation of an isometry. For example,

psvd(A, Axes(0,3), Axes(1,2), U, S, VT, n)
int psvd(const Tensor< Matrix, C > &a, std::vector< double > &s, const size_t target_rank)
Partial SVD for rank-2 tensor (matrix)
Definition: tensor_impl.hpp:1866

calculates the following decomposition.

\[ A_{abcd} \simeq \sum_{i=0}^{n-1} U_{adi} S_i (V^\dagger)_{ibc}. \]

Parameters
[in]aA tensor to be decomposed.
[in]axes_rowAxes for left singular vectors.
[in]axes_colAxes for right singular vectors.
[out]uTensor \( U \) corresponds to left singular vectors.
[out]sSingular values.
[out]vtTensor \( V^\dagger \) corresponds to right singular vectors.
[in]target_rankThe number of singular values to be calculated.
Returns
Information from linear-algebra library.
Attention
Computational cost is the same as full SVD.

◆ psvd() [3/4]

template<template< typename > class Matrix, typename C >
int mptensor::psvd ( const Tensor< Matrix, C > &  a,
std::vector< double > &  s,
const size_t  target_rank 
)

Partial SVD for rank-2 tensor (matrix)

This function only caluculates singular values.

Parameters
[in]aRank-2 tensor to be decomposed.
[out]sSingular values.
[in]target_rankThe number of singular values to be calculated.
Returns
Information from linear-algebra library.
Attention
Computational cost is the same as full SVD.

◆ psvd() [4/4]

template<template< typename > class Matrix, typename C >
int mptensor::psvd ( const Tensor< Matrix, C > &  a,
Tensor< Matrix, C > &  u,
std::vector< double > &  s,
Tensor< Matrix, C > &  vt,
const size_t  target_rank 
)

Partial SVD for rank-2 tensor (matrix)

\[ a_{ij} = \sum_k u_{ik} s_k (v^{\dagger})_{kj} \]

Parameters
[in]aRank-2 tensor to be decomposed.
[out]uTensor correspond to left singular vectors.
[out]sSingular values.
[out]vtTensor correspond to right singular vectors.
[in]target_rankThe number of singular values to be calculated.
Returns
Information from linear-algebra library.
Attention
Computational cost is the same as full SVD.

◆ qr() [1/2]

template<template< typename > class Matrix, typename C >
int mptensor::qr ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
Tensor< Matrix, C > &  q,
Tensor< Matrix, C > &  r 
)

QR decomposition of tensor.

This may be useful for creation of an isometry. For example,

qr(A, Axes(0,3), Axes(1,2), Q, R)
int qr(const Tensor< Matrix, C > &a, Tensor< Matrix, C > &q, Tensor< Matrix, C > &r)
QR decomposition of matrix (rank-2 tensor)
Definition: tensor_impl.hpp:1974

calculates the following decomposition,

\[ A_{abcd} = \sum_i Q_{adi} R_{ibc}, \]

where the tensor Q satisfies

\[ \sum_{ab} Q_{abi} Q_{abj}^{*} = \delta_{ij}, \]

and the tensor R corresponds to the upper triangular matrix.

Parameters
[in]aA tensor to be decomposed.
[in]axes_rowAxes for the orthogonal matrix.
[in]axes_colAxes for the upper triangular matrix.
[out]qDecomposed tensor corresponds to orthogonal matrix
[out]rDecomposed tensor corresponds to upper triangular matrix
Returns
Information from linear-algebra library.

◆ qr() [2/2]

template<template< typename > class Matrix, typename C >
int mptensor::qr ( const Tensor< Matrix, C > &  a,
Tensor< Matrix, C > &  q,
Tensor< Matrix, C > &  r 
)

QR decomposition of matrix (rank-2 tensor)

The sizes of q and r are reduced. When a is \( m\times n \) matrix and \( k=\min\{m,n\} \), the dimensions of q and r is \( m\times k \) and \( k\times n \).

Parameters
[in]aA matrix (rank-2 tensor) to be decomposed
[out]qOrthogonal matrix
[out]rUpper triangular matrix
Returns
Information from linear-algebra library.

◆ svd() [1/4]

template<template< typename > class Matrix, typename C >
int mptensor::svd ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
std::vector< double > &  s 
)

Singular value decomposition for tensor.

This function only calculates singular values.

Parameters
[in]aTensor to be decomposed.
[in]axes_rowAxes for left singular vectors.
[in]axes_colAxes for right singular vectors.
[out]sSingular values.
Returns
Information from linear-algebra library.

◆ svd() [2/4]

template<template< typename > class Matrix, typename C >
int mptensor::svd ( const Tensor< Matrix, C > &  a,
const Axes axes_row,
const Axes axes_col,
Tensor< Matrix, C > &  u,
std::vector< double > &  s,
Tensor< Matrix, C > &  vt 
)

Singular value decomposition for tensor.

This may be useful for creation of an isometry. For example,

svd(A, Axes(0,3), Axes(1,2), U, S, VT)
int svd(const Tensor< Matrix, C > &a, std::vector< double > &s)
Singular value decomposition for rank-2 tensor (matrix)
Definition: tensor_impl.hpp:1722

calculates the following decomposition.

\[ A_{abcd} = \sum_i U_{adi} S_i (V^\dagger)_{ibc}. \]

Parameters
[in]aA tensor to be decomposed.
[in]axes_rowAxes for left singular vectors.
[in]axes_colAxes for right singular vectors.
[out]uTensor \( U \) corresponds to left singular vectors.
[out]sSingular values.
[out]vtTensor \( V^\dagger \) corresponds to right singular vectors.
Returns
Information from linear-algebra library.

◆ svd() [3/4]

template<template< typename > class Matrix, typename C >
int mptensor::svd ( const Tensor< Matrix, C > &  a,
std::vector< double > &  s 
)

Singular value decomposition for rank-2 tensor (matrix)

This function only caluculates singular values.

Parameters
[in]aRank-2 tensor to be decomposed.
[out]sSingular values.
Returns
Information from linear-algebra library.

◆ svd() [4/4]

template<template< typename > class Matrix, typename C >
int mptensor::svd ( const Tensor< Matrix, C > &  a,
Tensor< Matrix, C > &  u,
std::vector< double > &  s,
Tensor< Matrix, C > &  vt 
)

Singular value decomposition for rank-2 tensor (matrix)

\[ a_{ij} = \sum_k u_{ik} s_k (v^{\dagger})_{kj} \]

Parameters
[in]aRank-2 tensor to be decomposed.
[out]uTensor correspond to left singular vectors.
[out]sSingular values.
[out]vtTensor correspond to right singular vectors.
Returns
Information from linear-algebra library.