FluxQubit#
- class scqubits.core.flux_qubit.FluxQubit(EJ1, EJ2, EJ3, ECJ1, ECJ2, ECJ3, ECg1, ECg2, ng1, ng2, flux, ncut, truncated_dim=6, id_str=None)[source]#
Flux Qubit
[1] Orlando et al., Physical Review B, 60, 15398 (1999). https://link.aps.org/doi/10.1103/PhysRevB.60.15398The original flux qubit as defined in [1], where the junctions are allowed to have varying junction energies and capacitances to allow for junction asymmetry. Typically, one takes \(E_{J1}=E_{J2}=E_J\), and \(E_{J3}=\alpha E_J\) where \(0\le \alpha \le 1\). The same relations typically hold for the junction capacitances. The Hamiltonian is given by
\[\begin{split}H_\text{flux}=&(n_{i}-n_{gi})4(E_\text{C})_{ij}(n_{j}-n_{gj}) \\ -&E_{J}\cos\phi_{1}-E_{J}\cos\phi_{2}-\alpha E_{J}\cos(2\pi f + \phi_{1} - \phi_{2}),\end{split}\]where \(i,j\in\{1,2\}\) is represented in the charge basis for both degrees of freedom. Initialize with, for example:
EJ = 35.0 alpha = 0.6 flux_qubit = scq.FluxQubit(EJ1 = EJ, EJ2 = EJ, EJ3 = alpha*EJ, ECJ1 = 1.0, ECJ2 = 1.0, ECJ3 = 1.0/alpha, ECg1 = 50.0, ECg2 = 50.0, ng1 = 0.0, ng2 = 0.0, flux = 0.5, ncut = 10)
- Parameters
EJ1 (float) – Josephson energy of the ith junction EJ1 = EJ2, with EJ3 = alpha * EJ1 and alpha <= 1
EJ2 (float) – Josephson energy of the ith junction EJ1 = EJ2, with EJ3 = alpha * EJ1 and alpha <= 1
EJ3 (float) – Josephson energy of the ith junction EJ1 = EJ2, with EJ3 = alpha * EJ1 and alpha <= 1
ECJ1 (float) – charging energy associated with the ith junction
ECJ2 (float) – charging energy associated with the ith junction
ECJ3 (float) – charging energy associated with the ith junction
ECg1 (float) – charging energy associated with the capacitive coupling to ground for the two islands
ECg2 (float) – charging energy associated with the capacitive coupling to ground for the two islands
ng1 (float) – offset charge associated with island i
ng2 (float) – offset charge associated with island i
flux (float) – magnetic flux through the circuit loop, measured in units of the flux quantum
ncut (int) – charge number cutoff for the charge on both islands n, n = -ncut, …, ncut
truncated_dim (
int) – desired dimension of the truncated quantum system; expected: truncated_dim > 1id_str (
Optional[str]) – optional string by which this instance can be referred to in HilbertSpace and ParameterSweep. If not provided, an id is auto-generated.
- Return type
Methods
Return the charging energy matrix
FluxQubit.__init__(EJ1, EJ2, EJ3, ECJ1, ...)FluxQubit.broadcast(event, **kwargs)Request a broadcast from CENTRAL_DISPATCH reporting event.
Return operator \(\cos \phi_1\) in the charge basis
Return operator \(\cos \phi_2\) in the charge basis
Use ipywidgets to create a new class instance
FluxQubit.create_from_file(filename)Read initdata and spectral data from file, and use those to create a new SpectrumData object.
Returns operator representing a derivative of the Hamiltonian with respect to EJ1.
Returns operator representing a derivative of the Hamiltonian with respect to EJ2.
Returns operator representing a derivative of the Hamiltonian with respect to EJ3.
Return dictionary with default parameter values for initialization of class instance
FluxQubit.deserialize(io_data)Take the given IOData and return an instance of the described class, initialized with the data stored in io_data.
Return a list of noise channels that are used when calculating the effective noise (i.e.
FluxQubit.eigensys([evals_count, filename, ...])Calculates eigenvalues and corresponding eigenvectors using scipy.linalg.eigh.
FluxQubit.eigenvals([evals_count, filename, ...])Calculates eigenvalues using scipy.linalg.eigh, returns numpy array of eigenvalues.
Return operator \(e^{i\phi_1}\) in the charge basis.
Return operator \(e^{i\phi_2}\) in the charge basis.
FluxQubit.filewrite(filename)Convenience method bound to the class.
FluxQubit.get_dispersion_vs_paramvals(...[, ...])Calculates eigenvalues/eigenstates for a varying system parameter, given an array of parameter values.
Returns dict appropriate for creating/initializing a new Serializable object.
Calculates matrix elements for a varying system parameter, given an array of parameter values.
Returns a list of all operator names for the quantum system.
FluxQubit.get_spectrum_vs_paramvals(...[, ...])Calculates eigenvalues/eigenstates for a varying system parameter, given an array of parameter values.
Return Hamiltonian in basis obtained by employing charge basis for both degrees of freedom
Return Hilbert space dimension.
Return the kinetic energy matrix.
FluxQubit.matrixelement_table(operator[, ...])Returns table of matrix elements for operator with respect to the eigenstates of the qubit.
Return charge number operator conjugate to \(\phi_1\)
Return charge number operator conjugate to \(\phi_2\)
FluxQubit.plot_coherence_vs_paramvals(...[, ...])Show plots of coherence for various channels supported by the qubit as they vary as a function of a changing parameter.
Generates a simple plot of a set of curves representing the charge or flux dispersion of transition energies.
FluxQubit.plot_evals_vs_paramvals(...[, ...])Generates a simple plot of a set of eigenvalues as a function of one parameter.
FluxQubit.plot_matelem_vs_paramvals(...[, ...])Generates a simple plot of a set of eigenvalues as a function of one parameter.
FluxQubit.plot_matrixelements(operator[, ...])Plots matrix elements for operator, given as a string referring to a class method that returns an operator matrix.
FluxQubit.plot_potential([phi_grid, ...])Draw contour plot of the potential energy.
Plot effective \(T_1\) coherence time (rate) as a function of changing parameter.
Plot effective \(T_2\) coherence time (rate) as a function of changing parameter.
FluxQubit.plot_wavefunction([esys, which, ...])Plots 2d phase-basis wave function.
FluxQubit.potential(phi1, phi2)Return value of the potential energy at phi1 and phi2, disregarding constants.
Return the potential energy matrix for the potential.
FluxQubit.receive(event, sender, **kwargs)Receive a message from CENTRAL_DISPATCH and initiate action on it.
Convert the content of the current class instance into IOData format.
FluxQubit.set_and_return(attr_name, value)Allows to set an attribute after which self is returned. This is useful for doing something like example::.
FluxQubit.set_params(**kwargs)Set new parameters through the provided dictionary.
Return operator \(\sin \phi_1\) in the charge basis
Return operator \(\sin \phi_2\) in the charge basis
Return a list of supported noise channels
FluxQubit.t1(i, j, noise_op, spectral_density)Calculate the transition time (or rate) using Fermi's Golden Rule due to a noise channel with a spectral density spectral_density and system noise operator noise_op.
FluxQubit.t1_capacitive([i, j, Q_cap, T, ...])\(T_1\) due to dielectric dissipation in the Josephson junction capacitances.
FluxQubit.t1_charge_impedance([i, j, Z, T, ...])Noise due to charge coupling to an impedance (such as a transmission line).
FluxQubit.t1_effective([noise_channels, ...])Calculate the effective \(T_1\) time (or rate).
FluxQubit.t1_flux_bias_line([i, j, M, Z, T, ...])Noise due to a bias flux line.
FluxQubit.t1_inductive([i, j, Q_ind, T, ...])\(T_1\) due to inductive dissipation in a superinductor.
FluxQubit.t1_quasiparticle_tunneling([i, j, ...])Noise due to quasiparticle tunneling across a Josephson junction.
FluxQubit.t2_effective([noise_channels, ...])Calculate the effective \(T_2\) time (or rate).
FluxQubit.tphi_1_over_f(A_noise, i, j, noise_op)Calculate the 1/f dephasing time (or rate) due to arbitrary noise source.
FluxQubit.tphi_1_over_f_cc([A_noise, i, j, ...])Calculate the 1/f dephasing time (or rate) due to critical-current noise from all three Josephson junctions \(EJ1\), \(EJ2\) and \(EJ3\).
FluxQubit.tphi_1_over_f_cc1([A_noise, i, j, ...])Calculate the 1/f dephasing time (or rate) due to critical current noise of junction associated with Josephson energy \(EJ1\).
FluxQubit.tphi_1_over_f_cc2([A_noise, i, j, ...])Calculate the 1/f dephasing time (or rate) due to critical current noise of junction associated with Josephson energy \(EJ2\).
FluxQubit.tphi_1_over_f_cc3([A_noise, i, j, ...])Calculate the 1/f dephasing time (or rate) due to critical current noise of junction associated with Josephson energy \(EJ3\).
FluxQubit.tphi_1_over_f_flux([A_noise, i, ...])Calculate the 1/f dephasing time (or rate) due to flux noise.
FluxQubit.tphi_1_over_f_ng([A_noise, i, j, ...])Calculate the 1/f dephasing time (or rate) due to charge noise.
FluxQubit.wavefunction([esys, which, phi_grid])Return a flux qubit wave function in phi1, phi2 basis
FluxQubit.widget([params])Use ipywidgets to modify parameters of class instance
Attributes
ECJ1Descriptor class for properties that are to be monitored for changes.
ECJ2Descriptor class for properties that are to be monitored for changes.
ECJ3Descriptor class for properties that are to be monitored for changes.
ECg1Descriptor class for properties that are to be monitored for changes.
ECg2Descriptor class for properties that are to be monitored for changes.
EJ1Descriptor class for properties that are to be monitored for changes.
EJ2Descriptor class for properties that are to be monitored for changes.
EJ3Descriptor class for properties that are to be monitored for changes.
fluxDescriptor class for properties that are to be monitored for changes.
id_strncutDescriptor class for properties that are to be monitored for changes.
ng1Descriptor class for properties that are to be monitored for changes.
ng2Descriptor class for properties that are to be monitored for changes.
truncated_dimDescriptor class for properties that are to be monitored for changes.
- broadcast(event, **kwargs)#
Request a broadcast from CENTRAL_DISPATCH reporting event.
- Parameters
event (
str) – event name from EVENTS**kwargs –
- Return type
None
- cos_phi_1_operator()[source]#
Return operator \(\cos \phi_1\) in the charge basis
- Return type
ndarray
- cos_phi_2_operator()[source]#
Return operator \(\cos \phi_2\) in the charge basis
- Return type
ndarray
- classmethod create()#
Use ipywidgets to create a new class instance
- Return type
- classmethod create_from_file(filename)#
Read initdata and spectral data from file, and use those to create a new SpectrumData object.
- Returns
new SpectrumData object, initialized with data read from file
- Return type
- Parameters
filename (str) –
- d_hamiltonian_d_EJ1()[source]#
Returns operator representing a derivative of the Hamiltonian with respect to EJ1.
- Return type
ndarray
- d_hamiltonian_d_EJ2()[source]#
Returns operator representing a derivative of the Hamiltonian with respect to EJ2.
- Return type
ndarray
- d_hamiltonian_d_EJ3()[source]#
Returns operator representing a derivative of the Hamiltonian with respect to EJ3.
- Return type
ndarray
- static default_params()[source]#
Return dictionary with default parameter values for initialization of class instance
- Return type
Dict[str,Any]
- classmethod deserialize(io_data)#
Take the given IOData and return an instance of the described class, initialized with the data stored in io_data.
- Return type
- Parameters
io_data (IOData) –
- classmethod effective_noise_channels()#
Return a list of noise channels that are used when calculating the effective noise (i.e. via t1_effective and t2_effective.
- Return type
List[str]
- eigensys(evals_count=6, filename=None, return_spectrumdata=False)#
Calculates eigenvalues and corresponding eigenvectors using scipy.linalg.eigh. Returns two numpy arrays containing the eigenvalues and eigenvectors, respectively.
- Parameters
evals_count (
int) – number of desired eigenvalues/eigenstates (default value = 6)filename (
Optional[str]) – path and filename without suffix, if file output desired (default value = None)return_spectrumdata (
bool) – if set to true, the returned data is provided as a SpectrumData object (default value = False)
- Return type
Union[Tuple[ndarray,ndarray],SpectrumData]- Returns
eigenvalues, eigenvectors as numpy arrays or in form of a SpectrumData object
- eigenvals(evals_count=6, filename=None, return_spectrumdata=False)#
Calculates eigenvalues using scipy.linalg.eigh, returns numpy array of eigenvalues.
- Parameters
evals_count (
int) – number of desired eigenvalues/eigenstates (default value = 6)filename (
Optional[str]) – path and filename without suffix, if file output desired (default value = None)return_spectrumdata (
bool) – if set to true, the returned data is provided as a SpectrumData object (default value = False)
- Return type
Union[SpectrumData,ndarray]- Returns
eigenvalues as ndarray or in form of a SpectrumData object
- exp_i_phi_1_operator()[source]#
Return operator \(e^{i\phi_1}\) in the charge basis.
- Return type
ndarray
- exp_i_phi_2_operator()[source]#
Return operator \(e^{i\phi_2}\) in the charge basis.
- Return type
ndarray
- filewrite(filename)#
Convenience method bound to the class. Simply accesses the write function.
- Return type
None- Parameters
filename (str) –
- get_dispersion_vs_paramvals(dispersion_name, param_name, param_vals, ref_param=None, transitions=(0, 1), levels=None, point_count=50, num_cpus=None)#
Calculates eigenvalues/eigenstates for a varying system parameter, given an array of parameter values. Returns a SpectrumData object with energy_data[n] containing eigenvalues calculated for parameter value param_vals[n].
- Parameters
dispersion_name (
str) – parameter inducing the dispersion, typically ‘ng’ or ‘flux’ (will be scanned over range from 0 to 1)param_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged inref_param (
Optional[str]) – optional, name of parameter to use as reference for the parameter value; e.g., to compute charge dispersion vs. EJ/EC, use EJ as param_name and EC as ref_paramtransitions (
Union[Tuple[int,int],Tuple[Tuple[int,int],...]]) – integer tuple or tuples specifying for which transitions dispersion is to be calculated (default: = (0,1))levels (
Union[int,Tuple[int,...],None]) – tuple specifying levels (rather than transitions) for which dispersion should be plotted; overrides transitions parameter when givenpoint_count (
int) – number of points scanned for the dispersion parameter for determining min and max values of transition energies (default: 50)num_cpus (
Optional[int]) – number of cores to be used for computation (default value: settings.NUM_CPUS)
- Return type
- get_initdata()#
Returns dict appropriate for creating/initializing a new Serializable object.
- Return type
Dict[str,Any]
- get_matelements_vs_paramvals(operator, param_name, param_vals, evals_count=6, num_cpus=None)#
Calculates matrix elements for a varying system parameter, given an array of parameter values. Returns a SpectrumData object containing matrix element data, eigenvalue data, and eigenstate data..
- Parameters
operator (
str) – name of class method in string form, returning operator matrixparam_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged inevals_count (
int) – number of desired eigenvalues (sorted from smallest to largest) (default value = 6)num_cpus (
Optional[int]) – number of cores to be used for computation (default value: settings.NUM_CPUS)
- Return type
- classmethod get_operator_names()#
Returns a list of all operator names for the quantum system. Note that this list omits any operators that start with “_”.
- Parameters
subsys – Class instance of quantum system
- Return type
List[str]- Returns
list of operator names
- get_spectrum_vs_paramvals(param_name, param_vals, evals_count=6, subtract_ground=False, get_eigenstates=False, filename=None, num_cpus=None)#
Calculates eigenvalues/eigenstates for a varying system parameter, given an array of parameter values. Returns a SpectrumData object with energy_data[n] containing eigenvalues calculated for parameter value param_vals[n].
- Parameters
param_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged inevals_count (
int) – number of desired eigenvalues (sorted from smallest to largest) (default value = 6)subtract_ground (
bool) – if True, eigenvalues are returned relative to the ground state eigenvalue (default value = False)get_eigenstates (
bool) – return eigenstates along with eigenvalues (default value = False)filename (
Optional[str]) – file name if direct output to disk is wantednum_cpus (
Optional[int]) – number of cores to be used for computation (default value: settings.NUM_CPUS)
- Return type
- hamiltonian()[source]#
Return Hamiltonian in basis obtained by employing charge basis for both degrees of freedom
- Return type
ndarray
- matrixelement_table(operator, evecs=None, evals_count=6, filename=None, return_datastore=False)#
Returns table of matrix elements for operator with respect to the eigenstates of the qubit. The operator is given as a string matching a class method returning an operator matrix. E.g., for an instance trm of Transmon, the matrix element table for the charge operator is given by trm.op_matrixelement_table(‘n_operator’). When esys is set to None, the eigensystem is calculated on-the-fly.
- Parameters
operator (
str) – name of class method in string form, returning operator matrix in qubit-internal basis.evecs (
Optional[ndarray]) – if not provided, then the necessary eigenstates are calculated on the flyevals_count (
int) – number of desired matrix elements, starting with ground state (default value = 6)filename (
Optional[str]) – output file namereturn_datastore (
bool) – if set to true, the returned data is provided as a DataStore object (default value = False)
- Return type
Union[DataStore,ndarray]
- plot_coherence_vs_paramvals(param_name, param_vals, noise_channels=None, common_noise_options=None, spectrum_data=None, scale=1, num_cpus=None, **kwargs)#
Show plots of coherence for various channels supported by the qubit as they vary as a function of a changing parameter.
For example, assuming qubit is a qubit object with flux being one of its parameters, one can see how coherence due to various noise channels vary as the flux changes:
qubit.plot_coherence_vs_paramvals(param_name='flux', param_vals=np.linspace(-0.5, 0.5, 100), scale=1e-3, ylabel=r"$\mu s$");
- Parameters
param_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged innoise_channels (
Union[str,List[str],List[Tuple[str,Dict]],None]) – channels to be plotted, if None then noise channels given by supported_noise_channels are usedcommon_noise_options (
Optional[Dict]) – common options used when calculating coherence timesspectrum_data (
Optional[SpectrumData]) – spectral data used during noise calculationsscale (float) – a number that all data is multiplied by before being plotted
num_cpus (
Optional[int]) – number of cores to be used for computation
- Return type
Figure, Axes
- plot_dispersion_vs_paramvals(dispersion_name, param_name, param_vals, ref_param=None, transitions=(0, 1), levels=None, point_count=50, num_cpus=None, **kwargs)#
Generates a simple plot of a set of curves representing the charge or flux dispersion of transition energies.
- Parameters
dispersion_name (
str) – parameter inducing the dispersion, typically ‘ng’ or ‘flux’ (will be scanned over range from 0 to 1)param_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged inref_param (
Optional[str]) – optional, name of parameter to use as reference for the parameter value; e.g., to compute charge dispersion vs. EJ/EC, use EJ as param_name and EC as ref_paramtransitions (
Union[Tuple[int,int],Tuple[Tuple[int,int],...]]) – integer tuple or tuples specifying for which transitions dispersion is to be calculated (default: = (0,1))levels (
Union[int,Tuple[int,...],None]) – int or tuple specifying level(s) (rather than transitions) for which dispersion should be plotted; overrides transitions parameter when givenpoint_count (
int) – number of points scanned for the dispersion parameter for determining min and max values of transition energies (default: 50)num_cpus (
Optional[int]) – number of cores to be used for computation (default value: settings.NUM_CPUS)**kwargs – standard plotting option (see separate documentation)
- Return type
Tuple[Figure,Axes]
- plot_evals_vs_paramvals(param_name, param_vals, evals_count=6, subtract_ground=False, num_cpus=None, **kwargs)#
Generates a simple plot of a set of eigenvalues as a function of one parameter. The individual points correspond to the a provided array of parameter values.
- Parameters
param_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged inevals_count (
int) – number of desired eigenvalues (sorted from smallest to largest) (default value = 6)subtract_ground (
bool) – whether to subtract ground state energy from all eigenvalues (default value = False)num_cpus (
Optional[int]) – number of cores to be used for computation (default value: settings.NUM_CPUS)**kwargs – standard plotting option (see separate documentation)
- Return type
Tuple[Figure,Axes]
- plot_matelem_vs_paramvals(operator, param_name, param_vals, select_elems=4, mode='abs', num_cpus=None, **kwargs)#
Generates a simple plot of a set of eigenvalues as a function of one parameter. The individual points correspond to the a provided array of parameter values.
- Parameters
operator (
str) – name of class method in string form, returning operator matrixparam_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged inselect_elems (
Union[int,List[Tuple[int,int]]]) – either maximum index of desired matrix elements, or list [(i1, i2), (i3, i4), …] of index tuples for specific desired matrix elements (default value = 4)mode (
str) – idx_entry from MODE_FUNC_DICTIONARY, e.g., ‘abs’ for absolute value (default value = ‘abs’)num_cpus (
Optional[int]) – number of cores to be used for computation (default value: settings.NUM_CPUS)**kwargs – standard plotting option (see separate documentation)
- Return type
Tuple[Figure,Axes]
- plot_matrixelements(operator, evecs=None, evals_count=6, mode='abs', show_numbers=False, show3d=True, **kwargs)#
Plots matrix elements for operator, given as a string referring to a class method that returns an operator matrix. E.g., for instance trm of Transmon, the matrix element plot for the charge operator n is obtained by trm.plot_matrixelements(‘n’). When esys is set to None, the eigensystem with which eigenvectors is calculated.
- Parameters
operator (
str) – name of class method in string form, returning operator matrixevecs (
Optional[ndarray]) – eigensystem data of evals, evecs; eigensystem will be calculated if set to None (default value = None)evals_count (
int) – number of desired matrix elements, starting with ground state (default value = 6)mode (
str) – idx_entry from MODE_FUNC_DICTIONARY, e.g., ‘abs’ for absolute value (default)show_numbers (
bool) – determines whether matrix element values are printed on top of the plot (default: False)show3d (
bool) – whether to show a 3d skyscraper plot of the matrix alongside the 2d plot (default: True)**kwargs – standard plotting option (see separate documentation)
- Return type
Union[Tuple[Figure,Tuple[Axes,Axes]],Tuple[Figure,Axes]]
- plot_potential(phi_grid=None, contour_vals=None, **kwargs)[source]#
Draw contour plot of the potential energy.
- Parameters
phi_grid (
Optional[Grid1d]) – used for setting a custom grid for phi; if None use self._default_gridcontour_vals (
Optional[ndarray]) – specific contours to draw**kwargs – plot options
- Return type
Tuple[Figure,Axes]
- plot_t1_effective_vs_paramvals(param_name, param_vals, noise_channels=None, common_noise_options=None, spectrum_data=None, get_rate=False, scale=1, num_cpus=None, **kwargs)#
Plot effective \(T_1\) coherence time (rate) as a function of changing parameter.
The effective \(T_1\) is calculated by considering a variety of depolarizing noise channels, according to the formula:
\[\frac{1}{T_{1}^{\rm eff}} = \frac{1}{2} \sum_k \frac{1}{T_{1}^{k}}\]where \(k\) runs over the channels that can contribute to the effective noise. By default all the depolarizing noise channels given by the method effective_noise_channels are included.
For example, assuming qubit is a qubit object with flux being one of its parameters, one can see how the effective \(T_1\) varies as the flux changes:
qubit.plot_t1_effective_vs_paramvals(param_name='flux', param_vals=np.linspace(-0.5, 0.5, 100), );
- Parameters
param_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged innoise_channels (
Union[str,List[str],List[Tuple[str,Dict]],None]) – channels to be plotted, if None then noise channels given by supported_noise_channels are usedcommon_noise_options (
Optional[Dict]) – common options used when calculating coherence timesspectrum_data (
Optional[SpectrumData]) – spectral data used during noise calculationsget_rate (
bool) – determines if rate or time should be plottedscale (
float) – a number that all data is multiplied by before being plottednum_cpus (
Optional[int]) – number of cores to be used for computation
- Return type
Figure, Axes
- plot_t2_effective_vs_paramvals(param_name, param_vals, noise_channels=None, common_noise_options=None, spectrum_data=None, get_rate=False, scale=1, num_cpus=None, **kwargs)#
Plot effective \(T_2\) coherence time (rate) as a function of changing parameter.
The effective \(T_2\) is calculated from both pure dephasing channels, as well as depolarization channels, according to the formula:
\[\frac{1}{T_{2}^{\rm eff}} = \sum_k \frac{1}{T_{\phi}^{k}} + \frac{1}{2} \sum_j \frac{1}{T_{1}^{j}}\]where \(k\) (\(j\)) run over the relevant pure dephasing ( depolarization) channels that can contribute to the effective noise. By default all noise channels given by the method effective_noise_channels are included.
For example, assuming qubit is a qubit object with flux being one of its parameters, one can see how the effective \(T_2\) varies as the flux changes:
qubit.plot_t2_effective_vs_paramvals(param_name='flux', param_vals=np.linspace(-0.5, 0.5, 100), );
- Parameters
param_name (
str) – name of parameter to be variedparam_vals (
ndarray) – parameter values to be plugged innoise_channels (
Union[str,List[str],List[Tuple[str,Dict]],None]) – channels to be plotted, if None then noise channels given by supported_noise_channels are usedcommon_noise_options (
Optional[Dict]) – common options used when calculating coherence timesspectrum_data (
Optional[SpectrumData]) – spectral data used during noise calculationsget_rate (
bool) – determines if rate or time should be plottedscale (
float) – a number that all data is multiplied by before being plottednum_cpus (
Optional[int]) – number of cores to be used for computation
- Return type
Figure, Axes
- plot_wavefunction(esys=None, which=0, phi_grid=None, mode='abs', zero_calibrate=True, **kwargs)[source]#
Plots 2d phase-basis wave function.
- Parameters
esys (
Optional[Tuple[ndarray,ndarray]]) – eigenvalues, eigenvectors as obtained from .eigensystem()which (
int) – index of wave function to be plotted (default value = (0)phi_grid (
Optional[Grid1d]) – used for setting a custom grid for phi; if None use self._default_gridmode (
str) – choices as specified in constants.MODE_FUNC_DICT (default value = ‘abs_sqr’)zero_calibrate (
bool) – if True, colors are adjusted to use zero wavefunction amplitude as the neutral color in the palette**kwargs – plot options
- Return type
Tuple[Figure,Axes]
- potential(phi1, phi2)[source]#
Return value of the potential energy at phi1 and phi2, disregarding constants.
- Return type
ndarray- Parameters
phi1 (ndarray) –
phi2 (ndarray) –
- receive(event, sender, **kwargs)#
Receive a message from CENTRAL_DISPATCH and initiate action on it.
- Parameters
event (
str) – event name from EVENTSsender (
DispatchClient) – original sender reporting the event**kwargs –
- Return type
None
- serialize()#
Convert the content of the current class instance into IOData format.
- Return type
- set_and_return(attr_name, value)#
Allows to set an attribute after which self is returned. This is useful for doing something like example:
qubit.set_and_return('flux', 0.23).some_method()
instead of example:
qubit.flux=0.23 qubit.some_method()
- Parameters
attr_name (
str) – name of class attribute in string formvalue (
Any) – value that the attribute is to be set to
- Return type
- Returns
self
- set_params(**kwargs)#
Set new parameters through the provided dictionary.
- sin_phi_1_operator()[source]#
Return operator \(\sin \phi_1\) in the charge basis
- Return type
ndarray
- sin_phi_2_operator()[source]#
Return operator \(\sin \phi_2\) in the charge basis
- Return type
ndarray
- classmethod supported_noise_channels()[source]#
Return a list of supported noise channels
- Return type
List[str]
- t1(i, j, noise_op, spectral_density, total=True, esys=None, get_rate=False)#
Calculate the transition time (or rate) using Fermi’s Golden Rule due to a noise channel with a spectral density spectral_density and system noise operator noise_op. Mathematically, it reads:
\[\frac{1}{T_1} = \frac{1}{\hbar^2} |\langle i| A_{\rm noise} | j \rangle|^2 S(\omega)\]We assume that the qubit energies (or the passed in eigenspectrum) has units of frequency (and not angular frequency).
The spectral_density argument should be a callable object (typically a function) of one argument, which is assumed to be an angular frequency (in the units currently set as system units.
- Parameters
i (int >=0) – state index that along with j defines a transition (i->j)
j (int >=0) – state index that along with i defines a transition (i->j)
noise_op (
Union[ndarray,csc_matrix]) – noise operatorspectral_density (
Callable) – defines a spectral density, must take one argument: omega (assumed to be in units of 2 pi * <system units>)total (
bool) – if False return a time/rate associated with a transition from state i to state j. if True return a time/rate associated with both i to j and j to i transitionsesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- t1_capacitive(i=1, j=0, Q_cap=None, T=0.015, total=True, esys=None, get_rate=False)#
\(T_1\) due to dielectric dissipation in the Josephson junction capacitances.
References: Nguyen et al (2019), Smith et al (2020)
- Parameters
i (int >=0) – state index that along with j defines a transition (i->j)
j (int >=0) – state index that along with i defines a transition (i->j)
Q_cap (
Union[float,Callable,None]) – capacitive quality factor; a fixed value or function of omegaT (
float) – temperature in Kelvintotal (
bool) – if False return a time/rate associated with a transition from state i to state j. if True return a time/rate associated with both i to j and j to i transitionsesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate –
- decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate
in inverse units.
- Return type
float
- t1_charge_impedance(i=1, j=0, Z=50, T=0.015, total=True, esys=None, get_rate=False)#
Noise due to charge coupling to an impedance (such as a transmission line).
References: Schoelkopf et al (2003), Ithier et al (2005)
- Parameters
i (int >=0) – state index that along with j defines a transition (i->j)
j (int >=0) – state index that along with i defines a transition (i->j)
Z (
Union[float,Callable]) – impedance; a fixed value or function of omegaT (
float) – temperature in Kelvintotal (
bool) – if False return a time/rate associated with a transition from state i to state j. if True return a time/rate associated with both i to j and j to i transitionsesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- t1_effective(noise_channels=None, common_noise_options=None, esys=None, get_rate=False, **kwargs)#
Calculate the effective \(T_1\) time (or rate).
The effective \(T_1\) is calculated by considering a variety of depolarizing noise channels, according to the formula:
\[\frac{1}{T_{1}^{\rm eff}} = \frac{1}{2} \sum_k \frac{1}{T_{1}^{k}}\]where \(k\) runs over the channels that can contribute to the effective noise. By default all the depolarizing noise channels given by the method effective_noise_channels are included. Users can also provide specific noise channels, with selected options, to be included in the effective \(T_1\) calculation. For example, assuming qubit is a qubit object, can can execute:
tune_tmon.t1_effective(noise_channels=['t1_charge_impedance', 't1_flux_bias_line'], common_noise_options=dict(T=0.050))
- Parameters
noise_channels (
Union[str,List[str],List[Tuple[str,Dict]],None]) – channels to be plotted, if None then noise channels given by supported_noise_channels are usedcommon_noise_options (
Optional[Dict]) – common options used when calculating coherence timesesys (
Optional[Tuple[ndarray,ndarray]]) – spectral data used during noise calculationsget_rate (
bool) – get rate or time
- Return type
float- Returns
- decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate
in inverse units.
- t1_flux_bias_line(i=1, j=0, M=400, Z=50, T=0.015, total=True, esys=None, get_rate=False)#
Noise due to a bias flux line.
References: Koch et al (2007), Groszkowski et al (2018)
- Parameters
i (int >=0) – state index that along with j defines a transition (i->j)
j (int >=0) – state index that along with i defines a transition (i->j)
M (
float) – Inductance in units of Phi_0 / AmpereZ (
Union[complex,float,Callable]) – A complex impedance; a fixed value or function of omegaT (
float) – temperature in Kelvintotal (
bool) – if False return a time/rate associated with a transition from state i to state j. if True return a time/rate associated with both i to j and j to i transitionsesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- t1_inductive(i=1, j=0, Q_ind=None, T=0.015, total=True, esys=None, get_rate=False)#
\(T_1\) due to inductive dissipation in a superinductor.
References: Nguyen et al (2019), Smith et al (2020)
- Parameters
i (int >=0) – state index that along with j defines a transition (i->j)
j (int >=0) – state index that along with i defines a transition (i->j)
Q_ind (
Union[float,Callable,None]) – inductive quality factor; a fixed value or function of omegaT (
float) – temperature in Kelvintotal (
bool) – if False return a time/rate associated with a transition from state i to state j. if True return a time/rate associated with both i to j and j to i transitionsesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- t1_quasiparticle_tunneling(i=1, j=0, Y_qp=None, x_qp=3e-06, T=0.015, Delta=0.00034, total=True, esys=None, get_rate=False)#
Noise due to quasiparticle tunneling across a Josephson junction.
References: Smith et al (2020), Catelani et al (2011), Pop et al (2014).
- Parameters
i (int >=0) – state index that along with j defines a transition (i->j)
j (int >=0) – state index that along with i defines a transition (i->j)
Y_qp (
Union[float,Callable,None]) – complex admittance; a fixed value or function of omegax_qp (
float) – quasiparticle density (in units of eV)T (
float) – temperature in KelvinDelta (
float) – superconducting gap (in units of eV)total (
bool) – if False return a time/rate associated with a transition from state i to state j. if True return a time/rate associated with both i to j and j to i transitionsesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- t2_effective(noise_channels=None, common_noise_options=None, esys=None, get_rate=False)#
Calculate the effective \(T_2\) time (or rate).
The effective \(T_2\) is calculated by considering a variety of pure dephasing and depolarizing noise channels, according to the formula:
\[\frac{1}{T_{2}^{\rm eff}} = \sum_k \frac{1}{T_{\phi}^{k}} + \frac{1}{2} \sum_j \frac{1}{T_{1}^{j}},\]where \(k\) (\(j\)) run over the relevant pure dephasing ( depolarization) channels that can contribute to the effective noise. By default all the noise channels given by the method effective_noise_channels are included. Users can also provide specific noise channels, with selected options, to be included in the effective \(T_2\) calculation. For example, assuming qubit is a qubit object, can can execute:
qubit.t2_effective(noise_channels=['t1_flux_bias_line', 't1_capacitive', ('tphi_1_over_f_flux', dict(A_noise=3e-6))], common_noise_options=dict(T=0.050))
- Parameters
noise_channels (None or str or list(str) or list(tuple(str, dict))) – channels to be plotted, if None then noise channels given by supported_noise_channels are used
common_noise_options (dict) – common options used when calculating coherence times
esys (tuple(evals, evecs)) – spectral data used during noise calculations
get_rate (bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- tphi_1_over_f(A_noise, i, j, noise_op, esys=None, get_rate=False, **kwargs)#
Calculate the 1/f dephasing time (or rate) due to arbitrary noise source.
We assume that the qubit energies (or the passed in eigenspectrum) has units of frequency (and not angular frequency).
- Parameters
A_noise (
float) – noise strengthi (int >=0) – state index that along with j defines a qubit
j (int >=0) – state index that along with i defines a qubit
noise_op (
Union[ndarray,csc_matrix]) – noise operator, typically Hamiltonian derivative w.r.t. noisy parameteresys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- tphi_1_over_f_cc(A_noise=1e-07, i=0, j=1, esys=None, get_rate=False, **kwargs)#
Calculate the 1/f dephasing time (or rate) due to critical-current noise from all three Josephson junctions \(EJ1\), \(EJ2\) and \(EJ3\). The combined noise is calculated by summing the rates from the individual contributions.
- Parameters
A_noise (
float) – noise strengthi (
int) – state index that along with j defines a qubitj (
int) – state index that along with i defines a qubitesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Return type
float- Returns
decoherence time in units of \(2\pi\) (system units), or rate in inverse units.
- tphi_1_over_f_cc1(A_noise=1e-07, i=0, j=1, esys=None, get_rate=False, **kwargs)#
Calculate the 1/f dephasing time (or rate) due to critical current noise of junction associated with Josephson energy \(EJ1\).
- Parameters
A_noise (
float) – noise strengthi (
int) – state index that along with j defines a qubitj (
int) – state index that along with i defines a qubitesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Return type
float- Returns
- decoherence time in units of \(2\pi ({\rm system\,\,units})\),
or rate in inverse units.
- tphi_1_over_f_cc2(A_noise=1e-07, i=0, j=1, esys=None, get_rate=False, **kwargs)#
Calculate the 1/f dephasing time (or rate) due to critical current noise of junction associated with Josephson energy \(EJ2\).
- Parameters
A_noise (
float) – noise strengthi (
int) – state index that along with j defines a qubitj (
int) – state index that along with i defines a qubitesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Return type
float- Returns
\(T_{\phi}\) time or rate: decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- tphi_1_over_f_cc3(A_noise=1e-07, i=0, j=1, esys=None, get_rate=False, **kwargs)#
Calculate the 1/f dephasing time (or rate) due to critical current noise of junction associated with Josephson energy \(EJ3\).
- Parameters
A_noise (
float) – noise strengthi (
int) – state index that along with j defines a qubitj (
int) – state index that along with i defines a qubitesys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Return type
float- Returns
decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- tphi_1_over_f_flux(A_noise=1e-06, i=0, j=1, esys=None, get_rate=False, **kwargs)#
Calculate the 1/f dephasing time (or rate) due to flux noise.
- Parameters
A_noise (
float) – noise strengthi (int >=0) – state index that along with j defines a qubit
j (int >=0) – state index that along with i defines a qubit
esys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- tphi_1_over_f_ng(A_noise=0.0001, i=0, j=1, esys=None, get_rate=False, **kwargs)#
Calculate the 1/f dephasing time (or rate) due to charge noise.
- Parameters
A_noise (
float) – noise strengthi (int >=0) – state index that along with j defines a qubit
j (int >=0) – state index that along with i defines a qubit
esys (
Optional[Tuple[ndarray,ndarray]]) – evals, evecs tupleget_rate (
bool) – get rate or time
- Returns
time or rate – decoherence time in units of \(2\pi ({\rm system\,\,units})\), or rate in inverse units.
- Return type
float
- wavefunction(esys=None, which=0, phi_grid=None)[source]#
Return a flux qubit wave function in phi1, phi2 basis
- Parameters
esys (
Optional[Tuple[ndarray,ndarray]]) – eigenvalues, eigenvectorswhich (
int) – index of desired wave function (default value = 0)phi_grid (
Optional[Grid1d]) – used for setting a custom grid for phi; if None use self._default_grid
- Return type
- widget(params=None)#
Use ipywidgets to modify parameters of class instance
- Parameters
params (Optional[Dict[str, Any]]) –