TunableTransmon¶
- class scqubits.core.transmon.TunableTransmon(EJmax, EC, d, flux, ng, ncut, truncated_dim=6, id_str=None, evals_method=None, evals_method_options=None, esys_method=None, esys_method_options=None)[source]¶
Class for the flux-tunable transmon qubit. The Hamiltonian is represented in dense form in the number basis, \(H_\text{CPB}=4E_\text{C}(\hat{ n}-n_g)^2-\frac{\mathcal{E}_\text{J}(\Phi)}{2}(|n\rangle\langle n+1|+\text{ h.c.})\), Here, the effective Josephson energy is flux-tunable: \(\mathcal{ E}_J(\Phi) = E_{J,\text{max}} \sqrt{\cos^2(\pi\Phi/\Phi_0) + d^2 \sin^2( \pi\Phi/\Phi_0)}\) and \(d=(E_{J2}-E_{J1})(E_{J1}+E_{J2})\) parametrizes the junction asymmetry.
Initialize with, for example:
TunableTransmon(EJmax=1.0, d=0.1, EC=2.0, flux=0.3, ng=0.2, ncut=30)
- Parameters:
EJmax (
float
) – maximum effective Josephson energy (sum of the Josephson energies of the two junctions)d (
float
) – junction asymmetry parameterEC (
float
) – charging energyflux (
float
) – flux threading the SQUID loop, in units of the flux quantumng (
float
) – offset chargencut (
int
) – charge basis cutoff, n = -ncut, …, ncuttruncated_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.esys_method (
Optional
[str
]) – method for esys diagonalization, callable or string representationesys_method_options (
Optional
[dict
]) – dictionary with esys diagonalization optionsevals_method (
Optional
[str
]) – method for evals diagonalization, callable or string representationevals_method_options (
Optional
[dict
]) – dictionary with evals diagonalization options
Methods
Returns the qubit's fundamental energy splitting, E_1 - E_0.
TunableTransmon.__init__
(EJmax, EC, d, flux, ...)Returns the qubit's anharmonicity, (E_2 - E_1) - (E_1 - E_0).
TunableTransmon.broadcast
(event, **kwargs)Request a broadcast from CENTRAL_DISPATCH reporting event.
TunableTransmon.cos_phi_operator
([energy_esys])Returns operator \(\cos \varphi\) in the charge or eigenenergy basis.
Use ipywidgets to create a new class instance
TunableTransmon.create_from_file
(filename)Read initdata and spectral data from file, and use those to create a new SpectrumData object.
TunableTransmon.d_hamiltonian_d_EJ
([energy_esys])Returns operator representing a derivative of the Hamiltonian with respect to EJ in the charge or eigenenergy basis.
Returns operator representing a derivative of the Hamiltonian with respect to flux in the charge or eigenenergy basis.
TunableTransmon.d_hamiltonian_d_ng
([energy_esys])Returns operator representing a derivative of the Hamiltonian with respect to charge offset ng in the charge or eigenenergy basis.
Return dictionary with default parameter values for initialization of class instance
TunableTransmon.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 default list of channels used when calculating effective t1 and t2 noise.
TunableTransmon.eigensys
([evals_count, ...])Calculates eigenvalues and corresponding eigenvectors using scipy.linalg.eigh.
TunableTransmon.eigenvals
([evals_count, ...])Calculates eigenvalues using scipy.linalg.eigh, returns numpy array of eigenvalues.
TunableTransmon.exp_i_phi_operator
([energy_esys])Returns operator \(e^{i\varphi}\) in the charge or eigenenergy basis.
TunableTransmon.filewrite
(filename)Convenience method bound to the class.
TunableTransmon.find_EJ_EC
(E01, anharmonicity)Finds the EJ and EC values given a qubit splitting E01 and anharmonicity.
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.
Calculates eigenvalues/eigenstates for a varying system parameter, given an array of parameter values.
TunableTransmon.hamiltonian
([energy_esys])Returns Hamiltonian in the charge or eigenenergy basis.
Returns Hilbert space dimension
TunableTransmon.matrixelement_table
(operator)Returns table of matrix elements for operator with respect to the eigenstates of the qubit.
TunableTransmon.n_operator
([energy_esys])Returns charge operator n in the charge or eigenenergy basis.
Return the transmon wave function in number basis.
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.
Generates a simple plot of a set of eigenvalues as a function of one parameter.
Generates a simple plot of a set of eigenvalues as a function of one parameter.
TunableTransmon.plot_matrixelements
(operator)Plots matrix elements for operator, given as a string referring to a class method that returns an operator matrix.
TunableTransmon.plot_n_wavefunction
([esys, ...])Plots transmon wave function in charge basis
Alias for plot_wavefunction
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.
TunableTransmon.plot_wavefunction
([which, ...])Plot 1d phase-basis wave function(s).
Transmon phase-basis potential evaluated at phi.
TunableTransmon.process_hamiltonian
(...[, ...])Return qubit Hamiltonian in chosen basis: either return unchanged (i.e., in native basis) or transform into eigenenergy basis
TunableTransmon.process_op
(native_op[, ...])Processes the operator native_op: either hand back native_op unchanged, or transform it into the energy eigenbasis.
TunableTransmon.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.
TunableTransmon.set_and_return
(attr_name, value)Allows to set an attribute after which self is returned. This is useful for doing something like example::.
TunableTransmon.set_params
(**kwargs)Set new parameters through the provided dictionary.
Set new parameters through the provided dictionary.
TunableTransmon.sin_phi_operator
([energy_esys])Returns operator \(\sin \varphi\) in the charge or eigenenergy basis.
Return a list of supported noise channels
TunableTransmon.t1
(i, j, noise_op, ...[, T, ...])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.
TunableTransmon.t1_capacitive
([i, j, Q_cap, ...])\(T_1\) due to dielectric dissipation in the Josephson junction capacitances.
TunableTransmon.t1_charge_impedance
([i, j, ...])Noise due to charge coupling to an impedance (such as a transmission line).
TunableTransmon.t1_effective
([...])Calculate the effective \(T_1\) time (or rate).
TunableTransmon.t1_flux_bias_line
([i, j, M, ...])Noise due to a bias flux line.
TunableTransmon.t1_inductive
([i, j, Q_ind, ...])\(T_1\) due to inductive dissipation in a superinductor.
Noise due to quasiparticle tunneling across a Josephson junction.
TunableTransmon.t2_effective
([...])Calculate the effective \(T_2\) time (or rate).
TunableTransmon.tphi_1_over_f
(A_noise, i, j, ...)Calculate the 1/f dephasing time (or rate) due to arbitrary noise source.
TunableTransmon.tphi_1_over_f_cc
([A_noise, ...])Calculate the 1/f dephasing time (or rate) due to critical current noise.
Calculate the 1/f dephasing time (or rate) due to flux noise.
TunableTransmon.tphi_1_over_f_ng
([A_noise, ...])Calculate the 1/f dephasing time (or rate) due to charge noise.
TunableTransmon.wavefunction
([esys, which, ...])Return the transmon wave function in phase basis.
Plot defaults for plotting.wavefunction1d.
TunableTransmon.widget
([params])Use ipywidgets to modify parameters of class instance
Attributes
EC
Descriptor class for properties that are to be monitored for changes.
This is the effective, flux dependent Josephson energy, playing the role of EJ in the parent class Transmon
EJmax
Descriptor class for properties that are to be monitored for changes.
d
Descriptor class for properties that are to be monitored for changes.
flux
Descriptor class for properties that are to be monitored for changes.
id_str
ncut
Descriptor class for properties that are to be monitored for changes.
ng
Descriptor class for properties that are to be monitored for changes.
truncated_dim
Descriptor class for properties that are to be monitored for changes.
- E01()¶
Returns the qubit’s fundamental energy splitting, E_1 - E_0.
- Return type:
float
- property EJ: float¶
This is the effective, flux dependent Josephson energy, playing the role of EJ in the parent class Transmon
- anharmonicity()¶
Returns the qubit’s anharmonicity, (E_2 - E_1) - (E_1 - E_0).
- Return type:
float
- broadcast(event, **kwargs)¶
Request a broadcast from CENTRAL_DISPATCH reporting event.
- Parameters:
event (
str
) – event name from EVENTS**kwargs
- Return type:
None
- cos_phi_operator(energy_esys=False)¶
Returns operator \(\cos \varphi\) in the charge or eigenenergy basis.
- Parameters:
energy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns operator \(\cos \varphi\) in the charge basis. If True, the energy eigenspectrum is computed, returns operator \(\cos \varphi\) in the energy eigenbasis. If energy_esys = esys, where esys is a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns operator \(\cos \varphi\) in the energy eigenbasis, and does not have to recalculate eigenspectrum.- Return type:
ndarray
- Returns:
Operator \(\cos \varphi\) in chosen basis as ndarray. If the eigenenergy basis is chosen, unless energy_esys is specified, \(\cos \varphi\) has dimensions of truncated_dim x truncated_dim. Otherwise, if eigenenergy basis is chosen, \(\cos \varphi\) has dimensions of m x m, for m given eigenvectors.
- 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_EJ(energy_esys=False)¶
Returns operator representing a derivative of the Hamiltonian with respect to EJ in the charge or eigenenergy basis.
- Parameters:
energy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns operator in the charge basis. If True, the energy eigenspectrum is computed, returns operator in the energy eigenbasis. If energy_esys = esys, where esys is a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns operator in the energy eigenbasis, and does not have to recalculate eigenspectrum.- Return type:
ndarray
- Returns:
Operator in chosen basis as ndarray. If the eigenenergy basis is chosen, unless energy_esys is specified, operator has dimensions of truncated_dim x truncated_dim. Otherwise, if eigenenergy basis is chosen, operator has dimensions of m x m, for m given eigenvectors.
- d_hamiltonian_d_flux(energy_esys=False)[source]¶
Returns operator representing a derivative of the Hamiltonian with respect to flux in the charge or eigenenergy basis.
Here, the derivative is taken with respect to flux before the qubit’s \(\phi\) degree of freedom in the Hamiltonian is shifted by a flux-dependent quantity \(\varphi_{0}\) (see Eq. 2.17 and surrounding text in PRA 76, 042319 (2007)). Then only after the flux derivative is taken, both the Hamiltonian as well as its flux derivative are assumed to be shifted by \(\varphi_{0}\).
- Parameters:
energy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns operator in the charge basis. If True, the energy eigenspectrum is computed, returns operator in the energy eigenbasis. If energy_esys = esys, where esys is a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns operator in the energy eigenbasis, and does not have to recalculate eigenspectrum.- Return type:
ndarray
- Returns:
Operator in chosen basis as ndarray. If the eigenenergy basis is chosen, unless energy_esys is specified, operator has dimensions of truncated_dim x truncated_dim. Otherwise, if eigenenergy basis is chosen, operator has dimensions of m x m, for m given eigenvectors.
- d_hamiltonian_d_ng(energy_esys=False)¶
Returns operator representing a derivative of the Hamiltonian with respect to charge offset ng in the charge or eigenenergy basis.
- Parameters:
energy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns operator in the charge basis. If True, the energy eigenspectrum is computed, returns operator in the energy eigenbasis. If energy_esys = esys, where esys is a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns operator in the energy eigenbasis, and does not have to recalculate eigenspectrum.- Return type:
ndarray
- Returns:
Operator in chosen basis as ndarray. If the eigenenergy basis is chosen, unless energy_esys is specified, operator has dimensions of truncated_dim x truncated_dim. Otherwise, if eigenenergy basis is chosen, operator has dimensions of m x m, for m given eigenvectors.
- 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:
TypeVar
(SerializableType
, bound= Serializable)- Parameters:
io_data (IOData)
- classmethod effective_noise_channels()¶
Return a default list of channels used when calculating effective t1 and t2 noise.
- 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 (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_operator(energy_esys=False)¶
Returns operator \(e^{i\varphi}\) in the charge or eigenenergy basis.
- Parameters:
energy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns operator \(e^{i\varphi}\) in the charge basis. If True, the energy eigenspectrum is computed, returns operator \(e^{i\varphi}\) in the energy eigenbasis. If energy_esys = esys, where esys is a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns operator \(e^{i\varphi}\) in the energy eigenbasis, and does not have to recalculate eigenspectrum.- Return type:
ndarray
- Returns:
Operator \(e^{i\varphi}\) in chosen basis as ndarray. If the eigenenergy basis is chosen, unless energy_esys is specified, \(e^{i\varphi}\) has dimensions of truncated_dim x truncated_dim. Otherwise, if eigenenergy basis is chosen, \(e^{i\varphi}\) has dimensions of m x m, for m given eigenvectors.
- filewrite(filename)¶
Convenience method bound to the class. Simply accesses the write function.
- Return type:
None
- Parameters:
filename (str)
- static find_EJ_EC(E01, anharmonicity, ng=0, ncut=30)¶
Finds the EJ and EC values given a qubit splitting E01 and anharmonicity.
- Parameters:
E01 (
float
) – qubit transition energyanharmonicity (
float
) – absolute qubit anharmonicity, (E2-E1) - (E1-E0)ng – offset charge (default: 0)
ncut – charge number cutoff (default: 30)
- Return type:
Tuple
[float
,float
]- Returns:
A tuple of the EJ and EC values representing the best fit.
- 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 (
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(energy_esys=False)¶
Returns Hamiltonian in the charge or eigenenergy basis.
- Parameters:
energy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns Hamiltonian in the charge basis. If True, the energy eigenspectrum is computed; returns Hamiltonian in the energy eigenbasis. If energy_esys = esys, where esys is a tuple containing two ndarrays (eigenvalues and energy eigenvectors); then return the Hamiltonian in the energy eigenbasis, do not recalculate eigenspectrum.- Return type:
ndarray
- Returns:
Hamiltonian in chosen basis as ndarray. For energy_esys=False, the Hamiltonian has dimensions of truncated_dim x truncated_dim. For energy_sys=esys, the Hamiltonian has dimensions of m x m, for m given eigenvectors.
- hilbertdim()¶
Returns Hilbert space dimension
- Return type:
int
- 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 (
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 (
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
]
- n_operator(energy_esys=False)¶
Returns charge operator n in the charge or eigenenergy basis.
- Parameters:
energy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns charge operator n in the charge basis. If True, energy eigenspectrum is computed, returns charge operator n in the energy eigenbasis. If energy_esys = esys, where esys is a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns charge operator n in the energy eigenbasis, and does not have to recalculate the eigenspectrum.- Return type:
ndarray
- Returns:
Charge operator n in chosen basis as ndarray. For energy_esys=True, n has dimensions of truncated_dim x truncated_dim. If an actual eigensystem is handed to energy_sys, then n has dimensions of m x m, where m is the number of given eigenvectors.
- numberbasis_wavefunction(esys=None, which=0)¶
Return the transmon wave function in number basis. The specific index of the wave function to be returned is which.
- Parameters:
esys (
Tuple
[ndarray
,ndarray
]) – if None, the eigensystem is calculated on the fly; otherwise, the provided eigenvalue, eigenvector arrays as obtained from .eigensystem(), are used (default value = None)which (
int
) – eigenfunction index (default value = 0)
- Return type:
- 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
]]]) – channels to be plotted, if None then noise channels given by supported_noise_channels are usedcommon_noise_options (
Dict
) – common options used when calculating coherence timesspectrum_data (
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 (
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_n_wavefunction(esys=None, mode='real', which=0, nrange=None, **kwargs)¶
Plots transmon wave function in charge basis
- Parameters:
esys (
Tuple
[ndarray
,ndarray
]) – eigenvalues, eigenvectorsmode (
str
) – ‘abs_sqr’, ‘abs’, ‘real’, ‘imag’which (
int
) – index or indices of wave functions to plot (default value = 0)nrange (
Tuple
[int
,int
]) – range of n to be included on the x-axis (default value = (-5,6))**kwargs – plotting parameters
- Return type:
Tuple
[Figure
,Axes
]
- plot_phi_wavefunction(esys=None, which=0, phi_grid=None, mode='abs_sqr', scaling=None, **kwargs)¶
Alias for plot_wavefunction
- Return type:
Tuple
[Figure
,Axes
]- Parameters:
esys (Tuple[ndarray, ndarray])
which (int)
phi_grid (Grid1d)
mode (str)
scaling (float)
- 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
]]]) – channels to be plotted, if None then noise channels given by supported_noise_channels are usedcommon_noise_options (
Dict
) – common options used when calculating coherence timesspectrum_data (
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
]]]) – channels to be plotted, if None then noise channels given by supported_noise_channels are usedcommon_noise_options (
Dict
) – common options used when calculating coherence timesspectrum_data (
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(which=0, mode='real', esys=None, phi_grid=None, scaling=None, **kwargs)¶
Plot 1d phase-basis wave function(s). Must be overwritten by higher-dimensional qubits like FluxQubits and ZeroPi.
- Parameters:
which (
Union
[int
,Iterable
[int
]]) – single index or tuple/list of integers indexing the wave function(s) to be plotted. If which is -1, all wavefunctions up to the truncation limit are plotted.mode (
str
) – choices as specified in constants.MODE_FUNC_DICT (default value = ‘abs_sqr’)esys (
Tuple
[ndarray
,ndarray
]) – eigenvalues, eigenvectorsphi_grid (
Grid1d
) – used for setting a custom grid for phi; if None use self._default_gridscaling (
Optional
[float
]) – custom scaling of wave function amplitude/modulus**kwargs – standard plotting option (see separate documentation)
- Return type:
Tuple
[Figure
,Axes
]
- potential(phi)¶
Transmon phase-basis potential evaluated at phi.
- Parameters:
phi (
Union
[float
,ndarray
]) – phase variable value- Return type:
ndarray
- process_hamiltonian(native_hamiltonian, energy_esys=False)¶
Return qubit Hamiltonian in chosen basis: either return unchanged (i.e., in native basis) or transform into eigenenergy basis
- Parameters:
native_hamiltonian (
Union
[ndarray
,csc_matrix
]) – Hamiltonian in native basisenergy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns Hamiltonian in the native basis If True, the energy eigenspectrum is computed, returns Hamiltonian in the energy eigenbasis if energy_esys is the energy eigenspectrum, in the form of a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns Hamiltonian in the energy eigenbasis, and does not have to recalculate eigenspectrum.
- Return type:
Union
[ndarray
,csc_matrix
]- Returns:
Hamiltonian, either unchanged in native basis, or transformed into eigenenergy basis
- process_op(native_op, energy_esys=False)¶
Processes the operator native_op: either hand back native_op unchanged, or transform it into the energy eigenbasis. (Native basis refers to the basis used internally by each qubit, e.g., charge basis in the case of Transmon.
- Parameters:
native_op (
Union
[ndarray
,csc_matrix
]) – operator in native basisenergy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns operator in the native basis If True, the energy eigenspectrum is computed, returns operator in the energy eigenbasis if energy_esys is the energy eigenspectrum, in the form of a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns operator in the energy eigenbasis, and does not have to recalculate eigenspectrum.
- Return type:
Union
[ndarray
,csc_matrix
]- Returns:
native_op either unchanged or transformed into eigenenergy basis
- 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.
- set_params_from_gui(change)¶
Set new parameters through the provided dictionary.
- sin_phi_operator(energy_esys=False)¶
Returns operator \(\sin \varphi\) in the charge or eigenenergy basis.
- Parameters:
energy_esys (
Union
[bool
,Tuple
[ndarray
,ndarray
]]) – If False (default), returns operator \(\sin \varphi\) in the charge basis. If True, the energy eigenspectrum is computed, returns operator \(\sin \varphi\) in the energy eigenbasis. If energy_esys = esys, where esys is a tuple containing two ndarrays (eigenvalues and energy eigenvectors), returns operator \(\sin \varphi\) in the energy eigenbasis, and does not have to recalculate eigenspectrum.- Return type:
ndarray
- Returns:
Operator \(\sin \varphi\) in chosen basis as ndarray. If the eigenenergy basis is chosen, unless energy_esys is specified, \(\sin \varphi\) has dimensions of truncated_dim x truncated_dim. Otherwise, if eigenenergy basis is chosen, \(\sin \varphi\) has dimensions of m x m, for m given eigenvectors.
- classmethod supported_noise_channels()[source]¶
Return a list of supported noise channels
- Return type:
List
[str
]
- t1(i, j, noise_op, spectral_density, T=0.015, 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 operatorT (
float
) – Temperature defined in Kelvinspectral_density (
Callable
) – defines a spectral density, must take two arguments: omega and T (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 (
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, noise_op=None, branch_params=None)¶
\(T_1\) due to dielectric dissipation in the Josephson junction capacitances.
References: Smith et al (2020), see also Nguyen et al (2019).
- 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
]) – 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 (
Tuple
[ndarray
,ndarray
]) – evals, evecs tupleget_rate (
bool
) – get rate or timenoise_op (ndarray | csc_matrix | Qobj | None)
branch_params (dict | None)
- 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_op=None)¶
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 (
Tuple
[ndarray
,ndarray
]) – evals, evecs tupleget_rate (
bool
) – get rate or timenoise_op (ndarray | csc_matrix | Qobj | None)
- 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
]]]) – channels to be plotted, if None then noise channels given by supported_noise_channels are usedcommon_noise_options (
Dict
) – common options used when calculating coherence timesesys (
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_op_method=None)¶
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 (
Tuple
[ndarray
,ndarray
]) – evals, evecs tupleget_rate (
bool
) – get rate or timenoise_op_method (Callable | None)
- 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, noise_op=None, branch_params=None)¶
\(T_1\) due to inductive dissipation in a superinductor.
References: Smith et al (2020), see also Nguyen et al (2019).
- 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
]) – 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 (
Tuple
[ndarray
,ndarray
]) – evals, evecs tupleget_rate (
bool
) – get rate or timenoise_op (ndarray | csc_matrix | Qobj | None)
branch_params (dict | None)
- 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_op=None)¶
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
]) – 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 (
Tuple
[ndarray
,ndarray
]) – evals, evecs tupleget_rate (
bool
) – get rate or timenoise_op (ndarray | csc_matrix | Qobj | None)
- 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 (
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.
- 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 (
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_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 (
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 (
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)¶
Return the transmon wave function in phase basis. The specific index of the wavefunction is which. esys can be provided, but if set to None then it is calculated on the fly.
- Parameters:
esys (
Optional
[Tuple
[ndarray
,ndarray
]]) – if None, the eigensystem is calculated on the fly; otherwise, the provided eigenvalue, eigenvector arrays as obtained from .eigensystem() are usedwhich (
int
) – eigenfunction index (default value = 0)phi_grid (
Grid1d
) – used for setting a custom grid for phi; if None use self._default_grid
- Return type:
- wavefunction1d_defaults(mode, evals, wavefunc_count)¶
Plot defaults for plotting.wavefunction1d.
- Parameters:
mode (
str
) – amplitude modifier, needed to give the correct default y labelevals (
ndarray
) – eigenvalues to include in plotwavefunc_count (
int
) – number of wave functions to be plotted
- Return type:
Dict
[str
,Any
]
- widget(params=None)¶
Use ipywidgets to modify parameters of class instance
- Parameters:
params (Dict[str, Any] | None)