filter_functions.numeric module
This module defines the functions to calculate everything related to filter functions.
Functions
calculate_control_matrix_from_atomic()Calculate the control matrix from those of atomic pulse sequences
calculate_control_matrix_from_scratch()Calculate the control matrix from scratch
calculate_control_matrix_periodic()Calculate the control matrix for a periodic Hamiltonian
calculate_noise_operators_from_atomic()Calculate the interaction picture noise operators from atomic segments. Same calculation as
calculate_control_matrix_from_atomic()except in Hilbert space.calculate_noise_operators_from_scratch()Calculate the interaction picture noise operators from scratch. Same calculation as
calculate_control_matrix_from_scratch()except in Hilbert space.calculate_cumulant_function()Calculate the cumulant function for a given
PulseSequenceobject.calculate_decay_amplitudes()Calculate the decay amplitudes, corresponding to first order terms of the Magnus expansion
calculate_frequency_shifts()Calculate the frequency shifts, corresponding to second order terms of the Magnus expansion
calculate_filter_function()Calculate the filter function from the control matrix
calculate_second_order_filter_function()Calculate the second order filter function used to compute the frequency shifts.
calculate_pulse_correlation_filter_function()Calculate the pulse correlation filter function from the control matrix
diagonalize()Diagonalize a Hamiltonian
error_transfer_matrix()Calculate the error transfer matrix of a pulse up to a unitary rotation
infidelity()Function to compute the infidelity of a pulse defined by a
PulseSequenceinstance for a given noise spectral density and frequencies
- calculate_control_matrix_from_atomic(phases: ndarray, control_matrix_atomic: ndarray, propagators_liouville: ndarray, show_progressbar: bool = False, which: str = 'total') ndarray | Tuple[ndarray, ndarray][source]
Calculate the control matrix from the control matrices of atomic segments.
- Parameters:
- phases: array_like, shape (G-1, n_omega)
The phase factors for \(g\in\{1, 2, \dots, G-1\}\). For \(g=0\), they are unity.
- control_matrix_atomic: array_like, shape (G, n_nops, d**2, n_omega)
The pulse control matrices for \(g\in\{1, 2, \dots, G\}\).
- propagators_liouville: array_like, shape (G-1, d**2, d**2)
The transfer matrices of the cumulative propagators for \(g\in\{1, 2, \dots, G-1\}\). For \(g=0\) it is the identity.
- show_progressbar: bool, optional
Show a progress bar for the calculation.
- which: str, (‘total’, ‘correlations’)
Compute the total control matrix (the sum of all time steps) or the correlation control matrix (first axis holds each pulses’ contribution)
- Returns:
- control_matrix: ndarray, shape ([G,] n_nops, d**2, n_omega)
The control matrix \(\tilde{\mathcal{B}}(\omega)\).
See also
calculate_control_matrix_from_scratchControl matrix from scratch.
liouville_representationLiouville representation for a given basis.
Notes
The control matrix is calculated by evaluating the sum
\[\tilde{\mathcal{B}}(\omega) = \sum_{g=1}^G e^{i\omega t_{g-1}} \tilde{\mathcal{B}}^{(g)}(\omega)\mathcal{Q}^{(g-1)}.\]
- calculate_control_matrix_from_scratch(eigvals: ndarray, eigvecs: ndarray, propagators: ndarray, omega: Sequence[float], basis: Basis, n_opers: Sequence[ndarray | Qobj], n_coeffs: Sequence[Sequence[float]], dt: Sequence[float], t: Sequence[float] | None = None, show_progressbar: bool = False, cache_intermediates: bool = False, out: ndarray | None = None) ndarray | Tuple[ndarray, Dict[str, ndarray]][source]
Calculate the control matrix from scratch, i.e. without knowledge of the control matrices of more atomic pulse sequences.
- Parameters:
- eigvals: array_like, shape (n_dt, d)
Eigenvalue vectors for each time pulse segment g with the first axis counting the pulse segment, i.e.
eigvals == array([D_0, D_1, ...]).- eigvecs: array_like, shape (n_dt, d, d)
Eigenvector matrices for each time pulse segment g with the first axis counting the pulse segment, i.e.
eigvecs == array([V_0, V_1, ...]).- propagators: array_like, shape (n_dt+1, d, d)
The propagators \(Q_g = P_g P_{g-1}\cdots P_0\) as a (d, d) array with d the dimension of the Hilbert space.
- omega: array_like, shape (n_omega,)
Frequencies at which the pulse control matrix is to be evaluated.
- basis: Basis, shape (d**2, d, d)
The basis elements in which the pulse control matrix will be expanded.
- n_opers: array_like, shape (n_nops, d, d)
Noise operators \(B_\alpha\).
- n_coeffs: array_like, shape (n_nops, n_dt)
The sensitivities of the system to the noise operators given by n_opers at the given time step.
- dt: array_like, shape (n_dt,)
Sequence duration, i.e. for the \(g\)-th pulse \(t_g - t_{g-1}\).
- t: array_like, shape (n_dt+1,), optional
The absolute times of the different segments. Can also be computed from dt.
- show_progressbar: bool, optional
Show a progress bar for the calculation.
- cache_intermediates: bool, optional
Keep and return intermediate terms of the calculation that can be reused in other computations (second order or gradients). Otherwise the sum is performed in-place. The default is False.
- out: ndarray, optional
A location into which the result is stored. See
numpy.ufunc().
- Returns:
- control_matrix: ndarray, shape (n_nops, d**2, n_omega)
The control matrix \(\tilde{\mathcal{B}}(\omega)\)
- intermediates: dict[str, ndarray]
Intermediate results of the calculation. Only if cache_intermediates is True.
See also
calculate_control_matrix_from_atomicControl matrix from concatenation.
calculate_control_matrix_periodicControl matrix for periodic system.
Notes
The control matrix is calculated according to
\[\tilde{\mathcal{B}}_{\alpha k}(\omega) = \sum_{g=1}^G e^{i\omega t_{g-1}} s_\alpha^{(g)}\mathrm{tr}\left( [\bar{B}_\alpha^{(g)}\circ I(\omega)] \bar{C}_k^{(g)} \right)\]where
\[\begin{split}I^{(g)}_{nm}(\omega) &= \int_0^{t_l - t_{g-1}}\mathrm{d}t\, e^{i(\omega+\omega_n-\omega_m)t} \\ &= \frac{e^{i(\omega+\omega_n-\omega_m) (t_l - t_{g-1})} - 1} {i(\omega+\omega_n-\omega_m)}, \\ \bar{B}_\alpha^{(g)} &= V^{(g)\dagger} B_\alpha V^{(g)}, \\ \bar{C}_k^{(g)} &= V^{(g)\dagger} Q_{g-1} C_k Q_{g-1}^\dagger V^{(g)},\end{split}\]and \(V^{(g)}\) is the matrix of eigenvectors that diagonalizes \(\tilde{\mathcal{H}}_n^{(g)}\), \(B_\alpha\) the \(\alpha\)-th noise operator \(s_\alpha^{(g)}\) the noise sensitivity during interval \(g\), and \(C_k\) the \(k\)-th basis element.
- calculate_control_matrix_periodic(phases: ndarray, control_matrix: ndarray, total_propagator_liouville: ndarray, repeats: int, check_invertible: bool = True) ndarray[source]
Calculate the control matrix of a periodic pulse given the phase factors, control matrix and transfer matrix of the total propagator, total_propagator_liouville, of the atomic pulse.
- Parameters:
- phases: ndarray, shape (n_omega,)
The phase factors \(e^{i\omega T}\) of the atomic pulse.
- control_matrix: ndarray, shape (n_nops, d**2, n_omega)
The control matrix \(\tilde{\mathcal{B}}^{(1)}(\omega)\) of the atomic pulse.
- total_propagator_liouville: ndarray, shape (d**2, d**2)
The transfer matrix \(\mathcal{Q}^{(1)}\) of the atomic pulse.
- repeats: int
The number of repetitions.
- check_invertible: bool, optional
Check for all frequencies if the inverse \(\mathbb{I} - e^{i\omega T} \mathcal{Q}^{(1)}\) exists by calculating the determinant. If it does not exist, the sum is evaluated explicitly for those points. The default is True.
- Returns:
- control_matrix: ndarray, shape (n_nops, d**2, n_omega)
The control matrix \(\tilde{\mathcal{B}}(\omega)\) of the repeated pulse.
Notes
The control matrix is computed as
\[\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 \(G\) the number of repetitions.
- calculate_cumulant_function(pulse: PulseSequence, spectrum: ndarray | None = None, omega: Sequence[float] | None = None, n_oper_identifiers: Sequence[str] | None = None, which: str = 'total', second_order: bool = False, decay_amplitudes: ndarray | None = None, frequency_shifts: ndarray | None = None, show_progressbar: bool = False, memory_parsimonious: bool = False, cache_intermediates: bool | None = None) ndarray[source]
Calculate the cumulant function \(\mathcal{K}(\tau)\).
The error transfer matrix is obtained from the cumulant function by exponentiation, \(\langle\tilde{\mathcal{U}}\rangle = \exp\mathcal{K}(\tau)\).
- Parameters:
- pulse: PulseSequence
The
PulseSequenceinstance for which to compute the cumulant function.- spectrum: array_like, shape ([[n_nops,] n_nops,] n_omega), optional
The noise power spectral density in units of inverse frequencies as an array of shape (n_omega,), (n_nops, n_omega), or (n_nops, n_nops, n_omega). In the first case, the same spectrum is taken for all noise operators, in the second, it is assumed that there are no correlations between different noise sources and thus there is one spectrum for each noise operator. In the third and most general case, there may be a spectrum for each pair of noise operators corresponding to the correlations between them. n_nops is the number of noise operators considered and should be equal to
len(n_oper_identifiers).- omega: array_like, shape (n_omega,), optional
The frequencies at which to evaluate the filter functions.
- n_oper_identifiers: array_like, optional
The identifiers of the noise operators for which to evaluate the cumulant function. The default is all.
- which: str, optional
Which decay amplitudes should be calculated, may be either ‘total’ (default) or ‘correlations’. See
infidelity()and Notes. Note that the latter is not available for the second order terms.- second_order: bool, optional
Also take into account the frequency shifts \(\Delta\) that correspond to second order Magnus expansion and constitute unitary terms. Default
False.- decay_amplitudes, array_like, shape ([[G, G,] n_nops,] n_nops, d**2, d**2), optional
A precomputed cumulant function. If given, spectrum, omega are not required.
- frequency_shifts, array_like, shape ([[G, G,] n_nops,] n_nops, d**2, d**2), optional
A precomputed frequency shift. If given, spectrum, omega are not required for second order terms.
- show_progressbar: bool, optional
Show a progress bar for the calculation of the control matrix.
- memory_parsimonious: bool, optional
Trade memory footprint for performance. See
calculate_decay_amplitudes(). The default isFalse.- cache_intermediates: bool, optional
Keep and return intermediate terms of the calculation of the control matrix that can be reused in other computations (second order or gradients). Otherwise the sum is performed in-place. Default is True if second_order=True, else False.
- Returns:
- cumulant_function: ndarray, shape ([[G, G,] n_nops,] n_nops, d**2, d**2)
The cumulant function. The individual noise operator contributions chosen by
n_oper_identifiersare on the third to last axis / axes, depending on whether the noise is cross-correlated or not. Ifwhich == 'correlations', the first two axes correspond to the contributions of the pulses in the sequence.
See also
calculate_decay_amplitudesCalculate the \(\Gamma_{\alpha\beta,kl}\)
error_transfer_matrixCalculate the error transfer matrix \(\exp\mathcal{K}\).
infidelityCalculate only infidelity of a pulse.
pulse_sequence.concatenateConcatenate
PulseSequenceobjects.calculate_pulse_correlation_filter_function
Notes
The cumulant function is given by
\[\begin{split}K_{\alpha\beta,ij}(\tau) = -\frac{1}{2} \sum_{kl}\biggl( &\Delta_{\alpha\beta,kl}\left( T_{klji} - T_{lkji} - T_{klij} + T_{lkij} \right) \\ + &\Gamma_{\alpha\beta,kl}\left( T_{klji} - T_{kjli} - T_{kilj} + T_{kijl} \right) \biggr)\end{split}\]Here, \(T_{ijkl} = \mathrm{tr}(C_i C_j C_k C_l)\) is a trivial function of the basis elements \(C_i\), and \(\Gamma_{\alpha\beta,kl}\) and \(\Delta_{\alpha\beta,kl}\) are the decay amplitudes and frequency shifts which correspond to first and second order in the Magnus expansion, respectively. Since the latter induce coherent errors, we can approximately neglect them if we assume that the pulse has been experimentally calibrated.
For a single qubit and represented in the Pauli basis, the above reduces to
\[\begin{split}K_{\alpha\beta,ij}(\tau) = \begin{cases} - \sum_{k\neq i}\Gamma_{\alpha\beta,kk} &\quad\mathrm{if}\: i = j, \\ - \Delta_{\alpha\beta,ij} + \Delta_{\alpha\beta,ji} + \Gamma_{\alpha\beta,ij} &\quad\mathrm{if}\: i\neq j, \end{cases}\end{split}\]for \(i\in\{1,2,3\}\) and \(K_{0j} = K_{i0} = 0\).
Lastly, the pulse correlation cumulant function resolves correlations in the cumulant function of a sequence of pulses \(g = 1,\dotsc,G\) such that the following holds:
\[K_{\alpha\beta,ij}(\tau) = \sum_{g,g'=1}^G K_{\alpha\beta,ij}^{(gg')}(\tau).\]
- calculate_decay_amplitudes(pulse: PulseSequence, spectrum: ndarray, omega: Sequence[float], n_oper_identifiers: Sequence[str] | None = None, which: str = 'total', show_progressbar: bool = False, cache_intermediates: bool = False, memory_parsimonious: bool = False) ndarray[source]
Get the decay amplitudes \(\Gamma_{\alpha\beta, kl}\) for noise sources \(\alpha,\beta\) and basis elements \(k,l\).
- Parameters:
- pulse: PulseSequence
The
PulseSequenceinstance for which to compute the decay amplitudes.- spectrum: array_like, shape ([[n_nops,] n_nops,] n_omega)
The noise power spectral density. If 1-d, the same spectrum is used for all noise operators. If 2-d, one (self-) spectrum for each noise operator is expected. If 3-d, should be a matrix of cross-spectra such that
spectrum[i, j] == spectrum[j, i].conj().- omega: array_like,
The frequencies at which to calculate the filter functions.
- n_oper_identifiers: array_like, optional
The identifiers of the noise operators for which to calculate the decay amplitudes. The default is all.
- which: str, optional
Which decay amplitudes should be calculated, may be either ‘total’ (default) or ‘correlations’. See
infidelity()and Notes.- show_progressbar: bool, optional
Show a progress bar for the calculation.
- cache_intermediates: bool, optional
Keep and return intermediate terms of the calculation that are useful in other places (if control matrix not already cached).
- memory_parsimonious: bool, optional
For large dimensions, the integrand
\[\tilde{\mathcal{B}}^\ast_{\alpha k}(\omega) S_{\alpha\beta}(\omega)\tilde{\mathcal{B}}_{\beta l}(\omega)\]can consume quite a large amount of memory if set up for all \(\alpha,\beta,k,l\) at once. If
True, it is only set up and integrated for a single \(k\) at a time and looped over. This is slower but requires much less memory. The default isFalse.
- Returns:
- decay_amplitudes: ndarray, shape ([[G, G,] n_nops,] n_nops, d**2, d**2)
The decay amplitudes.
- Raises:
- ValueError
If spectrum has incompatible shape.
See also
infidelityCompute the infidelity directly.
pulse_sequence.concatenateConcatenate
PulseSequenceobjects.calculate_frequency_shiftsSecond order (unitary) terms.
calculate_pulse_correlation_filter_function
Notes
The total decay amplitudes are given by
\[\Gamma_{\alpha\beta, kl} = \int\frac{\mathrm{d}\omega}{2\pi} \tilde{\mathcal{B}}^\ast_{\alpha k}(\omega) S_{\alpha\beta}(\omega)\tilde{\mathcal{B}}_{\beta l}(\omega).\]If pulse correlations are taken into account, they are given by
\[\Gamma_{\alpha\beta, kl}^{(gg')} = \int \frac{\mathrm{d}\omega}{2\pi} S_{\alpha\beta}(\omega) F_{\alpha\beta, kl}^{(gg')}(\omega).\]
- calculate_filter_function(control_matrix: ndarray, which: str = 'fidelity') ndarray[source]
Compute the filter function from the control matrix.
- Parameters:
- control_matrix: array_like, shape (n_nops, d**2, n_omega)
The control matrix.
- whichstr, optional
Which filter function to return. Either ‘fidelity’ (default) or ‘generalized’ (see Notes).
- Returns:
- filter_function: ndarray, shape (n_nops, n_nops, [d**2, d**2,] n_omega)
The filter functions for each noise operator correlation. The diagonal corresponds to the filter functions for uncorrelated noise sources.
See also
calculate_control_matrix_from_scratchControl matrix from scratch.
calculate_control_matrix_from_atomicControl matrix from concatenation.
calculate_pulse_correlation_filter_functionPulse correlations.
Notes
The 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).\]
- calculate_frequency_shifts(pulse: PulseSequence, spectrum: ndarray, omega: Sequence[float], n_oper_identifiers: Sequence[str] | None = None, show_progressbar: bool = False) ndarray[source]
Get the frequency shifts \(\Delta_{\alpha\beta, kl}\) for noise sources \(\alpha,\beta\) and basis elements \(k,l\).
- Parameters:
- pulse: PulseSequence
The
PulseSequenceinstance for which to compute the frequency shifts.- spectrum: array_like, shape ([[n_nops,] n_nops,] n_omega)
The two-sided noise power spectral density. If 1-d, the same spectrum is used for all noise operators. If 2-d, one (self-) spectrum for each noise operator is expected. If 3-d, should be a matrix of cross-spectra such that
spectrum[i, j] == spectrum[j, i].conj().- omega: array_like,
The frequencies. Note that the frequencies are assumed to be symmetric about zero.
- n_oper_identifiers: array_like, optional
The identifiers of the noise operators for which to calculate the frequency shifts. The default is all.
- show_progressbar: bool, optional
Show a progress bar for the calculation.
- Returns:
- Delta: ndarray, shape ([n_nops,] n_nops, d**2, d**2)
The frequency shifts.
- Raises:
- ValueError
If spectrum has incompatible shape.
See also
calculate_second_order_filter_functionCorresponding filter function.
calculate_decay_amplitudesFirst order (dissipative) terms.
infidelityCompute the infidelity directly.
pulse_sequence.concatenateConcatenate
PulseSequenceobjects.calculate_pulse_correlation_filter_function
Notes
The total frequency shifts are given by
\[\Delta_{\alpha\beta, kl} = \int_{-\infty}^\infty \frac{\mathrm{d}{\omega}}{2\pi} S_{\alpha\beta}(\omega) F_{\alpha\beta,kl}^{(2)}(\omega)\]with \(F_{\alpha\beta,kl}^{(2)}(\omega)\) the second order filter function.
- calculate_noise_operators_from_atomic(phases: ndarray, noise_operators_atomic: ndarray, propagators: ndarray, show_progressbar: bool = False) ndarray[source]
Calculate the interaction picutre noise operators from atomic segments.
- Parameters:
- phases: array_like, shape (G-1, n_omega)
The phase factors for \(g\in\{1, 2, \dots, G-1\}\). For \(g=0\) they are unity.
- noise_operators_atomic: array_like, shape (G, n_nops, d, d, n_omega)
The noise operators in the interaction picture of the g-th pulse, i.e. for \(g\in\{1, 2, \dots, G\}\).
- propagators: array_like, shape (G-1, d, d)
The cumulative propagators of the pulses \(g\in\{1, 2, \dots, G-1\}\). For \(g=0\) it is the identity.
- show_progressbar: bool, optional
Show a progress bar for the calculation.
- Returns:
- noise_operators: ndarray, shape (n_omega, n_nops, d, d)
The interaction picture noise operators \(\tilde{B}_\alpha(\omega)\).
See also
calculate_noise_operators_from_scratchCompute the operators from scratch.
calculate_control_matrix_from_atomicSame calculation in Liouville space.
Notes
The noise operators are calculated by evaluating the sum
\[\tilde{B}_\alpha(\omega) = \sum_{g=1}^G e^{i\omega t_{g-1}} Q_{g-1}^\dagger\tilde{B}_\alpha^{(g)}(\omega) Q_{g-1}.\]The control matrix then corresponds to the coefficients of expansion in an operator basis \(\{C_k\}_k\):
\[\tilde{\mathcal{B}}_{k\alpha}(\omega) = \mathrm{tr}(\tilde{B}_\alpha(\omega) C_k).\]Due to differences in implementation (for performance reasons), the axes of the result are transposed compared to the control matrix:
>>> ctrlmat = calculate_control_matrix_from_atomic(...) >>> ctrlmat.shape (n_nops, d**2, n_omega) >>> noiseops = calculate_noise_operators_from_atomic(...) >>> noiseops.shape (n_omega, n_nops, d, d) >>> ctrlmat_from_noiseops = basis.expand(noiseops) >>> np.allclose(ctrlmat, ctrlmat_from_noiseops.transpose(1, 2, 0)) True
- calculate_noise_operators_from_scratch(eigvals: ndarray, eigvecs: ndarray, propagators: ndarray, omega: Sequence[float], n_opers: Sequence[ndarray | Qobj], n_coeffs: Sequence[Sequence[float]], dt: Sequence[float], t: Sequence[float] | None = None, show_progressbar: bool = False, cache_intermediates: bool = False) ndarray | Tuple[ndarray, Dict[str, ndarray]][source]
Calculate the noise operators in interaction picture from scratch.
- Parameters:
- eigvals: array_like, shape (n_dt, d)
Eigenvalue vectors for each time pulse segment g with the first axis counting the pulse segment, i.e.
eigvals == array([D_0, D_1, ...]).- eigvecs: array_like, shape (n_dt, d, d)
Eigenvector matrices for each time pulse segment g with the first axis counting the pulse segment, i.e.
eigvecs == array([V_0, V_1, ...]).- propagators: array_like, shape (n_dt+1, d, d)
The propagators \(Q_g = P_g P_{g-1}\cdots P_0\) as a (d, d) array with d the dimension of the Hilbert space.
- omega: array_like, shape (n_omega,)
Frequencies at which the pulse control matrix is to be evaluated.
- n_opers: array_like, shape (n_nops, d, d)
Noise operators \(B_\alpha\).
- n_coeffs: array_like, shape (n_nops, n_dt)
The sensitivities of the system to the noise operators given by n_opers at the given time step.
- dt: array_like, shape (n_dt,)
Sequence duration, i.e. for the \(g\)-th pulse \(t_g - t_{g-1}\).
- t: array_like, shape (n_dt+1,), optional
The absolute times of the different segments. Can also be computed from dt.
- show_progressbar: bool, optional
Show a progress bar for the calculation.
- cache_intermediates: bool, optional
Keep and return intermediate terms of the calculation that can be reused in other computations (second order or gradients). Otherwise the sum is performed in-place. The default is False.
- Returns:
- noise_operators: ndarray, shape (n_omega, n_nops, d, d)
The interaction picture noise operators \(\tilde{B}_\alpha(\omega)\).
- intermediates: dict[str, ndarray]
Intermediate results of the calculation. Only if cache_intermediates is True.
See also
calculate_noise_operators_from_atomicCompute the operators from atomic segments.
calculate_control_matrix_from_scratchSame calculation in Liouville space.
Notes
The interaction picture noise operators are calculated according to
\[\tilde{B}_\alpha(\omega) = \sum_{g=1}^G e^{i\omega t_{g-1}} s_\alpha^{(g)} P^{(g)\dagger}\left[ \bar{B}^{(g)}_\alpha \circ I^{(g)}(\omega) \right] P^{(g)}\]where
\[\begin{split}I^{(g)}_{nm}(\omega) &= \int_0^{t_g - t_{g-1}}\mathrm{d}t\, e^{i(\omega+\omega_n-\omega_m)t} \\ &= \frac{e^{i(\omega+\omega_n-\omega_m) (t_g - t_{g-1})} - 1} {i(\omega+\omega_n-\omega_m)}, \\ \bar{B}_\alpha^{(g)} &= V^{(g)\dagger} B_\alpha V^{(g)}, \\ P^{(g)} &= V^{(g)\dagger} Q_{g-1},\end{split}\]and \(V^{(g)}\) is the matrix of eigenvectors that diagonalizes \(\tilde{\mathcal{H}}_n^{(g)}\), \(B_\alpha\) the \(\alpha\)-th noise operator, and \(s_\alpha^{(g)}\) the noise sensitivity during interval \(g\).
The control matrix then corresponds to the coefficients of expansion in an operator basis \(\{C_k\}_k\):
\[\tilde{\mathcal{B}}_{k\alpha}(\omega) = \mathrm{tr}(\tilde{B}_\alpha(\omega) C_k).\]Due to differences in implementation (for performance reasons), the axes of the result are transposed compared to the control matrix:
>>> ctrlmat = calculate_control_matrix_from_scratch(...) >>> ctrlmat.shape (n_nops, d**2, n_omega) >>> noiseops = calculate_noise_operators_from_scratch(...) >>> noiseops.shape (n_omega, n_nops, d, d) >>> ctrlmat_from_noiseops = basis.expand(noiseops) >>> np.allclose(ctrlmat, ctrlmat_from_noiseops.transpose(1, 2, 0)) True
- calculate_pulse_correlation_filter_function(control_matrix: ndarray, which: str = 'fidelity') ndarray[source]
Compute pulse correlation filter function from control matrix.
- Parameters:
- control_matrix: array_like, shape (n_pulses, n_nops, d**2, n_omega)
The control matrix.
- whichstr, optional
Which filter function to return. Either ‘fidelity’ (default) or ‘generalized’ (see Notes).
- Returns:
- filter_function_pc: ndarray, shape (G, G, n_nops, n_nops, [d**2, d**2,] n_omega)
The pulse correlation filter functions for each pulse and noise operator correlations. The first two axes hold the pulse correlations, the second two the noise correlations.
- whichstr, optional
Which filter function to return. Either ‘fidelity’ (default) or ‘generalized’ (see Notes).
See also
calculate_control_matrix_from_scratchControl matrix from scratch.
calculate_control_matrix_from_atomicControl matrix from concatenation.
calculate_filter_functionRegular filter function.
Notes
The generalized pulse correlation filter function is given by
\[F_{\alpha\beta,kl}^{(gg')}(\omega) = \bigl[ \mathcal{Q}^{(g'-1)\dagger} \tilde{\mathcal{B}}^{(g')\dagger}(\omega) \bigr]_{k\alpha} \bigl[ \tilde{\mathcal{B}}^{(g)}(\omega)\mathcal{Q}^{(g-1)} \bigr]_{\beta l} e^{i\omega(t_{g-1} - t_{g'-1})},\]with \(\tilde{\mathcal{B}}^{(g)}\) the control matrix of the \(g\)-th pulse. The fidelity pulse correlation function is obtained by tracing out the basis indices,
\[F_{\alpha\beta}^{(gg')}(\omega) = \sum_{k} F_{\alpha\beta,kk}^{(gg')}(\omega)\]
- calculate_second_order_filter_function_from_atomic(basis: Basis, filter_function_atomic: ndarray, control_matrix_atomic: ndarray, control_matrix_atomic_step: ndarray, control_matrix_atomic_cumulative: ndarray, propagators: ndarray, propagators_liouville: ndarray, intermediates: Sequence[Mapping[str, ndarray]], show_progressbar: bool = False)[source]
Calculate the second-order filter function from those of atomic segments.
- Parameters:
- basis: Basis, shape (d**2, d, d)
The operator basis for the filter function.
- filter_function_atomic: ndarray, shape (n_nops, n_nops, d**2, d**2, n_omega)
The filter function for the first segment, \(g = 1\).
- control_matrix_atomic: ndarray, shape (G, n_nops, d**2, n_omega)
The pulse control matrices for \(g\in\{1, 2, \dots, G\}\).
- control_matrix_atomic_step: ndarray, shape (G-1, n_nops, d**2, n_omega)
The pulse correlation control matrix, i.e., the summands of the first-order concatenated control matrix of the pulse sequence, \(\mathcal{G}^{(g)}(\omega)\).
- control_matrix_atomic_cumulative: ndarray, shape (G, n_nops, d**2, n_omega)
- propagators: ndarray, shape (G-1, d, d)
The cumulative propagators of the G pulses.
- propagators_liouville: ndarray, shape (G-1, d**2, d**2)
The transfer matrices of the cumulative propagators for \(g\in\{1, 2, \dots, G-1\}\). For \(g=0\) it is the identity.
- intermediates: Sequence[Dict[str, ndarray]}
Intermediate terms of the calculation of the control matrix and second-order filter function that can be reused here. Each entry of the sequence of length G should be a dictionary with the following required keys:
‘eigvecs_propagated’
‘n_opers_transformed’
‘second_order_integral’
‘second_order_complete_steps’
The first two are populated by
calculate_control_matrix_from_scratch(), the last two bycalculate_second_order_filter_function_from_scratch().- show_progressbar: bool, optional
Show a progress bar for the calculation.
- Returns:
- second_order_FF: ndarray, shape (n_nops, n_nops, d**2, d**2, n_omega)
The second-order filter function \(\mathcal{F}^{(2)}(\omega)\).
See also
calculate_control_matrix_from_scratchControl matrix from scratch.
calculate_control_matrix_from_atomicSimilar function for the first order FF.
calculate_second_order_filter_function_from_scratchSecond-order FF from scratch.
- calculate_second_order_filter_function_from_scratch(eigvals: ndarray, eigvecs: ndarray, propagators: ndarray, omega: Sequence[float], basis: Basis, n_opers: Sequence[ndarray | Qobj], n_coeffs: Sequence[Sequence[float]], dt: Sequence[float], intermediates: Dict[str, ndarray] | None = None, show_progressbar: bool = False, cache_intermediates: bool = False, cache_cumulative: bool = False) ndarray | tuple[ndarray, Dict[str, ndarray]][source]
Calculate the second order filter function for frequency shifts.
- Parameters:
- eigvals: array_like, shape (n_dt, d)
Eigenvalue vectors for each time pulse segment l with the first axis counting the pulse segment, i.e.
eigvals == array([D_0, D_1, ...]).- eigvecs: array_like, shape (n_dt, d, d)
Eigenvector matrices for each time pulse segment l with the first axis counting the pulse segment, i.e.
eigvecs == array([V_0, V_1, ...]).- propagators: array_like, shape (n_dt+1, d, d)
The propagators \(Q_l = P_l P_{l-1}\cdots P_0\) as a (d, d) array with d the dimension of the Hilbert space.
- omega: array_like, shape (n_omega,)
Frequencies at which the pulse control matrix is to be evaluated.
- basis: Basis, shape (d**2, d, d)
The basis elements in which the pulse control matrix will be expanded.
- n_opers: array_like, shape (n_nops, d, d)
Noise operators \(B_\alpha\).
- n_coeffs: array_like, shape (n_nops, n_dt)
The sensitivities of the system to the noise operators given by n_opers at the given time step.
- dt: array_like, shape (n_dt,)
Sequence duration, i.e. for the \(l\)-th pulse \(t_l - t_{l-1}\).
- intermediates: Dict[str, ndarray], optional
Intermediate terms of the calculation of the control matrix that can be reused here. If None (default), they are computed from scratch.
- show_progressbar: bool, optional
Show a progress bar for the calculation.
- cache_intermediates: bool, optional
Return a cache with intermediate results that can be recycled during concatenation.
- cache_cumulativebool, optional
Additionally cache the accumulated filter function for each step \(g\), that is, an array that holds the filter function of the sequence up to position g in element g of its first axis. Only if cache_intermediates is True.
Note
This might be an extremely large array.
- Returns:
- second_order_filter_functionndarray, shape (n_nops, n_nops, d**2, d**2, n_omega)
The second order filter function.
- intermediatesdict[str, ndarray]
Intermediate results of the calculation. Only if cache_intermediates is True.
See also
calculate_frequency_shiftsIntegrate over filter function times spectrum.
calculate_decay_amplitudesFirst order (dissipative) terms.
infidelityCompute the infidelity directly.
pulse_sequence.concatenateConcatenate
PulseSequenceobjects.calculate_pulse_correlation_filter_functioncalculate_second_order_filter_function_from_atomicSecond-order FF from atomic.
Notes
The second order filter function is given by
\[F_{\alpha\beta, kl}^{(2)} = \sum_{g=1}^G\left[ \mathcal{G}_{\alpha k}^{(g)\ast}(\omega) \sum_{g'=1}^{g-1}\mathcal{G}_{\beta l}^{(g')}(\omega) + s_\alpha^{(g)}\bar{B}_{\alpha,ij}^{(g)}\bar{C}_{k,ji}^{(g)} I_{ijmn}^{(g)}(\omega)\bar{C}_{l,nm}^{(g) \bar{B}_{\beta,mn}^{(g)}s_\beta^{(g)}} \right]\]with
\[\begin{split}\mathcal{G}^{(g)}(\omega) &= e^{i\omega t_{g-1}}\mathcal{B}^{(g)}(\omega) \mathcal{Q}^{(g-1)}, \\ I_{ijmn}^{(g)}(\omega) &= \int_{t_{g-1}}^{t_g}\mathrm{d}{t} e^{i\Omega_{ij}^{(g)}(t - t_{g-1}) - i\omega t} \int_{t_{g-1}}^{t}\mathrm{d}{t'} e^{i\Omega_{mn}^{(g)}(t' - t_{g-1}) + i\omega t'}.\end{split}\]
- diagonalize(hamiltonian: ndarray, dt: Sequence[float]) Tuple[ndarray, ndarray, ndarray][source]
Diagonalize a Hamiltonian.
Diagonalize the Hamiltonian which is piecewise constant during the times given by dt and return eigenvalues, eigenvectors, and the cumulative propagators \(Q_l\). Note that we calculate in units where \(\hbar\equiv 1\) so that
\[U(t, t_0) = \mathcal{T}\exp\left( -i\int_{t_0}^t\mathrm{d}t'\mathcal{H}(t') \right).\]- Parameters:
- hamiltonian: array_like, shape (n_dt, d, d)
Hamiltonian of shape (n_dt, d, d) with d the dimensionality of the system.
- dt: array_like, shape (n_dt,)
The time differences.
- Returns:
- eigvals: ndarray
Array of eigenvalues of shape (n_dt, d)
- eigvecs: ndarray
Array of eigenvectors of shape (n_dt, d, d)
- propagators: ndarray
Array of cumulative propagators of shape (n_dt+1, d, d)
- error_transfer_matrix(pulse: PulseSequence | None = None, spectrum: ndarray | None = None, omega: Sequence[float] | None = None, n_oper_identifiers: Sequence[str] | None = None, second_order: bool = False, cumulant_function: ndarray | None = None, show_progressbar: bool = False, memory_parsimonious: bool = False, cache_intermediates: bool = False) ndarray[source]
Compute the error transfer matrix up to unitary rotations.
- Parameters:
- pulse: PulseSequence
The
PulseSequenceinstance for which to compute the error transfer matrix.- spectrum: array_like, shape ([[n_nops,] n_nops,] n_omega)
The two-sided noise power spectral density in units of inverse frequencies as an array of shape (n_omega,), (n_nops, n_omega), or (n_nops, n_nops, n_omega). In the first case, the same spectrum is taken for all noise operators, in the second, it is assumed that there are no correlations between different noise sources and thus there is one spectrum for each noise operator. In the third and most general case, there may be a spectrum for each pair of noise operators corresponding to the correlations between them. n_nops is the number of noise operators considered and should be equal to
len(n_oper_identifiers).- omega: array_like, shape (n_omega,)
The frequencies at which to calculate the filter functions.
- n_oper_identifiers: array_like, optional
The identifiers of the noise operators for which to evaluate the error transfer matrix. The default is all. Note that, since in general contributions from different noise operators won’t commute, not selecting all noise operators results in neglecting terms of order \(\xi^4\).
- second_order: bool, optional
Also take into account the frequency shifts \(\Delta\) that correspond to second order Magnus expansion and constitute unitary terms. Default
False.- cumulant_function: ndarray, shape ([[G, G,] n_nops,] n_nops, d**2, d**2)
A precomputed cumulant function. If given, pulse, spectrum, omega are not required.
- show_progressbar: bool, optional
Show a progress bar for the calculation of the control matrix.
- memory_parsimonious: bool, optional
Trade memory footprint for performance. See
calculate_decay_amplitudes(). The default isFalse.- cache_intermediates: bool, optional
Keep and return intermediate terms of the calculation of the control matrix (if it is not already cached) that can be reused for second order or gradients. Can consume large amount of memory, but speed up the calculation.
- Returns:
- error_transfer_matrix: ndarray, shape (d**2, d**2)
The error transfer matrix. The individual noise operator contributions are summed up before exponentiating as they might not commute.
See also
calculate_cumulant_functionCalculate the cumulant function \(\mathcal{K}\)
calculate_decay_amplitudesCalculate the \(\Gamma_{\alpha\beta,kl}\)
infidelityCalculate only infidelity of a pulse.
Notes
The error transfer matrix is given by
\[\tilde{\mathcal{U}} = \exp\mathcal{K}(\tau)\]with \(\mathcal{K}(\tau)\) the cumulant function (see
calculate_cumulant_function()). For Gaussian noise this expression is exact when taking into account the decay amplitudes \(\Gamma\) and frequency shifts \(\Delta\). As the latter effects coherent errors it can be neglected if we assume that the experimenter has calibrated their pulse.For non-Gaussian noise the expression above is perturbative and includes noise up to order \(\xi^2\) and hence \(\tilde{\mathcal{U}} = \mathbb{1} + \mathcal{K}(\tau) + \mathcal{O}(\xi^4)\) (although it is evaluated as a matrix exponential in any case).
Given the above expression of the error transfer matrix, the entanglement fidelity is given by
\[\mathcal{F}_\mathrm{e} = \frac{1}{d^2}\mathrm{tr}\,\tilde{\mathcal{U}}.\]
- infidelity(pulse: PulseSequence, spectrum: Sequence[float] | Callable, omega: Sequence[float] | Dict[str, int | str], n_oper_identifiers: Sequence[str] | None = None, which: str = 'total', show_progressbar: bool = False, cache_intermediates: bool = False, return_smallness: bool = False, test_convergence: bool = False) ndarray | Any[source]
Calculate the leading order entanglement infidelity.
This function calculates the infidelity approximately from the leading peturbation (see Notes). To compute it exactly for Gaussian noise and vanishing coherent errors (second order Magnus terms), use
error_transfer_matrix()to obtain it from the full process matrix.- Parameters:
- pulse: PulseSequence
The
PulseSequenceinstance for which to calculate the infidelity for.- spectrum: array_like, shape ([[n_nops,] n_nops,] omega) or callable
The noise power spectral density in units of inverse frequencies as an array of shape (n_omega,), (n_nops, n_omega), or (n_nops, n_nops, n_omega). In the first case, the same spectrum is taken for all noise operators, in the second, it is assumed that there are no correlations between different noise sources and thus there is one spectrum for each noise operator. In the third and most general case, there may be a spectrum for each pair of noise operators corresponding to the correlations between them. n_nops is the number of noise operators considered and should be equal to
len(n_oper_identifiers).See Notes for a discussion on one- and two-sided power spectral densities.
If test_convergence is
True, a function handle to compute the power spectral density from a sequence of frequencies is expected.- omega: array_like or dict
The frequencies at which the integration is to be carried out. If test_convergence is
True, a dict with possible keys (‘omega_IR’, ‘omega_UV’, ‘spacing’, ‘n_min’, ‘n_max’, ‘n_points’), where all entries are integers except forspacingwhich should be a string, either ‘linear’ or ‘log’. ‘n_points’ controls how many steps are taken.- n_oper_identifiers: array_like, optional
The identifiers of the noise operators for which to calculate the infidelity contribution. If given, the infidelities for each noise operator will be returned. Otherwise, all noise operators will be taken into account.
- which: str, optional
Which infidelities should be calculated, may be either ‘total’ (default) or ‘correlations’. In the former case, one value is returned for each noise operator, corresponding to the total infidelity of the pulse (or pulse sequence). In the latter, an array of infidelities is returned where element (i,j) corresponds to the infidelity contribution of the correlations between pulses i and j (see Notes). Note that this option is only available if the pulse correlation filter functions have been computed during concatenation (see
calculate_pulse_correlation_filter_function()andconcatenate()).- show_progressbar: bool, optional
Show a progressbar for the calculation of the control matrix.
- cache_intermediates: bool, optional
Keep and return intermediate terms of the calculation of the control matrix (if it is not already cached) that can be reused for second order or gradients. The default is False.
- return_smallness: bool, optional
Return the smallness parameter \(\xi\) for the given spectrum.
- test_convergence: bool, optional
Test the convergence of the integral with respect to the number of frequency samples. Returns the number of frequency samples and the corresponding fidelities. See spectrum and omega for more information.
- Returns:
- infid: ndarray, shape ([[G, G,], n_nops,] n_nops)
Array with the infidelity contributions for each spectrum spectrum on the last axis or axes, depending on the shape of spectrum and which. If
whichiscorrelations, the first two axes are the individual pulse contributions. If spectrum is 2-d (3-d), the last axis (two axes) are the individual spectral contributions. Only if test_convergence isFalse.- n_samples: array_like
Array with number of frequency samples used for convergence test. Only if test_convergence is
True.- convergence_infids: array_like
Array with infidelities calculated in convergence test. Only if test_convergence is
True.
See also
calculate_decay_amplitudesCalculate the full matrix of first order terms.
error_transfer_matrixCalculate the full process matrix.
plotting.plot_infidelity_convergenceConvenience function to plot results.
pulse_sequence.concatenateConcatenate
PulseSequenceobjects.calculate_pulse_correlation_filter_function.
Notes
The infidelity is given by
\[\begin{split}\mathcal{I}_{\alpha\beta} &= 1 - \frac{1}{d^2}\mathrm{tr}\:\tilde{\mathcal{U}} \\ &= \frac{1}{d}\int_{-\infty}^{\infty} \frac{\mathrm{d}\omega}{2\pi}S_{\alpha\beta}(\omega) F_{\alpha\beta}(\omega) + \mathcal{O}\big(\xi^4\big) \\ &= \sum_{g,g'=1}^G \mathcal{I}_{\alpha\beta}^{(gg')}\end{split}\]with \(S_{\alpha\beta}(\omega)\) the two-sided noise spectral density and \(F_{\alpha\beta}(\omega)\) the first-order filter function for noise sources \(\alpha,\beta\). The noise spectrum may include correlated noise sources, that is, its entry at \((\alpha,\beta)\) corresponds to the correlations between sources \(\alpha\) and \(\beta\). \(\mathcal{I}_{\alpha\beta}^{(gg')}\) are the correlation infidelities that can be computed by setting
which='correlations'.One- and two-sided spectral densities
Since the real (imaginary) part of filter function \(F(\omega)\) is even (odd), it does not matter for integral whether \(S(\omega)\) is taken to be the one- or two-sided spectral density. However, care should be taken that, if it is one or the other, the frequencies \(\omega\) are positive or symmetric about zero, respectively.
To convert between one- and two-sided PSDs, use the following relationship:
\[S_\mathrm{onesided}(\omega) = 2 S_\mathrm{twosided}(\omega).\]Conversion to the Average Gate Infidelity (AGI)
To convert the entanglement infidelity to the average gate infidelity, use the following relation given by Horodecki et al. [Hor99] and Nielsen [Nie02]:
\[\mathcal{I}_\mathrm{avg} = \frac{d}{d+1}\mathcal{I}.\]Goodness of approximation
The smallness parameter is given by
\[\xi^2 = \sum_\alpha\left[ \lvert\lvert B_\alpha\rvert\rvert^2 \int_{-\infty}^\infty\frac{\mathrm{d}\omega}{2\pi} S_\alpha(\omega)\left(\sum_gs_\alpha^{(g)}\Delta t_g \right)^2 \right].\]Note that in practice, the integral is only evaluated on the interval \(\omega\in[\omega_\mathrm{min},\omega_\mathrm{max}]\).
References
[Hor99]Horodecki, M., Horodecki, P., & Horodecki, R. (1999). General teleportation channel, singlet fraction, and quasidistillation. Physical Review A - Atomic, Molecular, and Optical Physics, 60(3), 1888–1898. https://doi.org/10.1103/PhysRevA.60.1888
[Nie02]Nielsen, M. A. (2002). A simple formula for the average gate fidelity of a quantum dynamical operation. Physics Letters, Section A: General, Atomic and Solid State Physics, 303(4), 249–252. https://doi.org/10.1016/S0375-9601(02)01272-0