quspin.basis.tensor_basis

class quspin.basis.tensor_basis(*basis_list)[source]

Bases: basis

Constructs basis in tensor product Hilbert space.

The tensor_basis class combines two basis objects basis1 and basis2 together into a new basis object which can be then used, e.g., to create the Hamiltonian over the tensor product Hilbert space:

\[\mathcal{H}=\mathcal{H}_1\otimes\mathcal{H}_2\]

Notes

  • The tensor_basis operator strings are separated by a pipe symbol, “|”. However, the index array has NO pipe symbol.

  • If two fermion basis constructors are used, tensor_basis will assume that the two fermion species are distinguishable, i.e. their operators will commute (instead of anti-commute). For anticommuting fermion species, use the spinful_fermion_basis_* constructors.

The tensor_basis class does not allow one to make use of symmetries, save for particle conservation.

Examples

The following code shows how to construct the Fermi-Hubbard Hamiltonian by tensoring two spinless_fermion_basis_1d objects. This model can also be set up using the spinful_fermion_basis_1d class), which also allows the implementation of symmetries.

Notice that the operator strings for constructing Hamiltonians with a tensor_basis object are separated by a pipe symbol, ‘|’, while the index array has no splitting pipe character.

The code snippet below initiates the class, and is required to run the example codes for the function methods.

 1from quspin.basis import tensor_basis, spinless_fermion_basis_1d  # Hilbert spaces
 2from quspin.tools.measurements import obs_vs_time  # calculating dynamics
 3import numpy as np  # general math functions
 4
 5#
 6##### setting parameters for simulation
 7# physical parameters
 8L = 8  # system size
 9N = L // 2  # number of particles
10N_up = N // 2 + N % 2  # number of fermions with spin up
11N_down = N // 2  # number of fermions with spin down
12J = 1.0  # hopping strength
13U = 5.0  # interaction strength
14#
15###### create the basis
16# build the two bases to tensor together to spinful fermions
17basis_up = spinless_fermion_basis_1d(L, Nf=N_up)  # up basis
18basis_down = spinless_fermion_basis_1d(L, Nf=N_down)  # down basis
19basis = tensor_basis(basis_up, basis_down)  # spinful fermions
20print(basis)
21#
22##### create model
23# define site-coupling lists
24hop_right = [[-J, i, i + 1] for i in range(L - 1)]  # hopping to the right OBC
25hop_left = [[J, i, i + 1] for i in range(L - 1)]  # hopping to the left OBC
26int_list = [[U, i, i] for i in range(L)]  # onsite interaction
27# create static lists
28static = [
29    ["+-|", hop_left],  # up hop left
30    ["-+|", hop_right],  # up hop right
31    ["|+-", hop_left],  # down hop left
32    ["|-+", hop_right],  # down hop right
33    ["n|n", int_list],  # onsite interaction
34]
35#
36###### setting up operators
37# set up hamiltonian dictionary and observable (imbalance I)
38no_checks = dict(check_pcon=False, check_symm=False, check_herm=False)
39H = hamiltonian(static, [], basis=basis, **no_checks)
__init__(*basis_list)[source]

Initialises the tensor_basis object (basis for tensor product Hilbert spaces).

Parameters:
basis_listlist[basis]

List of basis objects to tensor together. Required minimum number is two.

Methods

Op(opstr, indx, J, dtype)

Constructs operator from a site-coupling list and an operator string in the tensor basis.

__init__(*basis_list)

Initialises the tensor_basis object (basis for tensor product Hilbert spaces).

check_hermitian(static, dynamic)

Checks operator string lists for hermiticity of the combined operator.

check_pcon(static, dynamic)

Checks operator string lists for particle number (magnetisation) conservartion of the combined operator.

check_symm(static, dynamic)

Checks operator string lists for the required symmetries of the combined operator.

ent_entropy(state[, sub_sys_A, return_rdm, ...])

Calculates entanglement entropy of subsystem A and the corresponding reduced density matrix

expanded_form([static, dynamic])

Splits up operator strings containing "x" and "y" into operator combinations of "+" and "-".

get_proj(dtype[, full_left, full_right])

Calculates transformation/projector from symmetry-reduced basis to full (symmetry-free) basis.

get_vec(v0[, sparse, full_left, full_right])

DEPRECATED (cf project_from).

index(*states)

Finds the index of user-defined Fock state in tensor basis.

inplace_Op(v_in, op_list, dtype[, ...])

Calculates the action of an operator on a state.

partial_trace(state[, sub_sys_A, ...])

Calculates reduced density matrix, through a partial trace of a quantum state in tensor_basis.

project_from(v0[, sparse, full_left, full_right])

Transforms state from symmetry-reduced basis to full (symmetry-free) basis.

Attributes

N

the value of N attribute from all the basis objects tensored together in a tuple ordered according to the input basis list.

Ns

number of states in the Hilbert space.

basis_left

first basis constructor out of the basis objects list to be tensored.

basis_right

all others basis constructors except for the first one of the basis objects list to be tensored.

blocks

contains the quantum numbers (blocks) for the symmetry sectors.

operators

set of available operator strings.

sps

number of states per site (ie, the on-site Hilbert space dimension).

property N

the value of N attribute from all the basis objects tensored together in a tuple ordered according to the input basis list.

Type:

tuple

property Ns

number of states in the Hilbert space.

Type:

int

Op(opstr, indx, J, dtype)[source]

Constructs operator from a site-coupling list and an operator string in the tensor basis.

Parameters:
opstrstr

Operator string in the tensor basis format. For instance: >>> opstr = “z|z”

indxlist(int)

List of integers to designate the sites the tensor basis operator is defined on. For instance: >>> indx = [1,5]

Jscalar

Coupling strength.

dtype‘type’

Data type (e.g. numpy.float64) to construct the operator with.

Returns:
tuple
(ME,row,col), where
  • numpy.ndarray(scalar): ME: matrix elements of type dtype.

  • numpy.ndarray(int): row: row indices of matrix representing the operator in the tensor basis,

    such that row[i] is the row index of ME[i].

  • numpy.ndarray(int): col: column index of matrix representing the operator in the tensor basis,

    such that col[i] is the column index of ME[i].

Examples

>>> J = 1.41
>>> indx = [1,5]
>>> opstr = "z|z"
>>> dtype = np.float64
>>> ME, row, col = Op(opstr,indx,J,dtype)
property basis_left

first basis constructor out of the basis objects list to be tensored.

Type:

basis

property basis_right

all others basis constructors except for the first one of the basis objects list to be tensored.

Type:

basis

property blocks

contains the quantum numbers (blocks) for the symmetry sectors.

Type:

dict

check_hermitian(static, dynamic)

Checks operator string lists for hermiticity of the combined operator.

Parameters:
static: list

Static operators formatted to be passed into the static argument of the hamiltonian class.

dynamic: list

Dynamic operators formatted to be passed into the dynamic argument of the hamiltonian class.

check_pcon(static, dynamic)

Checks operator string lists for particle number (magnetisation) conservartion of the combined operator.

Parameters:
static: list

Static operators formatted to be passed into the static argument of the hamiltonian class.

dynamic: list

Dynamic operators formatted to be passed into the dynamic argument of the hamiltonian class.

check_symm(static, dynamic)

Checks operator string lists for the required symmetries of the combined operator.

Parameters:
static: list

Static operators formatted to be passed into the static argument of the hamiltonian class.

dynamic: list

Dynamic operators formatted to be passed into the dynamic argument of the hamiltonian class.

ent_entropy(state, sub_sys_A='left', return_rdm=None, enforce_pure=False, return_rdm_EVs=False, sparse=False, alpha=1.0, sparse_diag=True, maxiter=None)[source]

Calculates entanglement entropy of subsystem A and the corresponding reduced density matrix

\[S_\mathrm{ent}(\alpha) = \frac{1}{1-\alpha}\log \mathrm{tr}_{A} \left( \mathrm{tr}_{A^c} \vert\psi\rangle\langle\psi\vert \right)^\alpha\]

Note: The logarithm used is the natural logarithm (base e).

Parameters:
stateobj

State of the quantum system. Can be either one of:

  • numpy.ndarray [shape (Ns,)]: pure state (default).

  • numpy.ndarray [shape (Ns,Ns)]: density matrix (DM).

sub_sys_Astr, optional

Defines subsystem A. Can be either one of:

  • “left”: refers to basis_left (Default).

  • “right”: refers to basis_right.

  • “both”: for initial mixed states the Renyi entropy of subsystem A and its complement

    B need not be the same. This option automatically sets return_rdm=both.

return_rdmstr, optional

Toggles returning the reduced DM. Can be tierh one of:

  • “A”: returns reduced DM of subsystem A.

  • “B”: returns reduced DM of subsystem B (complement of A).

  • “both”: returns reduced DM of both subsystems A and B.

enforce_purebool, optional

Whether or not to assume state is a collection of pure states or a mixed density matrix, if it is a square array. Default is False.

sparsebool, optional

Whether or not to return a sparse DM. Default is False.

return_rdm_EVsbool, optional

Whether or not to return the eigenvalues of rthe educed DM. If return_rdm is specified, the eigenvalues of the corresponding DM are returned. If return_rdm is NOT specified, the spectrum of rdm_A is returned by default. Default is False.

alphafloat, optional

Renyi \(\alpha\) parameter for the entanglement entropy. Default is \(\alpha=1\).

sparse_diagbool, optional

When sparse=True, this flag enforces the use of scipy.sparse.linalg.eigsh() to calculate the eigenvaues of the reduced DM.

maxiterint, optional

Specifies the number of iterations for Lanczos diagonalisation. Look up documentation for scipy.sparse.linalg.eigsh().

Returns:
dict
Dictionary with following keys, depending on input parameters:
  • “Sent_A”: entanglement entropy of subsystem A (default).

  • “Sent_B”: entanglement entropy of subsystem B.

  • “p_A”: singular values of reduced DM of subsystem A (default).

  • “p_B”: singular values of reduced DM of subsystem B.

  • “rdm_A”: reduced DM of subsystem A.

  • “rdm_B”: reduced DM of subsystem B.

Notes

Algorithm is based on both partial tracing and sigular value decomposition (SVD), optimised for speed.

Examples

>>> ent_entropy(state,sub_sys_A="left",return_rdm="A",enforce_pure=False,return_rdm_EVs=False,
>>>                             sparse=False,alpha=1.0,sparse_diag=True)
expanded_form(static=[], dynamic=[])

Splits up operator strings containing “x” and “y” into operator combinations of “+” and “-”. This function is useful for higher spin hamiltonians where “x” and “y” operators are not appropriate operators.

Parameters:
static: list

Static operators formatted to be passed into the static argument of the hamiltonian class.

dynamic: list

Dynamic operators formatted to be passed into the dynamic argument of the hamiltonian class.

Returns:
tuple
(static, dynamic), where
  • list: static: operator strings with “x” and “y” expanded into “+” and “-”, formatted to

    be passed into the static argument of the hamiltonian class.

  • list: dynamic: operator strings with “x” and “y” expanded into “+” and “-”, formatted to

    be passed into the dynamic argument of the hamiltonian class.

Notes

This function works with the tensor_basis and other basis which use the “|” symbol in the opstr.

Examples

>>> static = [["xx",[[1.0,0,1]]],["yy",[[1.0,0,1]]]]
>>> dynamic = [["y",[[1.0,0]],lambda t: t,[]]]
>>> expanded_form(static,dynamic)
get_proj(dtype, full_left=True, full_right=True)[source]

Calculates transformation/projector from symmetry-reduced basis to full (symmetry-free) basis.

Parameters:
dtype‘type’

Data type (e.g. numpy.float64) to construct the projector with.

full_leftbool, optional

Whether or not to transform the state to the full state in basis_left. Default is True.

full_rightbool, optional

Whether or not to transform the state to the full state in basis_right. Default is True.

Returns:
numpy.ndarray

Transformation/projector between the symmetry-reduced and the full basis.

Notes

Particularly useful when a given operation canot be carried away in the symmetry-reduced basis in a straightforward manner.

Examples

>>> P = get_proj(np.float64)
>>> print(P.shape)
get_vec(v0, sparse=True, full_left=True, full_right=True)[source]

DEPRECATED (cf project_from). Transforms state from symmetry-reduced basis to full (symmetry-free) basis.

Notes

This function is deprecated. Use project_from() instead (the inverse function, project_to(), is currently available in the basis_general classes only).

index(*states)[source]

Finds the index of user-defined Fock state in tensor basis.

Parameters:
stateslist(str)

List of strings which separately define the Fock state in each of the basis used to construct the tensor_basis object.

Returns:
int

Position of tensor Fock state in the tensor_basis.

Notes

Particularly useful for defining initial Fock states through a unit vector in the direction specified by index().

Examples

>>> s_1 = "".join("1" for i in range(2)) + "".join("0" for i in range(2))
>>> s_2 = "".join("1" for i in range(4))
>>> print( basis.index(s_1,s_2) )
inplace_Op(v_in, op_list, dtype, transposed=False, conjugated=False, a=1.0, v_out=None)

Calculates the action of an operator on a state.

Parameters:
v_inarray_like

state (or states stored in columns) to act on with the operator.

op_listlist

Operator string list which defines the operator to apply. Follows the format [[“z”,[i],Jz[i]] for i in range(L)], [“x”,[i],Jx[j]] for j in range(L)],…].

dtype‘type’

Data type (e.g. numpy.float64) to construct the operator with.

transposedbool, optional

if True this function will act with the trasposed operator.

conjugatedbool, optional

if True this function will act with the conjugated operator.

ascalar, optional

value to rescale resulting vector after performing the action of the operators. Same as rescaling all couplings by value a.

v_outarray_like

output array, must be the same shape as v_in and must match the type of the output.

Returns:
numpy.ndarray
  • if v_out is not None, this function modifies v_out inplace and returns it.

Notes

This function works with the tensor_basis and other basis which use the “|” symbol in the opstr.

Examples

>>> J = 1.41
>>> indx = [2,3]
>>> opstr = "zz"
>>> dtype = np.float64
>>> op_list=[[opstr,indx,J]]
>>> ME, row, col = inplace_Op(op_list,dtype)
property operators

set of available operator strings.

Type:

set

partial_trace(state, sub_sys_A='left', return_rdm='A', enforce_pure=False, sparse=False)[source]

Calculates reduced density matrix, through a partial trace of a quantum state in tensor_basis.

Parameters:
stateobj

State of the quantum system. Can be either one of:

  • numpy.ndarray [shape (Ns,)]: pure state (default).

  • numpy.ndarray [shape (Ns,Ns)]: density matrix (DM).

sub_sys_Astr, optional

Defines subsystem A. Can be either one of:

  • “left”: refers to basis_left (Default).

  • “right”: refers to basis_right.

  • “both”: for initial mixed states the Renyi entropy of subsystem A and its complement

    B need not be the same. This option automatically sets return_rdm=both.

return_rdmstr, required

Toggles returning the reduced DM. Can be either one of:

  • “A”: returns reduced DM of subsystem A.

  • “B”: returns reduced DM of subsystem B (complement of A).

  • “both”: returns reduced DM of both subsystems A and B.

enforce_purebool, optional

Whether or not to assume state is a collection of pure states or a mixed density matrix, if it is a square array. Default is False.

sparsebool, optional

Whether or not to return a sparse DM. Default is False.

Returns:
numpy.ndarray

Density matrix associated with state. Depends on optional arguments.

Examples

>>> partial_trace(state,sub_sys_A=None,return_rdm="A",enforce_pure=False,sparse=False)
project_from(v0, sparse=True, full_left=True, full_right=True)[source]

Transforms state from symmetry-reduced basis to full (symmetry-free) basis.

Parameters:
v0numpy.ndarray

Contains in its columns the states in the symmetry-reduced basis.

sparsebool, optional

Whether or not the output should be in sparse format. Default is True.

full_leftbool, optional

Whether or not to transform the state to the full state in basis_left. Default is True.

full_rightbool, optional

Whether or not to transform the state to the full state in basis_right. Default is True.

Returns:
numpy.ndarray

Array containing the state v0 in the full basis.

Notes

Particularly useful when a given operation canot be carried away in the symmetry-reduced basis in a straightforward manner.

Supports parallelisation to multiple states listed in the columns.

Examples

>>> v_full = project_from(v0)
>>> print(v_full.shape, v0.shape)
property sps

number of states per site (ie, the on-site Hilbert space dimension).

Type:

int