VMPS++
|
#include "tensors/DmrgContractions.h"
#include "VUMPS/VumpsTransferMatrix.h"
#include "ArnoldiSolver.h"
#include "Mpo.h"
Go to the source code of this file.
Functions | |
template<typename Symmetry , typename Scalar > | |
Eigenstate< Biped< Symmetry, Matrix< complex< Scalar >, Dynamic, Dynamic > > > | calc_LReigen (VMPS::DIRECTION::OPTION DIR, const vector< vector< Biped< Symmetry, Matrix< Scalar, Dynamic, Dynamic > > > > &Abra, const vector< vector< Biped< Symmetry, Matrix< Scalar, Dynamic, Dynamic > > > > &Aket, const Qbasis< Symmetry > &basisBra, const Qbasis< Symmetry > &basisKet, const vector< vector< qarray< Symmetry::Nq > > > &qloc, size_t dimK=100ul, double tol_input=1e-12, Biped< Symmetry, Matrix< Scalar, Dynamic, Dynamic > > *LReigenTop=NULL) |
template<typename Symmetry , typename MatrixType , typename MpoScalar > | |
Biped< Symmetry, MatrixType > | make_hL (const boost::multi_array< MpoScalar, 4 > &H2site, const vector< Biped< Symmetry, MatrixType > > &AL, const vector< qarray< Symmetry::Nq > > &qloc) |
template<typename Symmetry , typename MatrixType , typename MpoScalar > | |
Biped< Symmetry, MatrixType > | make_hR (const boost::multi_array< MpoScalar, 4 > &H2site, const vector< Biped< Symmetry, MatrixType > > &AR, const vector< qarray< Symmetry::Nq > > &qloc) |
template<typename Symmetry , typename MatrixType , typename MpoScalar > | |
Tripod< Symmetry, MatrixType > | make_YL (size_t b, const Tripod< Symmetry, MatrixType > &Lold, const vector< vector< Biped< Symmetry, MatrixType > > > &Abra, const vector< vector< vector< vector< Biped< Symmetry, SparseMatrix< MpoScalar > > > > > > &W, const vector< vector< Biped< Symmetry, MatrixType > > > &Aket, const vector< vector< qarray< Symmetry::Nq > > > &qloc, const vector< vector< qarray< Symmetry::Nq > > > &qOp, const std::unordered_map< pair< qarray< Symmetry::Nq >, size_t >, size_t > &basis_order_map) |
template<typename Symmetry , typename MatrixType , typename MpoScalar > | |
Tripod< Symmetry, MatrixType > | make_YR (size_t a, const Tripod< Symmetry, MatrixType > &Rold, const vector< vector< Biped< Symmetry, MatrixType > > > &Abra, const vector< vector< vector< vector< Biped< Symmetry, SparseMatrix< MpoScalar > > > > > > &W, const vector< vector< Biped< Symmetry, MatrixType > > > &Aket, const vector< vector< qarray< Symmetry::Nq > > > &qloc, const vector< vector< qarray< Symmetry::Nq > > > &qOp, const std::unordered_map< pair< qarray< Symmetry::Nq >, size_t >, size_t > &basis_order_map) |
template<typename Symmetry , typename MpoScalar > | |
void | contract_WW (const vector< vector< vector< SparseMatrix< MpoScalar > > > > &W12, const vector< qarray< Symmetry::Nq > > &qloc12, const vector< qarray< Symmetry::Nq > > &qOp12, const vector< vector< vector< SparseMatrix< MpoScalar > > > > &W34, const vector< qarray< Symmetry::Nq > > &qloc34, const vector< qarray< Symmetry::Nq > > &qOp34, vector< vector< vector< SparseMatrix< MpoScalar > > > > &W, vector< qarray< Symmetry::Nq > > &qloc, vector< qarray< Symmetry::Nq > > &qOp) |
Eigenstate< Biped< Symmetry, Matrix< complex< Scalar >, Dynamic, Dynamic > > > calc_LReigen | ( | VMPS::DIRECTION::OPTION | DIR, |
const vector< vector< Biped< Symmetry, Matrix< Scalar, Dynamic, Dynamic > > > > & | Abra, | ||
const vector< vector< Biped< Symmetry, Matrix< Scalar, Dynamic, Dynamic > > > > & | Aket, | ||
const Qbasis< Symmetry > & | basisBra, | ||
const Qbasis< Symmetry > & | basisKet, | ||
const vector< vector< qarray< Symmetry::Nq > > > & | qloc, | ||
size_t | dimK = 100ul , |
||
double | tol_input = 1e-12 , |
||
Biped< Symmetry, Matrix< Scalar, Dynamic, Dynamic > > * | LReigenTop = NULL |
||
) |
Definition at line 17 of file VumpsContractions.h.
void contract_WW | ( | const vector< vector< vector< SparseMatrix< MpoScalar > > > > & | W12, |
const vector< qarray< Symmetry::Nq > > & | qloc12, | ||
const vector< qarray< Symmetry::Nq > > & | qOp12, | ||
const vector< vector< vector< SparseMatrix< MpoScalar > > > > & | W34, | ||
const vector< qarray< Symmetry::Nq > > & | qloc34, | ||
const vector< qarray< Symmetry::Nq > > & | qOp34, | ||
vector< vector< vector< SparseMatrix< MpoScalar > > > > & | W, | ||
vector< qarray< Symmetry::Nq > > & | qloc, | ||
vector< qarray< Symmetry::Nq > > & | qOp | ||
) |
Definition at line 297 of file VumpsContractions.h.
Biped< Symmetry, MatrixType > make_hL | ( | const boost::multi_array< MpoScalar, 4 > & | H2site, |
const vector< Biped< Symmetry, MatrixType > > & | AL, | ||
const vector< qarray< Symmetry::Nq > > & | qloc | ||
) |
**Calculates the tensor (eq. (12)) from the explicit 4-legged 2-site Hamiltonian and .*/
Definition at line 72 of file VumpsContractions.h.
Biped< Symmetry, MatrixType > make_hR | ( | const boost::multi_array< MpoScalar, 4 > & | H2site, |
const vector< Biped< Symmetry, MatrixType > > & | AR, | ||
const vector< qarray< Symmetry::Nq > > & | qloc | ||
) |
Calculates the tensor (eq. (12)) from the explicit 4-legged 2-site Hamiltonian and .
Definition at line 150 of file VumpsContractions.h.
Tripod< Symmetry, MatrixType > make_YL | ( | size_t | b, |
const Tripod< Symmetry, MatrixType > & | Lold, | ||
const vector< vector< Biped< Symmetry, MatrixType > > > & | Abra, | ||
const vector< vector< vector< vector< Biped< Symmetry, SparseMatrix< MpoScalar > > > > > > & | W, | ||
const vector< vector< Biped< Symmetry, MatrixType > > > & | Aket, | ||
const vector< vector< qarray< Symmetry::Nq > > > & | qloc, | ||
const vector< vector< qarray< Symmetry::Nq > > > & | qOp, | ||
const std::unordered_map< pair< qarray< Symmetry::Nq >, size_t >, size_t > & | basis_order_map | ||
) |
Calculates the tensor (eq. (C17)) from the MPO tensor W
, the left transfer matrix L
and .
Definition at line 228 of file VumpsContractions.h.
Tripod< Symmetry, MatrixType > make_YR | ( | size_t | a, |
const Tripod< Symmetry, MatrixType > & | Rold, | ||
const vector< vector< Biped< Symmetry, MatrixType > > > & | Abra, | ||
const vector< vector< vector< vector< Biped< Symmetry, SparseMatrix< MpoScalar > > > > > > & | W, | ||
const vector< vector< Biped< Symmetry, MatrixType > > > & | Aket, | ||
const vector< vector< qarray< Symmetry::Nq > > > & | qloc, | ||
const vector< vector< qarray< Symmetry::Nq > > > & | qOp, | ||
const std::unordered_map< pair< qarray< Symmetry::Nq >, size_t >, size_t > & | basis_order_map | ||
) |
Calculates the tensor (eq. (C18)) from the MPO tensor W
, the left transfer matrix R
and .
Definition at line 263 of file VumpsContractions.h.