qtt.algorithms package
Methods for analysis of quantom dots and spin-qubits
Submodules
qtt.algorithms.allxy module
- qtt.algorithms.allxy.allxy_model(indices: float | ndarray, offset0: float, slope0: float, offset1: float, slope1: float, offset2: float, slope2: float) float | ndarray [source]
Model for AllXY experiment
The model consists of three linear segments. The segments correspond to the pairs of gates that result in fraction 0, 0.5 and 1 in the AllXY experiment.
- Parameters:
index – Indices of the allxy pairs or a single index
offset0 – Offset of first segment
slope0 – Slope of first segment
offset1 – Offset of second segment
slope1 – Slope of second segment
offset2 – Offset of last segment
slope2 – Slope of last segment
- Returns:
Fractions for the allxy pairs
- qtt.algorithms.allxy.fit_allxy(dataset: DataSet, initial_parameters: ndarray | None = None) Dict[str, Any] [source]
Fit AllXY measurement to piecewise linear model
- Parameters:
dataset – Dataset containing the AllXY measurement
initial_parameters – Optional set of initialization parameters
- Returns:
Dictionary with the fitting results
- qtt.algorithms.allxy.generate_allxy_combinations() List[Any] [source]
Generate all combinations of the AllXY sequence from Reed 2013
- qtt.algorithms.allxy.plot_allxy(dataset: DataSet, result: Dict[str, Any], fig: int | Axes | None = 1, plot_initial_estimate: bool = False)[source]
Plot the results of an AllXY fit
- Parameters:
dataset – Dataset containing the measurement data
result – Fitting result of fit_allxy
fig – Figure or axis handle. Is passed on to get_axis
plot_initial_guess – If True, then plot the initial estimate of the model
qtt.algorithms.anticrossing module
Functions to analyse anti-crossings in charge stability diagrams.
Created on Wed Jul 12 08:08:57 2017
@author: diepencjv / eendebakpt
- qtt.algorithms.anticrossing.fit_anticrossing(dataset, width_guess=None, angles_guess=None, psi=None, w=2.5, diff_dir='dx', plot=False, verbose=1, param={})[source]
Fits an anti-crossing model to a 2D scan.
The model fitted is without tunnel coupling.
- Parameters:
dataset (qcodes dataset) – the 2D scan measurement dataset
width_guess (None or float) – Initial estimate for width of anti-crossing in the same unit as the dataset axis.
angles_guess (None or float) – Initial estimate for angles of anti-crossing in radians
psi (None or float) – angle for cross section of anti-crossing
w (float) – Width of lines in anti-crossing model
diff_dir (str) – differentiation direction
plot (int) – If non-zero, then plot into specified matplotlib window
verbose (int) – verbosity level
param (dict) – additional parameters
- Returns:
- (dict) the parameters describing the fitted anti-cross
optionally the cost for fitting and the figure number
- Return type:
anticross_fit
The main fields in the result structure are:
centre (2x1 array): position of the fitted anti-crossing
fit_params (array): fitting parameters of the model. See
fitModel()
- fitpoints (dict): a dictionary with fitted points. centre is the centre point; left_point and right_point are the corners
of the anti-crossing. the four lines outwards are define by inner_points and outer_points
qtt.algorithms.awg_to_plunger module
Functionality to determine the awg_to_plunger ratio.
- qtt.algorithms.awg_to_plunger.analyse_awg_to_plunger(result: dict, method: str = 'hough', fig: int | None = None) dict [source]
- Determine the awg_to_plunger conversion factor from a 2D scan, two possible methods: ‘hough’ it fits the slope
of the addition line and calculates the correction to the awg_to_plunger conversion factor from there. if this doesn’t work for some reason, method ‘click’ can be used to find the addition lines by hand/eye.
- Parameters:
result – result dictionary of the function measure_awg_to_plunger, shape: result = {‘type’: ‘awg_to_plunger’, ‘awg_to_plunger’: None, ‘dataset’: ds.location}.
method – either ‘hough’ or ‘click’.
fig – determines of the analysis staps and the result is plotted.
- Returns:
- including to following entries:
angle (float): angle in radians. angle_degrees (float): angle in degrees. correction of awg_to_plunger (float): correction factor. dataset (str): location where the dataset is stored. type(str): type of calibration, ‘awg_to_plunger’.
- Return type:
result
- qtt.algorithms.awg_to_plunger.click_line(fig: int | None, show_points: bool = False) Tuple[float, float] [source]
Define a line through two points. The points are choosen by clicking on a position in a plot, two times. :param fig: number of figure to plot the points in, when None it will plot a figure. :param showpoints: If True, then plot the selected points in the figure
- Returns:
offset of the line. slope (float): slope of the line.
- Return type:
offset (float)
- qtt.algorithms.awg_to_plunger.measure_awg_to_plunger(station, gate, minstrument, scanrange=30, step=0.5)[source]
Performing a scan2Dfast measurement, same gate on both axis, where the one axis is sweeped with the awg and one axis is stepped with the DAC’s. Measurement should be centred around an addition line. From the slope of the addition line the awg to plunger conversion factor can be checked with the function analyse_awg_to_plunger.
- Parameters:
station (QCoDeS station) – measurement setup.
gate (str) – gate for which the awg to plunger conversion.
minstrument (str, int) – list with the name of the measurement instrument (str), and the channel number (int).
scanrange (float) – sweep- and steprange (mV), making a square 2d measurement.
step (float) – stepsize (mV).
- Returns:
- resultresult (dic): result dictionary of the function measure_awg_to_plunger,
shape: result = {‘type’: ‘awg_to_plunger’, ‘awg_to_plunger’: None, ‘dataset’: ds.location}.
- Return type:
result (dict)
qtt.algorithms.bias_triangles module
Functionality to analyse bias triangles
@author: amjzwerver
- qtt.algorithms.bias_triangles.E_charging(lev_arm, results, fig=None)[source]
Calculates the charging energy of a dot by using charge stability diagrams. Uses currently active figure.
- Parameters:
lev_arm (float) – lever arm for the gate to the dot.
results – dictionary returned from the function perpLineIntersect containing three points, the intersection point between a line through 1,2 and the third point and the length from points 3 to the intersection (horz/vert).
- Returns:
the charging energy for the dot.
- Return type:
E_charging (float)
- qtt.algorithms.bias_triangles.lever_arm(bias, results, fig=None)[source]
Calculates the lever arm of a dot by using bias triangles in charge sensing. Uses currently active figure.
- Parameters:
bias (float) – bias in uV between source and drain while taking the bias triangles.
results (dict) – dictionary returned from the function perpLineIntersect containing three points, the intersection point between a line through 1,2 and the third point and the length from points 3 to the intersection (horz/vert).
fig (bool) – adds lever arm to title of already existing figure with points.
- Returns:
the lever arm of the assigned dot in uV/mV.
- Return type:
lev_arm (float)
- qtt.algorithms.bias_triangles.perpLineIntersect(ds, description, vertical=True, points=None, fig=588, diff_dir='xy')[source]
Takes three points in a graph and calculates the length of a linepiece between a line through points 1,2 and a vertical/horizontal line through the third point. Uses the currently active figure.
- Parameters:
ds (dataset) – dataset with charge stability diagram and gate voltage in mV.
description – TODO
vertical (bool) – find intersection of point with line vertically (True) or horizontally (False).
points (None or array) – if None, then let the user select points.
fig (int) – figure number.
diff_dir (None or 'xy') – specification of differentiation direction.
- Returns:
- ‘intersection_point’ = intersection point
’distance’ = length of line from 3rd clicked point to line through clicked points 1 and 2 ‘clicked_points’ = coordinates of the three clicked points
- Return type:
(dict)
- qtt.algorithms.bias_triangles.plotAnalysedLines(clicked_pts, linePoints1_2, linePt3_vert, linePt3_horz, linePt3_ints, intersect_point)[source]
Plots lines based on three points clicked.
- Parameters:
clicked_pts (array) – lines between the three points (1-2), (2-3).
linePoints1_2 (array) – line fitted through points 1 and 2.
linePt3_vert (array) – vertical line through point 3.
linePt3_horz (array) – horizontal line through point 3.
linePt3_ints (array) – line through point 3 and its vert/horz intersection with the line through point 1,2.
intersect_point (array) – intersection point point 3, line1_2.
qtt.algorithms.chargesensor module
Analyse effect of sensing dot on fitting of tunnel barrier
@author: eendebakpt
- class qtt.algorithms.chargesensor.DataLinearizer(xsr, ysr)[source]
Bases:
object
- backward_curve(y)[source]
Evaluate a B-spline, see documentation scipy.interpolate.splev
- Parameters:
y (array) – to be fitted curve
- Returns:
backward array representing the evaluated spline function
- Return type:
x (array)
- forward(x)[source]
Apply fitted linear transformation to data
- Parameters:
x (array) – data to be fit
- Returns:
fitted data
- Return type:
y (array)
- forward_curve(x)[source]
Evaluate a B-spline, see documentation scipy.interpolate.splev
- Parameters:
x (array) – to be fitted curve
- Returns:
array representing the evaluated spline function
- Return type:
y (array)
- qtt.algorithms.chargesensor.correctChargeSensor(xscan, yscan, xs, ys, fig=None)[source]
Calculate correction for a non-linear charge sensor
- Parameters:
xscan – data of scan to be corrected
yscan – data of scan to be corrected
xs – scan of charge sensor response
ys – scan of charge sensor response
- Returns:
dictionary with intermediate results
- Return type:
results (dict)
qtt.algorithms.coulomb module
Functions to fit and analyse Coulomb peaks
- qtt.algorithms.coulomb.analyseCoulombPeaks(all_data, fig=None, verbose=1, parameters=None)[source]
Find Coulomb peaks in a 1D dataset
- qtt.algorithms.coulomb.analyseCoulombPeaksArray(x_data, y_data, fig=None, verbose=1, parameters=None)[source]
- Find Coulomb peaks in arrays of data. This is very similar to analyseCoulombPeaks,
but takes arrays of data as input. Hence the y_data can for example be either the I, Q or any combination of both obtained with RF reflectometry.
- qtt.algorithms.coulomb.analysePeaks(x, y, peaks, verbose=1, doplot=0, typicalhalfwidth=None, parameters=None, istep=None)[source]
Analyse Coulomb peaks
- Parameters:
x – arrays with data data in mV
y – arrays with data data in mV
peaks – list of detected peaks to be analysed
typicalhalfwidth – float typical width of peak (half side) in mV (mV ??)
- qtt.algorithms.coulomb.coulombPeaks(x_data, y_data, verbose=1, fig=None, plothalf=False, sampling_rate=None, parameters=None)[source]
Detects the Coulumb peaks in a 1D scan.
- qtt.algorithms.coulomb.filterOverlappingPeaks(goodpeaks, threshold=0.6, verbose=0)[source]
Filter peaks based on overlap
- qtt.algorithms.coulomb.filterPeaks(x, y, peaks, verbose=1, minheight=None)[source]
Filter the detected peaks
- Parameters:
x (array) – independent variable data
y (array) – signal data
peaks (list) – list of peaks
- Returns:
selected good peaks
- Return type:
Filtering criteria are:
minimal height
overlapping of peaks
- qtt.algorithms.coulomb.findBestSlope(x, y, minimal_derivative=None, fig=None, verbose=1)[source]
Find good slopes to use in sensing dot
- Parameters:
x (array) – independent variable data
y (array) – dependent variable data
minimal_derivative (None or float) – minimal derivative
- Returns:
slopes (…) results (object): additional data
- qtt.algorithms.coulomb.findSensingDotPosition(x, y, verbose=1, fig=None, plotLabels=True, plotScore=True, plothalf=False, useslopes=True, sampling_rate=None, parameters=None)[source]
Find best position for sensing dot
- Parameters:
x (array) – data
y (array) – data
verbose (int): output level
- Returns:
list of detected positions
- Return type:
goodpeaks (list)
- qtt.algorithms.coulomb.fitCoulombPeaks(x_data, y_data, lowvalue=None, verbose=1, fig=None, sampling_rate=1)[source]
Fit Coulumb peaks in a measurement series.
- qtt.algorithms.coulomb.fitPeaks(XX, YY, points, fig=None, verbose=0)[source]
Fit Gaussian model on local maxima
- qtt.algorithms.coulomb.peakFindBottom(x, y, peaks, fig=None, verbose=1)[source]
Find the left bottom of a detected peak
- qtt.algorithms.coulomb.peakScores(peaksx, x, y, hwtypical=10, verbose=1, fig=None)[source]
Calculate scores for list of peaks
- Parameters:
x (array) – independent variable data
y (array) – dependent variable data
peaksx (list) – list with detected peaks
- qtt.algorithms.coulomb.peakdataOrientation(x, y)[source]
For measured 1D scans order the data such that the independent variable is ordered
- Parameters:
x (array) – independent variable data
y (array) – dependent variable data
- Returns:
reordered data
- Return type:
x,y (array)
qtt.algorithms.fitting module
Fitting of various models.
- qtt.algorithms.fitting.extract_lmfit_parameters(lmfit_model: Model, lmfit_result: ModelResult) Dict[str, Any] [source]
Convert lmfit results to a dictionary
- Parameters:
lmfit_model – Model that was fitted
lmfit_result – Fitting results of lmfit
- Returns:
Dictionary with fitting parameters
- qtt.algorithms.fitting.fitFermiLinear(x_data, y_data, verbose=0, fig=None, lever_arm=1.16, l=None, use_lmfit=False)[source]
Fit data to a Fermi-Linear function
- Parameters:
x_data (array) – independent variable data
y_data (array) – dependent variable data
verbose (int) – verbosity (0 == silent). Not used
fig (int) – figure number
lever_arm (float) – leverarm passed to FermiLinear function
l (Any) – Deprecated parameter. Use lever_arm instead
use_lmfit (bool) – If True use lmfit for optimization, otherwise use scipy
- Returns:
fitted function parameters results (dict): extra fitting data
- Return type:
p (array)
See also
FermiLinear
- qtt.algorithms.fitting.fit_addition_line(dataset, trimborder=True)[source]
Fits a FermiLinear function to the addition line and finds the middle of the step.
- Parameters:
dataset (qcodes dataset) – The 1d measured data of addition line.
trimborder (bool) – determines if the edges of the data are taken into account for the fit.
- Returns:
x value of the middle of the addition line result_dict (dict): dictionary with the following results:
fit parameters (array): fit parameters of the Fermi Linear function parameters initial guess (array): parameters of initial guess dataset fit (qcodes dataset): dataset with fitted Fermi Linear function dataset initial guess (qcodes dataset):dataset with guessed Fermi Linear function
- Return type:
m_addition_line (float)
See also: FermiLinear and fitFermiLinear
- qtt.algorithms.fitting.fit_addition_line_array(x_data, y_data, trimborder=True)[source]
Fits a FermiLinear function to the addition line and finds the middle of the step.
Note: Similar to fit_addition_line but directly works with arrays of data.
- Parameters:
x_data (array) – independent variable data
y_data (array) – dependent variable data
trimborder (bool) – determines if the edges of the data are taken into account for the fit
- Returns:
x value of the middle of the addition line pfit (array): fit parameters of the Fermi Linear function pguess (array): parameters of initial guess
- Return type:
m_addition_line (float)
- qtt.algorithms.fitting.fit_double_gaussian(x_data, y_data, maxiter=None, maxfun=None, verbose=1, initial_params=None)[source]
Fitting of double gaussian
Fitting the Gaussians and finding the split between the up and the down state, separation between the max of the two gaussians measured in the sum of the std.
- Parameters:
x_data (array) – x values of the data
y_data (array) – y values of the data
maxiter (int) – Legacy argument, not used any more
maxfun (int) – Legacy argument, not used any more
verbose (int) – set to >0 to print convergence messages
initial_params (None or array) – optional, initial guess for the fit parameters: [A_dn, A_up, sigma_dn, sigma_up, mean_dn, mean_up]
- Returns:
fit parameters of the double gaussian: [A_dn, A_up, sigma_dn, sigma_up, mean_dn, mean_up] result_dict (dict): dictionary with results of the fit. Fields guaranteed in the dictionary:
parameters (array): Fitted parameters parameters initial guess (array): initial guess for the fit parameters, either the ones give to the
function, or generated by the function: [A_dn, A_up, sigma_dn, sigma_up, mean_dn, mean_up]
reduced_chi_squared (float): Reduced chi squared value of the fit separation (float): separation between the max of the two gaussians measured in the sum of the std split (float): value that separates the up and the down level left (array), right (array): Parameters of the left and right fitted Gaussian
- Return type:
par_fit (array)
- qtt.algorithms.fitting.fit_gaussian(x_data, y_data, maxiter=None, maxfun=None, verbose=0, initial_parameters=None, initial_params=None, estimate_offset=True)[source]
Fitting of a gaussian, see function ‘gaussian’ for the model that is fitted
The final optimization of the fit is performed with lmfit <https://lmfit.github.io/lmfit-py/> using the least_squares method.
- Parameters:
x_data (array) – x values of the data
y_data (array) – y values of the data
verbose (int) – set positive for verbose fit
initial_parameters (None or array) – optional, initial guess for the fit parameters: [mean, s, amplitude, offset]
estimate_offset (bool) – If True then include offset in the Gaussian parameters
maxiter (int) – Legacy argument, not used
maxfun (int) – Legacy argument, not used
- Returns:
fit parameters of the gaussian: [mean, s, amplitude, offset] result_dict (dict): result dictonary containging the fitparameters and the initial guess parameters
- Return type:
par_fit (array)
- qtt.algorithms.fitting.fit_sine(x_data: ndarray, y_data: ndarray, initial_parameters=None, positive_amplitude=True) Tuple[Dict[str, Any], Dict[str, Any]] [source]
Fit a sine wave for the inputted data; see sine function in functions.py for model
- Parameters:
x_data – x data points
y_data – data to be fitted
initial_parameters – list of 4 floats with initial guesses for: amplitude, frequency, phase and offset
positive_amplitude – If True, then enforce the amplitude to be positive
- Returns:
result_dict
- qtt.algorithms.fitting.initFermiLinear(x_data, y_data, fig=None)[source]
Initialization of fitting a FermiLinear function.
First the linear part is estimated, then the Fermi part of the function.
- Parameters:
x_data (array) – data for independent variable
y_data (array) – dependent variable
fig (int) – figure number
- Returns:
linear_part (array) fermi_part (array)
- qtt.algorithms.fitting.plot_FermiLinear(x_data, y_data, results, fig=10)[source]
Plot results for fitFermiLinear
qtt.algorithms.functions module
Mathematical functions and models
- qtt.algorithms.functions.Fermi(x, cc, A, T, kb=1)[source]
Fermi distribution
- Parameters:
- Returns:
value of the function
- Return type:
y (numpy array)
\[y = A*(1/ (1+\exp( (x-cc)/(kb*T) ) ) )\]
- qtt.algorithms.functions.FermiLinear(x, a, b, cc, A, T, l=1.16)[source]
Fermi distribution with linear function added
- Parameters:
- The default value of the leverarm is
(100 ueV/mV)/kb = (100*1e-6*scipy.constants.eV )/kb = 1.16.
For this value the input variable x should be in mV, the temperature T in K. We input the leverarm divided by kb for numerical stability.
- Returns:
value of the function
- Return type:
y (numpy array)
\[y = a*x + b + A*(1/ (1+\exp( l* (x-cc)/(T) ) ) )\]
- qtt.algorithms.functions.cost_exp_decay(x_data, y_data, params, threshold=None)[source]
Cost function for exponential decay.
- Parameters:
x_data (array) – the data for the input variable
y_data (array) – the data for the measured variable
params (array) – parameters of the exponential decay function, [A,B, gamma]
threshold (float or None or 'auto') –
- if the difference between data and model is larger then the threshold,
then the cost penalty is reduced.
If None use normal cost function. If ‘auto’ use automatic detection (at 95th percentile)
- Returns:
value which indicates the difference between the data and the fit
- Return type:
cost (float)
- qtt.algorithms.functions.cost_gauss_ramsey(x_data, y_data, params, weight_power=0)[source]
Cost function for gauss_ramsey.
- Parameters:
x_data (array) – the data for the input variable
y_data (array) – the data for the measured variable
params (array) – parameters of the gauss_ramsey function, [A,C,ramseyfreq,angle,B]
weight_power (float) –
- Returns:
value which indicates the difference between the data and the fit
- Return type:
cost (float)
- qtt.algorithms.functions.double_gaussian(x_data, params)[source]
A model for the sum of two Gaussian distributions.
- Parameters:
x_data (array) – x values of the data
params (array) – parameters of the two gaussians, [A_dn, A_up, sigma_dn, sigma_up, mean_dn, mean_up] amplitude of first (second) gaussian = A_dn (A_up) standard deviation of first (second) gaussian = sigma_dn (sigma_up) average value of the first (second) gaussian = mean_dn (mean_up)
- Returns:
model of a double gaussian
- Return type:
double_gauss (np.array)
- qtt.algorithms.functions.estimate_dominant_frequency(signal, sample_rate=1, remove_dc=True, fig=None)[source]
Estimate dominant frequency in a signal
- qtt.algorithms.functions.estimate_parameters_damped_sine_wave(x_data, y_data, exponent=2)[source]
Estimate initial parameters of a damped sine wave
The damped sine wave is described in https://en.wikipedia.org/wiki/Damped_sine_wave. This is a first estimate of the parameters, no numerical optimization is performed.
The amplitude is estimated from the minimum and maximum values of the data. The osciallation frequency using the dominant frequency in the FFT of the signal. The phase of the signal is calculated based on the first datapoint in the sequences and the other parameter estimates. Finally, the decay factor of the damped sine wave is determined by a heuristic rule.
Example
>>> estimate_parameters_damped_sine_wave(np.arange(10), np.sin(np.arange(10)))
- qtt.algorithms.functions.exp_function(x, a, b, c)[source]
Model for exponential function
$$y = a + b * np.exp(-c * x)$$
- Parameters:
x (array) – x values of the data
offset (a =) –
value (b = starting) –
time (c = 1/typical decay) –
- Returns:
model for exponantial decay
- Return type:
y (array)
- qtt.algorithms.functions.fit_double_gaussian(x_data, y_data, maxiter=None, maxfun=5000, verbose=1, initial_params=None)[source]
- qtt.algorithms.functions.fit_exp_decay(x_data, y_data, maxiter=None, maxfun=5000, verbose=1, initial_params=None, threshold=None, offset_parameter=None)[source]
Fit a exponential decay.
- Parameters:
x_data (array) – the data for the input variable
y_data (array) – the data for the measured variable
maxiter (int) – maximum number of iterations to perform
maxfun (int) – maximum number of function evaluations to make
verbose (int) – set to >0 to print convergence messages
initial_params (None or array) – optional, initial guess for the fit parameters: [A,B, gamma]
threshold (float or None) – threshold for the cost function. If the difference between data and model is larger then the threshold, these data are not taken into account for the fit. If None use automatic detection (at 95th percentile)
offset_parameter (None or float) – if None, then estimate the offset, otherwise fix the offset to the specified value
- Returns:
fit parameters of the exponential decay, [A, B, gamma]
- Return type:
fitted_parameters (array)
See:
exp_function()
- qtt.algorithms.functions.fit_gauss_ramsey(x_data, y_data, weight_power=None, maxiter=None, maxfun=5000, verbose=1, initial_params=None)[source]
Fit a gauss_ramsey. The function gauss_ramsey gives a model for the measurement result of a pulse Ramsey sequence while varying the free evolution time, the phase of the second pulse is made dependent on the free evolution time. This results in a gaussian decay multiplied by a sinus. Function as used by T.F. Watson et all., see function ‘gauss_ramsey’ and example in qtt/docs/notebooks/example_fit_ramsey.ipynb
- Parameters:
x_data (array) – the data for the independent variable
y_data (array) – the data for the measured variable
weight_power (float or None) – If a float then weight all the residual errors with a scale factor
maxiter (int) – maximum number of iterations to perform
maxfun (int) – maximum number of function evaluations to make
verbose (int) – set to >0 to print convergence messages
initial_params (None or array) – optional, initial guess for the fit parameters: [A,C,ramseyfreq,angle,B]
- Returns:
array with the fit parameters: [A,t2s,ramseyfreq,angle,B] result_dict (dict): dictionary containing a description, the par_fit and initial_params
- Return type:
par_fit (array)
- qtt.algorithms.functions.fit_gaussian(x_data, y_data, maxiter=None, maxfun=None, verbose=0, initial_parameters=None, initial_params=None, estimate_offset=True)[source]
- qtt.algorithms.functions.gauss_ramsey(x_data, params)[source]
Model for the measurement result of a pulse Ramsey sequence while varying the free evolution time, the phase of the second pulse is made dependent on the free evolution time. This results in a gaussian decay multiplied by a sinus. Function as used by T.F. Watson et all., example in qtt/docs/notebooks/example_fit_ramsey.ipynb
$$ gauss_ramsey = A * exp(-(x_data/t2s)**2) * sin(2pi*ramseyfreq * x_data - angle) +B $$
- Parameters:
x_data (array) – the data for the input variable
params (array) – parameters of the gauss_ramsey function, [A,t2s,ramseyfreq,angle,B]
- Result:
gauss_ramsey (array): model for the gauss_ramsey
- qtt.algorithms.functions.gaussian(x, mean, std, amplitude=1, offset=0)[source]
Model for Gaussian function
$$y = offset + amplitude * np.exp(-(1/2)*(x-mean)^2/s^2)$$
- Parameters:
x (array) – data points
mean – parameters
std – parameters
amplitude – parameters
offset – parameters
- Returns:
y (array)
- qtt.algorithms.functions.logistic(x, x0=0, alpha=1)[source]
Logistic function
Defines the logistic function
\[y = 1 / (1 + \exp(-2 * alpha * (x - x0)))\]- Parameters:
Example
y = logistic(0, 1, alpha=1)
- qtt.algorithms.functions.plot_gauss_ramsey_fit(x_data, y_data, fit_parameters, fig: int | Axes | None)[source]
Plot Gauss Ramsey fit
- Parameters:
x_data – Input array with time variable (in seconds)
y_data – Input array with signal
fit_parameters – Result of fit_gauss_ramsey (fitting units in seconds)
fig – Figure or axis handle. Is passed to get_axis
- qtt.algorithms.functions.raised_cosine(t: ndarray, roll_off_factor: float, symbol_period: float) ndarray [source]
Raised cosine impulse response function
See https://en.wikipedia.org/wiki/Raised-cosine_filter
- Parameters:
t – Independent variable
roll_off_factor – Roll-off factor
symbol_period – Symbol period
- Returns:
Calculated values of the raised cosine
- qtt.algorithms.functions.raised_cosine_frequency_domain(frequency: ndarray, roll_off_factor: float, symbol_period: float) ndarray [source]
Raised cosine frequency domain function
See https://en.wikipedia.org/wiki/Raised-cosine_filter
- Parameters:
frequency – Frequency
roll_off_factor – Roll-off factor
symbol_period – Symbol period
- Returns:
Calculated values of the raised cosine filter in frequency domain
- qtt.algorithms.functions.sine(x: float | ndarray, amplitude: float, frequency: float, phase: float, offset: float) float | ndarray [source]
Model for sine function
y = offset + amplitude * np.sin(2 * np.pi * frequency * x + phase)
- Parameters:
x – Independent data points
amplitude – Arguments for the sine model
frequency – Arguments for the sine model
phase – Arguments for the sine model
offset – Arguments for the sine model
- Returns:
Calculated data points
qtt.algorithms.gatesweep module
Functionality to analyse pinch-off scans
qtt.algorithms.generic module
Various functions
- qtt.algorithms.generic.boxcar_filter(signal: ndarray, kernel_size: ndarray | Tuple[int]) ndarray [source]
Perform boxcar filtering on an array. At the edges, the edge value is replicated beyond the edge as needed by the size of the kernel. This is the ‘nearest’ mode of scipy.ndimage.convolve. For details, see https://docs.scipy.org/doc/scipy/reference/generated/scipy.ndimage.convolve.html?highlight=mode
- Parameters:
signal – An array containing the signal to be filtered.
kernel_size – Multidimensional size of the filter box. Must have the same number of dimensions as the signal.
- Returns:
Array containing the filtered signal.
- qtt.algorithms.generic.detect_blobs_binary(bim)[source]
Simple blob detection in binary image
- Parameters:
bim (numpy array) – binary input image
- Output:
xx (numpy array): detected blob centres
Alternative implementation would be through findContours
- qtt.algorithms.generic.extent2fullextent(extent0, im)[source]
Convert extent to include half pixel border
- qtt.algorithms.generic.findCoulombDirection(im, ptx, step, widthmv=8, fig=None, verbose=1)[source]
Find direction of Coulomb peak using second order derivative
- qtt.algorithms.generic.flowField(im, fig=None, blocksize=11, ksize=3, resizefactor=1, eigenvec=1)[source]
Calculate flowfield of an image
- Parameters:
im (numpy array) – input image
fig (integer or None) – number of visualization window
- Returns:
flow ll (?): ??
- Return type:
flow (numpy array)
- qtt.algorithms.generic.getValuePixel(imx, pt)[source]
Return interpolated pixel value in an image
- Parameters:
imx (numpy array) – input image
pt (numpy array) – list of points
- Returns:
interpolated value
- Return type:
vv (numpy array)
- qtt.algorithms.generic.localMaxima(arr, radius=1, thr=None)[source]
Calculate local maxima in a 2D array
- qtt.algorithms.generic.makeCoulombFilter(theta0=-0.7853981633974483, step=1, ne=0, dphi=0.7853981633974483, widthmv=10, lengthmv=50.0, verbose=0, fig=None)[source]
Create filters to detect Coulomb peaks
- qtt.algorithms.generic.rescaleImage(im, imextent, mvx=None, mvy=None, verbose=0, interpolation=None, fig=None)[source]
Scale image to make pixels at specified resolution
- Parameters:
- Returns:
transformed image H (array): transformation matrix from units to pixels. H is the homogeneous transform from original to
scaled image
mvx (float): internal data mvy (float): internal data fx (float): internal data dy (float): internal data
- Return type:
ims (array)
- qtt.algorithms.generic.showCoulombDirection(ptx, ww, im=None, dd=None, fig=100)[source]
Show direction of Coulomb peak in an image
- qtt.algorithms.generic.showFlowField(im, flow, ll=None, ff=None, d=12, fig=-1)[source]
Show flow field
- Parameters:
im (numpy array) – input image
flow (numpy array) – input image
fig (integer or None) – number of visualization window
- qtt.algorithms.generic.smoothImage(im, k=3)[source]
Super simple image smoothing
- Parameters:
im (array) – input image
k (int) – kernel size
- Returns:
smoothed image
- Return type:
ims (array)
Example
ims = smoothImage(np.random.rand( 30,40))
- qtt.algorithms.generic.subpixelmax(A, mpos, verbose=0)[source]
Calculate maximum position with subpixel accuracy
For each position specified by mpos this method fits a parabola through 3 points and calculates the maximum position of the parabola.
- Parameters:
A (1D array) – Input data
mpos (array with integer indices) – Positions in the array A with maxima
verbose (int) – Verbosity level
- Returns:
Array with subpixel positions subval (array): Values at maximum positions
- Return type:
subpos (array)
qtt.algorithms.images module
Functionality to analyse and pre-process images
@author: eendebakpt
qtt.algorithms.markov_chain module
Class to generate signals with continous-time Markov chains
@author: pieter.eendebak@gmail.com
- class qtt.algorithms.markov_chain.ChoiceGenerator(number_of_states: int, cum_weights: ndarray, block_size: int = 5000)[source]
Bases:
object
Class to generate random elements with weighted selection
This is a replacement for random.choices that is efficient when a large number of choices has to be generated.
- Parameters:
number_of_states – number of choices that has to be generated
cum_weights (array[float]) – cumulative probabilities of the choices
block_size – size of blocks of choices to generate
- generate_choice() int [source]
Generate a single choice
- Returns:
Integer in the range 0 to the number of states
- class qtt.algorithms.markov_chain.ContinuousTimeMarkovModel(states: List[str], holding_parameters: List[float] | ndarray, jump_chain: ndarray)[source]
Bases:
object
- generate_sequence(length: int, delta_time: float, initial_state: None | int | ndarray = None, generators: List[ChoiceGenerator] | None = None) ndarray [source]
Generate a random sequence with the model
- Parameters:
length – number of elements in the sequence
delta_time – time step to be used. This is equal to one over the samplerate.
initial_state – This parameter determines how the first element of the generated sequence is chosen. If an int, then use that state is initial state. If None then take a random state weighted by the stationary distribution. If the initial_state is a list, then the list is interpreted as a probability distribution and the first element is sampled from all possible states according to the distribution specified.
generators – Optional list of generators to use
- Returns:
Array with generated sequence
- generate_sequences(length: int, delta_time: float = 1, initial_state: None | int | ndarray = None, number_of_sequences: int = 1) ndarray [source]
Generate multiple random sequences with the model
- Parameters:
length – number of elements in the sequence
delta_time – time step to be used
initial_state – This parameter determines how the first element of the generated sequences are chosen. The parameter is passed to the
generate_sequence()
method.number_of_sequences – Specified the number of sequences to generate
- Returns:
Array with generated sequences
- stationary_distribution() ndarray [source]
Return the stationary distribution of the model
The calculation method is taken from: https://vknight.org/unpeudemath/code/2015/08/01/simulating_continuous_markov_chains.html
- stationary_distribution_direct() ndarray [source]
Return the stationary distribution of the model
The calculation method is taken from: https://www.probabilitycourse.com/chapter11/11_3_2_stationary_and_limiting_distributions.php, Theorem 11.3
- static stationary_distribution_discrete(jump_chain) ndarray [source]
Return the stationary distrubution for a Markov chain
- qtt.algorithms.markov_chain.generate_traces(markov_model: ContinuousTimeMarkovModel, std_gaussian_noise: float = 1, state_mapping: ndarray | None = None, *args, **kwargs)[source]
Generate traces for a continuous-time Markov model with added noise
- Parameters:
markov_model – model to use for generation of traces
std_gaussian_noise – standard deviation of Gaussian noise to add to the output signal
state_mapping – If not None, replace each state in the generated trace by the corresponding element in the array
*args – passed to the generate_sequences function of the model
**kwargs –
passed to the generate_sequences function of the model
The traces are generated by the generate_sequences method from the model.
qtt.algorithms.misc module
Misc algorithms
Created on Wed Aug 31 16:19:30 2016
@author: eendebakpt
- qtt.algorithms.misc.fillPoly(im, poly_verts, color=None)[source]
Fill a polygon in an image with the specified color
Replacement for OpenCV function cv2.fillConvexPoly
- Parameters:
im (array) – array to plot into
poly_verts (kx2 array) – polygon vertices
color (array or float) – color to fill the polygon
- Returns:
resulting array
- Return type:
grid (array)
- qtt.algorithms.misc.point_in_poly(x, y, poly)[source]
Return true if a point is contained in a polygon
- qtt.algorithms.misc.points_in_poly(points, poly_verts)[source]
Determine whether points are contained in a polygon or not
- Parameters:
points (kx2 array) –
poly_verts (array) –
- qtt.algorithms.misc.polyarea(p)[source]
Return signed area of polygon
- Parameters:
p (2xN array or list of vertices) – vertices of polygon
- Returns:
area – area of polygon
- Return type:
>>> polyarea( [ [0,0], [1,0], [1,1], [0,2]] ) 1.5
qtt.algorithms.ohmic module
Functionality to fit scans of ohmic contacts
qtt.algorithms.onedot module
Functionality for analysis of single quantum dots
For more details see https://arxiv.org/abs/1603.02274
- qtt.algorithms.onedot.costscoreOD(a, b, pt, ww, verbose=0, output=False)[source]
Cost function for simple fit of one-dot open area
- qtt.algorithms.onedot.onedotGetBalance(dataset, verbose=1, fig=None, drawpoly=False, polylinewidth=2, linecolor='c', full_output=False, od=None)[source]
Determine tuning point from a 2D scan of a 1-dot
This function performs a simple fitting of the open (conducting region).
- Parameters:
od (one-dot structure or None) – data for one-dot
dd (2D dataset) – data containing charge stability diagram
- Returns:
dictionary with fitting results od (obj): modified one-dot object
- Return type:
fitresults (dict)
- qtt.algorithms.onedot.onedotGetBalanceFine(impixel=None, dd=None, verbose=1, fig=None, baseangle=-0.7853981633974483, units=None, full_output=False)[source]
Determine central position of Coulomb peak in 2D scan
The position is determined by scanning with Gabor filters and then performing blob detection
The image should be in pixel coordinates
- Returns:
detected point results (dict): dictionary with all results
- Return type:
pt (array)
qtt.algorithms.pat_fitting module
Functionality to fit PAT models
For more details see: https://arxiv.org/abs/1803.10352
@author: diepencjv / eendebakpt
- qtt.algorithms.pat_fitting.detect_peaks(x_data, y_data, imx, sigmamv=0.25, fig=400, period=0.001, model='one_ele')[source]
Detect peaks in sensor signal, e.g. from a pat scan.
- Parameters:
x_data (array) – detuning (mV)
y_data (array) – frequencies (Hz)
imx (array) – sensor signal of PAT scan, background is usually already subtracted
- Returns:
coordinates of detected peaks results (dict): additional fitting data
- Return type:
detected_peaks (array)
- qtt.algorithms.pat_fitting.fit_pat(x_data, y_data, z_data, background, trans='one_ele', period=0.001, even_branches=[True, True, True], par_guess=None, xoffset=None, verbose=1)[source]
Wrapper for fitting the energy transitions in a PAT scan.
For more details see: https://arxiv.org/abs/1803.10352
- Parameters:
- Returns:
fitted xoffset (mV), leverarm (ueV/mV) and t (ueV) results (dict): contains keys par_guess (array), imq (array) re-scaled and re-centered sensor signal, imextent (array), xd, yd, ydf
- Return type:
pp (array)
- qtt.algorithms.pat_fitting.fit_pat_to_peaks(pp, xd, yd, trans='one_ele', even_branches=[True, True, True], weights=None, xoffset=None, verbose=1, branch_reduction=None)[source]
Core fitting function for PAT measurements, based on detected resonance peaks (see detect_peaks).
- Parameters:
pp (array) – initial guess of fit parameters
xd (array) – x coordinates of peaks in sensor signal (mV)
yd (array) – y coordinates of peaks in sensor signal (Hz)
trans (string) – ‘one_ele’ or ‘two_ele’
xoffset (float) – the offset from zero detuning in voltage. If this has been determined before, then fixing this parameter reduces the fitting time.
- qtt.algorithms.pat_fitting.one_ele_pat_model(x_data, pp)[source]
Model for one electron pat
This is \(\phi=\sqrt{ { ( leverarm * (x-x_0) ) }^2 + 4 t^2 } \mathrm{ueV2Hz}\)
- Parameters:
x_data (array) – detuning (mV)
pp (array) – xoffset (mV), leverarm (ueV/mV) and t (ueV)
For more details see: https://arxiv.org/abs/1803.10352
- class qtt.algorithms.pat_fitting.pat_score(even_branches=[True, True, True], branch_reduction=None)[source]
Bases:
object
- qtt.algorithms.pat_fitting.plot_pat_fit(x_data, y_data, z_data, pp, trans='one_ele', fig=400, title='Fitted model', label='model')[source]
Plot the fitted model of the PAT transition(s)
- Parameters:
x_data (array) – detuning in millivolts
y_data (array) – frequencies
z_data (array) – sensor signal of PAT scan
pp (array) – xoffset (mV), leverarm (ueV/mV) and t (ueV)
model (function) – model describing the PAT transitions
- qtt.algorithms.pat_fitting.pre_process_pat(x_data, y_data, background, z_data, fig=None)[source]
Pre-process a pair of background and sensor signal from a pat scan.
- Parameters:
x_data (array) – detuning (mV)
y_data (array) – frequency (Hz)
background (array) – e.g. sensor signal of POL scan
z_data (array) – sensor signal of PAT scan
fig (None or int) –
- Returns:
imx (array) imq (array) backgr_sm (array)
qtt.algorithms.random_telegraph_signal module
Functionality to analyse random telegraph signals
Created on Wed Feb 28 10:20:46 2018
@author: riggelenfv /eendebakpt
- exception qtt.algorithms.random_telegraph_signal.FittingException[source]
Bases:
Exception
Fitting exception in RTS code
- qtt.algorithms.random_telegraph_signal.generate_RTS_signal(number_of_samples: int = 100000, std_gaussian_noise: float = 0.1, uniform_noise: float = 0.05, rate_up: float = 10000.0, rate_down: float = 15000.0, samplerate: float = 1000000.0) ndarray [source]
Generate a RTS signal
- Parameters:
number_of_samples – Length the the trace to be generated
std_normal_noise – std of Gaussian noise added to the signal
uniform_noise – uniform noise in the range +- uniform_noise/2 is added to the signal
rate_up – rate from down to up
rate_down – rate from up to down
samplerate – The samplerate of the signal to be generated
- Returns:
Array with generated signal (0 is down, 1 is up)
- qtt.algorithms.random_telegraph_signal.plot_two_level_threshold(results: dict, fig: int = 100, plot_initial_estimate: bool = False)[source]
- qtt.algorithms.random_telegraph_signal.rts2tunnel_ratio(binary_signal: ndarray) float [source]
Calculate ratio between tunnelrate down and up
From the mean and standard deviation of the RTS data we can determine the ratio between the two tunnel rates. See equations on https://en.wikipedia.org/wiki/Telegraph_process
- Parameters:
binary_signal – RTS signal with two levels 0 and 1
- Returns:
Ratio of tunnelrate up to down (l2) and down to up (l1)
- qtt.algorithms.random_telegraph_signal.transitions_durations(data: ndarray, split: float, add_start: bool = False, add_end: bool = False) Tuple[ndarray, ndarray] [source]
For data of a two level system (up and down) determine durations of segments
This function determines which datapoints belong to which level and finds the transitions, in order to determines how long the system stays in these levels.
- Parameters:
data (then include the segments at the end of the) – data from the two level system
split – value that separates the up and down level
add_start – If True, then include the segments at the start of the data
add_end: – If True:
data –
- Returns:
array of the durations (unit: data points) in the down level duration_up: array of durations (unit: data points) in the up level
- Return type:
duration_dn
- qtt.algorithms.random_telegraph_signal.tunnelrates_RTS(data: ndarray | DataSet, samplerate: float | None = None, min_sep: float = 2.0, max_sep: float = 7.0, min_duration: int = 5, num_bins: int | None = None, fig: int | None = None, ppt=None, verbose: int = 0, offset_parameter: float | None = None) Tuple[float | None, float | None, dict] [source]
This function takes an RTS dataset, fits a double gaussian, finds the split between the two levels, determines the durations in these two levels, fits a decaying exponential on two arrays of durations, which gives the tunneling frequency for both the levels. If the number of datapoints is too low to get enough points per bin for the exponential fit (for either the up or the down level), this analysis step is passed over. tunnelrate_dn and tunnelrate_up are returned as None, but similar information can be substracted from parameters[‘down_segments’] and parameters[‘up_segments’].
- Parameters:
data – qcodes DataSet (or 1d data array) with the RTS data
samplerate – sampling rate of the acquisition device, optional if given in the metadata of the measured data
min_sep – if the separation found for the fit of the double gaussian is less then this value, the fit probably failed and a FittingException is raised
max_sep – if the separation found for the fit of the double gaussian is more then this value, the fit probably failed and a FittingException is raised
min_duration – minimal number of datapoints a duration should last to be taking into account for the analysis
num_bins – number of bins for the histogram of signal values. If None, then determine based on the size of the data
fig – shows figures and sends them to the ppt when is not None
ppt – determines if the figures are send to a powerpoint presentation
verbose – prints info to the console when > 0
offset_parameter – Offset parameter for fitting of exponential decay
- Returns:
- tunneling rate of the down level to the up level (kHz) or None in case of
not enough datapoints
- tunnelrate_up: tunneling rate of the up level to the down level (kHz) or None in case of
not enough datapoints
- parameters: dictionary with relevent (fit) parameters. this includes:
tunnelrate_down (float): tunnel rate in Hz tunnelrate_up (float): tunnel rate up in Hz
- Return type:
tunnelrate_dn
- qtt.algorithms.random_telegraph_signal.two_level_threshold(data: ndarray, number_of_bins: int = 40) dict [source]
Determine threshold for separation of two-level signal
Typical examples of such a signal are an RTS signal or Elzerman readout.
- Parameters:
data – Two dimensional array with single traces
number_of_bins – Number of bins to use for calculation of double histogram
- Returns:
Dictionary with results. The key readout_threshold contains the calculated threshold
qtt.algorithms.tunneling module
Functionality for analysing inter-dot tunnel frequencies.
@author: diepencjv
- qtt.algorithms.tunneling.data_to_exc_ch(x_data, y_data, pol_fit)[source]
Convert y_data to units of excess charge.
Note: also re-centers to zero detuning in x-direction.
- Parameters:
x_data (1 x N array) – chemical potential difference in ueV.
y_data (1 x N array) – sensor data, e.g. from a sensing dot or QPC.
pol_fit (1 x 6 array) – fit parameters, see
polmod_all_2slopes()
.
- qtt.algorithms.tunneling.fit_pol_all(x_data, y_data, kT, model='one_ele', maxiter=None, maxfun=5000, verbose=1, par_guess=None, method='fmin')[source]
Polarization line fitting.
The default value for the maxiter argument of scipy.optimize.fmin is N*200 the number of variables, i.e. 1200 in our case. :param x_data: chemical potential difference in ueV. :type x_data: 1 x N array :param y_data: sensor data, e.g. from a sensing dot or QPC. :type y_data: 1 x N array :param kT: temperature in ueV. :type kT: float
- Returns:
fitted parameters, see
polmod_all_2slopes()
. par_guess (1 x 6 array): initial guess of parameters for fitting, seepolmod_all_2slopes()
. results (dictionary): dictionary with fitting results.- Return type:
par_fit (1 x 6 array)
- qtt.algorithms.tunneling.fit_pol_all_2(x_data, y_data, kT, model='one_ele', maxiter=None, maxfun=5000, verbose=1, par_guess=None, method='fmin', returnextra=False)[source]
- qtt.algorithms.tunneling.plot_polarization_fit(detuning, signal, results, fig, verbose=1)[source]
Plot the results of a polarization line fit.
- qtt.algorithms.tunneling.pol_mod_two_ele_boltz(x_data, par, kT)[source]
Model of the inter-dot transition with two electron spin states, also taking into account thermal occupation of the triplets.
- qtt.algorithms.tunneling.polmod_all_2slopes(x_data, par, kT, model=None)[source]
Polarization line model.
This model is based on [DiCarlo2004, Hensgens2017]. For an example see: https://github.com/VandersypenQutech/qtt/blob/master/examples/example_polFitting.ipynb
- Parameters:
x_data (1 x N array) – chemical potential difference in ueV.
par (1 x 6 array) – parameters for the model - par[0]: tunnel coupling in ueV - par[1]: offset in x_data for center of transition - par[2]: offset in background signal - par[3]: slope of sensor signal on left side - par[4]: slope of sensor signal on right side - par[5]: height of transition, i.e. sensitivity for electron transition.
kT (float) – temperature in ueV.
() (model) – Not used.
- Returns:
sensor data, e.g. from a sensing dot or QPC.
- Return type:
y_data (array)
- qtt.algorithms.tunneling.polweight_all_2slopes(x_data, y_data, par, kT, model='one_ele')[source]
Cost function for polarization fitting. :param x_data: chemical potential difference in ueV. :type x_data: 1 x N array :param y_data: sensor data, e.g. from a sensing dot or QPC. :type y_data: 1 x N array :param par: see polmod_all_2slopes. :type par: 1 x 6 array :param kT: temperature in ueV. :type kT: float
- Returns:
sum of residues.
- Return type:
total (float)
- qtt.algorithms.tunneling.polweight_all_2slopes_2(x_data, y_data, par, kT, model='one_ele')[source]
Cost function for polarization fitting. :param x_data: chemical potential difference in ueV. :type x_data: 1 x N array :param y_data: sensor data, e.g. from a sensing dot or QPC. :type y_data: 1 x N array :param par: see polmod_all_2slopes. :type par: 1 x 6 array :param kT: temperature in ueV. :type kT: float
- Returns:
sum of residues.
- Return type:
total (float)