VMPS++
Loading...
Searching...
No Matches
Mps< Symmetry, Scalar > Class Template Reference

Detailed Description

template<typename Symmetry, typename Scalar = double>
class Mps< Symmetry, Scalar >

Matrix Product State with conserved quantum numbers (Abelian and non abelian symmetries).

Template Parameters
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>
Note
We define the quantum number flow for the $A$-tensors as follows: left auxiliary leg $i$ is combined with the physical index $\sigma$ to obtain the right auxiliary leg $j$. For U(1) this means that $i+\sigma=j $. For SU(2) this means that the $A$-tensor decompose with the CGC $C^{i,\sigma\rightarrow j}_{m_i,m_\sigma\rightarrow m_j}$.

Definition at line 38 of file Mps.h.

#include <Mps.h>

Inheritance diagram for Mps< Symmetry, Scalar >:

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< NqQtarget () 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< NqQtot = 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
 

Member Typedef Documentation

◆ MatrixType

template<typename Symmetry , typename Scalar = double>
typedef Matrix<Scalar,Dynamic,Dynamic> Mps< Symmetry, Scalar >::MatrixType
private

Definition at line 40 of file Mps.h.

◆ qType

template<typename Symmetry , typename Scalar = double>
typedef Symmetry::qType Mps< Symmetry, Scalar >::qType
private

Definition at line 42 of file Mps.h.

Constructor & Destructor Documentation

◆ Mps() [1/4]

template<typename Symmetry , typename Scalar = double>
Mps< Symmetry, Scalar >::Mps ( )

Does nothing.

◆ Mps() [2/4]

template<typename Symmetry , typename Scalar = double>
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.

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() [3/4]

template<typename Symmetry , typename Scalar = double>
template<typename Hamiltonian >
Mps< Symmetry, Scalar >::Mps ( const Hamiltonian &  H,
size_t  Mmax,
qarray< Nq Qtot_input,
int  Nqmax_input 
)

Construct by pulling info from an Mpo.

Parameters
H: chain length and local basis will be retrieved from this Mpo
Mmax: size cutoff
Qtot_input: target quantum number
Nqmax_input: maximal initial number of symmetry blocks per site in the Mps

◆ Mps() [4/4]

template<typename Symmetry , typename Scalar = double>
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.

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

Member Function Documentation

◆ A_at() [1/2]

template<typename Symmetry , typename Scalar = double>
vector< Biped< Symmetry, MatrixType > > & Mps< Symmetry, Scalar >::A_at ( size_t  loc)
inline

Reference to the A-tensor at site loc.

Definition at line 459 of file Mps.h.

◆ A_at() [2/2]

template<typename Symmetry , typename Scalar = double>
const vector< Biped< Symmetry, MatrixType > > & Mps< Symmetry, Scalar >::A_at ( size_t  loc) const
inline

Const reference to the A-tensor at site loc.

Definition at line 456 of file Mps.h.

◆ absorb()

template<typename Symmetry , typename Scalar = double>
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.

Parameters
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

◆ add_site()

template<typename Symmetry , typename Scalar = double>
template<typename OtherScalar >
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.

◆ addScale()

template<typename Symmetry , typename Scalar = double>
template<typename OtherScalar >
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 $ \mathrel{+}= \alpha \cdot V_{in}$.

Parameters
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).
Warning
Not implemented for non-abelian symmetries.

◆ applyGate()

template<typename Symmetry , typename Scalar = double>
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.

◆ Asizes()

template<typename Symmetry , typename Scalar = double>
string Mps< Symmetry, Scalar >::Asizes ( ) const

Prints the sizes of the matrices for testing purposes.

◆ calc_Dmax()

template<typename Symmetry , typename Scalar = double>
size_t Mps< Symmetry, Scalar >::calc_Dmax ( ) const

Determines the maximal amount of rows or columns per site and per subspace.

◆ calc_fullMmax()

template<typename Symmetry , typename Scalar = double>
size_t Mps< Symmetry, Scalar >::calc_fullMmax ( ) const

For SU(2) symmetries, determines the equivalent U(1) bond dimension.

◆ calc_Mmax()

template<typename Symmetry , typename Scalar = double>
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).

◆ calc_N()

template<typename Symmetry , typename Scalar = double>
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.

Parameters
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.
Note
The nullspace is obtained by a full QR decomposition.
The nullspace is used for error estimation as suggested here: arXiv:1711.01104

◆ calc_Nqavg()

template<typename Symmetry , typename Scalar = double>
double Mps< Symmetry, Scalar >::calc_Nqavg ( ) const

Determines the average amount of subspaces per site.

◆ calc_Nqmax()

template<typename Symmetry , typename Scalar = double>
size_t Mps< Symmetry, Scalar >::calc_Nqmax ( ) const

Determines the maximal amount of subspaces per site.

◆ calc_Qlimits()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::calc_Qlimits ( )

Calculate quantum number bounds.

◆ canonize()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::canonize ( DMRG::DIRECTION::OPTION  DIR = DMRG::DIRECTION::LEFT)

Sweeps through the chain with DMRG::BROOM::QR, creating a canonical Mps.

Parameters
DIR: If DMRG::DIRECTION::LEFT, the result is left-canonical. If DMRG::DIRECTION::RIGHT, the result is right-canonical.

◆ cast()

template<typename Symmetry , typename Scalar = double>
template<typename OtherScalar >
Mps< Symmetry, OtherScalar > Mps< Symmetry, Scalar >::cast ( ) const

Casts the matrices from Scalar to OtherScalar.

◆ collapse()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::collapse ( )

For METTS.

Warning
Not tested and soon abandoned.

◆ dot()

template<typename Symmetry , typename Scalar = double>
Scalar Mps< Symmetry, Scalar >::dot ( const Mps< Symmetry, Scalar > &  Vket) const

Calculates the scalar product with another Mps.

◆ elongate_hetero()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::elongate_hetero ( size_t  Nleft = 0,
size_t  Nright = 0 
)
inline

Definition at line 506 of file Mps.h.

◆ enrich_left()

template<typename Symmetry , typename Scalar = double>
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.

◆ enrich_right()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::enrich_right ( size_t  loc,
PivotMatrix1< Symmetry, Scalar, Scalar > *  H 
)

◆ entanglementSpectrum()

template<typename Symmetry , typename Scalar = double>
vector< map< qarray< Nq >, ArrayXd > > Mps< Symmetry, Scalar >::entanglementSpectrum ( ) const
inline

Return the full entanglement spectrum, resolved by subspace quantum number.

Definition at line 471 of file Mps.h.

◆ entanglementSpectrumLoc()

template<typename Symmetry , typename Scalar = double>
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).

◆ entropy()

template<typename Symmetry , typename Scalar = double>
ArrayXd Mps< Symmetry, Scalar >::entropy ( ) const
inline

Returns the entropy for all bonds.

Definition at line 468 of file Mps.h.

◆ get_boundaryTensor()

template<typename Symmetry , typename Scalar = double>
Tripod< Symmetry, Matrix< Scalar, Dynamic, Dynamic > > Mps< Symmetry, Scalar >::get_boundaryTensor ( DMRG::DIRECTION::OPTION  DIR,
size_t  usePower = 1ul 
) const
inline

Definition at line 488 of file Mps.h.

◆ get_controlParams()

template<typename Symmetry , typename Scalar = double>
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.

◆ get_Min()

template<typename Symmetry , typename Scalar = double>
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.

Warning
Watch for overflow in large chains where one gets exponentially large values when multiplying all the matrices! The safer way is to randomize while sweeping, using Mps::setRandom(size_t loc).

◆ get_Mout()

template<typename Symmetry , typename Scalar = double>
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.

Warning
Watch for overflow in large chains where one gets exponentially large values when multiplying all the matrices! The safer way is to randomize while sweeping, using Mps::setRandom(size_t loc).

◆ get_pivot()

template<typename Symmetry , typename Scalar = double>
int Mps< Symmetry, Scalar >::get_pivot ( ) const
inline

Returns the pivot position.

Definition at line 462 of file Mps.h.

◆ get_truncWeight()

template<typename Symmetry , typename Scalar = double>
ArrayXd Mps< Symmetry, Scalar >::get_truncWeight ( ) const
inline

Returns the truncated weight for all sites.

Definition at line 465 of file Mps.h.

◆ graph()

template<typename Symmetry , typename Scalar = double>
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.

Parameters
filename: gets a ".dot" extension automatically

◆ inBasis() [1/2]

template<typename Symmetry , typename Scalar = double>
vector< Qbasis< Symmetry > > Mps< Symmetry, Scalar >::inBasis ( ) const
inline

Returns the target quantum number.

Definition at line 449 of file Mps.h.

◆ inBasis() [2/2]

template<typename Symmetry , typename Scalar = double>
Qbasis< Symmetry > Mps< Symmetry, Scalar >::inBasis ( size_t  loc) const
inline

Returns the auxiliary ingoing basis.

Definition at line 448 of file Mps.h.

◆ info()

template<typename Symmetry , typename Scalar = double>
string Mps< Symmetry, Scalar >::info ( ) const

Prints some info about this class.

◆ innerResize()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::innerResize ( size_t  Mmax)

Resizes the block matrices.

Parameters
Mmax: size cutoff
Note
The edges will have the analytically exact size, which may be smaller than Mmax.

◆ leftSplitStep()

template<typename Symmetry , typename Scalar = double>
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.

◆ leftSweepStep()

template<typename Symmetry , typename Scalar = double>
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.

Parameters
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

◆ locAvg()

template<typename Symmetry , typename Scalar = double>
template<typename MpoScalar >
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.

Parameters
O: local Mpo acting on the pivot side.
distance: distance to the end of the support of O.

◆ locBasis() [1/2]

template<typename Symmetry , typename Scalar = double>
vector< vector< qarray< Nq > > > Mps< Symmetry, Scalar >::locBasis ( ) const
inline

Returns the target quantum number.

Definition at line 445 of file Mps.h.

◆ locBasis() [2/2]

template<typename Symmetry , typename Scalar = double>
vector< qarray< Nq > > Mps< Symmetry, Scalar >::locBasis ( size_t  loc) const
inline

Returns the local basis.

Definition at line 444 of file Mps.h.

◆ memory()

template<typename Symmetry , typename Scalar = double>
double Mps< Symmetry, Scalar >::memory ( MEMUNIT  memunit = GB) const

Calculates the (theoretically) allocated memory (note: by default in GB).

◆ mend()

template<typename Symmetry , typename Scalar = double>
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.

Warning
: deprecated!

◆ Nqout()

template<typename Symmetry , typename Scalar = double>
size_t Mps< Symmetry, Scalar >::Nqout ( size_t  l)
inline

Sets all matrices to random using boost's uniform distribution from -1 to 1.

Warning
Watch for overflow in large chains where one gets exponentially large values when multiplying all the matrices! The safer way is to randomize while sweeping, using Mps::setRandom(size_t loc).

Definition at line 273 of file Mps.h.

◆ operator*=()

template<typename Symmetry , typename Scalar = double>
template<typename OtherScalar >
Mps< Symmetry, Scalar > & Mps< Symmetry, Scalar >::operator*= ( const OtherScalar &  alpha)

Performs $ \mathrel{*}= \alpha$. Applies it to the pivot site (if non-orthogonal, to the first site).

◆ operator+=()

template<typename Symmetry , typename Scalar = double>
Mps< Symmetry, Scalar > & Mps< Symmetry, Scalar >::operator+= ( const Mps< Symmetry, Scalar > &  Vin)

Performs Mps::addScale with alpha = +1.

Warning
Not implemented for non-abelian symmetries.

◆ operator-=()

template<typename Symmetry , typename Scalar = double>
Mps< Symmetry, Scalar > & Mps< Symmetry, Scalar >::operator-= ( const Mps< Symmetry, Scalar > &  Vin)

Performs Mps::addScale with alpha = -1.

Warning
Not implemented for non-abelian symmetries.

◆ operator/=()

template<typename Symmetry , typename Scalar = double>
template<typename OtherScalar >
Mps< Symmetry, Scalar > & Mps< Symmetry, Scalar >::operator/= ( const OtherScalar &  alpha)

Performs $ \mathrel{/}= \alpha$. Applies it to the pivot site (if non-orthogonal, to the first site).

◆ outBasis() [1/2]

template<typename Symmetry , typename Scalar = double>
vector< Qbasis< Symmetry > > Mps< Symmetry, Scalar >::outBasis ( ) const
inline

Returns the target quantum number.

Definition at line 453 of file Mps.h.

◆ outBasis() [2/2]

template<typename Symmetry , typename Scalar = double>
Qbasis< Symmetry > Mps< Symmetry, Scalar >::outBasis ( size_t  loc) const
inline

Returns the auxiliary outgoing basis.

Definition at line 452 of file Mps.h.

◆ outerResize() [1/3]

template<typename Symmetry , typename Scalar = double>
template<typename Hamiltonian >
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.

Parameters
H: chain length and local basis will be retrieved from this Mpo
Qtot_input: target quantum number
Nqmax_input: maximum initial number of symmetry blocks in the Mps per site

◆ outerResize() [2/3]

template<typename Symmetry , typename Scalar = double>
template<typename OtherMatrixType >
void Mps< Symmetry, Scalar >::outerResize ( const Mps< Symmetry, OtherMatrixType > &  V)

Determines all subspace quantum numbers and resizes the containers for the blocks. Memory for the matrices remains uninitiated. Pulls info from another Mps.

Parameters
V: chain length, local basis and target quantum number will be equal to this Mps

◆ outerResize() [3/3]

template<typename Symmetry , typename Scalar = double>
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.

Parameters
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

◆ outerResizeAll()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::outerResizeAll ( size_t  L_input,
vector< vector< qarray< Nq > > >  qloc_input,
qarray< Nq Qtot_input 
)

Resizes with all possible blocks.

Parameters
L_input: chain length
qloc_input: local basis
Qtot_input: target quantum number

◆ outerResizeNoSymm()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::outerResizeNoSymm ( )

◆ Qmultitarget()

template<typename Symmetry , typename Scalar = double>
vector< qarray< Nq > > Mps< Symmetry, Scalar >::Qmultitarget ( ) const
inline

Returns the multi-target quantum number for spectral functions.

Definition at line 441 of file Mps.h.

◆ Qtarget()

template<typename Symmetry , typename Scalar = double>
qarray< Nq > Mps< Symmetry, Scalar >::Qtarget ( ) const
inline

Returns the target quantum number.

Definition at line 438 of file Mps.h.

◆ resize_arrays()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::resize_arrays ( )

Shorthand to resize all the relevant arrays: A, inbase, outbase, truncWeight, S.

◆ rightSplitStep()

template<typename Symmetry , typename Scalar = double>
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.

◆ rightSweepStep()

template<typename Symmetry , typename Scalar = double>
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.

Parameters
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

◆ set_A_from_C()

template<typename Symmetry , typename Scalar = double>
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.

Warning
Not implemented for non-abelian symmetries, deprecated!

◆ set_open_bc()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::set_open_bc ( )
inline

Definition at line 483 of file Mps.h.

◆ set_Qlimits_to_inf()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::set_Qlimits_to_inf ( )

Set quantum number bounds to +-infinity.

◆ set_Qmultitarget()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::set_Qmultitarget ( const vector< qarray< Nq > > &  Qmulti_input)
inline

Sets all matrices to random using boost's uniform distribution from -1 to 1.

Warning
Watch for overflow in large chains where one gets exponentially large values when multiplying all the matrices! The safer way is to randomize while sweeping, using Mps::setRandom(size_t loc).

Definition at line 206 of file Mps.h.

◆ setProductState()

template<typename Symmetry , typename Scalar = double>
template<typename Hamiltonian >
void Mps< Symmetry, Scalar >::setProductState ( const Hamiltonian &  H,
const vector< qarray< Nq > > &  config 
)

Sets the Mps from a product state configuration.

Parameters
H: Hamiltonian, needed for Mps::outerResize
config: classical configuration, a vector of qarray, e.g. (+1,-1,+1,-1,...) for the Néel state

◆ setRandom() [1/2]

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::setRandom ( )

Sets all matrices to random using boost's uniform distribution from -1 to 1.

Warning
Watch for overflow in large chains where one gets exponentially large values when multiplying all the matrices! The safer way is to randomize while sweeping, using Mps::setRandom(size_t loc).

◆ setRandom() [2/2]

template<typename Symmetry , typename Scalar = double>
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.

◆ setZero()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::setZero ( )

Sets all matrices to zero.

◆ shift_hetero()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::shift_hetero ( int  Nshift = 0)
inline

Definition at line 551 of file Mps.h.

◆ squaredNorm()

template<typename Symmetry , typename Scalar = double>
double Mps< Symmetry, Scalar >::squaredNorm ( ) const

Calculates the squared norm. Exploits the canonical form if possible, calculates the dot product with itself otherwise.

◆ swap()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::swap ( Mps< Symmetry, Scalar > &  V)

Swaps with another Mps.

◆ sweepStep2()

template<typename Symmetry , typename Scalar = double>
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.

Parameters
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_SVif true, the singular value matrix is discarded (iseful for IDMRG)

◆ test_ortho()

template<typename Symmetry , typename Scalar = double>
string Mps< Symmetry, Scalar >::test_ortho ( double  tol = 1e-8) const

Tests the orthogonality of the Mps. Returns a string with "A"=left-canonical ( $\sum_s {A^s}^\dag A^s=I$), "B"=right-canonical ( $\sum_s B^s {B^s}^\dag=I$), "X"=both, "M"=neither; with the pivot site underlined.

Parameters
tol: The check is $\|\sum_s {A^s}^\dag A^s-I\|_{\infty} < tol$

◆ transform_base()

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::transform_base ( qarray< Symmetry::Nq >  Qtot,
int  L,
bool  PRINT = false 
)
inline

Definition at line 593 of file Mps.h.

◆ update_inbase() [1/2]

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::update_inbase ( )
inline

Definition at line 652 of file Mps.h.

◆ update_inbase() [2/2]

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::update_inbase ( size_t  loc)

Update the bases in case new blocks have appeared or old ones have disappeared

◆ update_outbase() [1/2]

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::update_outbase ( )
inline

Definition at line 653 of file Mps.h.

◆ update_outbase() [2/2]

template<typename Symmetry , typename Scalar = double>
void Mps< Symmetry, Scalar >::update_outbase ( size_t  loc)

◆ validate()

template<typename Symmetry , typename Scalar = double>
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.

Parameters
name: how to call the Mps in the output
Returns
: string with info

Friends And Related Function Documentation

◆ DmrgSolver

template<typename Symmetry , typename Scalar = double>
template<typename Symmetry_ , typename MpHamiltonian , typename Scalar_ >
friend class DmrgSolver
friend

Definition at line 45 of file Mps.h.

◆ HxV

template<typename Symmetry , typename Scalar = double>
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 
)
friend

◆ Mps

template<typename Symmetry , typename Scalar = double>
template<typename Symmetry_ , typename S_ >
friend class Mps
friend

Definition at line 55 of file Mps.h.

◆ MpsCompressor

template<typename Symmetry , typename Scalar = double>
template<typename Symmetry_ , typename S1 , typename S2 >
friend class MpsCompressor
friend

Definition at line 46 of file Mps.h.

◆ OxV

template<typename Symmetry , typename Scalar = double>
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 
)
friend

◆ OxV_exact

template<typename Symmetry , typename Scalar = double>
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 
)
friend

◆ TDVPPropagator

template<typename Symmetry , typename Scalar = double>
template<typename H , typename Symmetry_ , typename S1 , typename S2 , typename V >
friend class TDVPPropagator
friend

Definition at line 47 of file Mps.h.

Member Data Documentation

◆ A

template<typename Symmetry , typename Scalar = double>
vector<vector<Biped<Symmetry,MatrixType> > > Mps< Symmetry, Scalar >::A

A-tensor

Definition at line 624 of file Mps.h.

◆ Boundaries

template<typename Symmetry , typename Scalar = double>
MpsBoundaries<Symmetry,Scalar> Mps< Symmetry, Scalar >::Boundaries

Definition at line 481 of file Mps.h.

◆ inbase

template<typename Symmetry , typename Scalar = double>
vector<Qbasis<Symmetry> > Mps< Symmetry, Scalar >::inbase

bases on all ingoing and outgoing legs of the Mps

Definition at line 635 of file Mps.h.

◆ N_phys

template<typename Symmetry , typename Scalar = double>
size_t Mps< Symmetry, Scalar >::N_phys

volume of the system (normally (chain length) * (chain width))

Definition at line 612 of file Mps.h.

◆ Nq

template<typename Symmetry , typename Scalar = double>
constexpr size_t Mps< Symmetry, Scalar >::Nq = Symmetry::Nq
staticconstexprprivate

Definition at line 41 of file Mps.h.

◆ outbase

template<typename Symmetry , typename Scalar = double>
vector<Qbasis<Symmetry> > Mps< Symmetry, Scalar >::outbase

Definition at line 636 of file Mps.h.

◆ QinBot

template<typename Symmetry , typename Scalar = double>
vector<qarray<Nq> > Mps< Symmetry, Scalar >::QinBot

Definition at line 640 of file Mps.h.

◆ QinTop

template<typename Symmetry , typename Scalar = double>
vector<qarray<Nq> > Mps< Symmetry, Scalar >::QinTop

pre-calculated bounds for the quantum numbers that result from the finite system

Definition at line 639 of file Mps.h.

◆ qloc

template<typename Symmetry , typename Scalar = double>
vector<vector<qarray<Nq> > > Mps< Symmetry, Scalar >::qloc

local basis

Definition at line 615 of file Mps.h.

◆ Qmulti

template<typename Symmetry , typename Scalar = double>
vector<qarray<Nq> > Mps< Symmetry, Scalar >::Qmulti

multi-target quantum number for spectral functions

Definition at line 621 of file Mps.h.

◆ QoutBot

template<typename Symmetry , typename Scalar = double>
vector<qarray<Nq> > Mps< Symmetry, Scalar >::QoutBot

Definition at line 642 of file Mps.h.

◆ QoutTop

template<typename Symmetry , typename Scalar = double>
vector<qarray<Nq> > Mps< Symmetry, Scalar >::QoutTop

Definition at line 641 of file Mps.h.

◆ Qtot

template<typename Symmetry , typename Scalar = double>
qarray<Nq> Mps< Symmetry, Scalar >::Qtot = Symmetry::qvacuum()

total quantum number

Definition at line 618 of file Mps.h.

◆ S

template<typename Symmetry , typename Scalar = double>
ArrayXd Mps< Symmetry, Scalar >::S

entropy

Definition at line 630 of file Mps.h.

◆ SVspec

template<typename Symmetry , typename Scalar = double>
vector<map<qarray<Nq>,ArrayXd> > Mps< Symmetry, Scalar >::SVspec

Definition at line 632 of file Mps.h.

◆ truncWeight

template<typename Symmetry , typename Scalar = double>
ArrayXd Mps< Symmetry, Scalar >::truncWeight

truncated weight

Definition at line 627 of file Mps.h.


The documentation for this class was generated from the following files: