filter_functions.pulse_sequence module
This module defines the PulseSequence class, the package’s central object.
Classes
PulseSequenceThe pulse sequence defined by a Hamiltonian
Functions
concatenate()Function to concatenate different
PulseSequenceinstances and efficiently compute their joint filter functionconcatenate_periodic()Function to more efficiently concatenate many versions of the same
PulseSequenceinstances and compute their joint filter functionextend()Function to map several
PulseSequenceinstances to different qubits, efficiently scaling up cached attributes.
- class PulseSequence(*args, **kwargs)[source]
Bases:
objectA class for pulse sequences and their filter functions.
The Hamiltonian is separated into a control and a noise part with
\[\begin{split}\mathcal{H}_c &= \sum_i a_i(t) A_i \\ \mathcal{H}_n &= \sum_j s_j(t) b_j(t) B_j\end{split}\]where \(A_i\) and \(B_j\) are hermitian operators and \(b_j(t)\) are classically fluctuating noise variables captured in a power spectral density and not needed at instantiation of a
PulseSequence.- Parameters:
- H_c: list of lists
A nested list of n_cops nested lists as taken by QuTiP functions (see for example
qutip.propagator.propagator()) describing the control part of the Hamiltonian. The i-th entry of the list should be a list consisting of the i-th operator \(A_i\) making up the control Hamiltonian and a list or array \(a_i(t)\) describing the magnitude of that operator during the time intervals dt. Optionally, the list may also include operator identifiers. That is, H_c should look something like this:H = [[c_oper1, c_coeff1, c_oper_identifier1], [c_oper2, c_coeff2, c_oper_identifier2], ...]
The operators may be given either as NumPy arrays or QuTiP Qobjs and each coefficient array should have the same number of elements as dt, and should be given in units of \(\hbar\). If not every sublist (read: operator) was given a identifier, they are automatically filled up with ‘A_i’ where i is the position of the operator.
- H_n: list of lists
A nested list of n_nops nested lists as taken by QuTiP functions (see for example
qutip.propagator.propagator()) describing the noise part of the Hamiltonian. The j-th entry of the list should be a list consisting of the j-th operator \(B_j\) making up the noise Hamiltonian and a list or array describing the sensitivity \(s_j(t)\) of the system to the noise operator during the time intervals dt. Optionally, the list may also include operator identifiers. That is, H_n should look something like this:H = [[n_oper1, n_coeff1, n_oper_identifier1], [n_oper2, n_coeff2, n_oper_identifier2], ...]
The operators may be given either as NumPy arrays or QuTiP Qobjs and each coefficient array should have the same number of elements as dt, and should be given in units of \(\hbar\). If not every sublist (read: operator) was given a identifier, they are automatically filled up with ‘B_i’ where i is the position of the operator.
- dt: array_like, shape (n_dt,)
The segment durations of the Hamiltonian (i.e. durations of constant control). Internally, the control operation is taken to start at \(t_0\equiv 0\), i.e. the edges of the constant control segments are at times
t = [0, *np.cumsum(dt)].- basis: Basis, shape (d**2, d, d), optional
The operator basis in which to calculate. If a Generalized Gell-Mann basis (see
ggm()) is chosen, some calculations will be faster for large dimensions due to a simpler basis expansion. However, when extending the pulse sequence to larger qubit registers, cached filter functions cannot be retained since the GGM basis does not factor into tensor products. In this case a Pauli basis is preferable.
- Attributes:
- c_opers: ndarray, shape (n_cops, d, d)
Control operators. Note that they are stored sorted by their corresponding identifiers.
- n_opers: ndarray, shape (n_nops, d, d)
Noise operators. Note that they are stored sorted by their corresponding identifiers.
- c_oper_identifers: sequence of str
Identifiers for the control operators of the system. Stored sorted.
- n_oper_identifers: sequence of str
Identifiers for the noise operators of the system. Stored sorted.
- c_coeffs: ndarray, shape (n_cops, n_dt)
Control parameters in units of \(\hbar\). Note that they are stored sorted by their corresponding identifiers.
- n_coeffs: ndarray, shape (n_nops, n_dt)
Noise sensitivities in units of \(\hbar\). Note that they are stored sorted by their corresponding identifiers.
- dt: ndarray, shape (n_dt,)
Time steps
- t: ndarray, shape (n_dt + 1,)
Absolute times taken to start at \(t_0\equiv 0\)
- tau: float
Total duration. Equal to t[-1].
- d: int
Dimension of the Hamiltonian
- basis: Basis, shape (d**2, d, d)
The operator basis used for calculation
- nbytes: int
An estimate of the memory consumed by the PulseSequence instance and its attributes
- If the Hamiltonian was diagonalized, the eigenvalues and -vectors as
- well as the cumulative propagators are cached:
- eigvals: ndarray, shape (n_dt, d)
Eigenvalues \(D^{(g)}\)
- eigvecs: ndarray, shape (n_dt, d, d)
Eigenvectors \(V^{(g)}\)
- propagators: ndarray, shape (n_dt+1, d, d)
Cumulative propagators \(Q_g\)
- total_propagator: ndarray, shape (d, d)
The total propagator \(Q\) of the pulse alone. That is, \(|\psi(\tau)\rangle = Q|\psi(0)\rangle\).
- total_propagator_liouville: array_like, shape (d**2, d**2)
The transfer matrix for the total propagator of the pulse. Given by
liouville_representation(pulse.total_propagator, pulse.basis).- Furthermore, when the filter function is calculated, the frequencies
- are cached as well as other relevant quantities.
Methods
get_filter_function_derivative(omega[, ...])Calculate the pulse sequence's filter function derivative.
cleanup(method=’conservative’)
Delete cached attributes
is_cached(attr)
Checks if a given attribute of the
PulseSequenceis cacheddiagonalize()
Diagonalize the Hamiltonian of the pulse sequence, computing eigenvalues and -vectors as well as cumulative propagators
get_control_matrix(omega, show_progressbar=False)
Calculate the control matrix for frequencies omega
get_filter_function(omega, which=’fidelity’, show_progressbar=False)
Calculate the filter function for frequencies omega
get_pulse_correlation_filter_function(which=’fidelity’)
Get the pulse correlation filter function (only possible if computed during concatenation)
propagator_at_arb_t(t)
Calculate the propagator at arbitrary times
Notes
Due to the heavy use of NumPy’s
einsum()function, results have a floating point error of ~1e-13.Examples
A rotation by \(\pi\) around the axis between x and y preceeded and followed by a period of free evolution with the system subject to dephasing noise.
>>> import qutip as qt; import numpy as np >>> H_c = [[qt.sigmax(), [0, np.pi, 0]], [qt.sigmay(), [0, np.pi, 0]]] >>> # Equivalent pulse: >>> # H_c = [[qt.sigmax() + qt.sigmay(), [0, np.pi, 0]]] >>> # The noise sensitivity is constant >>> H_n = [[qt.sigmaz()/np.sqrt(2), [1, 1, 1], 'Z']] >>> dt = [1, 1, 1] >>> # Free evolution between t=0 and t=1, rotation between t=1 and t=2, >>> # and free evolution again from t=2 to t=3. >>> pulse = PulseSequence(H_c, H_n, dt) >>> pulse.c_oper_identifiers ['A_0', 'A_1'] >>> pulse.n_oper_identifiers ['Z'] >>> omega = np.logspace(-1, 2, 500) >>> F = pulse.get_filter_function(omega) # shape (1, 500) >>> # Plot the resulting filter function: >>> from filter_functions import plotting >>> fig, ax, leg = plotting.plot_filter_function(pulse)
- cache_control_matrix(omega: Sequence[float], control_matrix: ndarray | None = None, show_progressbar: bool = False, cache_intermediates: bool = False) None[source]
Cache the control matrix and the frequencies it was calculated for.
- Parameters:
- omega: array_like, shape (n_omega,)
The frequencies for which to cache the filter function.
- control_matrix: array_like, shape ([G,] n_nops, d**2, n_omega), optional
The control matrix for the frequencies omega. If
None, it is computed.- show_progressbar: bool
Show a progress bar for the calculation of the control matrix.
- cache_intermediates: bool, optional
Keep intermediate terms of the calculation that are also required by other computations. Only applies if control_matrix is not supplied.
- cache_filter_function(omega: Sequence[float], control_matrix: ndarray | None = None, filter_function: ndarray | None = None, which: str = 'fidelity', order: int = 1, show_progressbar: bool = False, cache_intermediates: bool = False, cache_second_order_cumulative: bool = False) None[source]
Cache the filter function. If control_matrix.ndim == 4, it is taken to be the ‘pulse correlation control matrix’ and summed along the first axis. In that case, also the pulse correlation filter function is calculated and cached. Total phase factors and transfer matrices of the the cumulative propagator are also cached so they can be reused during concatenation.
- Parameters:
- omega: array_like, shape (n_omega,)
The frequencies for which to cache the filter function.
- control_matrix: array_like, shape ([G,] n_nops, d**2, n_omega), optional
The control matrix for the frequencies omega. If
None, it is computed and the filter function derived from it.- filter_function: array_like, shape (n_nops, n_nops, [d**2, d**2,] n_omega), optional
The filter function for the frequencies omega. If
None, it is computed from R in the caseorder == 1and from scratch else.- which: str, optional
Which filter function to cache. Either ‘fidelity’ (default) or ‘generalized’.
- order: int, optional
First or second order filter function.
- show_progressbar: bool, optional
Show a progress bar for the calculation of the control matrix.
- cache_intermediates: bool, optional
Keep intermediate terms of the calculation that are also required by other computations.
- cache_second_order_cumulative: bool, optional
Also cache the accumulated filter function for each time step. Only if order is 2.
See also
PulseSequence.get_filter_functionGetter method
- cache_total_phases(omega: Sequence[float], total_phases: ndarray | None = None) None[source]
Cache the total phase factors for this pulse and omega.
- Parameters:
- omega: array_like, shape (n_omega,)
The frequencies for which to cache the phase factors.
- total_phases: array_like, shape (n_omega,), optional
The total phase factors for the frequencies omega. If
None, they are computed.
- cleanup(method: str = 'conservative') None[source]
Delete cached byproducts of the calculation of the filter function that are not necessarily needed anymore in order to free up memory.
- Parameters:
- method
If set to ‘conservative’ (the default), only the following attributes are deleted:
data[‘eigvals’]
data[‘eigvecs’]
data[‘propagators’]
If set to ‘greedy’, all of the above as well as the following attributes are deleted:
data[‘total_propagator’]
data[‘total_propagator_liouville’]
frequency_data[‘total_phases’]
frequency_data[‘control_matrix’]
frequency_data[‘control_matrix_pc’]
If set to ‘all’, all cached data is deleted.
If set to ‘frequency dependent’ only attributes that are functions of frequency are initalized to
None.Note that if this
PulseSequenceis concatenated with another one, some of the attributes might need to be calculated again, resulting in slower execution of the concatenation.
- property data: MappingProxyType
!! processed by numpydoc !!
- property frequency_data: MappingProxyType
!! processed by numpydoc !!
- classmethod from_arrays(c_opers: ndarray, c_oper_identifiers: ndarray, c_coeffs: ndarray, n_opers: ndarray, n_oper_identifiers: ndarray, n_coeffs: ndarray, dt: ndarray, basis: Basis | None = None)[source]
Instantiate a
PulseSequencefrom NumPy arrays.This is an alternative constructor to the nested list of lists structure that is used for initialization. See the class docstring for more information on the parameters.
- Parameters:
- c_opersndarray, shape (n_cops, d, d)
- c_oper_identifiersndarray, shape (n_cops,)
- c_coeffsndarray, shape (n_cops, n_dt)
- n_opersndarray, shape (n_nops, d, d)
- n_oper_identifiersndarray, shape (n_nops,)
- n_coeffsndarray, shape (n_nops, n_dt)
- dtndarray, shape (n_dt,)
- basisBasis, optional
- get_control_matrix(omega: Sequence[float], show_progressbar: bool = False, cache_intermediates: bool = False) ndarray[source]
Get the control matrix for the frequencies omega. If it has been cached for the same frequencies, the cached version is returned, otherwise it is calculated from scratch.
- Parameters:
- omega: array_like, shape (n_omega,)
The frequencies at which to evaluate the control matrix.
- show_progressbar: bool
Show a progress bar for the calculation of the control matrix.
- cache_intermediates: bool, optional
Keep intermediate terms of the calculation that are also required by other computations.
- Returns:
- control_matrix: ndarray, shape (n_nops, d**2, n_omega)
The control matrix for the noise operators.
- get_filter_function(omega: Sequence[float], which: str = 'fidelity', order: int = 1, show_progressbar: bool = False, cache_intermediates: bool = False, cache_second_order_cumulative: bool = False) ndarray[source]
Get the first or second order filter function.
The filter function is cached so it doesn’t need to be calculated twice for the same frequencies.
- Parameters:
- omega: array_like, shape (n_omega,)
The frequencies at which to evaluate the filter function.
- which: str, optional
Which filter function to return. Either ‘fidelity’ (default) or ‘generalized’ (see Notes). Only if
order == 1.- order: int, optional
First or second order filter function.
- show_progressbar: bool, optional
Show a progress bar for the calculation of the control matrix.
- cache_intermediates: bool, optional
Keep intermediate terms of the calculation that are also required by other computations.
- cache_second_order_cumulative: bool, optional
Also cache the accumulated filter function for each time step. Only if order is 2. This can take up a lot of memory, but is useful when interested in the time evolution. In that case, compute the filter function for the entire pulse and slice it afterwards:
pulse.cache_filter_function( omega, order=2, cache_second_order_cumulative=True ) pulses_t = [] for i in range(len(pulse)): pulses_t.append(pulse[:i])
pulses_twill contain the filter functions for each time step.
- Returns:
- filter_function: ndarray, shape (n_nops, n_nops, [d**2, d**2,] n_omega)
The filter function for each combination of noise operators as a function of omega.
Notes
The first-order generalized filter function is given by
\[F_{\alpha\beta,kl}(\omega) = \tilde{\mathcal{B}}_{\alpha k}^\ast(\omega) \tilde{\mathcal{B}}_{\beta l}(\omega),\]where \(\alpha,\beta\) are indices counting the noise operators \(B_\alpha\) and \(k,l\) indices counting the basis elements \(C_k\).
The fidelity filter function is obtained by tracing over the basis indices:
\[F_{\alpha\beta}(\omega) = \sum_{k} F_{\alpha\beta,kk}(\omega).\]
- get_filter_function_derivative(omega: Sequence[float], control_identifiers: Sequence[str] | None = None, n_oper_identifiers: Sequence[str] | None = None, n_coeffs_deriv: Sequence[Sequence[float]] | None = None) ndarray[source]
Calculate the pulse sequence’s filter function derivative.
- Parameters:
- omega: array_like, shape (n_omega,)
Frequencies at which the pulse control matrix is to be evaluated.
- control_identifiers: Sequence[str], shape (n_ctrl,)
Sequence of strings with the control identifiers to distinguish between control and drift Hamiltonian. The default is None, in which case the derivative is computed for all known non-noise operators.
- n_oper_identifiers: Sequence[str], shape (n_nops,)
Sequence of strings with the noise identifiers for which to compute the derivative contribution. The default is None, in which case it is computed for all known noise operators.
- n_coeffs_deriv: array_like, shape (n_nops, n_ctrl, n_dt)
The derivatives of the noise susceptibilities by the control amplitudes. The rows and columns should be in the same order as the corresponding identifiers above. Defaults to None, in which case the coefficients are assumed to be constant and hence their derivative vanishing.
Warning
Internally, control and noise terms of the Hamiltonian are stored alphanumerically sorted by their identifiers. If the noise and/or control identifiers above are not explicitly given, the rows and/or columns of this parameter need to be sorted in the same fashion.
- Returns:
- filter_function_deriv: ndarray, shape (n_nops, n_t, n_ctrl, n_omega)
The regular filter functions’ derivatives for variation in each control contribution. Sorted in the same fashion as n_coeffs_deriv or, if not given, alphanumerically by the identifiers.
- get_pulse_correlation_control_matrix() ndarray[source]
Get the pulse correlation control matrix if it was cached.
- get_pulse_correlation_filter_function(which: str = 'fidelity') ndarray[source]
Get the pulse correlation filter function given by
\[F_{\alpha\beta}^{(gg')}(\omega) = e^{i\omega(t_{g-1} - t_{g'-1})} \tilde{\mathcal{B}}^{(g)}(\omega)\mathcal{Q}^{(g-1)} \mathcal{Q}^{(g'-1)\dagger} \tilde{\mathcal{B}}^{(g')\dagger}(\omega),\]where \(g,g'\) index the pulse in the sequence and \(\alpha,\beta\) index the noise operators, if it was computed during concatenation. Since the calculation requires the individual pulse’s control matrices and phase factors, which are not retained after concatenation, the pulse correlation filter function cannot be computed afterwards.
Note that the frequencies for which the filter function was calculated are not stored.
- Returns:
- filter_function_pc: ndarray, shape (n_pls, n_pls, n_nops, n_nops, n_omega)
The pulse correlation filter function for each noise operator as a function of omega. The first two axes correspond to the pulses in the sequence, i.e. if the concatenated pulse sequence is \(C\circ B\circ A\), the first two axes are arranged like
\[\begin{split}F_{\alpha\beta}^{(gg')} &= \begin{pmatrix} F_{\alpha\beta}^{(AA)} & F_{\alpha\beta}^{(AB)} & F_{\alpha\beta}^{(AC)} \\ F_{\alpha\beta}^{(BA)} & F_{\alpha\beta}^{(BB)} & F_{\alpha\beta}^{(BC)} \\ F_{\alpha\beta}^{(CA)} & F_{\alpha\beta}^{(CB)} & F_{\alpha\beta}^{(CC)} \end{pmatrix}\end{split}\]for \(g,g'\in\{A, B, C\}\).
- get_total_phases(omega: Sequence[float]) ndarray[source]
Get the (cached) total phase factors for this pulse and omega.
- property intermediates: MappingProxyType
!! processed by numpydoc !!
- property nbytes: int
Return an estimate of the amount of memory consumed by this object (or, more precisely, the array attributes of this object).
- concatenate(pulses: Iterable[PulseSequence], calc_pulse_correlation_FF: bool = False, calc_filter_function: bool | None = None, calc_second_order_FF: bool | None = None, which: str = 'fidelity', omega: Sequence[float] | None = None, show_progressbar: bool = False) PulseSequence[source]
Concatenate an arbitrary number of pulses. Note that pulses are concatenated left-to-right, that is,
\[\mathtt{concatenate((A, B))} \equiv B \circ A\]so that \(A\) is executed before \(B\) when applying the concatenated pulse.
- Parameters:
- pulses: sequence of PulseSequences
The PulseSequence instances to be concatenated. If any of the instances have a cached filter function, the filter function for the composite pulse will also be calculated in order to make use of the speedup gained from concatenating the filter functions. If omega is given, calculation of the composite filter function is forced.
- calc_pulse_correlation_FF: bool, optional
Switch to control whether the pulse correlation filter function (see
PulseSequence.get_pulse_correlation_filter_function()) is calculated. If omega is not given, the cached frequencies of all pulses need to be equal.- calc_filter_function: bool, optional
Switch to force the calculation of the filter function to be carried out or not. Overrides the automatic behavior of calculating it if at least one pulse has a cached control matrix. If
Trueand no pulse has a cached control matrix, a list of frequencies must be supplied as omega. Overridden if either calc_pulse_correlation_FF or calc_second_order_FF are true.- calc_second_order_FFbool, optional
Compute the second-order filter function. Requires atomic pulses to have retained intermediates during the calculation of the control matrix.
Warning
This is an experimental feature and might have unexpected bugs.
- which: str, optional
Which filter function to compute. Either ‘fidelity’ (default) or ‘generalized’ (see
PulseSequence.get_filter_function()andPulseSequence.get_pulse_correlation_filter_function()).- omega: array_like, optional
Frequencies at which to evaluate the (pulse correlation) filter functions. If
None, an attempt is made to use cached frequencies.- show_progressbar: bool
Show a progress bar for the calculation of the control matrix.
- Returns:
- pulse: PulseSequence
The concatenated pulse.
- concatenate_periodic(pulse: PulseSequence, repeats: int, check_invertible: bool = True) PulseSequence[source]
Concatenate a pulse sequence pulse whose Hamiltonian is periodic repeats times. Although performing the same task, this function is much faster for concatenating many identical pulses with filter functions than
concatenate().Note that for large dimensions, the calculation of the control matrix using this function might be very memory intensive.
- Parameters:
- pulse: PulseSequence
The
PulseSequenceinstance to be repeated. If it has a cached filter function, the filter function for the new pulse will also be computed.- repeats: int
The number of repetitions
- check_invertiblebool
Test if the matrix inversion was successful and, where not, calculate the result ‘on foot’.
- Returns:
- newpulse: PulseSequence
The concatenated
PulseSequence
See also
concatenateConcatenate arbitrary PulseSequences.
Notes
The total control matrix is given by
\[\begin{split}\tilde{\mathcal{B}}(\omega) &= \tilde{\mathcal{B}}^{(1)}(\omega) \sum_{g=0}^{G-1} \left(e^{i\omega T}\right)^g \\ &= \tilde{\mathcal{B}}^{(1)}(\omega)\bigl( \mathbb{I} - e^{i\omega T} \mathcal{Q}^{(1)}\bigr)^{-1}\bigl( \mathbb{I} - \bigl(e^{i\omega T} \mathcal{Q}^{(1)}\bigr)^G\bigr).\end{split}\]with \(T\) the period of the control Hamiltonian and \(G\) the number of periods. The last equality is valid only if \(\mathbb{I} - e^{i\omega T}\mathcal{Q}^{(1)}\) is invertible.
- extend(pulse_to_qubit_mapping: Sequence[Sequence[PulseSequence | Sequence[int] | int | Mapping[str, str] | None]], N: int | None = None, d_per_qubit: int = 2, additional_noise_Hamiltonian: Sequence[Sequence[ndarray | Qobj | Sequence[float]]] | None = None, cache_diagonalization: bool | None = None, cache_filter_function: bool | None = None, omega: Sequence[float] | None = None, show_progressbar: bool = False) PulseSequence[source]
Map one or more pulse sequences to different qubits.
- Parameters:
- pulse_to_qubit_mapping: sequence of mapping tuples
A sequence of tuples with the first entry a
PulseSequenceinstance and the second anintor tuple ofints indicating the qubits that thePulseSequenceshould be mapped to. A mapping of operator identifiers may optionally be given as a third element of each tuple. By default, the index of the qubit the operator is mapped to is appended to its identifier.Pulse sequences defined for multiple qubits may also be extended to non-neighboring qubits. Note that for multi-qubit pulses the order of the qubits is respected, i.e. mapping a pulse to (1, 0) is different from mapping it to (0, 1).
- N: int
The total number of qubits the new
PulseSequenceshould be defined for. By default, this is inferred frompulse_to_qubit_mapping.- d_per_qubit: int
The size of the Hilbert space a single qubit requires.
- additional_noise_Hamiltonian: list of lists
Additional noise operators and corresponding sensitivities for the new pulse sequence.
- cache_diagonalization: bool
Force diagonalizing the new pulse sequence. By default, diagonalization is cached if all pulses in
pulse_to_qubit_mappinghave been diagonalized since it is much cheaper to get the relevant quantities as tensor products from the mapped pulses instead of diagonalizing the new pulse.- cache_filter_function: bool
Force computing the filter functions for the new pulse sequence. Noise operators of individual pulses will be extended to the new Hilbert space. By default, this is done if all pulses in
pulse_to_qubit_mappinghave their filter functions cached.Note that extending the filter functions is only possible if they the mapped pulses are using a separable basis like the Pauli basis.
- omega: array_like
Frequencies for which to compute the filter functions if
cache_filter_function == True. Defaults toNone, in which case the cached frequencies of the individual pulses need to be the same.- show_progressbar: bool
Show a progress bar for the calculation of the control matrix.
- Returns:
- newpulse: PulseSequence
The new pulse sequence on the larger qubit register. The noise operators (and possibly filter functions) are stored in the following order: first those of the multi-qubit pulses in the order they appeared in
pulse_to_qubit_mapping, then those of the single-qubit pulses, and lastly any additional ones that may be given byadditional_noise_Hamiltonian.
See also
remapMap PulseSequence to a different qubit.
concatenateConcatenate PulseSequences (in time).
concatenate_periodicPeriodically concatenate a PulseSequence.
Examples
>>> import filter_functions as ff >>> I, X, Y, Z = ff.util.paulis >>> X_pulse = ff.PulseSequence([[X, [np.pi/2], 'X']], ... [[X, [1], 'X'], [Z, [1], 'Z']], ... [1], basis=ff.Basis.pauli(1)) >>> XX_pulse = ff.extend([(X_pulse, 0), (X_pulse, 1)]) >>> XX_pulse.d 4 >>> XIX_pulse_1 = ff.extend([(X_pulse, 0), (X_pulse, 2)]) >>> XIX_pulse_1.d 8 >>> XXI_pulse = ff.extend([(X_pulse, 0), (X_pulse, 1)], N=3) >>> XXI_pulse.d 8
Filter functions are automatically cached if they are for mapped pulses:
>>> omega = ff.util.get_sample_frequencies(X_pulse) >>> X_pulse.cache_filter_function(omega) >>> XX_pulse = ff.extend([(X_pulse, 0), (X_pulse, 1)]) >>> XX_pulse.is_cached('filter_function') True
This behavior can also be overriden manually:
>>> XX_pulse = ff.extend([(X_pulse, 0), (X_pulse, 1)], ... cache_filter_function=False) >>> XX_pulse.is_cached('filter_function') False
Mapping pulses to non-neighboring qubits is also possible:
>>> Y_pulse = ff.PulseSequence([[Y, [np.pi/2], 'Y']], ... [[Y, [1], 'Y'], [Z, [1], 'Z']], ... [1], basis=ff.Basis.pauli(1)) >>> XXY_pulse = ff.extend([(XX_pulse, (0, 1)), (Y_pulse, 2)]) >>> XYX_pulse = ff.extend([(XX_pulse, (0, 2)), (Y_pulse, 1)])
Additionally, pulses can have the order of the qubits they are defined for permuted (see
remap()):>>> Z_pulse = ff.PulseSequence([[Z, [np.pi/2], 'Z']], [[Z, [1], 'Z']], ... [1], basis=ff.Basis.pauli(1)) >>> XY_pulse = ff.extend([(X_pulse, 0), (Y_pulse, 1)]) >>> YZX_pulse = ff.extend([(XY_pulse, (2, 0)), (Z_pulse, 1)])
Control and noise operator identifiers can be mapped according to a specified mapping:
>>> YX_pulse = ff.extend([(X_pulse, 1, {'X': 'IX', 'Z': 'IZ'}), ... (Y_pulse, 0, {'Y': 'YI', 'Z': 'ZI'})]) >>> YX_pulse.c_oper_identifiers array(['IX', 'YI'], dtype='<U2') >>> YX_pulse.n_oper_identifiers array(['IX', 'IZ', 'YI', 'ZI'], dtype='<U2')
We can also add an additional noise Hamiltonian:
>>> H_n = [[ff.util.tensor(Z, Z, Z), [1], 'ZZZ']] >>> XYX_pulse = ff.extend([(XX_pulse, (0, 2)), (Y_pulse, 1)], ... additional_noise_Hamiltonian=H_n) >>> 'ZZZ' in XYX_pulse.n_oper_identifiers True
- remap(pulse: PulseSequence, order: Sequence[int], d_per_qubit: int = 2, oper_identifier_mapping: Mapping[str, str] = None) PulseSequence[source]
Remap a PulseSequence by changing the order of qubits in the register. Cached attributes are automatically attempted to be retained.
Caution
This function simply permutes the order of the tensor product elements of control and noise operators. Thus, the resultant pulse will have its filter functions defined for different noise operators than the original one.
- Parameters:
- pulse: PulseSequence
The pulse whose qubit order should be permuted.
- order: sequence of ints
A list of permutation indices. E.g., if pulse is defined for two qubits,
order == [1, 0]will reverse the order of qubits.- d_per_qubit: int (default: 2)
The size of the Hilbert space a single qubit inhabitates.
- oper_identifier_mapping: dict_like
A mapping that maps operator identifiers from the old pulse to the remapped pulse. The default is the identity mapping.
- Returns:
- remapped_pulse: PulseSequence
A new
PulseSequenceinstance with the order of the qubits permuted according to order.
See also
extendMap PulseSequences to composite Hilbert spaces.
util.tensor_transposeTranspose the order of a tensor product.
Examples
>>> X, Y = util.paulis[1:3] >>> XY, YX = util.tensor(X, Y), util.tensor(Y, X) >>> pulse = PulseSequence([[XY, [np.pi/2], 'XY']], [[YX, [1], 'YX']], [1], ... Basis.pauli(2)) >>> mapping = {'XY': 'YX', 'YX': 'XY'} >>> remapped_pulse = remap(pulse, (1, 0), oper_identifier_mapping=mapping) >>> target_pulse = PulseSequence([[YX, [np.pi/2], 'YX']], ... [[XY, [1], 'XY']], [1], Basis.pauli(2)) >>> remapped_pulse == target_pulse True
Caching of attributes is automatically handled
>>> remapped_pulse.is_cached('filter_function') False >>> pulse.cache_filter_function(util.get_sample_frequencies(pulse)) >>> remapped_pulse = remap(pulse, (1, 0)) >>> remapped_pulse.is_cached('filter_function') True