SpinChainNNN

Inheritance Diagram

Inheritance diagram of tenpy.models.spins_nnn.SpinChainNNN

Methods

SpinChainNNN.__init__(model_params)

Initialize self.

SpinChainNNN.add_coupling(strength, u1, op1, …)

Add twosite coupling terms to the Hamiltonian, summing over lattice sites.

SpinChainNNN.add_coupling_term(strength, i, …)

Add a two-site coupling term on given MPS sites.

SpinChainNNN.add_local_term(strength, term)

Add a single term to self.

SpinChainNNN.add_onsite(strength, u, opname)

Add onsite terms to onsite_terms.

SpinChainNNN.add_onsite_term(strength, i, op)

Add an onsite term on a given MPS site.

SpinChainNNN.all_coupling_terms()

Sum of all coupling_terms.

SpinChainNNN.all_onsite_terms()

Sum of all onsite_terms.

SpinChainNNN.bond_energies(psi)

Calculate bond energies <psi|H_bond|psi>.

SpinChainNNN.calc_H_MPO([tol_zero])

Calculate MPO representation of the Hamiltonian.

SpinChainNNN.calc_H_MPO_from_bond([tol_zero])

Calculate the MPO Hamiltonian from the bond Hamiltonian.

SpinChainNNN.calc_H_bond([tol_zero])

calculate H_bond from coupling_terms and onsite_terms.

SpinChainNNN.calc_H_bond_from_MPO([tol_zero])

Calculate the bond Hamiltonian from the MPO Hamiltonian.

SpinChainNNN.calc_H_onsite([tol_zero])

Calculate H_onsite from self.onsite_terms.

SpinChainNNN.coupling_strength_add_ext_flux(…)

Add an external flux to the coupling strength.

SpinChainNNN.enlarge_mps_unit_cell([factor])

Repeat the unit cell for infinite MPS boundary conditions; in place.

SpinChainNNN.from_MPOModel(mpo_model)

Initialize a NearestNeighborModel from a model class defining an MPO.

SpinChainNNN.from_hdf5(hdf5_loader, h5gr, …)

Load instance from a HDF5 file.

SpinChainNNN.group_sites([n, grouped_sites])

Modify self in place to group sites.

SpinChainNNN.init_lattice(model_params)

Initialize a lattice for the given model parameters.

SpinChainNNN.init_sites(model_params)

Define the local Hilbert space and operators; needs to be implemented in subclasses.

SpinChainNNN.init_terms(model_params)

Add the onsite and coupling terms to the model; subclasses should implement this.

SpinChainNNN.save_hdf5(hdf5_saver, h5gr, subpath)

Export self into a HDF5 file.

SpinChainNNN.test_sanity()

Sanity check, raises ValueErrors, if something is wrong.

SpinChainNNN.trivial_like_NNModel()

Return a NearestNeighborModel with same lattice, but trivial (H=0) bonds.

class tenpy.models.spins_nnn.SpinChainNNN(model_params)[source]

Bases: tenpy.models.model.CouplingMPOModel, tenpy.models.model.NearestNeighborModel

Spin-S sites coupled by (next-)nearest neighbour interactions on a GroupedSite.

The Hamiltonian reads:

\[\begin{split}H = \sum_{\langle i,j \rangle, i < j} \mathtt{Jx} S^x_i S^x_j + \mathtt{Jy} S^y_i S^y_j + \mathtt{Jz} S^z_i S^z_j \\ + \sum_{\langle \langle i,j \rangle \rangle, i< j} \mathtt{Jxp} S^x_i S^x_j + \mathtt{Jyp} S^y_i S^y_j + \mathtt{Jzp} S^z_i S^z_j \\ - \sum_i \mathtt{hx} S^x_i + \mathtt{hy} S^y_i + \mathtt{hz} S^z_i\end{split}\]

Here, \(\langle i,j \rangle, i< j\) denotes nearest neighbors and \(\langle \langle i,j \rangle \rangle, i < j\) denotes next nearest neighbors. All parameters are collected in a single dictionary model_params, which is turned into a Config object.

Parameters

model_params (Config) – Parameters for the model. See SpinChainNNN below.

Options

config SpinChainNNN
option summary

bc_MPS

MPS boundary conditions. Coupling boundary conditions are chosen appropriately.

bc_x (from CouplingMPOModel) in FermionChain.init_lattice

``"open" | "periodic"``. [...]

bc_y (from CouplingMPOModel) in FermionChain.init_lattice

``"cylinder" | "ladder"``; only read out for 2D lattices. [...]

conserve

What should be conserved. See :class:`~tenpy.networks.Site.SpinSite`.

explicit_plus_hc (from CouplingMPOModel) in CouplingMPOModel

Whether the Hermitian conjugate of the MPO is computed at runtime, [...]

hx

Coupling as defined for the Hamiltonian above.

hy

Coupling as defined for the Hamiltonian above.

hz

Coupling as defined for the Hamiltonian above.

Jx

Coupling as defined for the Hamiltonian above.

Jxp

Coupling as defined for the Hamiltonian above.

Jy

Coupling as defined for the Hamiltonian above.

Jyp

Coupling as defined for the Hamiltonian above.

Jz

Coupling as defined for the Hamiltonian above.

Jzp

Coupling as defined for the Hamiltonian above.

L

Length of the chain in terms of :class:`~tenpy.networks.site.GroupedSite`, [...]

lattice (from CouplingMPOModel) in FermionChain.init_lattice

The name of a lattice pre-defined in TeNPy to be initialized. [...]

Lx (from CouplingMPOModel) in FermionChain.init_lattice

The length in x- and y-direction; only read out for 2D lattices. [...]

Ly (from CouplingMPOModel) in FermionChain.init_lattice

The length in x- and y-direction; only read out for 2D lattices. [...]

order (from CouplingMPOModel) in FermionChain.init_lattice

The order of sites within the lattice for non-trivial lattices, [...]

S

The 2S+1 local states range from m = -S, -S+1, ... +S.

sort_mpo_legs (from CouplingMPOModel) in CouplingMPOModel

Whether the virtual legs of the MPO should be sorted by charges, [...]

verbose (from Config) in Config

How much to print what's being done; higher means print more. [...]

option L: int

Length of the chain in terms of GroupedSite, i.e. we have 2*L spin sites.

option S: {0.5, 1, 1.5, 2, ...}

The 2S+1 local states range from m = -S, -S+1, … +S.

option conserve: 'best' | 'Sz' | 'parity' | None

What should be conserved. See SpinSite.

option Jx: float | array
option Jy: float | array
option Jz: float | array
option Jxp: float | array
option Jyp: float | array
option Jzp: float | array
option hx: float | array
option hy: float | array
option hz: float | array

Coupling as defined for the Hamiltonian above.

option bc_MPS: {'finite' | 'infinte'}

MPS boundary conditions. Coupling boundary conditions are chosen appropriately.

init_sites(model_params)[source]

Define the local Hilbert space and operators; needs to be implemented in subclasses.

This function gets called by init_lattice() to get the Site for the lattice unit cell.

Note

Initializing the sites requires to define the conserved quantum numbers. All pre-defined sites accept conserve=None to disable using quantum numbers. Many models in TeNPy read out the conserve model parameter, which can be set to "best" to indicate the optimal parameters.

Parameters

model_params (dict) – The model parameters given to __init__.

Returns

sites – The local sites of the lattice, defining the local basis states and operators.

Return type

(tuple of) Site

init_terms(model_params)[source]

Add the onsite and coupling terms to the model; subclasses should implement this.

add_coupling(strength, u1, op1, u2, op2, dx, op_string=None, str_on_first=True, raise_op2_left=False, category=None, plus_hc=False)[source]

Add twosite coupling terms to the Hamiltonian, summing over lattice sites.

Represents couplings of the form \(\sum_{x_0, ..., x_{dim-1}} strength[shift(\vec{x})] * OP0 * OP1\), where OP0 := lat.unit_cell[u0].get_op(op0) acts on the site (x_0, ..., x_{dim-1}, u1), and OP1 := lat.unit_cell[u1].get_op(op1) acts on the site (x_0+dx[0], ..., x_{dim-1}+dx[dim-1], u1). Possible combinations x_0, ..., x_{dim-1} are determined from the boundary conditions in possible_couplings().

The coupling strength may vary spatially if the given strength is a numpy array. The correct shape of this array is the coupling_shape returned by tenpy.models.lattice.possible_couplings() and depends on the boundary conditions. The shift(...) depends on dx, and is chosen such that the first entry strength[0, 0, ...] of strength is the prefactor for the first possible coupling fitting into the lattice if you imagine open boundary conditions.

The necessary terms are just added to coupling_terms; this function does not rebuild the MPO.

Deprecated since version 0.4.0: The arguments str_on_first and raise_op2_left will be removed in version 1.0.0.

Parameters
  • strength (scalar | array) – Prefactor of the coupling. May vary spatially (see above). If an array of smaller size is provided, it gets tiled to the required shape.

  • u1 (int) – Picks the site lat.unit_cell[u1] for OP1.

  • op1 (str) – Valid operator name of an onsite operator in lat.unit_cell[u1] for OP1.

  • u2 (int) – Picks the site lat.unit_cell[u2] for OP2.

  • op2 (str) – Valid operator name of an onsite operator in lat.unit_cell[u2] for OP2.

  • dx (iterable of int) – Translation vector (of the unit cell) between OP1 and OP2. For a 1D lattice, a single int is also fine.

  • op_string (str | None) – Name of an operator to be used between the OP1 and OP2 sites. Typical use case is the phase for a Jordan-Wigner transformation. The operator should be defined on all sites in the unit cell. If None, auto-determine whether a Jordan-Wigner string is needed, using op_needs_JW().

  • str_on_first (bool) – Whether the provided op_string should also act on the first site. This option should be chosen as True for Jordan-Wigner strings. When handling Jordan-Wigner strings we need to extend the op_string to also act on the ‘left’, first site (in the sense of the MPS ordering of the sites given by the lattice). In this case, there is a well-defined ordering of the operators in the physical sense (i.e. which of op1 or op2 acts first on a given state). We follow the convention that op2 acts first (in the physical sense), independent of the MPS ordering. Deprecated.

  • raise_op2_left (bool) – Raise an error when op2 appears left of op1 (in the sense of the MPS ordering given by the lattice). Deprecated.

  • category (str) – Descriptive name used as key for coupling_terms. Defaults to a string of the form "{op1}_i {op2}_j".

  • plus_hc (bool) – If True, the hermitian conjugate of the terms is added automatically.

Examples

When initializing a model, you can add a term \(J \sum_{<i,j>} S^z_i S^z_j\) on all nearest-neighbor bonds of the lattice like this:

>>> J = 1.  # the strength
>>> for u1, u2, dx in self.lat.pairs['nearest_neighbors']:
...     self.add_coupling(J, u1, 'Sz', u2, 'Sz', dx)

The strength can be an array, which gets tiled to the correct shape. For example, in a 1D Chain with an even number of sites and periodic (or infinite) boundary conditions, you can add alternating strong and weak couplings with a line like:

>>> self.add_coupling([1.5, 1.], 0, 'Sz', 0, 'Sz', dx)

Make sure to use the plus_hc argument if necessary, e.g. for hoppings:

>>> for u1, u2, dx in self.lat.pairs['nearest_neighbors']:
...     self.add_coupling(t, u1, 'Cd', u2, 'C', dx, plus_hc=True)

Alternatively, you can add the hermitian conjugate terms explictly. The correct way is to complex conjugate the strength, take the hermitian conjugate of the operators and swap the order (including a swap u1 <-> u2), and use the opposite direction -dx, i.e. the h.c. of add_coupling(t, u1, 'A', u2, 'B', dx)` is ``add_coupling(np.conj(t), u2, hc('B'), u1, hc('A'), -dx), where hc takes the hermitian conjugate of the operator names, see get_hc_op_name(). For spin-less fermions (FermionSite), this would be

>>> t = 1.  # hopping strength
>>> for u1, u2, dx in self.lat.pairs['nearest_neighbors']:
...     self.add_coupling(t, u1, 'Cd', u2, 'C', dx)
...     self.add_coupling(np.conj(t), u2, 'Cd', u1, 'C', -dx)  # h.c.

With spin-full fermions (SpinHalfFermions), it could be:

>>> for u1, u2, dx in self.lat.pairs['nearest_neighbors']:
...     self.add_coupling(t, u1, 'Cdu', u2, 'Cd', dx)  # Cdagger_up C_down
...     self.add_coupling(np.conj(t), u2, 'Cdd', u1, 'Cu', -dx)  # h.c. Cdagger_down C_up

Note that the Jordan-Wigner strings for the fermions are added automatically!

See also

add_onsite()

Add terms acting on one site only.

MultiCouplingModel.add_multi_coupling_term()

for terms on more than two sites.

add_coupling_term()

Add a single term without summing over \(vec{x}\).

add_coupling_term(strength, i, j, op_i, op_j, op_string='Id', category=None, plus_hc=False)[source]

Add a two-site coupling term on given MPS sites.

Wrapper for self.coupling_terms[category].add_coupling_term(...).

Warning

This function does not handle Jordan-Wigner strings! You might want to use add_local_term() instead.

Parameters
  • strength (float) – The strength of the coupling term.

  • j (i,) – The MPS indices of the two sites on which the operator acts. We require 0 <= i < N_sites and i < j, i.e., op_i acts “left” of op_j. If j >= N_sites, it indicates couplings between unit cells of an infinite MPS.

  • op2 (op1,) – Names of the involved operators.

  • op_string (str) – The operator to be inserted between i and j.

  • category (str) – Descriptive name used as key for coupling_terms. Defaults to a string of the form "{op1}_i {op2}_j".

  • plus_hc (bool) – If True, the hermitian conjugate of the term is added automatically.

add_local_term(strength, term, category=None, plus_hc=False)[source]

Add a single term to self.

The repesented term is strength times the product of the operators given in terms. Each operator is specified by the name and the site it acts on; the latter given by a lattice index, see Lattice.

Depending on the length of term, it can add an onsite term or a coupling term to onsite_terms or coupling_terms, respectively.

Parameters
  • strength (float/complex) – The prefactor of the term.

  • term (list of (str, array_like)) – List of tuples (opname, lat_idx) where opname is a string describing the operator acting on the site given by the lattice index lat_idx. Here, lat_idx is for example [x, y, u] for a 2D lattice, with u being the index within the unit cell.

  • category – Descriptive name used as key for onsite_terms or coupling_terms.

  • plus_hc (bool) – If True, the hermitian conjugate of the terms is added automatically.

add_onsite(strength, u, opname, category=None, plus_hc=False)[source]

Add onsite terms to onsite_terms.

Adds \(\sum_{\vec{x}} strength[\vec{x}] * OP\) to the represented Hamiltonian, where the operator OP=lat.unit_cell[u].get_op(opname) acts on the site given by a lattice index (x_0, ..., x_{dim-1}, u),

The necessary terms are just added to onsite_terms; doesn’t rebuild the MPO.

Parameters
  • strength (scalar | array) – Prefactor of the onsite term. May vary spatially. If an array of smaller size is provided, it gets tiled to the required shape.

  • u (int) – Picks a Site lat.unit_cell[u] out of the unit cell.

  • opname (str) – valid operator name of an onsite operator in lat.unit_cell[u].

  • category (str) – Descriptive name used as key for onsite_terms. Defaults to opname.

  • plus_hc (bool) – If True, the hermitian conjugate of the terms is added automatically.

See also

add_coupling()

Add a terms acting on two sites.

add_onsite_term()

Add a single term without summing over \(vec{x}\).

add_onsite_term(strength, i, op, category=None, plus_hc=False)[source]

Add an onsite term on a given MPS site.

Wrapper for self.onsite_terms[category].add_onsite_term(...).

Parameters
  • strength (float) – The strength of the term.

  • i (int) – The MPS index of the site on which the operator acts. We require 0 <= i < L.

  • op (str) – Name of the involved operator.

  • category (str) – Descriptive name used as key for onsite_terms. Defaults to op.

  • plus_hc (bool) – If True, the hermitian conjugate of the term is added automatically.

all_coupling_terms()[source]

Sum of all coupling_terms.

all_onsite_terms()[source]

Sum of all onsite_terms.

bond_energies(psi)[source]

Calculate bond energies <psi|H_bond|psi>.

Parameters

psi (MPS) – The MPS for which the bond energies should be calculated.

Returns

E_bond – List of bond energies: for finite bc, E_Bond[i] is the energy of bond i, i+1. (i.e. we omit bond 0 between sites L-1 and 0); for infinite bc E_bond[i] is the energy of bond i-1, i.

Return type

1D ndarray

calc_H_MPO(tol_zero=1e-15)[source]

Calculate MPO representation of the Hamiltonian.

Uses onsite_terms and coupling_terms to build an MPO graph (and then an MPO).

Parameters

tol_zero (float) – Prefactors with abs(strength) < tol_zero are considered to be zero.

Returns

H_MPO – MPO representation of the Hamiltonian.

Return type

MPO

calc_H_MPO_from_bond(tol_zero=1e-15)[source]

Calculate the MPO Hamiltonian from the bond Hamiltonian.

Parameters

tol_zero (float) – Arrays with norm < tol_zero are considered to be zero.

Returns

H_MPO – MPO representation of the Hamiltonian.

Return type

MPO

calc_H_bond(tol_zero=1e-15)[source]

calculate H_bond from coupling_terms and onsite_terms.

Parameters

tol_zero (float) – prefactors with abs(strength) < tol_zero are considered to be zero.

Returns

H_bond – Bond terms as required by the constructor of NearestNeighborModel. Legs are ['p0', 'p0*', 'p1', 'p1*']

Return type

list of Array

:raises ValueError : if the Hamiltonian contains longer-range terms.:

calc_H_bond_from_MPO(tol_zero=1e-15)[source]

Calculate the bond Hamiltonian from the MPO Hamiltonian.

Parameters

tol_zero (float) – Arrays with norm < tol_zero are considered to be zero.

Returns

H_bond – Bond terms as required by the constructor of NearestNeighborModel. Legs are ['p0', 'p0*', 'p1', 'p1*']

Return type

list of Array

:raises ValueError : if the Hamiltonian contains longer-range terms.:

calc_H_onsite(tol_zero=1e-15)[source]

Calculate H_onsite from self.onsite_terms.

Deprecated since version 0.4.0: This function will be removed in 1.0.0. Replace calls to this function by self.all_onsite_terms().remove_zeros(tol_zero).to_Arrays(self.lat.mps_sites()). You might also want to take explicit_plus_hc into account.

Parameters

tol_zero (float) – prefactors with abs(strength) < tol_zero are considered to be zero.

Returns

  • H_onsite (list of npc.Array)

  • onsite terms of the Hamiltonian. If explicit_plus_hc is True, – Hermitian conjugates of the onsite terms will be included.

coupling_strength_add_ext_flux(strength, dx, phase)[source]

Add an external flux to the coupling strength.

When performing DMRG on a “cylinder” geometry, it might be useful to put an “external flux” through the cylinder. This means that a particle hopping around the cylinder should pick up a phase given by the external flux [Resta1997]. This is also called “twisted boundary conditions” in literature. This function adds a complex phase to the strength array on some bonds, such that particles hopping in positive direction around the cylinder pick up exp(+i phase).

Warning

For the sign of phase it is important that you consistently use the creation operator as op1 and the annihilation operator as op2 in add_coupling().

Parameters
  • strength (scalar | array) – The strength to be used in add_coupling(), when no external flux would be present.

  • dx (iterable of int) – Translation vector (of the unit cell) between op1 and op2 in add_coupling().

  • phase (iterable of float) – The phase of the external flux for hopping in each direction of the lattice. E.g., if you want flux through the cylinder on which you have an infinite MPS, you should give phase=[0, phi] souch that particles pick up a phase phi when hopping around the cylinder.

Returns

strength – The strength array to be used as strength in add_coupling() with the given dx.

Return type

complex array

Examples

Let’s say you have an infinite MPS on a cylinder, and want to add nearest-neighbor hopping of fermions with the FermionSite. The cylinder axis is the x-direction of the lattice, so to put a flux through the cylinder, you want particles hopping around the cylinder to pick up a phase phi given by the external flux.

>>> strength = 1. # hopping strength without external flux
>>> phi = np.pi/4 # determines the external flux strength
>>> strength_with_flux = self.coupling_strength_add_ext_flux(strength, dx, [0, phi])
>>> for u1, u2, dx in self.lat.pairs['nearest_neighbors']:
...     self.add_coupling(strength_with_flux, u1, 'Cd', u2, 'C', dx)
...     self.add_coupling(np.conj(strength_with_flux), u2, 'Cd', u1, 'C', -dx)
enlarge_mps_unit_cell(factor=2)[source]

Repeat the unit cell for infinite MPS boundary conditions; in place.

This has to be done after finishing initialization and can not be reverted.

Parameters

factor (int) – The new number of sites in the MPS unit cell will be increased from N_sites to factor*N_sites_per_ring. Since MPS unit cells are repeated in the x-direction in our convetion, the lattice shape goes from (Lx, Ly, ..., Lu) to (Lx*factor, Ly, ..., Lu).

classmethod from_MPOModel(mpo_model)[source]

Initialize a NearestNeighborModel from a model class defining an MPO.

This is especially usefull in combination with MPOModel.group_sites().

Parameters

mpo_model (MPOModel) – A model instance implementing the MPO. Does not need to be a NearestNeighborModel, but should only have nearest-neighbor couplings.

Examples

The SpinChainNNN2 has next-nearest-neighbor couplings and thus only implements an MPO:

>>> from tenpy.models.spins_nnn import SpinChainNNN2
>>> nnn_chain = SpinChainNNN2({'L': 20})
parameter 'L'=20 for SpinChainNNN2
>>> print(isinstance(nnn_chain, NearestNeighborModel))
False
>>> print("range before grouping:", nnn_chain.H_MPO.max_range)
range before grouping: 2

By grouping each two neighboring sites, we can bring it down to nearest neighbors.

>>> nnn_chain.group_sites(2)
>>> print("range after grouping:", nnn_chain.H_MPO.max_range)
range after grouping: 1

Yet, TEBD will not yet work, as the model doesn’t define H_bond. However, we can initialize a NearestNeighborModel from the MPO:

>>> nnn_chain_for_tebd = NearestNeighborModel.from_MPOModel(nnn_chain)
classmethod from_hdf5(hdf5_loader, h5gr, subpath)[source]

Load instance from a HDF5 file.

This method reconstructs a class instance from the data saved with save_hdf5().

Parameters
  • hdf5_loader (Hdf5Loader) – Instance of the loading engine.

  • h5gr (Group) – HDF5 group which is represent the object to be constructed.

  • subpath (str) – The name of h5gr with a '/' in the end.

Returns

obj – Newly generated class instance containing the required data.

Return type

cls

group_sites(n=2, grouped_sites=None)[source]

Modify self in place to group sites.

Group each n sites together using the GroupedSite. This might allow to do TEBD with a Trotter decomposition, or help the convergence of DMRG (in case of too long range interactions).

This has to be done after finishing initialization and can not be reverted.

Parameters
  • n (int) – Number of sites to be grouped together.

  • grouped_sites (None | list of GroupedSite) – The sites grouped together.

Returns

grouped_sites – The sites grouped together.

Return type

list of GroupedSite

init_lattice(model_params)[source]

Initialize a lattice for the given model parameters.

This function reads out the model parameter lattice. This can be a full Lattice instance, in which case it is just returned without further action. Alternatively, the lattice parameter can be a string giving the name of one of the predefined lattices, which then gets initialized. Depending on the dimensionality of the lattice, this requires different model parameters.

Parameters

model_params (dict) – The model parameters given to __init__.

Returns

lat – An initialized lattice.

Return type

Lattice

Options

option CouplingMPOModel.lattice: str | Lattice

The name of a lattice pre-defined in TeNPy to be initialized. Alternatively, a (possibly self-defined) Lattice instance. In the latter case, no further parameters are read out.

option CouplingMPOModel.bc_MPS: str

Boundary conditions for the MPS.

option CouplingMPOModel.order: str

The order of sites within the lattice for non-trivial lattices, e.g, 'default', 'snake', see ordering(). Only used if lattice is a string.

option CouplingMPOModel.L: int

The length in x-direction; only read out for 1D lattices. For an infinite system the length of the unit cell.

option CouplingMPOModel.Lx: int
option CouplingMPOModel.Ly: int

The length in x- and y-direction; only read out for 2D lattices. For "infinite" bc_MPS, the system is infinite in x-direction and Lx is the number of “rings” in the infinite MPS unit cell, while Ly gives the circumference around the cylinder or width of th the rung for a ladder (depending on bc_y).

option CouplingMPOModel.bc_y: str

"cylinder" | "ladder"; only read out for 2D lattices. The boundary conditions in y-direction.

option CouplingMPOModel.bc_x: str

"open" | "periodic". Can be used to force “periodic” boundaries for the lattice, i.e., for the couplings in the Hamiltonian, even if the MPS is finite. Defaults to "open" for bc_MPS="finite" and "periodic" for bc_MPS="infinite. If you are not aware of the consequences, you should probably not use “periodic” boundary conditions. (The MPS is still “open”, so this will introduce long-range couplings between the first and last sites of the MPS!)

save_hdf5(hdf5_saver, h5gr, subpath)[source]

Export self into a HDF5 file.

This method saves all the data it needs to reconstruct self with from_hdf5().

This implementation saves the content of __dict__ with save_dict_content(), storing the format under the attribute 'format'.

Parameters
  • hdf5_saver (Hdf5Saver) – Instance of the saving engine.

  • h5gr (:class`Group`) – HDF5 group which is supposed to represent self.

  • subpath (str) – The name of h5gr with a '/' in the end.

test_sanity()[source]

Sanity check, raises ValueErrors, if something is wrong.

trivial_like_NNModel()[source]

Return a NearestNeighborModel with same lattice, but trivial (H=0) bonds.