quspin.basis.photon_basis

class quspin.basis.photon_basis(basis_constructor, *constructor_args, **blocks)[source]

Bases: tensor_basis

Constructs basis in photon-particle Hilbert space.

The photon_basis is child class of tensor_basis which allows the user to define a basis which couples a lattice particle to a single photon mode.

There are two types of photon_basis objects that one can create:
  • conserving basis: lattice particle basis and photon number conserved separately. This object is constructed using the tensor_basis class.

  • non-conserving basis: only total sum of lattice particles and photons conserved.

The operator strings for the photon and lattice sectors are the same as for the lattice bases, respectively. The photon_basis operator string uses the pipe character ‘|’ to distinguish between lattice operators (left) and photon operators (right).

\[\begin{array}{cccc} \texttt{basis}/\texttt{opstr} & \texttt{"I"} & \texttt{"+"} & \texttt{"-"} & \texttt{"n"} & \texttt{"z"} & \texttt{"x"} & \texttt{"y"} \newline \texttt{spin_basis_*} & \hat{1} & \hat S^+(\hat\sigma^+) & \hat S^-(\hat\sigma^-) & - & \hat S^z(\hat\sigma^z) & \hat S^x(\hat\sigma^x) & \hat S^y(\hat\sigma^y) \ \newline \texttt{boson_basis_*}& \hat{1} & \hat b^\dagger & \hat b & \hat b^\dagger \hat b & \hat b^\dagger\hat b - \frac{\mathrm{sps}-1}{2} & - & - \newline \texttt{*_fermion_basis_*}& \hat{1} & \hat c^\dagger & \hat c & \hat c^\dagger \hat c & \hat c^\dagger\hat c - \frac{1}{2} & \hat c + \hat c^\dagger & -i\left( \hat c - \hat c^\dagger\right) \newline \end{array}\]

Examples

For the conserving basis, one can specify the total number of quanta using the the Ntot keyword argument:

>>> p_basis = photon_basis(basis_class,*basis_args,Ntot=...,**symmetry_blocks)

For the non-conserving basis, one must specify the total number of photon (a.k.a. harmonic oscillator) states using the Nph argument:

>>> p_basis = photon_basis(basis_class,*basis_args,Nph=...,**symmetry_blocks)

The code snippet below shows how to use the photon_basis class to construct the Jaynes-Cummings Hamiltonian. As an initial state, we choose a coherent state in the photon sector and the ground state of the two-level system (atom).

 1from quspin.operators import hamiltonian  # Hamiltonian and observables
 2from quspin.basis.photon import coherent_state  # HO coherent state
 3import numpy as np  # generic math functions
 4
 5#
 6##### define model parameters #####
 7Nph_tot = 60  # maximum photon occupation
 8Nph = Nph_tot / 2  # mean number of photons in initial coherent state
 9Omega = 3.5  # drive frequency
10A = 0.8  # spin-photon coupling strength (drive amplitude)
11Delta = 1.0  # difference between atom energy levels
12#
13##### set up photon-atom Hamiltonian #####
14# define operator site-coupling lists
15ph_energy = [[Omega]]  # photon energy
16at_energy = [[Delta, 0]]  # atom energy
17absorb = [[A / (2.0 * np.sqrt(Nph)), 0]]  # absorption term
18emit = [[A / (2.0 * np.sqrt(Nph)), 0]]  # emission term
19# define static and dynamics lists
20static = [["|n", ph_energy], ["x|-", absorb], ["x|+", emit], ["z|", at_energy]]
21dynamic = []
22# compute atom-photon basis
23basis = photon_basis(spin_basis_1d, L=1, Nph=Nph_tot)
24print(basis)
25# compute atom-photon Hamiltonian H
26H = hamiltonian(static, dynamic, dtype=np.float64, basis=basis)
27#
28##### define initial state #####
29# define atom ground state
30psi_at_i = np.array([1.0, 0.0])  # spin-down eigenstate of \sigma^z
31# define photon coherent state with mean photon number Nph
32psi_ph_i = coherent_state(np.sqrt(Nph), Nph_tot + 1)
33# compute atom-photon initial state as a tensor product
34psi_i = np.kron(psi_at_i, psi_ph_i)
__init__(basis_constructor, *constructor_args, **blocks)[source]

Initialises the photon_basis object.

Parameters:
basis_constructor: :obj:`basis`

basis constructor for the lattice part of the photon_basis.

constructor_args: obj

Required arguments required by the specific basis constructor.

blocks: obj

Optional keyword arguments for basis_constructor which include (but are not limited to):

Nph (int) - specify the dimension of photon (HO) Hilbert space.

Ntot (int) - specify total number of particles (photons + lattice).

anyblock (int) - specify any lattice symmetry blocks

Methods

Op(opstr, indx, J, dtype)

Constructs operator from a site-coupling list and anoperator string in the photon basis.

__init__(basis_constructor, ...)

Initialises the photon_basis object.

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

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

get_vec(v0[, sparse, Nph, full_part])

DEPRECATED (see 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 photon_basis.

project_from(v0[, sparse, Nph, full_part])

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_particle

lattice basis.

basis_photon

photon basis.

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.

chain_N

number of sites on lattice.

chain_Ns

number of states in the photon Hilbert space only.

operators

set of available operator strings.

particle_N

number of sites on lattice.

particle_Ns

number of states in the lattice Hilbert space only.

particle_sps

number of lattice states per site.

photon_Ns

number of states in the photon Hilbert space only.

photon_sps

number of photon states per site.

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 anoperator string in the photon basis.

Parameters:
opstr: str

Operator string in the photon basis format. Photon operators stand on the right. For instance: >>> opstr = “x|n”

indx: list(int)

List of integers to designate the sites the photon basis operator is defined on. The photon_basis assumes that the single photon couples globally to the lattice, hence the indx requires only the lattice site position. For instance: >>> indx = [5]

J: scalar

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 photon 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 photon basis,

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

Examples

>>> J = 1.41
>>> indx = [5]
>>> opstr = "x|n"
>>> 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_particle

lattice basis.

Type:

basis

property basis_photon

photon basis.

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

property chain_N

number of sites on lattice.

Type:

int

property chain_Ns

number of states in the photon Hilbert space only.

Type:

int

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='particles', 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:
state: obj

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_A: str, optional

Subsystem to calculate the density matrix of. Can be either one of:

  • “particles”: refers to lattice subsystem (Default).

  • “photons”: refers to photon subsystem.

return_rdm: str, 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_pure: bool, 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.

sparse: bool, optional

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

return_rdm_EVs: bool, 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.

alpha: float, optional

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

sparse_diag: bool, optional

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

maxiter: int, 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="photons",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, Nph=None, full_part=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.

Nph: int, optional

Nuber of photons. Required for conserving photon basis.

full_part: bool, optional

Whether or not to transform the state to the full state in the lattice basis in the case of a conserving photon_basis. 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, Nph=None, full_part=True)[source]

DEPRECATED (see 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)

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='particles', return_rdm=None, enforce_pure=False, sparse=False)[source]

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

Parameters:
state: obj

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_A: str, optional

Subsystem to calculate the density matrix of. Can be either one of:

  • “particles”: refers to lattice subsystem (Default).

  • “photons”: refers to photon subsystem.

return_rdm: str, 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_pure: bool, 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.

sparse: bool, 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)
property particle_N

number of sites on lattice.

Type:

int

property particle_Ns

number of states in the lattice Hilbert space only.

Type:

int

property particle_sps

number of lattice states per site.

Type:

int

property photon_Ns

number of states in the photon Hilbert space only.

Type:

int

property photon_sps

number of photon states per site.

Type:

int

project_from(v0, sparse=True, Nph=None, full_part=True)[source]

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

Parameters:
v0: numpy.ndarray

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

sparse: bool, optional

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

Nph: int, optional

Nuber of photons. Required for conserving photon basis.

full_part: bool, optional

Whether or not to transform the state to the full state in the lattice basis in the case of a conserving photon_basis. 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