quimb.calc¶
Functions for more advanced calculations of quantities and properties of quantum objects.
Attributes¶
Calculate the entropy of a pure states' subsystem, optionally switching |
|
Compute the trace sqrt of a subsystem, possibly using an approximate |
|
Decompose an operator into Paulis. |
|
Decompose an operator into bell-states. |
Functions¶
|
Fidelity between two quantum states. By default, the unsquared fidelity |
|
Take state rho and purify it into a wavefunction of squared |
|
Operate on a mixed state with a set of kraus operators: |
|
Compute the projector for the target |
|
Measure state |
|
Simulate measuring each qubit of |
|
Dephase |
|
Compute the (von Neumann) entropy. |
|
Find the mutual information for a bipartition of a state. |
|
Make sure all indices found in the tuples |
|
Calculate the mutual information of two subsystems of a pure state, |
|
Find the schmidt gap of the bipartition of |
|
Return the trace of the sqrt of a positive semidefinite operator. |
|
Partial transpose of a density operator. |
|
Compute the norm of the partial transpose for (log)-negativity, |
|
Compute logarithmic negativity between two subsytems. |
|
Compute the logarithmic negativity between two subsystems of a pure |
|
Compute negativity between two subsytems. |
|
Concurrence of two-qubit state. |
|
One way classical information for two qubit density operator. |
|
Quantum Discord for two qubit density operator. |
|
Trace distance between two states: |
|
Print a state in the computational basis. |
|
Decomposes an operator via the Hilbert-schmidt inner product. |
|
Calculate the correlation between two sites given two operators. |
|
Calculate the correlation between sites for a list of operator pairs |
|
Calculate the pair-wise function ent_fn between all sites or blocks |
|
|
|
Check if operator has any degenerate eigenvalues, determined relative |
|
Determines whether a vector is an eigenvector of an operator. |
|
Calculate the page entropy, i.e. expected entropy for a subsytem |
Get the analytic isotropic heisenberg chain ground energy for length L. |
Module Contents¶
- quimb.calc.fidelity(p1, p2, squared=False)[source]¶
Fidelity between two quantum states. By default, the unsquared fidelity is used, equivalent in the pure state case to
|<psi|phi>|
.
- quimb.calc.purify(rho)[source]¶
Take state rho and purify it into a wavefunction of squared dimension.
- Parameters:
rho (operator) – Density operator to purify.
- Returns:
The purified ket.
- Return type:
vector
- quimb.calc.kraus_op(rho, Ek, dims=None, where=None, check=False)[source]¶
Operate on a mixed state with a set of kraus operators:
\[\sigma = \sum_k E_k \rho E_k^{\dagger}\]- Parameters:
rho ((d, d) density matrix) – The density operator to perform operation on.
Ek ((K, d, d) array or sequence of K (d, d) arrays) – The Kraus operator(s).
dims (sequence of int, optional) – The subdimensions of
rho
.where (int or sequence of int, optional) – Which of the subsystems to apply the operation on.
check (bool, optional) – Where to check
sum_k(Ek.H @ Ek) == 1
.
- Returns:
sigma – The state after the kraus operation.
- Return type:
density matrix
Examples
The depolarising channel:
In [1]: import quimb as qu In [2]: rho = qu.rand_rho(2) In [3]: I, X, Y, Z = (qu.pauli(s) for s in 'IXYZ') In [4]: [qu.expec(rho, A) for A in (X, Y, Z)] Out[4]: [-0.652176185230622, -0.1301762132792484, 0.22362918368272583] In [5]: p = 0.1 In [6]: Ek = [ ...: (1 - p)**0.5 * I, ...: (p / 3)**0.5 * X, ...: (p / 3)**0.5 * Y, ...: (p / 3)**0.5 * Z ...: ] In [7]: sigma = qu.kraus_op(rho, Ek) In [8]: [qu.expec(sigma, A) for A in (X, Y, Z)] Out[8]: [-0.5652193605332058, -0.11281938484201527, 0.1938119591916957]
- quimb.calc.projector(A, eigenvalue=1.0, tol=1e-12, autoblock=False)[source]¶
Compute the projector for the target
eigenvalue
of operatorA
.- Parameters:
A (operator or tuple[array, operator].) – The hermitian observable. If a tuple is supplied, assume this is the eigendecomposition of the observable.
eigenvalue (float, optional) – The target eigenvalue to construct the projector for, default: 1.
tol (float, optional) – The tolerance with which to select all eigenvalues.
autoblock (bool, optional) – Whether to use automatic symmetry exploitation.
- Returns:
P – The projector onto the target eigenspace.
- Return type:
dense matrix
- quimb.calc.measure(p, A, eigenvalue=None, tol=1e-12)[source]¶
Measure state
p
with observableA
, collapsing the state in the process and returning the relevant eigenvalue.\[\left| \psi \right\rangle \rightarrow \frac{ P_i \left| \psi \right\rangle }{ \sqrt{\left\langle \psi \right| P_i \left| \psi \right\rangle} }\]and
\[\rho \rightarrow \frac{ P_i \rho P_i^{\dagger} }{ \text{Tr} \left[ P_i \rho \right] }\]along with the corresponding eigenvalue \(\lambda_i\).
- Parameters:
p (vector or matrix) – The quantum state to measure.
A (matrix or tuple[array, matrix]) – The hermitian operator corresponding to an observable. You can supply also a pre-diagonalized operator here as a tuple of eigenvalues and eigenvectors.
tol (float, optional) – The tolerance within which to group eigenspaces.
eigenvalue (float, optional) – If specified, deterministically collapse to this result. Otherwise randomly choose a result as in ‘real-life’.
- Returns:
result (float) – The result of the measurement.
p_after (vector or matrix) – The quantum state post measurement.
- quimb.calc.simulate_counts(p, C, phys_dim=2, seed=None)[source]¶
Simulate measuring each qubit of
p
in the computational basis, producing output like that ofqiskit
.- Parameters:
- Returns:
results – The counts for each bit string measured.
- Return type:
Examples
Simulate measuring the state of each qubit in a GHZ-state:
>>> import quimb as qu >>> psi = qu.ghz_state(3) >>> qu.simulate_counts(psi, 1024) {'000': 514, '111': 510}
- quimb.calc.dephase(rho, p, rand_rank=None)[source]¶
Dephase
rho
by amountp
, that is, mix it with the maximally mixed state:rho -> (1 - p) * rho + p * I / d
- Parameters:
- Returns:
rho_dephase – The dephased density operator.
- Return type:
operator
- quimb.calc.entropy(a, rank=None)[source]¶
Compute the (von Neumann) entropy.
- Parameters:
a (operator or 1d array) – Positive operator or list of positive eigenvalues.
rank (int (optional)) – If operator has known rank, then a partial decomposition can be used to accelerate the calculation.
- Returns:
The von Neumann entropy.
- Return type:
See also
mutinf
,entropy_subsys
,entropy_subsys_approx
- quimb.calc.entropy_subsys[source]¶
Calculate the entropy of a pure states’ subsystem, optionally switching to an approximate lanczos method when the subsystem is very large.
- Parameters:
psi_ab (vector) – Bipartite state.
dims (sequence of int) – The sub-dimensions of the state.
sysa (sequence of int) – The indices of which dimensions to calculate the entropy for.
approx_thresh (int, optional) – The size of sysa at which to switch to the approx method. Set to
None
to never use the approximation.**approx_opts – Supplied to
entropy_subsys_approx()
, if used.
- Returns:
The subsytem entropy.
- Return type:
See also
entropy
,entropy_subsys_approx
,mutinf_subsys
- quimb.calc.mutinf(p, dims=(2, 2), sysa=0, rank=None)[source]¶
Find the mutual information for a bipartition of a state.
That is,
H(A) + H(B) - H(AB)
, for von Neumann entropyH
, and two subsystems A and B.- Parameters:
p (vector or operator) – State, can be vector or operator.
sysa (int, optional) – Index of first subsystem, A.
sysb (int, optional) – Index of second subsystem, B.
rank (int, optional) – If known, the rank of rho_ab, to speed calculation of
H(AB)
up. For example, ifp
comes from tracing out three qubits from a system, then its rank is 2^3 = 8 etc.
- Return type:
See also
entropy
,mutinf_subsys
,entropy_subsys_approx
- quimb.calc.check_dims_and_indices(dims, *syss)[source]¶
Make sure all indices found in the tuples
syss
are inrange(len(dims))
.
- quimb.calc.mutinf_subsys(psi_abc, dims, sysa, sysb, approx_thresh=2**13, **approx_opts)[source]¶
Calculate the mutual information of two subsystems of a pure state, possibly using an approximate lanczos method for large subsytems.
- Parameters:
psi_abc (vector) – Tri-partite pure state.
dims (sequence of int) – The sub dimensions of the state.
sysa (sequence of int) – The index(es) of the subsystem(s) to consider part of ‘A’.
sysb (sequence of int) – The index(es) of the subsystem(s) to consider part of ‘B’.
approx_thresh (int, optional) – The size of subsystem at which to switch to the approximate lanczos method. Set to
None
to never use the approximation.approx_opts – Supplied to
entropy_subsys_approx()
, if used.
- Returns:
The mutual information.
- Return type:
See also
mutinf
,entropy_subsys
,entropy_subsys_approx
,logneg_subsys
- quimb.calc.schmidt_gap(psi_ab, dims, sysa)[source]¶
Find the schmidt gap of the bipartition of
psi_ab
. That is, the difference between the two largest eigenvalues of the reduced density operator.
- quimb.calc.tr_sqrt(A, rank=None)[source]¶
Return the trace of the sqrt of a positive semidefinite operator.
- quimb.calc.tr_sqrt_subsys[source]¶
Compute the trace sqrt of a subsystem, possibly using an approximate lanczos method when the subsytem is big.
- Parameters:
psi_ab (vector) – Bipartite state.
dims (sequence of int) – The sub-dimensions of the state.
sysa (sequence of int) – The indices of which dimensions to calculate the trace sqrt for.
approx_thresh (int, optional) – The size of sysa at which to switch to the approx method. Set to
None
to never use the approximation.**approx_opts – Supplied to
tr_sqrt_subsys_approx()
, if used.
- Returns:
The subsytem entropy.
- Return type:
See also
tr_sqrt
,tr_sqrt_subsys_approx
,partial_transpose_norm
- quimb.calc.partial_transpose(p, dims=(2, 2), sysa=0)[source]¶
Partial transpose of a density operator.
- Parameters:
- Return type:
operator
See also
- quimb.calc.partial_transpose_norm(p, dims, sysa)[source]¶
Compute the norm of the partial transpose for (log)-negativity, taking a shortcut (trace sqrt of reduced subsytem), when system is a vector.
- quimb.calc.logneg(p, dims=(2, 2), sysa=0)[source]¶
Compute logarithmic negativity between two subsytems. This is defined as log_2( | rho_{AB}^{T_B} | ). This only handles bipartitions (including pure states efficiently), and will not trace anything out.
- Parameters:
- Return type:
See also
negativity
,partial_transpose
,logneg_subsys_approx
- quimb.calc.logneg_subsys(psi_abc, dims, sysa, sysb, approx_thresh=2**13, **approx_opts)[source]¶
Compute the logarithmic negativity between two subsystems of a pure state, possibly using an approximate lanczos for large subsystems. Uses a special method if the two subsystems form a bipartition of the state.
- Parameters:
psi_abc (vector) – Tri-partite pure state.
dims (sequence of int) – The sub dimensions of the state.
sysa (sequence of int) – The index(es) of the subsystem(s) to consider part of ‘A’.
sysb (sequence of int) – The index(es) of the subsystem(s) to consider part of ‘B’.
approx_thresh (int, optional) – The size of subsystem at which to switch to the approximate lanczos method. Set to
None
to never use the approximation.approx_opts – Supplied to
logneg_subsys_approx()
, if used.
- Returns:
The logarithmic negativity.
- Return type:
See also
logneg
,mutinf_subsys
,logneg_subsys_approx
- quimb.calc.negativity(p, dims=(2, 2), sysa=0)[source]¶
Compute negativity between two subsytems.
This is defined as (| rho_{AB}^{T_B} | - 1) / 2. If
len(dims) > 2
, then the non-target dimensions will be traced out first.- Parameters:
- Return type:
See also
logneg
,partial_transpose
,negativity_subsys_approx
- quimb.calc.concurrence(p, dims=(2, 2), sysa=0, sysb=1)[source]¶
Concurrence of two-qubit state.
If
len(dims) > 2
, then the non-target dimensions will be traced out first.- Parameters:
- Return type:
- quimb.calc.one_way_classical_information(p_ab, prjs, precomp_func=False)[source]¶
One way classical information for two qubit density operator.
- Parameters:
p_ab (operator) – State of two qubits
prjs (sequence of matrices) – The POVMs.
precomp_func (bool, optional) – Whether to return a pre-computed function, closed over the actual state.
- Returns:
The one-way classical information or the function to compute it for the given state which takes a set of POVMs as its single argument.
- Return type:
float or callable
- quimb.calc.quantum_discord(p, dims=(2, 2), sysa=0, sysb=1, method='COBYLA', tol=1e-12, maxiter=2**14)[source]¶
Quantum Discord for two qubit density operator.
If
len(dims) > 2
, then the non-target dimensions will be traced out first.- Parameters:
- Return type:
- quimb.calc.trace_distance(p1, p2)[source]¶
Trace distance between two states:
\[\delta(\rho, \sigma) = \frac{1}{2} \left| \rho - \sigma \right|_\mathrm{tr}\]If two wavefunctions are supplied the trace distance will be computed via the more efficient expression:
\[\delta(|\psi\rangle\langle\psi|, |\phi\rangle\langle\phi|) = \sqrt{1 - \langle \psi | \phi \rangle^2}\]- Parameters:
p1 (ket or density operator) – The first state.
p2 (ket or density operator) – The second state.
- Return type:
- quimb.calc.cprint(psi, prec=6)[source]¶
Print a state in the computational basis.
- Parameters:
psi (vector) – The pure state.
prec (int, optional) – How many significant figures to print.
Examples
>>> cprint(rand_ket(2**2)) (-0.0751069-0.192635j) |00> (0.156837-0.562213j) |01> (-0.307381+0.0291168j) |10> (-0.14302+0.707661j) |11>
>>> cprint(w_state(4)) (0.5+0j) |0001> (0.5+0j) |0010> (0.5+0j) |0100> (0.5+0j) |1000>
- quimb.calc.decomp(a, fn, fn_args, fn_d, nmlz_func, mode='p', tol=0.001)[source]¶
Decomposes an operator via the Hilbert-schmidt inner product.
Can both print the decomposition or return it.
- Parameters:
a (ket or density operator) – Operator to decompose.
fn (callable) – Function to generate operator/state to decompose with.
fn_args – Sequence of args whose permutations will be supplied to
fn
.fn_d (int) – The dimension of the operators that fn produces.
nmlz_func (callable) – Function to produce a normlization coefficient given the
n
permutations of operators that will be produced.mode – String, include
'p'
to print the decomp and/or'c'
to return OrderedDict, sorted by size of contribution.tol – Print operators with contirbution above
tol
only.
- Returns:
Pauli operator name and expec with
a
.- Return type:
None or OrderedDict
See also
- quimb.calc.pauli_decomp¶
Decompose an operator into Paulis.
- quimb.calc.bell_decomp¶
Decompose an operator into bell-states.
- quimb.calc.correlation(p, A, B, sysa, sysb, dims=None, sparse=None, precomp_func=False)[source]¶
Calculate the correlation between two sites given two operators.
- Parameters:
p (ket or density operator) – State to compute correlations for, ignored if
precomp_func=True
.A (operator) – Operator to act on first subsystem.
B (operator) – Operator to act on second subsystem.
sysa (int) – Index of first subsystem.
sysb (int) – Index of second subsystem.
dims (tuple of int, optional) – Internal dimensions of
p
, will be assumed to be qubits if not given.sparse (bool, optional) – Whether to compute with sparse operators.
precomp_func (bool, optional) – Whether to return result or single arg function closed over precomputed operator.
- Returns:
The correlation, <ab> - <a><b>, or a function to compute for an arbitrary state.
- Return type:
float or callable
- quimb.calc.pauli_correlations(p, ss=('xx', 'yy', 'zz'), sysa=0, sysb=1, sum_abs=False, precomp_func=False)[source]¶
Calculate the correlation between sites for a list of operator pairs choisen from the pauli matrices.
- Parameters:
p (ket or density operator) – State to compute correlations for. Ignored if
precomp_func=True
.ss (tuple or str) – List of pairs specifiying pauli matrices.
sysa (int, optional) – Index of first site.
sysb (int, optional) – Index of second site.
sum_abs (bool, optional) – Whether to sum over the absolute values of each correlation
precomp_func (bool, optional) – whether to return the values or a single argument function closed over precomputed operators etc.
- Returns:
Either the value(s) of each correlation or the function(s) to compute the correlations for an arbitrary state, depending on
sum_abs
andprecomp_func
.- Return type:
- quimb.calc.ent_cross_matrix(p, sz_blc=1, ent_fn=logneg, calc_self_ent=True, upscale=False)[source]¶
Calculate the pair-wise function ent_fn between all sites or blocks of a state.
- Parameters:
p (ket or density operator) – State.
sz_blc (int) – Size of the blocks to partition the state into. If the number of individual sites is not a multiple of this then the final (smaller) block will be ignored.
ent_fn (callable) – Bipartite function, notionally entanglement
calc_self_ent (bool) – Whether to calculate the function for each site alone, purified. If the whole state is pure then this is the entanglement with the whole remaining system.
upscale (bool, optional) – Whether, if sz_blc != 1, to upscale the results so that the output array is the same size as if it was.
- Returns:
array of pairwise ent_fn results.
- Return type:
2D-array
- quimb.calc.qid(p, dims, inds, precomp_func=False, sparse_comp=True, norm_func=norm, power=2, coeff=1)[source]¶
- quimb.calc.is_degenerate(op, tol=1e-12)[source]¶
Check if operator has any degenerate eigenvalues, determined relative to mean spacing of all eigenvalues.
- quimb.calc.is_eigenvector(x, A, tol=1e-14)[source]¶
Determines whether a vector is an eigenvector of an operator.
- quimb.calc.page_entropy(sz_subsys, sz_total)[source]¶
Calculate the page entropy, i.e. expected entropy for a subsytem of a random state in Hilbert space.
- quimb.calc.heisenberg_energy(L)[source]¶
Get the analytic isotropic heisenberg chain ground energy for length L. Useful for testing. Assumes the heisenberg model is defined with spin operators not pauli matrices (overall factor of 2 smaller). Taken from [1].
[1] Nickel, Bernie. “Scaling corrections to the ground state energy of the spin-½ isotropic anti-ferromagnetic Heisenberg chain.” Journal of Physics Communications 1.5 (2017): 055021