GraceQ/MPS2
A high-performance matrix product state algorithms library based on GraceQ/tensor
|
A generic MPO generator. More...
#include <mpogen.h>
Public Member Functions | |
MPOGenerator (const SiteVec< TenElemT, QNT > &, const QNT &) | |
Create a MPO generator. More... | |
void | AddTerm (const TenElemT, const GQTensorVec &, const std::vector< size_t > &) |
The most generic API for adding a many-body term to the MPO generator. More... | |
void | AddTerm (const TenElemT, const GQTensorVec &, const std::vector< size_t > &, const GQTensorVec &, const std::vector< std::vector< size_t > > &inst_ops_idxs_set=kNullUintVecVec) |
void | AddTerm (const TenElemT, const GQTensorT &, const size_t, const GQTensorT &op2=GQTensorT(), const size_t op2_idx=0, const GQTensorT &inst_op=GQTensorT(), const std::vector< size_t > &inst_op_idxs=kNullUintVec) |
Add one-body or two-body interaction term. More... | |
FSM | GetFSM (void) |
MPO< GQTensorT > | Gen (void) |
A generic MPO generator.
A matrix-product operator (MPO) generator which can generate an efficient MPO for a quantum many-body system with any type of n-body interaction term.
TenElemType | Element type of the MPO tensors, can be GQTEN_Double or GQTEN_Complex. |
gqmps2::MPOGenerator< TenElemT, QNT >::MPOGenerator | ( | const SiteVec< TenElemT, QNT > & | site_vec, |
const QNT & | zero_div | ||
) |
Create a MPO generator.
Create a MPO generator using the sites of the system which is described by a SiteVec.
site_vec | The local Hilbert spaces of each sites of the system. |
zero_div | The zero value of the given quantum number type which is used to set the divergence of the MPO. |
void gqmps2::MPOGenerator< TenElemT, QNT >::AddTerm | ( | const TenElemT | coef, |
const GQTensorVec & | local_ops, | ||
const std::vector< size_t > & | local_ops_idxs | ||
) |
The most generic API for adding a many-body term to the MPO generator.
Notice that the indexes of the operators have to be ascending sorted.
coef | The coefficient of the term. |
local_ops | All the local (on-site) operators in the term. |
local_ops_idxs | The site indexes of these local operators. |
void gqmps2::MPOGenerator< TenElemT, QNT >::AddTerm | ( | const TenElemT | coef, |
const GQTensorT & | op1, | ||
const size_t | op1_idx, | ||
const GQTensorT & | op2 = GQTensorT() , |
||
const size_t | op2_idx = 0 , |
||
const GQTensorT & | inst_op = GQTensorT() , |
||
const std::vector< size_t > & | inst_op_idxs = kNullUintVec |
||
) |
Add one-body or two-body interaction term.
coef | The coefficient of the term. |
op1 | The first physical operator for the term. |
op1_idx | The site index of the first physical operator. |
op2 | The second physical operator for the term. |
op2_idx | The site index of the second physical operator. |
inst_op | The insertion operator for the two-body interaction term. |
inst_op_idxs | The explicit site indexes of the insertion operator. |