Generic abelian-symmetric iPEPS¶
Class describes abelian-symmetric iPEPS with rectangular unit cell, together with convenience functions for reading and writing such iPEPS from and into (JSON) files.
- class ipeps.ipeps_abelian.IPEPS_ABELIAN(settings, sites, vertexToSite=None, lX=None, lY=None, peps_args=<config.PEPSARGS object>, global_args=<config.GLOBALARGS object>)[source]¶
- Parameters:
settings (NamedTuple or SimpleNamespace (TODO link to definition)) – YAST configuration
sites (dict[tuple(int,int) : yast.Tensor]) – map from elementary unit cell to on-site tensors
vertexToSite (function(tuple(int,int))->tuple(int,int)) – function mapping arbitrary vertex of a square lattice into a vertex within elementary unit cell
lX (int) – length of the elementary unit cell in X direction
lY (int) – length of the elementary unit cell in Y direction
build_open_dl (bool) – build complementary
IPEPS_ABELIAN
with with open double-layer on-site tensorspeps_args (PEPSARGS) – ipeps configuration
global_args (GLOBALARGS) – global configuration
Member
sites
is a dictionary of non-equivalent on-site tensors indexed by tuple of coordinates (x,y) within the elementary unit cell. The index-position convetion for on-site tensors is defined as follows:(-1)u (-1)s |/ (-1)l--a--(+1)r <=> a[s,u,l,d,r] with reference symmetry signature [-1,-1,-1,1,1] | (+1)d
where s denotes physical index, and u,l,d,r label four principal directions up, left, down, right in anti-clockwise order starting from up. Member
vertexToSite
is a mapping function from any vertex (x,y) on a square lattice passed in as tuple(int,int) to a corresponding vertex within elementary unit cell.On-site tensor of an IPEPS_ABELIAN object
wfc
at vertex (x,y) is conveniently accessed through the member functionsite
, which internally usesvertexToSite
mapping:coord= (0,0) a_00= wfc.site(coord)
By combining the appropriate
vertexToSite
mapping function with elementary unit cell specified throughsites
, various tilings of a square lattice can be achieved:# Example 1: 1-site translational iPEPS sites={(0,0): a} def vertexToSite(coord): return (0,0) wfc= IPEPS_ABELIAN(sites,vertexToSite) # resulting tiling: # y\x -2 -1 0 1 2 # -2 a a a a a # -1 a a a a a # 0 a a a a a # 1 a a a a a # Example 2: 2-site bipartite iPEPS sites={(0,0): a, (1,0): b} def vertexToSite(coord): x = (coord[0] + abs(coord[0]) * 2) % 2 y = abs(coord[1]) return ((x + y) % 2, 0) wfc= IPEPS_ABELIAN(sites,vertexToSite) # resulting tiling: # y\x -2 -1 0 1 2 # -2 A b a b a # -1 B a b a b # 0 A b a b a # 1 B a b a b # Example 3: iPEPS with 3x2 unit cell with PBC sites={(0,0): a, (1,0): b, (2,0): c, (0,1): d, (1,1): e, (2,1): f} wfc= IPEPS_ABELIAN(sites,lX=3,lY=2) # resulting tiling: # y\x -2 -1 0 1 2 # -2 b c a b c # -1 e f d e f # 0 b c a b c # 1 e f d e f
where in the last example a default setting for
vertexToSite
is used, which maps square lattice into elementary unit cell of sizelX
xlY
assuming periodic boundary conditions (PBC) along both X and Y directions.Performance-wise, it is favourable to construct complementary iPEPS formed by open double-layer on-site tensors. Such tensors are repeatedly used as building blocks of reduced density matrices when computing observables. If build_open_dl=True they are pre-computed and accessible through member sites_dl_open or site_dl_open(coord) convenience function.
Note
in case of differentiation through reduced density matrix construction, the sites_dl_open computation must be a member of computation graph for correct gradients. It can be explicitly recomputed by invoking build_sites_dl_open().
- add_noise(noise=0, peps_args=<config.PEPSARGS object>)[source]¶
- Parameters:
noise (float) – magnitude of the noise
- Returns:
a copy of state with noisy on-site tensors. For default value of
noise
being zeroself
is returned.- Return type:
Create a new state by adding random uniform noise with magnitude
noise
to all copies of on-site tensors. The noise is added to all blocks making up the individual on-site tensors.
- build_sites_dl()[source]¶
If
config.PEPSARGS.build_dl_open
, build complementary on-site double-layer iPEPS.
- build_sites_dl_open(fusion_level='full')[source]¶
If
config.PEPSARGS.build_dl
, build complementary open on-site double-layer iPEPS.- Parameters:
fusion_level (str) – see _fused_open_dl_site()
- get_checkpoint()[source]¶
- Returns:
serializable representation of IPEPS_ABELIAN state
- Return type:
dict
Return dict containing serialized on-site (block-sparse) tensors. The individual blocks are serialized into Numpy ndarrays. This function is called by optimizer to create checkpoints during the optimization process.
- get_parameters()[source]¶
- Returns:
variational parameters of iPEPS
- Return type:
iterable
This function is called by optimizer to access variational parameters of the state.
- load_checkpoint(checkpoint_file)[source]¶
- Parameters:
checkpoint_file (str or file object) – path to checkpoint file
Initializes IPEPS_ABELIAN from checkpoint file.
Note
The vertexToSite mapping function is not a part of checkpoint and must be provided either when instantiating IPEPS_ABELIAN or afterwards.
- site(coord)[source]¶
- Parameters:
coord (tuple(int,int)) – tuple (x,y) specifying vertex on a square lattice
- Returns:
on-site tensor corresponding to the vertex (x,y)
- Return type:
yast.Tensor
- site_dl(coord)[source]¶
- Parameters:
coord (tuple(int,int)) – tuple (x,y) specifying vertex on a square lattice
- Returns:
double-layer on-site tensor corresponding to the vertex (x,y)
- Return type:
yast.Tensor
If
config.PEPSARGS.build_dl
, then precomputed double-layer on-site tensor is returned. Otherwise, Exception is raised.
- site_dl_open(coord)[source]¶
- Parameters:
coord (tuple(int,int)) – tuple (x,y) specifying vertex on a square lattice
- Returns:
open double-layer on-site tensor corresponding to the vertex (x,y)
- Return type:
yast.Tensor
If
config.PEPSARGS.build_dl_open
, then precomputed open double-layer on-site tensor is returned. Otherwise, Exception is raised.
- sync_precomputed()[source]¶
Force recomputation of double-layer and open double-layer on-site tensors if corresponding options
config.PEPSARGS.build_dl
andconfig.PEPSARGS.build_dl_open
areTrue
.Note
In active autograd regions, it might be necessary to force recomputation if the corresponding double-layer tensors are to be part of computational graph and hence differentiated.
- to(device)[source]¶
- Parameters:
device (str) – device identifier
- Returns:
returns a copy of the state on
device
. If the state already resides on device returnsself
.- Return type:
Move the entire state to
device
.
- to_dense(peps_args=<config.PEPSARGS object>, global_args=<config.GLOBALARGS object>)[source]¶
- Returns:
returns equivalent dense state with all on-site tensors in their dense representation on PyTorch backend.
- Return type:
Create an IPEPS state with all on-site tensors as dense possesing no explicit block structure (symmetry). This operations preserves gradients on returned dense state.
- write_to_file(outputfile, tol=None, normalize=False)[source]¶
Writes state to file. See
write_ipeps()
.
- ipeps.ipeps_abelian.read_ipeps(jsonfile, settings, vertexToSite=None, peps_args=<config.PEPSARGS object>, global_args=<config.GLOBALARGS object>)[source]¶
- Parameters:
jsonfile (str or Path object) – input file describing IPEPS_ABELIAN in json format
settings (NamedTuple or SimpleNamespace (TODO link to definition)) – YAST configuration
vertexToSite (function(tuple(int,int))->tuple(int,int)) – function mapping arbitrary vertex of a square lattice into a vertex within elementary unit cell
peps_args (PEPSARGS) – ipeps configuration
global_args (GLOBALARGS) – global configuration
- Returns:
wavefunction
- Return type:
Read state in JSON format from file.
A simple PBC
vertexToSite
function is used by default
- ipeps.ipeps_abelian.write_ipeps(state, outputfile, tol=None, normalize=False, peps_args=<config.PEPSARGS object>, global_args=<config.GLOBALARGS object>)[source]¶
- Parameters:
state (IPEPS_ABELIAN) – wavefunction to write out in json format
outputfile – target file
tol (float) – minimum magnitude of tensor elements which are written out
normalize (bool) – if True, on-site tensors are normalized before writing
Write state to file.
Generic abelian-symmetric iPEPS with weights¶
- class ipeps.ipeps_abelian.IPEPS_ABELIAN_WEIGHTED(settings, sites, weights, vertexToSite=None, lX=None, lY=None, peps_args=<config.PEPSARGS object>, global_args=<config.GLOBALARGS object>)[source]¶
- Parameters:
settings (NamedTuple or SimpleNamespace (TODO link to definition)) – YAST configuration
sites (dict[tuple(int,int) : yast.Tensor]) – map from elementary unit cell to on-site tensors
weights (dict[tuple(tuple(int,int), tuple(int,int)) : yast.Tensor]) – map from edges within unit cell to weight tensors
vertexToSite (function(tuple(int,int))->tuple(int,int)) – function mapping arbitrary vertex of a square lattice into a vertex within elementary unit cell
lX (int) – length of the elementary unit cell in X direction
lY (int) – length of the elementary unit cell in Y direction
peps_args (PEPSARGS) – ipeps configuration
global_args (GLOBALARGS) – global configuration
IPEPS_ABELIAN_WEIGHTED augments basic IPEPS_ABELIAN with a tensor on each bond within elementary unit cell. In case of diagonal and positive semi-definite tensors, these are called weights. Such augmented ansatz provides basic structure for iTEBD algorithms such as Simple Update.
The keys of weights dictionary index tensors by tuple of (coord, dxy) where coord specifies site within elementary unit cell and (dxy) is a directional vector specifying up, left, down, or right bond of that site as (0,-1), (-1,0), (0,1) or (1,0) respectively. Thus the weights is not injective dictionary, instead keys (coord,dxy) and (coord+dxy,-dxy) should index identical tensor.
- absorb_weights(peps_args=<config.PEPSARGS object>, global_args=<config.GLOBALARGS object>)[source]¶
- Parameters:
build_open_dl (bool) – see IPEPS_ABELIAN
- Returns:
regular IPEPS_ABELIAN obtained by symmetricaly absorbing weights of IPEPS_ABELIAN_WEIGHTED into its on-site tensors
- Return type:
Reduce weighted iPEPS to regular iPEPS by splitting its weights symmetrically as W = sqrt(W)sqrt(W) and absorbing them into on-site tensors:
\sqrt(W) |/s |/s \sqrt(W)--a--\sqrt(W) = --a'-- | | \sqrt(W)
Note
assumes weight tensors are diagonal and positive semi-definite
- weight(weight_id)[source]¶
- Parameters:
weight_id (tuple(tuple(int,int), tuple(int,int))) – tuple with (x,y) coords specifying vertex on a square lattice and tuple with (dx,dy) coords specifying on of the directions (0,-1), (-1,0), (0,1), (1,0) corresponding to up, left, down, and right respectively.
- Returns:
diagonal weight tensor
- Return type:
yast.Tensor
- ipeps.ipeps_abelian.get_weighted_ipeps(state, weights, peps_args=<config.PEPSARGS object>, global_args=<config.GLOBALARGS object>)[source]
- Parameters:
state (IPEPS_ABELIAN) – abelian-symmetric iPEPS
weights (dict[tuple(tuple(int,int), tuple(int,int)) : yamps.tensor.Tensor]) – map from edges within unit cell to weight tensors
peps_args (PEPSARGS) – ipeps configuration
global_args (GLOBALARGS) – global configuration
- Returns:
iPEPS wavefunction augmented with weights
- Return type:
Create IPEPS_ABELIAN_WEIGHTED from regular IPEPS_ABELIAN by mapping weight tensor to each bond in the elementary unit cell.