silq.analysis package¶
Submodules¶
silq.analysis.analysis module¶
-
silq.analysis.analysis.
find_high_low
(traces, plot=False, threshold_peak=0.02, attempts=8, threshold_method='config', min_voltage_difference='config', threshold_requires_high_low='config', min_SNR=None)[source]¶ Find high and low voltages of traces using histograms
This function determines the high and low voltages of traces by binning them into 30 bins, and trying to discern two peaks. Useful for determining the threshold value for measuring a blip.
If no two peaks can be discerned after all attempts, None is returned for each of the returned dict keys except DC_voltage.
- Parameters
traces (
ndarray
) – 2D array of acquisition tracesplot (
bool
) – Whether to plot the histogramsthreshold_peak (
float
) – Threshold for discerning a peak. Will be varied if too many/few peaks are foundattempts (
int
) – Maximum number of attempts for discerning two peaks. Each attempt the threshold_peak is decreased/increased depending on if too many/few peaks were foundthreshold_method (
str
) –Method used to determine the threshold voltage. Allowed methods are:
mean: average of high and low voltage.
{n}*std_low: n standard deviations above mean low voltage, where n is a float (ignore slash in raw docstring).
{n}*std_high: n standard deviations below mean high voltage, where n is a float (ignore slash in raw docstring).
config: Use threshold method provided in
config.analysis.threshold_method
(mean
if not specified)
min_voltage_difference (
Union
[float
,str
]) – minimum difference between high and low voltage. If not satisfied, all results are None. Will try to retrieve from config.analysis.min_voltage_difference, else defaults to 0.3Vthreshold_requires_high_low (
Union
[bool
,str
]) – Whether or not both a high and low voltage must be discerned before returning a threshold voltage. If set to False and threshold_method is notmean
, a threshold voltage is always determined, even if no two voltage peaks can be discerned. In this situation, there usually aren’t any blips, or the blips are too short-lived to have a proper high current. When the threshold_method isstd_low
(std_high
), the top (bottom) 20% of voltages are scrapped to ensure any short-lived blips with a high (low) current aren’t included. The average is then taken over the remaining 80% of voltages, which is then the average low (high) voltage. Default is True. Can be set by config.analysis.threshold_requires_high_lowmin_SNR (
Optional
[float
]) – Minimum SNR between high and low voltages required to determine a threshold voltage (default None).
- Returns
low (float): Mean low voltage,
None
if two peaks cannot be discernedhigh (float): Mean high voltage,
None
if no two peaks cannot be discernedthreshold_voltage (float): Threshold voltage for a blip. If SNR is below
min_SNR
or no two peaks can be discerned, returnsNone
.voltage_difference (float): Difference between low and high voltage. If not two peaks can be discerned, returns
None
.DC_voltage (float): Average voltage of traces.
- Return type
Dict[str, Any]
-
silq.analysis.analysis.
edge_voltage
(traces, edge, state, threshold_voltage=None, points=5, start_idx=0)[source]¶ Test traces for having a high/low voltage at begin/end
- Parameters
traces (
ndarray
) – 2D array of acquisition tracesedge (
str
) – which side of traces to test, eitherbegin
orend
state (
str
) – voltage that the edge must have, eitherlow
orhigh
threshold_voltage (
Optional
[float
]) – threshold voltage for ahigh
voltage (blip). If not specified,find_high_low
is used to determine thresholdpoints (
int
) – Number of data points to average over to determine statestart_idx (
int
) – index of first point to use. Useful if there is some capacitive voltage spike occuring at the start. Only used if edge isbegin
.
- Return type
ndarray
- Returns
1D boolean array, True if the trace has the correct state at the edge
-
silq.analysis.analysis.
find_up_proportion
(traces, threshold_voltage=None, return_array=False, start_idx=0, filter_window=0)[source]¶ Determine the up proportion of traces (traces that have blips)
- Parameters
traces (
ndarray
) – 2D array of acquisition tracesthreshold_voltage (
Optional
[float
]) – threshold voltage for ahigh
voltage (blip). If not specified,find_high_low
is used to determine thresholdreturn_array (
bool
) – whether to return the boolean array or the up proportionstart_idx (
int
) – index of first point to use. Useful if there is some capacitive voltage spike occuring at the start. Only used if edge isbegin
filter_window (
int
) – number of points of smoothing (0 means no smoothing)
- Return type
- Returns
- if return_array is False
(float) The proportion of traces with a blip
- else
Boolean array, True if the trace has a blip
-
silq.analysis.analysis.
count_blips
(traces, threshold_voltage, sample_rate, t_skip, ignore_final=False)[source]¶ Count number of blips and durations in high/low state.
- Parameters
- Returns
blips (float): Number of blips per trace.
blips_per_second (float): Number of blips per second.
low_blip_duration (np.ndarray): Durations in low-voltage state.
high_blip_duration (np.ndarray): Durations in high-voltage state.
mean_low_blip_duration (float): Average duration in low state.
mean_high_blip_duration (float): Average duration in high state.
- Return type
Dict[str, Any]
-
silq.analysis.analysis.
analyse_traces
(traces, sample_rate, filter=None, min_filter_proportion=0.5, t_skip=0, t_read=None, segment='begin', threshold_voltage=None, threshold_method='config', plot=False, plot_high_low=False)[source]¶ Analyse voltage, up proportions, and blips of acquisition traces
- Parameters
traces (
ndarray
) – 2D array of acquisition traces.sample_rate (
float
) – acquisition sample rate (per second).filter (
Optional
[str
]) – only use traces that begin in low or high voltage. Allowed values arelow
,high
orNone
(do not filter).min_filter_proportion (
float
) – minimum proportion of traces that satisfy filter. If below this value, up_proportion etc. are not calculated.t_skip (
float
) – initial time to skip for each trace (ms).t_read (
Optional
[float
]) – duration of each trace to use for calculating up_proportion etc. e.g. for a long trace, you want to compare up proportion of start and end segments.segment (
str
) – Use beginning or end of trace fort_read
. Allowed values arebegin
andend
.threshold_voltage (
Optional
[float
]) – threshold voltage for ahigh
voltage (blip). If not specified,find_high_low
is used to determine threshold.threshold_method (
str
) –Method used to determine the threshold voltage. Allowed methods are:
mean: average of high and low voltage.
{n}*std_low: n standard deviations above mean low voltage, where n is a float (ignore slash in raw docstring).
{n}*std_high: n standard deviations below mean high voltage, where n is a float (ignore slash in raw docstring).
config: Use threshold method provided in
config.analysis.threshold_method
(mean
if not specified)
plot (
Union
[bool
,Axis
]) – Whether to plot traces with results. If True, will create a MatPlot object and add results. Can also pass a MatPlot axis, in which case that will be used. Each trace is preceded by a block that can be green (measured blip during start), red (no blip measured), or white (trace was filtered out).
- Returns
up_proportion (float): proportion of traces that has a blip
end_high (float): proportion of traces that end with high voltage
end_low (float): proportion of traces that end with low voltage
num_traces (int): Number of traces that satisfy filter
filtered_traces_idx (np.ndarray): 1D bool array, True if that trace satisfies filter
voltage_difference (float): voltage difference between high and low voltages
average_voltage (float): average voltage over all traces
threshold_voltage (float): threshold voltage for counting a blip (high voltage). Is calculated if not provided as input arg.
blips (float): average blips per trace.
mean_low_blip_duration (float): average duration in low state
mean_high_blip_duration (float): average duration in high state
- Return type
Dict[str, Any]
Note
- If no threshold voltage is provided, and no two peaks can be discerned,
all results except average_voltage are set to an initial value (either 0 or undefined)
- If the filtered trace proportion is less than min_filter_proportion,
the results
up_proportion
,end_low
,end_high
are set to an initial value
-
silq.analysis.analysis.
analyse_EPR
(empty_traces, plunge_traces, read_traces, sample_rate, t_skip, t_read, min_filter_proportion=0.5, threshold_voltage=None, filter_traces=True, plot=False)[source]¶ Analyse an empty-plunge-read sequence
- Parameters
empty_traces (
ndarray
) – 2D array of acquisition traces in empty (ionized) stateplunge_traces (
ndarray
) – 2D array of acquisition traces in plunge (neutral) stateread_traces (
ndarray
) – 2D array of acquisition traces in read statesample_rate (
float
) – acquisition sample rate (per second).t_skip (
float
) – initial time to skip for each trace (ms).t_read (
float
) – duration of each trace to use for calculating up_proportion etc. e.g. for a long trace, you want to compare up proportion of start and end segments.min_filter_proportion (
float
) – minimum proportion of traces that satisfy filter. If below this value, up_proportion etc. are not calculated.threshold_voltage (
Optional
[float
]) – threshold voltage for ahigh
voltage (blip). If not specified,find_high_low
is used to determine threshold.
- Returns
fidelity_empty (float): proportion of empty traces that end ionized (high voltage). Traces are filtered out that do not start neutral (low voltage).
voltage_difference_empty (float): voltage difference between high and low state in empty traces
fidelity_load (float): proportion of plunge traces that end neutral (low voltage). Traces are filtered out that do not start ionized (high voltage).
voltage_difference_load (float): voltage difference between high and low state in plunge traces
up_proportion (float): proportion of read traces that have blips For each trace, only up to t_read is considered.
dark_counts (float): proportion of read traces that have dark counts. For each trace, only the final t_read is considered.
contrast (float): =up_proportion - dark_counts
voltage_difference_read (float): voltage difference between high and low state in read traces
blips (float): average blips per trace in read traces.
mean_low_blip_duration (float): average duration in low state.
mean_high_blip_duration (float): average duration in high state.
- Return type
-
silq.analysis.analysis.
analyse_flips
(up_proportions_arrs, threshold_up_proportion)[source]¶ Analyse flipping between NMR states
For each up_proportion array, it will count the number of flips (transition between high/low up proportion) for each sample.
If more than one up_proportion array is passed, combined flipping events between each pair of states is also considered (i.e. one goes from low to high up proportion, the other from high to low). Furthermore, filtering is also performed where flips are only counted if the nuclear state remains in the subspace spanned by the pair of states for all samples.
- Parameters
up_proportions_arrs (
List
[ndarray
]) – 3D Arrays of up proportions, calculated fromanalyse_traces
. Up proportion is the proportion of traces that contain blips. First and second dimensions are arbitrary (can be a singleton), third dimension is samples.threshold_up_proportion (
Union
[Sequence
,float
]) – Threshold for when an up_proportion is high enough to count the nucleus as being in that state (consistently able to flip the electron). Can also be a tuple with two elements, in which case the first is a lower threshold, below which we can say the nucleus is not in the state, and the second is an upper threshold . If any up proportion is not in that state, it is in an undefined state, and not counted for flipping. For any filtering on up_proportion pairs, the whole row is considered invalid (NaN).
- Returns
flips(_{idx}) (int): Flips between high/low up_proportion. If more than one up_proportion_arr is provided, a zero-based index is added to specify the up_proportion_array. If an up_proportion is between the lower and upper threshold, it is not counted.
flip_probability(_{idx}) (float): proportion of flips compared to maximum flips (samples - 1). If more than one up_proportion_arr is provided, a zero-based index is added to specify the up_proportion_array.
The following results are between neighbouring pairs of up_proportion_arrs (idx1-idx2== +-1), and are only returned if more than one up_proportion_arr is given:
combined_flips_{idx1}{idx2}: combined flipping events between up_proportion_arrs idx1 and idx2, where one of the up_proportions switches from high to low, and the other from low to high.
combined_flip_probability_{idx1}{idx2}: Flipping probability of the combined flipping events.
Additionally, each of the above results will have another result with the same name, but prepended with
filtered_
, and appended with_{idx1}{idx2}
if not already present. Here, all the values are filtered out where the corresponding pair of up_proportion samples do not have exactly one high and one low for each sample. The values that do not satisfy the filter are set to np.nan.filtered_scans_{idx1}{idx2}: 2D bool array, True if pair of up_proportion rows remain in subspace
- Return type
Dict[str, Any]
silq.analysis.fit_toolbox module¶
-
class
silq.analysis.fit_toolbox.
AMSineFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
-
find_initial_parameters
(xvals, ydata, initial_parameters={}, plot=False)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals – x-coordinates of data points
ydata – data points
initial_parameters – Fixed initial parameters to be skipped.
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(t, amplitude, frequency, phase, offset, amplitude_AM, frequency_AM, phase_AM)[source]¶ Amplitude-Modulated Sinusoidal fit using time as x-coordinates
-
sweep_parameter
= 't'¶
-
-
class
silq.analysis.fit_toolbox.
BayesianUpdateFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a ‘Bayesian update’ function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(t, tau_1, tau_2, prior, offset)[source]¶ Exponential function using time as x-coordinate
-
sweep_parameter
= 't'¶
-
-
class
silq.analysis.fit_toolbox.
DoubleExponentialFit
(*args, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a double exponential function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(t, tau_1, A_1, A_2, tau_2, offset)[source]¶ Exponential function using time as x-coordinate
- Parameters
- Return type
- Returns
exponential data points
-
sweep_parameter
= 't'¶
-
-
class
silq.analysis.fit_toolbox.
DoubleFermiFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a double Fermi-Dirac distribution function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(V, A_1, A_2, U_1, U_2, T, alpha, offset)[source]¶ Exponential function using time as x-coordinate
- Parameters
- Return type
- Returns
exponential data points
-
kB
= 8.61733034e-05¶
-
sweep_parameter
= 'V'¶
-
-
class
silq.analysis.fit_toolbox.
ExponentialFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for an exponential function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(t, tau, amplitude, offset)[source]¶ Exponential function using time as x-coordinate
-
sweep_parameter
= 't'¶
-
-
class
silq.analysis.fit_toolbox.
ExponentialSineFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
-
find_initial_parameters
(xvals, ydata, initial_parameters={}, plot=False)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals – x-coordinates of data points
ydata – data points
initial_parameters – Fixed initial parameters to be skipped.
- Returns
Parameters object containing initial parameters.
-
sweep_parameter
= 't'¶
-
-
class
silq.analysis.fit_toolbox.
FermiFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a Fermi-Dirac distribution function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
kB
= 8.61733034e-05¶
-
sweep_parameter
= 'V'¶
-
-
class
silq.analysis.fit_toolbox.
Fit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
object
Base fitting class.
To fit a specific function, this class should be subclassed.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
add_to_plot
(ax, N=201, xrange=None, xscale=1, yscale=1, **kwargs)[source]¶ Add fit to existing plot axis
- Parameters
ax (
Axis
) – Axis to add plot toN (
int
) – number of points to use as x values (to smoothe fit curve)xrange (
Optional
[Tuple
[float
]]) – Optional range for x values (min, max)xscale (
float
) – value to multiple x values by to rescale axisyscale (
float
) – value to multiple y values by to rescale axiskwargs – Additional plot kwargs. By default Fit.plot_kwargs are used
- Returns
plot_handle of fit curve
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
find_nearest_index
(array, value)[source]¶ Find index in array that is closest to a value
- Parameters
array – array from which to find the nearest index
value – target value for which we want the index of the nearest array element
- Returns
Index of element in array that is closest to target value
-
get_parameters
(xvals, ydata, initial_parameters={}, fixed_parameters={}, parameter_constraints={}, weights=None, **kwargs)[source]¶ Get parameters for fitting :param xvals: x-coordinates of data points :param ydata: Data points :param initial_parameters: {parameter: initial_value} combination,
to specify the initial value of certain parameters. The initial values of other parameters are found using
Fit.find_initial_parameters
.- Parameters
fixed_parameters – {parameter: fixed_value} combination, to specify parameters whose values should not be varied.
parameter_constraints – {parameter: {constraint : value, …}, …} combination to further constrain existing parameters. e.g. {‘frequency’ : {‘min’ : 0}} ensures only positive frequencies can be fit.
weights – Weights for data points, must have same shape as ydata
-
perform_fit
(parameters=None, print=False, plot=None, **kwargs)[source]¶ Perform fitting routine
- Return type
ModelResult
- Returns
ModelResult object containing fitting results
-
plot_kwargs
= {'color': 'cyan', 'linestyle': '--', 'lw': 3}¶
-
sweep_parameter
= None¶
-
-
class
silq.analysis.fit_toolbox.
GaussianFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a Gaussian function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(x, amplitude, mean, sigma, offset)[source]¶ Gaussian function using x as x-coordinate
-
sweep_parameter
= 'x'¶
-
-
class
silq.analysis.fit_toolbox.
LinearFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a linear function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
sweep_parameter
= 't'¶
-
-
class
silq.analysis.fit_toolbox.
LorentzianFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a Lorentzian function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(x, amplitude, mean, gamma, offset)[source]¶ Lorentzian function using x as x-coordinate
-
sweep_parameter
= 'x'¶
-
-
class
silq.analysis.fit_toolbox.
RabiFrequencyFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
-
find_initial_parameters
(xvals, ydata, initial_parameters={}, plot=False)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals – x-coordinates of data points
ydata – data points
initial_parameters – Fixed initial parameters to be skipped.
- Returns
Parameters object containing initial parameters.
-
sweep_parameter
= 'f'¶
-
-
class
silq.analysis.fit_toolbox.
SineFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
-
find_initial_parameters
(xvals, ydata, initial_parameters={}, plot=False)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals – x-coordinates of data points
ydata – data points
initial_parameters – Fixed initial parameters to be skipped.
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(t, amplitude, frequency, phase, offset)[source]¶ Sinusoidal fit using time as x-coordinates
-
sweep_parameter
= 't'¶
-
-
class
silq.analysis.fit_toolbox.
T1Fit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a 1/T1 vs magnetic field (B) function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
sweep_parameter
= 'x'¶
-
-
class
silq.analysis.fit_toolbox.
VoigtFit
(fit=True, print=False, plot=None, **kwargs)[source]¶ Bases:
silq.analysis.fit_toolbox.Fit
Fitting class for a Voigt function.
To fit data to a function, use the method
Fit.perform_fit
. This will find its initial parameters viaFit.find_initial_parameters
, after which it will fit the data toFit.fit_function
.Note
The fitting routine uses lmfit, a wrapper package around scipy.optimize.
-
find_initial_parameters
(xvals, ydata, initial_parameters)[source]¶ Estimate initial parameters from data.
This is needed to ensure that the fitting will converge.
- Parameters
xvals (
ndarray
) – x-coordinates of data pointsydata (
ndarray
) – data pointsinitial_parameters (
dict
) – Fixed initial parameters to be skipped.
- Return type
Parameters
- Returns
Parameters object containing initial parameters.
-
static
fit_function
(x, amplitude, mean, gamma, sigma, offset)[source]¶ Gaussian function using x as x-coordinate
-
sweep_parameter
= 'x'¶
-
silq.analysis.tunneling_times module¶
-
silq.analysis.tunneling_times.
analyse_tunnel_times
(tunnel_times_arr, bins=50, range=None, silent=True, double_exponential=True, plot=False, fixed_parameters={}, initial_parameters={}, plot_fit_kwargs={}, **kwargs)[source]¶