Documentation¶
Permutational Invariant Quantum Solver (PIQS)
This module calculates the Liouvillian for the dynamics of ensembles of identical two-level systems (TLS) in the presence of local and collective processes by exploiting permutational symmetry and using the Dicke basis.
-
class
piqs.dicke.
Dicke
(N, hamiltonian=None, emission=0.0, dephasing=0.0, pumping=0.0, collective_emission=0.0, collective_dephasing=0.0, collective_pumping=0.0) Bases:
object
The Dicke class which builds the Lindbladian and Liouvillian matrix.
Example
>>> from piqs import Dicke, jspin >>> N = 2 >>> jx, jy, jz = jspin(N) >>> jp = jspin(N, "+") >>> jm = jspin(N, "-") >>> ensemble = Dicke(N, emission=1.) >>> L = ensemble.liouvillian()
Parameters: - N (int) -- The number of two-level systems.
- hamiltonian --
A Hamiltonian in the Dicke basis.
The matrix dimensions are (nds, nds), with nds being the number of Dicke states. The Hamiltonian can be built with the operators given by the jspin functions.
- emission (float) -- Incoherent emission coefficient (also nonradiative emission). default: 0.0
- dephasing (float) -- Local dephasing coefficient. default: 0.0
- pumping (float) -- Incoherent pumping coefficient. default: 0.0
- collective_emission (float) -- Collective (superradiant) emmission coefficient. default: 0.0
- collective_pumping (float) -- Collective pumping coefficient. default: 0.0
- collective_dephasing (float) -- Collective dephasing coefficient. default: 0.0
-
N
¶ int -- The number of two-level systems.
-
hamiltonian
¶ :class: qutip.Qobj -- A Hamiltonian in the Dicke basis.
The matrix dimensions are (nds, nds), with nds being the number of Dicke states. The Hamiltonian can be built with the operators given by the jspin function in the "dicke" basis.
-
emission
¶ float -- Incoherent emission coefficient (also nonradiative emission). default: 0.0
-
dephasing
¶ float -- Local dephasing coefficient. default: 0.0
-
pumping
¶ float -- Incoherent pumping coefficient. default: 0.0
-
collective_emission
¶ float -- Collective (superradiant) emmission coefficient. default: 0.0
-
collective_dephasing
¶ float -- Collective dephasing coefficient. default: 0.0
-
collective_pumping
¶ float -- Collective pumping coefficient. default: 0.0
-
nds
¶ int -- The number of Dicke states.
-
dshape
¶ tuple -- The shape of the Hilbert space in the Dicke or uncoupled basis. default: (nds, nds).
-
__repr__
() Print the current parameters of the system.
-
c_ops
()¶ Build collapse operators in the full Hilbert space 2^N.
Returns: c_ops_list -- The list with the collapse operators in the 2^N Hilbert space. Return type: list
-
coefficient_matrix
()¶ Build coefficient matrix for ODE for a diagonal problem.
Returns: M -- The matrix M of the coefficients for the ODE dp/dt = M p. p is the vector of the diagonal matrix elements of the density matrix rho in the Dicke basis. Return type: ndarray
-
lindbladian
()¶ Build the Lindbladian superoperator of the dissipative dynamics.
Returns: lindbladian -- The Lindbladian matrix as a qutip.Qobj. Return type: class: qutip.Qobj
-
liouvillian
()¶ Build the total Liouvillian using the Dicke basis.
Returns: liouv -- The Liouvillian matrix for the system. Return type: class: qutip.Qobj
-
pisolve
(initial_state, tlist, options=None)¶ Solve for diagonal Hamiltonians and initial states faster.
Parameters: - initial_state -- An initial state specified as a density matrix of qutip.Qbj type
- tlist (ndarray) -- A 1D numpy array of list of timesteps to integrate
- options -- The options for the solver.
Returns: result -- A dictionary of the type qutip.solver.Result which holds the results of the evolution.
Return type:
-
prune_eigenstates
(liouvillian)¶ Remove spurious eigenvalues and eigenvectors of the Liouvillian.
Spurious means that the given eigenvector has elements outside of the block-diagonal matrix.
Parameters: liouvillian_eigenstates (list) -- A list with the eigenvalues and eigenvectors of the Liouvillian including spurious ones. Returns: correct_eigenstates -- The list with the correct eigenvalues and eigenvectors of the Liouvillian. Return type: list
-
class
piqs.dicke.
Pim
(N, emission=0.0, dephasing=0, pumping=0, collective_emission=0, collective_pumping=0, collective_dephasing=0)¶ Bases:
object
The Permutation Invariant Matrix class.
Initialize the class with the parameters for generating a Permutation Invariant matrix which evolves a given diagonal initial state p as:
dp/dt = MpParameters: - N (int) -- The number of two-level systems.
- emission (float) -- Incoherent emission coefficient (also nonradiative emission). default: 0.0
- dephasing (float) -- Local dephasing coefficient. default: 0.0
- pumping (float) -- Incoherent pumping coefficient. default: 0.0
- collective_emission (float) -- Collective (superradiant) emmission coefficient. default: 0.0
- collective_pumping (float) -- Collective pumping coefficient. default: 0.0
- collective_dephasing (float) -- Collective dephasing coefficient. default: 0.0
-
N
¶ int -- The number of two-level systems.
-
emission
¶ float -- Incoherent emission coefficient (also nonradiative emission). default: 0.0
-
dephasing
¶ float -- Local dephasing coefficient. default: 0.0
-
pumping
¶ float -- Incoherent pumping coefficient. default: 0.0
-
collective_emission
¶ float -- Collective (superradiant) emmission coefficient. default: 0.0
-
collective_dephasing
¶ float -- Collective dephasing coefficient. default: 0.0
-
collective_pumping
¶ float -- Collective pumping coefficient. default: 0.0
-
M
¶ dict -- A nested dictionary of the structure {row: {col: val}} which holds non zero elements of the matrix M
-
calculate_j_m
(dicke_row, dicke_col)¶ Get the value of j and m for the particular Dicke space element.
Parameters: Returns: j, m -- The j and m values.
Return type:
-
calculate_k
(dicke_row, dicke_col)¶ Get k value from the current row and column element in the Dicke space.
Parameters: Returns: k -- The row index for the matrix M for given Dicke space element
Return type:
-
coefficient_matrix
()¶ Generate the matrix M governing the dynamics for diagonal cases.
If the initial density matrix and the Hamiltonian is diagonal, the evolution of the system is given by the simple ODE: dp/dt = Mp.
-
isdicke
(dicke_row, dicke_col)¶ Check if an element in a matrix is a valid element in the Dicke space. Dicke row: j value index. Dicke column: m value index. The function returns True if the element exists in the Dicke space and False otherwise.
Parameters:
-
solve
(rho0, tlist, options=None)¶ Solve the ODE for the evolution of diagonal states and Hamiltonians.
-
tau1
(j, m)¶ Calculate the element of the coefficient matrix relative to p_jmm.
-
tau2
(j, m)¶ Calculate the element of the coefficient matrix relative to p_jm+1m+1.
-
tau3
(j, m)¶ Calculate the element of the coefficient matrix relative to p_j+1m+1m+1.
-
tau4
(j, m)¶ Calculate the element of the coefficient matrix relative to p_j-1m+1m+1.
-
tau5
(j, m)¶ Calculate the element of the coefficient matrix relative to p_j+1mm.
-
tau6
(j, m)¶ Calculate the element of the coefficient matrix relative to p_j-1mm.
-
tau7
(j, m)¶ Calculate the element of the coefficient matrix relative to p_j+1m-1m-1.
-
tau8
(j, m)¶ Calculate the element of the coefficient matrix relative to p_jm-1m-1.
-
tau9
(j, m)¶ Calculate the element of the coefficient matrix relative to p_j-1m-1m-1.
-
tau_valid
(dicke_row, dicke_col)¶ Find the Tau functions which are valid for this value of (dicke_row, dicke_col) given the number of TLS. This calculates the valid tau values and reurns a dictionary specifying the tau function name and the value.
Parameters: Returns: taus -- A dictionary of key, val as {tau: value} consisting of the valid taus for this row and column of the Dicke space element.
Return type:
-
piqs.dicke.
am
(j, m)¶ Calculate the coefficient am by applying J_- |j, m>.
The action of am is given by: \(J_{-}|j, m\rangle = A_{-}(j, m)|j, m-1\rangle\)
Parameters: Returns: a_minus -- The value of \(a_{-}\).
Return type:
-
piqs.dicke.
ap
(j, m)¶ Calculate the coefficient ap by applying J_+ |j, m>.
The action of ap is given by: \(J_{+}|j, m\rangle = A_{+}(j, m)|j, m+1\rangle\)
Parameters: Returns: a_plus -- The value of \(a_{+}\).
Return type:
-
piqs.dicke.
block_matrix
(N)¶ Construct the block-diagonal matrix for the Dicke basis.
Parameters: N (int) -- Number of two-level systems. Returns: block_matr -- A 2D block-diagonal matrix of ones with dimension (nds,nds), where nds is the number of Dicke states for N two-level systems. Return type: ndarray
-
piqs.dicke.
collapse_uncoupled
(N, emission=0.0, dephasing=0.0, pumping=0.0, collective_emission=0.0, collective_dephasing=0.0, collective_pumping=0.0)¶ Create the collapse operators (c_ops) of the Lindbladian in the uncoupled basis.
These operators are in the uncoupled basis of the two-level system (TLS) SU(2) Pauli matrices.
Notes
The collapse operator list can be given to qutip.mesolve. Notice that the operators are placed in a Hilbert space of dimension \(2^N\) Thus the method is suitable only for small N (of the order of 10).
Parameters: - N (int) -- The number of two-level systems.
- emission (float) -- Incoherent emission coefficient (also nonradiative emission). default: 0.0
- dephasing (float) -- Local dephasing coefficient. default: 0.0
- pumping (float) -- Incoherent pumping coefficient. default: 0.0
- collective_emission (float) -- Collective (superradiant) emmission coefficient. default: 0.0
- collective_pumping (float) -- Collective pumping coefficient. default: 0.0
- collective_dephasing (float) -- Collective dephasing coefficient. default: 0.0
Returns: c_ops -- The list of collapse operators as qutip.Qobj for the system.
Return type:
-
piqs.dicke.
css
(N, x=0.7071067811865475, y=0.7071067811865475, basis='dicke', coordinates='cartesian')¶ Generate the density matrix of the Coherent Spin State (CSS).
It can be defined as, \(|CSS \rangle = \prod_i^N(a|1\rangle_i + b|0\rangle_i)\) with \(a = sin(\frac{\theta}{2})\), \(b = e^{i \phi}\cos(\frac{\theta}{2})\). The default basis is that of Dicke space \(|j, m\rangle \langle j, m'|\). The default state is the symmetric CSS, \(|CSS\rangle = |+\rangle\).
Parameters: Returns: rho -- The CSS state density matrix.
Return type: class: qutip.Qobj
-
piqs.dicke.
dicke
(N, j, m)¶ Generate a Dicke state as a pure density matrix in the Dicke basis.
For instance, the superradiant state given by \(|j, m\rangle = |1, 0\rangle\) for N = 2, and the state is represented as a density matrix of size (nds, nds) or (4, 4), with the (1, 1) element set to 1.
Parameters: Returns: rho -- The density matrix.
Return type: class: qutip.Qobj
-
piqs.dicke.
dicke_basis
(N, jmm1=None)¶ Initialize the density matrix of a Dicke state for several (j, m, m1).
This function can be used to build arbitrary states in the Dicke basis \(|j, m\rangle \langle j, m^{\prime}|\). We create coefficients for each (j, m, m1) value in the dictionary jmm1. The mapping for the (i, k) index of the density matrix to the |j, m> values is given by the cythonized function jmm1_dictionary. A density matrix is created from the given dictionary of coefficients for each (j, m, m1).
Parameters: Returns: rho -- The density matrix in the Dicke basis.
Return type: class: qutip.Qobj
-
piqs.dicke.
energy_degeneracy
(N, m)¶ Calculate the number of Dicke states with same energy.
The use of the Decimals class allows to explore N > 1000, unlike the built-in function scipy.special.binom
Parameters: Returns: degeneracy -- The energy degeneracy
Return type:
-
piqs.dicke.
excited
(N, basis='dicke')¶ Generate the density matrix for the excited state.
This state is given by (N/2, N/2) in the default Dicke basis. If the argument basis is "uncoupled" then it generates the state in a \(2^N\) dim Hilbert space.
Parameters: Returns: state -- The excited state density matrix in the requested basis.
Return type: class: qutip.Qobj
-
piqs.dicke.
ghz
(N, basis='dicke')¶ Generate the density matrix of the GHZ state.
If the argument basis is "uncoupled" then it generates the state in a \(2^N\) dim Hilbert space.
Parameters: Returns: state -- The GHZ state density matrix in the requested basis.
Return type: class: qutip.Qobj
-
piqs.dicke.
ground
(N, basis='dicke')¶ Generate the density matrix of the ground state.
This state is given by (N/2, -N/2) in the Dicke basis. If the argument basis is "uncoupled" then it generates the state in a \(2^N\) dim Hilbert space.
Parameters: Returns: state -- The ground state density matrix in the requested basis.
Return type: class: qutip.Qobj
-
piqs.dicke.
identity_uncoupled
(N)¶ Generate the identity in a \(2^N\) dimensional Hilbert space.
The identity matrix is formed from the tensor product of N TLSs.
Parameters: N (int) -- The number of two-level systems. Returns: identity -- The identity matrix. Return type: class: qutip.Qobj
-
piqs.dicke.
isdiagonal
(mat)¶ Check if the input matrix is diagonal.
Parameters: mat (ndarray/Qobj) -- A 2D numpy array Returns: diag -- True/False depending on whether the input matrix is diagonal. Return type: bool
-
piqs.dicke.
jspin
(N, op=None, basis='dicke')¶ Calculate the list of collective operators of the total algebra.
The Dicke basis \(|j,m\rangle\langle j,m'|\) is used by default. Otherwise with "uncoupled" the operators are in a \(2^N\) space.
Parameters: Returns: j_alg -- A list of qutip.Qobj representing all the operators in the "dicke" or "uncoupled" basis or a single operator requested.
Return type: list or :class: qutip.Qobj
-
piqs.dicke.
m_degeneracy
(N, m)¶ Calculate the number of Dicke states \(|j, m\rangle\) with same energy.
Parameters: Returns: degeneracy -- The m-degeneracy.
Return type:
-
piqs.dicke.
num_dicke_ladders
(N)¶ Calculate the total number of ladders in the Dicke space.
For a collection of N two-level systems it counts how many different "j" exist or the number of blocks in the block-diagonal matrix.
Parameters: N (int) -- The number of two-level systems. Returns: Nj -- The number of Dicke ladders. Return type: int
-
piqs.dicke.
num_dicke_states
(N)¶ Calculate the number of Dicke states.
Parameters: N (int) -- The number of two-level systems. Returns: nds -- The number of Dicke states. Return type: int
-
piqs.dicke.
num_tls
(nds)¶ Calculate the number of two-level systems.
Parameters: nds (int) -- The number of Dicke states. Returns: N -- The number of two-level systems. Return type: int
-
piqs.dicke.
spin_algebra
(N, op=None)¶ Create the list [sx, sy, sz] with the spin operators.
The operators are constructed for a collection of N two-level systems (TLSs). Each element of the list, i.e., sx, is a vector of qutip.Qobj objects (spin matrices), as it cointains the list of the SU(2) Pauli matrices for the N TLSs. Each TLS operator sx[i], with i = 0, ..., (N-1), is placed in a \(2^N\)-dimensional Hilbert space.
Notes
sx[i] is \(\frac{\sigma_x}{2}\) in the composite Hilbert space.
Parameters: N (int) -- The number of two-level systems. Returns: spin_operators -- A list of qutip.Qobj operators - [sx, sy, sz] or the requested operator. Return type: list or :class: qutip.Qobj
-
piqs.dicke.
state_degeneracy
(N, j)¶ Calculate the degeneracy of the Dicke state.
Each state \(|j, m\rangle\) includes D(N,j) irreducible representations \(|j, m, \alpha\rangle\).
Uses Decimals to calculate higher numerator and denominators numbers.
Parameters: Returns: degeneracy -- The state degeneracy.
Return type:
-
piqs.dicke.
superradiant
(N, basis='dicke')¶ Generate the density matrix of the superradiant state.
This state is given by (N/2, 0) or (N/2, 0.5) in the Dicke basis. If the argument basis is "uncoupled" then it generates the state in a 2**N dim Hilbert space.
Parameters: Returns: state -- The superradiant state density matrix in the requested basis.
Return type: class: qutip.Qobj
-
piqs.dicke.
tau_column
(tau, k, j)¶ Determine the column index for the non-zero elements of the matrix for a particular row k and the value of j from the Dicke space.
Parameters:
Citation generator¶
Citation generator for PIQS