![]() |
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... | |
Functions to decompose a tensor into some tensors.
| 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.
| [in] | a | A tensor |
| [in] | axes_row | Axes for left singular vectors. |
| [in] | axes_col | Axes for right singular vectors. |
| [out] | w | Eigenvalues |
| 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,
calculates the following decomposition.
\[ A_{abcd} = A_{(ad)(bc)} = \sum_i U_{adi} w_i (U^\dagger)_{ibc}. \]
| [in] | a | A tensor |
| [in] | axes_row | Axes for left singular vectors. |
| [in] | axes_col | Axes for right singular vectors. |
| [out] | w | Eigenvalues |
| [out] | z | Tensor corresponds to Eigenvectors |
| int mptensor::eig | ( | const Tensor< Matrix, C > & | a, |
| std::vector< complex > & | w | ||
| ) |
Compute only the eigenvalues of a general matrix (rank-2 tensor)
| [in] | a | Rank-2 tensor |
| [out] | w | Eigenvalues |
| 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)
| [in] | a | Rank-2 tensor |
| [out] | w | Eigenvalues |
| [out] | z | Eigenvectors |
| 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,
calculates the following decomposition.
\[ A_{abcd} = A_{(ad)(bc)} = \sum_i U_{adi} w_i (U^\dagger)_{ibc}. \]
| [in] | a | A tensor |
| [in] | axes_row | Axes for left singular vectors. |
| [in] | axes_col | Axes for right singular vectors. |
| [out] | w | Eigenvalues |
| [out] | z | Tensor corresponds to Eigenvectors |
| 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.
| [in] | a | A tensor |
| [in] | axes_row | Axes for left singular vectors. |
| [in] | axes_col | Axes for right singular vectors. |
| [out] | w | Eigenvalues |
| 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
| [in] | a | A complex Hermitian or real symmetric tensor. |
| [in] | axes_row_a | Axes of tensor a to be row after matricization. |
| [in] | axes_col_a | Axes of tensor a to be column after matricization. |
| [in] | b | A complex Hermitian or real symmetric definitive positive tensor. |
| [in] | axes_row_b | Axes of tensor b to be row after matricization. |
| [in] | axes_col_b | Axes of tensor b to be column after matricization. |
| [out] | w | Eigenvalues |
| [out] | z | Tensor corresponds to Eigenvectors |
For example,
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)} \).
a and b should be square and have the same size.z is determined by axes_col_a. | 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)
| [in] | a | Rank-2 tensor |
| [out] | w | Eigenvalues |
| 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)
| [in] | a | Rank-2 tensor |
| [out] | w | Eigenvalues |
| [out] | z | Eigenvectors |
| 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.
| [in] | a | Tensor to be decomposed. |
| [in] | axes_row | Axes for left singular vectors. |
| [in] | axes_col | Axes for right singular vectors. |
| [out] | s | Singular values. |
| [in] | target_rank | The number of singular values to be calculated. |
| 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,
calculates the following decomposition.
\[ A_{abcd} \simeq \sum_{i=0}^{n-1} U_{adi} S_i (V^\dagger)_{ibc}. \]
| [in] | a | A tensor to be decomposed. |
| [in] | axes_row | Axes for left singular vectors. |
| [in] | axes_col | Axes for right singular vectors. |
| [out] | u | Tensor \( U \) corresponds to left singular vectors. |
| [out] | s | Singular values. |
| [out] | vt | Tensor \( V^\dagger \) corresponds to right singular vectors. |
| [in] | target_rank | The number of singular values to be calculated. |
| 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.
| [in] | a | Rank-2 tensor to be decomposed. |
| [out] | s | Singular values. |
| [in] | target_rank | The number of singular values to be calculated. |
| 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} \]
| [in] | a | Rank-2 tensor to be decomposed. |
| [out] | u | Tensor correspond to left singular vectors. |
| [out] | s | Singular values. |
| [out] | vt | Tensor correspond to right singular vectors. |
| [in] | target_rank | The number of singular values to be calculated. |
| 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,
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.
| [in] | a | A tensor to be decomposed. |
| [in] | axes_row | Axes for the orthogonal matrix. |
| [in] | axes_col | Axes for the upper triangular matrix. |
| [out] | q | Decomposed tensor corresponds to orthogonal matrix |
| [out] | r | Decomposed tensor corresponds to upper triangular matrix |
| 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 \).
| [in] | a | A matrix (rank-2 tensor) to be decomposed |
| [out] | q | Orthogonal matrix |
| [out] | r | Upper triangular matrix |
| 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.
| [in] | a | Tensor to be decomposed. |
| [in] | axes_row | Axes for left singular vectors. |
| [in] | axes_col | Axes for right singular vectors. |
| [out] | s | Singular values. |
| 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,
calculates the following decomposition.
\[ A_{abcd} = \sum_i U_{adi} S_i (V^\dagger)_{ibc}. \]
| [in] | a | A tensor to be decomposed. |
| [in] | axes_row | Axes for left singular vectors. |
| [in] | axes_col | Axes for right singular vectors. |
| [out] | u | Tensor \( U \) corresponds to left singular vectors. |
| [out] | s | Singular values. |
| [out] | vt | Tensor \( V^\dagger \) corresponds to right singular vectors. |
| 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.
| [in] | a | Rank-2 tensor to be decomposed. |
| [out] | s | Singular values. |
| 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} \]
| [in] | a | Rank-2 tensor to be decomposed. |
| [out] | u | Tensor correspond to left singular vectors. |
| [out] | s | Singular values. |
| [out] | vt | Tensor correspond to right singular vectors. |