2. Package API

2.2. Submodules

2.3. qmix.circuit module

This module contains classes and functions to describe the embedding circuit.

Description:

In experimental systems, SIS junctions are embedded within complex RF networks. These networks are referred to as the embedding circuit. Since all of the components in these embedding circuits are linear, the embedding circuit can be reduced to a Thevenin equivalent circuit, with one for each tone and harmonic.

To fully describe the embedding circuit, 3 bits of information are needed for each signal that is applied to the junction:

  1. the frequency of the applied signal,

  2. the Thevenin voltage of the embedding circuit at this freq., and

  3. the Thevenin impedance of the embedding circuit at this freq.

The main class in this module (EmbeddingCircuit) allows the user to build an embedding circuit in the proper format.

class qmix.circuit.EmbeddingCircuit(num_f=1, num_p=1, vb_min=0, vb_max=2, vb_npts=201, fgap=None, vgap=None, rn=None, name='')

Bases: object

Class for building and describing the embedding circuit.

This includes the frequencies, Thevenin voltages and Thevenin impedances of all signals applied to the junction.

Note

Unless specified otherwise, all input values are normalized. The voltages are normalized to the gap voltage, resistances are normalized to the normal-state resistance, currents are normalized to gap current, and frequencies are normalized to the gap frequency. Refer to the argument descriptions below to see if the value should be normalized or not.

Creating an instance of this class will set the sizes and data types of all of the class attributes, but the actual values will need to be set manually. In this way, this class is sort of like a fancy struct. The class attributes that have to be set manually are:

  • freq: frequencies (normalized to the gap frequency)

  • vt: Thevenin voltage (normalized to the gap voltage)

  • zt: Thevenin impedance (normalized to the normal resistance)

See the attribute descriptions below for information about how to do this.

Example

Here we will create an instance of the embedding circuit class with 2 tones and 3 harmonics:

>>> cct1 = EmbeddingCircuit(2, 3, name='Test 1')
>>> print(cct1)
Embedding circuit (Tones:2, Harmonics:3): Test 1

Once initialized, we can begin defining the properties of the embedding circuit. I normally start with the frequencies:

>>> cct1.freq[1] = 0.30  # first tone
>>> cct1.freq[2] = 0.32  # second tone

Then, we have to set the voltages and impedances for all of the different signals. For example, for the 1st harmonic of the 1st tone:

>>> cct1.vt[1,1] = 0.5             # Thevenin voltage
>>> cct1.zt[1,1] = 0.3 + 1j * 0.1  # Thevenin impedance

This has to be done for each signal (a total of 6 times in this case).

In order to use non-normalized values (e.g., set the available power of a signal in units [W]), we need to define the electrical properties of the junction during initialization. For example:

>>> cct2 = EmbeddingCircuit(1, 1, vgap=2.8e-3, rn=14.)

You can now set the frequencies. E.g.:

>>> cct2.set_freq(250, f=1, units='GHz')
>>> round(cct2.freq[1], 4)
0.3693

Once an impedance has been set, you can also set the power of the signal using absolute units. Here we will set the available power for the first harmonic of the first tone to 10 nW (10e-9 W).

>>> cct2.zt[1,1] = 0.5
>>> cct2.set_available_power(10, 1, 1, 'nW')

And we can then display this power in units [dBm].

>>> cct2.available_power(1, 1, 'dBm')
-50.0
Parameters
  • num_f (int, optional, default is 1) – Number of fundamental frequencies (tones) applied to the junction.

  • num_p (int, optional, default is 1) – Number of harmonics included for each tone.

  • vb_min (float, optional, default is 0) – Minimum bias voltage.

  • vb_max (float, optional, default is 2) – Maximum bias voltage.

  • vb_npts (int, optional, default is 201) – Number of points in bias voltage sweep.

  • fgap (float, optional) – Gap frequency of the junction in units [Hz]. This is equal to e*Vgap/h, where e is the charge of an electron, Vgap is the gap voltage, and h is the Planck constant. fgap is used to normalize and de-normalize frequency values (that’t it!).

  • vgap (float, optional) – Gap voltage of the junction in units [V]. This is the voltage where the sharp non-linearity in the DC I-V curve occurs (a.k.a., the transition voltage).

  • rn (float, optional) – Normal-state resistance of the junction in units [ohms]. This is the resistance of the junction at a temperature slight above the critical temperature. It is found by calculating the dynamic resistance of the DC I-V curve above the gap voltage.

  • name (str, optional) – Name used to describe this specific instance.

freq

An array for the frequencies normalized to the gap frequency. This is a 1-dimensional array of real numbers. It contains the frequencies of all of the fundamental tones that are applied to the junction. The normalized frequency is equivalent to the normalized photon voltage. The photon voltage is defined as hf/e, where h is the Planck constant, f is the frequency of the fundamental tone, and e is the charge of an electron. Note that this array is 1-based, meaning that the frequency of the 1st tone is found in .freq[1]. (The index represents the tone number.) This attribute must be set after initialization!

Type

numpy.ndarray

vt

An array for the Thevenin voltage normalized to the gap voltage. This is a 2-dimensional array of complex values. It contains the voltages for all of the Thevenin equivalent circuits (which describe the embedding circuit). In order, the indices are: .vt[f,p] for tone f and harmonic p. Note that this array is 1-based, meaning that the voltage for tone number 2 / harmonic number 3 is stored in vt[2,3]. This attribute must be set after initialization!

Type

numpy.ndarray

zt

An array for the Thevenin impedance array normalized to the normal-state resistance. This is a 2-dimensional array of complex values. It contains the impedances of all of the Thevenin equivalent circuits (which describe the embedding circuit). In order, the indices are: .zt[f,p], for tone f and harmonic p. Note that this array is 1-based, meaning that the impedance for tone number 2, and harmonic number 3 is stored in zt[2,3]. This attribute must be set after initialization!

Type

numpy.ndarray

num_f

Number of fundamental frequencies (tones) applied to the junction.

Type

int

num_p

Number of harmonics included for each tone.

Type

int

num_n

Total number of signals. This is equal to num_f*num_p.

Type

int

fgap

Gap frequency of the junction in units [Hz]. This is equal to e*Vgap/h, where e is the charge of an electron, Vgap is the gap voltage, and h is the Planck constant. Note that E=fgap*e is the energy required to break Cooper pairs, so at frequencies above the fgap the superconductors will begin to become lossy. Here, fgap is used to normalize and de-normalize frequency values (that’t it).

Type

float

vgap

Gap voltage of the junction in units [V]. This is the voltage where the sharp non-linearity in the DC I-V curve occurs (i.e., the transition voltage). This value is used to normalize and de-normalize voltages.

Type

float

igap

Gap current of the junction in units [A]. This is equal to vgap/rn. This value is used to normalize and de-normalize currents.

Type

float

rn

Normal-state resistance of the junction in units [ohms]. This is the resistance of the junction at a temperature slightly above the critical temperature. It is found by calculating the dynamic resistance of the DC I-V curve above the gap voltage. This value is used to normalize and de-normalize resistances.

Type

float

vb

Array for the DC bias voltage sweep. This value is normalized to the gap voltage.

Type

numpy.ndarray

vb_npts

Number of points in the bias voltage sweep.

Type

int

name

A name to describe this instance of the embedding circuit class.

Type

str

comment

A list of comments to describe the different signals. For example, to describe tone 1/harmonic 1 as the local-oscillator signal, you might use cct.comment[1][1] = "LO". This has to be set after the initialization of the EmbeddingCircuit class.

Type

list

available_power(f=1, p=1, units='W')

Return available power of tone f and harmonic p.

Note

Gap voltage and normal resistance must be set prior to using this method. If they are not, an error will be raised.

Parameters
  • f (int, optional, default is 1) – Tone index number.

  • p (int, optional, default is 1) – Harmonic index number.

  • units (str, optional, default is 'W') – Units for power. One of ‘W’, ‘mW’, ‘uW’, ‘nW’, ‘pW’, ‘fW’, ‘dBm’, or ‘dBW’.

Returns

Available power in specified units

Return type

float

initialize_vj()

Initialize junction voltage array.

Returns an empty matrix that is the shape that vj should be (the voltage across the junction). Strictly speaking, vj shouldn’t be saved within this class, but it is okay for this class to initialize vj since it has all the data about what the matrix sizes should be.

This function is useful when you want to set the voltage across the junction directly (skipping the harmonic balance procedure).

Returns

An empty matrix for the junction voltage

Return type

ndarray

lock()

Make all Numpy arrays contained within this class unwriteable.

This can be useful for debugging. An error will be raised if you try to change the values of the Numpy arrays while they are locked.

print_info()

Print information about the embedding circuit to the terminal.

save_info(filename='embedding-circuit.txt')

Save this embedding circuit to a text file.

This text file can then be read in by read_circuit in order to regenerate the embedding circuit.

Parameters

filename (string) – Filename for embedding circuit file

set_alpha(alpha, f=1, p=1, zj=0.66)

Set the drive level of tone f and harmonic p (approximate).

This method guesses what the Thevenin voltage should be in order to get the desired drive level, but you won’t actually know what the drive level is until you run the simulation.

Note

Frequency and Thevenin impedance must be set prior to using this method. Otherwise, an assertion error will be raised.

Parameters
  • alpha (float) – drive level, alpha = voltage / vph

  • f (int, optional, default is 1) – tone

  • p (int, optional, default is 1) – harmonic

  • zj (float, optional, default is 0.66) – the impedance to assume for the junction (normalized to the normal-state resistance). This value will depend on frequency and pump level.

set_available_power(power, f=1, p=1, units='W')

Set available power of tone f and harmonic p.

This method will set the Thevenin voltage in order to provide the correct power level.

Note

The gap voltage, normal resistance and Thevenin impedance must be set prior to using this method. Otherwise, an assertion error will be raised.

Parameters
  • power (float) – power, in given units

  • f (int, optional, default is 1) – tone

  • p (int, optional, default is 1) – harmonic

  • units (str, optional, default is 'W') – units for power. One of ‘W’,

  • 'mW'

  • 'uW'

  • 'nW'

  • 'pW'

  • 'fW'

  • 'dBm'

  • 'dBW'. (or) –

set_freq(value, f=1, units='Hz')

Set the frequency of tone f.

Normally, this can be done by setting the value of the attribute directly. E.g.:

>>> cct = EmbeddingCircuit(1, 1, vgap=2.8e-3, rn=14.)
>>> cct.freq[1] = 0.5

However, if you would instead like to use non-normalized units, this method can be very handy. E.g.:

>>> cct.set_freq(350, f=1, units='GHz')
>>> round(cct.freq[1], 2)
0.52

Note

The gap frequency or the gap voltage must be defined in order to use this method.

Parameters
  • value (float) – value to set using given units

  • f (int, optional, default is 1) – tone number

  • units (str, optional, default is 'Hz') – units for input value, ‘Hz’ for frequency in units Hz, ‘V’ for photon voltage in units V, and ‘norm’ for either normalized photon voltage or normalized frequency. SI prefixes can also be included: ‘MHz’, ‘GHz’, ‘THz’, and ‘mV’.

set_name(name, f=1, p=1)

Set a name for a given tone and harmonic combination.

This has no effect on the simulation. It’s just nice for keeping track of the different signals.

Parameters
  • name (str) – name of tone/harmonic

  • f (int, optional, default is 1) – frequency number to set

  • p (int, optional, default is 1) – harmonic number to set

unlock()

Make all Numpy arrays contained within this class writeable.

This can be useful for debugging.

property vph

Get photon voltage.

For backwards compatibility.

qmix.circuit.read_circuit(filename)

Build an embedding circuit from an embedding circuit file.

This function will build an instance of the EmbeddingCircuit class based on a file previously generated by the EmbeddingCircuit.save_info method.

Parameters

filename (str) – filename of the embedding circuit file

Returns

instance of the embedding circuit class

Return type

qmix.circuit.EmbeddingCircuit

2.4. qmix.harmonic_balance module

This module contains functions to perform harmonic balance of non-linear SIS mixer circuits.

Description

Each signal that is applied to an SIS junction can be represented by a Thevenin equivalent circuit (see qmix.circuit). This circuit will then induce a voltage across the SIS junction. The exact voltage that is induced depends on the impedance of the SIS junction. However, this impedance changes depending on the other signals that are applied to the junction.

Harmonic balance is a procedure to solve for the voltage across the SIS junction for each signal that is applied to the junction. This techniques uses Newton’s method to find the solution numerically.

qmix.harmonic_balance.check_hb_error(vj_check, cct, resp, num_b=15, stop_rerror=0.001)

Check the results from the harmonic_balance function.

Just to double check. Mostly for debugging purposes.

Parameters
  • vj_check (ndarray) – The voltage across junction to check

  • cct (qmix.circuit.EmbeddingCircuit) – Embedding circuit

  • resp (qmix.respfn.RespFn) – Response function

  • num_b (optional) – Summation limits for phase factor coefficients, default is 15

  • stop_rerror (float, optional) – Maximum acceptable relative error, default is 0.001

Raises

AssertionError – If the stop error is not met

qmix.harmonic_balance.harmonic_balance(cct, resp, num_b=15, max_it=10, stop_rerror=0.001, vj_initial=None, damp_coeff=1.0, mode='o', verbose=True, zj_guess=0.67)

Perform harmonic balance.

Determine the harmonic balance of the junction + embedding circuit system. Uses Newton’s method to find the solution. For more information, see Garrett (2018); Kittara (2002); Kittara, Winthington & Yassin (2007); or Withington, Kittara & Yassin (2003). [Full references in online docs.]

Parameters
  • cct (qmix.circuit.EmbeddingCircuit) – Embedding circuit

  • resp (qmix.respfn.RespFn) – Response function

  • num_b (int/tuple, optional) – Summation limits for phase factor coefficients, default is 15

  • max_it (int, optional) – Maximum number of iterations, default is 10

  • stop_rerror (float, optional) – Maximum acceptable relative error, default is 0.001

  • vj_initial (ndarray, optional) – Initial guess of junction voltage (vj), default is None

  • damp_coeff (float, optional) – Dampening coefficient for correction factor (0-1), default is 1

  • mode (string, optional) – output vj (‘o’), print (‘p’), output extra data (‘x’), default is “o”

  • verbose (bool, optional) – print info to terminal if true, default is True

Returns

junction voltage that satisfies the circuit

Return type

ndarray

2.5. qmix.qtcurrent module

This module contains functions to calculate the quasiparticle tunneling currents passing through an SIS junction.

Description:

Given the voltages applied across an SIS junction, the quasiparticle tunneling currents can be calculated using multi-tone spectral domain analysis (MTSDA; see references in online docs).

Note

This code is largely based on P. Kittara’s 2002 DPhil thesis (see references in online docs). I include some inline comments to refer to specific equations.

Also, all of the values in this module are normalized, i.e., voltages are normalized to the gap voltage, frequencies are normalized to the gap frequency, etc.

qmix.qtcurrent.calculate_phase_factor_coeff(vj, freq, num_f, num_p, num_b)

Calculate the overall phase factor spectrum coefficients.

Runs once per qtcurrent function call.

Eqns. 5.7 and 5.12 in Kittara’s thesis.

Parameters
  • vj (ndarray) – Voltage across the SIS junction

  • freq (ndarray) – Frequencies

  • num_f (int) – Number of non-harmonically related frequencies

  • num_p (int) – Number of harmonics

  • num_b (int) – Summation limits for phase factor coefficients

Returns

Phase factor spectrum coefficients (C_k(H) in Kittara)

Return type

ndarray

qmix.qtcurrent.interpolate_respfn(cct, resp, num_b)

Interpolate the response function at all necessary voltages.

I have included this as a stand-alone function because if you are going to be running qtcurrent over and over again with the same input signal frequencies, it can save time by pre-interpolating the response function.

Parameters
Returns

The interpolated response function as a matrix.

Return type

ndarray

qmix.qtcurrent.qtcurrent(vj, cct, resp, freq_list, num_b=15, verbose=True, resp_matrix=None)

Calculate the quasiparticle tunneling current.

This function uses multi-tone spectral domain analysis (MTSDA; see references in online docs). The current is calculated based on the voltage applied across the junction (vj).

Note

This function will return the tunneling current for all of the frequencies listed in freq_list (normalized to the gap frequency). E.g., to solve for the dc tunneling current and the ac tunneling current at 230 GHz, the freq_list would be [0, 230e9 / fgap] where fgap is the gap frequency.

Maximum of 4 non-harmonically related tones.

Parameters
  • vj (ndarray) – Voltage across the SIS junction

  • cct (qmix.circuit.EmbeddingCircuit) – Embedding circuit

  • resp (qmix.respfn.RespFn) – Response function

  • freq_list – Calculate the tunneling currents for these frequencies (normalized to the gap frequency)

  • num_b (float/tuple, optional) – Summation limits for phase factor coefficients, default is 15

  • verbose (bool, optional) – Print info to the terminal if true, default is True

  • resp_matrix (ndarray, optional) – The interpolated response function matrix, generated by interpolate_respfn(), default is None

Returns

Quasiparticle tunneling current

Return type

ndarray

2.6. qmix.respfn module

This module contains classes to represent the response function of the SIS junction.

There are several different types of response function classes:

  • Response functions generated directly from I-V data:

    • RespFn: This is the base class for all of the other response function classes. This class will generate a response function based on a DC I-V curve (i.e., DC voltage and current data). Note that this class assumes that you have “pre-processed” the data. This means that it will use the voltage and current data to generate the interpolation directly. Normally, you want to have more data points around curvier regions in order to minimize how much time the interpolation takes. If you haven’t done this, it is a good idea to use RespFnFromIVData instead.

    • RespFnFromIVData: This class will generate a response function based on a DC I-V curve (i.e., DC voltage and current data). Unlike RespFn, this class will resample the response function in order to optimize the interpolation.

  • Response functions generated from I-V curve models:

    • RespFnPerfect: This class will generate a response function based on an ideal DC I-V curve (i.e., the subgap current is exactly zero below the gap voltage and exactly equal to the bias voltage above the gap, assuming normalized values). This DC I-V curve has an infinitely sharp transition. Note however that you can smear the transition using the v_smear argument. This will convolve the ideal reasponse function with a Gaussian distribution, allowing you to control the sharpness of the transition.

    • RespFnPolynomial: This class will generate a response function based on the polynomial model from Kennedy (1999). The order of the polynomial controls the sharpness of the transition, so this class can be used to simulate the effect of the transition’s sharpness (e.g., how does the gain change when the gap is more or less sharp?).

    • RespFnExponential: This class will generate a response function based on the exponential model from Rashid et al. (2016). This model is very similar to RespFnPolynomial, except that you can include a leakage current (i.e., a finite subgap resistance).

Upon initialization, these classes will:

  1. Load or calculate the DC I-V curve (which is the imagainary part of response function).

  2. Calculate the Kramers-Kronig transform of the DC I-V curve (which is the real part of the response function).

  3. Setup the density of the data points to optimize interpolation.

    • The response function needs enough data points that it can be interpolated accurately, but at the same time, not so many points that the interpolation takes too long.

  4. Calculate cubic spline fits for the I-V curve and the KK transform.

    • Doing this once at the start allows the data to be interpolated very quickly later on.

Once initialized, the classes allow the user to interpolated the DC I-V curve, the KK transform, the derivative of the I-V curve, the derivative of the KK transform, and the response function (a complex value). These classes are optimized to interpolate very quickly.

Examples

For a quick example, we will generate a response function using the polynomial model, with polynomial order 50:

>>> resp = RespFnPolynomial(50, verbose=False)

You can then interpolate the DC I-V curve and the KK transform:

>>> bias_voltage = np.array([0.5, 1.0, 2.0])
>>> dc_current = resp.idc(bias_voltage)
>>> np.around(dc_current, 1)
array([0. , 0.5, 2. ])
>>> kk_current = resp.ikk(bias_voltage)
>>> np.around(kk_current, 1)
array([-0.5,  1.1,  0.1])

You can also interpolate the response function directly, which is a complex array:

>>> resp = resp(bias_voltage)
>>> np.around(resp, 1)
array([-0.5+0.j ,  1.1+0.5j,  0.1+2.j ])

Here, the real part is the KK transform (same as the previous kk_current results) and the imaginary part is the DC I-V curve (same as the previous dc_current results).

class qmix.respfn.RespFn(voltage, current, **params)

Bases: object

Generate the response function from pre-processed I-V data.

This class expects “pre-processed” I-V data, meaning that there are more data points around the curvier regions of the I-V curve and fewer points in the linear regions. Sampling in this way helps the interpolation run as quick as possible.

If you have not “pre-processed” your data, please use RespFnFromIVData.

Parameters
  • voltage (ndarray) – normalized DC bias voltage

  • current (ndarray) – normalized DC tunneling current

Keyword Arguments
  • verbose (bool, default is True) – print info to terminal?

  • max_npts_dc (int, default is 101) – maximum number of points in DC I-V curve

  • max_npts_kk (int, default is 151) – maximum number of points in KK transform

  • max_interp_error (float, default is 0.001) – maximum interpolation error (in units of normalized current)

  • check_error (bool, default is False) – check interpolation error?

  • v_smear (float, default is None) – smear DC I-V curve by convolving with a Gaussian dist. with this std. dev.

  • kk_n (int, default is 50) – padding for Hilbert transform (see qmix.mathfn.kktrans.kk_trans)

  • spline_order (int, default is 3) – spline order for interpolations

voltage

The DC bias voltage values.

Type

ndarray

current

The DC tunneling current values.

Type

ndarray

voltage_kk

The DC bias voltage values that correspond to current_kk.

Type

ndarray

current_kk

The values of the KK transform of the DC I-V curve.

Type

ndarray

didc(vbias)

Interpolate the derivative of the DC I-V curve at the given bias voltage.

This is defined as d(idc) / d(vb) where idc is the DC tunneling current and vb is the bias voltage.

Note

This method is not used directly by QMix, but it can be useful if you are calculating the tunneling currents using Tucker theory (see: Tucker and Feldman, 1985).

Parameters

vbias (ndarray) – Bias voltage (normalized)

Returns

derivative of the DC tunneling current

Return type

ndarray

dikk(vbias)

Interpolate the derivative of the Kramers-Kronig transform.

This is defined as d(ikk) / d(vb) where ikk is the Kramers- Kronig transform of the DC tunneling current and vb is the bias voltage.

Note

This method is not used directly by QMix, but it can be useful if you are calculating the tunneling currents using Tucker theory (see: Tucker and Feldman, 1985).

Parameters

vbias (ndarray) – Bias voltage (normalized)

Returns

derivative of the KK transform of the DC I-V curve

Return type

ndarray

idc(vbias)

Interpolate the DC I-V curve.

This is the imaginary component of the respones function, and it is used to calculate the quasiparticle tunneling currents in qmix.qtcurrent.qtcurrent.

Parameters

vbias (ndarray) – Bias voltage (normalized)

Returns

DC tunneling current

Return type

ndarray

ikk(vbias)

Interpolate the Kramers-Kronig transform of the DC I-V curve at the given bias voltage.

This is the real component of the response function, and it is used to calculate the quasiparticle tunneling currents in qmix.qtcurrent.qtcurrent.

Parameters

vbias (ndarray) – Bias voltage (normalized)

Returns

KK transform of the DC I-V curve

Return type

ndarray

plot(fig_name=None, ax=None)

Plot the response function.

This will plot the real and imaginary parts separately.

Note: If fig_name is provided, this method will save the plot to the specified folder and then close the plot. This means that the Matplotlib axis object will not be returned in this case. This is done to prevent too many plots from being open at the time.

Parameters
  • fig_name (str, default is None) – name of figure file name, if you wish to save

  • ax (matplotlib.axes.Axes, default is None) – figure axis, if you would like to add to an existing figure

Returns

figure axis (only if fig_name is

None)

Return type

matplotlib.axes.Axes

plot_interpolation(fig_name=None, ax=None)

Plot the interpolation of the response function.

This can be used to check the interpolation of the response function. (Mostly just a sanity check.)

Note: If fig_name is provided, this method will save the plot to the specified folder and then close the plot. This means that the Matplotlib axis object will not be returned in this case. This is done to prevent too many plots from being open at the time.

Parameters
  • fig_name (str, default is None) – name of figure file name, if you wish to save

  • ax (matplotlib.axes.Axes, default is None) – figure axis, if you would like to add to an existing figure

Returns

figure axis

Return type

matplotlib.axes.Axes

resp(vbias)

Interpolate the response function.

Parameters

vbias (ndarray) – Bias voltage (normalized)

Returns

Response function (a complex value)

Return type

ndarray

resp_conj(vbias)

Interpolate the complex conjugate of the response function.

Note

This method is not used directly by QMix, but it can be useful if you are calculating the tunneling currents using Tucker theory (see: Tucker and Feldman, 1985).

This method is included because it might be slightly faster than np.conj(resp(vb)) where resp is an instance of this class.

Parameters

vbias (ndarray) – Bias voltage (normalized)

Returns

Complex conjugate of the response function

Return type

ndarray

resp_swap(vbias)

Interpolate the response function, with the real and imaginary components swapped.

Note

This method is not used directly by QMix, but it can be useful if you are calculating the tunneling currents using Tucker theory (see: Tucker and Feldman, 1985).

This method is included because it might be slightly faster than 1j*np.conj(resp(vb)) where resp is an instance of this class. This is the normal way that you would swap the real and imaginary components.

Parameters

vbias (ndarray) – Bias voltage (normalized)

Returns

Response function with the real and imaginary components

swapped

Return type

ndarray

show_current(fig_name=None, ax=None)

Plot the interpolation of the response function.

This can be used to check the interpolation of the response function. (Mostly just a sanity check.)

Note: If fig_name is provided, this method will save the plot to the specified folder and then close the plot. This means that the Matplotlib axis object will not be returned in this case. This is done to prevent too many plots from being open at the time.

Warning

This function is deprecated. Please use plot_interpolation instead. I renamed this function to be more consistent across the QMix package.

Parameters
  • fig_name (string) – figure name if saved

  • ax – figure axis

class qmix.respfn.RespFnExponential(vgap=0.0028, rn=14, rsg=1000, agap=40000.0, **kwargs)

Bases: qmix.respfn.RespFn

Response function based on the exponential I-V curve model.

This model is from Rashid et al. (2016). Through this model you can set the sharpness of the non-linearity and the subgap resistance.

See qmix.mathfn.ivcurve_models.exponential for the model.

Parameters
  • vgap (float) – Gap voltage (un-normalized)

  • rsg (float) – Sub-gap resistance (un-normalized)

  • rn (float) – Normal resistance (un-normalized)

  • a (float) – Gap smearing parameter (4e4 is typical)

Keyword Arguments
  • verbose (bool, default is True) – print info to terminal?

  • max_npts_dc (int, default is 101) – maximum number of points in DC I-V curve

  • max_npts_kk (int, default is 151) – maximum number of points in KK transform

  • max_interp_error (float, default is 0.001) – maximum interpolation error (in units of normalized current)

  • check_error (bool, default is False) – check interpolation error?

  • v_smear (float, default is None) – smear DC I-V curve by convolving with a Gaussian dist. with this std. dev.

  • kk_n (int, default is 50) – padding for Hilbert transform (see qmix.mathfn.kktrans.kk_trans)

  • spline_order (int, default is 3) – spline order for interpolations

class qmix.respfn.RespFnFromIVData(voltage, current, **kwargs)

Bases: qmix.respfn.RespFn

Generate the response function from I-V data.

Unlike RespFn, this class will resample the I-V data to optimize the interpolation.

Note

This class expects normalized I-V data that extends from at least vb=0 to vb=vlimit, where vb is the bias voltage and vlimit is one of the keyword arguments.

Parameters
  • voltage (ndarray) – normalized DC bias voltage

  • current (ndarray) – normalized DC tunneling current

Keyword Arguments
  • verbose (bool, default is True) – print info to terminal?

  • max_npts_dc (int, default is 101) – maximum number of points in DC I-V curve

  • max_npts_kk (int, default is 151) – maximum number of points in KK transform

  • max_interp_error (float, default is 0.001) – maximum interpolation error (in units of normalized current)

  • check_error (bool, default is False) – check interpolation error?

  • v_smear (float, default is None) – smear DC I-V curve by convolving with a Gaussian dist. with this std. dev.

  • kk_n (int, default is 50) – padding for Hilbert transform (see qmix.mathfn.kktrans.kk_trans)

  • spline_order (int, default is 3) – spline order for interpolations

  • vlimit (float, default is 1.8) – import all DC I-V data from vb=0 to vb=vlimit, where vb is the bias voltage normalized to the gap voltage.

class qmix.respfn.RespFnPerfect(**params)

Bases: qmix.respfn.RespFn

Response function based on the perfect I-V curve model.

The perfect I-V curve has zero subgap current below the transition, and a current exactly equal to vb * Rn, where vb is the bias voltage and Rn is the normal resistance, above the transition.

See qmix.mathfn.ivcurve_models.perfect for the model.

Keyword Arguments
  • verbose (bool, default is True) – print info to terminal?

  • max_npts_dc (int, default is 101) – maximum number of points in DC I-V curve

  • max_npts_kk (int, default is 151) – maximum number of points in KK transform

  • max_interp_error (float, default is 0.001) – maximum interpolation error (in units of normalized current)

  • check_error (bool, default is False) – check interpolation error?

  • v_smear (float, default is None) – smear DC I-V curve by convolving with a Gaussian dist. with this std. dev.

  • kk_n (int, default is 50) – padding for Hilbert transform (see qmix.mathfn.kktrans.kk_trans)

  • spline_order (int, default is 3) – spline order for interpolations

class qmix.respfn.RespFnPolynomial(p_order=50, **kwargs)

Bases: qmix.respfn.RespFn

Response function based on the polynomial I-V curve model.

This model is from Kennedy (1999). The order of the polynomial (p_order) controls the sharpness of the non-linearity.

See qmix.mathfn.ivcurve_models.polynomial for the model.

Parameters

p_order (int) – Order of the polynomial

Keyword Arguments
  • verbose (bool, default is True) – print info to terminal?

  • max_npts_dc (int, default is 101) – maximum number of points in DC I-V curve

  • max_npts_kk (int, default is 151) – maximum number of points in KK transform

  • max_interp_error (float, default is 0.001) – maximum interpolation error (in units of normalized current)

  • check_error (bool, default is False) – check interpolation error?

  • v_smear (float, default is None) – smear DC I-V curve by convolving with a Gaussian dist. with this std. dev.

  • kk_n (int, default is 50) – padding for Hilbert transform (see qmix.mathfn.kktrans.kk_trans)

  • spline_order (int, default is 3) – spline order for interpolations

2.7. Module contents

Quantum Mixing Software

QMix is used to simulate Superconductor/Insulator/Superconductor (SIS) mixers. It uses multi-tone spectral domain analysis, which makes QMix ideal for simulating higher-order harmonics, power saturation and wide IF bandwidth devices.