quimb.tensor.circuit ==================== .. py:module:: quimb.tensor.circuit .. autoapi-nested-parse:: Tools for quantum circuit simulation using tensor networks. Attributes ---------- .. autoapisummary:: 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._MPS_METHODS Classes ------- .. autoapisummary:: quimb.tensor.circuit.Gate quimb.tensor.circuit.Circuit quimb.tensor.circuit.CircuitMPS quimb.tensor.circuit.CircuitPermMPS quimb.tensor.circuit.CircuitDense Functions --------- .. autoapisummary:: quimb.tensor.circuit.recursive_stack quimb.tensor.circuit._convert_ints_and_floats quimb.tensor.circuit._put_registers_last quimb.tensor.circuit.parse_qsim_str quimb.tensor.circuit.parse_qsim_file quimb.tensor.circuit.parse_qsim_url quimb.tensor.circuit.to_clean_list quimb.tensor.circuit.multi_replace quimb.tensor.circuit.get_openqasm2_regexes quimb.tensor.circuit.parse_openqasm2_str quimb.tensor.circuit.parse_openqasm2_file quimb.tensor.circuit.parse_openqasm2_url quimb.tensor.circuit.register_constant_gate quimb.tensor.circuit.register_param_gate quimb.tensor.circuit.register_special_gate quimb.tensor.circuit.rx_gate_param_gen quimb.tensor.circuit.ry_gate_param_gen quimb.tensor.circuit.rz_gate_param_gen quimb.tensor.circuit.u3_gate_param_gen quimb.tensor.circuit.u2_gate_param_gen quimb.tensor.circuit.u1_gate_param_gen quimb.tensor.circuit.cu3_param_gen quimb.tensor.circuit.cu2_param_gen quimb.tensor.circuit.cu1_param_gen quimb.tensor.circuit.crx_param_gen quimb.tensor.circuit.cry_param_gen quimb.tensor.circuit.crz_param_gen quimb.tensor.circuit.fsim_param_gen quimb.tensor.circuit.fsimg_param_gen quimb.tensor.circuit.givens_param_gen quimb.tensor.circuit.givens2_param_gen quimb.tensor.circuit.xx_plus_yy_param_gen quimb.tensor.circuit.xx_minus_yy_param_gen quimb.tensor.circuit.rxx_param_gen quimb.tensor.circuit.ryy_param_gen quimb.tensor.circuit.rzz_param_gen quimb.tensor.circuit.su4_gate_param_gen quimb.tensor.circuit.apply_swap quimb.tensor.circuit.build_controlled_gate_htn quimb.tensor.circuit._apply_controlled_gate_mps quimb.tensor.circuit._apply_controlled_gate_htn quimb.tensor.circuit.apply_controlled_gate quimb.tensor.circuit._cached_param_gate_build quimb.tensor.circuit.sample_bitstring_from_prob_ndarray quimb.tensor.circuit.rehearsal_dict quimb.tensor.circuit.parse_to_gate Module Contents --------------- .. py:function:: recursive_stack(x) .. py:function:: _convert_ints_and_floats(x) .. py:function:: _put_registers_last(x) .. py:function:: parse_qsim_str(contents) Parse a 'qsim' input format string into circuit information. The format is described here: https://quantumai.google/qsim/input_format. :param contents: The full string of the qsim file. :type contents: str :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. :rtype: dict .. py:function:: parse_qsim_file(fname, **kwargs) Parse a qsim file. .. py:function:: parse_qsim_url(url, **kwargs) Parse a qsim url. .. py:function:: to_clean_list(s, delimiter) Split, strip and filter a string by a given character into a list. .. py:function:: multi_replace(s, replacements) Replace multiple substrings in a string. .. py:function:: get_openqasm2_regexes() .. py:function:: parse_openqasm2_str(contents) 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. .. py:function:: parse_openqasm2_file(fname, **kwargs) Parse an OpenQASM 2.0 file. .. py:function:: parse_openqasm2_url(url, **kwargs) Parse an OpenQASM 2.0 url. .. py:data:: ALL_GATES .. py:data:: ONE_QUBIT_GATES .. py:data:: TWO_QUBIT_GATES .. py:data:: ALL_PARAM_GATES .. py:data:: ONE_QUBIT_PARAM_GATES .. py:data:: TWO_QUBIT_PARAM_GATES .. py:data:: GATE_TAGS .. py:data:: GATE_SIZE .. py:data:: CONSTANT_GATES .. py:data:: PARAM_GATES .. py:data:: SPECIAL_GATES .. py:function:: register_constant_gate(name, G, num_qubits, tag=None) .. py:function:: register_param_gate(name, param_fn, num_qubits, tag=None) .. py:function:: register_special_gate(name, fn, num_qubits, tag=None, array=None) .. py:function:: rx_gate_param_gen(params) .. py:function:: ry_gate_param_gen(params) .. py:function:: rz_gate_param_gen(params) .. py:function:: u3_gate_param_gen(params) .. py:function:: u2_gate_param_gen(params) .. py:function:: u1_gate_param_gen(params) .. py:function:: cu3_param_gen(params) .. py:function:: cu2_param_gen(params) .. py:function:: cu1_param_gen(params) .. py:function:: crx_param_gen(params) Parametrized controlled X-rotation. .. py:function:: cry_param_gen(params) Parametrized controlled Y-rotation. .. py:function:: crz_param_gen(params) Parametrized controlled Z-rotation. .. py:function:: fsim_param_gen(params) .. py:function:: fsimg_param_gen(params) .. py:function:: givens_param_gen(params) .. py:function:: givens2_param_gen(params) .. py:function:: xx_plus_yy_param_gen(params) .. py:function:: xx_minus_yy_param_gen(params) .. py:function:: rxx_param_gen(params) Parametrized two qubit XX-rotation. .. math:: \mathrm{RXX}(\theta) = \exp(-i \frac{\theta}{2} X_i X_j) .. py:function:: ryy_param_gen(params) Parametrized two qubit YY-rotation. .. math:: \mathrm{RYY}(\theta) = \exp(-i \frac{\theta}{2} Y_i Y_j) .. py:function:: rzz_param_gen(params) Parametrized two qubit ZZ-rotation. .. math:: \mathrm{RZZ}(\theta) = \exp(-i \frac{\theta}{2} Z_i Z_j) .. py:function:: su4_gate_param_gen(params) 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, .. py:data:: _MPS_METHODS .. py:function:: apply_swap(psi, i, j, **gate_opts) .. py:function:: build_controlled_gate_htn(ncontrol, gate, upper_inds, lower_inds, tags_each=None, tags_all=None, bond_ind=None) Build a low rank hyper tensor network (CP-decomp like) representation of a multi controlled gate. .. py:function:: _apply_controlled_gate_mps(psi, gate, tags=None, **gate_opts) Apply a multi-controlled gate to a state represented as an MPS. .. py:function:: _apply_controlled_gate_htn(psi, gate, tags=None, propagate_tags='register', **gate_opts) .. py:function:: apply_controlled_gate(psi, gate, tags=None, contract='auto-split-gate', propagate_tags='register', **gate_opts) .. py:function:: _cached_param_gate_build(fn, params) .. py:class:: Gate(label, params, qubits=None, controls=None, round=None, parametrize=False) A simple class for storing the details of a quantum circuit gate. :param label: The name or 'identifier' of the gate. :type label: str :param params: The parameters of the gate. :type params: Iterable[float] :param qubits: Which qubits the gate acts on. :type qubits: Iterable[int], optional :param controls: Which qubits are the controls. :type controls: Iterable[int], optional :param round: If given, which round or layer the gate is part of. :type round: int, optional :param parametrize: Whether the gate will correspond to a parametrized tensor. :type parametrize: bool, optional .. py:attribute:: __slots__ :value: ('_label', '_params', '_qubits', '_controls', '_round', '_parametrize', '_tag', '_special',... .. py:attribute:: _label .. py:attribute:: _params .. py:attribute:: _round :value: 0 .. py:attribute:: _parametrize :value: False .. py:attribute:: _tag .. py:attribute:: _special .. py:attribute:: _constant .. py:attribute:: _array :value: None .. py:method:: from_raw(U, qubits=None, controls=None, round=None) :classmethod: .. py:method:: copy() .. py:property:: label .. py:property:: params .. py:property:: qubits .. py:property:: total_qubit_count .. py:property:: controls .. py:property:: round .. py:property:: special .. py:property:: parametrize .. py:property:: tag .. py:method:: copy_with(**kwargs) Take a copy of this gate but with some attributes changed. .. py:method:: build_array() Build the array representation of the gate. For controlled gates this *excludes* the control qubits. .. py:property:: array .. py:method:: build_mpo(L=None, **kwargs) Build an MPO representation of this gate. .. py:method:: __repr__() .. py:function:: sample_bitstring_from_prob_ndarray(p, seed=None) Sample a bitstring from n-dimensional tensor ``p`` of probabilities. .. rubric:: 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' .. py:function:: rehearsal_dict(tn, tree) .. py:function:: parse_to_gate(gate_id, *gate_args, params=None, qubits=None, controls=None, gate_round=None, parametrize=None) Map all types of gate specification into a `Gate` object. .. py:class:: 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, gate_tag_id='GATE_{}', tag_gate_rounds=True, round_tag_id='ROUND_{}', tag_gate_labels=True, bra_site_ind_id='b{}', dtype=None, to_backend=None, convert_eager=False) Class for simulating quantum circuits using tensor networks. The class keeps a list of :class:`Gate` objects in sync with a tensor network representing the current state of the circuit. :param N: The number of qubits. :type N: int, optional :param psi0: The initial state, assumed to be ``|00000....0>`` if not given. The state is always copied and the tag ``PSI0`` added. :type psi0: TensorNetwork1DVector, optional :param gate_opts: Default keyword arguments to supply to each :func:`~quimb.tensor.tensor_1d.gate_TN_1D` call during the circuit. :type gate_opts: dict_like, optional :param gate_contract: Shortcut for setting the default `'contract'` option in `gate_opts`. :type gate_contract: str, optional :param gate_propagate_tags: Shortcut for setting the default `'propagate_tags'` option in `gate_opts`. :type gate_propagate_tags: str, optional :param tags: 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``). :type tags: str or sequence of str, optional :param psi0_dtype: Ensure the initial state has this dtype. :type psi0_dtype: str, optional :param psi0_tag: Ensure the initial state has this tag. :type psi0_tag: str, optional :param tag_gate_numbers: Whether to tag each gate tensor with its number in the circuit, like ``"GATE_{g}"``. This is required for updating the circuit parameters. :type tag_gate_numbers: bool, optional :param gate_tag_id: The format string for tagging each gate tensor, by default e.g. ``"GATE_{g}"``. :type gate_tag_id: str, optional :param tag_gate_rounds: Whether to tag each gate tensor with its number in the circuit, like ``"ROUND_{r}"``. :type tag_gate_rounds: bool, optional :param round_tag_id: The format string for tagging each round of gates, by default e.g. ``"ROUND_{r}"``. :type round_tag_id: str, optional :param tag_gate_labels: Whether to tag each gate tensor with its gate type label, e.g. ``{"X_1/2", "ISWAP", "CCX", ...}``.. :type tag_gate_labels: bool, optional :param bra_site_ind_id: Use this to label 'bra' site indices when creating certain (mostly internal) intermediate tensor networks. :type bra_site_ind_id: str, optional :param dtype: A default dtype to perform calculations in. Depending on `convert_eager`, this is enforced *after* circuit construction and simplification (the default for exact simulation), or eagerly to the initial state and as gates are applied (the default for MPS simulation). :type dtype: str, optional :param to_backend: If given, apply this function to both the initial state arrays and to every gate as it is applied. :type to_backend: callable, optional :param convert_eager: Whether to eagerly perform dtype casting and application of `to_backend` as gates are supplied, or wait until after the necessary TNs for a particular task such as sampling are formed and simplified. Deferred conversion (`convert_eager=False`) is the default mode for full contraction. :type convert_eager: bool, optional .. attribute:: psi The current circuit wavefunction as a tensor network. :type: TensorNetwork1DVector .. attribute:: uni The current circuit unitary operator as a tensor network. :type: TensorNetwork1DOperator .. attribute:: gates The gates in the circuit. :type: tuple[Gate] .. rubric:: 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 >>> 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 .. seealso:: :py:obj:`Gate` .. py:attribute:: tag_gate_numbers :value: True .. py:attribute:: tag_gate_rounds :value: True .. py:attribute:: tag_gate_labels :value: True .. py:attribute:: dtype :value: None .. py:attribute:: to_backend :value: None .. py:attribute:: convert_eager :value: False .. py:attribute:: _backend_gate_cache .. py:attribute:: gate_opts .. py:attribute:: _gates :value: [] .. py:attribute:: _ket_site_ind_id .. py:attribute:: _bra_site_ind_id :value: 'b{}' .. py:attribute:: _gate_tag_id :value: 'GATE_{}' .. py:attribute:: _round_tag_id :value: 'ROUND_{}' .. py:attribute:: _sample_n_gates :value: -1 .. py:attribute:: _storage .. py:attribute:: _sampled_conditionals .. py:method:: copy() Copy the circuit and its state. .. py:method:: _maybe_convert(obj, dtype=None) .. py:method:: apply_to_arrays(fn) Apply a function to all the arrays in the circuit. .. py:method:: get_params() Get a pytree - in this case a dict - of all the parameters in the circuit. :returns: A dictionary mapping gate numbers to their parameters. :rtype: dict[int, tuple] .. py:method:: set_params(params) Set the parameters of the circuit. :param params: A dictionary mapping gate numbers to the new parameters. :type params: dict` .. py:method:: from_qsim_str(contents, progbar=False, **circuit_opts) :classmethod: Generate a ``Circuit`` instance from a 'qsim' string. .. py:method:: from_qsim_file(fname, progbar=False, **circuit_opts) :classmethod: Generate a ``Circuit`` instance from a 'qsim' file. The qsim file format is described here: https://quantumai.google/qsim/input_format. .. py:method:: from_qsim_url(url, progbar=False, **circuit_opts) :classmethod: Generate a ``Circuit`` instance from a 'qsim' url. .. py:attribute:: from_qasm .. py:attribute:: from_qasm_file .. py:attribute:: from_qasm_url .. py:method:: from_openqasm2_str(contents, progbar=False, **circuit_opts) :classmethod: Generate a ``Circuit`` instance from an OpenQASM 2.0 string. .. py:method:: from_openqasm2_file(fname, progbar=False, **circuit_opts) :classmethod: Generate a ``Circuit`` instance from an OpenQASM 2.0 file. .. py:method:: from_openqasm2_url(url, progbar=False, **circuit_opts) :classmethod: Generate a ``Circuit`` instance from an OpenQASM 2.0 url. .. py:method:: from_gates(gates, N=None, progbar=False, **kwargs) :classmethod: Generate a ``Circuit`` instance from a sequence of gates. :param gates: The sequence of gates to apply. :type gates: sequence[Gate] or sequence[tuple] :param N: The number of qubits. If not given, will be inferred from the gates. :type N: int, optional :param progbar: Whether to show a progress bar. :type progbar: bool, optional :param kwargs: Supplied to the ``Circuit`` constructor. .. py:property:: gates .. py:property:: num_gates .. py:method:: ket_site_ind(i) Get the site index for the given qubit. .. py:method:: bra_site_ind(i) Get the 'bra' site index for the given qubit, if forming an operator. .. py:method:: gate_tag(g) Get the tag for the given gate, indexed linearly. .. py:method:: round_tag(r) Get the tag for the given round (/layer). .. py:method:: _init_state(N, dtype='complex128') .. py:method:: _apply_gate(gate, tags=None, **gate_opts) Apply a ``Gate`` to this ``Circuit``. This is the main method that all calls to apply a gate should go through. :param gate: The gate to apply. :type gate: Gate :param tags: Tags to add to the gate tensor(s). :type tags: str or sequence of str, optional .. py:method:: apply_gate(gate_id, *gate_args, params=None, qubits=None, controls=None, gate_round=None, parametrize=None, **gate_opts) 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) :param gate_id: 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 case ``gate_args`` should be supplied with ``(*params, *qubits)``. - A raw array, in which case ``gate_args`` should be supplied with ``(*qubits,)``. :type gate_id: Gate, str, or array_like :param gate_args: The arguments to supply to it. :type gate_args: list[str] :param gate_round: The gate round. If ``gate_id`` is integer-like, will also be taken from here, with then ``gate_id, gate_args = gate_args[0], gate_args[1:]``. :type gate_round: int, optional :param gate_opts: Supplied to the gate function, options here will override the default ``gate_opts``. .. py:method:: apply_gate_raw(U, where, controls=None, gate_round=None, **gate_opts) Apply the raw array ``U`` as a gate on qubits in ``where``. It will be assumed to be unitary for the sake of computing reverse lightcones. .. py:method:: apply_gates(gates, progbar=False, **gate_opts) Apply a sequence of gates to this tensor network quantum circuit. :param gates: The sequence of gates to apply. :type gates: Sequence[Gate] or Sequence[Tuple] :param gate_opts: Supplied to :meth:`~quimb.tensor.circuit.Circuit.apply_gate`. .. py:method:: h(i, gate_round=None, **kwargs) .. py:method:: x(i, gate_round=None, **kwargs) .. py:method:: y(i, gate_round=None, **kwargs) .. py:method:: z(i, gate_round=None, **kwargs) .. py:method:: s(i, gate_round=None, **kwargs) .. py:method:: sdg(i, gate_round=None, **kwargs) .. py:method:: t(i, gate_round=None, **kwargs) .. py:method:: tdg(i, gate_round=None, **kwargs) .. py:method:: sx(i, gate_round=None, **kwargs) .. py:method:: sxdg(i, gate_round=None, **kwargs) .. py:method:: x_1_2(i, gate_round=None, **kwargs) .. py:method:: y_1_2(i, gate_round=None, **kwargs) .. py:method:: z_1_2(i, gate_round=None, **kwargs) .. py:method:: w_1_2(i, gate_round=None, **kwargs) .. py:method:: hz_1_2(i, gate_round=None, **kwargs) .. py:method:: cnot(i, j, gate_round=None, **kwargs) .. py:method:: cx(i, j, gate_round=None, **kwargs) .. py:method:: cy(i, j, gate_round=None, **kwargs) .. py:method:: cz(i, j, gate_round=None, **kwargs) .. py:method:: iswap(i, j, gate_round=None, **kwargs) .. py:method:: iden(i, gate_round=None) .. py:method:: swap(i, j, gate_round=None, **kwargs) .. py:method:: rx(theta, i, gate_round=None, parametrize=False, **kwargs) .. py:method:: ry(theta, i, gate_round=None, parametrize=False, **kwargs) .. py:method:: rz(theta, i, gate_round=None, parametrize=False, **kwargs) .. py:method:: u3(theta, phi, lamda, i, gate_round=None, parametrize=False, **kwargs) .. py:method:: u2(phi, lamda, i, gate_round=None, parametrize=False, **kwargs) .. py:method:: u1(lamda, i, gate_round=None, parametrize=False, **kwargs) .. py:method:: phase(lamda, i, gate_round=None, parametrize=False, **kwargs) .. py:method:: cu3(theta, phi, lamda, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: cu2(phi, lamda, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: cu1(lamda, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: cphase(lamda, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: fsim(theta, phi, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: fsimg(theta, zeta, chi, gamma, phi, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: givens(theta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: givens2(theta, phi, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: xx_plus_yy(theta, beta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: xx_minus_yy(theta, beta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: rxx(theta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: ryy(theta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: rzz(theta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: crx(theta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: cry(theta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: crz(theta, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: su4(theta1, phi1, lamda1, theta2, phi2, lamda2, theta3, phi3, lamda3, theta4, phi4, lamda4, t1, t2, t3, i, j, gate_round=None, parametrize=False, **kwargs) .. py:method:: ccx(i, j, k, gate_round=None, **kwargs) .. py:method:: ccnot(i, j, k, gate_round=None, **kwargs) .. py:method:: toffoli(i, j, k, gate_round=None, **kwargs) .. py:method:: ccy(i, j, k, gate_round=None, **kwargs) .. py:method:: ccz(i, j, k, gate_round=None, **kwargs) .. py:method:: cswap(i, j, k, gate_round=None, **kwargs) .. py:method:: fredkin(i, j, k, gate_round=None, **kwargs) .. py:property:: psi Tensor network representation of the wavefunction. .. py:method:: get_uni(transposed=False) Tensor network representation of the unitary operator (i.e. with the initial state removed). .. py:property:: uni .. py:method:: get_reverse_lightcone_tags(where) Get the tags of gates in this circuit corresponding to the 'reverse' lightcone propagating backwards from registers in ``where``. :param where: The register or register to get the reverse lightcone of. :type where: int or sequence of int :returns: The sequence of gate tags (``GATE_{i}``, ...) corresponding to the lightcone. :rtype: tuple[str] .. py:method:: get_psi_reverse_lightcone(where, keep_psi0=False) Get just the bit of the wavefunction in the reverse lightcone of sites in ``where`` - i.e. causally linked. :param where: The sites to propagate the the lightcone back from, supplied to :meth:`~quimb.tensor.circuit.Circuit.get_reverse_lightcone_tags`. :type where: int, or sequence of int :param keep_psi0: Keep the tensors corresponding to the initial wavefunction regardless of whether they are outside of the lightcone. :type keep_psi0: bool, optional :returns: **psi_lc** :rtype: TensorNetwork1DVector .. py:method:: clear_storage() Clear all cached data. .. py:method:: _maybe_init_storage() .. py:method:: get_psi_simplified(seq='ADCRS', atol=1e-12, equalize_norms=False) Get the full wavefunction post local tensor network simplification. :param seq: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type seq: str, optional :param atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type atol: float, optional :param equalize_norms: Actively renormalize tensor norms during simplification. :type equalize_norms: bool, optional :returns: **psi** :rtype: TensorNetwork1DVector .. py:method:: get_rdm_lightcone_simplified(where, seq='ADCRS', atol=1e-12, equalize_norms=False) Get a simplified TN of the norm of the wavefunction, with gates outside reverse lightcone of ``where`` cancelled, and physical indices within ``where`` preserved so that they can be fixed (sliced) or used as output indices. :param where: The region assumed to be the target density matrix essentially. Supplied to :meth:`~quimb.tensor.circuit.Circuit.get_reverse_lightcone_tags`. :type where: int or sequence of int :param seq: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type seq: str, optional :param atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type atol: float, optional :param equalize_norms: Actively renormalize tensor norms during simplification. :type equalize_norms: bool, optional :rtype: TensorNetwork .. py:method:: amplitude(b, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-12, simplify_equalize_norms=True, backend=None, dtype=None, rehearse=False) Get the amplitude coefficient of bitstring ``b``. .. math:: c_b = \langle b | \psi \rangle :param b: The bitstring to compute the transition amplitude for. :type b: str or sequence of int :param optimize: 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). :type optimize: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param backend: Backend to perform the contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param rehearse: 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. :type rehearse: bool or "tn", optional .. py:method:: amplitude_rehearse(b='random', simplify_sequence='ADCRS', simplify_atol=1e-12, simplify_equalize_norms=True, optimize='auto-hq', dtype=None, rehearse=True) 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. :param b: The bitstring to rehearse computing the transition amplitude for, if ``'random'`` (the default) a random bitstring will be used. :type b: 'random', str or sequence of int :param optimize: 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). :type optimize: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param backend: Backend to perform the marginal contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :rtype: dict .. py:attribute:: amplitude_tn .. py:method:: partial_trace(keep, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-12, simplify_equalize_norms=True, backend=None, dtype=None, rehearse=False) Perform the partial trace on the circuit wavefunction, retaining only qubits in ``keep``, and making use of reverse lightcone cancellation: .. math:: \rho_{\bar{q}} = Tr_{\bar{p}} |\psi_{\bar{q}} \rangle \langle \psi_{\bar{q}}| Where :math:`\bar{q}` is the set of qubits to keep, :math:`\psi_{\bar{q}}` is the circuit wavefunction only with gates in the causal cone of this set, and :math:`\bar{p}` is the remaining qubits. :param keep: The qubit(s) to keep as we trace out the rest. :type keep: int or sequence of int :param optimize: 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). :type optimize: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param backend: Backend to perform the marginal contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param rehearse: 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. :type rehearse: bool or "tn", optional :rtype: array or dict .. py:attribute:: partial_trace_rehearse .. py:attribute:: partial_trace_tn .. py:method:: local_expectation(G, where, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-12, simplify_equalize_norms=True, backend=None, dtype=None, rehearse=False) Compute the a single expectation value of operator ``G``, acting on sites ``where``, making use of reverse lightcone cancellation. .. math:: \langle \psi_{\bar{q}} | G_{\bar{q}} | \psi_{\bar{q}} \rangle where :math:`\bar{q}` is the set of qubits :math:`G` acts one and :math:`\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. :param G: The raw operator(s) to find the expectation of. :type G: array or sequence[array] :param where: Which qubits the operator acts on. :type where: int or sequence of int :param optimize: 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). :type optimize: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param backend: Backend to perform the marginal contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param gate_opts: Options to use when applying ``G`` to the wavefunction. :type gate_opts: None or dict_like :param rehearse: 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. :type rehearse: bool or "tn", optional :rtype: scalar, tuple[scalar] or dict .. py:attribute:: local_expectation_rehearse .. py:attribute:: local_expectation_tn .. py:method:: 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) Compute the probability tensor of qubits in ``where``, given possibly fixed qubits in ``fix`` and tracing everything else having removed redundant unitary gates. :param where: The qubits to compute the marginal probability distribution of. :type where: sequence of int :param fix: Measurement results on other qubits to fix. :type fix: None or dict[int, str], optional :param optimize: 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). :type optimize: str, optional :param backend: Backend to perform the marginal contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param rehearse: Whether to perform the marginal contraction or just return the associated TN and contraction tree. :type rehearse: bool or "tn", optional .. py:attribute:: compute_marginal_rehearse .. py:attribute:: compute_marginal_tn .. py:method:: calc_qubit_ordering(qubits=None, method='greedy-lightcone') Get a order to measure ``qubits`` in, by greedily choosing whichever has the smallest reverse lightcone followed by whichever expands this lightcone *least*. :param qubits: The qubits to generate a lightcone ordering for, if ``None``, assume all qubits. :type qubits: None or sequence of int :returns: The order to 'measure' qubits in. :rtype: tuple[int] .. py:method:: _parse_qubits_order(qubits=None, order=None) Simply initializes the default of measuring all qubits, and the default order, or checks that ``order`` is a permutation of ``qubits``. .. py:method:: _group_order(order, group_size=1) Take the qubit ordering ``order`` and batch it in groups of size ``group_size``, sorting the qubits (for caching reasons) within each group. .. py:method:: get_qubit_distances(method='dijkstra', alpha=2) Get a nested dictionary of qubit distances. This is computed from a graph representing qubit interactions. The graph has an edge between qubits if they are acted on by the same gate, and the distance-weight of the edge is exponentially small in the number of gates between them. :param method: The method to use to compute the qubit distances. See :func:`networkx.all_pairs_dijkstra_path_length` and :func:`networkx.resistance_distance`. :type method: {'dijkstra', 'resistance'}, optional :param alpha: The distance weight between qubits is ``alpha**(num_gates - 1 )``. :type alpha: float, optional :returns: The distance between each pair of qubits, accessed like ``distances[q1][q2]``. If two qubits are not connected, the distance is missing. :rtype: dict[int, dict[int, float]] .. py:method:: reordered_gates_dfs_clustered() Get the gates reordered by a depth first search traversal of the multi-qubit gate graph that greedily selects successive gates which are 'close' in graph distance, and shifts single qubit gates to be adjacent to multi-qubit gates where possible. .. py:method:: sample(C, qubits=None, order=None, group_size=10, max_marginal_storage=2**20, seed=None, optimize='auto-hq', backend=None, dtype='complex64', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=True) 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, :math:`\{q_0, q_1, q_2, q_3, \ldots\}` we first compute: .. math:: 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 :math:`\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 :math:`r_0` we can then move on to the next marginal: .. math:: 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 :math:`r_1`. Then we compute: .. math:: 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 .. math:: |\langle r_0 r_1 r_2 r_3 \ldots | \psi \rangle|^2 since there is nothing left to trace out. :param C: The number of times to sample. :type C: int :param qubits: Which qubits to measure, defaults (``None``) to all qubits. :type qubits: None or sequence of int, optional :param order: 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 of ``qubits``. :type order: None or sequence of int, optional :param group_size: 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``. :type group_size: int, optional :param max_marginal_storage: The total cumulative number of marginal probabilites to cache, once this is exceeded caching will be turned off. :type max_marginal_storage: int, optional :param seed: A random seed, passed to ``numpy.random.seed`` if given. :type seed: None or int, optional :param optimize: Contraction path optimizer to use for the marginals, shouldn't be a non-reusable path optimizer as called on many different TNs. Passed to :func:`cotengra.array_contract_tree`. :type optimize: str, optional :param backend: Backend to perform the marginal contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :Yields: **bitstrings** (*sequence of str*) .. py:method:: sample_rehearse(qubits=None, order=None, group_size=10, result=None, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=True, rehearse=True, progbar=False) Perform the preparations and contraction tree findings for :meth:`~quimb.tensor.circuit.Circuit.sample`, caching various intermedidate objects, but don't perform the main contractions. :param qubits: Which qubits to measure, defaults (``None``) to all qubits. :type qubits: None or sequence of int, optional :param order: Which order to measure the qubits in, defaults (``None``) to an order based on greedily expanding the smallest reverse lightcone. :type order: None or sequence of int, optional :param group_size: 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``. :type group_size: int, optional :param result: Explicitly check the computational cost of this result, assumed to be all zeros if not given. :type result: None or dict[int, str], optional :param optimize: Contraction path optimizer to use for the marginals, shouldn't be a non-reusable path optimizer as called on many different TNs. Passed to :func:`cotengra.array_contract_tree`. :type optimize: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param progbar: Whether to show the progress of finding each contraction tree. :type progbar: bool, optional :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'). :rtype: dict[tuple[int], dict] .. py:attribute:: sample_tns .. py:method:: sample_chaotic(C, marginal_qubits, fix=None, max_marginal_storage=2**20, seed=None, optimize='auto-hq', backend=None, dtype='complex64', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=True) 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 :math:`r_0 r_1 r_2 \ldots r_{N - 6}` for the remaining :math:`N - 5` qubits will be chosen, then the final marginal will be computed as .. math:: 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. :param C: The number of times to sample. :type C: int :param marginal_qubits: 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]``. :type marginal_qubits: int or sequence of int :param fix: Measurement results on other qubits to fix. These will be randomly sampled if ``fix`` is not given or a qubit is missing. :type fix: None or dict[int, str], optional :param seed: A random seed, passed to ``numpy.random.seed`` if given. :type seed: None or int, optional :param optimize: 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). :type optimize: str, optional :param backend: Backend to perform the marginal contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :Yields: *str* .. py:method:: sample_chaotic_rehearse(marginal_qubits, result=None, optimize='auto-hq', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=True, dtype='complex64', rehearse=True) Rehearse chaotic sampling (perform just the TN simplifications and contraction tree finding). :param marginal_qubits: 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]``. :type marginal_qubits: int or sequence of int :param result: Explicitly check the computational cost of this result, assumed to be all zeros if not given. :type result: None or dict[int, str], optional :param optimize: 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). :type optimize: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :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. :rtype: dict[tuple[int], dict] .. py:attribute:: sample_chaotic_tn .. py:method:: get_gate_by_gate_circuits(group_size=10) Get a sequence of circuits by partitioning the gates into groups such circuit `i + 1` acts on at most ``group_size`` new qubits compared to circuit `i`. :param group_size: The maximum number of new qubits that can be acted on by a circuit compared to its predecessor. :type group_size: int, optional :returns: A sequence of dicts, each with keys ``'circuit'`` and ``'where'``, where the former is a :class:`~quimb.tensor.circuit.Circuit` and the latter the tuple of new qubits that it acts on comparaed to the previous circuit. :rtype: Sequence[dict] .. py:method:: sample_gate_by_gate(C, group_size=10, seed=None, max_marginal_storage=2**20, optimize='auto-hq', backend=None, dtype='complex64', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=True) Sample this circuit using the gate-by-gate method, where we 'evolve' a result bitstring by sequentially including more and more gates, at each step updating the result by computing a full conditional marginal. See "How to simulate quantum measurement without computing marginals" by Sergey Bravyi, David Gosset, Yinchen Liu (https://arxiv.org/abs/2112.08499). The overall complexity of this is guaranteed to be similar to that of computing a single amplitude which can be much better than the naive "qubit-by-qubit" (`.sample`) method. However, it requires evaluting a number of tensor networks that scales linearly with the number of gates which can offset any practical advantages for shallow circuits for example. :param C: The number of samples to generate. :type C: int :param group_size: The maximum number of qubits that can be acted on by a circuit compared to its predecessor. This will be the dimension of the marginal computed at each step. :type group_size: int, optional :param seed: A random seed, passed to ``numpy.random.seed`` if given. :type seed: None or int, optional :param max_marginal_storage: The total cumulative number of marginal probabilites to cache, once this is exceeded caching will be turned off. :type max_marginal_storage: int, optional :param optimize: Contraction path optimizer to use for the marginals, shouldn't be a non-reusable path optimizer as called on many different TNs. Passed to :func:`cotengra.array_contract_tree`. :type optimize: str, optional :param backend: Backend to perform the marginal contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param rehearse: 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. :type rehearse: bool, optional :Yields: *str* .. py:method:: sample_gate_by_gate_rehearse(group_size=10, optimize='auto-hq', dtype='complex64', simplify_sequence='ADCRS', simplify_atol=1e-06, simplify_equalize_norms=True, rehearse=True, progbar=False) Perform the preparations and contraction tree findings for :meth:`~quimb.tensor.circuit.Circuit.sample_gate_by_gate`, caching various intermedidate objects, but don't perform the main contractions. :param group_size: The maximum number of qubits that can be acted on by a circuit compared to its predecessor. This will be the dimension of the marginal computed at each step. :type group_size: int, optional :param optimize: Contraction path optimizer to use for the marginals, shouldn't be a non-reusable path optimizer as called on many different TNs. Passed to :func:`cotengra.array_contract_tree`. :type optimize: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param rehearse: If ``True``, generate and cache the simplified tensor network and contraction tree but don't actually perform the contraction. If "tn", only generate the simplified tensor networks. :type rehearse: True or "tn", optional :rtype: Sequence[dict] or Sequence[TensorNetwork] .. py:attribute:: sample_gate_by_gate_tns .. py:method:: to_dense(reverse=False, optimize='auto-hq', simplify_sequence='R', simplify_atol=1e-12, simplify_equalize_norms=True, backend=None, dtype=None, rehearse=False) Generate the dense representation of the final wavefunction. :param reverse: Whether to reverse the order of the subsystems, to match the convention of qiskit for example. :type reverse: bool, optional :param optimize: 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). :type optimize: str, optional :param dtype: If given, convert the tensors to this dtype prior to contraction. :type dtype: dtype or str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param simplify_equalize_norms: Actively renormalize tensor norms during simplification. :type simplify_equalize_norms: bool, optional :param backend: Backend to perform the contraction with, e.g. ``'numpy'``, ``'cupy'`` or ``'jax'``. Passed to ``cotengra``. :type backend: str, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param rehearse: 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. :type rehearse: bool, optional :returns: **psi** -- The densely represented wavefunction with ``dtype`` data. :rtype: qarray .. py:attribute:: to_dense_rehearse .. py:attribute:: to_dense_tn .. py:method:: simulate_counts(C, seed=None, reverse=False, **to_dense_opts) Simulate measuring all qubits in the computational basis many times. Unlike :meth:`~quimb.tensor.circuit.Circuit.sample`, this generates all the samples simultaneously using the full wavefunction constructed from :meth:`~quimb.tensor.circuit.Circuit.to_dense`, then calling :func:`~quimb.calc.simulate_counts`. .. warning:: Because this constructs the full wavefunction it always requires exponential memory in the number of qubits, regardless of circuit depth and structure. :param C: The number of 'experimental runs', i.e. total counts. :type C: int :param seed: A seed for reproducibility. :type seed: int, optional :param reverse: Whether to reverse the order of the subsystems, to match the convention of qiskit for example. :type reverse: bool, optional :param to_dense_opts: Suppled to :meth:`~quimb.tensor.circuit.Circuit.to_dense`. :returns: **results** -- The number of recorded counts for each :rtype: dict[str, int] .. py:method:: schrodinger_contract(*args, **contract_opts) .. py:method:: xeb(samples_or_counts, cache=None, cache_maxsize=2**20, progbar=False, **amplitude_opts) Compute the linear cross entropy benchmark (XEB) for samples or counts, amplitude per amplitude. :param samples_or_counts: Either the raw bitstring samples or a dict mapping bitstrings to the number of counts observed. :type samples_or_counts: Iterable[str] or Dict[str, int] :param cache: A dictionary to store the probabilities in, if not supplied ``quimb.utils.LRU(cache_maxsize)`` will be used. :type cache: dict, optional :param cache_maxsize: The maximum size of the cache to be used. :param optional: The maximum size of the cache to be used. :param progbar: Whether to show progress as the bitstrings are iterated over. :param optional: Whether to show progress as the bitstrings are iterated over. :param amplitude_opts: Supplied to :meth:`~quimb.tensor.circuit.Circuit.amplitude`. .. py:method:: xeb_ex(optimize='auto-hq', simplify_sequence='R', simplify_atol=1e-12, simplify_equalize_norms=True, dtype=None, backend=None, autojit=False, progbar=False, **contract_opts) 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. :param optimize: Contraction path optimizer. :type optimize: str or PathOptimizer, optional :param simplify_sequence: Simplifications to apply to tensor network prior to contraction. :type simplify_sequence: str, optional :param simplify_sequence: Which local tensor network simplifications to perform and in which order, see :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_sequence: str, optional :param simplify_atol: The tolerance with which to compare to zero when applying :meth:`~quimb.tensor.tensor_core.TensorNetwork.full_simplify`. :type simplify_atol: float, optional :param dtype: Data type to cast the TN to before contraction. :type dtype: str, optional :param backend: Convert tensors to, and then use contractions from, this library. :type backend: str, optional :param autojit: Apply ``autoray.autojit`` to the contraciton and map-reduce. :type autojit: bool, optional :param progbar: Show progress in terms of number of wavefunction chunks processed. :type progbar: bool, optional .. py:method:: update_params_from(tn) Assuming ``tn`` is a tensor network with tensors tagged ``GATE_{i}`` corresponding to this circuit (e.g. from ``circ.psi`` or ``circ.uni``) but with updated parameters, update the current circuit parameters and tensors with those values. This is an inplace modification of the ``Circuit``. :param tn: The tensor network to find the updated parameters from. :type tn: TensorNetwork .. py:method:: draw(figsize=None, radius=1 / 3, drawcolor=(0.5, 0.5, 0.5), linewidth=1) Draw a simple linear schematic of the circuit. :param figsize: The size of the figure, if not given will be set based on the number of gates and qubits. :type figsize: tuple, optional :param radius: The radius of the gates. :type radius: float, optional :param drawcolor: The color of the wires. :type drawcolor: tuple, optional :param linewidth: The linewidth of the wires. :type linewidth: float, optional :returns: * **fig** (*matplotlib.Figure*) -- The figure object. * **ax** (*matplotlib.Axes*) -- The axis object. .. py:method:: __repr__() .. py:class:: CircuitMPS(N=None, *, psi0=None, max_bond=None, cutoff=1e-10, gate_opts=None, gate_contract='auto-mps', dtype=None, to_backend=None, convert_eager=True, **circuit_opts) Bases: :py:obj:`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. :param N: The number of qubits in the circuit. :type N: int, optional :param psi0: The initial state, assumed to be ``|00000....0>`` if not given. The state is always copied and the tag ``PSI0`` added. :type psi0: TensorNetwork1DVector, optional :param max_bond: The maximum bond dimension to truncate to when applying gates, if any. This is simply a shortcut for setting ``gate_opts['max_bond']``. :type max_bond: int, optional :param cutoff: The singular value cutoff to use when truncating the state. This is simply a shortcut for setting ``gate_opts['cutoff']``. :type cutoff: float, optional :param gate_opts: Default options to pass to each gate, for example, "max_bond" and "cutoff" etc. :type gate_opts: dict, optional :param gate_contract: 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. See :func:`tensor_network_1d_compress`. :type gate_contract: str, optional :param dtype: The data type to use for the state tensor. :type dtype: str, optional :param to_backend: A function to convert tensor data to a particular backend. :type to_backend: callable, optional :param convert_eager: Whether to eagerly perform dtype casting and application of `to_backend` as gates are supplied, or wait until after the necessary TNs for a particular task such as sampling are formed and simplified. Eager conversion (`convert_eager=True`) is the default mode for MPS simulation, unlike full contraction. :type convert_eager: bool, optional :param circuit_opts: Supplied to :class:`~quimb.tensor.circuit.Circuit`. .. attribute:: psi The current state of the circuit, always in MPS form. :type: MatrixProductState .. rubric:: 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, ) ) .. py:method:: _init_state(N, dtype='complex128') .. py:method:: apply_gates(gates, progbar=False, **gate_opts) Apply a sequence of gates to this tensor network quantum circuit. :param gates: The sequence of gates to apply. :type gates: Sequence[Gate] or Sequence[Tuple] :param gate_opts: Supplied to :meth:`~quimb.tensor.circuit.Circuit.apply_gate`. .. py:property:: psi Tensor network representation of the wavefunction. .. py:property:: uni .. py:method:: calc_qubit_ordering(qubits=None) MPS already has a natural ordering. .. py:method:: get_psi_reverse_lightcone(where, keep_psi0=False) Override ``get_psi_reverse_lightcone`` as for an MPS the lightcone is not meaningful. .. py:method:: sample(C, seed=None, dtype=None, *, qubits=None, order=None, group_size=None, max_marginal_storage=None, optimize=None, backend=None, simplify_sequence=None, simplify_atol=None, simplify_equalize_norms=None) Sample the MPS circuit ``C`` times. :param C: The number of samples to generate. :type C: int :param seed: A random seed or generator to use for reproducibility. :type seed: None, int, or generator, optional .. py:method:: fidelity_estimate() Estimate the fidelity of the current state based on its norm, which tracks how much the state has been truncated: .. math:: \tilde{F} = \left| \langle \psi | \psi \rangle \right|^2 \approx \left|\langle \psi_\mathrm{ideal} | \psi \rangle\right|^2 .. seealso:: :py:obj:`error_estimate` .. py:method:: error_estimate() Estimate the error in the current state based on the norm of the discarded part of the state: .. math:: \epsilon = 1 - \tilde{F} .. seealso:: :py:obj:`fidelity_estimate` .. py:method:: local_expectation(G, where, normalized=False, dtype=None, *, simplify_sequence=None, simplify_atol=None, simplify_equalize_norms=None, backend=None, rehearse=None, **contract_opts) Compute the local expectation value of a local operator at ``where`` (via forming the reduced density matrix). Note this moves the orthogonality around inplace, and records it in `info`. :param G: The local operator tensor. :type G: Tensor :param where: The qubit to compute the expectation value at. :type where: int :param normalized: Whether to normalize the expectation value by the norm of the state. :type normalized: bool, optional :param dtype: If given, ensure the TN is cast to this dtype before contracting. :type dtype: dtype, optional :rtype: float .. py:class:: CircuitPermMPS(N=None, psi0=None, gate_opts=None, gate_contract='swap+split', **circuit_opts) Bases: :py:obj:`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``. The ``psi`` 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. .. py:attribute:: qubits .. py:method:: _apply_gate(gate, tags=None, **gate_opts) Apply a ``Gate`` to this ``Circuit``. This is the main method that all calls to apply a gate should go through. :param gate: The gate to apply. :type gate: Gate :param tags: Tags to add to the gate tensor(s). :type tags: str or sequence of str, optional .. py:method:: calc_qubit_ordering(qubits=None) Given by the current qubit permutation. .. py:method:: get_psi_unordered() Return the MPS representing the state but without reordering the sites. .. py:method:: sample(C, seed=None, dtype=None) Sample the PermMPS circuit ``C`` times. :param C: The number of samples to generate. :type C: int :param seed: A random seed or generator to use for reproducibility. :type seed: None, int, or generator, optional :Yields: *str* -- The next sample bitstring. .. py:property:: psi Tensor network representation of the wavefunction. .. py:class:: CircuitDense(N=None, psi0=None, gate_opts=None, gate_contract=True, tags=None, convert_eager=True, **circuit_opts) Bases: :py:obj:`Circuit` Quantum circuit simulation keeping the state in full dense form. .. py:property:: psi Tensor network representation of the wavefunction. .. py:property:: uni .. py:method:: calc_qubit_ordering(qubits=None) Qubit ordering doesn't matter for a dense wavefunction. .. py:method:: get_psi_reverse_lightcone(where, keep_psi0=False) Override ``get_psi_reverse_lightcone`` as for a dense wavefunction the lightcone is not meaningful.