Permutational Invariance

Permutational Invariant Quantum Solver (PIQS)

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), \(\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.

Consider a system evolving according to the equation

\[ \begin{align}\begin{aligned}\dot{\rho} = \mathcal{D}_\text{TLS}(\rho)=-\frac{i}{\hbar}\lbrack H,\rho \rbrack +\frac{\gamma_\text{CE}}{2}\mathcal{L}_{J_{-}}[\rho] +\frac{\gamma_\text{CD}}{2}\mathcal{L}_{J_{z}}[\rho] +\frac{\gamma_\text{CP}}{2}\mathcal{L}_{J_{+}}[\rho]\\+\sum_{n=1}^{N}\left( \frac{\gamma_\text{E}}{2}\mathcal{L}_{J_{-,n}}[\rho] +\frac{\gamma_\text{D}}{2}\mathcal{L}_{J_{z,n}}[\rho] +\frac{\gamma_\text{P}}{2}\mathcal{L}_{J_{+,n}}[\rho]\right)\end{aligned}\end{align} \]

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) emission

  • collective_dephasing defines \(\gamma_\text{CD}\), collective dephasing

  • collective_pumping defines \(\gamma_\text{CP}\), collective pumping.

  • emission defines \(\gamma_\text{E}\), incoherent emission (losses)

  • dephasing defines \(\gamma_\text{D}\), local dephasing

  • pumping defines \(\gamma_\text{P}\), incoherent pumping.

Then the system.lindbladian() creates the total TLS Lindbladian 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)

For more example of use, see the “Permutational Invariant Lindblad Dynamics” section in the tutorials section of the website, http://qutip.org/tutorials.html.

Useful PIQS functions.

Operators

Command

Description

Collective spin algebra \(J_x,\ J_y,\ J_z\)

jspin(N)

The collective spin algebra \(J_x,\ J_y,\ J_z\) for \(N\) TLSs

Collective spin \(J_x\)

jspin(N, "x")

The collective spin operator \(Jx\). Requires \(N\) number of TLSs

Collective spin \(J_y\)

jspin(N, "y")

The collective spin operator \(J_y\). Requires \(N\) number of TLSs

Collective spin \(J_z\)

jspin(N, "z")

The collective spin operator \(J_z\). 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 \(J_z\) in uncoupled basis

jspin(N, "z", basis='uncoupled')

The collective spin operator \(J_z\) in the uncoupled basis of dimension \(2^N\).

Dicke state \(|j,m\rangle\) density matrix

dicke(N, j, m)

The density matrix for the Dicke state given by \(|j,m\rangle\)

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 qutip.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 API documentation. Attributes to the qutip.piqs.Dicke and qutip.piqs.Pim class can also be found there.