quspin.basis.user_basis

class quspin.basis.user_basis(basis_dtype, N, op_dict, sps=2, pcon_dict=None, pre_check_state=None, allowed_ops=None, parallel=False, Ns_block_est=None, _make_basis=True, block_order=None, noncommuting_bits=[], _Np=None, **blocks)[source]

Bases: basis_general

Constructs basis for USER-DEFINED functionality of a basis object.

The user_basis unveils the inner workings of QuSpin. This is the most advanced usage of the package, and requires some understanding of python, the numba package used to interface QuSpin’s underlying cpp code with python, and some experience with bitwise operations to manipulate integers.

Since we believe that the users will benefit from a more detailed discussion on how the user_basis is intended to work, we also provide a detailed tutorial: user_basis tutorial, which covers the general concepts and provides six complete examples of various complexity.

Examples

The following example shows how to use the user_basis class to construct the Hamiltonian

\[H = \sum_j P_{j-1}\sigma^x_j P_{j+1},\quad P_j = |\downarrow_j\rangle\langle\downarrow_j|\]

using translation and reflection symmetry. The projector operator \(P_j\), which only allows a spin-up state in the basis to be preceded and succeeded by a spin-down, is incorporated by constructing the corresponding user_basis object. One can then just build the Hamiltonian \(H=\sum_j\sigma^x_j\) in the constrained Hilbert space.

More examples (including explanations of the class methods and attributes) can be found at: user_basis tutorial.

  1#
  2from quspin.operators import hamiltonian
  3from quspin.basis.user import user_basis  # Hilbert space user basis
  4from quspin.basis.user import (
  5    next_state_sig_32,
  6    pre_check_state_sig_32,
  7    op_sig_32,
  8    map_sig_32,
  9)  # user_basis dtypes
 10from numba import carray, cfunc  # numba helper functions
 11from numba import uint32, int32  # numba data types
 12import numpy as np
 13
 14#
 15N = 14  # lattice sites
 16
 17
 18#
 19######  function to call when applying operators
 20@cfunc(op_sig_32, locals=dict(s=int32, b=uint32))
 21def op(op_struct_ptr, op_str, ind, N, args):
 22    # using struct pointer to pass op_struct_ptr back to C++ see numba Records
 23    op_struct = carray(op_struct_ptr, 1)[0]
 24    err = 0
 25    ind = N - ind - 1  # convention for QuSpin for mapping from bits to sites.
 26    s = (((op_struct.state >> ind) & 1) << 1) - 1
 27    b = 1 << ind
 28    #
 29    if op_str == 120:  # "x" is integer value 120 (check with ord("x"))
 30        op_struct.state ^= b
 31    elif op_str == 121:  # "y" is integer value 120 (check with ord("y"))
 32        op_struct.state ^= b
 33        op_struct.matrix_ele *= 1.0j * s
 34    elif op_str == 122:  # "z" is integer value 120 (check with ord("z"))
 35        op_struct.matrix_ele *= s
 36    else:
 37        op_struct.matrix_ele = 0
 38        err = -1
 39    #
 40    return err
 41
 42
 43#
 44op_args = np.array([], dtype=np.uint32)
 45
 46
 47#
 48######  function to filter states/project states out of the basis
 49#
 50@cfunc(
 51    pre_check_state_sig_32,
 52    locals=dict(s_shift_left=uint32, s_shift_right=uint32),
 53)
 54def pre_check_state(s, N, args):
 55    """imposes that that a bit with 1 must be preceded and followed by 0,
 56    i.e. a particle on a given site must have empty neighboring sites.
 57    #
 58    Works only for lattices of up to N=32 sites (otherwise, change mask)
 59    #
 60    """
 61    mask = 0xFFFFFFFF >> (32 - N)  # works for lattices of up to 32 sites
 62    # cycle bits left by 1 periodically
 63    s_shift_left = ((s << 1) & mask) | ((s >> (N - 1)) & mask)
 64    #
 65    # cycle bits right by 1 periodically
 66    s_shift_right = ((s >> 1) & mask) | ((s << (N - 1)) & mask)
 67    #
 68    return (((s_shift_right | s_shift_left) & s)) == 0
 69
 70
 71#
 72pre_check_state_args = None
 73
 74
 75#
 76######  define symmetry maps
 77#
 78@cfunc(
 79    map_sig_32,
 80    locals=dict(
 81        shift=uint32,
 82        xmax=uint32,
 83        x1=uint32,
 84        x2=uint32,
 85        period=int32,
 86        l=int32,
 87    ),
 88)
 89def translation(x, N, sign_ptr, args):
 90    """works for all system sizes N."""
 91    shift = args[0]  # translate state by shift sites
 92    period = N  # periodicity/cyclicity of translation
 93    xmax = args[1]
 94    #
 95    l = (shift + period) % period
 96    x1 = x >> (period - l)
 97    x2 = (x << l) & xmax
 98    #
 99    return x2 | x1
100
101
102T_args = np.array([1, (1 << N) - 1], dtype=np.uint32)
103
104
105#
106@cfunc(
107    map_sig_32,
108    locals=dict(
109        out=uint32,
110        s=int32,
111    ),
112)
113def parity(x, N, sign_ptr, args):
114    """works for all system sizes N."""
115    out = 0
116    s = args[0]  # N-1
117    #
118    out ^= x & 1
119    x >>= 1
120    while x:
121        out <<= 1
122        out ^= x & 1
123        x >>= 1
124        s -= 1
125    #
126    out <<= s
127    return out
128
129
130P_args = np.array([N - 1], dtype=np.uint32)
131#
132######  construct user_basis
133# define maps dict
134maps = dict(
135    T_block=(translation, N, 0, T_args),
136    P_block=(parity, 2, 0, P_args),
137)
138# define particle conservation and op dicts
139op_dict = dict(op=op, op_args=op_args)
140# define pre_check_state
141pre_check_state = (
142    pre_check_state,
143    pre_check_state_args,
144)  # None gives a null pinter to args
145# create user basis
146basis = user_basis(
147    np.uint32,
148    N,
149    op_dict,
150    allowed_ops=set("xyz"),
151    sps=2,
152    pre_check_state=pre_check_state,
153    Ns_block_est=300000,
154    **maps,
155)
156# print basis
157print(basis)
158#
159###### construct Hamiltonian
160# site-coupling lists
161h_list = [[1.0, i] for i in range(N)]
162# operator string lists
163static = [
164    ["x", h_list],
165]
166# compute Hamiltonian, no checks have been implemented
167no_checks = dict(check_symm=False, check_pcon=False, check_herm=False)
168H = hamiltonian(static, [], basis=basis, dtype=np.float64, **no_checks)
__init__(basis_dtype, N, op_dict, sps=2, pcon_dict=None, pre_check_state=None, allowed_ops=None, parallel=False, Ns_block_est=None, _make_basis=True, block_order=None, noncommuting_bits=[], _Np=None, **blocks)[source]

Intializes the user_basis_general object (basis for user defined ED calculations).

Parameters:
basis_dtype: numpy.dtype object

the data type used to represent the states in the basis: must be either uint32 or uint64.

N: int

Number of sites.

op_dict: dict
used to define the basis.Op function; the dictionary contais the following items:
  • op(op_struct_ptr,op_str,site_ind,N,args): numba.CFunc object

    This is a numba-compiled function (CFunc) which calculates the matrix elements \(\mathrm{me}\) given a state \(|s\rangle\) together with a character to represent the operator, an integer site_ind specifying the site of that local operator, the total number of sites N, and a set of optional uint-dtype arguments args. See the above example for how one would use this for spin-1/2 system.

  • op_args: np.ndarray[basis_dtype]

    used to pass the arguments args to the CFunc op(…,args). The corresponding key must be a np.ndarray[basis_dtype].

pcon_dict: dict, optional
This dictionary contains the following items which are required to use particle conservation in this basis:
minimum requirements:
  • Np: tuple/int, list(tuple/int)

    specifies the particle sector(s).

  • next_state(s,counter,N,args): numba.CFunc object

    given a quantum state \(|s\rangle\) in the integer-representation s, this CFunc generates the next lexicographically ordered particle conservation state. counter is an intrinsic variable which increments by unity every time the function is called, N is the total number of lattice sites, and args holds any optional arguments stored in a np.ndarray[basis_dtype].

  • next_state_args: np.ndarray(basis_dtype)

    optional arguments for next_state(…,args).

  • get_Ns_pcon(N,Np): python function

    when called as get_Ns_pcon(N,Np), this python function returns the size of the symmetery-free particle conservation basis, given the N lattice sites and Np (see above).

  • get_s0_pcon(N,Np): python function

    when called as get_s0_pcon(N,Np), this python function returns the starting state to generate the whole particle conservation basis by repeatedly calling next_state().

advanced requirements to access basis.Op_bra_ket() functionality (on top of the minimum requirements):
  • n_sectors: int, list(int)

    number of integers which parameterize the particle sectors, e.g. with spinful fermions there is a particle number for both the up and the down sectors, and hence n_sectors=2.

  • count_particles(s,p_number_ptr,args): numba.CFunc object

    For a quantum state s in the integer representation, this CFunc counts the number of particles in each particle sector and places them into a pointer p_number_ptr (count_particles does not return any output). The pointer provided will have n_sector slots of memory allocated. The components of the pointer p_number_ptr must correspond to the ordering of Np. The integer s cannot be changed.

  • count_particles_args: np.ndarray(int)

    compulsory arguments for count_particles(…,args) (whenever used).

pre_check_state(s,N,args): numba.CFunc object or tuple(numba.CFunc object,ndarray(C-contiguous,dtype=basis_dtype)), optional

This CFunc allows the user to specify a boolean criterion used to discard/filter states from the basis. In the low-level code, this function is applied before checking if a given state is representative state (i.e. belogs to a given symmetry sector) or not. This allows the user to, e.g., enforce a local Hilbert-space constraint (e.g. for a spinful fermion basis to never have a doubly occupied site). One can pass additional arguments args using a np.ndarray[basis_dtype].

allowed_ops: list/set, optional

A list of allowed characters, each of which is to be passed in to the op in the form of op_str (see above).

parallel: bool, optional

turns on parallel code when constructing the basis even when no symmetries are implemented. Useful when constructing highly constrained Hilbert spaces with pre_check_state.

sps: int, optional

The number of states per site (i.e. the local on-site Hilbert space dimension).

Ns_full: int, optional

Total number of states in the Hilbert space without any symmetries. For a single species this value is sps**N

Ns_block_est: int, optional

An estimate for the size of the symmetry reduced block, QuSpin does a simple estimate which is not always correct.

block_order: tuple/list, optional

A list of strings containing the names of the symmetry blocks which specifies the order in which the symmetries will be applied to the state when calculating the basis. The first element in the list is applied to the state first followed by the second element, etc. If the list is not specificed the ordering is such that the symmetry with the largest cycle is the first, followed by the second largest, etc.

noncommuting_bits: list, optional

A list of tuples specifying if bits belong to a group of sites that do not commute. The first element in each tuple represents the group of sites, and the second element represents the phase-factor that is given during the exchange.

**blocks: optional

keyword arguments which pass the symmetry generator arrays. For instance:

>>> basis(...,kxblock=(CFunc,m_Q,q,args),...)

The key names of the symmetry sector, e.g. kxblock, can be defined arbitrarily by the user. The values are tuples where the first entry contains the numba-CFunc which generates the symmetry transformation \(Q\) acting on the state (see class example), the second entry is an integer \(m_Q\) which gives the periodicity of the symmetry sector (\(Q^{m_Q} = 1\)), and \(q\) is the quantum number for the given sector. Optional arguments can be passed using the`args` argument which is a np.ndarray[basis_dtype]. Note that if the periodicity is wrong the basis will give undefined behavior.

Methods

Op(opstr, indx, J, dtype)

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

Op_bra_ket(opstr, indx, J, dtype, ket_states)

Finds bra states which connect given ket states by operator from a site-coupling list and an operator string.

Op_shift_sector(other_basis, op_list, v_in)

Applies symmetry non-conserving operator to state in symmetry-reduced basis.

__init__(basis_dtype, N, op_dict[, sps, ...])

Intializes the user_basis_general object (basis for user defined ED calculations).

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_amp(states[, out, amps, mode])

Computes the rescale factor of state amplitudes between the symmetry-reduced and full basis.

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.

make([Ns_block_est, N_p])

Creates the entire basis by calling the basis constructor.

make_basis_blocks([N_p])

Creates/modifies the bounds for representatives based on prefix tages.

normalization(states[, out])

Computes normalization of basis states.

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.

project_to(v0[, sparse, pcon])

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

representative(states[, out, return_g, ...])

Maps states to their representatives under the basis symmetries.

state_to_int(state)

Finds integer representation of a state defined in string format.

Attributes

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 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)
Op_bra_ket(opstr, indx, J, dtype, ket_states, reduce_output=True)

Finds bra states which connect given ket states by operator from a site-coupling list and an operator string.

Given a set of ket states \(|s\rangle\), the function returns the bra states \(\langle s'|\) which connect to them through an operator, together with the corresponding matrix elements.

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 matrix elements with.

ket_statesnumpy.ndarray(int)

Ket states in integer representation. Must be of same data type as basis.

reduce_output: bool, optional

If set to False, the returned arrays have the same size as ket_states; If set to True zeros are purged.

Returns:
tuple
(ME,bra,ket), where
  • numpy.ndarray(scalar): ME: matrix elements of type dtype, which connects the ket and bra states.

  • numpy.ndarray(int): bra: bra states, obtained by applying the matrix representing the operator in the lattice basis,

    to the ket states, such that bra[i] corresponds to ME[i] and connects to ket[i].

  • numpy.ndarray(int): ket: ket states, such that ket[i] corresponds to ME[i] and connects to bra[i].

Notes

  • Similar to Op but instead of returning the matrix indices (row,col), it returns the states (bra,ket) in integer representation.

  • Does NOT require the full basis (see basis optional argument make_basis).

  • If a state from ket_states does not have a non-zero matrix element, it is removed from the returned list. See otional argument reduce_output.

Examples

>>> J = 1.41
>>> indx = [2,3]
>>> opstr = "zz"
>>> dtype = np.float64
>>> ME, bra, ket = Op_bra_ket(opstr,indx,J,dtype,ket_states)
Op_shift_sector(other_basis, op_list, v_in, v_out=None, dtype=None)

Applies symmetry non-conserving operator to state in symmetry-reduced basis.

An operator, which does not conserve a symmetry, induces a change in the quantum number of a state defined in the corresponding symmetry sector. Hence, when the operator is applied on a quantum state, the state shifts the symmetry sector. Op_shift_sector() handles this automatically.

NOTE: One has to make sure that (i) the operator moves the state between the two sectors, and (ii) the two bases objects have the same symmetries. This function will not give the correct results otherwise.

Formally equivalent to:

>>> P1 = basis_sector_1.get_proj(np.complex128) # projector between full and initial basis
>>> P2 = basis_sector_2.get_proj(np.complex128) # projector between full and target basis
>>> v_in_full = P1.dot(v_in) # go from initial basis to to full basis
>>> v_out_full = basis_full.inplace_Op(v_in_full,op_list,np.complex128) # apply Op
>>> v_out = P2.T.conj().dot(v_out_full) # project to target basis
Parameters:
other_basisbasis object

basis_general object for the initial symmetry sector. Must be the same basis class type as the basis whose instance is Op_shift_sector() (i.e. the basis in basis.Op_shift_sector()).

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

v_inarray_like, (other_basis.Ns,…)

Initial state to apply the symmetry non-conserving operator on. Must have the same length as other_basis.Ns.

v_outarray_like, (basis.Ns,…), optional

Optional array to write the result for the final/target state in.

dtypenumpy dtype for matrix elements, optional

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

Returns:
(basis.Ns, ) numpy.ndarray

Array containing the state v_out in the current basis, i.e. the basis in basis.Op_shift_sector().

Notes

  • particularly useful when computing correlation functions.

  • supports parallelization to multiple states listed in the columns of v_in.

  • the user is strongly advised to use the code under “Formally equivalent” above to check the results of this function for small system sizes.

Examples

>>> v_out = basis.Op_shift_sector(initial_basis, op_list, v_in)
>>> print(v_out.shape, basis.Ns, v_in.shape, initial_basis.Ns)
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_amp(states, out=None, amps=None, mode='representative')

Computes the rescale factor of state amplitudes between the symmetry-reduced and full basis.

Given a quantum state \(s\) and a state amplitude in the full basis \(\psi_s\), its representative (under the symemtries) \(r(s)\) with a corresponding amplitude \(\psi^\text{sym}_r\), the function computes the ratio \(C\), defined as

\[\psi_s = C\psi_r^\text{sym}\]
Parameters:
statesarray_like(int)

Fock-basis (z-basis) states to find the amplitude rescale factor \(C\) of. States are stored in integer representations.

outnumpy.ndarray(float), optional

variable to store the rescale factors \(C\) of the states in. Must be a real or complex-valued numpy.ndarray of the same shape as states.

ampsnumpy.ndarray(float), optional

array of amplitudes to rescale by the amplitude factor \(C\) (see mode). Updated in-place. Must be a real or complex-valued numpy.ndarray of the same shape as states.

modestring, optional
  • if mode=’representative’ (default), then the function assumes that
    1. states already contains representatives (i.e. states in the symmetry-reduced basis);

    2. amps (if passed) are amplitudes in the symmetry-reduced basis (\(\psi_r^\text{symm}\)). The function will update amps in-place to \(\psi_s\).

  • if mode=’full_basis’, then the function assumes that
    1. states contains full-basis states (the funciton will compute the corresponding representatives);

    2. amps (if passed) are amplitudes in the full basis (\(\psi_s\)). The function will update amps in-place to \(\psi_r^\text{symm}\);

      Note: the function will also update the variable states in place with the corresponding representatives.

Returns:
array_like(float)

amplitude rescale factor \(C\) (see expression above).

Notes

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

  • To transform an entire state from a symmetry-reduced basis to the full (symmetry-free) basis, use the basis.get_vec() function.

  • Returns zero, if the state passed to the function is not part of the symmetry-reduced basis.

  • If amps is passed, the user has to make sure that the input data in amps correspond to the states.

  • The function assumes that states comply with the particle conservation symmetry the basis was constructed with.

Examples

>>> C = get_amp(states,out=None,amps=None,mode='representative')
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.

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 out in the symmetry-reduced basis in a straightforward manner.

  • see also Op_shift_sector().

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; see also the inverse function project_to().

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)
make(Ns_block_est=None, N_p=None)

Creates the entire basis by calling the basis constructor.

Parameters:
Ns_block_est: int, optional

Overwrites the internal estimate of the size of the reduced Hilbert space for the given symmetries. This can be used to help conserve memory if the exact size of the H-space is known ahead of time.

N_p: int, optional

number of bits to use in the prefix label used to generate blocks for searching positions of representatives.

Returns:
int

Total number of states in the (symmetry-reduced) Hilbert space.

Notes

The memory stored in the basis grows exponentially as exactly \(2^{N_p+1}\). The default behavior is to use N_p such that the size of the stored information for the representative bounds is approximately as large as the basis. This is not as effective for basis which small particle numbers as the blocks have very uneven sizes. To not use the blocks just set N_p=0.

Examples

>>> N, Nup = 8, 4
>>> basis=spin_basis_general(N,Nup=Nup,make_basis=False)
>>> print(basis)
>>> basis.make()
>>> print(basis)
make_basis_blocks(N_p=None)

Creates/modifies the bounds for representatives based on prefix tages.

Parameters:
N_p: int, optional

number of bits to use in the prefix label used to generate blocks for searching positions of representatives.

Notes

The memory stored in the basis grows exponentially as exactly \(2^{N_p+1}\). The default behavior is to use N_p such that the size of the stored information for the representative bounds is approximately as large as the basis. This is not as effective for basis which small particle numbers as the blocks have very uneven sizes. To not use the blocks just set N_p=0.

Examples

>>> N, Nup = 8, 4
>>> basis=spin_basis_general(N,Nup=Nup,make_basis=False)
>>> print(basis)
>>> basis.make()
>>> print(basis)
property noncommuting_bits

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

Type:

list

normalization(states, out=None)

Computes normalization of basis states.

Parameters:
statesarray_like(int)

Fock-basis (z-basis) states to find the normalizations of. States are stored in integer representations.

outnumpy.ndarray(unsigned int), optional

variable to store the normalizations of the states in. Must be a numpy.ndarray of datatype unsigned int (e.g. numpy.uint16), and same shape as states.

Returns:
array_like(int)

normalizations of states for the given (symmetry-reduced) basis.

Notes

  • Returns zero, if the state is not part of the symmetry-reduced basis.

  • The normalizations can be used to compute matrix elements in the symmetry-reduced basis.

Examples

>>> basis=spin_basis_general(N,Nup=Nup,make_basis=False)
>>> s = 17
>>> norm_s = basis.normalization(s)
>>> print(s,norm_s)
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 cannot be carried out in the symmetry-reduced basis in a straightforward manner.

  • supports parallelisation to multiple states listed in the columns.

  • inverse function to project_to.

Examples

>>> v_full = project_from(v0)
>>> print(v_full.shape, v0.shape)
project_to(v0, sparse=True, pcon=False)

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

Parameters:
v0numpy.ndarray

Contains in its columns the states in the full (symmetry-free) 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 symmetry-reduced basis.

Notes

  • particularly useful when a given operation cannot be carried out in the full basis.

  • supports parallelisation to multiple states listed in the columns.

  • inverse function to project_from.

Examples

>>> v_symm = project_to(v0)
>>> print(v_symm.shape, v0.shape)
representative(states, out=None, return_g=False, return_sign=False)

Maps states to their representatives under the basis symmetries.

Parameters:
statesarray_like(int)

Fock-basis (z-basis) states to find the representatives of. States are stored in integer representations.

outnumpy.ndarray(int), optional

variable to store the representative states in. Must be a numpy.ndarray of same datatype as basis, and same shape as states.

return_gbool, optional

if set to True, the function also returns the integer g corresponding to the number of times each basis symmetry needs to be applied to a given state to obtain its representative.

return_signbool, optional

if set to True, the function returns the sign of the representative relative to the original state (nontrivial only for fermionic bases).

Returns:
tuple

( representatives, g_array, sign_array ) * array_like(int): representatives: Representatives under basis symmetries, corresponding to states. * array_like(int): g_array of size (number of states, number of symmetries). Requires return_g=True. Contains integers corresponding to the number of times each basis symmetry needs to be applied to a given state to obtain its representative. * array_like(int): sign_array of size (number of states,). Requires return_sign=True. Contains sign of the representative relative to the original state (nontrivial only for fermionic bases).

Examples

>>> basis=spin_basis_general(N,Nup=Nup,make_basis=False)
>>> s = 17
>>> r = basis.representative(s)
>>> print(s,r)
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)