|
Cytnx v1.0.0
|
linear algebra related functions. More...
Functions | |
| cytnx::UniTensor | Add (const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt) |
| The addition function between two UniTensor. | |
| template<class T > | |
| cytnx::UniTensor | Add (const T &lc, const cytnx::UniTensor &Rt) |
| The addition function between a template type and a UniTensor. | |
| template<class T > | |
| cytnx::UniTensor | Add (const cytnx::UniTensor &Lt, const T &rc) |
| The addition function between a UniTensor and a template type. | |
| cytnx::UniTensor | Sub (const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt) |
| The subtraction function between two UniTensor. | |
| template<class T > | |
| cytnx::UniTensor | Sub (const T &lc, const cytnx::UniTensor &Rt) |
| The subtraction function between a UniTensor and a template type. | |
| template<class T > | |
| cytnx::UniTensor | Sub (const cytnx::UniTensor &Lt, const T &rc) |
| The subtraction function between a UniTensor and a template type. | |
| cytnx::UniTensor | Mul (const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt) |
| The multiplication function between two UniTensor. | |
| template<class T > | |
| cytnx::UniTensor | Mul (const T &lc, const cytnx::UniTensor &Rt) |
| The multiplication function between a UniTensor and a template type. | |
| template<class T > | |
| cytnx::UniTensor | Mul (const cytnx::UniTensor &Lt, const T &rc) |
| The multiplication function between a UniTensor and a template type. | |
| cytnx::UniTensor | Div (const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt) |
| The division function between two UniTensor. | |
| template<class T > | |
| cytnx::UniTensor | Div (const T &lc, const cytnx::UniTensor &Rt) |
| The division function between a UniTensor and a template type. | |
| template<class T > | |
| cytnx::UniTensor | Div (const cytnx::UniTensor &Lt, const T &rc) |
| The division function between a UniTensor and a template type. | |
| cytnx::UniTensor | Mod (const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt) |
| element-wise modulo | |
| template<class T > | |
| cytnx::UniTensor | Mod (const T &lc, const cytnx::UniTensor &Rt) |
| The modulo function between a UniTensor and a template type. | |
| template<class T > | |
| cytnx::UniTensor | Mod (const cytnx::UniTensor &Lt, const T &rc) |
| The modulo function between a UniTensor and a template type. | |
| std::vector< cytnx::UniTensor > | Svd (const cytnx::UniTensor &Tin, const bool &is_UvT=true) |
| Perform Singular-Value decomposition on a UniTensor using divide-and-conquer method. | |
| std::vector< cytnx::UniTensor > | Rsvd (const cytnx::UniTensor &Tin, cytnx_uint64 keepdim, bool is_U=true, bool is_vT=true, cytnx_uint64 power_iteration=2, unsigned int seed=random::__static_random_device()) |
| Perform randomized Singular-Value decomposition on a UniTensor using the ?gesvd method. | |
| std::vector< cytnx::UniTensor > | Gesvd (const cytnx::UniTensor &Tin, const bool &is_U=true, const bool &is_vT=true) |
| Perform Singular-Value decomposition on a UniTensor using ?gesvd method. | |
| std::vector< cytnx::UniTensor > | Svd_truncate (const cytnx::UniTensor &Tin, const cytnx_uint64 &keepdim, const double &err=0., const bool &is_UvT=true, const unsigned int &return_err=0, const cytnx_uint64 &mindim=1) |
| Perform Singular-Value decomposition on a UniTensor with truncation. | |
| std::vector< cytnx::UniTensor > | Svd_truncate (const cytnx::UniTensor &Tin, const cytnx_uint64 &keepdim, std::vector< cytnx_uint64 > min_blockdim, const double &err=0., const bool &is_UvT=true, const unsigned int &return_err=0, const cytnx_uint64 &mindim=1) |
Perform Singular-Value decomposition on a UniTensor with truncation and keep at most min_blockdim singular values in each block. | |
| std::vector< cytnx::UniTensor > | Gesvd_truncate (const cytnx::UniTensor &Tin, const cytnx_uint64 &keepdim, const double &err=0., const bool &is_U=true, const bool &is_vT=true, const unsigned int &return_err=0, const cytnx_uint64 &mindim=1) |
| Perform Singular-Value decomposition on a UniTensor with truncation. | |
| std::vector< cytnx::UniTensor > | Gesvd_truncate (const cytnx::UniTensor &Tin, const cytnx_uint64 &keepdim, std::vector< cytnx_uint64 > min_blockdim, const double &err=0., const bool &is_U=true, const bool &is_vT=true, const unsigned int &return_err=0, const cytnx_uint64 &mindim=1) |
Perform Singular-Value decomposition on a UniTensor with truncation and keep at most min_blockdim singular values in each block. | |
| std::vector< cytnx::UniTensor > | Rsvd_truncate (const cytnx::UniTensor &Tin, cytnx_uint64 keepdim, double err=0., bool is_U=true, bool is_vT=true, unsigned int return_err=0, cytnx_uint64 mindim=1, cytnx_uint64 oversampling_summand=10, double oversampling_factor=1., cytnx_uint64 power_iteration=0, unsigned int seed=random::__static_random_device()) |
| Perform a truncated Singular-Value decomposition of a UniTensor. | |
| std::vector< cytnx::UniTensor > | Hosvd (const cytnx::UniTensor &Tin, const std::vector< cytnx_uint64 > &mode, const bool &is_core=true, const bool &is_Ls=false, const std::vector< cytnx_int64 > &trucate_dim=std::vector< cytnx_int64 >()) |
| template<typename T > | |
| cytnx::UniTensor | ExpH (const cytnx::UniTensor &Tin, const T &a, const T &b=0) |
| Perform the exponential function on a UniTensor, which the blocks are Hermitian matrix. | |
| template<typename T > | |
| cytnx::UniTensor | ExpM (const cytnx::UniTensor &Tin, const T &a, const T &b=0) |
| Perform the exponential function on a UniTensor. | |
| cytnx::UniTensor | ExpH (const cytnx::UniTensor &Tin) |
| Perform the exponential function on a UniTensor, which the blocks are Hermitian matrix. | |
| cytnx::UniTensor | ExpM (const cytnx::UniTensor &Tin) |
| Perform the exponential function on a UniTensor. | |
| cytnx::UniTensor | Trace (const cytnx::UniTensor &Tin, const cytnx_int64 &a=0, const cytnx_int64 &b=1) |
| cytnx::UniTensor | Trace (const cytnx::UniTensor &Tin, const std::string &a, const std::string &b) |
| Perform trace over two legs of a UniTensor. | |
| std::vector< cytnx::UniTensor > | Qr (const cytnx::UniTensor &Tin, const bool &is_tau=false) |
| Perform the QR decomposition on a UniTensor. | |
| std::vector< cytnx::UniTensor > | Qdr (const cytnx::UniTensor &Tin, const bool &is_tau=false) |
| Perform the QDR decomposition on a UniTensor. | |
| UniTensor | Pow (const cytnx::UniTensor &Tin, const double &p) |
| take power p on all the elements in UniTensor. | |
| void | Pow_ (UniTensor &Tin, const double &p) |
| Take power p on all the elements in UniTensor, inplacely. | |
| cytnx::UniTensor | Conj (const cytnx::UniTensor &UT) |
| Elementwise conjugate of the UniTensor. | |
| void | Conj_ (cytnx::UniTensor &UT) |
| Inplace elementwise conjugate of the UniTensor. | |
| Tensor | Add (const Tensor &Lt, const Tensor &Rt) |
| The addition function for Tensor. | |
| template<class T > | |
| Tensor | Add (const T &lc, const Tensor &Rt) |
| The addition function for Tensor. | |
| template<class T > | |
| Tensor | Add (const Tensor &Lt, const T &rc) |
| The addition function for Tensor. | |
| void | iAdd (Tensor &Lt, const Tensor &Rt) |
| The addition function for Tensor, inplacely. | |
| Tensor | Sub (const Tensor &Lt, const Tensor &Rt) |
| The subtraction function for Tensor. | |
| template<class T > | |
| Tensor | Sub (const T &lc, const Tensor &Rt) |
| The subtraction function for Tensor. | |
| template<class T > | |
| Tensor | Sub (const Tensor &Lt, const T &rc) |
| The subtraction function for Tensor. | |
| void | iSub (Tensor &Lt, const Tensor &Rt) |
| The subtraction function for Tensot, inplscely. | |
| Tensor | Mul (const Tensor &Lt, const Tensor &Rt) |
| The multiplication function for Tensor. | |
| template<class T > | |
| Tensor | Mul (const T &lc, const Tensor &Rt) |
| The multiplication function for Tensor. | |
| template<class T > | |
| Tensor | Mul (const Tensor &Lt, const T &rc) |
| The multiplication function for Tensor. | |
| void | iMul (Tensor &Lt, const Tensor &Rt) |
| The multiplication function for Tensor, inplacely. | |
| Tensor | Div (const Tensor &Lt, const Tensor &Rt) |
| The division function for Tensor. | |
| template<class T > | |
| Tensor | Div (const T &lc, const Tensor &Rt) |
| The division function for Tensor. | |
| template<class T > | |
| Tensor | Div (const Tensor &Lt, const T &rc) |
| The division function for Tensor. | |
| void | iDiv (Tensor &Lt, const Tensor &Rt) |
| The inplace division function for Tensor, inplacely. | |
| Tensor | Mod (const Tensor &Lt, const Tensor &Rt) |
| The mod function for Tensor. | |
| template<class T > | |
| Tensor | Mod (const T &lc, const Tensor &Rt) |
| The mod function for Tensor. | |
| template<class T > | |
| Tensor | Mod (const Tensor &Lt, const T &rc) |
| The mod function for Tensor. | |
| Tensor | Cpr (const Tensor &Lt, const Tensor &Rt) |
| The comparison function for Tensor. | |
| template<class T > | |
| Tensor | Cpr (const T &lc, const Tensor &Rt) |
| The comparison function for Tensor. | |
| template<class T > | |
| Tensor | Cpr (const Tensor &Lt, const T &rc) |
| The comparison function for Tensor. | |
| Tensor | Norm (const Tensor &Tl) |
| Calculate the norm of a tensor. | |
| Tensor | Norm (const cytnx::UniTensor &uTl) |
| Calculate the norm of an UniTensor. | |
| Tensor | Det (const Tensor &Tl) |
| Calculate the determinant of a tensor. | |
| Tensor | Rand_isometry (const Tensor &Tin, const cytnx_uint64 &keepdim, const cytnx_uint64 &power_iteration=2, const unsigned int &seed=random::__static_random_device()) |
| Generate an isometrized left isometry for a rank-2 Tensor (a matrix), such that Q * Qdag * Tin approximates Tin (Qdag is the hermitean conjugate of Q) | |
| std::vector< Tensor > | Svd (const Tensor &Tin, const bool &is_UvT=true) |
| Perform Singular-Value decomposition on a rank-2 Tensor (a matrix). | |
| std::vector< Tensor > | Gesvd (const Tensor &Tin, const bool &is_U=true, const bool &is_vT=true) |
| Perform Singular-Value decomposition on a rank-2 Tensor (a matrix). | |
| std::vector< Tensor > | Rsvd (const Tensor &Tin, cytnx_uint64 keepdim, bool is_U=true, bool is_vT=true, cytnx_uint64 power_iteration=2, unsigned int seed=random::__static_random_device()) |
| Perform randomized Singular-Value decomposition on a rank-2 Tensor (a matrix) using the ?gesvd method. | |
| std::vector< Tensor > | Svd_truncate (const Tensor &Tin, const cytnx_uint64 &keepdim, const double &err=0., const bool &is_UvT=true, const unsigned int &return_err=0, const cytnx_uint64 &mindim=1) |
| Perform a truncated Singular-Value decomposition of a rank-2 Tensor (a matrix). | |
| std::vector< Tensor > | Rsvd_truncate (const Tensor &Tin, cytnx_uint64 keepdim, double err=0., bool is_U=true, bool is_vT=true, unsigned int return_err=0, cytnx_uint64 mindim=1, cytnx_uint64 oversampling_summand=10, double oversampling_factor=1., cytnx_uint64 power_iteration=0, unsigned int seed=random::__static_random_device()) |
| Perform a truncated Singular-Value decomposition of a rank-2 Tensor (a matrix). | |
| std::vector< Tensor > | Gesvd_truncate (const Tensor &Tin, const cytnx_uint64 &keepdim, const double &err=0., const bool &is_U=true, const bool &is_vT=true, const unsigned int &return_err=0, const cytnx_uint64 &mindim=1) |
| Perform a truncated Singular-Value decomposition of a rank-2 Tensor (a matrix) by a randomized SVD. | |
| std::vector< Tensor > | Hosvd (const Tensor &Tin, const std::vector< cytnx_uint64 > &mode, const bool &is_core=true, const bool &is_Ls=false, const std::vector< cytnx_int64 > &trucate_dim=std::vector< cytnx_int64 >()) |
| std::vector< Tensor > | Qr (const Tensor &Tin, const bool &is_tau=false) |
| Perform QR decomposition on a rank-2 Tensor. | |
| std::vector< Tensor > | Qdr (const Tensor &Tin, const bool &is_tau=false) |
| Perform QDR decomposition on a rank-2 Tensor. | |
| std::vector< Tensor > | Eigh (const Tensor &Tin, const bool &is_V=true, const bool &row_v=false) |
| eigen-value decomposition for Hermitian matrix | |
| std::vector< UniTensor > | Eigh (const cytnx::UniTensor &Tin, const bool &is_V=true, const bool &row_v=false) |
| std::vector< Tensor > | Eig (const Tensor &Tin, const bool &is_V=true, const bool &row_v=false) |
| eigen-value decomposition for generic square matrix | |
| std::vector< UniTensor > | Eig (const cytnx::UniTensor &Tin, const bool &is_V=true, const bool &row_v=false) |
| Tensor | Trace (const Tensor &Tn, const cytnx_uint64 &axisA=0, const cytnx_uint64 &axisB=1) |
| perform trace over index. | |
| Tensor | Min (const Tensor &Tn) |
| get the minimum element. | |
| Tensor | Max (const Tensor &Tn) |
| get the maximum element. | |
| Tensor | Sum (const Tensor &Tn) |
| get the sum of all the elements. | |
| Tensor | Matmul (const Tensor &TL, const Tensor &TR) |
| perform matrix multiplication on two tensors. | |
| Tensor | Matmul_dg (const Tensor &Tl, const Tensor &Tr) |
| perform matrix multiplication on two Tensors with one rank-1 and the other rank-2 where the rank-1 represent the diagonal elements of the specific tensor. | |
| Tensor | InvM (const Tensor &Tin) |
| Matrix inverse. | |
| UniTensor | InvM (const cytnx::UniTensor &Tin) |
| void | InvM_ (Tensor &Tin) |
| inplace matrix inverse. | |
| void | InvM_ (UniTensor &Tin) |
| Tensor | Inv (const Tensor &Tin, const double &clip) |
| Element-wise inverse with clip. | |
| void | Inv_ (Tensor &Tin, const double &clip) |
| inplace perform Element-wise inverse with clip. | |
| Tensor | Conj (const Tensor &Tin) |
| Conjugate all the element in Tensor. | |
| void | Conj_ (Tensor &Tin) |
| inplace perform Conjugate on all the element in Tensor. | |
| Tensor | Exp (const Tensor &Tin) |
| Exponential all the element in Tensor. | |
| Tensor | Expf (const Tensor &Tin) |
| Exponential all the element in Tensor. | |
| void | Exp_ (Tensor &Tin) |
| inplace perform Exponential on all the element in Tensor. | |
| void | Expf_ (Tensor &Tin) |
| inplace perform Exponential on all the element in Tensor. | |
| Tensor | Pow (const Tensor &Tin, const double &p) |
| take power p on all the elements in Tensor. | |
| void | Pow_ (Tensor &Tin, const double &p) |
| inplace perform power on all the elements in Tensor. | |
| Tensor | Abs (const Tensor &Tin) |
| Elementwise absolute value. | |
| void | Abs_ (Tensor &Tin) |
inplace perform elementwiase absolute value. @This is just a inplace version of Abs. The input Tensor Tin will be modified. | |
| Tensor | Diag (const Tensor &Tin) |
| return a diagonal tensor with diagonal elements provided as Tin. | |
| Tensor | Tensordot (const Tensor &Tl, const Tensor &Tr, const std::vector< cytnx_uint64 > &idxl, const std::vector< cytnx_uint64 > &idxr, const bool &cacheL=false, const bool &cacheR=false) |
| perform tensor dot by sum out the indices assigned of two Tensors. | |
| Tensor | Tensordot_dg (const Tensor &Tl, const Tensor &Tr, const std::vector< cytnx_uint64 > &idxl, const std::vector< cytnx_uint64 > &idxr, const bool &diag_L) |
| perform tensor dot by sum out the indices assigned of two Tensors, with either one of them to be a rank-2 diagonal tensor represented by a rank-2 tensor. | |
| Tensor | Outer (const Tensor &Tl, const Tensor &Tr) |
| perform outer produces of two rank-1 Tensor. | |
| Tensor | Kron (const Tensor &Tl, const Tensor &Tr, const bool &Tl_pad_left=false, const bool &Tr_pad_left=false) |
| perform kronecker produces of two Tensor. | |
| Tensor | Directsum (const Tensor &T1, const Tensor &T2, const std::vector< cytnx_uint64 > &shared_axes) |
| perform directsum of two Tensor. | |
| Tensor | Vectordot (const Tensor &Tl, const Tensor &Tr, const bool &is_conj=false) |
| perform inner product of vectors | |
| Tensor | Dot (const Tensor &Tl, const Tensor &Tr) |
| dot product of two arrays. | |
| std::vector< Tensor > | Tridiag (const Tensor &Diag, const Tensor &Sub_diag, const bool &is_V=true, const bool &is_row=false, bool throw_excp=false) |
| perform diagonalization of symmetric tri-diagnoal matrix. | |
| template<typename T > | |
| Tensor | ExpH (const Tensor &in, const T &a, const T &b=0) |
| perform matrix exponential for Hermitian matrix | |
| Tensor | ExpH (const Tensor &in) |
| perform matrix exponential for Hermitian matrix | |
| template<typename T > | |
| Tensor | ExpM (const Tensor &in, const T &a, const T &b=0) |
| perform matrix exponential for generic matrix | |
| Tensor | ExpM (const Tensor &in) |
| perform matrix exponential for generic matrix | |
| std::vector< Tensor > | Arnoldi (LinOp *Hop, const Tensor &Tin=Tensor(), const std::string which="LM", const cytnx_uint64 &maxiter=10000, const cytnx_double &cvg_crit=0, const cytnx_uint64 &k=1, const bool &is_V=true, const cytnx_int32 &ncv=0, const bool &verbose=false) |
| Performs Arnoldi iteration for matrices or linear functions. | |
| std::vector< UniTensor > | Arnoldi (LinOp *Hop, const cytnx::UniTensor &Tin, const std::string which="LM", const cytnx_uint64 &maxiter=10000, const cytnx_double &cvg_crit=1.0e-9, const cytnx_uint64 &k=1, const bool &is_V=true, const cytnx_int32 &ncv=0, const bool &verbose=false) |
| Performs Arnoldi iteration for matrices or linear functions. | |
| std::vector< Tensor > | Lanczos (LinOp *Hop, const Tensor &Tin=Tensor(), const std::string method="Gnd", const double &CvgCrit=1.0e-14, const unsigned int &Maxiter=10000, const cytnx_uint64 &k=1, const bool &is_V=true, const bool &is_row=false, const cytnx_uint32 &max_krydim=0, const bool &verbose=false) |
| perform Lanczos for hermitian/symmetric matrices or linear function. | |
| std::vector< Tensor > | Lanczos (LinOp *Hop, const Tensor &Tin=Tensor(), const std::string which="SA", const cytnx_uint64 &maxiter=10000, const cytnx_double &cvg_crit=0, const cytnx_uint64 &k=1, const bool &is_V=true, const cytnx_int32 &ncv=0, const bool &verbose=false) |
| Performs Lanczos iteration for matrices or linear functions. | |
| std::vector< UniTensor > | Lanczos (LinOp *Hop, const cytnx::UniTensor &Tin, const std::string which="SA", const cytnx_uint64 &maxiter=10000, const cytnx_double &cvg_crit=0, const cytnx_uint64 &k=1, const bool &is_V=true, const cytnx_int32 &ncv=0, const bool &verbose=false) |
| Performs Lanczos iteration for matrices or linear functions. | |
| std::vector< UniTensor > | Lanczos (LinOp *Hop, const cytnx::UniTensor &Tin=UniTensor(), const std::string method="Gnd", const double &CvgCrit=1.0e-14, const unsigned int &Maxiter=10000, const cytnx_uint64 &k=1, const bool &is_V=true, const bool &is_row=false, const cytnx_uint32 &max_krydim=4, const bool &verbose=false) |
| perform Lanczos for hermitian/symmetric matrices or linear function. | |
| std::vector< Tensor > | Lanczos_ER (LinOp *Hop, const cytnx_uint64 &k=1, const bool &is_V=true, const cytnx_uint64 &maxiter=10000, const double &CvgCrit=1.0e-14, const bool &is_row=false, const Tensor &Tin=Tensor(), const cytnx_uint32 &max_krydim=4, const bool &verbose=false) |
| perform Lanczos for hermitian/symmetric matrices or linear function. | |
| std::vector< Tensor > | Lanczos_Gnd (LinOp *Hop, const double &CvgCrit=1.0e-14, const bool &is_V=true, const Tensor &Tin=Tensor(), const bool &verbose=false, const unsigned int &Maxiter=100000) |
| perform Lanczos for hermitian/symmetric matrices or linear function to get ground state and lowest eigen value | |
| std::vector< UniTensor > | Lanczos_Gnd_Ut (LinOp *Hop, const cytnx::UniTensor &Tin, const double &CvgCrit=1.0e-14, const bool &is_V=true, const bool &verbose=false, const unsigned int &Maxiter=100000) |
| perform Lanczos for hermitian/symmetric matrices or linear function to get ground state and lowest eigen value | |
| UniTensor | Lanczos_Exp (LinOp *Hop, const cytnx::UniTensor &v, const Scalar &tau, const double &CvgCrit=1.0e-10, const unsigned int &Maxiter=100000, const bool &verbose=false) |
| Perform the Lanczos algorithm for hermitian operator \(H\) to approximate \(e^{H\tau}v\). | |
| std::vector< Tensor > | Lstsq (const Tensor &A, const Tensor &b, const float &rcond=-1) |
| Return the least-squares solution to a linear matrix equation. | |
| Tensor | Axpy (const Scalar &a, const Tensor &x, const Tensor &y=Tensor()) |
| Blas Axpy, performing \( a\textbf{x} + \textbf{y} \), inplacely. | |
| void | Axpy_ (const Scalar &a, const Tensor &x, Tensor &y) |
| Blas Axpy, performing \( \textbf{y} = a\textbf{x} + \textbf{y} \), inplacely. | |
| Tensor | Ger (const Tensor &x, const Tensor &y, const Scalar &a=Scalar()) |
| Blas Ger, performing return = a*vec(x)*vec(y)^T. | |
| void | Gemm_ (const Scalar &a, const Tensor &x, const Tensor &y, const Scalar &b, Tensor &c) |
| Blas Gemm, performing \( \textbf{c} = a\textbf{x}\textbf{y} + b\textbf{c} \), inplacely. | |
| Tensor | Gemm (const Scalar &a, const Tensor &x, const Tensor &y) |
| Blas Gemm, performing \( a\textbf{x}\textbf{y} -> \) return. | |
| void | Gemm_Batch (const std::vector< cytnx_int64 > &m_array, const std::vector< cytnx_int64 > &n_array, const std::vector< cytnx_int64 > &k_array, const std::vector< Scalar > &alpha_array, const std::vector< Tensor > &a_tensors, const std::vector< Tensor > &b_tensors, const std::vector< Scalar > &beta_array, std::vector< Tensor > &c_tensors, const cytnx_int64 group_count, const std::vector< cytnx_int64 > &group_size) |
| Blas Gemm_Batch, performing many(batch) \( \textbf{c} = \alpha\textbf{a}\textbf{b} +
\beta\textbf{c} \), inplacely. You do not need to consider the row-major or column-major, just provide the correct shape. | |
linear algebra related functions.
This namespace contains all the linear algebra related functions. For example, the matrix multiplication, the singular value decomposition, etc. If the linear algebra can only perfom on a matrix, then in most cases,
| void cytnx::linalg::Abs_ | ( | Tensor & | Tin | ) |
inplace perform elementwiase absolute value. @This is just a inplace version of Abs. The input Tensor Tin will be modified.
| [in] | Tin,the | input Tensor. |
| cytnx::UniTensor cytnx::linalg::Add | ( | const cytnx::UniTensor & | Lt, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The addition function between two UniTensor.
This is the addition function for UniTensor. It will perform the element-wise addition. That means if the left UniTensor Lt is given as \( T_L \) and the right UniTensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] + T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the UniTensor \( T_L \) and \( T_R \). It will perform the element-wise addition and note that it will return a new UniTensor object.
Lt and Rt must have the same shape. | cytnx::UniTensor cytnx::linalg::Add | ( | const cytnx::UniTensor & | Lt, |
| const T & | rc | ||
| ) |
The addition function between a UniTensor and a template type.
This is the addition function for UniTensor. It will add the UniTensor and a template type together and add every element in the UniTensor with the template type. That means if the UniTensor Lt is given as \( T_i \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_i[i] + c, \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | Lt | The left UniTensor. |
| [in] | rc | The right template type. supported type: Scalar, cytnx::cytnx_complex128, cytnx::cytnx_complex64, cytnx::cytnx_double, cytnx::cytnx_float, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
rc will be casted to the same type as the UniTensor Lt. | cytnx::UniTensor cytnx::linalg::Add | ( | const T & | lc, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The addition function between a template type and a UniTensor.
This is the addition function for UniTensor. It will add the UniTensor and a template type together and add every element in the UniTensor with the template type. That means if the template type lc is given as \( c \) and the UniTensor Rt is given as \( T_i \), then the result will be:
\[ T_o[i] = c + T_i[i], \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | lc | The left template type. supported type: Scalar, cytnx::cytnx_complex128, cytnx::cytnx_complex64, cytnx::cytnx_double, cytnx::cytnx_float, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
| [in] | Rt | The right UniTensor. |
lc will be casted to the same type as the UniTensor Rt. The addition function for Tensor.
This is the addition function between a Tensor and a template type. It will perform the element-wise addition. That means if the template type lc is given as \( c \) and the Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = c + T_L[i], \]
where \( T_R[i] \) is the elements in the Tensor \( T_R \). It will perform the element-wise addition and note that it will return a new Tensor object.
| [in] | lc | The left template type. |
| [in] | Rt | The right Tensor. |
The addition function for Tensor.
This is the addition function between a Tensor and a template type. It will perform the element-wise addition. That means if the Tensor Lt is given as \( T_L \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_L[i] + c, \]
where \( T_L[i] \) is the elements in the Tensor \( T_L \). It will perform the element-wise addition and note that it will return a new Tensor object.
| [in] | Lt | The left Tensor. |
| [in] | rc | The right template type. |
The addition function for Tensor.
This is the addition function between two Tensor. It will perform the element-wise addition. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] + T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise addition and note that it will return a new Tensor object.
Lt and Rt must be the same. | std::vector< UniTensor > cytnx::linalg::Arnoldi | ( | LinOp * | Hop, |
| const cytnx::UniTensor & | Tin, | ||
| const std::string | which = "LM", |
||
| const cytnx_uint64 & | maxiter = 10000, |
||
| const cytnx_double & | cvg_crit = 1.0e-9, |
||
| const cytnx_uint64 & | k = 1, |
||
| const bool & | is_V = true, |
||
| const cytnx_int32 & | ncv = 0, |
||
| const bool & | verbose = false |
||
| ) |
Performs Arnoldi iteration for matrices or linear functions.
This function calculates the eigenvalue problem using the Arnoldi algorithm. It calls ARPACK routines.
| [in] | Hop | The Linear Operator defined by the LinOp class or its inheritance (see LinOp). |
| [in] | Tin | the initial UniTensor. |
| [in] | which | Specifies the order in which eigenvalues and corresponding eigenvectors should be found. The supported options are:
|
| [in] | maxiter | The maximum number of iteration steps. |
| [in] | cvg_crit | The convergence criterion of the energy. If set to 0, it reaches machine precision. |
| [in] | k | The number of eigenvalues to compute. |
| [in] | is_V | If set to true, the eigenvectors will be returned. |
| [in] | ncv | Related to number of Arnoldi Vectors generated (the Krylov subspace). If 0 (default), the value is typically determined automatically as min(2*k+10, dim), where dim is the dimension of the linear operator (e.g., nx in LinOp). |
| [in] | verbose | Prints out iteration information. |
k <= nx. Where nx is the dimension of the vector.Tin should be a rank-1 Tensor, and the dimension need to same as LinOp->nx().ncv \( \neq 0\), then k+2<= ncv<= dim, where dim = LinOp->nx(). which equal to SI or LI, and the data type of LinOp is real rather than complex, the output of eigenvectors are same as ARPACK. For example, if the full spectrum is \([1.0+0i, 0.7+0.1i,
0.7-0.1i, 0.2+0.2i, 0.2-0.2i]\), and you choose which="LI", k=2, then the output eigenvalues are \([0.2+0.2i, 0.2-0.2i]\) rather than \([0.2+0.2i, 0.7+0.1i]\) (as ARPACK prioritizes returning complete complex conjugate pairs for real operators). | std::vector< Tensor > cytnx::linalg::Arnoldi | ( | LinOp * | Hop, |
| const Tensor & | Tin = Tensor(), |
||
| const std::string | which = "LM", |
||
| const cytnx_uint64 & | maxiter = 10000, |
||
| const cytnx_double & | cvg_crit = 0, |
||
| const cytnx_uint64 & | k = 1, |
||
| const bool & | is_V = true, |
||
| const cytnx_int32 & | ncv = 0, |
||
| const bool & | verbose = false |
||
| ) |
Performs Arnoldi iteration for matrices or linear functions.
This function calculates the eigenvalue problem using the Arnoldi algorithm. It calls ARPACK routines.
| [in] | Hop | The Linear Operator defined by the LinOp class or its inheritance (see LinOp). |
| [in] | Tin | The initial vector; this should be a rank-1 Tensor. |
| [in] | which | Specifies the order in which eigenvalues and corresponding eigenvectors should be found. The supported options are:
|
| [in] | maxiter | The maximum number of iteration steps. |
| [in] | cvg_crit | The convergence criterion of the energy. If set to 0, it reaches machine precision. |
| [in] | k | The number of eigenvalues to compute. |
| [in] | is_V | If set to true, the eigenvectors will be returned. |
| [in] | ncv | Related to number of Arnoldi Vectors generated (the Krylov subspace). If 0 (default), the value is typically determined automatically as min(2*k+10, dim), where dim is the dimension of the linear operator (e.g., nx in LinOp). |
| [in] | verbose | Prints out iteration information. |
std::vector<Tensor> containing: The first Tensor holds the computed eigenvalues. If is_V is true, subsequent Tensor objects contain the corresponding eigenvectors. which equal to SI or LI, and the data type of LinOp is real rather than complex, the output of eigenvectors are same as ARPACK. For example, if the full spectrum is \([1.0+0i,
0.7+0.1i, 0.7-0.1i, 0.2+0.2i, 0.2-0.2i]\), and you choose which="LI", k=2, then the output eigenvalues are \([0.2+0.2i, 0.2-0.2i]\) rather than \([0.2+0.2i, 0.7+0.1i]\) (as ARPACK prioritizes returning complete complex conjugate pairs for real operators). Blas Axpy, performing \( a\textbf{x} + \textbf{y} \), inplacely.
This function performs
\[ a\textbf{x} + \textbf{y}, \]
where \( \textbf{x},\textbf{y} \) are Tensor and \( a \) is a Scalar. The dtype of return Tensor will be the strongest among x, y and a.
Blas Axpy, performing \( \textbf{y} = a\textbf{x} + \textbf{y} \), inplacely.
This function performs
\[ \textbf{y} = a\textbf{x} + \textbf{y}, \]
where \( \textbf{x},\textbf{y} \) are Tensor and a is a Scalar. The dtype of return Tensor will be the strongest among x, y and a.
| cytnx::UniTensor cytnx::linalg::Conj | ( | const cytnx::UniTensor & | UT | ) |
Elementwise conjugate of the UniTensor.
| [in] | UT | The input UniTensor. |
| void cytnx::linalg::Conj_ | ( | cytnx::UniTensor & | UT | ) |
Inplace elementwise conjugate of the UniTensor.
| [in] | UT | The input UniTensor. |
| void cytnx::linalg::Conj_ | ( | Tensor & | Tin | ) |
The comparison function for Tensor.
This is the comparison function between a Tensor and a template type. It will perform the element-wise comparison. That means if the left template type lc is given as \( c \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = \left\{\begin{array}{ll} 1, & \text{if } c=T_R[i]\\ 0, & \text{else} \end{array}\right. \]
where \( T_o[i]\), \( c \) and \( T_R[i] \) are the elements in the Tensor \( T_o\), \( c \) and \( T_R \) and \( T_o[i] \) is the output Tensor which is a boolean type (see cytnx::Type). It will perform the element-wise comparison and note that it will return a new Tensor object.
| [in] | lc | The left template type. |
| [in] | Rt | The right Tensor. |
The comparison function for Tensor.
This is the comparison function between a Tensor and a template type. It will perform the element-wise comparison. That means if the left Tensor Lt is given as \( T_L \) and the right template type rc is given as \( c \), then the result will be:
\[ T_o[i] = \left\{\begin{array}{ll} 1, & \text{if } T_L[i]=c\\ 0, & \text{else} \end{array}\right. \]
where \( T_o[i]\), \( T_L[i] \) and \( c \) are the elements in the Tensor \( T_o\), \( T_L \) and \( c \) and \( T_o[i] \) is the output Tensor which is a boolean type (see cytnx::Type). It will perform the element-wise comparison and note that it will return a new Tensor object.
| [in] | Lt | The left Tensor. |
| [in] | rc | The right template type. |
The comparison function for Tensor.
This is the comparison function between two Tensor. It will perform the element-wise comparison. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = \left\{\begin{array}{ll} 1, & \text{if } T_L[i]=T_R[i]\\ 0, & \text{else} \end{array}\right. \]
where \( T_o[i]\), \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_o\), \( T_L \) and \( T_R \) and \( T_o[i] \) is the output Tensor which is a boolean type (see cytnx::Type). It will perform the element-wise comparison and note that it will return a new Tensor object.
Lt and Rt should have the same shape. Calculate the determinant of a tensor.
This function will calculate the determinant of a square matrix . That means the input tensor should be a rank-2 tensor with shape (N,N).
| [in] | Tl | input a Tensor with shape (N,N) |
return a diagonal tensor with diagonal elements provided as Tin.
This function will return a diagonal tensor with diagonal elements provided as Tin. Furthermore, the return Tensor will be rank-2, with shape=(L, L); where L is the number of elements in Tin.
Tin should be a rank-2 Tensor. | Tensor cytnx::linalg::Directsum | ( | const Tensor & | T1, |
| const Tensor & | T2, | ||
| const std::vector< cytnx_uint64 > & | shared_axes | ||
| ) |
perform directsum of two Tensor.
The function assume two tensor has the same rank, and axes indicated in <shared_axes> are the same for both T1 and T2. The output tensors will have same rank as T1 and T2, with the dimension of rest of the axes being the sum of dimensions of T1 and T2. e.g., the out put shape = (i1+i2,j1+j2, share_axis_1, k1+k2, share_axis_2, ...); where T1.shape = (i1,j1,share_axis_1,k1,share_axis_2 ...) and T2.shape = (i2,j2,share_axis_1,k2,share_axis_2 ...)
| [in] | T1 | rank-n Tensor #1 |
| [in] | T2 | rank-n Tensor #2 |
| [in] | shared_axes | The axes that are shared by two tensors |
| cytnx::UniTensor cytnx::linalg::Div | ( | const cytnx::UniTensor & | Lt, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The division function between two UniTensor.
This is the division function for UniTensor. It will divide the left UniTensor and the right UniTensor together. It will divide every element in the left UniTensor with the right UniTensor. That means if the left UniTensor Lt is given as \( T_L \) and the right UniTensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] \div T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the UniTensor \( T_L \) and \( T_R \). It will perform the element-wise division and note that it will return a new UniTensor object.
Lt and Rt must have the same shape. | cytnx::UniTensor cytnx::linalg::Div | ( | const cytnx::UniTensor & | Lt, |
| const T & | rc | ||
| ) |
The division function between a UniTensor and a template type.
This is the division function for UniTensor. It will divide the UniTensor and a template type together. It will divide every element in the UniTensor with the template type. That means if the UniTensor Lt is given as \( T_i \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_i[i] \div c, \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | Lt | The left UniTensor. |
| [in] | rc | The right template type. supported type: Scalar, cytnx::cytnx_complex128, cytnx::cytnx_complex64, cytnx::cytnx_double, cytnx::cytnx_float, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
rc will be casted to the same type as the UniTensor Lt.| cytnx::UniTensor cytnx::linalg::Div | ( | const T & | lc, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The division function between a UniTensor and a template type.
This is the division function for UniTensor. It will divide the UniTensor and a template type together. It will divide every element in the UniTensor with the template type. That means if the template type lc is given as \( c \) and the UniTensor Rt is given as \( T_i \), then the result will be:
\[ T_o[i] = c \div T_i[i], \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | Lt | The left UniTensor. |
| [in] | rc | The right template type. supported type: Scalar, cytnx::cytnx_complex128, cytnx::cytnx_complex64, cytnx::cytnx_double, cytnx::cytnx_float, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
lc will be casted to the same type as the UniTensor Rt.The division function for Tensor.
This is the division function between a Tensor and a template type. It will perform the element-wise division. That means if the left template type lc is given as \( c \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = c / T_R[i] \]
where \( T_L[i] \) is the elements in the Tensor \( T_L \). It will perform the element-wise division and note that it will return a new Tensor object.
| [in] | lc | The left template type. |
| [in] | Rt | The right Tensor. |
Rt should not contain any zero element. The division function for Tensor.
This is the division function between a Tensor and a template type. It will perform the element-wise division. That means if the left Tensor Lt is given as \( T_L \) and the right template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_L[i] / c \]
where \( T_L[i] \) is the elements in the Tensor \( T_L \). It will perform the element-wise division and note that it will return a new Tensor object.
| [in] | Lt | The left Tensor. |
| [in] | rc | The right template type. |
rc should not be zero. The division function for Tensor.
This is the division function between two Tensor. It will perform the element-wise division. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] / T_R[i] \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise division and note that it will return a new Tensor object.
Rt should not contain any zero element. dot product of two arrays.
| std::vector< UniTensor > cytnx::linalg::Eig | ( | const cytnx::UniTensor & | Tin, |
| const bool & | is_V = true, |
||
| const bool & | row_v = false |
||
| ) |
| std::vector< Tensor > cytnx::linalg::Eig | ( | const Tensor & | Tin, |
| const bool & | is_V = true, |
||
| const bool & | row_v = false |
||
| ) |
eigen-value decomposition for generic square matrix
This function will perform eigen-value decomposition for generic square matrix. Given a matrix Tin as \( M \), then the result will be:
\[ M = V D V^{-1}, \]
where \( V \) is a invertible matrix contains the eigen vectors, and \( D \) is a diagonal matrix contains the eigen values.
| [in] | Tin | A square matrix (a rank-2 Tensor) |
| [in] | is_V | whether need to return the eigen vectors. |
| [in] | row_V | if set to ture, the return eigen vectors will be row form. |
is_V = true. Furthermore, if row_V = true, then the eigen vectors will be row form. Otherwise, the eigen vectors will be column form. Tin should be a square matrix. | std::vector< UniTensor > cytnx::linalg::Eigh | ( | const cytnx::UniTensor & | Tin, |
| const bool & | is_V = true, |
||
| const bool & | row_v = false |
||
| ) |
| std::vector< Tensor > cytnx::linalg::Eigh | ( | const Tensor & | Tin, |
| const bool & | is_V = true, |
||
| const bool & | row_v = false |
||
| ) |
eigen-value decomposition for Hermitian matrix
This function will perform eigen-value decomposition for Hermitian matrix. Given a matrix Tin as \( M \), then the result will be:
\[ M = V D V^\dagger, \]
where \( V \) is a unitary matrix contains the eigen vectors, and \( D \) is a diagonal matrix contains the eigen values.
| [in] | Tin | a Tensor , it should be a rank-2 tensor (matrix) |
| [in] | is_V | whether need to return the eigen vectors. |
| [in] | row_V | if set to ture, the return eigen vectors will be row form. |
is_V = true. Furthermore, if row_V = true, then the eigen vectors will be row form. Otherwise, the eigen vectors will be column form. Tin should be a Hermitian matrix. Tin is not a Hermitian matrix, only the lower triangular matrix will be used. (This is strongly not recommended, please use Eig(const Tensor &Tin, const bool &is_V, const bool &row_v) instead. | void cytnx::linalg::Exp_ | ( | Tensor & | Tin | ) |
inplace perform Exponential on all the element in Tensor.
This function will perform Exponential on all the elements in Tensor Tin. Furthermore,
| [in] | Tin,the | input Tensor. |
| void cytnx::linalg::Expf_ | ( | Tensor & | Tin | ) |
inplace perform Exponential on all the element in Tensor.
This function will perform Exponential on all the elements in Tensor Tin. Furthermore,
| [in] | Tin,the | input Tensor. |
| cytnx::UniTensor cytnx::linalg::ExpH | ( | const cytnx::UniTensor & | Tin | ) |
Perform the exponential function on a UniTensor, which the blocks are Hermitian matrix.
This function performs the exponential function on a UniTensor Tin, which the blocks are Hermitian matrix. For more details, please refer to the documentation of the function ExpH(const Tensor &Tin)
| cytnx::UniTensor cytnx::linalg::ExpH | ( | const cytnx::UniTensor & | Tin, |
| const T & | a, | ||
| const T & | b = 0 |
||
| ) |
Perform the exponential function on a UniTensor, which the blocks are Hermitian matrix.
This function performs the exponential function on a UniTensor Tin, which the blocks are Hermitian matrix. For more details, please refer to the documentation of the function ExpH(const Tensor &Tin, const T &a, const T &b).
perform matrix exponential for Hermitian matrix
This function perform matrix exponential for Hermitian matrix, That is,
\[ O = \exp{M} \]
in should be a Hermitian matrix. in is not a Hermitian matrix, only the lower triangular matrix will be used. (This is strongly not recommended, please use ExpM(const Tensor &in) instead).| Tensor cytnx::linalg::ExpH | ( | const Tensor & | in, |
| const T & | a, | ||
| const T & | b = 0 |
||
| ) |
perform matrix exponential for Hermitian matrix
This function perform matrix exponential for Hermitian matrix, That is,
\[ O = \exp{(aM + b)} \]
| [in] | in | input Tensor, should be Hermitian |
| [in] | a | rescale factor |
| [in] | b | bias |
in should be a Hermitian matrix. in is not a Hermitian matrix, only the lower triangular matrix will be used. (This is strongly not recommended, please use ExpM(const Tensor &in, const T &a, const T &b) instead). | cytnx::UniTensor cytnx::linalg::ExpM | ( | const cytnx::UniTensor & | Tin | ) |
Perform the exponential function on a UniTensor.
This function performs the exponential function on a UniTensor Tin. For more details, please refer to the documentation of the function ExpM(const Tensor &Tin)
| cytnx::UniTensor cytnx::linalg::ExpM | ( | const cytnx::UniTensor & | Tin, |
| const T & | a, | ||
| const T & | b = 0 |
||
| ) |
Perform the exponential function on a UniTensor.
This function performs the exponential function on a UniTensor Tin. For more details, please refer to the documentation of the function ExpM(const Tensor &Tin, const T &a, const T &b).
| void cytnx::linalg::Gemm_ | ( | const Scalar & | a, |
| const Tensor & | x, | ||
| const Tensor & | y, | ||
| const Scalar & | b, | ||
| Tensor & | c | ||
| ) |
Blas Gemm, performing \( \textbf{c} = a\textbf{x}\textbf{y} + b\textbf{c} \), inplacely.
This function performs
\[ \textbf{c} = a\textbf{x}\textbf{y} + b\textbf{c}, \]
where \( \textbf{x},\textbf{y},\textbf{c} \) are rank-2 Tensor and a, b are Scalars. The dtype of return Tensor will be the strongest among x, y and a b c.
| void cytnx::linalg::Gemm_Batch | ( | const std::vector< cytnx_int64 > & | m_array, |
| const std::vector< cytnx_int64 > & | n_array, | ||
| const std::vector< cytnx_int64 > & | k_array, | ||
| const std::vector< Scalar > & | alpha_array, | ||
| const std::vector< Tensor > & | a_tensors, | ||
| const std::vector< Tensor > & | b_tensors, | ||
| const std::vector< Scalar > & | beta_array, | ||
| std::vector< Tensor > & | c_tensors, | ||
| const cytnx_int64 | group_count, | ||
| const std::vector< cytnx_int64 > & | group_size | ||
| ) |
Blas Gemm_Batch, performing many(batch) \( \textbf{c} = \alpha\textbf{a}\textbf{b} + \beta\textbf{c} \), inplacely. You do not need to consider the row-major or column-major, just provide the correct shape.
see https://www.intel.com/content/www/us/en/developer/articles/technical/introducing-batch-gemm-operations.html This function performs tensor type check and type conversion, then call the corresponding blas function.
| [in] | m_array | array of cytnx_int64, each element is the number of rows of a_tensors |
| [in] | n_array | array of cytnx_int64, each element is the number of columns of b_tensors |
| [in] | k_array | array of cytnx_int64, each element is the number of columns of a_tensors and the number of rows of b_tensors |
| [in] | alpha_array | array of Scalar, each element is the scalar alpha |
| [in] | a_tensors | array of Tensor, each element is a rank-2 Tensor with shape (m_array[i],k_array[i]) |
| [in] | b_tensors | array of Tensor, each element is a rank-2 Tensor with shape (k_array[i],n_array[i]) |
| [in] | beta_array | array of Scalar, each element is the scalar beta |
| [in,out] | c_tensors | array of Tensor, each element is a rank-2 Tensor with shape (m_array[i],n_array[i]), {must be properly initialized with the correct shape}. |
| [in] | group_count | cytnx_int64, the number of groups |
| [in] | group_size | array of cytnx_int64, each element is the number of matrices in each group |
Blas Ger, performing return = a*vec(x)*vec(y)^T.
This function performs a*x*y^T where x,y are rank-1 Tensor with dimension nx and ny respectively; and a is a Scalar. The dtype of return Tensor will be the strongest among x,y and a.
| [in] | x | Tensor, rank-1 with size nx |
| [in] | y | Tensor, rank-1 with size ny |
| [in] | a | Scalar, if not provided a = 1. |
| std::vector< cytnx::UniTensor > cytnx::linalg::Gesvd | ( | const cytnx::UniTensor & | Tin, |
| const bool & | is_U = true, |
||
| const bool & | is_vT = true |
||
| ) |
Perform Singular-Value decomposition on a UniTensor using ?gesvd method.
This function performs the Singular-Value decomposition on a UniTensor Tin. The result will depend on the rowrank of the UniTensor Tin. For more details, please refer to the documentation of the functions Gesvd(const Tensor &Tin, const bool &is_U, const
bool &is_vT) and Svd(const Tensor &Tin, const bool &is_UvT).
| std::vector< Tensor > cytnx::linalg::Gesvd | ( | const Tensor & | Tin, |
| const bool & | is_U = true, |
||
| const bool & | is_vT = true |
||
| ) |
Perform Singular-Value decomposition on a rank-2 Tensor (a matrix).
This function will perform Singular-Value decomposition on a matrix (a rank-2 Tensor). That means given a matrix Tin as \( M \), then the result will be:
\[ M = U S V^\dagger, \]
where \( U \) is a left uniform matrix, \( S \) is a diagonal matrix with singular values, and \( V^\dagger \) is the conjugate transpose of the right uniform matrix \( V \). Furthermore, \( U \) and \( V \) are unitary matrices, and \( S \) is a non-negative diagonal matrix.
| [in] | Tin | a Tensor, it should be a rank-2 tensor (matrix) |
| [in] | is_U | whether need to return left unitary matrix. |
| [in] | is_vT | whether need to return right unitary matrix |
[std::vector<Tensors>]
is_U is true, then the tensor \( U \) will be pushed back to the vector, and if is_vT is true, \( V^\dagger \) will be pushed back to the vector. | std::vector< cytnx::UniTensor > cytnx::linalg::Gesvd_truncate | ( | const cytnx::UniTensor & | Tin, |
| const cytnx_uint64 & | keepdim, | ||
| const double & | err = 0., |
||
| const bool & | is_U = true, |
||
| const bool & | is_vT = true, |
||
| const unsigned int & | return_err = 0, |
||
| const cytnx_uint64 & | mindim = 1 |
||
| ) |
Perform Singular-Value decomposition on a UniTensor with truncation.
This function performs the Singular-Value decomposition of a UniTensor Tin and truncates the singular values. The result will depend on the rowrank of the UniTensor Tin. This version uses the ?gesvd method. See references below for more details.
| std::vector< cytnx::UniTensor > cytnx::linalg::Gesvd_truncate | ( | const cytnx::UniTensor & | Tin, |
| const cytnx_uint64 & | keepdim, | ||
| std::vector< cytnx_uint64 > | min_blockdim, | ||
| const double & | err = 0., |
||
| const bool & | is_U = true, |
||
| const bool & | is_vT = true, |
||
| const unsigned int & | return_err = 0, |
||
| const cytnx_uint64 & | mindim = 1 |
||
| ) |
Perform Singular-Value decomposition on a UniTensor with truncation and keep at most min_blockdim singular values in each block.
This function performs the Singular-Value decomposition of a UniTensor Tin and truncates the singular values. This version uses the ?gesvd method. See references below for more details.
| std::vector< Tensor > cytnx::linalg::Gesvd_truncate | ( | const Tensor & | Tin, |
| const cytnx_uint64 & | keepdim, | ||
| const double & | err = 0., |
||
| const bool & | is_U = true, |
||
| const bool & | is_vT = true, |
||
| const unsigned int & | return_err = 0, |
||
| const cytnx_uint64 & | mindim = 1 |
||
| ) |
Perform a truncated Singular-Value decomposition of a rank-2 Tensor (a matrix) by a randomized SVD.
| [in] | Tin | a Tensor, it should be a rank-2 tensor (matrix) |
| [in] | keepdim | the number (at most) of singular values to keep. |
| [in] | err | the cutoff error (the singular values smaller than err will be truncated.) |
| [in] | is_UvT | if true, the left- and right- unitary matrices (isometries) are returned. |
| [in] | return_err | whether the error shall be returned. If return_err is true, then the largest error will be pushed back to the vector (The smallest singular value in the return singular values matrix \( S \).) If return_err is a positive int, then the full list of truncated singular values will be returned. |
[std::vector<Tensors>]
is_UvT is true, then the tensors \( U \) and \( V^\dagger \) will be pushed back to the vector.return_err is true, then the error will be pushed back to the vector. keepdim for degenerate singular values: if the largest \( n \) truncated singular values would be exactly equal to the smallest kept singular value, then the bond dimension is enlarged to keepdim \( + n \). Example: if the singular values are (1 2 2 2 2 3) and keepdim = 3, then the bond dimension will be 5 in order to keep all the degenerate singular values. | std::vector< cytnx::UniTensor > cytnx::linalg::Hosvd | ( | const cytnx::UniTensor & | Tin, |
| const std::vector< cytnx_uint64 > & | mode, | ||
| const bool & | is_core = true, |
||
| const bool & | is_Ls = false, |
||
| const std::vector< cytnx_int64 > & | trucate_dim = std::vector< cytnx_int64 >() |
||
| ) |
| std::vector< Tensor > cytnx::linalg::Hosvd | ( | const Tensor & | Tin, |
| const std::vector< cytnx_uint64 > & | mode, | ||
| const bool & | is_core = true, |
||
| const bool & | is_Ls = false, |
||
| const std::vector< cytnx_int64 > & | trucate_dim = std::vector< cytnx_int64 >() |
||
| ) |
The addition function for Tensor, inplacely.
This is the inplace addition function between two Tensor. It will perform the element-wise addition. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_L[i] = T_L[i] + T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise addition and note that it will modify the left Tensor Lt.
Lt and Rt must be the same. Lt. The inplace division function for Tensor, inplacely.
This is the inplace division function between two Tensor. It will perform the element-wise division. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] / T_R[i] \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise division and note that it will modify the left Tensor Lt.
Rt should not contain any zero element. Lt. The multiplication function for Tensor, inplacely.
This is the multiplication function between two Tensor. It will perform the element-wise multiplication. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] * T_R[i] \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise multiplication and note that it will modify the left Tensor Lt.
Lt. Element-wise inverse with clip.
This function will perform Element-wise inverse with clip. If A[i] < clip, then 1/A[i] = 0 will be set. That is, the out put will be:
\[ A_{out} = \left\{ \begin{array}{ll} 1/A[i] & \mathrm{if} \ A[i] \geq \mathrm{clip} \\ 0 & \mathrm{otherwise} \end{array} \right. \]
| [in] | Tin | a Tensor |
| [in] | clip | the clip value. |
| void cytnx::linalg::Inv_ | ( | Tensor & | Tin, |
| const double & | clip | ||
| ) |
inplace perform Element-wise inverse with clip.
This function will perform Element-wise inverse with clip. This function is just as same as Inv, but it will modify the input Tensor inplace.
| [in] | Tin | a Tensor |
| [in] | clip | the clip value. |
| UniTensor cytnx::linalg::InvM | ( | const cytnx::UniTensor & | Tin | ) |
| void cytnx::linalg::InvM_ | ( | Tensor & | Tin | ) |
inplace matrix inverse.
This function will perform matrix inverse on the input matrix Tin, inplacely.
| void cytnx::linalg::InvM_ | ( | UniTensor & | Tin | ) |
The subtraction function for Tensot, inplscely.
This is the subtraction function between two Tensor. It will perform the element-wise subtraction. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_L[i] = T_L[i] - T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise subtraction and note that it will modify the left Tensor Lt.
Lt and Rt must be the same. Lt. | Tensor cytnx::linalg::Kron | ( | const Tensor & | Tl, |
| const Tensor & | Tr, | ||
| const bool & | Tl_pad_left = false, |
||
| const bool & | Tr_pad_left = false |
||
| ) |
perform kronecker produces of two Tensor.
This function will perform kronecker produces of two Tensor Tl and Tr. Furthermore, the function assume two tensor has the same rank. In case where two tensors have different ranks, the small one will be extend by adding redundant dimension to the beginning of axis (T<x>_pad_right=true) or by adding redundant dim to the last axis (if T<x>_pad_left=false [default]). if the Tensor #1 has shape=(i1,j1,k1,l1...), and Tensor #2 has shape=(i2,j2,k2,l2...); then the return Tensor will have shape=(i1*i2,j1*j2,k1*k2...)
| [in] | Tl | rank-n Tensor #1 |
| [in] | Tr | rank-m Tensor #2 |
| [in] | Tl_pad_left | The padding scheme for Tl if Tl.rank != Tr.rank |
| [in] | Tr_pad_left | The padding scheme for Tr if Tl.rank != Tr.rank |
| std::vector< UniTensor > cytnx::linalg::Lanczos | ( | LinOp * | Hop, |
| const cytnx::UniTensor & | Tin, | ||
| const std::string | which = "SA", |
||
| const cytnx_uint64 & | maxiter = 10000, |
||
| const cytnx_double & | cvg_crit = 0, |
||
| const cytnx_uint64 & | k = 1, |
||
| const bool & | is_V = true, |
||
| const cytnx_int32 & | ncv = 0, |
||
| const bool & | verbose = false |
||
| ) |
Performs Lanczos iteration for matrices or linear functions.
This function calculates the eigenvalue problem using the Lanczos algorithm. It calls ARPACK routines.
| [in] | Hop | The Linear Operator defined by the LinOp class or its inheritance (see LinOp). |
| [in] | Tin | the initial UniTensor. |
| [in] | which | Specifies the order in which eigenvalues and corresponding eigenvectors should be found. The supported options are:
|
| [in] | maxiter | The maximum number of iteration steps. |
| [in] | cvg_crit | The convergence criterion of the energy. If set to 0, it reaches machine precision. |
| [in] | k | The number of eigenvalues to compute. |
| [in] | is_V | If set to true, the eigenvectors will be returned. |
| [in] | ncv | Related to number of Arnoldi Vectors generated (the Krylov subspace). If 0 (default), the value is typically determined automatically as min(2*k+10, dim), where dim is the dimension of the linear operator (e.g., nx in LinOp). |
| [in] | verbose | Prints out iteration information. |
k <= nx. Where nx is the dimension of the vector.Tin should be a rank-1 Tensor, and the dimension need to same as LinOp->nx().ncv \( \neq 0\), then k+2<= ncv<= dim, where dim = LinOp->nx(). | std::vector< UniTensor > cytnx::linalg::Lanczos | ( | LinOp * | Hop, |
| const cytnx::UniTensor & | Tin = UniTensor(), |
||
| const std::string | method = "Gnd", |
||
| const double & | CvgCrit = 1.0e-14, |
||
| const unsigned int & | Maxiter = 10000, |
||
| const cytnx_uint64 & | k = 1, |
||
| const bool & | is_V = true, |
||
| const bool & | is_row = false, |
||
| const cytnx_uint32 & | max_krydim = 4, |
||
| const bool & | verbose = false |
||
| ) |
perform Lanczos for hermitian/symmetric matrices or linear function.
This function calculate the eigen value problem using explicitly restarted Lanczos. #Performance tune: For small linear dimension, try to reduce max_krydim.
| [in] | Hop | the Linear Operator defined by LinOp class or it's inheritance (see LinOp). |
| [in] | Tin | the initial vector, this should be rank-1. |
| [in] | method | the desired Lanczos method to use, the supported options are: 'ER' : explicitly restarted Lanczos 'Gnd' : naive Lanczos |
| [in] | CvgCrit | the convergence criterion of the energy. |
| [in] | maxiter | the maximum interation steps for each k. |
| [in] | k | the number of lowest k eigen values. |
| [in] | is_V | if set to true, the eigen vectors will be returned. |
| [in] | is_row | whether the return eigen vectors should be in row-major form. |
| [in] | max_krydim | the maximum krylov subspace dimension for each iteration. |
| [in] | verbose | print out iteration info. |
| std::vector< Tensor > cytnx::linalg::Lanczos | ( | LinOp * | Hop, |
| const Tensor & | Tin = Tensor(), |
||
| const std::string | method = "Gnd", |
||
| const double & | CvgCrit = 1.0e-14, |
||
| const unsigned int & | Maxiter = 10000, |
||
| const cytnx_uint64 & | k = 1, |
||
| const bool & | is_V = true, |
||
| const bool & | is_row = false, |
||
| const cytnx_uint32 & | max_krydim = 0, |
||
| const bool & | verbose = false |
||
| ) |
perform Lanczos for hermitian/symmetric matrices or linear function.
This function calculate the eigen value problem using explicitly restarted Lanczos. #Performance tune: For small linear dimension, try to reduce max_krydim.
| [in] | Hop | the Linear Operator defined by LinOp class or it's inheritance (see LinOp). |
| [in] | Tin | the initial vector, this should be rank-1. |
| [in] | method | the desired Lanczos method to use, the supported options are: 'ER' : explicitly restarted Lanczos 'Gnd' : naive Lanczos |
| [in] | CvgCrit | the convergence criterion of the energy. |
| [in] | maxiter | the maximum interation steps for each k. |
| [in] | k | the number of lowest k eigen values. |
| [in] | is_V | if set to true, the eigen vectors will be returned. |
| [in] | is_row | whether the return eigen vectors should be in row-major form. |
| [in] | max_krydim | the maximum krylov subspace dimension for each iteration. |
| [in] | verbose | print out iteration info. |
| std::vector< Tensor > cytnx::linalg::Lanczos | ( | LinOp * | Hop, |
| const Tensor & | Tin = Tensor(), |
||
| const std::string | which = "SA", |
||
| const cytnx_uint64 & | maxiter = 10000, |
||
| const cytnx_double & | cvg_crit = 0, |
||
| const cytnx_uint64 & | k = 1, |
||
| const bool & | is_V = true, |
||
| const cytnx_int32 & | ncv = 0, |
||
| const bool & | verbose = false |
||
| ) |
Performs Lanczos iteration for matrices or linear functions.
This function calculates the eigenvalue problem using the Lanczos algorithm. It calls ARPACK routines.
| [in] | Hop | The Linear Operator defined by the LinOp class or its inheritance (see LinOp). |
| [in] | Tin | the initial vector, this should be rank-1. |
| [in] | which | Specifies the order in which eigenvalues and corresponding eigenvectors should be found. The supported options are:
|
| [in] | maxiter | The maximum number of iteration steps. |
| [in] | cvg_crit | The convergence criterion of the energy. If set to 0, it reaches machine precision. |
| [in] | k | The number of eigenvalues to compute. |
| [in] | is_V | If set to true, the eigenvectors will be returned. |
| [in] | ncv | Related to number of Arnoldi Vectors generated (the Krylov subspace). If 0 (default), the value is typically determined automatically as min(2*k+10, dim), where dim is the dimension of the linear operator (e.g., nx in LinOp). |
| [in] | verbose | Prints out iteration information. |
| std::vector< Tensor > cytnx::linalg::Lanczos_ER | ( | LinOp * | Hop, |
| const cytnx_uint64 & | k = 1, |
||
| const bool & | is_V = true, |
||
| const cytnx_uint64 & | maxiter = 10000, |
||
| const double & | CvgCrit = 1.0e-14, |
||
| const bool & | is_row = false, |
||
| const Tensor & | Tin = Tensor(), |
||
| const cytnx_uint32 & | max_krydim = 4, |
||
| const bool & | verbose = false |
||
| ) |
perform Lanczos for hermitian/symmetric matrices or linear function.
This function calculate the eigen value problem using explicitly restarted Lanczos. #Performance tune: For small linear dimension, try to reduce max_krydim.
| [in] | Hop | the Linear Operator defined by LinOp class or it's inheritance (see LinOp). |
| [in] | k | the number of lowest k eigen values. |
| [in] | is_V | if set to true, the eigen vectors will be returned. |
| [in] | maxiter | the maximum interation steps for each k. |
| [in] | CvgCrit | the convergence criterion of the energy. |
| [in] | is_row | whether the return eigen vectors should be in row-major form. |
| [in] | Tin | the initial vector, this should be rank-1 |
| [in] | max_krydim | the maximum krylov subspace dimension for each iteration. |
| [in] | verbose | print out iteration info. |
| UniTensor cytnx::linalg::Lanczos_Exp | ( | LinOp * | Hop, |
| const cytnx::UniTensor & | v, | ||
| const Scalar & | tau, | ||
| const double & | CvgCrit = 1.0e-10, |
||
| const unsigned int & | Maxiter = 100000, |
||
| const bool & | verbose = false |
||
| ) |
Perform the Lanczos algorithm for hermitian operator \(H\) to approximate \(e^{H\tau}v\).
This function perform the Lanczos-like algorithm for hermitian linear operator \(H\) to approximate
\[ e^{H\tau}v \]
and return the state \(w\) such that
\[ |\exp(H\tau)v - w| < \delta. \]
Here \(v\) is a given vector or a state.
| [in] | Hop | the Linear Operator defined by LinOp class or it's inheritance (see LinOp). The operation method \(Hv\) need to be defined in it. |
| [in] | v | The input vector (or state). The norm \(|v|\) should be equal to 1. |
| [in] | tau | A scalar, it can be complex number. |
| [in] | CvgCrit | \(\delta\), the convergence criterion. |
| [in] | Maxiter | the maximum interation steps for each k. |
| [in] | verbose | print out iteration info. |
| std::vector< Tensor > cytnx::linalg::Lanczos_Gnd | ( | LinOp * | Hop, |
| const double & | CvgCrit = 1.0e-14, |
||
| const bool & | is_V = true, |
||
| const Tensor & | Tin = Tensor(), |
||
| const bool & | verbose = false, |
||
| const unsigned int & | Maxiter = 100000 |
||
| ) |
perform Lanczos for hermitian/symmetric matrices or linear function to get ground state and lowest eigen value
This function calculate the eigen value problem using naive Lanczos to get ground state and lowest eigen value.
| [in] | Hop | the Linear Operator defined by LinOp class or it's inheritance (see LinOp). |
| [in] | CvgCrit | the convergence criterion of the energy. |
| [in] | is_V | if set to true, the eigen vectors will be returned. |
| [in] | Tin | the initial vector, this should be rank-1 |
| [in] | verbose | print out iteration info. |
| [in] | maxiter | the maximum interation steps for each k. |
| std::vector< UniTensor > cytnx::linalg::Lanczos_Gnd_Ut | ( | LinOp * | Hop, |
| const cytnx::UniTensor & | Tin, | ||
| const double & | CvgCrit = 1.0e-14, |
||
| const bool & | is_V = true, |
||
| const bool & | verbose = false, |
||
| const unsigned int & | Maxiter = 100000 |
||
| ) |
perform Lanczos for hermitian/symmetric matrices or linear function to get ground state and lowest eigen value
This function calculate the eigen value problem using naive Lanczos to get ground state and lowest eigen value.
| [in] | Hop | the Linear Operator defined by LinOp class or it's inheritance (see LinOp). |
| [in] | CvgCrit | the convergence criterion of the energy. |
| [in] | is_V | if set to true, the eigen vectors will be returned. |
| [in] | Tin | the initial vector, this should be a UniTensor. |
| [in] | verbose | print out iteration info. |
| [in] | maxiter | the maximum interation steps for each k. |
| std::vector< Tensor > cytnx::linalg::Lstsq | ( | const Tensor & | A, |
| const Tensor & | b, | ||
| const float & | rcond = -1 |
||
| ) |
Return the least-squares solution to a linear matrix equation.
Computes the vector x that approximatively solves the equation A @ x = b. The equation may be under-, well-, or over-determined independent columns. If a is square and of full rank, then x (but for round-off error) is the "exact" solution of the equation. Else, x minimizes the Euclidean 2-norm \( || b - a x ||_2 \).
| [in] | A | "Coefficient" matrix, must be two-dimensional. |
| [in] | b | Ordinate or "dependent variable" values, must be two-dimensional, the least-squares solution is calculated for each of the K columns of b. |
| [in] | 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, If it is negative, the machine precision is used. |
get the maximum element.
| [in] | Tn | a cytnx::Tensor |
get the minimum element.
| [in] | Tn | a cytnx::Tensor |
| cytnx::UniTensor cytnx::linalg::Mod | ( | const cytnx::UniTensor & | Lt, |
| const cytnx::UniTensor & | Rt | ||
| ) |
element-wise modulo
The modulo function between two UniTensor.
This is the modulo function for UniTensor. It will perform the element-wise modulo between two UniTensor. That means if the UniTensor Lt is given as \( T_L \) and the UniTensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] \mod T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the UniTensor \( T_L \) and \( T_R \).
Lt and Rt must have the same shape.Lt and Rt need to be integer type. | cytnx::UniTensor cytnx::linalg::Mod | ( | const cytnx::UniTensor & | Lt, |
| const T & | rc | ||
| ) |
The modulo function between a UniTensor and a template type.
This is the modulo function for UniTensor. It will perform the element-wise modulo between a UniTensor and a template type. That means if the UniTensor Lt is given as \( T_i \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_i[i] \mod c, \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | Lt | The left UniTensor. |
| [in] | rc | The right template type. supported type: Scalar, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
Lt and rc need to be integer type. rc will be casted to the same type as the UniTensor Lt. | cytnx::UniTensor cytnx::linalg::Mod | ( | const T & | lc, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The modulo function between a UniTensor and a template type.
This is the modulo function for UniTensor. It will perform the element-wise modulo between a UniTensor and a template type. That means if the template type lc is given as \( c \) and the UniTensor Rt is given as \( T_i \), then the result will be:
\[ T_o[i] = c \mod T_i[i], \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| lc | The left template type. The right template type. supported type: Scalar, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
| Rt | The right UniTensor. |
lc and Rt need to be integer type. lc will be casted to the same type as the UniTensor Rt. The mod function for Tensor.
This is the mod function between a Tensor and a template type. It will perform the element-wise mod. That means if the left template type lc is given as \( c \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = c % T_R[i] \]
where \( T_R[i] \) is the elements in the Tensor \( T_R \). It will perform the element-wise mod and note that it will return a new Tensor object.
| [in] | lc | The left template type. |
| [in] | Rt | The right Tensor. |
rc should be integer type. The mod function for Tensor.
This is the mod function between a Tensor and a template type. It will perform the element-wise mod. That means if the left Tensor Lt is given as \( T_L \) and the right template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_L[i] % c \]
where \( T_L[i] \) is the elements in the Tensor \( T_L \). It will perform the element-wise mod and note that it will return a new Tensor object.
| [in] | Lt | The left Tensor. |
| [in] | rc | The right template type. |
rc should be integer type. The mod function for Tensor.
This is the mod function between two Tensor. It will perform the element-wise mod. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] % T_R[i] \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise mod and note that it will return a new Tensor object.
Lt and Rt should have the same shape and need to be integer type. | cytnx::UniTensor cytnx::linalg::Mul | ( | const cytnx::UniTensor & | Lt, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The multiplication function between two UniTensor.
This is the multiplication function for UniTensor. It will multiply the two UniTensor together. It will multiply every element in the UniTensor Lt with the corresponding element in the UniTensor Rt. That means if the UniTensor Lt is given as \( T_L \) and the UniTensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] \times T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the UniTensor \( T_L \) and \( T_R \). It will perform the element-wise multiplication and note that it will return a new UniTensor object.
Lt and Rt must have the same shape. | cytnx::UniTensor cytnx::linalg::Mul | ( | const cytnx::UniTensor & | Lt, |
| const T & | rc | ||
| ) |
The multiplication function between a UniTensor and a template type.
This is the multiplication function for UniTensor. It will multiply the UniTensor and a template type together. It will multiply every element in the UniTensor with the template type. That means if the UniTensor Lt is given as \( T_i \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_i[i] \times c, \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | Lt | The left UniTensor. |
| [in] | rc | The right template type. supported type: Scalar, cytnx::cytnx_complex128, cytnx::cytnx_complex64, cytnx::cytnx_double, cytnx::cytnx_float, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
rc will be casted to the same type as the UniTensor Lt. | cytnx::UniTensor cytnx::linalg::Mul | ( | const T & | lc, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The multiplication function between a UniTensor and a template type.
This is the multiplication function for UniTensor. It will multiply the UniTensor and a template type together. It will multiply every element in the UniTensor with the template type. That means if the template type lc is given as \( c \) and the UniTensor Rt is given as \( T_i \), then the result will be:
\[ T_o[i] = c \times T_i[i], \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | Lt | The left UniTensor. |
| [in] | rc | The right template type. supported type: Scalar, cytnx::cytnx_complex128, cytnx::cytnx_complex64, cytnx::cytnx_double, cytnx::cytnx_float, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
lc will be casted to the same type as the UniTensor Rt. The multiplication function for Tensor.
This is the multiplication function between a Tensor and a template type. It will perform the element-wise multiplication. That means if the left Tensor Lt is given as \( T_L \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_L[i] * c, \]
where \( T_L[i] \) is the elements in the Tensor \( T_L \). It will perform the element-wise multiplication and note that it will return a new Tensor object.
| [in] | Lt | The left Tensor. |
| [in] | rc | The right template type. |
The multiplication function for Tensor.
This is the multiplication function between a Tensor and a template type. It will perform the element-wise multiplication. That means if the left Tensor Lt is given as \( T_L \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_L[i] * c, \]
where \( T_L[i] \) is the elements in the Tensor \( T_L \). It will perform the element-wise multiplication and note that it will return a new Tensor object.
| [in] | Lt | The left Tensor. |
| [in] | rc | The right template type. |
The multiplication function for Tensor.
This is the multiplication function between two Tensor. It will perform the element-wise multiplication. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] * T_R[i] \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise multiplication and note that it will return a new Tensor object.
| Tensor cytnx::linalg::Norm | ( | const cytnx::UniTensor & | uTl | ) |
Calculate the norm of an UniTensor.
This function will calculate the norm of an UniTensor. If the input UniTensor is rank-1, then the frobenius norm will be calculated. If the input UniTensor is rank-N with N>=2, then the blocks will be flatten (see flatten()) to 1d first, and then calculate the frobenius norm.
| [in] | uTl | input UniTensor |
Calculate the norm of a tensor.
This function will calculate the norm of a tensor. If the input tensor is rank-1, then the frobenius norm will be calculated. If the input tensor is rank-N with N>=2, then the tensor will be flatten (see flatten()) to 1d first, and then calculate the frobenius norm.
| [in] | Tl | input Tensor |
perform outer produces of two rank-1 Tensor.
This function will perform outer produces of two rank-1 Tensor Tl and Tr. Furthermore, if the Tensor Tl has [shape_1], and Tensor Tr has [shape_2]; then the return Tensor will have shape: concate(shape_1,shape_2)
| UniTensor cytnx::linalg::Pow | ( | const cytnx::UniTensor & | Tin, |
| const double & | p | ||
| ) |
take power p on all the elements in UniTensor.
This function will take power p on all the elements in UniTensor.
| [in] | Tin | the input UniTensor |
| [in] | p | the power |
Tin is a real UniTensor and containt negative elements, then p must be an integer. | void cytnx::linalg::Pow_ | ( | Tensor & | Tin, |
| const double & | p | ||
| ) |
| void cytnx::linalg::Pow_ | ( | UniTensor & | Tin, |
| const double & | p | ||
| ) |
Take power p on all the elements in UniTensor, inplacely.
This function will take power p on all the elements in UniTensor, inplacely.
| [in,out] | Tin | the input UniTensor |
| [in] | p | the power |
Tin is a real UniTensor and containt negative elements, then p must be an integer. | std::vector< cytnx::UniTensor > cytnx::linalg::Qdr | ( | const cytnx::UniTensor & | Tin, |
| const bool & | is_tau = false |
||
| ) |
Perform the QDR decomposition on a UniTensor.
This function performs the QDR decomposition on a UniTensor Tin. The result will depend on the rowrank of the UniTensor Tin. For more details, please refer to the documentation of the function Qdr(const Tensor &Tin, const bool &is_tau).
Perform QDR decomposition on a rank-2 Tensor.
| [in] | Tin | a cytnx::Tensor, it should be a rank-2 tensor (matrix) |
| [in] | is_tau | if return the tau that contains the Householder reflectors that generate q along with r. The tau array contains scaling factors for the reflectors |
[std::vector<Tensors>]
is_tau = true. | std::vector< cytnx::UniTensor > cytnx::linalg::Qr | ( | const cytnx::UniTensor & | Tin, |
| const bool & | is_tau = false |
||
| ) |
Perform the QR decomposition on a UniTensor.
This function performs the QR decomposition on a UniTensor Tin. The result will depend on the rowrank of the UniTensor Tin. For more details, please refer to the documentation of the function Qr(const Tensor &Tin, const bool &is_tau).
Perform QR decomposition on a rank-2 Tensor.
This function will perform QR decomposition on a matrix (a rank-2 Tensor). That means given a matrix Tin as \( M \), then the result will be:
\[ M = Q R, \]
where \( Q \) is a orthogonal matrix, and \( R \) is a right-upper triangular matrix.
| [in] | Tin | a Tensor, it should be a rank-2 tensor (a matrix) |
| [in] | is_tau | if return the tau that contains the Householder reflectors that generate q along with r. The tau array contains scaling factors for the reflectors |
[std::vector<Tensors>]
is_tau = true. | Tensor cytnx::linalg::Rand_isometry | ( | const Tensor & | Tin, |
| const cytnx_uint64 & | keepdim, | ||
| const cytnx_uint64 & | power_iteration = 2, |
||
| const unsigned int & | seed = random::__static_random_device() |
||
| ) |
Generate an isometrized left isometry for a rank-2 Tensor (a matrix), such that Q * Qdag * Tin approximates Tin (Qdag is the hermitean conjugate of Q)
The power iteration algorithm with Gaussian random matrices is used: 1) Generate Gaussian random matrix X 2) Y = Tin * X 3) QR factorization: Y = Q * R 4) [U', S, vT] = SVD(Qdag * Tom) where Qdag is the hermitian conjugate of Q 5) U = Q * U' If power_iteration > 0, additional steps are inserted between 3) and 4): 3.1) Y = Tdag * Q where Tdag is the hermitian conjugate of Q 3.2) QR factorization: Y = Q * R 3.3) Y = Tin * Q 3.4) QR factorization: Y = Q * R steps 3.1) to 3.4) are repeated power_iteration times; this power method results in an exponential convergence of the singular values and is especially suitable if the singular values decay slowly. The result will depend on the rowrank of the UniTensor Tin. For more details, please refer to the documentation of the function Gesvd(const Tensor &Tin, const bool &is_UvT).
| [in] | Tin | a Tensor, it should be a rank-2 tensor (matrix) |
| [in] | keepdim | the number of singular values to keep. Note that this is the size of the subspace, and the smallest singular values are less reliable. Use truncation with oversampling. |
| [in] | power_iteration | number of iterations for the power method: Y = (A * Adag)^power_iteration * A * Tin |
| [in] | seed | the seed for the random generator. [Default] Using device entropy. |
| std::vector< cytnx::UniTensor > cytnx::linalg::Rsvd | ( | const cytnx::UniTensor & | Tin, |
| cytnx_uint64 | keepdim, | ||
| bool | is_U = true, |
||
| bool | is_vT = true, |
||
| cytnx_uint64 | power_iteration = 2, |
||
| unsigned int | seed = random::__static_random_device() |
||
| ) |
Perform randomized Singular-Value decomposition on a UniTensor using the ?gesvd method.
A randomized SVD of a UniTensor Tin; uses randomized isometries. The result will depend on the rowrank of the UniTensor Tin.
| [in] | Tin | a UniTensor, with the correct rowrank set to interpret it as a matrix |
| [in] | keepdim | the number of singular values to keep. Note that this is the size of the subspace, and the smallest singular values are less reliable. Use truncation with oversampling. |
| [in] | is_U | if true, the left-unitary UniTensor U (isometry) is returned. |
| [in] | is_vT | if true, the right-unitary UniTensor vT (isometry) is returned. |
| [in] | power_iteration | number of iterations for the power method: Y = (A * Adag)^power_iteration * A * Tin |
| [in] | seed | the seed for the random generator. [Default] Using device entropy. |
| std::vector< Tensor > cytnx::linalg::Rsvd | ( | const Tensor & | Tin, |
| cytnx_uint64 | keepdim, | ||
| bool | is_U = true, |
||
| bool | is_vT = true, |
||
| cytnx_uint64 | power_iteration = 2, |
||
| unsigned int | seed = random::__static_random_device() |
||
| ) |
Perform randomized Singular-Value decomposition on a rank-2 Tensor (a matrix) using the ?gesvd method.
A randomized SVD of a matrix Tin; uses randomized isometries.
| [in] | keepdim | the number of singular values to keep. Note that this is the size of the subspace, and the smallest singular values are less reliable. Use truncation with oversampling. |
| [in] | is_U | if true, the left-unitary UniTensor U (isometry) is returned. |
| [in] | is_vT | if true, the right-unitary UniTensor vT (isometry) is returned. |
| [in] | power_iteration | number of iterations for the power method: Y = (A * Adag)^power_iteration * A * Tin |
| [in] | seed | the seed for the random generator. [Default] Using device entropy. |
| std::vector< cytnx::UniTensor > cytnx::linalg::Rsvd_truncate | ( | const cytnx::UniTensor & | Tin, |
| cytnx_uint64 | keepdim, | ||
| double | err = 0., |
||
| bool | is_U = true, |
||
| bool | is_vT = true, |
||
| unsigned int | return_err = 0, |
||
| cytnx_uint64 | mindim = 1, |
||
| cytnx_uint64 | oversampling_summand = 10, |
||
| double | oversampling_factor = 1., |
||
| cytnx_uint64 | power_iteration = 0, |
||
| unsigned int | seed = random::__static_random_device() |
||
| ) |
Perform a truncated Singular-Value decomposition of a UniTensor.
This function will perform a truncated Singular-Value decomposition of a UniTensor. It uses the ?gesvd method for the SVD. It will perform the randomized SVD first, and then truncate the singular values to the given cutoff err. That means, given a UniYrndot Tin as \( M \), then the result will be:
\[ M = U S V^\dagger, \]
where \( S
\) is a singular values matrix with the singular values truncated to the given cutoff err. The dimension of \( S \) is at most keepdim.
The truncation order is as following (later constraints might be violated by previous ones):
1) Keep at most keepdim singular values; there might be an exception in case of exact degeneracies, see note below
2) Keep at least mindim singular values;
3) Drop all singular values smaller than err (no normalization applied to the singular values)
| [in] | Tin | a UniTensor, it should be a rank-2 tensor (matrix) |
| [in] | keepdim | the number (at most) of singular values to keep. |
| [in] | err | the cutoff error (the singular values smaller than err will be truncated.) |
| [in] | is_UvT | if true, the left- and right- unitary matrices (isometries) are returned. |
| [in] | return_err | whether the error shall be returned. If return_err is true, then the largest error will be pushed back to the vector (The smallest singular value in the return singular values matrix \( S \).) If return_err is a positive int, then the full list of truncated singular values will be returned. |
| [in] | oversampling_summand | the randomized SVD computes (1 + oversampling_fact) * keepdim + oversampling_summand before further truncating to keepdim |
| [in] | oversampling_fact | the randomized SVD computes (1 + oversampling_fact) * keepdim + oversampling_summand before further truncating to keepdim |
| [in] | power_iteration | number of iterations for the power method: Y = (A * Adag)^power_iteration * A * Tin |
| [in] | seed | the seed for the random generator. [Default] Using device entropy. |
[std::vector<Tensors>]
is_UvT is true, then the tensors \( U \) and \( V^\dagger \) will be pushed back to the vector.return_err is true, then the error will be pushed back to the vector. | std::vector< Tensor > cytnx::linalg::Rsvd_truncate | ( | const Tensor & | Tin, |
| cytnx_uint64 | keepdim, | ||
| double | err = 0., |
||
| bool | is_U = true, |
||
| bool | is_vT = true, |
||
| unsigned int | return_err = 0, |
||
| cytnx_uint64 | mindim = 1, |
||
| cytnx_uint64 | oversampling_summand = 10, |
||
| double | oversampling_factor = 1., |
||
| cytnx_uint64 | power_iteration = 0, |
||
| unsigned int | seed = random::__static_random_device() |
||
| ) |
Perform a truncated Singular-Value decomposition of a rank-2 Tensor (a matrix).
This function will perform a truncated Singular-Value decomposition of a matrix (a rank-2 Tensor). It uses the ?gesvd method for the SVD. It will perform the full SVD first, and then truncate the singular values to the given cutoff err. That means, given a matrix Tin as \( M \), then the result will be:
\[ M = U S V^\dagger, \]
where \( S \) is a singular values matrix with the singular values truncated to the given cutoff err. The dimension of \( S \) is at most keepdim.
The truncation order is as following (later constraints might be violated by previous ones):
1) Keep at most keepdim singular values; there might be an exception in case of exact degeneracies, see note below
2) Keep at least mindim singular values;
3) Drop all singular values smaller than err (no normalization applied to the singular values)
| [in] | Tin | a Tensor, it should be a rank-2 tensor (matrix) |
| [in] | keepdim | the number (at most) of singular values to keep. |
| [in] | err | the cutoff error (the singular values smaller than err will be truncated.) |
| [in] | is_UvT | if true, the left- and right- unitary matrices (isometries) are returned. |
| [in] | return_err | whether the error shall be returned. If return_err is true, then the largest error will be pushed back to the vector (The smallest singular value in the return singular values matrix \( S \).) If return_err is a positive int, then the full list of truncated singular values will be returned. |
| [in] | oversampling_summand | the randomized SVD computes (1 + oversampling_fact) * keepdim + oversampling_summand before further truncating to keepdim |
| [in] | oversampling_fact | the randomized SVD computes (1 + oversampling_fact) * keepdim + oversampling_summand before further truncating to keepdim |
| [in] | power_iteration | number of iterations for the power method: Y = (A * Adag)^power_iteration * A * Tin |
| [in] | seed | the seed for the random generator. [Default] Using device entropy. |
[std::vector<Tensors>]
is_UvT is true, then the tensors \( U \) and \( V^\dagger \) will be pushed back to the vector.return_err is true, then the error will be pushed back to the vector. | cytnx::UniTensor cytnx::linalg::Sub | ( | const cytnx::UniTensor & | Lt, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The subtraction function between two UniTensor.
This is the subtraction function for UniTensor. It will subtract the UniTensor and a template type together. It will subtract every element in the UniTensor with the template type. That means if the UniTensor Lt is given as \( T_L \) and the UniTensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] - T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the UniTensor \( T_L \) and \( T_R \). It will perform the element-wise subtraction and note that it will return a new UniTensor object.
Lt and Rt must have the same shape. | cytnx::UniTensor cytnx::linalg::Sub | ( | const cytnx::UniTensor & | Lt, |
| const T & | rc | ||
| ) |
The subtraction function between a UniTensor and a template type.
This is the subtraction function for UniTensor. It will subtract the UniTensor and a template type together. It will subtract every element in the UniTensor with the template type. That means if the UniTensor Lt is given as \( T_i \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_i[i] - c, \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | Lt | The left UniTensor. |
| [in] | rc | The right template type. supported type: Scalar, cytnx::cytnx_complex128, cytnx::cytnx_complex64, cytnx::cytnx_double, cytnx::cytnx_float, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
rc will be casted to the same type as the UniTensor Lt. | cytnx::UniTensor cytnx::linalg::Sub | ( | const T & | lc, |
| const cytnx::UniTensor & | Rt | ||
| ) |
The subtraction function between a UniTensor and a template type.
This is the subtraction function for UniTensor. It will subtract the UniTensor and a template type together. It will subtract every element in the UniTensor with the template type. That means if the template type lc is given as \( c \) and the UniTensor Rt is given as \( T_i \), then the result will be:
\[ T_o[i] = c - T_i[i], \]
where \( T_i[i] \) and \( T_o[i] \) are the elements in the UniTensor \( T_i \) and \( T_o \).
| [in] | Lt | The left UniTensor. |
| [in] | rc | The right template type. supported type: Scalar, cytnx::cytnx_complex128, cytnx::cytnx_complex64, cytnx::cytnx_double, cytnx::cytnx_float, cytnx::cytnx_int64, cytnx::cytnx_uint64, cytnx::cytnx_int32, cytnx::cytnx_uint32, cytnx::cytnx_int16, cytnx::cytnx_uint16, cytnx::cytnx_bool. |
lc will be casted to the same type as the UniTensor Rt. The subtraction function for Tensor.
This is the subtraction function between a Tensor and a template type. It will perform the element-wise subtraction. That means if the template type lc is given as \( c \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = c - T_R[i], \]
where \( T_R[i] \) is the elements in the Tensor \( T_R \). It will perform the element-wise subtraction and note that it will return a new Tensor object.
| [in] | lc | The left template type. |
| [in] | Rt | The right Tensor. |
The subtraction function for Tensor.
This is the subtraction function between a Tensor and a template type. It will perform the element-wise subtraction. That means if the left Tensor Lt is given as \( T_L \) and the template type rc is given as \( c \), then the result will be:
\[ T_o[i] = T_L[i] - c, \]
where \( T_L[i] \) is the elements in the Tensor \( T_L \). It will perform the element-wise subtraction and note that it will return a new Tensor object.
| [in] | Lt | The left Tensor. |
| [in] | rc | The right template type. |
The subtraction function for Tensor.
This is the subtraction function between two Tensor. It will perform the element-wise subtraction. That means if the left Tensor Lt is given as \( T_L \) and the right Tensor Rt is given as \( T_R \), then the result will be:
\[ T_o[i] = T_L[i] - T_R[i], \]
where \( T_L[i] \) and \( T_R[i] \) are the elements in the Tensor \( T_L \) and \( T_R \). It will perform the element-wise subtraction and note that it will return a new Tensor object.
get the sum of all the elements.
| [in] | Tn | a cytnx::Tensor |
| std::vector< cytnx::UniTensor > cytnx::linalg::Svd | ( | const cytnx::UniTensor & | Tin, |
| const bool & | is_UvT = true |
||
| ) |
Perform Singular-Value decomposition on a UniTensor using divide-and-conquer method.
This function performs the Singular-Value decomposition on a UniTensor Tin. The result will depend on the rowrank of the UniTensor Tin. For more details, please refer to the documentation of the function Svd(const Tensor &Tin, const bool &is_UvT).
Perform Singular-Value decomposition on a rank-2 Tensor (a matrix).
This function will perform Singular-Value decomposition on a matrix (a rank-2 Tensor). That means given a matrix Tin as \( M \), then the result will be:
\[ M = U S V^\dagger, \]
where \( U \) is a left uniform matrix, \( S \) is a diagonal matrix with singular values, and \( V^\dagger \) is the conjugate transpose of the right uniform matrix \( V \). Furthermore, \( U \) and \( V \) are unitary matrices, and \( S \) is a non-negative diagonal matrix.
| [in] | Tin | a Tensor, it should be a rank-2 tensor (matrix) |
| [in] | is_UvT | whether need to return a left unitary matrix. |
[std::vector<Tensors>]
is_UvT is true, then the tensors \( U,V^\dagger \) will be pushed back to the vector. | std::vector< cytnx::UniTensor > cytnx::linalg::Svd_truncate | ( | const cytnx::UniTensor & | Tin, |
| const cytnx_uint64 & | keepdim, | ||
| const double & | err = 0., |
||
| const bool & | is_UvT = true, |
||
| const unsigned int & | return_err = 0, |
||
| const cytnx_uint64 & | mindim = 1 |
||
| ) |
Perform Singular-Value decomposition on a UniTensor with truncation.
This function performs the Singular-Value decomposition of a UniTensor Tin and truncates the singular values. The result will depend on the rowrank of the UniTensor Tin. For more details, please refer to the references below.
| std::vector< cytnx::UniTensor > cytnx::linalg::Svd_truncate | ( | const cytnx::UniTensor & | Tin, |
| const cytnx_uint64 & | keepdim, | ||
| std::vector< cytnx_uint64 > | min_blockdim, | ||
| const double & | err = 0., |
||
| const bool & | is_UvT = true, |
||
| const unsigned int & | return_err = 0, |
||
| const cytnx_uint64 & | mindim = 1 |
||
| ) |
Perform Singular-Value decomposition on a UniTensor with truncation and keep at most min_blockdim singular values in each block.
This function performs the Singular-Value decomposition of a UniTensor Tin and truncates the singular values. The result will depend on the rowrank of the UniTensor Tin. For each block, the minimum dimension can be chosen. This can be helpful to avoid loosing symmetry sectors in the truncated SVD. For more details, please refer to the documentation of the function Svd_truncate(const Tensor &Tin, const cytnx_uint64
&keepdim, const double &err, const bool &is_UvT, const unsigned int &return_err).
The truncation order is as following (later constraints might be violated by previous ones):
1) Keep the largest min_blockdim singular values in each block; reduce keepdim and mindim by the number of already kept singular values
2) Keep at most keepdim singular values; there might be an exception in case of exact degeneracies where more singular values are kept
3) Keep at least mindim singular values;
4) Drop all singular values smaller than err (no normalization applied to the singular values)
| [in] | Tin | a BlockUniTensor, with the correct rowrank set to interpret it as a matrix |
| [in] | keepdim | the number (at most) of singular values to keep. |
| [in] | min_blockdim | a vector containing the minimum dimension of each block; alternatively, a vector with only one element can be given to have the same min_blockdim for each block |
| [in] | err | the cutoff error (the singular values smaller than err will be truncated.) |
| [in] | is_UvT | if true, the left- and right- unitary UniTensors (isometries) are returned. |
| [in] | return_err | whether the error shall be returned. If return_err is true, then the largest error will be pushed back to the vector (The smallest singular value in the return singular values matrix \( S \).) If return_err is a positive int, then the full list of truncated singular values will be returned. |
[std::vector<Tensors>]
Tin. keepdim for degenerate singular values: if the largest \( n \) truncated singular values would be exactly equal to the smallest kept singular value, then the bond dimension is enlarged to keepdim \( + n \). Example: if the singular values are (1 2 2 2 2 3) and keepdim = 3, then the bond dimension will be 5 in order to keep all the degenerate singular values. | std::vector< Tensor > cytnx::linalg::Svd_truncate | ( | const Tensor & | Tin, |
| const cytnx_uint64 & | keepdim, | ||
| const double & | err = 0., |
||
| const bool & | is_UvT = true, |
||
| const unsigned int & | return_err = 0, |
||
| const cytnx_uint64 & | mindim = 1 |
||
| ) |
Perform a truncated Singular-Value decomposition of a rank-2 Tensor (a matrix).
This function will perform a truncated Singular-Value decomposition of a matrix (a rank-2 Tensor). It will perform the full SVD first, and then truncate the singular values to the given cutoff err. That means, given a matrix Tin as \( M \), then the result will be:
\[ M = U S V^\dagger, \]
where \( S \) is a singular values matrix with the singular values truncated to the given cutoff err. The dimension of \( S \) is at most keepdim.
The truncation order is as following (later constraints might be violated by previous ones):
1) Keep at most keepdim singular values; there might be an exception in case of exact degeneracies, see note below
2) Keep at least mindim singular values;
3) Drop all singular values smaller than err (no normalization applied to the singular values)
| [in] | Tin | a Tensor, it should be a rank-2 tensor (matrix) |
| [in] | keepdim | the number (at most) of singular values to keep. |
| [in] | err | the cutoff error (the singular values smaller than err will be truncated.) |
| [in] | is_UvT | if true, the left- and right- unitary matrices (isometries) are returned. |
| [in] | return_err | whether the error shall be returned. If return_err is true, then the largest error will be pushed back to the vector (The smallest singular value in the return singular values matrix \( S \).) If return_err is a positive int, then the full list of truncated singular values will be returned. |
[std::vector<Tensors>]
is_UvT is true, then the tensors \( U \) and \( V^\dagger \) will be pushed back to the vector.return_err is true, then the error will be pushed back to the vector. keepdim for degenerate singular values: if the largest \( n \) truncated singular values would be exactly equal to the smallest kept singular value, then the bond dimension is enlarged to keepdim \( + n \). Example: if the singular values are (1 2 2 2 2 3) and keepdim = 3, then the bond dimension will be 5 in order to keep all the degenerate singular values. | Tensor cytnx::linalg::Tensordot | ( | const Tensor & | Tl, |
| const Tensor & | Tr, | ||
| const std::vector< cytnx_uint64 > & | idxl, | ||
| const std::vector< cytnx_uint64 > & | idxr, | ||
| const bool & | cacheL = false, |
||
| const bool & | cacheR = false |
||
| ) |
perform tensor dot by sum out the indices assigned of two Tensors.
| [in] | Tl | Tensor #1 |
| [in] | Tr | Tensor #2 |
| [in] | idxl | the indices of rank of Tensor #1 that is going to sum with Tensor #2 |
| [in] | idxr | the indices of rank of Tensor #2 that is going to sum with Tensor #1 |
| [in] | cacheL | cache Tensor #1 (See user-guide for details) |
| [in] | cacheR | cache Tensor #2 (See user-guide for details) |
| Tensor cytnx::linalg::Tensordot_dg | ( | const Tensor & | Tl, |
| const Tensor & | Tr, | ||
| const std::vector< cytnx_uint64 > & | idxl, | ||
| const std::vector< cytnx_uint64 > & | idxr, | ||
| const bool & | diag_L | ||
| ) |
perform tensor dot by sum out the indices assigned of two Tensors, with either one of them to be a rank-2 diagonal tensor represented by a rank-2 tensor.
| [in] | Tl | Tensor #1 |
| [in] | Tr | Tensor #2 |
| [in] | idxl | the indices of rank of Tensor #1 that is going to sum with Tensor #2 |
| [in] | idxr | the indices of rank of Tensor #2 that is going to sum with Tensor #1 |
| [in] | diag_L | if Tl(true)/Tr(false) is a diagnal matrix, represented by a rank-1 tensor. |
| cytnx::UniTensor cytnx::linalg::Trace | ( | const cytnx::UniTensor & | Tin, |
| const cytnx_int64 & | a = 0, |
||
| const cytnx_int64 & | b = 1 |
||
| ) |
| cytnx::UniTensor cytnx::linalg::Trace | ( | const cytnx::UniTensor & | Tin, |
| const std::string & | a, | ||
| const std::string & | b | ||
| ) |
| Tensor cytnx::linalg::Trace | ( | const Tensor & | Tn, |
| const cytnx_uint64 & | axisA = 0, |
||
| const cytnx_uint64 & | axisB = 1 |
||
| ) |
perform trace over index.
This function will perform trace over index axisA and axisB. For example, if Tn is a rank-4 tensor \( T \), then axisA = 0 and axisB = 2, then the result will be:
\[ \mathrm{Tr}_{j,l}(T) = \sum_{i,k} T_{i,j,k,l} \]
| [in] | Tn | a Tensor |
| [in] | axisA | the first index to perform trace. |
| [in] | axisB | the second index to perform trace. |
Tn should be at-least rank-2 Tensor. | std::vector< Tensor > cytnx::linalg::Tridiag | ( | const Tensor & | Diag, |
| const Tensor & | Sub_diag, | ||
| const bool & | is_V = true, |
||
| const bool & | is_row = false, |
||
| bool | throw_excp = false |
||
| ) |
perform diagonalization of symmetric tri-diagnoal matrix.
| [in] | Diag | Tensor #1 |
| [in] | Sub_diag | Tensor #2 |
| [in] | is_V | if calculate the eigen value. |
| [in] | k | Return k lowest eigen vector if is_V=True |
| [in] | throw_excp | Whether to throw exception when error occurs in Tridiag internal function |
| Tensor cytnx::linalg::Vectordot | ( | const Tensor & | Tl, |
| const Tensor & | Tr, | ||
| const bool & | is_conj = false |
||
| ) |
perform inner product of vectors
| [in] | Tl | Tensor #1 |
| [in] | Tr | Tensor #2 |
| [in] | if | the Tl should be conjugated (only work for complex. For real Tensor, no function), default: false |