Permutational Invariance
- class 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,
The Dicke class which builds the Lindbladian and Liouvillian matrix.
- Parameters:
- N: int
The number of two-level systems.
- hamiltonian
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 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
- Attributes:
- N: int
The number of two-level systems.
- hamiltonian
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).
Examples
>>> from qutip.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()
- c_ops()[source]
Build collapse operators in the full Hilbert space 2^N.
- Returns:
- c_ops_list: list
The list with the collapse operators in the 2^N Hilbert space.
- coefficient_matrix()[source]
Build coefficient matrix for ODE for a diagonal problem.
- Returns:
- M: ndarray
The matrix M of the coefficients for the ODE dp/dt = Mp. p is the vector of the diagonal matrix elements of the density matrix rho in the Dicke basis.
- lindbladian()[source]
Build the Lindbladian superoperator of the dissipative dynamics.
- Returns:
- lindbladian
Qobj
The Lindbladian matrix as a qutip.Qobj.
- lindbladian
- liouvillian()[source]
Build the total Liouvillian using the Dicke basis.
- Returns:
- liouv
Qobj
The Liouvillian matrix for the system.
- liouv
- pisolve(initial_state, tlist)[source]
Solve for diagonal Hamiltonians and initial states faster.
- Parameters:
- initial_state
Qobj
An initial state specified as a density matrix of qutip.Qbj type.
- tlist: ndarray
A 1D numpy array of list of timesteps to integrate
- initial_state
- Returns:
- result: list
A dictionary of the type qutip.piqs.Result which holds the results of the evolution.
- class Pim(
- N,
- emission=0.0,
- dephasing=0,
- pumping=0,
- collective_emission=0,
- collective_pumping=0,
- collective_dephasing=0,
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 = Mp
- 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
- Attributes:
- 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)[source]
Get the value of j and m for the particular Dicke space element.
- Parameters:
- dicke_row, dicke_col: int
The row and column from the Dicke space matrix
- Returns:
- j, m: float
The j and m values.
- calculate_k(dicke_row, dicke_col)[source]
Get k value from the current row and column element in the Dicke space.
- Parameters:
- dicke_row, dicke_col: int
The row and column from the Dicke space matrix.
- Returns
- ——-
- k: int
The row index for the matrix M for given Dicke space element.
- coefficient_matrix()[source]
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)[source]
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:
- dicke_row, dicke_colint
Index of the element in Dicke space which needs to be checked
- tau_valid(dicke_row, dicke_col)[source]
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:
- dicke_row, dicke_colint
Index of the element in Dicke space which needs to be checked.
- Returns:
- taus: dict
A dictionary of key, val as {tau: value} consisting of the valid taus for this row and column of the Dicke space element.
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. It also allows to characterize nonlinear functions of the density matrix.
- am(j, m)[source]
Calculate the operator
am
used later.The action of
ap
is given by: \(J_{-}\lvert j,m\rangle = A_{-}(jm)\lvert j,m-1\rangle\)- Parameters:
- j: float
The value for j.
- m: float
The value for m.
- Returns:
- a_minus: float
The value of \(a_{-}\).
- ap(j, m)[source]
Calculate the coefficient
ap
by applying \(J_+\lvert j,m\rangle\).The action of ap is given by: \(J_{+}\lvert j, m\rangle = A_{+}(j, m) \lvert j, m+1\rangle\)
- Parameters:
- j, m: float
The value for j and m in the dicke basis \(\lvert j, m\rangle\).
- Returns:
- a_plus: float
The value of \(a_{+}\).
- block_matrix(N, elements='ones')[source]
Construct the block-diagonal matrix for the Dicke basis.
- Parameters:
- Nint
Number of two-level systems.
- elementsstr {‘ones’, ‘degeneracy’}, default: ‘ones’
- Returns:
- block_matrndarray
A 2D block-diagonal matrix with dimension (nds,nds), where nds is the number of Dicke states for N two-level systems. Filled with ones or the value of degeneracy at each matrix element.
- 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.
- Parameters:
- N: int
The number of two-level systems.
- emission: float, default: 0.0
Incoherent emission coefficient (also nonradiative emission).
- dephasing: float, default: 0.0
Local dephasing coefficient.
- pumping: float, default: 0.0
Incoherent pumping coefficient.
- collective_emission: float, default: 0.0
Collective (superradiant) emmission coefficient.
- collective_pumping: float, default: 0.0
Collective pumping coefficient.
- collective_dephasing: float, default: 0.0
Collective dephasing coefficient.
- Returns:
- c_ops: list
The list of collapse operators as
Qobj
for the system.
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).
- 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, \(\lvert CSS\rangle = \prod_i^N(a\lvert1\rangle_i+b\lvert0\rangle_i)\) with \(a = sin(\frac{\theta}{2})\), \(b = e^{i \phi}\cos(\frac{\theta}{2})\). The default basis is that of Dicke space \(\lvert j, m\rangle \langle j, m'\rvert\). The default state is the symmetric CSS, \(\lvert CSS\rangle = \lvert+\rangle\).
- Parameters:
- N: int
The number of two-level systems.
- x, y: float, default: sqrt(1/2)
The coefficients of the CSS state.
- basis: str {“dicke”, “uncoupled”}, default: “dicke”
The basis to use.
- coordinates: str {“cartesian”, “polar”}, default: “cartesian”
If polar then the coefficients are constructed as \(sin(x/2), cos(x/2)e^(iy)`\).
- Returns:
- rho:
Qobj
The CSS state density matrix.
- rho:
- dicke(N, j, m)[source]
Generate a Dicke state as a pure density matrix in the Dicke basis.
For instance, the superradiant state given by \(\lvert j, m\rangle = \lvert 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:
- N: int
The number of two-level systems.
- j: float
The eigenvalue j of the Dicke state (j, m).
- m: float
The eigenvalue m of the Dicke state (j, m).
- Returns:
- rho:
Qobj
The density matrix.
- rho:
- dicke_basis(N, jmm1)[source]
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 \(\lvert j, m\rangle\langle j, m'\rvert\). 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 \(\lvert j, m\rangle\) 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:
- N: int
The number of two-level systems.
- jmm1: dict
A dictionary of {(j, m, m1): p} that gives a density p for the (j, m, m1) matrix element.
- Returns:
- rho:
Qobj
The density matrix in the Dicke basis.
- rho:
- dicke_blocks(rho)[source]
Create the list of blocks for block-diagonal density matrix in the Dicke basis.
- Parameters:
- rho
Qobj
A 2D block-diagonal matrix of ones with dimension (nds,nds), where nds is the number of Dicke states for N two-level systems.
- rho
- Returns:
- square_blocks: list of np.ndarray
Give back the blocks list.
- dicke_blocks_full(rho)[source]
Give the full (2^N-dimensional) list of blocks for a Dicke-basis matrix.
- Parameters:
- rho
Qobj
A 2D block-diagonal matrix of ones with dimension (nds,nds), where nds is the number of Dicke states for N two-level systems.
- rho
- Returns:
- full_blockslist
The list of blocks expanded in the 2^N space for N qubits.
- dicke_function_trace(f, rho)[source]
Calculate the trace of a function on a Dicke density matrix. :param f: A Taylor-expandable function of rho. :type f: function :param rho: A density matrix in the Dicke basis. :type rho:
Qobj
- Returns:
- resfloat
Trace of a nonlinear function on rho.
- energy_degeneracy(N, m)[source]
Calculate the number of Dicke states with same energy.
The use of the
Decimals
class allows to explore N > 1000, unlike the built-in functionscipy.special.binom
.- Parameters:
- N: int
The number of two-level systems.
- m: float
Total spin z-axis projection eigenvalue. This is proportional to the total energy.
- Returns:
- degeneracy: int
The energy degeneracy
- entropy_vn_dicke(rho)[source]
Von Neumann Entropy of a Dicke-basis density matrix.
- Parameters:
- rho
Qobj
A 2D block-diagonal matrix of ones with dimension (nds, nds), where nds is the number of Dicke states for N two-level systems.
- rho
- Returns:
- entropy_dm: float
Entropy. Use degeneracy to multiply each block.
- excited(N, basis='dicke')[source]
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:
- N: int
The number of two-level systems.
- basis: str, {“dicke”, “uncoupled”}, default: “dicke”
The basis to use.
- Returns:
- state:
Qobj
The excited state density matrix in the requested basis.
- state:
- ghz(N, basis='dicke')[source]
Generate the density matrix of the GHZ state.
If the argument
basis
is “uncoupled” then it generates the state in a \(2^N\)-dimensional Hilbert space.- Parameters:
- N: int
The number of two-level systems.
- basis: str, {“dicke”, “uncoupled”}, default: “dicke”
The basis to use.
- Returns:
- state:
Qobj
The GHZ state density matrix in the requested basis.
- state:
- ground(N, basis='dicke')[source]
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\)-dimensional Hilbert space.- Parameters:
- N: int
The number of two-level systems.
- basis: str, {“dicke”, “uncoupled”}, default: “dicke”
The basis to use.
- Returns:
- state:
Qobj
The ground state density matrix in the requested basis.
- state:
- identity_uncoupled(N)[source]
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:
Qobj
The identity matrix.
- identity:
- isdiagonal(mat)[source]
Check if the input matrix is diagonal.
- Parameters:
- mat: ndarray/Qobj
A 2D numpy array
- Returns:
- diag: bool
True/False depending on whether the input matrix is diagonal.
- jspin(N, op=None, basis='dicke')[source]
Calculate the list of collective operators of the total algebra.
The Dicke basis \(\lvert j,m\rangle\langle j,m'\rvert\) is used by default. Otherwise with “uncoupled” the operators are in a \(2^N\) space.
- Parameters:
- N: int
Number of two-level systems.
- op: str {‘x’, ‘y’, ‘z’, ‘+’, ‘-‘}, optional
The operator to return ‘x’,’y’,’z’,’+’,’-‘. If no operator given, then output is the list of operators for [‘x’,’y’,’z’].
- basis: str {“dicke”, “uncoupled”}, default: “dicke”
The basis of the operators.
- Returns:
- j_alg: list or
Qobj
A list of qutip.Qobj representing all the operators in the “dicke” or “uncoupled” basis or a single operator requested.
- j_alg: list or
- m_degeneracy(N, m)[source]
Calculate the number of Dicke states \(\lvert j, m\rangle\) with same energy.
- Parameters:
- N: int
The number of two-level systems.
- m: float
Total spin z-axis projection eigenvalue (proportional to the total energy).
- Returns:
- degeneracy: int
The m-degeneracy.
- num_dicke_ladders(N)[source]
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: int
The number of Dicke ladders.
- num_dicke_states(N)[source]
Calculate the number of Dicke states.
- Parameters:
- N: int
The number of two-level systems.
- Returns:
- nds: int
The number of Dicke states.
- num_tls(nds)[source]
Calculate the number of two-level systems.
- Parameters:
- nds: int
The number of Dicke states.
- Returns:
- N: int
The number of two-level systems.
- purity_dicke(rho)[source]
Calculate purity of a density matrix in the Dicke basis. It accounts for the degenerate blocks in the density matrix.
- Parameters:
- rho
Qobj
Density matrix in the Dicke basis of qutip.piqs.jspin(N), for N spins.
- rho
- Returns:
- purityfloat
The purity of the quantum state. It’s 1 for pure states, 0<=purity<1 for mixed states.
- spin_algebra(N, op=None)[source]
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.
- Parameters:
- N: int
The number of two-level systems.
- Returns:
- spin_operators: list or
Qobj
A list of qutip.Qobj operators - [sx, sy, sz] or the requested operator.
- spin_operators: list or
Notes
sx[i] is \(\frac{\sigma_x}{2}\) in the composite Hilbert space.
- state_degeneracy(N, j)[source]
Calculate the degeneracy of the Dicke state.
Each state \(\lvert j, m\rangle\) includes D(N,j) irreducible representations \(\lvert j, m, \alpha\rangle\).
Uses Decimals to calculate higher numerator and denominators numbers.
- Parameters:
- N: int
The number of two-level systems.
- j: float
Total spin eigenvalue (cooperativity).
- Returns:
- degeneracy: int
The state degeneracy.
- superradiant(N, basis='dicke')[source]
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:
- N: int
The number of two-level systems.
- basis: str, {“dicke”, “uncoupled”}, default: “dicke”
The basis to use.
- Returns:
- state:
Qobj
The superradiant state density matrix in the requested basis.
- state:
- tau_column(tau, k, j)[source]
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:
- tau: str
The tau function to check for this k and j.
- k: int
The row of the matrix M for which the non zero elements have to be calculated.
- j: float
The value of j for this row.