quimb.tensor.circuit
¶
Tools for quantum circuit simulation using tensor networks.
TODO: - [ ] gate-by-gate sampling - [ ] sub-MPO apply for MPS simulation - [ ] multi qubit gates via MPO for MPS simulation
Module Contents¶
Classes¶
A simple class for storing the details of a quantum circuit gate. |
|
Class for simulating quantum circuits using tensor networks. |
|
Quantum circuit simulation keeping the state always in a MPS form. If |
|
Quantum circuit simulation keeping the state always in an MPS form, but |
|
Quantum circuit simulation keeping the state in full dense form. |
Functions¶
|
Parse a 'qsim' input format string into circuit information. |
|
Parse a qsim file. |
|
Parse a qsim url. |
|
Split, strip and filter a string by a given character into a list. |
|
Replace multiple substrings in a string. |
|
Parse the string contents of an OpenQASM 2.0 file. This parser does not |
|
Parse an OpenQASM 2.0 file. |
|
Parse an OpenQASM 2.0 url. |
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
Parametrized controlled X-rotation. |
|
Parametrized controlled Y-rotation. |
|
Parametrized controlled Z-rotation. |
|
|
|
|
|
|
|
Parametrized two qubit XX-rotation. |
|
Parametrized two qubit YY-rotation. |
|
Parametrized two qubit ZZ-rotation. |
|
See https://arxiv.org/abs/quant-ph/0308006 - Fig. 7. |
|
|
|
Build a low rank hyper tensor network (CP-decomp like) representation of |
|
Apply a multi-controlled gate to a state represented as an MPS. |
|
|
|
|
|
|
Sample a bitstring from n-dimensional tensor |
|
|
|
|
Map all types of gate specification into a Gate object. |
Attributes¶
- quimb.tensor.circuit.parse_qsim_str(contents)[source]¶
Parse a ‘qsim’ input format string into circuit information.
The format is described here: https://quantumai.google/qsim/input_format.
- Parameters:
contents (str) – The full string of the qsim file.
- Returns:
circuit_info – Information about the circuit:
circuit_info[‘n’]: the number of qubits
circuit_info[‘n_gates’]: the number of gates in total
circuit_info[‘gates’]: list[list[str]], list of gates, each of which is a list of strings read from a line of the qsim file.
- Return type:
- quimb.tensor.circuit.to_clean_list(s, delimiter)[source]¶
Split, strip and filter a string by a given character into a list.
- quimb.tensor.circuit.multi_replace(s, replacements)[source]¶
Replace multiple substrings in a string.
- quimb.tensor.circuit.parse_openqasm2_str(contents)[source]¶
Parse the string contents of an OpenQASM 2.0 file. This parser does not support classical control flow is not guaranteed to check the full openqasm grammar.
- quimb.tensor.circuit.ALL_GATES¶
- quimb.tensor.circuit.ONE_QUBIT_GATES¶
- quimb.tensor.circuit.TWO_QUBIT_GATES¶
- quimb.tensor.circuit.ALL_PARAM_GATES¶
- quimb.tensor.circuit.ONE_QUBIT_PARAM_GATES¶
- quimb.tensor.circuit.TWO_QUBIT_PARAM_GATES¶
- quimb.tensor.circuit.GATE_TAGS¶
- quimb.tensor.circuit.GATE_SIZE¶
- quimb.tensor.circuit.CONSTANT_GATES¶
- quimb.tensor.circuit.PARAM_GATES¶
- quimb.tensor.circuit.SPECIAL_GATES¶
- quimb.tensor.circuit.rxx_param_gen(params)[source]¶
Parametrized two qubit XX-rotation.
\[\mathrm{RXX}(\theta) = \exp(-i \frac{\theta}{2} X_i X_j)\]
- quimb.tensor.circuit.ryy_param_gen(params)[source]¶
Parametrized two qubit YY-rotation.
\[\mathrm{RYY}(\theta) = \exp(-i \frac{\theta}{2} Y_i Y_j)\]
- quimb.tensor.circuit.rzz_param_gen(params)[source]¶
Parametrized two qubit ZZ-rotation.
\[\mathrm{RZZ}(\theta) = \exp(-i \frac{\theta}{2} Z_i Z_j)\]
- quimb.tensor.circuit.su4_gate_param_gen(params)[source]¶
See https://arxiv.org/abs/quant-ph/0308006 - Fig. 7. params: # theta1, phi1, lamda1, # theta2, phi2, lamda2, # theta3, phi3, lamda3, # theta4, phi4, lamda4, # t1, t2, t3,
- quimb.tensor.circuit._MPS_METHODS¶
- quimb.tensor.circuit.build_controlled_gate_htn(ncontrol, gate, upper_inds, lower_inds, tags_each=None, tags_all=None, bond_ind=None)[source]¶
Build a low rank hyper tensor network (CP-decomp like) representation of a multi controlled gate.
- quimb.tensor.circuit._apply_controlled_gate_mps(psi, gate, tags=None, **gate_opts)[source]¶
Apply a multi-controlled gate to a state represented as an MPS.
- quimb.tensor.circuit._apply_controlled_gate_htn(psi, gate, tags=None, propagate_tags='register', **gate_opts)[source]¶
- quimb.tensor.circuit.apply_controlled_gate(psi, gate, tags=None, contract='auto-split-gate', propagate_tags='register', **gate_opts)[source]¶
- class quimb.tensor.circuit.Gate(label, params, qubits=None, controls=None, round=None, parametrize=False)[source]¶
A simple class for storing the details of a quantum circuit gate.
- Parameters:
label (str) – The name or ‘identifier’ of the gate.
params (Iterable[float]) – The parameters of the gate.
qubits (Iterable[int], optional) – Which qubits the gate acts on.
controls (Iterable[int], optional) – Which qubits are the controls.
round (int, optional) – If given, which round or layer the gate is part of.
parametrize (bool, optional) – Whether the gate will correspond to a parametrized tensor.
- property label¶
- property params¶
- property qubits¶
- property total_qubit_count¶
- property controls¶
- property round¶
- property special¶
- property parametrize¶
- property tag¶
- property array¶
- __slots__ = ('_label', '_params', '_qubits', '_controls', '_round', '_parametrize', '_tag', '_special',...¶
- quimb.tensor.circuit.sample_bitstring_from_prob_ndarray(p)[source]¶
Sample a bitstring from n-dimensional tensor
p
of probabilities.Examples
>>> import numpy as np >>> p = np.zeros(shape=(2, 2, 2, 2, 2)) >>> p[0, 1, 0, 1, 1] = 1.0 >>> sample_bitstring_from_prob_ndarray(p) '01011'
- quimb.tensor.circuit.parse_to_gate(gate_id, *gate_args, params=None, qubits=None, controls=None, gate_round=None, parametrize=None)[source]¶
Map all types of gate specification into a Gate object.
- class quimb.tensor.circuit.Circuit(N=None, psi0=None, gate_opts=None, gate_contract='auto-split-gate', gate_propagate_tags='register', tags=None, psi0_dtype='complex128', psi0_tag='PSI0', tag_gate_numbers=True, tag_gate_rounds=True, tag_gate_labels=True, bra_site_ind_id='b{}', to_backend=None)[source]¶
Class for simulating quantum circuits using tensor networks.
- Parameters:
N (int, optional) – The number of qubits.
psi0 (TensorNetwork1DVector, optional) – The initial state, assumed to be
|00000....0>
if not given. The state is always copied and the tagPSI0
added.gate_opts (dict_like, optional) – Default keyword arguments to supply to each
gate_TN_1D()
call during the circuit.gate_contract (str, optional) – Shortcut for setting the default ‘contract’ option in gate_opts.
gate_propagate_tags (str, optional) – Shortcut for setting the default ‘propagate_tags’ option in gate_opts.
tags (str or sequence of str, optional) – Tag(s) to add to the initial wavefunction tensors (whether these are propagated to the rest of the circuit’s tensors depends on
gate_opts
).psi0_dtype (str, optional) – Ensure the initial state has this dtype.
psi0_tag (str, optional) – Ensure the initial state has this tag.
tag_gate_numbers (bool, optional) – Whether to tag each gate tensor with its number in the circuit, like
"GATE_{g}"
.tag_gate_rounds (bool, optional) – Whether to tag each gate tensor with its number in the circuit, like
"ROUND_{r}"
.tag_gate_labels (bool, optional) – Whether to tag each gate tensor with its gate type label, e.g.
{"X_1/2", "ISWAP", "CCX", ...}
..bra_site_ind_id (str, optional) – Use this to label ‘bra’ site indices when creating certain (mostly internal) intermediate tensor networks.
- psi¶
The current wavefunction.
- Type:
Examples
Create 3-qubit GHZ-state:
>>> qc = qtn.Circuit(3) >>> gates = [ ('H', 0), ('H', 1), ('CNOT', 1, 2), ('CNOT', 0, 2), ('H', 0), ('H', 1), ('H', 2), ] >>> qc.apply_gates(gates) >>> qc.psi <TensorNetwork1DVector(tensors=12, indices=14, L=3, max_bond=2)>
>>> qc.psi.to_dense().round(4) qarray([[ 0.7071+0.j], [ 0. +0.j], [ 0. +0.j], [-0. +0.j], [-0. +0.j], [ 0. +0.j], [ 0. +0.j], [ 0.7071+0.j]])
>>> for b in qc.sample(10): ... print(b) 000 000 111 000 111 111 000 111 000 000
- property psi¶
Tensor network representation of the wavefunction.
- property uni¶
- property num_gates¶
- classmethod from_qsim_str(contents, **circuit_opts)[source]¶
Generate a
Circuit
instance from a ‘qsim’ string.
- classmethod from_qsim_file(fname, **circuit_opts)[source]¶
Generate a
Circuit
instance from a ‘qsim’ file.The qsim file format is described here: https://quantumai.google/qsim/input_format.
- classmethod from_qsim_url(url, **circuit_opts)[source]¶
Generate a
Circuit
instance from a ‘qsim’ url.
- classmethod from_openqasm2_str(contents, **circuit_opts)[source]¶
Generate a
Circuit
instance from an OpenQASM 2.0 string.
- classmethod from_openqasm2_file(fname, **circuit_opts)[source]¶
Generate a
Circuit
instance from an OpenQASM 2.0 file.
- classmethod from_openqasm2_url(url, **circuit_opts)[source]¶
Generate a
Circuit
instance from an OpenQASM 2.0 url.
- classmethod from_gates(gates, N=None, progbar=False, **kwargs)[source]¶
Generate a
Circuit
instance from a sequence of gates.
- _apply_gate(gate, tags=None, **gate_opts)[source]¶
Apply a
Gate
to thisCircuit
. This is the main method that all calls to apply a gate should go through.
- apply_gate(gate_id, *gate_args, params=None, qubits=None, controls=None, gate_round=None, parametrize=None, **gate_opts)[source]¶
Apply a single gate to this tensor network quantum circuit. If
gate_round
is supplied the tensor(s) added will be tagged with'ROUND_{gate_round}'
. Alternatively, putting an integer first like so:circuit.apply_gate(10, 'H', 7)
Is automatically translated to:
circuit.apply_gate('H', 7, gate_round=10)
- Parameters:
gate_id (Gate, str, or array_like) –
Which gate to apply. This can be:
A
Gate
instance, i.e. with parameters and qubits already specified.A string, e.g.
'H'
,'U3'
, etc. in which casegate_args
should be supplied with(*params, *qubits)
.A raw array, in which case
gate_args
should be supplied with(*qubits,)
.
gate_round (int, optional) – The gate round. If
gate_id
is integer-like, will also be taken from here, with thengate_id, gate_args = gate_args[0], gate_args[1:]
.gate_opts – Supplied to the gate function, options here will override the default
gate_opts
.
- apply_gate_raw(U, where, controls=None, gate_round=None, **gate_opts)[source]¶
Apply the raw array
U
as a gate on qubits inwhere
. It will be assumed to be unitary for the sake of computing reverse lightcones.
- apply_gates(gates, progbar=False, **gate_opts)[source]¶
Apply a sequence of gates to this tensor network quantum circuit.
- Parameters:
gates (Sequence[Gate] or Sequence[Tuple]) – The sequence of gates to apply.
gate_opts – Supplied to
apply_gate()
.
- su4(theta1, phi1, lamda1, theta2, phi2, lamda2, theta3, phi3, lamda3, theta4, phi4, lamda4, t1, t2, t3, i, j, gate_round=None, parametrize=False, **kwargs)[source]¶
- get_uni(transposed=False)[source]¶
Tensor network representation of the unitary operator (i.e. with the initial state removed).
- get_reverse_lightcone_tags(where)[source]¶
Get the tags of gates in this circuit corresponding to the ‘reverse’ lightcone propagating backwards from registers in
where
.
- get_psi_reverse_lightcone(where, keep_psi0=False)[source]¶
Get just the bit of the wavefunction in the reverse lightcone of sites in
where
- i.e. causally linked.- Parameters:
where (int, or sequence of int) – The sites to propagate the the lightcone back from, supplied to
get_reverse_lightcone_tags()
.keep_psi0 (bool, optional) – Keep the tensors corresponding to the initial wavefunction regardless of whether they are outside of the lightcone.
- Returns:
psi_lc
- Return type:
- get_psi_simplified(seq='ADCRS', atol=1e-12, equalize_norms=False)[source]¶
Get the full wavefunction post local tensor network simplification.
- Parameters:
seq (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
- Returns:
psi
- Return type:
- get_rdm_lightcone_simplified(where, seq='ADCRS', atol=1e-12, equalize_norms=False)[source]¶
Get a simplified TN of the norm of the wavefunction, with gates outside reverse lightcone of
where
cancelled, and physical indices withinwhere
preserved so that they can be fixed (sliced) or used as output indices.- Parameters:
where (int or sequence of int) – The region assumed to be the target density matrix essentially. Supplied to
get_reverse_lightcone_tags()
.seq (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
- Return type:
- amplitude(b, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-12, simplify_equalize_norms=False, backend=None, dtype='complex128', rehearse=False)[source]¶
Get the amplitude coefficient of bitstring
b
.\[c_b = \langle b | \psi \rangle\]- Parameters:
b (str or sequence of int) – The bitstring to compute the transition amplitude for.
optimize (str, optional) – Contraction path optimizer to use for the amplitude, can be a non-reusable path optimizer as only called once (though path won’t be cached for later use in that case).
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
backend (str, optional) – Backend to perform the contraction with, e.g.
'numpy'
,'cupy'
or'jax'
. Passed tocotengra
.dtype (str, optional) – Data type to cast the TN to before contraction.
rehearse (bool or "tn", optional) – If
True
, generate and cache the simplified tensor network and contraction tree but don’t actually perform the contraction. Returns a dict with keys"tn"
and'tree'
with the tensor network that will be contracted and the corresponding contraction tree if so.
- amplitude_rehearse(b='random', simplify_sequence='ADCRS', simplify_atol=1e-12, simplify_equalize_norms=False, optimize='auto-hq', dtype='complex128', rehearse=True)[source]¶
Perform just the tensor network simplifications and contraction tree finding associated with computing a single amplitude (caching the results) but don’t perform the actual contraction.
- Parameters:
b ('random', str or sequence of int) – The bitstring to rehearse computing the transition amplitude for, if
'random'
(the default) a random bitstring will be used.optimize (str, optional) – Contraction path optimizer to use for the marginal, can be a non-reusable path optimizer as only called once (though path won’t be cached for later use in that case).
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
backend (str, optional) – Backend to perform the marginal contraction with, e.g.
'numpy'
,'cupy'
or'jax'
. Passed tocotengra
.dtype (str, optional) – Data type to cast the TN to before contraction.
- Return type:
- partial_trace(keep, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-12, simplify_equalize_norms=False, backend=None, dtype='complex128', rehearse=False)[source]¶
Perform the partial trace on the circuit wavefunction, retaining only qubits in
keep
, and making use of reverse lightcone cancellation:\[\rho_{\bar{q}} = Tr_{\bar{p}} |\psi_{\bar{q}} \rangle \langle \psi_{\bar{q}}|\]Where \(\bar{q}\) is the set of qubits to keep, \(\psi_{\bar{q}}\) is the circuit wavefunction only with gates in the causal cone of this set, and \(\bar{p}\) is the remaining qubits.
- Parameters:
keep (int or sequence of int) – The qubit(s) to keep as we trace out the rest.
optimize (str, optional) – Contraction path optimizer to use for the reduced density matrix, can be a non-reusable path optimizer as only called once (though path won’t be cached for later use in that case).
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
backend (str, optional) – Backend to perform the marginal contraction with, e.g.
'numpy'
,'cupy'
or'jax'
. Passed tocotengra
.dtype (str, optional) – Data type to cast the TN to before contraction.
rehearse (bool or "tn", optional) – If
True
, generate and cache the simplified tensor network and contraction tree but don’t actually perform the contraction. Returns a dict with keys"tn"
and'tree'
with the tensor network that will be contracted and the corresponding contraction tree if so.
- Return type:
array or dict
- local_expectation(G, where, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-12, simplify_equalize_norms=False, backend=None, dtype='complex128', rehearse=False)[source]¶
Compute the a single expectation value of operator
G
, acting on siteswhere
, making use of reverse lightcone cancellation.\[\langle \psi_{\bar{q}} | G_{\bar{q}} | \psi_{\bar{q}} \rangle\]where \(\bar{q}\) is the set of qubits \(G\) acts one and \(\psi_{\bar{q}}\) is the circuit wavefunction only with gates in the causal cone of this set. If you supply a tuple or list of gates then the expectations will be computed simultaneously.
- Parameters:
G (array or sequence[array]) – The raw operator(s) to find the expectation of.
where (int or sequence of int) – Which qubits the operator acts on.
optimize (str, optional) – Contraction path optimizer to use for the local expectation, can be a non-reusable path optimizer as only called once (though path won’t be cached for later use in that case).
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
backend (str, optional) – Backend to perform the marginal contraction with, e.g.
'numpy'
,'cupy'
or'jax'
. Passed tocotengra
.dtype (str, optional) – Data type to cast the TN to before contraction.
gate_opts (None or dict_like) – Options to use when applying
G
to the wavefunction.rehearse (bool or "tn", optional) – If
True
, generate and cache the simplified tensor network and contraction tree but don’t actually perform the contraction. Returns a dict with keys'tn'
and'tree'
with the tensor network that will be contracted and the corresponding contraction tree if so.
- Return type:
- compute_marginal(where, fix=None, optimize='auto-hq', backend=None, dtype='complex64', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=True, rehearse=False)[source]¶
Compute the probability tensor of qubits in
where
, given possibly fixed qubits infix
and tracing everything else having removed redundant unitary gates.- Parameters:
where (sequence of int) – The qubits to compute the marginal probability distribution of.
fix (None or dict[int, str], optional) – Measurement results on other qubits to fix.
optimize (str, optional) – Contraction path optimizer to use for the marginal, can be a non-reusable path optimizer as only called once (though path won’t be cached for later use in that case).
backend (str, optional) – Backend to perform the marginal contraction with, e.g.
'numpy'
,'cupy'
or'jax'
. Passed tocotengra
.dtype (str, optional) – Data type to cast the TN to before contraction.
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
rehearse (bool or "tn", optional) – Whether to perform the marginal contraction or just return the associated TN and contraction tree.
- calc_qubit_ordering(qubits=None, method='greedy-lightcone')[source]¶
Get a order to measure
qubits
in, by greedily choosing whichever has the smallest reverse lightcone followed by whichever expands this lightcone least.
- _parse_qubits_order(qubits=None, order=None)[source]¶
Simply initializes the default of measuring all qubits, and the default order, or checks that
order
is a permutation ofqubits
.
- _group_order(order, group_size=1)[source]¶
Take the qubit ordering
order
and batch it in groups of sizegroup_size
, sorting the qubits (for caching reasons) within each group.
- sample(C, qubits=None, order=None, group_size=1, max_marginal_storage=2**20, seed=None, optimize='auto-hq', backend=None, dtype='complex64', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=False)[source]¶
Sample the circuit given by
gates
,C
times, using lightcone cancelling and caching marginal distribution results. This is a generator. This proceeds as a chain of marginal computations.Assuming we have
group_size=1
, and some ordering of the qubits, \(\{q_0, q_1, q_2, q_3, \ldots\}\) we first compute:\[p(q_0) = \mathrm{diag} \mathrm{Tr}_{1, 2, 3,\ldots} | \psi_{0} \rangle \langle \psi_{0} |\]I.e. simply the probability distribution on a single qubit, conditioned on nothing. The subscript on \(\psi\) refers to the fact that we only need gates from the causal cone of qubit 0. From this we can sample an outcome, either 0 or 1, if we call this \(r_0\) we can then move on to the next marginal:
\[p(q_1 | r_0) = \mathrm{diag} \mathrm{Tr}_{2, 3,\ldots} \langle r_0 | \psi_{0, 1} \rangle \langle \psi_{0, 1} | r_0 \rangle\]I.e. the probability distribution of the next qubit, given our prior result. We can sample from this to get \(r_1\). Then we compute:
\[p(q_2 | r_0 r_1) = \mathrm{diag} \mathrm{Tr}_{3,\ldots} \langle r_0 r_1 | \psi_{0, 1, 2} \rangle \langle \psi_{0, 1, 2} | r_0 r_1 \rangle\]Eventually we will reach the ‘final marginal’, which we can compute as
\[|\langle r_0 r_1 r_2 r_3 \ldots | \psi \rangle|^2\]since there is nothing left to trace out.
- Parameters:
C (int) – The number of times to sample.
qubits (None or sequence of int, optional) – Which qubits to measure, defaults (
None
) to all qubits.order (None or sequence of int, optional) – Which order to measure the qubits in, defaults (
None
) to an order based on greedily expanding the smallest reverse lightcone. If specified it should be a permutation ofqubits
.group_size (int, optional) – How many qubits to group together into marginals, the larger this is the fewer marginals need to be computed, which can be faster at the cost of higher memory. The marginal themselves will each be of size
2**group_size
.max_marginal_storage (int, optional) – The total cumulative number of marginal probabilites to cache, once this is exceeded caching will be turned off.
seed (None or int, optional) – A random seed, passed to
numpy.random.seed
if given.optimize (str, optional) – Contraction path optimizer to use for the marginals, shouldn’t be a non-reusable path optimizer as called on many different TNs. Passed to
cotengra.array_contract_tree()
.backend (str, optional) – Backend to perform the marginal contraction with, e.g.
'numpy'
,'cupy'
or'jax'
. Passed tocotengra
.dtype (str, optional) – Data type to cast the TN to before contraction.
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
- Yields:
bitstrings (sequence of str)
- sample_rehearse(qubits=None, order=None, group_size=1, result=None, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=False, rehearse=True, progbar=False)[source]¶
Perform the preparations and contraction tree findings for
sample()
, caching various intermedidate objects, but don’t perform the main contractions.- Parameters:
qubits (None or sequence of int, optional) – Which qubits to measure, defaults (
None
) to all qubits.order (None or sequence of int, optional) – Which order to measure the qubits in, defaults (
None
) to an order based on greedily expanding the smallest reverse lightcone.group_size (int, optional) – How many qubits to group together into marginals, the larger this is the fewer marginals need to be computed, which can be faster at the cost of higher memory. The marginal’s size itself is exponential in
group_size
.result (None or dict[int, str], optional) – Explicitly check the computational cost of this result, assumed to be all zeros if not given.
optimize (str, optional) – Contraction path optimizer to use for the marginals, shouldn’t be a non-reusable path optimizer as called on many different TNs. Passed to
cotengra.array_contract_tree()
.simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
progbar (bool, optional) – Whether to show the progress of finding each contraction tree.
- Returns:
One contraction tree object per grouped marginal computation. The keys of the dict are the qubits the marginal is computed for, the values are a dict containing a representative simplified tensor network (key: ‘tn’) and the main contraction tree (key: ‘tree’).
- Return type:
- sample_chaotic(C, marginal_qubits, max_marginal_storage=2**20, seed=None, optimize='auto-hq', backend=None, dtype='complex64', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=False)[source]¶
Sample from this circuit, assuming it to be chaotic. Which is to say, only compute and sample correctly from the final marginal, assuming that the distribution on the other qubits is uniform. Given
marginal_qubits=5
for instance, for each sample a random bit-string \(r_0 r_1 r_2 \ldots r_{N - 6}\) for the remaining \(N - 5\) qubits will be chosen, then the final marginal will be computed as\[p(q_{N-5}q_{N-4}q_{N-3}q_{N-2}q_{N-1} | r_0 r_1 r_2 \ldots r_{N-6}) = |\langle r_0 r_1 r_2 \ldots r_{N - 6} | \psi \rangle|^2\]and then sampled from. Note the expression on the right hand side has 5 open indices here and so is a tensor, however if
marginal_qubits
is not too big then the cost of contracting this is very similar to a single amplitude.Note
This method assumes the circuit is chaotic, if its not, then the samples produced will not be an accurate representation of the probability distribution.
- Parameters:
C (int) – The number of times to sample.
marginal_qubits (int or sequence of int) – The number of qubits to treat as marginal, or the actual qubits. If an int is given then the qubits treated as marginal will be
circuit.calc_qubit_ordering()[:marginal_qubits]
.seed (None or int, optional) – A random seed, passed to
numpy.random.seed
if given.optimize (str, optional) – Contraction path optimizer to use for the marginal, can be a non-reusable path optimizer as only called once (though path won’t be cached for later use in that case).
backend (str, optional) – Backend to perform the marginal contraction with, e.g.
'numpy'
,'cupy'
or'jax'
. Passed tocotengra
.dtype (str, optional) – Data type to cast the TN to before contraction.
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
- Yields:
str
- sample_chaotic_rehearse(marginal_qubits, result=None, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=False, dtype='complex64', rehearse=True)[source]¶
Rehearse chaotic sampling (perform just the TN simplifications and contraction tree finding).
- Parameters:
marginal_qubits (int or sequence of int) – The number of qubits to treat as marginal, or the actual qubits. If an int is given then the qubits treated as marginal will be
circuit.calc_qubit_ordering()[:marginal_qubits]
.result (None or dict[int, str], optional) – Explicitly check the computational cost of this result, assumed to be all zeros if not given.
optimize (str, optional) – Contraction path optimizer to use for the marginal, can be a non-reusable path optimizer as only called once (though path won’t be cached for later use in that case).
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
dtype (str, optional) – Data type to cast the TN to before contraction.
- Returns:
The contraction path information for the main computation, the key is the qubits that formed the final marginal. The value is itself a dict with keys
'tn'
- a representative tensor network - and'tree'
- the contraction tree.- Return type:
- to_dense(reverse=False, optimize='auto-hq', simplify_sequence='R', simplify_atol=1e-12, simplify_equalize_norms=False, backend=None, dtype=None, rehearse=False)[source]¶
Generate the dense representation of the final wavefunction.
- Parameters:
reverse (bool, optional) – Whether to reverse the order of the subsystems, to match the convention of qiskit for example.
optimize (str, optional) – Contraction path optimizer to use for the contraction, can be a non-reusable path optimizer as only called once (though path won’t be cached for later use in that case).
dtype (str, optional) – If given, convert the tensors to this dtype prior to contraction.
simplify_sequence (str, optional) – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.simplify_equalize_norms (bool, optional) – Actively renormalize tensor norms during simplification.
backend (str, optional) – Backend to perform the contraction with, e.g.
'numpy'
,'cupy'
or'jax'
. Passed tocotengra
.dtype – Data type to cast the TN to before contraction.
rehearse (bool, optional) – If
True
, generate and cache the simplified tensor network and contraction tree but don’t actually perform the contraction. Returns a dict with keys'tn'
and'tree'
with the tensor network that will be contracted and the corresponding contraction tree if so.
- Returns:
psi – The densely represented wavefunction with
dtype
data.- Return type:
- simulate_counts(C, seed=None, reverse=False, **to_dense_opts)[source]¶
Simulate measuring all qubits in the computational basis many times. Unlike
sample()
, this generates all the samples simultaneously using the full wavefunction constructed fromto_dense()
, then callingsimulate_counts()
.Warning
Because this constructs the full wavefunction it always requires exponential memory in the number of qubits, regardless of circuit depth and structure.
- Parameters:
C (int) – The number of ‘experimental runs’, i.e. total counts.
seed (int, optional) – A seed for reproducibility.
reverse (bool, optional) – Whether to reverse the order of the subsystems, to match the convention of qiskit for example.
to_dense_opts – Suppled to
to_dense()
.
- Returns:
results – The number of recorded counts for each
- Return type:
- xeb(samples_or_counts, cache=None, cache_maxsize=2**20, progbar=False, **amplitude_opts)[source]¶
Compute the linear cross entropy benchmark (XEB) for samples or counts, amplitude per amplitude.
- Parameters:
samples_or_counts (Iterable[str] or Dict[str, int]) – Either the raw bitstring samples or a dict mapping bitstrings to the number of counts observed.
cache (dict, optional) – A dictionary to store the probabilities in, if not supplied
quimb.utils.LRU(cache_maxsize)
will be used.cache_maxsize – The maximum size of the cache to be used.
optional – The maximum size of the cache to be used.
progbar – Whether to show progress as the bitstrings are iterated over.
optional – Whether to show progress as the bitstrings are iterated over.
amplitude_opts – Supplied to
amplitude()
.
- xeb_ex(optimize='auto-hq', simplify_sequence='R', simplify_atol=1e-12, simplify_equalize_norms=False, dtype=None, backend=None, autojit=False, progbar=False, **contract_opts)[source]¶
Compute the exactly expected XEB for this circuit. The main feature here is that if you supply a cotengra optimizer that searches for sliced indices then the XEB will be computed without constructing the full wavefunction.
- Parameters:
optimize (str or PathOptimizer, optional) – Contraction path optimizer.
simplify_sequence (str, optional) – Simplifications to apply to tensor network prior to contraction.
simplify_sequence – Which local tensor network simplifications to perform and in which order, see
full_simplify()
.simplify_atol (float, optional) – The tolerance with which to compare to zero when applying
full_simplify()
.dtype (str, optional) – Data type to cast the TN to before contraction.
backend (str, optional) – Convert tensors to, and then use contractions from, this library.
autojit (bool, optional) – Apply
autoray.autojit
to the contraciton and map-reduce.progbar (bool, optional) – Show progress in terms of number of wavefunction chunks processed.
- update_params_from(tn)[source]¶
Assuming
tn
is a tensor network with tensors taggedGATE_{i}
corresponding to this circuit (e.g. fromcirc.psi
orcirc.uni
) but with updated parameters, update the current circuit parameters and tensors with those values.This is an inplace modification of the
Circuit
.- Parameters:
tn (TensorNetwork) – The tensor network to find the updated parameters from.
- class quimb.tensor.circuit.CircuitMPS(N=None, *, psi0=None, max_bond=None, cutoff=1e-10, gate_opts=None, gate_contract='auto-mps', **circuit_opts)[source]¶
Bases:
Circuit
Quantum circuit simulation keeping the state always in a MPS form. If you think the circuit will not build up much entanglement, or you just want to keep a rigorous handle on how much entanglement is present, this can be useful.
- Parameters:
N (int, optional) – The number of qubits in the circuit.
psi0 (TensorNetwork1DVector, optional) – The initial state, assumed to be
|00000....0>
if not given. The state is always copied and the tagPSI0
added.max_bond (int, optional) – The maximum bond dimension to truncate to when applying gates, if any. This is simply a shortcut for setting
gate_opts['max_bond']
.cutoff (float, optional) – The singular value cutoff to use when truncating the state. This is simply a shortcut for setting
gate_opts['cutoff']
.gate_opts (dict, optional) – Default options to pass to each gate, for example, “max_bond” and “cutoff” etc.
gate_contract (str, optional) –
The default method for applying gates. Relevant MPS options are:
'auto-mps'
: automatically choose a method that maintains the MPS form (default). This uses'swap+split'
for 2-qubit gates and'nonlocal'
for 3+ qubit gates.'swap+split'
: swap nonlocal qubits to be next to each other, before applying the gate, then swapping them back'nonlocal'
: turn the gate into a potentially nonlocal (sub) MPO and apply it directly. Seetensor_network_1d_compress()
.
circuit_opts – Supplied to
Circuit
.
- psi¶
The current state of the circuit, always in MPS form.
- Type:
Examples
Create a circuit object that always uses the “nonlocal” method for contracting in gates, and the “dm” compression method within that, using a large cutoff and maximum bond dimension:
circ = qtn.CircuitMPS( N=56, gate_opts=dict( contract="nonlocal", method="dm", max_bond=1024, cutoff=1e-3, ) )
- property psi¶
Tensor network representation of the wavefunction.
- property uni¶
- apply_gates(gates, progbar=False, **gate_opts)[source]¶
Apply a sequence of gates to this tensor network quantum circuit.
- Parameters:
gates (Sequence[Gate] or Sequence[Tuple]) – The sequence of gates to apply.
gate_opts – Supplied to
apply_gate()
.
- get_psi_reverse_lightcone(where, keep_psi0=False)[source]¶
Override
get_psi_reverse_lightcone
as for an MPS the lightcone is not meaningful.
- fidelity_estimate()[source]¶
Estimate the fidelity of the current state based on its norm, which tracks how much the state has been truncated:
\[\tilde{F} = \left| \langle \psi | \psi \rangle \right|^2 \approx \left|\langle \psi_\mathrm{ideal} | \psi \rangle\right|^2\]See also
- class quimb.tensor.circuit.CircuitPermMPS(N=None, psi0=None, gate_opts=None, gate_contract='swap+split', **circuit_opts)[source]¶
Bases:
CircuitMPS
Quantum circuit simulation keeping the state always in an MPS form, but lazily tracking the qubit ordering rather than ‘swapping back’ qubits after applying non-local gates. This can be useful for circuits with no expectation of locality. The qubit ordering is always tracked in the attribute
qubits
. Thepsi
attribute returns the TN with the sites reindexed and retagged according to the current qubit ordering, meaning it is no longer an MPS. Use circ.get_psi_unordered() to get the unpermuted MPS and use circ.qubits to get the current qubit ordering if you prefer.- property psi¶
Tensor network representation of the wavefunction.