Cytnx v1.0.0
Loading...
Searching...
No Matches
Functions
cytnx::linalg Namespace Reference

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::UniTensorSvd (const cytnx::UniTensor &Tin, const bool &is_UvT=true)
 Perform Singular-Value decomposition on a UniTensor using divide-and-conquer method.
 
std::vector< cytnx::UniTensorGesvd (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::UniTensorSvd_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::UniTensorSvd_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::UniTensorGesvd_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::UniTensorGesvd_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::UniTensorHosvd (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::UniTensorQr (const cytnx::UniTensor &Tin, const bool &is_tau=false)
 Perform the QR decomposition on a UniTensor.
 
std::vector< cytnx::UniTensorQdr (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.
 
std::vector< TensorSvd (const Tensor &Tin, const bool &is_UvT=true)
 Perform Singular-Value decomposition on a rank-2 Tensor (a matrix).
 
std::vector< TensorGesvd (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< TensorSvd_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< TensorGesvd_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).
 
std::vector< TensorHosvd (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< TensorQr (const Tensor &Tin, const bool &is_tau=false)
 Perform QR decomposition on a rank-2 Tensor.
 
std::vector< TensorQdr (const Tensor &Tin, const bool &is_tau=false)
 Perform QDR decomposition on a rank-2 Tensor.
 
std::vector< TensorEigh (const Tensor &Tin, const bool &is_V=true, const bool &row_v=false)
 eigen-value decomposition for Hermitian matrix
 
std::vector< UniTensorEigh (const cytnx::UniTensor &Tin, const bool &is_V=true, const bool &row_v=false)
 
std::vector< TensorEig (const Tensor &Tin, const bool &is_V=true, const bool &row_v=false)
 eigen-value decomposition for generic square matrix
 
std::vector< UniTensorEig (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< TensorTridiag (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< TensorArnoldi (LinOp *Hop, const Tensor &Tin=Tensor(), 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 bool &verbose=false)
 perform Arnoldi for matrices or linear function.
 
std::vector< UniTensorArnoldi (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 bool &verbose=false)
 perform Arnoldi for matrices or linear function.
 
std::vector< TensorLanczos (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< UniTensorLanczos (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< TensorLanczos_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< TensorLanczos_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< UniTensorLanczos_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< TensorLstsq (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.
 

Detailed Description

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,

  1. If the object is Tensor, then it need to be rank-2.
  2. If the object is UniTensor, then the result will depend on the UniTensor's rowrank.

Function Documentation

◆ Abs()

Tensor cytnx::linalg::Abs ( const Tensor Tin)

Elementwise absolute value.

This function will perform Elementwise absolute value on all the elements in Tensor Tin. That is, the output will be:

\[ T_{o}[i] = |T_{i}[i]| \]

Parameters
[in]Tintensor.
Returns
[Tensor]

◆ Abs_()

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.

Parameters
[in]Tin,theinput Tensor.
Note
on return, the elements in Tin will be modified to it's absolute value. Note that if the input tensor is complex, it will be modified to real type.

◆ Add() [1/6]

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.

Parameters
[in]LtThe left UniTensor.
[in]RtThe right UniTensor.
Returns
The result UniTensor.
Precondition
Lt and Rt must have the same shape.
See also
UniTensor::Add(const cytnx::UniTensor &Rt) const, operator+(const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt)

◆ Add() [2/6]

template<class T >
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 \).

Parameters
[in]LtThe 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.

Returns
The result UniTensor.
Precondition
The supported template type shown above.
Note
The inpute template type rc will be casted to the same type as the UniTensor Lt.
See also
operator+(const cytnx::UniTensor &Lt, const T &rc), Add(const T &lc, const cytnx::UniTensor &Rt)

◆ Add() [3/6]

template<class T >
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 \).

Parameters
[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]RtThe right UniTensor.
Returns
The result UniTensor.
Precondition
The supported template type shown above.
Note
The inpute template type lc will be casted to the same type as the UniTensor Rt.
See also
operator+(const T &lc, const cytnx::UniTensor &Rt), Add(const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt)

◆ Add() [4/6]

template<class T >
Tensor cytnx::linalg::Add ( const T &  lc,
const Tensor 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.

Parameters
[in]lcThe left template type.
[in]RtThe right Tensor.
Returns
The result Tensor.
See also
Add(const Tensor &Lt, const Tensor &Rt), Add(const Tensor &Lt, const T &rc), iAdd(Tensor &Lt, const Tensor &Rt), operator+(const Tensor &Lt, const Tensor &Rt)

◆ Add() [5/6]

template<class T >
Tensor cytnx::linalg::Add ( const Tensor Lt,
const T &  rc 
)

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.

Parameters
[in]LtThe left Tensor.
[in]rcThe right template type.
Returns
The result Tensor.
See also
Add(const Tensor &Lt, const Tensor &Rt), Add(const T &lc, const Tensor &Rt), iAdd(Tensor &Lt, const Tensor &Rt), operator+(const Tensor &Lt, const Tensor &Rt)

◆ Add() [6/6]

Tensor cytnx::linalg::Add ( const Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in]LtThe left Tensor.
[in]RtThe right Tensor.
Returns
The result Tensor.
Precondition
The shape of Lt and Rt must be the same.
See also
Add(const T &lc, const Tensor &Rt), Add(const Tensor &Lt, const T &rc), iAdd(Tensor &Lt, const Tensor &Rt), operator+(const Tensor &Lt, const Tensor &Rt)

◆ Arnoldi() [1/2]

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 bool &  verbose = false 
)

perform Arnoldi for matrices or linear function.

This function calculate the eigen value problem using Arnoldi algorithm.

Parameters
[in]Hopthe Linear Operator defined by LinOp class or it's inheritance (see LinOp).
[in]Tinthe initial UniTensor.
[in]which

which order eigenvlues and corresponding eigenvectors should be find, the supported options are:

'LM' : largest magnitude 'LR' : largest real part 'LI' : largest imaginary part 'SR' : smallest real part 'SI' : smallest imaginary part

[in]maxiterthe maximum interation steps for each k.
[in]cvg_critthe convergence criterion of the energy.
[in]kthe number of lowest k eigen values.
[in]is_Vif set to true, the eigen vectors will be returned.
[in]verboseprint out iteration info.
Returns
[eigvals (UniTensor), eigvec_1, eivec_2, ..., eigvec_k]. The first UniTensor contains eigenvalues.
Note
To use, define a linear operator with LinOp class either by assign a custom function or create a class that inherit LinOp (see LinOp for further details)
Precondition
  1. The initial UniTensor cannot be empty.
  2. The UniTensor version of the Arnoldi not support which = 'SM'.

◆ Arnoldi() [2/2]

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 = 1.0e-9,
const cytnx_uint64 &  k = 1,
const bool &  is_V = true,
const bool &  verbose = false 
)

perform Arnoldi for matrices or linear function.

This function calculate the eigen value problem using Arnoldi algorithm.

Parameters
[in]Hopthe Linear Operator defined by LinOp class or it's inheritance (see LinOp).
[in]Tinthe initial vector, this should be rank-1.
[in]which

which order eigenvlues and corresponding eigenvectors should be find, the supported options are:

'LM' : largest magnitude 'LR' : largest real part 'LI' : largest imaginary part 'SM' : smallest magnitude 'SR' : smallest real part 'SI' : smallest imaginary part

[in]maxiterthe maximum interation steps for each k.
[in]cvg_critthe convergence criterion of the energy.
[in]kthe number of lowest k eigen values.
[in]is_Vif set to true, the eigen vectors will be returned.
[in]verboseprint out iteration info.
Returns
[eigvals (Tensor), eigvecs (Tensor)(option)]
Note
To use, define a linear operator with LinOp class either by assign a custom function or create a class that inherit LinOp (see LinOp for further details)

◆ Axpy()

Tensor cytnx::linalg::Axpy ( const Scalar &  a,
const Tensor x,
const Tensor y = Tensor() 
)

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.

Parameters
[in]aScalar.
[in]xTensor, can be any rank
[in]yTensor, can be any rank
Returns
[Tensor] If \( \textbf{y} \) is not specify, then it performs \( a\textbf{x} \) -> return
Note
This will return a new tensor.

◆ Axpy_()

void cytnx::linalg::Axpy_ ( const Scalar &  a,
const Tensor x,
Tensor y 
)

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.

Parameters
[in]aScalar.
[in]xTensor, can be any rank
[in]yTensor, can be any rank
Returns
[Tensor] If \( \textbf{y} \) is not specify, then it performs \( a\textbf{x} \) -> return
Note
Compared to Axpy(const Scalar &a, const Tensor &x, const Tensor &y = Tensor()), this function will perform inplacely.

◆ Conj() [1/2]

cytnx::UniTensor cytnx::linalg::Conj ( const cytnx::UniTensor UT)

Elementwise conjugate of the UniTensor.

Parameters
[in]UTThe input UniTensor.
Returns
[UniTensor] The UniTensor with all element being conjugated
See also
See UniTensor.Conj() for further details

◆ Conj() [2/2]

Tensor cytnx::linalg::Conj ( const Tensor Tin)

Conjugate all the element in Tensor.

This function take the complex conjugate of all the elements in Tensor Tin. That is, the output will be:

\[ T_{o}[i] = T_{in}[i]^* \]

Furthermore,

  1. if the input Tensor is complex, then return a new Tensor with all the elements are conjugated.
  2. if the input Tensor is real, then return a copy of input Tensor.
    Parameters
    [in]Tina Tensor
    Returns
    [Tensor]

◆ Conj_() [1/2]

void cytnx::linalg::Conj_ ( cytnx::UniTensor UT)

Inplace elementwise conjugate of the UniTensor.

Parameters
[in]UTThe input UniTensor.
See also
See UniTensor.Conj_() for further details

◆ Conj_() [2/2]

void cytnx::linalg::Conj_ ( Tensor Tin)

inplace perform Conjugate on all the element in Tensor.

This function take the complex conjugate of all the elements in Tensor Tin. This function is just a inplace version of Conj.

Note
  1. if the input Tensor is complex, the elements of input Tensor will all be conjugated.
  2. if the input Tensor is real, then nothing act.

◆ Cpr() [1/3]

template<class T >
Tensor cytnx::linalg::Cpr ( const T &  lc,
const Tensor Rt 
)

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.

Parameters
[in]lcThe left template type.
[in]RtThe right Tensor.
Returns
The result Tensor.
See also
Cpr(const Tensor &Lt, const Tensor &Rt), Cpr(const Tensor &Lt, const T &rc)

◆ Cpr() [2/3]

template<class T >
Tensor cytnx::linalg::Cpr ( const Tensor Lt,
const T &  rc 
)

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.

Parameters
[in]LtThe left Tensor.
[in]rcThe right template type.
Returns
The result Tensor.
See also
Cpr(const Tensor &Lt, const Tensor &Rt), Cpr(const T &lc, const Tensor &Rt)

◆ Cpr() [3/3]

Tensor cytnx::linalg::Cpr ( const Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in]LtThe left Tensor.
[in]RtThe right Tensor.
Returns
The result Tensor.
Precondition
The input tensors Lt and Rt should have the same shape.
See also
Cpr(const T &lc, const Tensor &Rt), Cpr(const Tensor &Lt, const T &rc)

◆ Det()

Tensor cytnx::linalg::Det ( const Tensor Tl)

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).

Parameters
[in]Tlinput a Tensor with shape (N,N)
Returns
Tensor
Precondition
the input tensor should be a rank-2 tensor with shape (N,N). (a square matrix)

◆ Diag()

Tensor cytnx::linalg::Diag ( const Tensor Tin)

return a diagonal tensor with diagonal elements provided as Tin.

Returns
[Tensor]

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.

Precondition
Tin should be a rank-2 Tensor.

◆ Directsum()

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 out put 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 ...)

Parameters
[in]T1rank-n Tensor #1
[in]T2rank-n Tensor #2
[in]shared_axesThe axes that are shared by two tensors
Returns
[Tensor]
Precondition
two tensor should on same device.

◆ Div() [1/6]

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.

Parameters
[in]LtThe left UniTensor.
[in]RtThe right UniTensor.
Returns
The result UniTensor.
Precondition
Lt and Rt must have the same shape.
See also
UniTensor::Div(const cytnx::UniTensor &Rt) const, operator/(const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt)

◆ Div() [2/6]

template<class T >
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 \).

Parameters
[in]LtThe 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.

Returns
The result UniTensor.
Precondition
The supported template type shown above.
Note
  1. The inpute template type rc will be casted to the same type as the UniTensor Lt.
  2. The division by zero is not allowed.
See also
operator/(const cytnx::UniTensor &Lt, const T &rc), Div(const cytnx::UniTensor &Lt, const T &rc)

◆ Div() [3/6]

template<class T >
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 \).

Parameters
[in]LtThe 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.

Returns
The result UniTensor.
Precondition
The supported template type shown above.
Note
  1. The inpute template type lc will be casted to the same type as the UniTensor Rt.
  2. The division by zero is not allowed.
See also
operator/(const T &lc, const cytnx::UniTensor &Rt), Div(const T &lc, const cytnx::UniTensor &Rt)

◆ Div() [4/6]

template<class T >
Tensor cytnx::linalg::Div ( const T &  lc,
const Tensor 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.

Parameters
[in]lcThe left template type.
[in]RtThe right Tensor.
Returns
The result Tensor.
Precondition
the right tensor Rt should not contain any zero element.
See also
Div(const Tensor &Lt, const Tensor &Rt), Div(const Tensor &Lt, const T &rc), iDiv(Tensor &Lt, const Tensor &Rt), operator/(const Tensor &Lt, const Tensor &Rt)

◆ Div() [5/6]

template<class T >
Tensor cytnx::linalg::Div ( const Tensor Lt,
const T &  rc 
)

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.

Parameters
[in]LtThe left Tensor.
[in]rcThe right template type.
Returns
The result Tensor.
Precondition
the right template type rc should not be zero.
See also
Div(const Tensor &Lt, const Tensor &Rt), Div(const T &lc, const Tensor &Rt), iDiv(Tensor &Lt, const Tensor &Rt), operator/(const Tensor &Lt, const Tensor &Rt)

◆ Div() [6/6]

Tensor cytnx::linalg::Div ( const Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in]LtThe left Tensor.
[in]RtThe right Tensor.
Returns
The result Tensor.
Precondition
the right Tensor Rt should not contain any zero element.
See also
Div(const T &lc, const Tensor &Rt), Div(const Tensor &Lt, const T &rc), iDiv(Tensor &Lt, const Tensor &Rt), operator/(const Tensor &Lt, const Tensor &Rt)

◆ Dot()

Tensor cytnx::linalg::Dot ( const Tensor Tl,
const Tensor Tr 
)

dot product of two arrays.

  1. if both Tl and Tr are 1d arrays, it is inner product of vectors (no complex conj), it calls linalg.Vectordot with is_conj=false.
  2. if both Tl and Tr are 2d arrays, it calls linalg.Matmul to compute the matrix multiplication
  3. if Tl is Nd array (with N>=2, and Tr is 1-D array, it is sum product over the last axis of a with b
    Parameters
    [in]TlTensor #1
    [in]TrTensor #2
    Returns
    [Tensor]
    Note
    performance tune: This function have better performance when two arrays with same types, and are one of following type: cytnx_double, cytnx_float, cytnx_complex64 or cytnx_complex128. [Python] In Python API, operator@ is overloaded as a shorthand of linalg::Dot.

◆ Eig() [1/2]

std::vector< UniTensor > cytnx::linalg::Eig ( const cytnx::UniTensor Tin,
const bool &  is_V = true,
const bool &  row_v = false 
)

◆ Eig() [2/2]

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.

Parameters
[in]TinA square matrix (a rank-2 Tensor)
[in]is_Vwhether need to return the eigen vectors.
[in]row_Vif set to ture, the return eigen vectors will be row form.
Returns
[std::vector<Tensors>], the first tensor is the eigen values, a 1-d tensor (vector). The second tensor is the eigen vectors, a 2-d tensor (matrix). This tensor will only return when is_V = true. Furthermore, if row_V = true, then the eigen vectors will be row form. Otherwise, the eigen vectors will be column form.
Precondition
the Tin should be a square matrix.

◆ Eigh() [1/2]

std::vector< UniTensor > cytnx::linalg::Eigh ( const cytnx::UniTensor Tin,
const bool &  is_V = true,
const bool &  row_v = false 
)

◆ Eigh() [2/2]

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.

Parameters
[in]Tina Tensor , it should be a rank-2 tensor (matrix)
[in]is_Vwhether need to return the eigen vectors.
[in]row_Vif set to ture, the return eigen vectors will be row form.
Returns
[std::vector<Tensors>], the first tensor is the eigen values, a 1-d tensor (vector). The second tensor is the eigen vectors, a 2-d tensor (matrix). This tensor will only return when is_V = true. Furthermore, if row_V = true, then the eigen vectors will be row form. Otherwise, the eigen vectors will be column form.
Precondition
the Tin should be a Hermitian matrix.
Warning
If 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.

◆ Exp()

Tensor cytnx::linalg::Exp ( const Tensor Tin)

Exponential all the element in Tensor.

This function will perform Exponential on all the elements in Tensor Tin. That is, the output will be:

\[ T_{o}[i] = e^{T_{i}[i]} \]

Note that it will cast to Double type or ComplexDouble type.

Parameters
[in]Tina Tensor
Returns
[Double Tensor] or [ComplexDouble Tensor]

◆ Exp_()

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,

  1. on return, the elements in Tin will be modified to it's exponetial value.
  2. For Real, if the type is not Double, change the type of the input tensor to Double.
  3. For Complex, if input is ComplexFloat, promote to ComplexDouble.
    Parameters
    [in]Tin,theinput Tensor.

◆ Expf()

Tensor cytnx::linalg::Expf ( const Tensor Tin)

Exponential all the element in Tensor.

This function will perform Exponential on all the elements in Tensor Tin. That is, the output will be:

\[ T_{o}[i] = e^{T_{i}[i]} \]

Note that it will cast to Float type or ComplexFloat type.

Parameters
[in]Tina Tensor
Returns
[Float Tensor] or [ComplexFloat Tensor]

◆ Expf_()

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,

  1. on return, the elements in Tin will be modified to it's exponetial value.
  2. For Real, if the type is not Float, change the type of the input tensor to Float.
  3. For Complex, if input is ComplexDouble, promote to ComplexFloat.
    Parameters
    [in]Tin,theinput Tensor.

◆ ExpH() [1/4]

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)

See also
ExpH(const Tensor &Tin)

◆ ExpH() [2/4]

template<typename T >
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).

See also
ExpH(const Tensor &Tin, const T &a, const T &b)

◆ ExpH() [3/4]

Tensor cytnx::linalg::ExpH ( const Tensor in)

perform matrix exponential for Hermitian matrix

This function perform matrix exponential for Hermitian matrix, That is,

\[ O = \exp{M} \]

Precondition
the in should be a Hermitian matrix.
Warning
If 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).
See also
ExpH(const Tensor &in, const T &a, const T &b)

◆ ExpH() [4/4]

template<typename T >
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)} \]

Parameters
[in]ininput Tensor, should be Hermitian
[in]arescale factor
[in]bbias
Returns
[Tensor]
Precondition
the in should be a Hermitian matrix.
Warning
If 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).

◆ ExpM() [1/4]

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)

See also
ExpM(const Tensor &Tin)

◆ ExpM() [2/4]

template<typename T >
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).

See also
ExpM(const Tensor &Tin, const T &a, const T &b)

◆ ExpM() [3/4]

Tensor cytnx::linalg::ExpM ( const Tensor in)

perform matrix exponential for generic matrix

This function perform matrix exponential for generic matrix, That is,

\[ O = \exp{M} \]

Parameters
[in]ininput Tensor, should be a square rank-2.
Returns
[Tensor]
See also
ExpM(const Tensor &in, const T &a, const T &b)

◆ ExpM() [4/4]

template<typename T >
Tensor cytnx::linalg::ExpM ( const Tensor in,
const T &  a,
const T &  b = 0 
)

perform matrix exponential for generic matrix

This function perform matrix exponential for generic matrix, That is,

\[ O = \exp{(aM + b)} \]

Parameters
[in]ininput Tensor, should be a square rank-2.
[in]arescale factor
[in]bbias
Precondition
the in should be a square matrix.
Returns
[Tensor]

◆ Gemm()

Tensor cytnx::linalg::Gemm ( const Scalar &  a,
const Tensor x,
const Tensor y 
)

Blas Gemm, performing \( a\textbf{x}\textbf{y} -> \) return.

This function performs

\[ return = a\textbf{x}\textbf{y}, \]

where \( \textbf{x},\textbf{y} \) are rank-2 Tensor and a is Scalar. The dtype of return Tensor will be the strongest among x, y and a.

Parameters
[in]aScalar.
[in]xTensor, rank-2 with shape (M,N)
[in]yTensor, rank-2 with shape (N,K)
Returns
[Tensor] with shape (M,K)

◆ Gemm_()

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.

Parameters
[in]aScalar.
[in]xTensor, rank-2 with shape (M,N)
[in]yTensor, rank-2 with shape (N,K)
[in]bScalar.
[in,out]cTensor, rank-2 with shape (M,K), must be properly initialized with the correct shape.

◆ Gemm_Batch()

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.

Parameters
[in]m_arrayarray of cytnx_int64, each element is the number of rows of a_tensors
[in]n_arrayarray of cytnx_int64, each element is the number of columns of b_tensors
[in]k_arrayarray of cytnx_int64, each element is the number of columns of a_tensors and the number of rows of b_tensors
[in]alpha_arrayarray of Scalar, each element is the scalar alpha
[in]a_tensorsarray of Tensor, each element is a rank-2 Tensor with shape (m_array[i],k_array[i])
[in]b_tensorsarray of Tensor, each element is a rank-2 Tensor with shape (k_array[i],n_array[i])
[in]beta_arrayarray of Scalar, each element is the scalar beta
[in,out]c_tensorsarray 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_countcytnx_int64, the number of groups
[in]group_sizearray of cytnx_int64, each element is the number of matrices in each group

◆ Ger()

Tensor cytnx::linalg::Ger ( const Tensor x,
const Tensor y,
const Scalar &  a = Scalar() 
)

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.

Parameters
[in]xTensor, rank-1 with size nx
[in]yTensor, rank-1 with size ny
[in]aScalar, if not provided a = 1.
Returns
[Tensor with shape (nx,ny)]
Note
This will return a new tensor.

◆ Gesvd() [1/2]

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).

See also
Gesvd(const Tensor &Tin, const bool &is_U, const bool &is_vT), Svd(const Tensor &Tin, const bool &is_UvT).

◆ Gesvd() [2/2]

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.

Parameters
[in]Tina Tensor, it should be a rank-2 tensor (matrix)
[in]is_Uwhether need to return left unitary matrix.
[in]is_vTwhether need to return right unitary matrix
Returns

[std::vector<Tensors>]

  1. The first tensor is a 1-d tensor contaning the singular values
  2. If 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.
Precondition
The input tensor should be a rank-2 tensor (matrix).
See also
Gesvd_truncate(const Tensor &Tin, const cytnx_uint64 &keepdim, const double &err, const bool &is_U, const bool &is_vT, const unsigned int &return_err, const cytnx_uint64& mindim)

◆ Gesvd_truncate() [1/3]

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.

See also
Svd_truncate(const cytnx::UniTensor &Tin, const cytnx_uint64 &keepdim, const std::vector<cytnx_uint64> min_blockdim, const double &err, const bool &is_UvT, const unsigned int &return_err, const cytnx_uint64 &mindim), Gesvd(const cytnx::UniTensor &Tin, const bool &is_U, const bool &is_vT)

◆ Gesvd_truncate() [2/3]

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.

See also
Svd_truncate(const cytnx::UniTensor &Tin, const cytnx_uint64 &keepdim, const std::vector<cytnx_uint64> min_blockdim, const double &err, const bool &is_UvT, const unsigned int &return_err, const cytnx_uint64 &mindim), Gesvd(const cytnx::UniTensor &Tin, const bool &is_U, const bool &is_vT)

◆ Gesvd_truncate() [3/3]

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).

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)

Parameters
[in]Tina Tensor, it should be a rank-2 tensor (matrix)
[in]keepdimthe number (at most) of singular values to keep.
[in]errthe cutoff error (the singular values smaller than err will be truncated.)
[in]is_UvTif true, the left- and right- unitary matrices (isometries) are returned.
[in]return_errwhether 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.
Returns

[std::vector<Tensors>]

  1. The first tensor is a 1-d tensor containing the singular values
  2. If is_UvT is true, then the tensors \( U \) and \( V^\dagger \) will be pushed back to the vector.
  3. If return_err is true, then the error will be pushed back to the vector.
Precondition
The input tensor should be a rank-2 tensor (matrix).
See also
Gesvd(const Tensor &Tin, const bool &is_U, const bool &is_vT)
Note
The truncated bond dimension can be larger than 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.

◆ Hosvd() [1/2]

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 >() 
)

◆ Hosvd() [2/2]

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 >() 
)

◆ iAdd()

void cytnx::linalg::iAdd ( Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in,out]LtThe left Tensor.
[in]RtThe right Tensor.
Precondition
The shape of Lt and Rt must be the same.
Note
Compare to the function Add(const Tensor &Lt, const Tensor &Rt), this is a inplace function and it will modify the left Tensor Lt.
See also
Add(const Tensor &Lt, const Tensor &Rt), Add(const T &lc, const Tensor &Rt), Add(const Tensor &Lt, const T &rc), operator+(const Tensor &Lt, const Tensor &Rt)

◆ iDiv()

void cytnx::linalg::iDiv ( Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in,out]LtThe left Tensor.
[in]RtThe right Tensor.
Precondition
the right Tensor Rt should not contain any zero element.
Note
compare to the Div(const Tensor &Lt, const Tensor &Rt) function, this is a inplace function, which will modify the left Tensor Lt.
See also
Div(const Tensor &Lt, const Tensor &Rt), Div(const T &lc, const Tensor &Rt), Div(const Tensor &Lt, const T &rc), operator/(const Tensor &Lt, const Tensor &Rt)

◆ iMul()

void cytnx::linalg::iMul ( Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in,out]LtThe left Tensor.
[in]RtThe right Tensor.
Note
Compare to Mul(const Tensor &Lt, const Tensor &Rt), this is inplace function and will modify the left Tensor Lt.
See also
Mul(const Tensor &Lt, const Tensor &Rt), Mul(const T &lc, const Tensor &Rt), Mul(const Tensor &Lt, const T &rc), operator*(const Tensor &Lt, const Tensor &Rt)

◆ Inv()

Tensor cytnx::linalg::Inv ( const Tensor Tin,
const double &  clip 
)

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. \]

Parameters
[in]Tina Tensor
[in]clipthe clip value.
Returns
[Tensor]
Note
For complex type Tensors, the square norm is used to determine the clip.

◆ Inv_()

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.

Parameters
[in]Tina Tensor
[in]clipthe clip value.
Returns
[Tensor]
Note
  1. For complex type Tensors, the square norm is used to determine the clip.
  2. on return, all the elements will be modified to it's inverse. if Tin is integer type, it will automatically promote to Type.Double.

◆ InvM() [1/2]

UniTensor cytnx::linalg::InvM ( const cytnx::UniTensor Tin)

◆ InvM() [2/2]

Tensor cytnx::linalg::InvM ( const Tensor Tin)

Matrix inverse.

This function will perform matrix inverse on the input matrix Tin.

Returns
[Tensor] the inversion of the input matrix.
Precondition
Tin should be a rank-2 Tensor.

◆ InvM_() [1/2]

void cytnx::linalg::InvM_ ( Tensor Tin)

inplace matrix inverse.

This function will perform matrix inverse on the input matrix Tin, inplacely.

Note
Compare to InvM, this is inlpace function. The input matrix will be modified to it's inverse.
Precondition
the Tin should be a rank-2 Tensor.

◆ InvM_() [2/2]

void cytnx::linalg::InvM_ ( UniTensor Tin)

◆ iSub()

void cytnx::linalg::iSub ( Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in,out]LtThe left Tensor.
[in]RtThe right Tensor.
Precondition
The shape of Lt and Rt must be the same.
Note
Compare to the function Sub(const Tensor &Lt, const Tensor &Rt), this is a inplace function and it will modify the left Tensor Lt.
See also
Sub(const Tensor &Lt, const Tensor &Rt), Sub(const T &lc, const Tensor &Rt), Sub(const Tensor &Lt, const T &rc), operator-(const Tensor &Lt, const Tensor &Rt)

◆ Kron()

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...)

Parameters
[in]Tlrank-n Tensor #1
[in]Trrank-m Tensor #2
[in]Tl_pad_leftThe padding scheme for Tl if Tl.rank != Tr.rank
[in]Tr_pad_leftThe padding scheme for Tr if Tl.rank != Tr.rank
Returns
[Tensor]
Precondition
two tensor should on same device.

◆ Lanczos() [1/2]

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.

Parameters
[in]Hopthe Linear Operator defined by LinOp class or it's inheritance (see LinOp).
[in]Tinthe 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

Parameters
[in]CvgCritthe convergence criterion of the energy.
[in]maxiterthe maximum interation steps for each k.
[in]kthe number of lowest k eigen values.
[in]is_Vif set to true, the eigen vectors will be returned.
[in]is_rowwhether the return eigen vectors should be in row-major form.
[in]max_krydimthe maximum krylov subspace dimension for each iteration.
[in]verboseprint out iteration info.
Returns
[eigvals (UniTensor), eigvecs (UniTensor)(option)]
Note
To use, define a linear operator with LinOp class either by assign a custom function or create a class that inherit LinOp (see LinOp for further details)

◆ Lanczos() [2/2]

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.

Parameters
[in]Hopthe Linear Operator defined by LinOp class or it's inheritance (see LinOp).
[in]Tinthe 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

Parameters
[in]CvgCritthe convergence criterion of the energy.
[in]maxiterthe maximum interation steps for each k.
[in]kthe number of lowest k eigen values.
[in]is_Vif set to true, the eigen vectors will be returned.
[in]is_rowwhether the return eigen vectors should be in row-major form.
[in]max_krydimthe maximum krylov subspace dimension for each iteration.
[in]verboseprint out iteration info.
Returns
[eigvals (Tensor), eigvecs (Tensor)(option)]
Note
To use, define a linear operator with LinOp class either by assign a custom function or create a class that inherit LinOp (see LinOp for further details)

◆ Lanczos_ER()

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.

Parameters
[in]Hopthe Linear Operator defined by LinOp class or it's inheritance (see LinOp).
[in]kthe number of lowest k eigen values.
[in]is_Vif set to true, the eigen vectors will be returned.
[in]maxiterthe maximum interation steps for each k.
[in]CvgCritthe convergence criterion of the energy.
[in]is_rowwhether the return eigen vectors should be in row-major form.
[in]Tinthe initial vector, this should be rank-1
[in]max_krydimthe maximum krylov subspace dimension for each iteration.
[in]verboseprint out iteration info.
Returns
[eigvals (Tensor), eigvecs (Tensor)(option)]
Note
To use, define a linear operator with LinOp class either by assign a custom function or create a class that inherit LinOp (see LinOp for further details)

◆ Lanczos_Exp()

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.

Parameters
[in]Hopthe Linear Operator defined by LinOp class or it's inheritance (see LinOp). The operation method \(Hv\) need to be defined in it.
[in]vThe input vector (or state). The norm \(|v|\) should be equal to 1.
[in]tauA scalar, it can be complex number.
[in]CvgCrit\(\delta\), the convergence criterion.
[in]Maxiterthe maximum interation steps for each k.
[in]verboseprint out iteration info.
Returns
UniTensor \(w\)
Note
To use, define a linear operator with LinOp class either by assign a custom function or create a class that inherit LinOp (see LinOp for further details)
Warning
User need to guarantee that the input operator \(H\) is Hermitian , and the exponetiate \(e^{-H\tau}\) will converged. Ohterwise, the function will return the wrong results without any warning.

◆ Lanczos_Gnd()

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.

Parameters
[in]Hopthe Linear Operator defined by LinOp class or it's inheritance (see LinOp).
[in]CvgCritthe convergence criterion of the energy.
[in]is_Vif set to true, the eigen vectors will be returned.
[in]Tinthe initial vector, this should be rank-1
[in]verboseprint out iteration info.
[in]maxiterthe maximum interation steps for each k.
Returns
[eigvals (Tensor), eigvecs (Tensor)(option)]
Note
To use, define a linear operator with LinOp class either by assign a custom function or create a class that inherit LinOp (see LinOp for further details)

◆ Lanczos_Gnd_Ut()

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.

Parameters
[in]Hopthe Linear Operator defined by LinOp class or it's inheritance (see LinOp).
[in]CvgCritthe convergence criterion of the energy.
[in]is_Vif set to true, the eigen vectors will be returned.
[in]Tinthe initial vector, this should be a UniTensor.
[in]verboseprint out iteration info.
[in]maxiterthe maximum interation steps for each k.
Returns
[eigvals (UniTensor::Dense), eigvecs (UniTensor)(option)]
Note
To use, define a linear operator with LinOp class either by assign a custom function or create a class that inherit LinOp (see LinOp for further details)

◆ Lstsq()

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 \).

Parameters
[in]A“Coefficient” matrix, must be two-dimensional.
[in]bOrdinate or “dependent variable” values, must be two-dimensional, the least-squares solution is calculated for each of the K columns of b.
[in]rcondCut-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.
Returns
[std::vector<Tensors>]
  1. the first tensor is least-squares solutions in the K columns.
  2. the second tensor is the sums of squared residuals: Squared Euclidean 2-norm for each column in b - a @ x. If the rank of a is < N or M <= N, this is a zero Tensor.
  3. the third tensor is the rank of matrix A.
  4. the forth tensor is singular values of A.
Author
Ke

◆ Matmul()

Tensor cytnx::linalg::Matmul ( const Tensor TL,
const Tensor TR 
)

perform matrix multiplication on two tensors.

This function will perform matrix multiplication on two matrices (2-rank Tensor) TL and TR. The result will be:

\[ T = T_L T_R \]

Parameters
[in]TLa left Tensor
[in]TRa right Tensor
Precondition
the TL and TR should be rank-2 Tensor.

◆ Matmul_dg()

Tensor cytnx::linalg::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.

Note
one of TL and TR should be rank-1 Tensor and the other should be rank-2 Tensor.

◆ Max()

Tensor cytnx::linalg::Max ( const Tensor Tn)

get the maximum element.

Parameters
[in]Tna cytnx::Tensor
Note
For complex TN, only real part is compared.

◆ Min()

Tensor cytnx::linalg::Min ( const Tensor Tn)

get the minimum element.

Parameters
[in]Tna cytnx::Tensor
Note
For complex TN, only real part is compared.

◆ Mod() [1/6]

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 \).

Parameters
[in]LtThe left UniTensor.
[in]RtThe right UniTensor.
Returns
The result UniTensor.
Precondition
  1. Lt and Rt must have the same shape.
  2. The input UniTensor Lt and Rt need to be integer type.
See also
Tensor::Mod(const T &rhs), operator%(const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt)

◆ Mod() [2/6]

template<class T >
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 \).

Parameters
[in]LtThe 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.

Returns
The result UniTensor.
Precondition
The input Lt and rc need to be integer type.
Note
The inpute template type rc will be casted to the same type as the UniTensor Lt.
See also
operator%(const cytnx::UniTensor &Lt, const T &rc), Mod(const cytnx::UniTensor &Lt, const T &rc)

◆ Mod() [3/6]

template<class T >
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 \).

Parameters
lcThe 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.

RtThe right UniTensor.
Returns
The result UniTensor.
Precondition
The input lc and Rt need to be integer type.
Note
The inpute template type lc will be casted to the same type as the UniTensor Rt.
See also
operator%(const cytnx::UniTensor &Lt, const T &rc), Mod(const cytnx::UniTensor &Lt, const T &rc)

◆ Mod() [4/6]

template<class T >
Tensor cytnx::linalg::Mod ( const T &  lc,
const Tensor 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.

Parameters
[in]lcThe left template type.
[in]RtThe right Tensor.
Returns
The result Tensor.
Precondition
the right template type rc should be integer type.
See also
Mod(const Tensor &Lt, const Tensor &Rt), Mod(const Tensor &Lt, const T &rc)

◆ Mod() [5/6]

template<class T >
Tensor cytnx::linalg::Mod ( const Tensor Lt,
const T &  rc 
)

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.

Parameters
[in]LtThe left Tensor.
[in]rcThe right template type.
Returns
The result Tensor.
Precondition
the right template type rc should be integer type.
See also
Mod(const Tensor &Lt, const Tensor &Rt), Mod(const T &lc, const Tensor &Rt)

◆ Mod() [6/6]

Tensor cytnx::linalg::Mod ( const Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in]LtThe left Tensor.
[in]RtThe right Tensor.
Returns
The result Tensor.
Precondition
The input tensors Lt and Rt should have the same shape and need to be integer type.
See also
Mod(const T &lc, const Tensor &Rt), Mod(const Tensor &Lt, const T &rc)

◆ Mul() [1/6]

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.

Parameters
[in]LtThe left UniTensor.
[in]RtThe right UniTensor.
Returns
The result UniTensor.
Precondition
Lt and Rt must have the same shape.
See also
UniTensor::Mul(const cytnx::UniTensor &Rt) const, operator*(const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt)

◆ Mul() [2/6]

template<class T >
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 \).

Parameters
[in]LtThe 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.

Returns
The result UniTensor.
Precondition
The supported template type shown above.
Note
The inpute template type rc will be casted to the same type as the UniTensor Lt.
See also
operator*(const cytnx::UniTensor &Lt, const T &rc), Mul(const cytnx::UniTensor &Lt, const T &rc)

◆ Mul() [3/6]

template<class T >
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 \).

Parameters
[in]LtThe 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.

Returns
The result UniTensor.
Precondition
The supported template type shown above.
Note
The inpute template type lc will be casted to the same type as the UniTensor Rt.
See also
operator*(const T &lc, const cytnx::UniTensor &Rt), Mul(const T &lc, const cytnx::UniTensor &Rt)

◆ Mul() [4/6]

template<class T >
Tensor cytnx::linalg::Mul ( const T &  lc,
const Tensor 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.

Parameters
[in]LtThe left Tensor.
[in]rcThe right template type.
Returns
The result Tensor.
See also
Mul(const Tensor &Lt, const Tensor &Rt), Mul(const T &lc, const Tensor &Rt), iMul(Tensor &Lt, const Tensor &Rt), operator*(const Tensor &Lt, const Tensor &Rt)

◆ Mul() [5/6]

template<class T >
Tensor cytnx::linalg::Mul ( const Tensor Lt,
const T &  rc 
)

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.

Parameters
[in]LtThe left Tensor.
[in]rcThe right template type.
Returns
The result Tensor.
See also
Mul(const Tensor &Lt, const Tensor &Rt), Mul(const T &lc, const Tensor &Rt), iMul(Tensor &Lt, const Tensor &Rt), operator*(const Tensor &Lt, const Tensor &Rt)

◆ Mul() [6/6]

Tensor cytnx::linalg::Mul ( const Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in]LtThe left Tensor.
[in]RtThe right Tensor.
Returns
The result Tensor.
See also
Mul(const T &lc, const Tensor &Rt), Mul(const Tensor &Lt, const T &rc), iMul(Tensor &Lt, const Tensor &Rt), operator*(const Tensor &Lt, const Tensor &Rt)

◆ Norm() [1/2]

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.

Parameters
[in]uTlinput UniTensor
Returns
Tensor

◆ Norm() [2/2]

Tensor cytnx::linalg::Norm ( const Tensor Tl)

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.

Parameters
[in]Tlinput Tensor
Returns
Tensor

◆ Outer()

Tensor cytnx::linalg::Outer ( const Tensor Tl,
const Tensor Tr 
)

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)

Parameters
[in]Tlrank-1 Tensor #1
[in]Trrank-1 Tensor #2
Returns
[Tensor]
Precondition
two tensor should on same device.

◆ Pow() [1/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.

Parameters
[in]Tinthe input UniTensor
[in]pthe power
Precondition
If Tin is a real UniTensor and containt negative elements, then p must be an integer.
Returns
UniTensor with the same shape as Tin, but with the elements are the power of Tin.
Note
Compare to the Pow_(UniTensor &Tin, const double &p) function, this function will not modify the input UniTensor and return a new UniTensor.
See also
Pow_(UniTensor &Tin, const double &p)

◆ Pow() [2/2]

Tensor cytnx::linalg::Pow ( const Tensor Tin,
const double &  p 
)

take power p on all the elements in Tensor.

This function will perform power p on all the elements in Tensor Tin. That is, the output will be:

\[ T_{o}[i] = T_{i}[i]^{p} \]

Parameters
[in]p,thepower
Returns
[Tensor]

◆ Pow_() [1/2]

void cytnx::linalg::Pow_ ( Tensor Tin,
const double &  p 
)

inplace perform power on all the elements in Tensor.

this is just a inplace version of Pow. The input Tensor Tin will be modified.

Parameters
[in]Tin,theinput Tensor.
[in]p,thepower.

◆ Pow_() [2/2]

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.

Parameters
[in,out]Tinthe input UniTensor
[in]pthe power
Precondition
If Tin is a real UniTensor and containt negative elements, then p must be an integer.
Note
Compare to the Pow function, this is an inplacely function, which will modify the input UniTensor.
See also
Pow(const cytnx::UniTensor &Tin, const double &p)

◆ Qdr() [1/2]

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).

See also
Qdr(const Tensor &Tin, const bool &is_tau)

◆ Qdr() [2/2]

std::vector< Tensor > cytnx::linalg::Qdr ( const Tensor Tin,
const bool &  is_tau = false 
)

Perform QDR decomposition on a rank-2 Tensor.

Parameters
[in]Tina cytnx::Tensor, it should be a rank-2 tensor (matrix)
[in]is_tauif return the tau that contains the Householder reflectors that generate q along with r. The tau array contains scaling factors for the reflectors
Returns

[std::vector<Tensors>]

  1. the first tensor is the orthomormal matrix \( Q \), a 2-d tensor (matrix)
  2. the second tensor is the diagonal matrix \( D \), a 1-d tensor (diagonal matrix)
  3. the third tensor is the right-upper triangular matrix \( R \), a 2-d tensor (matrix)
  4. the forth tensor is the Householder reflectors \( H \), a 1-d tensor (vector). This tensor will only return when is_tau = true.
Precondition
The input tensor should be a rank-2 tensor (matrix).
See also
Qr(const Tensor &Tin, const bool &is_tau)

◆ Qr() [1/2]

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).

See also
Qr(const Tensor &Tin, const bool &is_tau)

◆ Qr() [2/2]

std::vector< Tensor > cytnx::linalg::Qr ( const Tensor Tin,
const bool &  is_tau = false 
)

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.

Parameters
[in]Tina Tensor, it should be a rank-2 tensor (a matrix)
[in]is_tauif return the tau that contains the Householder reflectors that generate q along with r. The tau array contains scaling factors for the reflectors
Returns

[std::vector<Tensors>]

  1. the first tensor is the orthomormal matrix \( Q \), a 2-d tensor (matrix)
  2. the second tensor is the right-upper triangular matrix \( R \), a 2-d tensor (matrix)
  3. the third tensor is the Householder reflectors \( H \), a 1-d tensor (vector). This tensor will only return when is_tau = true.
Precondition
The input tensor should be a rank-2 tensor (matrix).
See also
Qdr(const Tensor &Tin, const bool &is_tau)

◆ Sub() [1/6]

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.

Parameters
[in]LtThe left UniTensor.
[in]RtThe right UniTensor.
Returns
The result UniTensor.
Precondition
Lt and Rt must have the same shape.
See also
UniTensor::Sub(const cytnx::UniTensor &Rt) const, operator-(const cytnx::UniTensor &Lt, const cytnx::UniTensor &Rt)

◆ Sub() [2/6]

template<class T >
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 \).

Parameters
[in]LtThe 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.

Returns
The result UniTensor.
Precondition
The supported template type shown above.
Note
The inpute template type rc will be casted to the same type as the UniTensor Lt.
See also
operator-(const cytnx::UniTensor &Lt, const T &rc), Sub(const cytnx::UniTensor &Lt, const T &rc)

◆ Sub() [3/6]

template<class T >
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 \).

Parameters
[in]LtThe 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.

Returns
The result UniTensor.
Precondition
The supported template type shown above.
Note
The inpute template type lc will be casted to the same type as the UniTensor Rt.
See also
operator-(const T &lc, const cytnx::UniTensor &Rt), Sub(const T &lc, const cytnx::UniTensor &Rt)

◆ Sub() [4/6]

template<class T >
Tensor cytnx::linalg::Sub ( const T &  lc,
const Tensor 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.

Parameters
[in]lcThe left template type.
[in]RtThe right Tensor.
Returns
The result Tensor.
See also
Sub(const Tensor &Lt, const Tensor &Rt), Sub(const Tensor &Lt, const T &rc), iSub(Tensor &Lt, const Tensor &Rt), operator-(const Tensor &Lt, const Tensor &Rt)

◆ Sub() [5/6]

template<class T >
Tensor cytnx::linalg::Sub ( const Tensor Lt,
const T &  rc 
)

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.

Parameters
[in]LtThe left Tensor.
[in]rcThe right template type.
Returns
The result Tensor.
See also
Sub(const Tensor &Lt, const Tensor &Rt), Sub(const T &lc, const Tensor &Rt), iSub(Tensor &Lt, const Tensor &Rt), operator-(const Tensor &Lt, const Tensor &Rt)

◆ Sub() [6/6]

Tensor cytnx::linalg::Sub ( const Tensor Lt,
const Tensor Rt 
)

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.

Parameters
[in]LtThe left Tensor.
[in]RtThe right Tensor.
Returns
The result Tensor.
See also
Sub(const T &lc, const Tensor &Rt), Sub(const Tensor &Lt, const T &rc), iSub(Tensor &Lt, const Tensor &Rt), operator-(const Tensor &Lt, const Tensor &Rt)

◆ Sum()

Tensor cytnx::linalg::Sum ( const Tensor Tn)

get the sum of all the elements.

Parameters
[in]Tna cytnx::Tensor

◆ Svd() [1/2]

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).

See also
Svd(const Tensor &Tin, const bool &is_UvT)

◆ Svd() [2/2]

std::vector< Tensor > cytnx::linalg::Svd ( const Tensor Tin,
const bool &  is_UvT = 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.

Parameters
[in]Tina Tensor, it should be a rank-2 tensor (matrix)
[in]is_UvTwhether need to return a left unitary matrix.
Returns

[std::vector<Tensors>]

  1. The first tensor is a 1-d tensor contanin the singular values
  2. If is_UvT is true, then the tensors \( U,V^\dagger \) will be pushed back to the vector.
Precondition
The input tensor should be a rank-2 tensor (matrix).
See also
Svd_truncate(const Tensor &Tin, const cytnx_uint64 &keepdim, const double &err, const bool &is_UvT, const unsigned int &return_err, const cytnx_uint64& mindim)

◆ Svd_truncate() [1/3]

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.

See also
Svd_truncate(const cytnx::UniTensor &Tin, const cytnx_uint64 &keepdim, const std::vector<cytnx_uint64> min_blockdim, const double &err, const bool &is_UvT, const unsigned int &return_err, const cytnx_uint64 &mindim), Svd_truncate(const Tensor &Tin, const cytnx_uint64 &keepdim, const double &err, const bool &is_UvT, const unsigned int &return_err, const cytnx_uint64& mindim)

◆ Svd_truncate() [2/3]

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)

Parameters
[in]Tina BlockUniTensor, with the correct rowrank set to interpret it as a matrix
[in]keepdimthe number (at most) of singular values to keep.
[in]min_blockdima 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]errthe cutoff error (the singular values smaller than err will be truncated.)
[in]is_UvTif true, the left- and right- unitary UniTensors (isometries) are returned.
[in]return_errwhether 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.
Returns

[std::vector<Tensors>]

  1. The first UniTensor is a diagonal UniTensor containing the singular values
  2. If is_UvT is true, then the UniTensor \( U \) and \( V^\dagger \) will be pushed back to the vector.
  3. If return_err is true, then the error will be pushed back to the vector.
Precondition
This function assumes a BlockUniTensor as input for Tin.
See also
Svd_truncate(const Tensor &Tin, const cytnx_uint64 &keepdim, const double &err, const bool &is_UvT, const unsigned int &return_err, const cytnx_uint64& mindim)
Note
The truncated bond dimension can be larger than 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.

◆ Svd_truncate() [3/3]

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)

Parameters
[in]Tina Tensor, it should be a rank-2 tensor (matrix)
[in]keepdimthe number (at most) of singular values to keep.
[in]errthe cutoff error (the singular values smaller than err will be truncated.)
[in]is_UvTif true, the left- and right- unitary matrices (isometries) are returned.
[in]return_errwhether 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.
Returns

[std::vector<Tensors>]

  1. The first tensor is a 1-d tensor containing the singular values
  2. If is_UvT is true, then the tensors \( U \) and \( V^\dagger \) will be pushed back to the vector.
  3. If return_err is true, then the error will be pushed back to the vector.
Precondition
The input tensor should be a rank-2 tensor (matrix).
See also
Svd(const Tensor &Tin, const bool &is_UvT)
Note
The truncated bond dimension can be larger than 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.

◆ Tensordot()

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.

Parameters
[in]TlTensor #1
[in]TrTensor #2
[in]idxlthe indices of rank of Tensor #1 that is going to sum with Tensor #2
[in]idxrthe indices of rank of Tensor #2 that is going to sum with Tensor #1
[in]cacheLcache Tensor #1 (See user-guide for details)
[in]cacheRcache Tensor #2 (See user-guide for details)
Returns
[Tensor]
Note
  1. the elements in idxl and idxr have one to one correspondence.
  2. two tensors should on same device.

◆ Tensordot_dg()

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.

Parameters
[in]TlTensor #1
[in]TrTensor #2
[in]idxlthe indices of rank of Tensor #1 that is going to sum with Tensor #2
[in]idxrthe indices of rank of Tensor #2 that is going to sum with Tensor #1
[in]diag_Lif Tl(true)/Tr(false) is a diagnal matrix, represented by a rank-1 tensor.
Returns
[Tensor]
Note
  1. the elements in idxl and idxr have one to one correspondence.
  2. two tensors should on same device.
  3. if diag_L=true, Tl should be a rank-1 tensor as the diagonal elements of a diagonal matrix. if false, Tr should be a rank-1 tensor

◆ Trace() [1/3]

cytnx::UniTensor cytnx::linalg::Trace ( const cytnx::UniTensor Tin,
const cytnx_int64 &  a = 0,
const cytnx_int64 &  b = 1 
)

◆ Trace() [2/3]

cytnx::UniTensor cytnx::linalg::Trace ( const cytnx::UniTensor Tin,
const std::string &  a,
const std::string &  b 
)

Perform trace over two legs of a UniTensor.

This function performs trace over two legs of a UniTensor Tin. The two legs are specified by a and b. For more details, please refer to the documentation of the function Trace(const Tensor &Tin, const cytnx_int64 &a, const cytnx_int64 &b).

See also
Trace(const Tensor &Tin, const cytnx_uint64 &a, const cytnx_uint64 &b)

◆ Trace() [3/3]

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} \]

Parameters
[in]Tna Tensor
[in]axisAthe first index to perform trace.
[in]axisBthe second index to perform trace.
Precondition
the Tn should be at-least rank-2 Tensor.

◆ Tridiag()

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.

Parameters
[in]DiagTensor #1
[in]Sub_diagTensor #2
[in]is_Vif calculate the eigen value.
[in]kReturn k lowest eigen vector if is_V=True
[in]throw_excpWhether to throw exception when error occurs in Tridiag internal function
Returns
[vector<Tensor>] if is_V = True, the first tensor is the eigen value, and second tensor is eigenvector of shape [k,L].
Precondition
two Tensors must be Rank-1, with length of Diag = L and Sub_diag length = L-1.
Note
performance tune: This function have better performance when two vectors with same types, and are one of following type: cytnx_double, cytnx_float. In general all real type can be use as input, which will be promote to floating point type for calculation.

◆ Vectordot()

Tensor cytnx::linalg::Vectordot ( const Tensor Tl,
const Tensor Tr,
const bool &  is_conj = false 
)

perform inner product of vectors

Parameters
[in]TlTensor #1
[in]TrTensor #2
[in]ifthe Tl should be conjugated (only work for complex. For real Tensor, no function), default: false
Returns
[Tensor] Rank-0
Precondition
two Tensors must be Rank-1, with same length.
Note
performance tune: This function have better performance when two vectors with same types, and are one of following type: cytnx_double, cytnx_float, cytnx_complex64 or cytnx_complex128.