quspin.basis.spin_basis_1d

class quspin.basis.spin_basis_1d(L, Nup=None, m=None, S='1/2', pauli=1, **blocks)[source]

Bases: basis_1d

Constructs basis for spin operators in a specified 1-d symmetry sector.

The supported operator strings for spin_basis_1d are:

\[\begin{array}{cccc} \texttt{basis}/\texttt{opstr} & \texttt{"I"} & \texttt{"+"} & \texttt{"-"} & \texttt{"z"} & \texttt{"x"} & \texttt{"y"} \newline \texttt{spin_basis_1d} & \hat{1} & \hat\sigma^+ & \hat\sigma^- & \hat\sigma^z & (\hat\sigma^x) & (\hat\sigma^y) \ \newline \end{array}\]
Note:
  • The relation between spin and Pauli matrices is \(\vec S = \vec \sigma/2\).

  • The default operators for spin-1/2 are the Pauli matrices, NOT the spin operators. To change this, see the argument pauli of the spin_basis class. Higher spins can only be defined using the spin operators, and do NOT support the operator strings “x” and “y”.

Examples

The code snippet below shows how to use the spin_boson_1d class to construct the basis in the zero momentum sector of positive parity for the spin Hamiltonian.

\[H(t) = \sum_j J\sigma^z_{j+1}\sigma^z_j + h\sigma^z_j + g\cos\Omega t\sigma^x_j\]
 1from quspin.basis import spin_basis_1d  # Hilbert space spin basis
 2import numpy as np  # generic math functions
 3
 4#
 5##### define model parameters #####
 6L = 4  # system size
 7J = 1.0  # spin interaction
 8g = 0.809  # transverse field
 9h = 0.9045  # longitudinal field
10##### define periodic drive #####
11Omega = 4.5  # drive frequency
12
13
14def drive(t, Omega):
15    return np.cos(Omega * t)
16
17
18drive_args = [Omega]
19#
20##### construct basis in the 0-total momentum and +1-parity sector
21basis = spin_basis_1d(L=L, a=1, kblock=0, pblock=1)
22print(basis)
23# define PBC site-coupling lists for operators
24x_field = [[g, i] for i in range(L)]
25z_field = [[h, i] for i in range(L)]
26J_nn = [[J, i, (i + 1) % L] for i in range(L)]  # PBC
27# static and dynamic lists
28static = [["zz", J_nn], ["z", z_field]]
29dynamic = [["x", x_field, drive, drive_args]]
30###### construct Hamiltonian
31H = hamiltonian(static, dynamic, dtype=np.float64, basis=basis)
__init__(L, Nup=None, m=None, S='1/2', pauli=1, **blocks)[source]

Intializes the spin_basis_1d object (basis for spin operators).

Parameters:
L: int

Length of chain/number of sites.

Nup: {int,list}, optional

Total magnetisation, \(\sum_j S^z_j\), projection. Can be integer or list to specify one or more particle sectors. Negative values are taken to be subtracted from the fully polarized up state as: Nup_max + Nup + 1. e.g. to get the Nup_max state use Nup = -1, for Nup_max-1 state use Nup = -2, etc.

m: float, optional

Density of spin up in chain (spin up per site).

S: str, optional

Size of local spin degrees of freedom. Can be any (half-)integer from: “1/2”,”1”,”3/2”,…,”9999/2”,”5000”.

pauli: int, optional (requires `S=”1/2”`)
  • for pauli=0 the code uses spin-1/2 operators:

\[\begin{split}S^x = \frac{1}{2}\begin{pmatrix}0 & 1\\ 1 & 0\end{pmatrix},\quad S^y = \frac{1}{2}\begin{pmatrix}0 & -i\\ i & 0\end{pmatrix},\quad S^z = \frac{1}{2}\begin{pmatrix}1 & 0\\ 0 & -1\end{pmatrix},\quad S^+ = \begin{pmatrix}0 & 1\\ 0 & 0\end{pmatrix},\quad S^- = \begin{pmatrix}0 & 0\\ 1 & 0\end{pmatrix}\end{split}\]
  • for pauli=1 the code uses Pauli matrices with:

\[\begin{split}\sigma^x = \begin{pmatrix}0 & 1\\ 1 & 0\end{pmatrix},\quad \sigma^y = \begin{pmatrix}0 & -i\\ i & 0\end{pmatrix},\quad \sigma^z = \begin{pmatrix}1 & 0\\ 0 & -1\end{pmatrix},\quad \sigma^+ = \begin{pmatrix}0 & 2\\ 0 & 0\end{pmatrix},\quad \sigma^- = \begin{pmatrix}0 & 0\\ 2 & 0\end{pmatrix}\end{split}\]
  • for pauli=-1 the code uses Pauli matrices with:

\[\begin{split}\sigma^x = \begin{pmatrix}0 & 1\\ 1 & 0\end{pmatrix},\quad \sigma^y = \begin{pmatrix}0 & -i\\ i & 0\end{pmatrix},\quad \sigma^z = \begin{pmatrix}1 & 0\\ 0 & -1\end{pmatrix},\quad \sigma^+ = \begin{pmatrix}0 & 1\\ 0 & 0\end{pmatrix},\quad \sigma^- = \begin{pmatrix}0 & 0\\ 1 & 0\end{pmatrix}\end{split}\]
**blocks: optional

extra keyword arguments which include:

a (int) - specifies unit cell size for translation.

kblock (int) - specifies momentum block. The physical manifestation of this symmetry transformation is translation by a lattice sites.

pblock (int) - specifies parity block. The physical manifestation of this symmetry transformation is reflection about the middle of the chain.

zblock (int) - specifies spin inversion symmetry block. The physical manifestation of this symmetry transformation is to flip the sign of the spin-z component.

pzblock (int) - specifies parity followed by spin inversion symmetry block. The physical manifestation of this symmetry transformation is reflection about the middle of the chain and a simultaneous flip of the sign of the spin-z component.

zAblock (int) - specifies spin inversion symmetry block for sublattice A (defined as all even lattice sites). The physical manifestation of this symmetry transformation is to flip the sign of the spin-z component on all even sites.

zBblock (int) - specifies spin inversion symmetry block for sublattice B (defined as all odd lattice sites). The physical manifestation of this symmetry transformation is to flip the sign of the spin-z component on all odd sites.

Methods

Op(opstr, indx, J, dtype)

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

__init__(L[, Nup, m, S, pauli])

Intializes the spin_basis_1d object (basis for spin operators).

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, density, ...])

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[, pcon])

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

get_vec(v0[, sparse, pcon])

DEPRECATED (cf project_from).

index(s)

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

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

Calculates the action of an operator on a state.

int_to_state(state[, bracket_notation])

Finds string representation of a state defined in integer representation.

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

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

project_from(v0[, sparse, pcon])

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

state_to_int(state)

Finds integer representation of a state defined in string format.

Attributes

L

length of lattice.

N

number of sites the basis is constructed with.

Ns

number of states in the Hilbert space.

blocks

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

description

information about basis object.

dtype

data type of basis state integers.

noncommuting_bits

list of bits that represent sites that do not commute along with the phase required from commuting sites

operators

set of available operator strings.

sps

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

states

basis states stored in their integer representation.

property L

length of lattice.

Type:

int

property N

number of sites the basis is constructed with.

Type:

int

property Ns

number of states in the Hilbert space.

Type:

int

Op(opstr, indx, J, dtype)

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

Parameters:
opstrstr

Operator string in the lattice basis format. For instance:

>>> opstr = "zz"
indxlist(int)

List of integers to designate the sites the lattice basis operator is defined on. For instance:

>>> indx = [2,3]
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 lattice 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 lattice basis,

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

Examples

>>> J = 1.41
>>> indx = [2,3]
>>> opstr = "zz"
>>> dtype = np.float64
>>> ME, row, col = Op(opstr,indx,J,dtype)
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.

property description

information about basis object.

Type:

str

property dtype

data type of basis state integers.

Type:

numpy.dtype

ent_entropy(state, sub_sys_A=None, density=True, subsys_ordering=True, return_rdm=None, enforce_pure=False, return_rdm_EVs=False, sparse=False, alpha=1.0, sparse_diag=True, maxiter=None, svd_solver=None, svd_kwargs=None)

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

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

where the normalization \(N_A\) can be switched on and off using the optional argument density. This expression reduces to the familiar von Neumann entropy in the limit \(\alpha=1\).

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_Atuple/list, optional

Defines the sites contained in subsystem A [by python convention the first site of the chain is labelled j=0]. Default is tuple(range(N//2)) with N the number of lattice sites.

densitybool, optional

Toggles whether to return entanglement entropy normalized by the number of sites in the subsystem.

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.

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

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.

subsys_orderingbool, optional

Whether or not to reorder the sites in sub_sys_A in ascending order. Default is True.

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().

svd_solverobject, optional

Specifies the svd solver to be used, e.g. numpy.linalg.svd or scipy.linalg.svd, or a custom solver. Effective when enforce_pure=True or sparse=False.

svd_kwargsdict, optional

Specifies additional arguments for svd_solver.

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=[0,3,4,7],return_rdm="A",enforce_pure=False,return_rdm_EVs=False,
>>>                             sparse=False,alpha=1.0,sparse_diag=True,subsys_ordering=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, pcon=False)

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.

sparsebool, optional

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

pconbool, optional

Whether or not to return the projector to the particle number (magnetisation) conserving basis (useful in bosonic/single particle systems). Default is pcon=False.

Returns:
scipy.sparse.csc_matrix

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,pcon=False)
>>> print(P.shape)
get_vec(v0, sparse=True, pcon=False)

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(s)

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

Parameters:
s{str, int}

Defines the Fock state with number of particles (spins) per site in underlying lattice basis.

Returns:
int

Position of the Fock state in the lattice basis.

Notes

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

Examples

>>> i0 = index("111000") # pick state from basis set
>>> print(basis)
>>> print(i0)
>>> psi = np.zeros(basis.Ns,dtype=np.float64)
>>> psi[i0] = 1.0 # define state corresponding to the string "111000"
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)
int_to_state(state, bracket_notation=True)

Finds string representation of a state defined in integer representation.

Parameters:
stateint

Defines the Fock state in integer representation in underlying lattice basis.

bracket_notationbool, optional

Toggles whether to return the state in |str> notation.

Returns:
str

String corresponding to the Fock state in the lattice basis.

Notes

This function is the inverse of state_to_int.

Examples

>>> s = basis[0] # pick state from basis set
>>> s_str = basis.int_to_state(s)
>>> print(s_str)
property noncommuting_bits

list of bits that represent sites that do not commute along with the phase required from commuting sites

Type:

list

property operators

set of available operator strings.

Type:

set

partial_trace(state, sub_sys_A=None, subsys_ordering=True, return_rdm='A', enforce_pure=False, sparse=False)

Calculates reduced density matrix, through a partial trace of a quantum state in a lattice 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_Atuple/list, optional

Defines the sites contained in subsystem A [by python convention the first site of the chain is labelled j=0]. Default is tuple(range(N//2)) with N the number of lattice sites.

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.

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

subsys_orderingbool, optional

Whether or not to reorder the sites in sub_sys_A in ascending order. Default is True.

enforce_purebool, optional

Whether or not to assume state is a colelction 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.

Notes

This function can also be applied to trace out operators/observables defined by the input state, in which case one has to additionally normalize the final output by the Hilbert space dimension of the traced-out space. However, if an operator is defined in a symmetry-reduced basis, there is a caveat. In such a case, one has to:
  1. use the basis.get_proj() function to lift the operator to the full basis;

  2. apply basis.partial_trace();

  3. repeat this procedure for all symmetry sectors, and sum up the resulting reduced operators [this is becauce one has to add in the information about how the operator acts on the full Hilbert space].

Examples

>>> partial_trace(state,sub_sys_A=tuple(range(basis.N//2),return_rdm="A",enforce_pure=False,sparse=False,subsys_ordering=True)
project_from(v0, sparse=True, pcon=False)

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.

pconbool, optional

Whether or not to return the output in the particle number (magnetisation) conserving basis (useful in bosonic/single particle systems). Default is pcon=False.

Returns:
numpy.ndarray

Array containing the state v0 in the full basis.

Notes

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

Supports parallelisation to multiple states listed in the columns.

Examples

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

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

Type:

int

state_to_int(state)

Finds integer representation of a state defined in string format.

Parameters:
statestr

Defines the Fock state with number of particles (spins) per site in underlying lattice basis.

Returns:
int

Integer corresponding to the Fock state in the lattice basis.

Notes

This function is the einverse of int_to_state.

Examples

>>> s_str = "111000" # pick state from basis set
>>> s = basis.state_to_int(s_str)
>>> print(s)
property states

basis states stored in their integer representation.

Type:

numpy.ndarray(int)