User Guide¶
The Permutational Invariant Quantum Solver (PIQS) is an open-source Python solver to study the exact Lindbladian dynamics of open quantum systems consisting of identical qubits. It is integrated in QuTiP and can be imported as as a model.
Using this library, the Liouvillian of an ensemble of \(N\) qubits, or two-level systems (TLSs), \(\mathcal{D}_{TLS}(\rho)\), 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 [1].
Consider a system evolving according to the equation
where \(J_{\alpha,n}=\frac{1}{2}\sigma_{\alpha,n}\) are SU(2) Pauli spin operators, with \({\alpha=x,y,z}\) and \(J_{\pm,n}=\sigma_{\pm,n}\). The collective spin operators are \(J_{\alpha} = \sum_{n}J_{\alpha,n}\) . The Lindblad super-operators are \(\mathcal{L}_{A} = 2A\rho A^\dagger - A^\dagger A \rho - \rho A^\dagger A\).
The inclusion of local processes in the dynamics lead to using a Liouvillian space of dimension \(4^N\). By exploiting the permutational invariance of identical particles [2-8], the Liouvillian \(\mathcal{D}_\text{TLS}(\rho)\) can be built as a block-diagonal matrix in the basis of Dicke states \(|j, m \rangle\).
The system under study is defined by creating an object of the
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 \(\gamma_\text{CE}\), collective (superradiant) emissioncollective_dephasing
defines \(\gamma_\text{CD}\), collective dephasingcollective_pumping
defines \(\gamma_\text{CP}\), collective pumping.emission
defines \(\gamma_\text{E}\), incoherent emission (losses)dephasing
defines \(\gamma_\text{D}\), local dephasingpumping
defines \(\gamma_\text{P}\), incoherent pumping.
Then the system.lindbladian()
creates the total TLS Linbladian superoperator matrix. Similarly, system.hamiltonian
defines the TLS hamiltonian of the system \(H_\text{TLS}\).
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 piqs import Dicke
from qutip import steadystate
N = 10
system = Dicke(N, emission = 1, pumping = 2)
L = system.liouvillian()
steady = steadystate(L)
Operators |
Command |
Description |
---|---|---|
Collective spin Jx |
|
The collective spin operator Jx. Requires N number of TLS |
Collective spin J+ |
|
The collective spin operator J+. |
Collective spin J- |
|
The collective spin operator Jz. |
Collective spin Jz in uncoupled basis |
|
The collective spin operator Jz in the uncoupled basis |
Dicke state |j, m> |
|
A Dicke state given by |j, m> |
Excited state in uncoupled basis |
|
The excited state in the uncoupled basis |
GHZ state in the Dicke basis |
|
The GHZ state in the Dicke (default) basis for N number of TLS |
Collapse operators of the ensemble |
|
The collapse operators for the ensemble can be called by the c_ops method of the dicke class. |