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

Detailed Description

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

Uniform Matrix Product State. Analogue of the Mps class.

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>

Definition at line 41 of file Umps.h.

#include <Umps.h>

Public Member Functions

 Umps ()
 
template<typename Hamiltonian >
 Umps (const Hamiltonian &H, qarray< Nq > Qtot_input, size_t L_input, size_t Mmax, size_t Nqmax, bool INIT_TO_HALF_INTEGER_SPIN)
 
 Umps (const vector< qarray< Symmetry::Nq > > &qloc_input, qarray< Nq > Qtot_input, size_t L_input, size_t Mmax, size_t Nqmax, bool INIT_TO_HALF_INTEGER_SPIN)
 
 Umps (const vector< vector< qarray< Symmetry::Nq > > > &qloc_input, qarray< Nq > Qtot_input, size_t L_input, size_t Mmax, size_t Nqmax, bool INIT_TO_HALF_INTEGER_SPIN)
 
string info () const
 
void graph (string filename) const
 
string test_ortho (double tol=1e-6) const
 
void setRandom ()
 
void normalize_C ()
 
void resize (size_t Mmax_input, size_t Nqmax_input, bool INIT_TO_HALF_INTEGER_SPIN)
 
template<typename OtherMatrixType >
void resize (const Umps< Symmetry, OtherMatrixType > &V)
 
void resize_arrays ()
 
void svdDecompose (size_t loc, GAUGE::OPTION gauge=GAUGE::C)
 
void polarDecompose (size_t loc, GAUGE::OPTION gauge=GAUGE::C)
 
VectorXd entropy () const
 
vector< map< qarray< Nq >, tuple< ArrayXd, int > > > entanglementSpectrum () const
 
std::pair< vector< qarray< Symmetry::Nq > >, ArrayXd > entanglementSpectrumLoc (size_t loc) const
 
template<typename OtherScalar >
Umps< Symmetry, OtherScalar > cast () const
 
Umps< Symmetry, double > real () const
 
vector< qarray< Symmetry::Nq > > locBasis (size_t loc) const
 
vector< vector< qarray< Symmetry::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
 
size_t get_frst_rows () const
 
size_t get_last_cols () const
 
size_t length () const
 
Scalar calc_epsLRsq (GAUGE::OPTION gauge, size_t loc) const
 
size_t calc_Mmax () const
 
size_t calc_fullMmax () const
 
size_t calc_Dmax () const
 
size_t calc_Nqmax () const
 
double memory (MEMUNIT memunit) const
 
double dot (const Umps< Symmetry, Scalar > &Vket) const
 
const vector< Biped< Symmetry, MatrixType > > & A_at (GAUGE::OPTION g, size_t loc) const
 
qarray< Symmetry::Nq > Qtop (size_t loc) const
 
qarray< Symmetry::Nq > Qbot (size_t loc) const
 
size_t minus1modL (size_t l) const
 
qarray< NqQtarget () const
 
void calc_N (DMRG::DIRECTION::OPTION DIR, size_t loc, vector< Biped< Symmetry, MatrixType > > &N) const
 
void truncate (bool SET_AC_RANDOM=true)
 
void orthogonalize_left (GAUGE::OPTION g, vector< Biped< Symmetry, MatrixType > > &G_L)
 
void orthogonalize_right (GAUGE::OPTION g, vector< Biped< Symmetry, MatrixType > > &G_R)
 
std::pair< complex< double >, Biped< Symmetry, Matrix< complex< double >, Dynamic, Dynamic > > > calc_dominant_1symm (GAUGE::OPTION g, DMRG::DIRECTION::OPTION DIR, const Mpo< Symmetry, complex< double > > &R, bool TRANSPOSE, bool CONJUGATE) const
 
std::pair< complex< double >, Biped< Symmetry, Matrix< complex< double >, Dynamic, Dynamic > > > calc_dominant_2symm (GAUGE::OPTION g, DMRG::DIRECTION::OPTION DIR, const Mpo< Symmetry, complex< double > > &R1, const Mpo< Symmetry, complex< double > > &R2) const
 
vector< std::pair< complex< double >, Biped< Symmetry, Matrix< complex< double >, Dynamic, Dynamic > > > > calc_dominant (GAUGE::OPTION g=GAUGE::R, DMRG::DIRECTION::OPTION DIR=DMRG::DIRECTION::RIGHT, int N=2, double tol=1e-15, int dimK=-1, qarray< Symmetry::Nq > Qtot=Symmetry::qvacuum(), string label="") const
 
template<typename MpoScalar >
vector< std::pair< complex< double >, Tripod< Symmetry, Matrix< complex< double >, Dynamic, Dynamic > > > > calc_dominant_Q (const Mpo< Symmetry, MpoScalar > &O, GAUGE::OPTION g=GAUGE::R, DMRG::DIRECTION::OPTION DIR=DMRG::DIRECTION::RIGHT, int N=2, double tol=1e-15, int dimK=-1, string label="") const
 
void adjustQN (const size_t number_cells)
 
void sort_A (size_t loc, GAUGE::OPTION g, bool SORT_ALL_GAUGES=false)
 
void updateC (size_t loc)
 
void updateAC (size_t loc, GAUGE::OPTION g)
 
void enrich (size_t loc, GAUGE::OPTION g, const vector< Biped< Symmetry, MatrixType > > &P)
 
template<typename MpoScalar >
ArrayXXcd intercellSF (const Mpo< Symmetry, MpoScalar > &Oalfa, const Mpo< Symmetry, MpoScalar > &Obeta, int Lx, double kmin=0., double kmax=2.*M_PI, int kpoints=51, DMRG::VERBOSITY::OPTION VERB=DMRG::VERBOSITY::ON_EXIT, double tol=1e-12)
 
template<typename MpoScalar >
complex< Scalar > intercellSFpoint (const Mpo< Symmetry, MpoScalar > &Oalfa, const Mpo< Symmetry, MpoScalar > &Obeta, int Lx, double kval, DMRG::VERBOSITY::OPTION VERB=DMRG::VERBOSITY::ON_EXIT)
 
template<typename MpoScalar >
ArrayXXcd SF (const ArrayXXcd &cellAvg, const vector< Mpo< Symmetry, MpoScalar > > &Oalfa, const vector< Mpo< Symmetry, MpoScalar > > &Obeta, int Lx, double kmin, double kmax, int kpoints, DMRG::VERBOSITY::OPTION VERB=DMRG::VERBOSITY::ON_EXIT, double tol=1e-12)
 
template<typename MpoScalar >
complex< Scalar > SFpoint (const ArrayXXcd &cellAvg, const vector< Mpo< Symmetry, MpoScalar > > &Oalfa, const vector< Mpo< Symmetry, MpoScalar > > &Obeta, int Lx, double kval, DMRG::VERBOSITY::OPTION VERB=DMRG::VERBOSITY::ON_EXIT)
 
void calc_entropy (size_t loc, bool PRINT=false)
 
void calc_entropy (bool PRINT=false)
 
void update_inbase (size_t loc, GAUGE::OPTION g=GAUGE::C)
 
void update_outbase (size_t loc, GAUGE::OPTION g=GAUGE::C)
 
void update_inbase (GAUGE::OPTION g=GAUGE::C)
 
void update_outbase (GAUGE::OPTION g=GAUGE::C)
 

Public Attributes

size_t N_sites
 
size_t Mmax
 
size_t Nqmax
 
double eps_svd = 1e-13
 
double eps_truncWeight = 1e-14
 
size_t max_Nsv =100000ul
 
size_t min_Nsv =1ul
 
int max_Nrich
 
qarray< NqQtot
 
ArrayXd truncWeight
 
vector< vector< qarray< Symmetry::Nq > > > qloc
 
std::array< vector< vector< Biped< Symmetry, MatrixType > > >, 3 > A
 
vector< Biped< Symmetry, MatrixType > > C
 
VectorXd S
 
vector< map< qarray< Nq >, tuple< ArrayXd, int > > > SVspec
 
vector< Qbasis< Symmetry > > inbase
 
vector< Qbasis< Symmetry > > outbase
 

Private Types

typedef Matrix< Scalar, Dynamic, Dynamic > MatrixType
 

Static Private Attributes

static constexpr size_t Nq = Symmetry::Nq
 

Friends

template<typename Symmetry_ , typename MpHamiltonian , typename Scalar_ >
class VumpsSolver
 
template<typename Symmetry_ , typename S1 , typename S2 >
class MpsCompressor
 

Member Typedef Documentation

◆ MatrixType

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

Definition at line 43 of file Umps.h.

Constructor & Destructor Documentation

◆ Umps() [1/4]

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

Does nothing.

Definition at line 52 of file Umps.h.

◆ Umps() [2/4]

template<typename Symmetry , typename Scalar >
template<typename Hamiltonian >
Umps< Symmetry, Scalar >::Umps ( const Hamiltonian &  H,
qarray< Nq Qtot_input,
size_t  L_input,
size_t  Mmax,
size_t  Nqmax,
bool  INIT_TO_HALF_INTEGER_SPIN 
)

Constructs a Umps with fixed bond dimension with the info from the Hamiltonian.

Definition at line 425 of file Umps.h.

◆ Umps() [3/4]

template<typename Symmetry , typename Scalar >
Umps< Symmetry, Scalar >::Umps ( const vector< qarray< Symmetry::Nq > > &  qloc_input,
qarray< Nq Qtot_input,
size_t  L_input,
size_t  Mmax,
size_t  Nqmax,
bool  INIT_TO_HALF_INTEGER_SPIN 
)

Constructs a Umps with fixed bond dimension with a uniform given basis.

Definition at line 434 of file Umps.h.

◆ Umps() [4/4]

template<typename Symmetry , typename Scalar >
Umps< Symmetry, Scalar >::Umps ( const vector< vector< qarray< Symmetry::Nq > > > &  qloc_input,
qarray< Nq Qtot_input,
size_t  L_input,
size_t  Mmax,
size_t  Nqmax,
bool  INIT_TO_HALF_INTEGER_SPIN 
)

Constructs a Umps with fixed bond dimension with a given basis.

Definition at line 445 of file Umps.h.

Member Function Documentation

◆ A_at()

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

Returns $A_L$, $A_R$ or $A_C$ at site loc as const ref.

Definition at line 201 of file Umps.h.

◆ adjustQN()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::adjustQN ( const size_t  number_cells)

This functions transforms all quantum numbers in the Umps (Umps::qloc and QN in Umps::A) by $q \rightarrow q * N_{cells}$. It is used for avg(Umps V, Mpo O, Umps V) in VumpsLinearAlgebra.h when O.length() > V.length(). In this case the quantum numbers in the Umps are transformed in correspondence with V.length() and this is incompatible with the quantum numbers in O.length() which are transformed in correspondence to O.length().

Parameters
number_cells: $N_{cells}$

Definition at line 2024 of file Umps.h.

◆ calc_Dmax()

template<typename Symmetry , typename Scalar >
size_t Umps< Symmetry, Scalar >::calc_Dmax

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

Definition at line 486 of file Umps.h.

◆ calc_dominant()

template<typename Symmetry , typename Scalar >
vector< std::pair< complex< double >, Biped< Symmetry, Matrix< complex< double >, Dynamic, Dynamic > > > > Umps< Symmetry, Scalar >::calc_dominant ( GAUGE::OPTION  g = GAUGE::R,
DMRG::DIRECTION::OPTION  DIR = DMRG::DIRECTION::RIGHT,
int  N = 2,
double  tol = 1e-15,
int  dimK = -1,
qarray< Symmetry::Nq >  Qtot = Symmetry::qvacuum(),
string  label = "" 
) const

Definition at line 2763 of file Umps.h.

◆ calc_dominant_1symm()

template<typename Symmetry , typename Scalar >
std::pair< complex< double >, Biped< Symmetry, Matrix< complex< double >, Dynamic, Dynamic > > > Umps< Symmetry, Scalar >::calc_dominant_1symm ( GAUGE::OPTION  g,
DMRG::DIRECTION::OPTION  DIR,
const Mpo< Symmetry, complex< double > > &  R,
bool  TRANSPOSE,
bool  CONJUGATE 
) const

Calculates either the right or the left fixed point of the transfer-matrix build up with A-tensors in GAUGE g.

Parameters
g: GAUGE to orthogonalize.
DIR: LEFT or RIGHT fixed point.
Note
The return values are of type complex<double>. Can we choose them sometimes to be real?

Definition at line 2933 of file Umps.h.

◆ calc_dominant_2symm()

template<typename Symmetry , typename Scalar >
std::pair< complex< double >, Biped< Symmetry, Matrix< complex< double >, Dynamic, Dynamic > > > Umps< Symmetry, Scalar >::calc_dominant_2symm ( GAUGE::OPTION  g,
DMRG::DIRECTION::OPTION  DIR,
const Mpo< Symmetry, complex< double > > &  R1,
const Mpo< Symmetry, complex< double > > &  R2 
) const

Definition at line 2985 of file Umps.h.

◆ calc_dominant_Q()

template<typename Symmetry , typename Scalar >
template<typename MpoScalar >
vector< std::pair< complex< double >, Tripod< Symmetry, Matrix< complex< double >, Dynamic, Dynamic > > > > Umps< Symmetry, Scalar >::calc_dominant_Q ( const Mpo< Symmetry, MpoScalar > &  O,
GAUGE::OPTION  g = GAUGE::R,
DMRG::DIRECTION::OPTION  DIR = DMRG::DIRECTION::RIGHT,
int  N = 2,
double  tol = 1e-15,
int  dimK = -1,
string  label = "" 
) const

Definition at line 2844 of file Umps.h.

◆ calc_entropy() [1/2]

template<typename Symmetry , typename Scalar = double>
void Umps< Symmetry, Scalar >::calc_entropy ( bool  PRINT = false)
inline

Calculate entropy for all sites.

Definition at line 359 of file Umps.h.

◆ calc_entropy() [2/2]

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::calc_entropy ( size_t  loc,
bool  PRINT = false 
)

Calculate entropy at site loc.

Definition at line 1156 of file Umps.h.

◆ calc_epsLRsq()

template<typename Symmetry , typename Scalar >
Scalar Umps< Symmetry, Scalar >::calc_epsLRsq ( GAUGE::OPTION  gauge,
size_t  loc 
) const

Calculates the left and right decomposition error as $\epsilon_L=\big|A_C-A_LC\big|^2$ and $\epsilon_R=\big|A_C-CA_R\big|^2$ (eq. (18)).

Definition at line 1223 of file Umps.h.

◆ calc_fullMmax()

template<typename Symmetry , typename Scalar >
std::size_t Umps< Symmetry, Scalar >::calc_fullMmax

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

Definition at line 512 of file Umps.h.

◆ calc_Mmax()

template<typename Symmetry , typename Scalar >
size_t Umps< Symmetry, Scalar >::calc_Mmax

Determines the maximal bond dimension per site (sum of A.rows or A.cols over all subspaces).

Definition at line 499 of file Umps.h.

◆ calc_N()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::calc_N ( DMRG::DIRECTION::OPTION  DIR,
size_t  loc,
vector< Biped< Symmetry, MatrixType > > &  N 
) const

Definition at line 1781 of file Umps.h.

◆ calc_Nqmax()

template<typename Symmetry , typename Scalar >
size_t Umps< Symmetry, Scalar >::calc_Nqmax

Determines the maximal amount of subspaces per site.

Definition at line 473 of file Umps.h.

◆ cast()

template<typename Symmetry , typename Scalar >
template<typename OtherScalar >
Umps< Symmetry, OtherScalar > Umps< Symmetry, Scalar >::cast

Casts the matrices from Scalar to OtherScalar.

Definition at line 1972 of file Umps.h.

◆ dot()

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

Calculates the scalar product with another Umps by finding the dominant eigenvalue of the transfer matrix. See arXiv:0804.2509 and Phys. Rev. B 78, 155117.

Definition at line 1121 of file Umps.h.

◆ enrich()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::enrich ( size_t  loc,
GAUGE::OPTION  g,
const vector< Biped< Symmetry, MatrixType > > &  P 
)

Enlarges the tensors of the Umps with an enrichment tensor P and resizes everything necessary with zeros. The tensor P needs to be calculated in advance. This is done directly in the VumpsSolver.

Parameters
loc: location of the site to enrich.
g: The gauge to enrich. L means, we need to update site tensor at loc+1 accordingly. R means updating site loc-1 with zeros.
P: the tensor with the enrichment. It is calculated after Eq. (A31).

Definition at line 2374 of file Umps.h.

◆ entanglementSpectrum()

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

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

Definition at line 112 of file Umps.h.

◆ entanglementSpectrumLoc()

template<typename Symmetry , typename Scalar >
std::pair< vector< qarray< Symmetry::Nq > >, ArrayXd > Umps< Symmetry, Scalar >::entanglementSpectrumLoc ( size_t  loc) const

Return the entanglement spectrum at the site loc (values all subspaces merged and sorted).

Definition at line 3060 of file Umps.h.

◆ entropy()

template<typename Symmetry , typename Scalar = double>
VectorXd Umps< Symmetry, Scalar >::entropy ( ) const
inline

Returns the entropy for all sites.

Definition at line 109 of file Umps.h.

◆ get_frst_rows()

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::get_frst_rows ( ) const
inline

Returns the amount of rows of first tensor. Useful for environment tensors in contractions.

Definition at line 141 of file Umps.h.

◆ get_last_cols()

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::get_last_cols ( ) const
inline

Returns the amount of columns of last tensor. Useful for environment tensors in contractions.

Definition at line 144 of file Umps.h.

◆ graph()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::graph ( string  filename) const

Prints a graph. See Mps::graph.

Definition at line 899 of file Umps.h.

◆ inBasis() [1/2]

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

Definition at line 134 of file Umps.h.

◆ inBasis() [2/2]

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

Returns the ingoing basis.

Definition at line 133 of file Umps.h.

◆ info()

template<typename Symmetry , typename Scalar >
string Umps< Symmetry, Scalar >::info

Prints some info about this class.

Definition at line 392 of file Umps.h.

◆ intercellSF()

template<typename Symmetry , typename Scalar >
template<typename MpoScalar >
ArrayXXcd Umps< Symmetry, Scalar >::intercellSF ( const Mpo< Symmetry, MpoScalar > &  Oalfa,
const Mpo< Symmetry, MpoScalar > &  Obeta,
int  Lx,
double  kmin = 0.,
double  kmax = 2.*M_PI,
int  kpoints = 51,
DMRG::VERBOSITY::OPTION  VERB = DMRG::VERBOSITY::ON_EXIT,
double  tol = 1e-12 
)

Calculates the static structure factor between cells according to "Tangent-space methods for uniform matrix product states" (2018), chapter 2.5.

Note
For unit cells, the Fourier transform is done between cells only, the sublattice indices are fixed by the Mpos Oalfa, Obeta.
The result has the k-values as the first column and the SSF as the second column.
Parameters
Oalfa: first operator of correlation
Obeta: second operator of correlation
Lx: length of the unit cell in x-direction
kmin: start with this k-value
kmax: end with this k-value (include it)
kpoints: number of equidistant points in interval
VERB: how much information to print

Definition at line 3083 of file Umps.h.

◆ intercellSFpoint()

template<typename Symmetry , typename Scalar >
template<typename MpoScalar >
complex< Scalar > Umps< Symmetry, Scalar >::intercellSFpoint ( const Mpo< Symmetry, MpoScalar > &  Oalfa,
const Mpo< Symmetry, MpoScalar > &  Obeta,
int  Lx,
double  kval,
DMRG::VERBOSITY::OPTION  VERB = DMRG::VERBOSITY::ON_EXIT 
)

Calculates the static structure factor between cells for one k-point only. See the more general function above.

Definition at line 3241 of file Umps.h.

◆ length()

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::length ( ) const
inline

Returns the amount of sites, i.e. the size of the unit cell.

Definition at line 147 of file Umps.h.

◆ locBasis() [1/2]

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

Definition at line 130 of file Umps.h.

◆ locBasis() [2/2]

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

Returns the local basis.

Definition at line 129 of file Umps.h.

◆ memory()

template<typename Symmetry , typename Scalar >
double Umps< Symmetry, Scalar >::memory ( MEMUNIT  memunit) const

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

Definition at line 456 of file Umps.h.

◆ minus1modL()

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::minus1modL ( size_t  l) const
inline

Safely calculates $l-1 mod L$ without overflow for size_t.

Definition at line 208 of file Umps.h.

◆ normalize_C()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::normalize_C

Normalizes the state, so that $Tr C^{\dagger} C = 1$

Definition at line 851 of file Umps.h.

◆ orthogonalize_left()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::orthogonalize_left ( GAUGE::OPTION  g,
vector< Biped< Symmetry, MatrixType > > &  G_L 
)

Orthogonalize the tensor with GAUGE g to be left-orthonormal using algorithm 2 from https://arxiv.org/abs/1810.07006.

Parameters
g: GAUGE to orthogonalize.
G_L: Gauge-Transformation which performs the orthogonalization: $A[g] \rightarrow G_L*A[g]*G_L^{-1}$
Note
: Call this function with GAUGE::L to reorthogonalize AL.

Definition at line 2639 of file Umps.h.

◆ orthogonalize_right()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::orthogonalize_right ( GAUGE::OPTION  g,
vector< Biped< Symmetry, MatrixType > > &  G_R 
)

Orthogonalize the tensor with GAUGE g to be right-orthonormal using the analogue to algorithm 2 from https://arxiv.org/abs/1810.07006.

Parameters
g: GAUGE to orthogonalize.
G_R: Gauge-Transformation which performs the orthogonalization: $A[g] \rightarrow G_R^{-1}*A[g]*G_R$
Note
: Call this function with GAUGE::R to reorthogonalize AR.

Definition at line 2574 of file Umps.h.

◆ outBasis() [1/2]

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

Definition at line 138 of file Umps.h.

◆ outBasis() [2/2]

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

Returns the outgoing basis.

Definition at line 137 of file Umps.h.

◆ polarDecompose()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::polarDecompose ( size_t  loc,
GAUGE::OPTION  gauge = GAUGE::C 
)

Calculates $A_L$ and $A_R$ from $A_C$ and $C$ at site loc using the polar decomposition (eq. (21),(22)). This is supposed to be non-optimal, but accurate.

Definition at line 1357 of file Umps.h.

◆ Qbot()

template<typename Symmetry , typename Scalar >
qarray< Symmetry::Nq > Umps< Symmetry, Scalar >::Qbot ( size_t  loc) const

Definition at line 548 of file Umps.h.

◆ Qtarget()

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

Returns the total quantum number of the Umps.

Definition at line 211 of file Umps.h.

◆ Qtop()

template<typename Symmetry , typename Scalar >
qarray< Symmetry::Nq > Umps< Symmetry, Scalar >::Qtop ( size_t  loc) const

quantum number bounds for compatibility in contract_AA

Definition at line 541 of file Umps.h.

◆ real()

template<typename Symmetry , typename Scalar >
Umps< Symmetry, double > Umps< Symmetry, Scalar >::real

Returns a real Umps containing the real part of this.

Warning
Does not check, whether the imaginary part is zero.

Definition at line 1998 of file Umps.h.

◆ resize() [1/2]

template<typename Symmetry , typename Scalar >
template<typename OtherMatrixType >
void Umps< Symmetry, Scalar >::resize ( const Umps< 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 umps

Definition at line 582 of file Umps.h.

◆ resize() [2/2]

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::resize ( size_t  Mmax_input,
size_t  Nqmax_input,
bool  INIT_TO_HALF_INTEGER_SPIN 
)

Resizes the bond dimension to Dmax and sets Nqmax blocks per site.

Definition at line 622 of file Umps.h.

◆ resize_arrays()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::resize_arrays

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

Definition at line 555 of file Umps.h.

◆ setRandom()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::setRandom

Sets all matrices $A_L$, $A_R$, $A_C$, $C$) to random using boost's uniform distribution from -1 to 1.

Definition at line 862 of file Umps.h.

◆ SF()

template<typename Symmetry , typename Scalar >
template<typename MpoScalar >
ArrayXXcd Umps< Symmetry, Scalar >::SF ( const ArrayXXcd &  cellAvg,
const vector< Mpo< Symmetry, MpoScalar > > &  Oalfa,
const vector< Mpo< Symmetry, MpoScalar > > &  Obeta,
int  Lx,
double  kmin,
double  kmax,
int  kpoints,
DMRG::VERBOSITY::OPTION  VERB = DMRG::VERBOSITY::ON_EXIT,
double  tol = 1e-12 
)

Calculates the full static structure factor for a range of k-points.

Parameters
cellAvg: all expectation values $<Oalfa_i Obeta_j>$ within unit cell
Oalfa: first operator of correlation at each cell point
Obeta: second operator of correlation at each cell point
Lx: length of the unit cell in x-direction
kmin: start with this k-value
kmax: end with this k-value (include it)
kpoints: number of equidistant points in interval
VERB: how much information to print

Definition at line 3282 of file Umps.h.

◆ SFpoint()

template<typename Symmetry , typename Scalar >
template<typename MpoScalar >
complex< Scalar > Umps< Symmetry, Scalar >::SFpoint ( const ArrayXXcd &  cellAvg,
const vector< Mpo< Symmetry, MpoScalar > > &  Oalfa,
const vector< Mpo< Symmetry, MpoScalar > > &  Obeta,
int  Lx,
double  kval,
DMRG::VERBOSITY::OPTION  VERB = DMRG::VERBOSITY::ON_EXIT 
)

Calculates the full static structure factor between cells for one k-point only. See the more general function above.

Definition at line 3250 of file Umps.h.

◆ sort_A()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::sort_A ( size_t  loc,
GAUGE::OPTION  g,
bool  SORT_ALL_GAUGES = false 
)

Sorts the A tensors of a specific gauge. If SORT_ALL_GAUGES is true, then obviously all A tensors get sorted.

Definition at line 2295 of file Umps.h.

◆ svdDecompose()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::svdDecompose ( size_t  loc,
GAUGE::OPTION  gauge = GAUGE::C 
)

Calculates $A_L$ and $A_R$ from $A_C$ and $C$ at site loc using SVD (eqs. (19),(20)). This is supposed to be optimal, but not accurate.

Definition at line 1512 of file Umps.h.

◆ test_ortho()

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

Tests the orthogonality of the Umps.

Definition at line 971 of file Umps.h.

◆ truncate()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::truncate ( bool  SET_AC_RANDOM = true)

Performs a truncation of an Umps by the singular values of the center-matrix C. Updates AL and AR with the truncated isometries from the SVD and reorthogonalize them afterwards.

Parameters
SET_AC_RANDOM: bool to decide, whether to set AC to random or to C*AR. Truncation of a converged Umps can be done with SET_AC_RANDOM=false to evaluate e.g. observables after the truncation. Truncation during the variational optimization should be done with SET_AC_RANDOM=true because otherwise this would be a bias.

Definition at line 2702 of file Umps.h.

◆ update_inbase() [1/2]

template<typename Symmetry , typename Scalar = double>
void Umps< Symmetry, Scalar >::update_inbase ( GAUGE::OPTION  g = GAUGE::C)
inline

Definition at line 387 of file Umps.h.

◆ update_inbase() [2/2]

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::update_inbase ( size_t  loc,
GAUGE::OPTION  g = GAUGE::C 
)

update basis

Definition at line 525 of file Umps.h.

◆ update_outbase() [1/2]

template<typename Symmetry , typename Scalar = double>
void Umps< Symmetry, Scalar >::update_outbase ( GAUGE::OPTION  g = GAUGE::C)
inline

Definition at line 388 of file Umps.h.

◆ update_outbase() [2/2]

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::update_outbase ( size_t  loc,
GAUGE::OPTION  g = GAUGE::C 
)

Definition at line 533 of file Umps.h.

◆ updateAC()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::updateAC ( size_t  loc,
GAUGE::OPTION  g 
)

Updates the tensor AC with zeros if the auxiallary basis has changed, e.g. after an enrichment process

Parameters
loc: location of the C tensor for the update.
g: Pull information about changed dimension from either A[GAUGE::L] or A[GAUGE::R].
Warning
Do not insert g = GAUGE::C here.

Definition at line 2345 of file Umps.h.

◆ updateC()

template<typename Symmetry , typename Scalar >
void Umps< Symmetry, Scalar >::updateC ( size_t  loc)

Updates the tensor C with zeros if the auxiallary basis has changed, e.g. after an enrichment process

Parameters
loc: location of the C tensor for the update.

Definition at line 2316 of file Umps.h.

Friends And Related Function Documentation

◆ MpsCompressor

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

Definition at line 47 of file Umps.h.

◆ VumpsSolver

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

Definition at line 46 of file Umps.h.

Member Data Documentation

◆ A

template<typename Symmetry , typename Scalar = double>
std::array<vector<vector<Biped<Symmetry,MatrixType> > >,3> Umps< Symmetry, Scalar >::A

A-tensors in the three gauges L, R, C

Definition at line 368 of file Umps.h.

◆ C

template<typename Symmetry , typename Scalar = double>
vector<Biped<Symmetry,MatrixType> > Umps< Symmetry, Scalar >::C

center matrix C

Definition at line 371 of file Umps.h.

◆ eps_svd

template<typename Symmetry , typename Scalar = double>
double Umps< Symmetry, Scalar >::eps_svd = 1e-13

Definition at line 348 of file Umps.h.

◆ eps_truncWeight

template<typename Symmetry , typename Scalar = double>
double Umps< Symmetry, Scalar >::eps_truncWeight = 1e-14

Definition at line 349 of file Umps.h.

◆ inbase

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

bases on all ingoing and outgoing legs of the Umps

Definition at line 381 of file Umps.h.

◆ max_Nrich

template<typename Symmetry , typename Scalar = double>
int Umps< Symmetry, Scalar >::max_Nrich

Definition at line 351 of file Umps.h.

◆ max_Nsv

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::max_Nsv =100000ul

Definition at line 350 of file Umps.h.

◆ min_Nsv

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::min_Nsv =1ul

Definition at line 350 of file Umps.h.

◆ Mmax

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::Mmax

Definition at line 347 of file Umps.h.

◆ N_sites

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::N_sites

parameter

Definition at line 346 of file Umps.h.

◆ Nq

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

Definition at line 44 of file Umps.h.

◆ Nqmax

template<typename Symmetry , typename Scalar = double>
size_t Umps< Symmetry, Scalar >::Nqmax

Definition at line 347 of file Umps.h.

◆ outbase

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

Definition at line 382 of file Umps.h.

◆ qloc

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

local basis

Definition at line 365 of file Umps.h.

◆ Qtot

template<typename Symmetry , typename Scalar = double>
qarray<Nq> Umps< Symmetry, Scalar >::Qtot

Definition at line 353 of file Umps.h.

◆ S

template<typename Symmetry , typename Scalar = double>
VectorXd Umps< Symmetry, Scalar >::S

null space (see eq. (25) and surrounding text)

Definition at line 376 of file Umps.h.

◆ SVspec

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

Definition at line 378 of file Umps.h.

◆ truncWeight

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

truncated weight

Definition at line 362 of file Umps.h.


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