VMPS++
|
Matrix Product State with conserved quantum numbers (Abelian and non abelian symmetries).
Symmetry | : A class providing all relevant functions and infos that are determined by the Symmetry. Click here for more information. |
Scalar | : double or complex<double> |
#include <Mps.h>
Public Member Functions | |
Mps () | |
Mps (size_t L_input, vector< vector< qarray< Nq > > > qloc_input, qarray< Nq > Qtot_input, size_t N_phys_input, int Nqmax_input, bool TRIVIAL_BOUNDARIES=true) | |
template<typename Hamiltonian > | |
Mps (const Hamiltonian &H, size_t Mmax, qarray< Nq > Qtot_input, int Nqmax_input) | |
Mps (size_t L_input, const vector< vector< Biped< Symmetry, MatrixType > > > &As, const vector< vector< qarray< Nq > > > &qloc_input, qarray< Nq > Qtot_input, size_t N_phys_input) | |
void | set_open_bc () |
Tripod< Symmetry, Matrix< Scalar, Dynamic, Dynamic > > | get_boundaryTensor (DMRG::DIRECTION::OPTION DIR, size_t usePower=1ul) const |
void | elongate_hetero (size_t Nleft=0, size_t Nright=0) |
void | shift_hetero (int Nshift=0) |
void | transform_base (qarray< Symmetry::Nq > Qtot, int L, bool PRINT=false) |
void | calc_Qlimits () |
void | set_Qlimits_to_inf () |
void | update_inbase (size_t loc) |
void | update_outbase (size_t loc) |
void | update_inbase () |
void | update_outbase () |
void | resize_arrays () |
void | outerResizeNoSymm () |
template<typename OtherScalar > | |
void | add_site (size_t loc, OtherScalar alpha, const Mps< Symmetry, Scalar > &Vin) |
void | enrich_left (size_t loc, PivotMatrix1< Symmetry, Scalar, Scalar > *H) |
void | enrich_right (size_t loc, PivotMatrix1< Symmetry, Scalar, Scalar > *H) |
void | setRandom () |
void | setRandom (size_t loc) |
void | setZero () |
void | canonize (DMRG::DIRECTION::OPTION DIR=DMRG::DIRECTION::LEFT) |
void | outerResize (size_t L_input, vector< vector< qarray< Nq > > > qloc_input, qarray< Nq > Qtot_input, int Nqmax_input=500) |
void | outerResizeAll (size_t L_input, vector< vector< qarray< Nq > > > qloc_input, qarray< Nq > Qtot_input) |
template<typename Hamiltonian > | |
void | outerResize (const Hamiltonian &H, qarray< Nq > Qtot_input, int Nqmax_input=500) |
template<typename OtherMatrixType > | |
void | outerResize (const Mps< Symmetry, OtherMatrixType > &V) |
void | innerResize (size_t Mmax) |
template<typename Hamiltonian > | |
void | setProductState (const Hamiltonian &H, const vector< qarray< Nq > > &config) |
void | mend () |
void | set_A_from_C (size_t loc, const vector< Tripod< Symmetry, MatrixType > > &C, DMRG::BROOM::OPTION TOOL=DMRG::BROOM::SVD) |
void | set_Qmultitarget (const vector< qarray< Nq > > &Qmulti_input) |
string | test_ortho (double tol=1e-8) const |
string | info () const |
string | Asizes () const |
string | validate (string name="Mps") const |
void | graph (string filename) const |
size_t | calc_Mmax () const |
size_t | get_Min (size_t loc) const |
size_t | get_Mout (size_t loc) const |
size_t | calc_fullMmax () const |
size_t | calc_Dmax () const |
size_t | calc_Nqmax () const |
size_t | Nqout (size_t l) |
double | calc_Nqavg () const |
double | memory (MEMUNIT memunit=GB) const |
template<typename OtherScalar > | |
void | addScale (OtherScalar alpha, const Mps< Symmetry, Scalar > &Vin, bool SVD_COMPRESS=false) |
Mps< Symmetry, Scalar > & | operator+= (const Mps< Symmetry, Scalar > &Vin) |
Mps< Symmetry, Scalar > & | operator-= (const Mps< Symmetry, Scalar > &Vin) |
template<typename OtherScalar > | |
Mps< Symmetry, Scalar > & | operator*= (const OtherScalar &alpha) |
template<typename OtherScalar > | |
Mps< Symmetry, Scalar > & | operator/= (const OtherScalar &alpha) |
void | applyGate (const TwoSiteGate< Symmetry, Scalar > &gate, size_t l, DMRG::DIRECTION::OPTION DIR) |
template<typename OtherScalar > | |
Mps< Symmetry, OtherScalar > | cast () const |
Scalar | dot (const Mps< Symmetry, Scalar > &Vket) const |
double | squaredNorm () const |
template<typename MpoScalar > | |
Scalar | locAvg (const Mpo< Symmetry, MpoScalar > &O, size_t distance=0) const |
void | swap (Mps< Symmetry, Scalar > &V) |
void | get_controlParams (const Mps< Symmetry, Scalar > &V) |
void | collapse () |
void | rightSweepStep (size_t loc, DMRG::BROOM::OPTION BROOM, PivotMatrix1< Symmetry, Scalar, Scalar > *H=NULL, bool DISCARD_V=false) |
void | leftSweepStep (size_t loc, DMRG::BROOM::OPTION BROOM, PivotMatrix1< Symmetry, Scalar, Scalar > *H=NULL, bool DISCARD_U=false) |
void | calc_N (DMRG::DIRECTION::OPTION DIR, size_t loc, vector< Biped< Symmetry, MatrixType > > &N) |
void | sweepStep2 (DMRG::DIRECTION::OPTION DIR, size_t loc, const vector< Biped< Symmetry, MatrixType > > &Apair, bool SEPARATE_SV=false) |
void | leftSplitStep (size_t loc, Biped< Symmetry, MatrixType > &C) |
void | rightSplitStep (size_t loc, Biped< Symmetry, MatrixType > &C) |
void | absorb (size_t loc, DMRG::DIRECTION::OPTION DIR, const Biped< Symmetry, MatrixType > &C) |
qarray< Nq > | Qtarget () const |
vector< qarray< Nq > > | Qmultitarget () const |
vector< qarray< Nq > > | locBasis (size_t loc) const |
vector< vector< qarray< Nq > > > | locBasis () const |
Qbasis< Symmetry > | inBasis (size_t loc) const |
vector< Qbasis< Symmetry > > | inBasis () const |
Qbasis< Symmetry > | outBasis (size_t loc) const |
vector< Qbasis< Symmetry > > | outBasis () const |
const vector< Biped< Symmetry, MatrixType > > & | A_at (size_t loc) const |
vector< Biped< Symmetry, MatrixType > > & | A_at (size_t loc) |
int | get_pivot () const |
ArrayXd | get_truncWeight () const |
ArrayXd | entropy () const |
vector< map< qarray< Nq >, ArrayXd > > | entanglementSpectrum () const |
std::pair< vector< qarray< Symmetry::Nq > >, ArrayXd > | entanglementSpectrumLoc (size_t loc) const |
Public Member Functions inherited from DmrgJanitor< PivotMatrixType > | |
DmrgJanitor () | |
DmrgJanitor (size_t L_input) | |
void | set_defaultCutoffs () |
void | set_pivot (int pivot_input) |
void | sweep (size_t loc, DMRG::BROOM::OPTION TOOL, PivotMatrixType *H=NULL) |
void | skim (DMRG::DIRECTION::OPTION DIR, DMRG::BROOM::OPTION TOOL, PivotMatrixType *H=NULL) |
void | skim (DMRG::BROOM::OPTION TOOL, PivotMatrixType *H=NULL) |
void | entropy_skim () |
void | sweepStep (DMRG::DIRECTION::OPTION DIR, size_t loc, DMRG::BROOM::OPTION TOOL, PivotMatrixType *H=NULL, bool DISCARD_U_or_V=false) |
virtual void | rightSweepStep (size_t loc, DMRG::BROOM::OPTION TOOL, PivotMatrixType *H=NULL, bool DISCARD_V=false) |
virtual void | leftSweepStep (size_t loc, DMRG::BROOM::OPTION TOOL, PivotMatrixType *H=NULL, bool DISCARD_U=false) |
size_t | length () const |
Public Attributes | |
MpsBoundaries< Symmetry, Scalar > | Boundaries |
size_t | N_phys |
vector< vector< qarray< Nq > > > | qloc |
qarray< Nq > | Qtot = Symmetry::qvacuum() |
vector< qarray< Nq > > | Qmulti |
vector< vector< Biped< Symmetry, MatrixType > > > | A |
ArrayXd | truncWeight |
ArrayXd | S |
vector< map< qarray< Nq >, ArrayXd > > | SVspec |
vector< Qbasis< Symmetry > > | inbase |
vector< Qbasis< Symmetry > > | outbase |
vector< qarray< Nq > > | QinTop |
vector< qarray< Nq > > | QinBot |
vector< qarray< Nq > > | QoutTop |
vector< qarray< Nq > > | QoutBot |
Public Attributes inherited from DmrgJanitor< PivotMatrixType > | |
double | eps_svd |
double | eps_truncWeight |
double | alpha_rsvd |
size_t | max_Nsv |
size_t | min_Nsv |
int | max_Nrich |
Private Types | |
typedef Matrix< Scalar, Dynamic, Dynamic > | MatrixType |
typedef Symmetry::qType | qType |
Static Private Attributes | |
static constexpr size_t | Nq = Symmetry::Nq |
Friends | |
template<typename Symmetry_ , typename MpHamiltonian , typename Scalar_ > | |
class | DmrgSolver |
template<typename Symmetry_ , typename S1 , typename S2 > | |
class | MpsCompressor |
template<typename H , typename Symmetry_ , typename S1 , typename S2 , typename V > | |
class | TDVPPropagator |
template<typename Symmetry_ , typename S_ > | |
class | Mps |
template<typename Symmetry_ , typename S1 , typename S2 > | |
void | HxV (const Mpo< Symmetry_, S1 > &H, const Mps< Symmetry_, S2 > &Vin, Mps< Symmetry_, S2 > &Vout, DMRG::VERBOSITY::OPTION VERBOSITY) |
template<typename Symmetry_ , typename S1 , typename S2 > | |
void | OxV (const Mpo< Symmetry_, S1 > &H, const Mps< Symmetry_, S2 > &Vin, Mps< Symmetry_, S2 > &Vout, DMRG::BROOM::OPTION TOOL) |
template<typename Symmetry_ , typename S1 , typename S2 > | |
void | OxV_exact (const Mpo< Symmetry_, S1 > &H, const Mps< Symmetry_, S2 > &Vin, Mps< Symmetry_, S2 > &Vout, double tol_compr, DMRG::VERBOSITY::OPTION VERBOSITY, int max_halfsweeps, int min_halfsweeps, int Minit) |
Additional Inherited Members | |
Protected Attributes inherited from DmrgJanitor< PivotMatrixType > | |
int | pivot = -1 |
size_t | N_sites |
|
private |
|
private |
Does nothing.
Mps< Symmetry, Scalar >::Mps | ( | size_t | L_input, |
vector< vector< qarray< Nq > > > | qloc_input, | ||
qarray< Nq > | Qtot_input, | ||
size_t | N_phys_input, | ||
int | Nqmax_input, | ||
bool | TRIVIAL_BOUNDARIES = true |
||
) |
Construct by setting all the relevant parameters.
L_input | : chain length |
qloc_input | : local basis |
Qtot_input | : target quantum number |
N_phys_input | : the volume of the system (normally (chain length) * (chain width)) |
Nqmax_input | : maximal initial number of symmetry blocks per site in the Mps |
Mps< Symmetry, Scalar >::Mps | ( | size_t | L_input, |
const vector< vector< Biped< Symmetry, MatrixType > > > & | As, | ||
const vector< vector< qarray< Nq > > > & | qloc_input, | ||
qarray< Nq > | Qtot_input, | ||
size_t | N_phys_input | ||
) |
Construct by explicitly providing the A-matrices. Basically only for testing purposes.
L_input | : chain length |
As | : vector of vectors of A matrices ([site][local basis index]) |
qloc_input | : vector of local bases for all sites |
Qtot_input | : target quantum number |
N_phys_input | : the volume of the system (normally (chain length) * (chain width)) |
|
inline |
|
inline |
void Mps< Symmetry, Scalar >::absorb | ( | size_t | loc, |
DMRG::DIRECTION::OPTION | DIR, | ||
const Biped< Symmetry, MatrixType > & | C | ||
) |
Absorbs the zero-site tensor C
(as obtained after an SVD split) into the Mps. Used in TDVPPropagator.
loc | : site to do the absorption |
DIR | : specifies whether the absorption is on the left-sweep or right-sweep |
C | : the zero-site tensor to be absorbed |
void Mps< Symmetry, Scalar >::add_site | ( | size_t | loc, |
OtherScalar | alpha, | ||
const Mps< Symmetry, Scalar > & | Vin | ||
) |
Adds one site at a time in addScale, conserving memory.
void Mps< Symmetry, Scalar >::addScale | ( | OtherScalar | alpha, |
const Mps< Symmetry, Scalar > & | Vin, | ||
bool | SVD_COMPRESS = false |
||
) |
Adds another Mps to the given one and scales by alpha
, i.e. performs .
alpha | : scalar for scaling |
Vin | : Mps to be added |
SVD_COMPRESS | : If true , the resulting Mps is compressed using SVD. If false , the summation is exact (direct sum of the matrices). |
void Mps< Symmetry, Scalar >::applyGate | ( | const TwoSiteGate< Symmetry, Scalar > & | gate, |
size_t | l, | ||
DMRG::DIRECTION::OPTION | DIR | ||
) |
Act with a two-site gate on the MPS at sites l
and l+1
.
string Mps< Symmetry, Scalar >::Asizes | ( | ) | const |
Prints the sizes of the matrices for testing purposes.
size_t Mps< Symmetry, Scalar >::calc_Dmax | ( | ) | const |
Determines the maximal amount of rows or columns per site and per subspace.
size_t Mps< Symmetry, Scalar >::calc_fullMmax | ( | ) | const |
For SU(2) symmetries, determines the equivalent U(1) bond dimension.
size_t Mps< Symmetry, Scalar >::calc_Mmax | ( | ) | const |
Determines the maximal bond dimension per site (sum of A.rows
or A.cols
over all subspaces).
void Mps< Symmetry, Scalar >::calc_N | ( | DMRG::DIRECTION::OPTION | DIR, |
size_t | loc, | ||
vector< Biped< Symmetry, MatrixType > > & | N | ||
) |
Calculates the nullspace of the site tensor on site loc
when blocked into direction DIR
.
DIR | : direction of the sweep, either LEFT or RIGHT. |
loc | : site to perform the sweep on; afterwards the pivot is shifted to loc-1 or loc+1 |
N | : tensor to write the nullsapce to. |
double Mps< Symmetry, Scalar >::calc_Nqavg | ( | ) | const |
Determines the average amount of subspaces per site.
size_t Mps< Symmetry, Scalar >::calc_Nqmax | ( | ) | const |
Determines the maximal amount of subspaces per site.
void Mps< Symmetry, Scalar >::calc_Qlimits | ( | ) |
Calculate quantum number bounds.
void Mps< Symmetry, Scalar >::canonize | ( | DMRG::DIRECTION::OPTION | DIR = DMRG::DIRECTION::LEFT | ) |
Sweeps through the chain with DMRG::BROOM::QR, creating a canonical Mps.
DIR | : If DMRG::DIRECTION::LEFT, the result is left-canonical. If DMRG::DIRECTION::RIGHT, the result is right-canonical. |
Mps< Symmetry, OtherScalar > Mps< Symmetry, Scalar >::cast | ( | ) | const |
Casts the matrices from Scalar
to OtherScalar
.
void Mps< Symmetry, Scalar >::collapse | ( | ) |
For METTS.
Scalar Mps< Symmetry, Scalar >::dot | ( | const Mps< Symmetry, Scalar > & | Vket | ) | const |
Calculates the scalar product with another Mps.
|
inline |
void Mps< Symmetry, Scalar >::enrich_left | ( | size_t | loc, |
PivotMatrix1< Symmetry, Scalar, Scalar > * | H | ||
) |
Enriches the search space in order to dynamically find the right bond dimension. Used in sweeps with the option RICH_SVD.
void Mps< Symmetry, Scalar >::enrich_right | ( | size_t | loc, |
PivotMatrix1< Symmetry, Scalar, Scalar > * | H | ||
) |
std::pair< vector< qarray< Symmetry::Nq > >, ArrayXd > Mps< Symmetry, Scalar >::entanglementSpectrumLoc | ( | size_t | loc | ) | const |
Return the entanglement spectrum at the site loc
(values all subspaces merged and sorted).
|
inline |
|
inline |
void Mps< Symmetry, Scalar >::get_controlParams | ( | const Mps< Symmetry, Scalar > & | V | ) |
Copies the control parameters from another Mps, i.e. all the cutoff tolerances specified in DmrgJanitor.
size_t Mps< Symmetry, Scalar >::get_Min | ( | size_t | loc | ) | const |
Sets all matrices to random using boost's uniform distribution from -1 to 1.
size_t Mps< Symmetry, Scalar >::get_Mout | ( | size_t | loc | ) | const |
Sets all matrices to random using boost's uniform distribution from -1 to 1.
|
inline |
|
inline |
void Mps< Symmetry, Scalar >::graph | ( | string | filename | ) | const |
Writes the subspace connections as a directed graph into a file. Run
dot filename.dot -Tpdf -o filename.pdf
to create a shiny pdf.
filename | : gets a ".dot" extension automatically |
string Mps< Symmetry, Scalar >::info | ( | ) | const |
Prints some info about this class.
void Mps< Symmetry, Scalar >::innerResize | ( | size_t | Mmax | ) |
Resizes the block matrices.
Mmax | : size cutoff |
Mmax
. void Mps< Symmetry, Scalar >::leftSplitStep | ( | size_t | loc, |
Biped< Symmetry, MatrixType > & | C | ||
) |
Performs an SVD split to the left and writes the zero-site tensor to C
. Used in TDVPPropagator.
void Mps< Symmetry, Scalar >::leftSweepStep | ( | size_t | loc, |
DMRG::BROOM::OPTION | BROOM, | ||
PivotMatrix1< Symmetry, Scalar, Scalar > * | H = NULL , |
||
bool | DISCARD_U = false |
||
) |
Performs a sweep step to the left.
loc | : site to perform the sweep on; afterwards the pivot is shifted to loc-1 |
BROOM | : choice of decomposition |
H | : non-local information from transfer matrices is provided here when BROOM is DMRG::BROOM::RICH_SVD |
DISCARD_U | : if true , don't multiply the U-matrix onto the next site |
Scalar Mps< Symmetry, Scalar >::locAvg | ( | const Mpo< Symmetry, MpoScalar > & | O, |
size_t | distance = 0 |
||
) | const |
Calculates the expectation value with a local operator at the pivot site.
O | : local Mpo acting on the pivot side. |
distance | : distance to the end of the support of O . |
double Mps< Symmetry, Scalar >::memory | ( | MEMUNIT | memunit = GB | ) | const |
Calculates the (theoretically) allocated memory (note: by default in GB).
void Mps< Symmetry, Scalar >::mend | ( | ) |
Finds broken paths through the quantum number subspaces and mends them by resizing with appropriate zeros. This is needed when applying an Mpo which changes quantum numbers, making some paths impossible. For example, one can add a particle at the beginning or end of the chain with the same target particle number, but if an annihilator is applied in the middle, only the first path survives.
|
inline |
Sets all matrices to random using boost's uniform distribution from -1 to 1.
Mps< Symmetry, Scalar > & Mps< Symmetry, Scalar >::operator*= | ( | const OtherScalar & | alpha | ) |
Performs . Applies it to the pivot site (if non-orthogonal, to the first site).
Mps< Symmetry, Scalar > & Mps< Symmetry, Scalar >::operator+= | ( | const Mps< Symmetry, Scalar > & | Vin | ) |
Performs Mps::addScale with alpha
= +1.
Mps< Symmetry, Scalar > & Mps< Symmetry, Scalar >::operator-= | ( | const Mps< Symmetry, Scalar > & | Vin | ) |
Performs Mps::addScale with alpha
= -1.
Mps< Symmetry, Scalar > & Mps< Symmetry, Scalar >::operator/= | ( | const OtherScalar & | alpha | ) |
Performs . Applies it to the pivot site (if non-orthogonal, to the first site).
void Mps< Symmetry, Scalar >::outerResize | ( | const Hamiltonian & | H, |
qarray< Nq > | Qtot_input, | ||
int | Nqmax_input = 500 |
||
) |
Determines all subspace quantum numbers and resizes the containers for the blocks. Memory for the matrices remains uninitiated. Pulls info from an Mpo.
void Mps< Symmetry, Scalar >::outerResize | ( | size_t | L_input, |
vector< vector< qarray< Nq > > > | qloc_input, | ||
qarray< Nq > | Qtot_input, | ||
int | Nqmax_input = 500 |
||
) |
Determines all subspace quantum numbers and resizes the containers for the blocks. Memory for the matrices remains uninitialized.
L_input | : chain length |
qloc_input | : local basis |
Qtot_input | : target quantum number |
Nqmax_input | : maximum initial number of symmetry blocks in the Mps per site |
void Mps< Symmetry, Scalar >::outerResizeAll | ( | size_t | L_input, |
vector< vector< qarray< Nq > > > | qloc_input, | ||
qarray< Nq > | Qtot_input | ||
) |
Resizes with all possible blocks.
L_input | : chain length |
qloc_input | : local basis |
Qtot_input | : target quantum number |
void Mps< Symmetry, Scalar >::outerResizeNoSymm | ( | ) |
void Mps< Symmetry, Scalar >::resize_arrays | ( | ) |
Shorthand to resize all the relevant arrays: A
, inbase
, outbase
, truncWeight
, S
.
void Mps< Symmetry, Scalar >::rightSplitStep | ( | size_t | loc, |
Biped< Symmetry, MatrixType > & | C | ||
) |
Performs an SVD split to the right and writes the zero-site tensor to C
. Used in TDVPPropagator.
void Mps< Symmetry, Scalar >::rightSweepStep | ( | size_t | loc, |
DMRG::BROOM::OPTION | BROOM, | ||
PivotMatrix1< Symmetry, Scalar, Scalar > * | H = NULL , |
||
bool | DISCARD_V = false |
||
) |
Performs a sweep step to the right.
loc | : site to perform the sweep on; afterwards the pivot is shifted to loc+1 |
BROOM | : choice of decomposition |
H | : non-local information from transfer matrices is provided here when BROOM is DMRG::BROOM::RICH_SVD |
DISCARD_V | : if true , don't multiply the V-matrix onto the next site |
void Mps< Symmetry, Scalar >::set_A_from_C | ( | size_t | loc, |
const vector< Tripod< Symmetry, MatrixType > > & | C, | ||
DMRG::BROOM::OPTION | TOOL = DMRG::BROOM::SVD |
||
) |
Sets the A-matrix at a given site by performing SVD on the C-tensor.
|
inline |
void Mps< Symmetry, Scalar >::set_Qlimits_to_inf | ( | ) |
Set quantum number bounds to +-infinity.
|
inline |
Sets all matrices to random using boost's uniform distribution from -1 to 1.
void Mps< Symmetry, Scalar >::setProductState | ( | const Hamiltonian & | H, |
const vector< qarray< Nq > > & | config | ||
) |
Sets the Mps from a product state configuration.
H | : Hamiltonian, needed for Mps::outerResize |
config | : classical configuration, a vector of qarray , e.g. (+1,-1,+1,-1,...) for the Néel state |
void Mps< Symmetry, Scalar >::setRandom | ( | ) |
Sets all matrices to random using boost's uniform distribution from -1 to 1.
void Mps< Symmetry, Scalar >::setRandom | ( | size_t | loc | ) |
Sets all matrices at site loc
to random using boost's uniform distribution from -1 to 1.
void Mps< Symmetry, Scalar >::setZero | ( | ) |
Sets all matrices to zero.
|
inline |
double Mps< Symmetry, Scalar >::squaredNorm | ( | ) | const |
Calculates the squared norm. Exploits the canonical form if possible, calculates the dot product with itself otherwise.
void Mps< Symmetry, Scalar >::swap | ( | Mps< Symmetry, Scalar > & | V | ) |
Swaps with another Mps.
void Mps< Symmetry, Scalar >::sweepStep2 | ( | DMRG::DIRECTION::OPTION | DIR, |
size_t | loc, | ||
const vector< Biped< Symmetry, MatrixType > > & | Apair, | ||
bool | SEPARATE_SV = false |
||
) |
Performs a two-site sweep.
DIR | : direction of the sweep, either LEFT or RIGHT. |
loc | : site to perform the sweep on; afterwards the pivot is shifted to loc-1 or loc+1 |
Apair | : pair of two Mps site tensors which are split via an SVD |
SEPARATE_SV | if true , the singular value matrix is discarded (iseful for IDMRG) |
string Mps< Symmetry, Scalar >::test_ortho | ( | double | tol = 1e-8 | ) | const |
Tests the orthogonality of the Mps. Returns a string with "A"=left-canonical ( ), "B"=right-canonical ( ), "X"=both, "M"=neither; with the pivot site underlined.
tol | : The check is |
|
inline |
void Mps< Symmetry, Scalar >::update_inbase | ( | size_t | loc | ) |
Update the bases in case new blocks have appeared or old ones have disappeared
|
inline |
void Mps< Symmetry, Scalar >::update_outbase | ( | size_t | loc | ) |
string Mps< Symmetry, Scalar >::validate | ( | string | name = "Mps< Symmetry, Scalar >" | ) | const |
Checks if the sizes of the block matrices are consistent, so that they can be multiplied.
name | : how to call the Mps in the output |
|
friend |
|
friend |
|
friend |
|
friend |
|
friend |
|
friend |
|
friend |
vector<vector<Biped<Symmetry,MatrixType> > > Mps< Symmetry, Scalar >::A |
MpsBoundaries<Symmetry,Scalar> Mps< Symmetry, Scalar >::Boundaries |
size_t Mps< Symmetry, Scalar >::N_phys |
|
staticconstexprprivate |
ArrayXd Mps< Symmetry, Scalar >::S |
ArrayXd Mps< Symmetry, Scalar >::truncWeight |