This module implements functions that return unitaries for transforming
density matrices between different bases. It currently includes transformations
related to the generalized Gell-Mann basis and the element-wise (computational
or canonical) basis, with potential for future extensions to other transformations.
Returns the unitary matrix that transforms from the element-wise
(computational) basis to the coefficients in the generalized
Gell-Mann normalized basis.
>>> frombase_change_unitariesimportel_gm>>> dim=3>>> U_el_gm=gm_el(dim)>>> # Define a density matrix in the computational basis>>> rho=np.diag([0.6,0.3,0.1])>>> # Transform the density matrix to the Gell-Mann basis>>> rho_gm=U_el_gm@rho.reshape(dim**2)>>> print("rho in Gell-Mann basis:")>>> print(rho_gm)[0.38723499+0.j 0. +0.j 0.21213203+0.j 0. +0.j 0.38723499+0.j 0. +0.j 0.21213203+0.j 0. +0.j 0.2647605 +0.j]
Returns the unitary matrix that transforms from the coefficients in the
generalized Gell-Mann normalized basis to the element-wise (computational) basis.
Gell-Mann matrices are a generalization of the Pauli matrices used to
describe higher-dimensional systems in quantum mechanics [1]. They form a
complete, orthogonal basis for the space of Hermitian matrices.
The d>3 generalization is taken as in [2] (Eq. 8).
References
Examples
>>> frombase_change_unitariesimportgm_el>>> U_el_gm=el_gm(2)# Transforms from the normalized Pauli basis>>> # Define a density matrix in the Pauli basis>>> rho_pauli=np.array([1/2,1/2,0,0])>>> # Transform the density matrix to the element-wise basis>>> rho=(U_el_gm@rho_pauli).reshape((2,2))>>> print("rho in element-wise basis:")>>> print(rho)[[ 0.5 0.5 ] [ 0.5 0.5 ]]
This module implements functions that return the transition matrices or
Kraus representation describing different families of quantum channels.
These functions cover a variety of channels including classical permutations,
dephasing, depolarizing, embedding classical channels, initializers,
POVMs, probabilistic damping, and probabilistic unitaries.
It describes the coherent damping from state |x> to |y>.
Reduces to the usual qubit amplitude damping for dim=0, x=1, y=0.
Parameters:
dim (int) – The dimension of the Hilbert space.
lamb (float) – The damping parameter ranging from 0 to 1.
x (np.ndarray | int, optional) – The origin state for the damping, by default 1. If int, it is interpreted as a basis state.
y (np.ndarray | int, optional) – The target state for the damping, by default 0. If int, it is interpreted as a basis state.
kraus_representation (bool, optional) – If true, the function returns the Kraus representation of the channel with shape
(Kraus operator index, dim2, dim1), defaults to False.
Returns:
The transition matrix or Kraus operators representing the amplitude damping channel.
Return type:
np.ndarray
Notes
For the general d-dimensional formulation we considered Definition 1 in [3].
References
Raises:
ValueError – If lambda not in [0, 1]; x or y are int and are not in [0, dim), or are np.ndarray and shape not compatible with (dim,1).
Returns the transition matrix representing a classical permutation
channel.
Parameters:
dim (int) – The dimension of the Hilbert space.
perm (list) – The array defining the permutation that maps i to perm[i].
kraus_representation (bool, optional) – Defaults to False. WARNING: This map does not admit a Kraus representation, if it
is set to True it will rise a ValueError.
Returns:
The transition matrix representing the classical permutation channel.
Returns the transition matrix for a dephasing channel.
Parameters:
dim (int) – The dimension of the Hilbert space.
g (float | np.ndarray) – Dephasing strength. If float, it is used as an uniform damping for all
off-diagonal terms; if np.ndarray, g[i,j] describes the damping of the
element (i,j).
The validity of g to define an appropriate dephasing channel is not verified.
u (np.ndarray, optional) – The unitary matrix defining the basis in which dephasing occurs.
Defaults to the identity matrix.
kraus_representation (bool, optional) – If true, the function returns the Kraus representation of the channel with shape
(Kraus operator index, dim2, dim1), defaults to False.
Returns:
The transition matrix or Kraus operators for the dephasing channel.
Return type:
np.ndarray
Notes
This function implements dephasing by dampening the off-diagonal terms
in a given basis. If g is a matrix, it should be real, symmetric, positive
semi-definite and all diagonal terms equal to one.
For other generalized dephasing, one could define the corresponding
PVMs and include an identity with some dampening, then pass them to
the povm function.
Returns the transition matrix for a depolarizing channel.
Parameters:
dim (int) – The dimension of the Hilbert space.
p (float) – Depolarizing probability.
r (np.ndarray, optional) – The stationary state of the depolarizing channel. Defaults to the
completely mixed state.
kraus_representation (bool, optional) – If true, the function returns the Kraus representation of the channel with shape
(Kraus operator index, dim2, dim1), defaults to False.
kraus_atol (float, optional) – The tolerance to neglect a Kraus operator, defaults to 1e-7.
Returns:
The transition matrix or Kraus operators for the depolarizing channel.
Embeds a classical stochastic matrix into a quantum transition matrix.
The resulting process sets off-diagonal terms to zero and acts on the diagonal
terms following the classical stochastic process.
Parameters:
dim (int) – The dimension of the input Hilbert space.
stoch_mat (np.ndarray) – The classical column-stochastic matrix to be embedded.
kraus_representation (bool, optional) – If true, the function returns the Kraus representation of the channel with shape
(Kraus operator index, dim2, dim1), defaults to False.
kraus_atol (float, optional) – The tolerance to neglect a Kraus operator, defaults to 1e-7.
Returns:
The transition matrix or Kraus operators representing the embedded classical channel.
Returns a transition matrix for initializing a quantum state.
The initializer channel takes a classical probability distribution as
input and produces a density matrix.
The mode argument defines the format of the input and output states.
The initial state can be a classical array (c) or a density matrix whose
diagonal elements are used as the input probability distribution (q),
discarding off-diagonal terms. The output state can store only the
initialized states (q) or also the input probability distribution (qc)
Parameters:
dim (int) – The dimension of the output Hilbert space.
states (np.ndarray) – Array of states defining the choices for initialization.
The expected shape is either (number of states, dim, dim) or
(number of states, dim) , such that states[i] represents the
state prepared when the classical variable takes the value i (described
as a density matrix or as a vector).
mode (str, optional) – Mode of initialization (‘c-q’ for classical to quantum,
‘c-qc’ for classical to quantum and classical, ‘q-q’ for
quantum to quantum, ‘q-qc’ for quantum to quantum and
classical). Defaults to ‘c-q’.
kraus_representation (bool, optional) – If true, the function returns the Kraus representation of the channel with shape
(Kraus operator index, dim2, dim1), defaults to False.
Some modes do not admit a Kraus representation.
kraus_atol (float, optional) – The tolerance to neglect a Kraus operator, defaults to 1e-7.
Returns:
The transition matrix or Kraus operators for the initializer.
Returns a transition matrix for a quantum channel defined by a
positive operator-valued measure (POVM).
Parameters:
dim (int) – The dimension of the Hilbert space.
pos (np.ndarray) – The POVM elements defining the channel. The expected shape is
(number of positive operators, dim, dim), such that pos[i] is
the ith positive operator.
mode (str, optional) – Mode of the channel (‘q-q’ for quantum to quantum,
‘q-c’ for quantum to classical, ‘q-qc’ for quantum to
quantum and classical). Defaults to ‘q-q’.
kraus_representation (bool, optional) – If true, the function returns the Kraus representation of the channel with shape
(Kraus operator index, dim2, dim1), defaults to False.
Some modes do not admit a Kraus representation.
Returns:
The transition matrix or Kraus operators for the POVM-based quantum channel.
Returns a transition matrix for a probabilistic damping channel.
Parameters:
dim (int) – The dimension of the Hilbert space.
p (float) – Damping probability.
kraus_representation (bool, optional) – If true, the function returns the Kraus representation of the channel with shape
(Kraus operator index, dim2, dim1), defaults to False.
Returns:
The transition matrix or Kraus operators for the probabilistic damping channel.
Returns a transition matrix for a channel that applies unitaries
probabilistically.
Parameters:
dim (int) – The dimension of the Hilbert space.
p_arr (np.ndarray) – Probabilities for each unitary operation.
u_arr (np.ndarray) – Array of unitary matrices. The expected shape is
(number of unitaries, dim, dim), such that u_arr[i] is
the ith possible unitary.
kraus_representation (bool, optional) – If true, the function returns the Kraus representation of the channel with shape
(Kraus operator index, dim2, dim1), defaults to False.
Returns:
The transition matrix or Kraus operators for the probabilistic unitary channel.
Return type:
np.ndarray
Examples
This example demonstrates the implementation of a bit-flip channel.
Returns a transition matrix for the transposition operator. This is positive but
not completely positive, thus it is not a quantum channel.
Parameters:
dim (int) – The dimension of the Hilbert space.
u (np.ndarray or None) – Basis in which the transposition is performed, if None the canonical basis is used.
Defaults to None.
kraus_representation (bool, optional) – Defaults to False. WARNING: This map does not admit a Kraus representation, if it
is set to True it will rise a ValueError.
This module implements various operations on quantum channels. It includes
functions for obtaining the Choi state of a channel, performing tensor
products of channel representations, and finding fixed points of a channel.
Computes the Doeblin coefficient of a quantum channel.
The Doeblin coefficient is the maximum erasure probability for which an erasure channel
can be degraded into the given channel. This implementation uses the SDP construction
described in reference [1].
The function dependends on the optional package cvxpy. Please install it with
pip install cvxpy before using the function.
Parameters:
channel (np.ndarray) – The transition matrix representing the quantum channel.
transpose (bool, optional) – Whether to check transpose-degradability instead of the usual degradability.
Defaults to False.
subspace_projection (np.ndarray, optional) – An orthogonal projector into a subspace of dimension at least 2. If provided,
it restricts the search for the coefficient to that subspace. Defaults to None.
Returns:
The Doeblin coefficient of the quantum channel. Returns None if the problem is
infeasible or unbounded.
Return type:
float
Notes
This function implements the SDP (Semidefinite Programming) construction given
in reference [1] to find the Doeblin coefficient.
Returns the Kraus operators of a given quantum channel.
Parameters:
t (np.ndarray) – The transition matrix representing the quantum channel,
with shape (dim2**2, dim1**2)
atol (float, optional) – The absolute tolerance for the norm of a Kraus operator to be
considered as a non-trivial Kraus operator, defaults to 1e-6.
Returns:
The Kraus operators corresponding to the quantum channel,
with shape (kraus_rank, dim2, dim1).
Return type:
np.ndarray
Notes
Currently, the method uses the spectral decomposition of the Choi matrix and discards
eigenvalues that are smaller than atol. This is equivalent to discarding Kraus operators
with np.trace(K.T @ K) <= atol * dim1.
Using the matrix square root of the Choi matrix is a valid alternative if
the canonical form is not needed. It seems more efficient in cases where the
matrix is sparse.
Example 1: Tensor product of a list of matrices
Demonstrates how to construct an inhomogeneous local depolarizer that
acts on two qubits, but only adds noise to the state of the first qubit.
Generates a random quantum channel by generating a random isometry
(truncated Haar-random unitary) and tracing out the environment.
Parameters:
dim_in (int) – Input dimension of the quantum channel.
dim_out (int, optional) – Output dimension of the quantum channel. If None, defaults to dim_in.
kraus_rank (int, optional) – Rank of the Kraus operators representing the quantum channel. If None, defaults to dim_in * dim_out.
random_state (int | np.random.Generator, optional) – Used for drawing random variates. If seed is None, the RandomState singleton is used.
If seed is an int, a new RandomState instance is used, seeded with seed. If seed is
already a RandomState or Generator instance, then that object is used. Default is None.
Returns:
The transition matrix representing the random quantum channel.
Generates a random quantum state by sampling Haar random
pure states on dimension dim * rank and computing the
partiual trace on rank degrees of freedom.
Parameters:
dim (int) – Dimension of the quantum state.
rank (int, optional) – Rank of the state. If None, defaults to dim.
random_state (int | np.random.Generator, optional) – Used for drawing random variates. If seed is None, the RandomState singleton is used.
If seed is an int, a new RandomState instance is used, seeded with seed. If seed is
already a RandomState or Generator instance, then that object is used. Default is None.
Generates random quantum states by sampling their spectrum from the homogeneous Dirichlet distribution
with concentration parameter alpha. The basis are rotated with Haar random unitary matrices.
Parameters:
dim (int) – Dimension of the quantum state.
alpha (float) – Concentration parameter of the Dirichlet distribution.
random_state (int | np.random.Generator, optional) – Used for drawing random variates. If seed is None, the RandomState singleton is used.
If seed is an int, a new RandomState instance is used, seeded with seed. If seed is
already a RandomState or Generator instance, then that object is used. Default is None.
Returns:
The random quantum state.
Return type:
np.ndarray
Note
As alpha -> 0 generated states are more pure, alpha = 1 is
equivalent to a Hilbert-Schmidt uniform measure, and as
alpha -> infty it concentrates around the maximally mixed
state.
random_state (int | np.random.Generator, optional) – Used for drawing random variates. If seed is None, the RandomState singleton is used.
If seed is an int, a new RandomState instance is used, seeded with seed. If seed is
already a RandomState or Generator instance, then that object is used. Default is None.
Returns:
The transition matrix representing the random unitary quantum channel.
This module contains functions to directly transform density matrices. Such as,
subsystem permutations or applying a quantum channel in a subsystem without
having to build the transfer matrix for the whole system.
Applies the channel to the state in the chosen sites.
Parameters:
state (np.ndarray) – Density matrix of the state.
system (tuple[int]) – The system structure represented by a tuple with the local
dimension of the constituent subsystems.
active_sites (tuple[int] | int) – The sites in which the channel acts on. If the channel’s input and output dimensions
are equal, the order of the sites is used to arrange the sites before applying the channel
and recover them back. If the channel’s output dimension differs from the input active_sites collapse
to a single site in the relative position corresponding to active_sites[0].
channel (np.ndarray) – The transition matrix of the quantum channel.
Returns:
Density matrix of the state resulting from applying the channel.
Generates a low rank representation of a density matrix with shape (effective rank, dim), or the reverse transformation.
Parameters:
state (np.ndarray) – Density matrix representation of the input state, or a low rank representation if reverse is True.
rank_atol (float, optional) – Absolute tolerance of accumulated probability in the discarded eigenvalues, by default 1e-7. It is overwritten by output_rank.
output_rank (int, optional) – Set the rank of the representation to a fixed value, by default None and rank_atol is used.
reverse (bool, optional) – If true, implements the reverse transformation, by default False.
Returns:
Low rank representation of the state, or the density matrix representation if reverse is True.