The Permutational Invariant Quantum Solver (PIQS) is a QuTiP module that allows to study the dynamics of an open quantum system consisting of an ensemble of identical qubits that can dissipate through local and collective baths according to a Lindblad master equation.
The Liouvillian of an ensemble of N qubits, or two-level systems (TLSs), 𝒟TLS(ρ), can be built using only polynomial – instead of exponential – resources. This has many applications for the study of realistic quantum optics models of many TLSs and in general as a tool in cavity QED.
Consider a system evolving according to the equation
where
The inclusion of local processes in the dynamics lead to using a Liouvillian space of dimension 4N. By exploiting the permutational invariance of identical particles [2-8], the Liouvillian 𝒟TLS(ρ) can be built as a block-diagonal matrix in the basis of Dicke states |j, m⟩.
The system under study is defined by creating an object of the ~qutip.piqs.piqs.Dicke
class, e.g. simply named system
, whose first attribute is
system.N
, the number of TLSs of the system N.
The rates for collective and local processes are simply defined as
collective_emission
defines γCE, collective (superradiant) emissioncollective_dephasing
defines γCD, collective dephasingcollective_pumping
defines γCP, collective pumping.emission
defines γE, incoherent emission (losses)dephasing
defines γD, local dephasingpumping
defines γP, incoherent pumping.
Then the system.lindbladian()
creates the total TLS Lindbladian superoperator matrix. Similarly, system.hamiltonian
defines the TLS hamiltonian of the system HTLS.
The system's Liouvillian can be built using system.liouvillian()
. The properties of a Piqs object can be visualized by simply calling system
. We give two basic examples on the use of PIQS. In the first example the incoherent emission of N driven TLSs is considered.
from qutip import piqs
N = 10
system = piqs.Dicke(N, emission = 1, pumping = 2)
L = system.liouvillian()
steady = steadystate(L)
For more example of use, see the "Permutational Invariant Lindblad Dynamics" section in the tutorials section of the website, https://qutip.org/tutorials.html.
Useful PIQS functions.Operators | Command | Description |
---|---|---|
Collective spin algebra Jx, Jy, Jz | jspin(N) |
The collective spin algebra Jx, Jy, Jz for N TLSs |
Collective spin Jx | jspin(N, "x") |
The collective spin operator Jx. Requires N number of TLSs |
Collective spin Jy | jspin(N, "y") |
The collective spin operator Jy. Requires N number of TLSs |
Collective spin Jz | jspin(N, "z") |
The collective spin operator Jz. Requires N number of TLSs |
Collective spin J+ | jspin(N, "+") |
The collective spin operator J+. |
Collective spin J− | jspin(N, "-") |
The collective spin operator J−. |
Collective spin Jz in uncoupled basis | jspin(N, "z", basis='uncoupled') |
The collective spin operator Jz in the uncoupled basis of dimension 2N. |
Dicke state |j, m⟩ density matrix | dicke(N, j, m) |
The density matrix for the Dicke state given by |j, m⟩ |
Excited-state density matrix in Dicke basis | excited(N) |
The excited state in the Dicke basis |
Excited-state density matrix in uncoupled basis | excited(N, basis="uncoupled") |
The excited state in the uncoupled basis |
Ground-state density matrix in Dicke basis | ground(N) |
The ground state in the Dicke basis |
GHZ-state density matrix in the Dicke basis | ghz(N) |
The GHZ-state density matrix in the Dicke (default) basis for N number of TLS |
Collapse operators of the ensemble | Dicke.c_ops() |
The collapse operators for the ensemble can be called by the c_ops method of the Dicke class. |
Note that the mathematical object representing the density matrix of the full system that is manipulated (or obtained from steadystate) in the Dicke-basis formalism used here is a representative of the density matrix. This representative object is of linear size N^2, whereas the full density matrix is defined over a 2^N Hilbert space. In order to calculate nonlinear functions of such density matrix, such as the Von Neumann entropy or the purity, it is necessary to take into account the degeneracy of each block of such block-diagonal density matrix. Note that as long as one calculates expected values of operators, being Tr[A*rho] a linear function of rho, the representative density matrix give straightforwardly the correct result. When a nonlinear function of the density matrix needs to be calculated, one needs to weigh each degenerate block correctly; this is taken care by the dicke_function_trace in .piqs
, and the user can use it to define general nonlinear functions that can be described as the trace of a Taylor expandable function. Two nonlinear functions that use dicke_function_trace and are already implemented are purity_dicke, to calculate the purity of a density matrix in the Dicke basis, and entropy_vn_dicke, which can be used to calculate the Von Neumann entropy.
More functions relative to the qutip.piqs
module can be found at apidoc
. Attributes to the .piqs.Dicke
and .piqs.Pim
class can also be found there.