GraceQ/MPS2
A high-performance matrix product state algorithms library based on GraceQ/tensor
Functions
two_site_update_finite_vmps_impl.h File Reference

Implementation details for two-site finite variational MPS algorithm. More...

#include "gqmps2/algorithm/vmps/two_site_update_finite_vmps.h"
#include "gqmps2/one_dim_tn/mpo/mpo.h"
#include "gqmps2/one_dim_tn/mps/finite_mps/finite_mps.h"
#include "gqmps2/utilities.h"
#include "gqmps2/one_dim_tn/framework/ten_vec.h"
#include "gqmps2/consts.h"
#include "gqten/gqten.h"
#include "gqten/utility/timer.h"
#include <iostream>
#include <iomanip>
#include <vector>
#include <string>
#include <stdio.h>
#include <assert.h>

Go to the source code of this file.

Functions

template<typename DTenT >
double gqmps2::MeasureEE (const DTenT &s, const size_t sdim)
 
std::string gqmps2::GenEnvTenName (const std::string &dir, const long blk_len, const std::string temp_path)
 
void gqmps2::RemoveFile (const std::string &file)
 
template<typename TenElemT , typename QNT >
GQTEN_Double gqmps2::TwoSiteFiniteVMPS (FiniteMPS< TenElemT, QNT > &mps, const MPO< GQTensor< TenElemT, QNT >> &mpo, const SweepParams &sweep_params)
 Function to perform two-site update finite vMPS algorithm. More...
 
template<typename TenElemT , typename QNT >
void gqmps2::InitEnvs (FiniteMPS< TenElemT, QNT > &mps, const MPO< GQTensor< TenElemT, QNT >> &mpo, const SweepParams &sweep_params)
 
template<typename TenElemT , typename QNT >
double gqmps2::TwoSiteFiniteVMPSSweep (FiniteMPS< TenElemT, QNT > &mps, const MPO< GQTensor< TenElemT, QNT >> &mpo, const SweepParams &sweep_params)
 Function to perform a single two-site finite vMPS sweep. More...
 
template<typename TenElemT , typename QNT >
double gqmps2::TwoSiteFiniteVMPSUpdate (FiniteMPS< TenElemT, QNT > &mps, TenVec< GQTensor< TenElemT, QNT >> &lenvs, TenVec< GQTensor< TenElemT, QNT >> &renvs, const MPO< GQTensor< TenElemT, QNT >> &mpo, const SweepParams &sweep_params, const char dir, const size_t target_site)
 
template<typename TenElemT , typename QNT >
void gqmps2::LoadRelatedTens (FiniteMPS< TenElemT, QNT > &mps, TenVec< GQTensor< TenElemT, QNT >> &lenvs, TenVec< GQTensor< TenElemT, QNT >> &renvs, const size_t target_site, const char dir, const SweepParams &sweep_params)
 
template<typename TenElemT , typename QNT >
void gqmps2::DumpRelatedTens (FiniteMPS< TenElemT, QNT > &mps, TenVec< GQTensor< TenElemT, QNT >> &lenvs, TenVec< GQTensor< TenElemT, QNT >> &renvs, const size_t target_site, const char dir, const SweepParams &sweep_params)
 

Detailed Description

Implementation details for two-site finite variational MPS algorithm.

Function Documentation

◆ TwoSiteFiniteVMPS()

template<typename TenElemT , typename QNT >
GQTEN_Double gqmps2::TwoSiteFiniteVMPS ( FiniteMPS< TenElemT, QNT > &  mps,
const MPO< GQTensor< TenElemT, QNT >> &  mpo,
const SweepParams sweep_params 
)

Function to perform two-site update finite vMPS algorithm.

Note
The input MPS will be considered an empty one.

◆ TwoSiteFiniteVMPSSweep()

template<typename TenElemT , typename QNT >
double gqmps2::TwoSiteFiniteVMPSSweep ( FiniteMPS< TenElemT, QNT > &  mps,
const MPO< GQTensor< TenElemT, QNT >> &  mpo,
const SweepParams sweep_params 
)

Function to perform a single two-site finite vMPS sweep.

Note
Before the sweep and after the sweep, the MPS is empty.