quspin.basis.spinful_fermion_basis_1d

class quspin.basis.spinful_fermion_basis_1d(L, Nf=None, nf=None, double_occupancy=True, **blocks)[source]

Bases: spinless_fermion_basis_1d, basis_1d

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

The supported operator strings for spinful_fermion_basis_1d are:

\[\begin{array}{cccc} \texttt{basis}/\texttt{opstr} & \texttt{"I"} & \texttt{"+"} & \texttt{"-"} & \texttt{"n"} & \texttt{"z"} \newline \texttt{spinful_fermion_basis_1d}& \hat{1} & \hat c^\dagger & \hat c & \hat c^\dagger c & \hat c^\dagger\hat c - \frac{1}{2} \newline \end{array}\]

Notes

  • The spinful_fermion_basis operator strings are separated by a pipe symbol, “|”, to distinguish the spin-up from spin-down species. However, the index array has NO pipe symbol.

  • Particle-hole like symmetries for fermions can be programmed using the spinful_fermion_basis_general class.

Examples

The code snippet below shows how to use the spinful_fermion_basis_1d class to construct the basis in the zero momentum sector of positive fermion spin for the Fermi-Hubbard Hamiltonian

\[H=-J\sum_{j,\sigma} c^\dagger_{j+1\sigma}c_{j\sigma} + \mathrm{h.c.} - \mu\sum_{j,\sigma} n_{j\sigma} + U \sum_j n_{j\uparrow} n_{j\downarrow}\]

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

 1from quspin.basis import spinful_fermion_basis_1d  # Hilbert space spinful fermion basis
 2import numpy as np  # generic math functions
 3
 4#
 5##### define model parameters #####
 6L = 6  # system size
 7J = 1.0  # hopping strength
 8U = np.sqrt(2)  # onsite interaction strength
 9#
10##### construct basis at half-filling in the 0-total momentum and +1-spin flip sector
11basis = spinful_fermion_basis_1d(L=L, Nf=(L // 2, L // 2), a=1, kblock=0, sblock=1)
12print(basis)
13#
14##### define PBC site-coupling lists for operators
15# define site-coupling lists
16hop_right = [[-J, i, (i + 1) % L] for i in range(L)]  # hopping to the right PBC
17hop_left = [[J, i, (i + 1) % L] for i in range(L)]  # hopping to the left PBC
18int_list = [[U, i, i] for i in range(L)]  # onsite interaction
19# static and dynamic lists
20static = [
21    ["+-|", hop_left],  # up hop left
22    ["-+|", hop_right],  # up hop right
23    ["|+-", hop_left],  # down hop left
24    ["|-+", hop_right],  # down hop right
25    ["n|n", int_list],  # onsite interaction
26]
27dynamic = []
28###### construct Hamiltonian
29H = hamiltonian(static, dynamic, dtype=np.float64, basis=basis)
__init__(L, Nf=None, nf=None, double_occupancy=True, **blocks)[source]

Intializes the fermion_basis_1d object (basis for fermionic operators).

Parameters:
L: int

Length of chain/number of sites.

Nf: tuple(int,list), optional

Number of fermions in chain. First (left) entry refers to spin-up and second (right) entry refers to spin-down. Each of the two entries can be integer or list to specify one or more particle sectors.

nf: tuple(float), optional

Density of fermions in chain (fermions per site). First (left) entry refers to spin-up. Second (right) entry refers to spin-down.

double_occupancy: bool, optional

Boolean to toggle the presence of doubly-occupied sites (both a spin up and a spin-down fermion present on the same lattice site) in the basis. Default is double_occupancy=True, for which doubly-occupied states are present.

**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.

sblock (int) - specifies fermion spin inversion block. The physical manifestation of this symmetry transformation is the exchange of a spin-up and a spin-down fermion on a fixed lattice site.

psblock (int) - specifies parity followed by fermion spin inversion symmetry block. The physical manifestation of this symmetry transformation is reflection about the middle of the chain, and a simultaneous exchange of a spin-up and a spin-down fermion on a fixed lattice site.

Methods

Op(opstr, indx, J, dtype)

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

__init__(L[, Nf, nf, double_occupancy])

Intializes the fermion_basis_1d object (basis for fermionic 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(up_state, down_state)

Finds the index of user-defined Fock state in spinful fermion 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, density, ...])

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

Total number of sites (spin-up + spin-down) the basis is constructed with; N=2L.

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

Total number of sites (spin-up + spin-down) the basis is constructed with; N=2L.

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)[source]

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

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

  • dict(‘V_states’,V_states) [shape (Ns,Nvecs)]: collection of Nvecs states stored in the columns of V_states.

sub_sys_Atuple, 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),range(N//2)) with N the number of physical lattice sites (e.g. sites which both species of fermions can occupy). The format is (spin_up_subsys,spin_down_subsys) (see example below).

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 colelction 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\):

\[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).

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

>>> sub_sys_A_up=range(basis.L//2) # subsystem for spin-up fermions
>>> sub_sys_A_down=range(basis.L//2+1) # subsystem for spin-down fermions
>>> subsys_A=(sub_sys_A_up,sub_sys_A_down)
>>> state=1.0/np.sqrt(basis.Ns)*np.ones(basis.Ns) # infinite temperature state
>>> ent_entropy(state,sub_sys_A=subsys_A,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(up_state, down_state)[source]

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

Parameters:
up_statestr

string which define the Fock state for the spin up fermions.

down_statestr

string which define the Fock state for the spin down fermions.

Returns:
int

Position of the Fock state in the spinful_fermion_basis_1d.

Notes

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

Examples

>>> s_up = "".join("1" for i in range(2)) + "".join("0" for i in range(2))
>>> s_down = "".join("0" for i in range(2)) + "".join("1" for i in range(2))
>>> print( basis.index(s_up,s_down) )
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)[source]

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

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

  • dict(‘V_states’,V_states) [shape (Ns,Nvecs)]: collection of Nvecs states stored in the columns of V_states.

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),range(N//2)) with N the number of physical lattice sites (e.g. sites which both species of fermions can occupy). The format is (spin_up_subsys,spin_down_subsys) (see example below).

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.

Examples

>>> sub_sys_A_up=range(basis.L//2) # subsystem for spin-up fermions
>>> sub_sys_A_down=range(basis.L//2+1) # subsystem for spin-down fermions
>>> subsys_A=(sub_sys_A_up,sub_sys_A_down)
>>> state=1.0/np.sqrt(basis.Ns)*np.ones(basis.Ns) # infinite temperature state
>>> partial_trace(state,sub_sys_A=subsys_A,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)[source]

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)