qetpy.detcal package

Submodules

qetpy.detcal.cut module

qetpy.detcal.cut.symmetrizedist(vals)

Function to symmetrize a distribution about zero. Useful for if the distribution of some value centers around a nonzero value, but should center around zero. An example of this would be when most of the measured slopes are nonzero, but we want the slopes with zero values (e.g. lots of muon tails, which we want to cut out). To do this, the algorithm randomly chooses points in a histogram to cut out until the histogram is symmetric about zero.

Parameters:
vals : ndarray

A 1-d array of the values that will be symmetrized.

Returns:
czeromeanslope : ndarray

A boolean mask of the values that should be kept.

qetpy.detcal.cut.pileupcut(traces, fs=625000.0, outlieralgo='removeoutliers', nsig=2, removemeans=False)

Function to automatically cut out outliers of the optimum filter amplitudes of the inputted traces.

Parameters:
traces : ndarray

2-dimensional array of traces to do cuts on

fs : float, optional

Digitization rate that the data was taken at

outlieralgo : string, optional

Which outlier algorithm to use. If set to “removeoutliers”, uses the removeoutliers algorithm that removes data based on the skewness of the dataset. If set to “iterstat”, uses the iterstat algorithm to remove data based on being outside a certain number of standard deviations from the mean

nsig : float, optional

If outlieralgo is “iterstat”, this can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the optimum filter amplitudes. Default is 2.

removemeans : boolean, optional

Boolean flag on if the mean of each trace should be removed before doing the optimal filter (True) or if the means should not be removed (False). This is useful for dIdV traces, when we want to cut out pulses that have smaller amplitude than the dIdV overshoot. Default is False.

Returns:
cpileup : ndarray

Boolean array giving which indices to keep or throw out based on the outlier algorithm

qetpy.detcal.cut.slopecut(traces, fs=625000.0, outlieralgo='removeoutliers', nsig=2, is_didv=False, symmetrizeflag=False, sgfreq=100.0)

Function to automatically cut out outliers of the slopes of the inputted traces. Includes a routine that attempts to symmetrize the distribution of slopes around zero, which is useful when the majority of traces have a slope.

Parameters:
traces : ndarray

2-dimensional array of traces to do cuts on

fs : float, optional

Digitization rate that the data was taken at

outlieralgo : string, optional

Which outlier algorithm to use. If set to “removeoutliers”, uses the removeoutliers algorithm that removes data based on the skewness of the dataset. If set to “iterstat”, uses the iterstat algorithm to remove data based on being outside a certain number of standard deviations from the mean

nsig : float, optional

If outlieralgo is “iterstat”, this can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the slopes. Default is 2.

is_didv : bool, optional

Boolean flag on whether or not the trace is a dIdV curve

symmetrizeflag : bool, optional

Flag for whether or not the slopes should be forced to have an average value of zero. Should be used if most of the traces have a slope

sgfreq : float, optional

If is_didv is True, then the sgfreq is used to know where the flat parts of the traces should be

Returns:
cslope : ndarray

Boolean array giving which indices to keep or throw out based on the outlier algorithm

qetpy.detcal.cut.baselinecut(traces, fs=625000.0, outlieralgo='removeoutliers', nsig=2, is_didv=False, sgfreq=100.0)

Function to automatically cut out outliers of the baselines of the inputted traces.

Parameters:
traces : ndarray

2-dimensional array of traces to do cuts on

fs : float, optional

Digitization rate that the data was taken at

outlieralgo : string, optional

Which outlier algorithm to use. If set to “removeoutliers”, uses the removeoutliers algorithm that removes data based on the skewness of the dataset. If set to “iterstat”, uses the iterstat algorithm to remove data based on being outside a certain number of standard deviations from the mean

nsig : float, optional

If outlieralgo is “iterstat”, this can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the baselines. Default is 2.

is_didv : bool, optional

Boolean flag on whether or not the trace is a dIdV curve

sgfreq : float, optional

If is_didv is True, then the sgfreq is used to know where the flat parts of the traces should be

Returns:
cbaseline : ndarray

Boolean array giving which indices to keep or throw out based on the outlier algorithm

qetpy.detcal.cut.chi2cut(traces, fs=625000.0, outlieralgo='iterstat', nsig=2)

Function to automatically cut out outliers of the baselines of the inputted traces.

Parameters:
traces : ndarray

2-dimensional array of traces to do cuts on

fs : float, optional

Digitization rate that the data was taken at

outlieralgo : string, optional

Which outlier algorithm to use. If set to “removeoutliers”, uses the removeoutliers algorithm that removes data based on the skewness of the dataset. If set to “iterstat”, uses the iterstat algorithm to remove data based on being outside a certain number of standard deviations from the mean

nsig : float, optional

If outlieralgo is “iterstat”, this can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the Chi2s. Default is 2.

Returns:
cchi2 : ndarray

Boolean array giving which indices to keep or throw out based on the outlier algorithm

qetpy.detcal.cut.autocuts(traces, fs=625000.0, is_didv=False, sgfreq=200.0, symmetrizeflag=False, outlieralgo='removeoutliers', lgcpileup1=True, lgcslope=True, lgcbaseline=True, lgcpileup2=True, lgcchi2=True, nsigpileup1=2, nsigslope=2, nsigbaseline=2, nsigpileup2=2, nsigchi2=3)

Function to automatically cut out bad traces based on the optimum filter amplitude, slope, baseline, and chi^2 of the traces.

Parameters:
traces : ndarray

2-dimensional array of traces to do cuts on

fs : float, optional

Sample rate that the data was taken at

is_didv : bool, optional

Boolean flag on whether or not the trace is a dIdV curve

sgfreq : float, optional

If is_didv is True, then the sgfreq is used to know where the flat parts of the traces should be

symmetrizeflag : bool, optional

Flag for whether or not the slopes should be forced to have an average value of zero. Should be used if most of the traces have a slope

outlieralgo : string, optional

Which outlier algorithm to use. If set to “removeoutliers”, uses the removeoutliers algorithm that removes data based on the skewness of the dataset. If set to “iterstat”, uses the iterstat algorithm to remove data based on being outside a certain number of standard deviations from the mean

lgcpileup1 : boolean, optional

Boolean value on whether or not do the pileup1 cut (this is the initial pileup cut that is always done whether or not we have dIdV data). Default is True.

lgcslope : boolean, optional

Boolean value on whether or not do the slope cut. Default is True.

lgcbaseline : boolean, optional

Boolean value on whether or not do the baseline cut. Default is True.

lgcpileup2 : boolean, optional

Boolean value on whether or not do the pileup2 cut (this cut is only done when is_didv is also True). Default is True.

lgcchi2 : boolean, optional

Boolean value on whether or not do the chi2 cut. Default is True.

nsigpileup1 : float, optional

If outlieralgo is “iterstat”, this can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the optimum filter amplitudes. Default is 2.

nsigslope : float, optional

If outlieralgo is “iterstat”, this can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the slopes. Default is 2.

nsigbaseline : float, optional

If outlieralgo is “iterstat”, this can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the baselines. Default is 2.

nsigpileup2 : float, optional

If outlieralgo is “iterstat”, this can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the optimum filter amplitudes after the mean has been subtracted. (only used if is_didv is True). Default is 2.

nsigchi2 : float, optional

This can be used to tune the number of standard deviations from the mean to cut outliers from the data when using iterstat on the chi^2 values. Default is 3. This is always used, as iterstat is always used for the chi^2 cut.

Returns:
ctot : ndarray

Boolean array giving which indices to keep or throw out based on the autocuts algorithm

qetpy.detcal.didv module

qetpy.detcal.didv.didvinitfromdata(tmean, didvmean, didvstd, offset, offset_err, fs, sgfreq, sgamp, rshunt, r0=0.3, r0_err=0.001, rload=0.01, rload_err=0.001, priors=None, invpriorscov=None, add180phase=False, dt0=1e-05)

Function to initialize and process a dIdV dataset without having all of the traces, but just the parameters that are required for fitting. After running, this returns a DIDV class object that is ready for fitting.

Parameters:
tmean : ndarray

The average trace in time domain, units of Amps

didvstd : ndarray

The complex standard deviation of the didv in frequency space for each frequency

didvmean : ndarray

The average trace converted to didv

offset : float

The offset (i.e. baseline value) of the didv trace, in Amps

offset_err : float

The error in the offset of the didv trace, in Amps

fs : float

Sample rate of the data taken, in Hz

sgfreq : float

Frequency of the signal generator, in Hz

sgamp : float

Amplitude of the signal generator, in Amps (equivalent to jitter in the QET bias)

rshunt : float

Shunt resistance in the circuit, Ohms

r0 : float, optional

Resistance of the TES in Ohms

r0_err : float, optional

Error in the resistance of the TES (Ohms)

rload : float, optional

Load resistance of the circuit (rload = rshunt + rparasitic), Ohms

rload_err : float, optional

Error in the load resistance, Ohms

priors : ndarray, optional

Prior known values of Irwin’s TES parameters for the trace. Should be in the order of (rload,r0,beta,l,L,tau0,dt)

invpriorscov : ndarray, optional

Inverse of the covariance matrix of the prior known values of Irwin’s TES parameters for the trace (any values that are set to zero mean that we have no knowledge of that parameter)

add180phase : boolean, optional

If the signal generator is out of phase (i.e. if it looks like –__ instead of __–), then this should be set to True. Adds half a period of the signal generator to the dt0 attribute

dt0 : float, optional

The value of the starting guess for the time offset of the didv when fitting. The best way to use this value if it isn’t converging well is to run the fit multiple times, setting dt0 equal to the fit’s next value, and seeing where the dt0 value converges. The fit can have a difficult time finding the value on the first run if it the initial value is far from the actual value, so a solution is to do this iteratively.

Returns:
didvobj : Object

A DIDV class object that can be used to fit the dIdV and return the fit parameters.

qetpy.detcal.didv.onepoleimpedance(freq, A, tau2)

Function to calculate the impedance (dvdi) of a TES with the 1-pole fit

Parameters:
freq : array_like, float

The frequencies for which to calculate the admittance (in Hz)

A : float

The fit parameter A in the complex impedance (in Ohms), superconducting: A=rload, normal: A = rload+Rn

tau2 : float

The fit parameter tau2 in the complex impedance (in s), superconducting: tau2=L/rload, normal: tau2=L/(rload+Rn)

Returns:
dvdi : array_like, float

The complex impedance of the TES with the 1-pole fit

qetpy.detcal.didv.onepoleadmittance(freq, A, tau2)

Function to calculate the admittance (didv) of a TES with the 1-pole fit

Parameters:
freq : array_like, float

The frequencies for which to calculate the admittance (in Hz)

A : float

The fit parameter A in the complex impedance (in Ohms), superconducting: A=rload, normal: A = rload+Rn

tau2 : float

The fit parameter tau2 in the complex impedance (in s), superconducting: tau2=L/rload, normal: tau2=L/(rload+Rn)

Returns:
1.0/dvdi : array_like, float

The complex admittance of the TES with the 1-pole fit

qetpy.detcal.didv.twopoleimpedance(freq, A, B, tau1, tau2)

Function to calculate the impedance (dvdi) of a TES with the 2-pole fit

Parameters:
freq : array_like, float

The frequencies for which to calculate the admittance (in Hz)

A : float

The fit parameter A in the complex impedance (in Ohms), A = rload + r0*(1+beta)

B : float

The fit parameter B in the complex impedance (in Ohms), B = r0*l*(2+beta)/(1-l) (where l is Irwin’s loop gain)

tau1 : float

The fit parameter tau1 in the complex impedance (in s), tau1=tau0/(1-l)

tau2 : float

The fit parameter tau2 in the complex impedance (in s), tau2=L/(rload+r0*(1+beta))

Returns:
dvdi : array_like, float

The complex impedance of the TES with the 2-pole fit

qetpy.detcal.didv.twopoleadmittance(freq, A, B, tau1, tau2)

Function to calculate the admittance (didv) of a TES with the 2-pole fit

Parameters:
freq : array_like, float

The frequencies for which to calculate the admittance (in Hz)

A : float

The fit parameter A in the complex impedance (in Ohms), A = rload + r0*(1+beta)

B : float

The fit parameter B in the complex impedance (in Ohms), B = r0*l*(2+beta)/(1-l) (where l is Irwin’s loop gain)

tau1 : float

The fit parameter tau1 in the complex impedance (in s), tau1=tau0/(1-l)

tau2 : float

The fit parameter tau2 in the complex impedance (in s), tau2=L/(rload+r0*(1+beta))

Returns:
1.0/dvdi : array_like, float

The complex admittance of the TES with the 2-pole fit

qetpy.detcal.didv.threepoleimpedance(freq, A, B, C, tau1, tau2, tau3)

Function to calculate the impedance (dvdi) of a TES with the 3-pole fit

Parameters:
freq : array_like, float

The frequencies for which to calculate the admittance (in Hz)

A : float

The fit parameter A in the complex impedance (in Ohms)

B : float

The fit parameter B in the complex impedance (in Ohms)

C : float

The fit parameter C in the complex impedance

tau1 : float

The fit parameter tau1 in the complex impedance (in s)

tau2 : float

The fit parameter tau2 in the complex impedance (in s)

tau3 : float

The fit parameter tau3 in the complex impedance (in s)

Returns:
dvdi : array_like, float

The complex impedance of the TES with the 3-pole fit

qetpy.detcal.didv.threepoleadmittance(freq, A, B, C, tau1, tau2, tau3)

Function to calculate the admittance (didv) of a TES with the 3-pole fit

Parameters:
freq : array_like, float

The frequencies for which to calculate the admittance (in Hz)

A : float

The fit parameter A in the complex impedance (in Ohms)

B : float

The fit parameter B in the complex impedance (in Ohms)

C : float

The fit parameter C in the complex impedance

tau1 : float

The fit parameter tau1 in the complex impedance (in s)

tau2 : float

The fit parameter tau2 in the complex impedance (in s)

tau3 : float

The fit parameter tau3 in the complex impedance (in s)

Returns:
1.0/dvdi : array_like, float

The complex admittance of the TES with the 3-pole fit

qetpy.detcal.didv.twopoleimpedancepriors(freq, rload, r0, beta, l, L, tau0)

Function to calculate the impedance (dvdi) of a TES with the 2-pole fit from Irwin’s TES parameters

Parameters:
freq : array_like, float

The frequencies for which to calculate the admittance (in Hz)

rload : float

The load resistance of the TES (in Ohms)

r0 : float

The resistance of the TES (in Ohms)

beta : float

The current sensitivity of the TES

l : float

Irwin’s loop gain

L : float

The inductance in the TES circuit (in Henrys)

tau0 : float

The thermal time constant of the TES (in s)

Returns:
dvdi : array_like, float

The complex impedance of the TES with the 2-pole fit from Irwin’s TES parameters

qetpy.detcal.didv.twopoleadmittancepriors(freq, rload, r0, beta, l, L, tau0)

Function to calculate the admittance (didv) of a TES with the 2-pole fit from Irwin’s TES parameters

Parameters:
freq : array_like, float

The frequencies for which to calculate the admittance (in Hz)

rload : float

The load resistance of the TES (in Ohms)

r0 : float

The resistance of the TES (in Ohms)

beta : float

The current sensitivity of the TES, beta=d(log R)/d(log I)

l : float

Irwin’s loop gain, l = P0*alpha/(G*Tc)

L : float

The inductance in the TES circuit (in Henrys)

tau0 : float

The thermal time constant of the TES (in s), tau0=C/G

Returns:
1.0/dvdi : array_like, float

The complex admittance of the TES with the 2-pole fit from Irwin’s TES parameters

qetpy.detcal.didv.convolvedidv(x, A, B, C, tau1, tau2, tau3, sgamp, rshunt, sgfreq, dutycycle)

Function to convert the fitted TES parameters for the complex impedance to a TES response to a square wave jitter in time domain.

Parameters:
x : array_like

Time values for the trace (in s)

A : float

The fit parameter A in the complex impedance (in Ohms)

B : float

The fit parameter B in the complex impedance (in Ohms)

C : float

The fit parameter C in the complex impedance

tau1 : float

The fit parameter tau1 in the complex impedance (in s)

tau2 : float

The fit parameter tau2 in the complex impedance (in s)

tau3 : float

The fit parameter tau3 in the complex impedance (in s)

sgamp : float

The peak-to-peak size of the square wave jitter (in Amps)

rshunt : float

The shunt resistance of the TES electronics (in Ohms)

sgfreq : float

The frequency of the square wave jitter (in Hz)

dutycycle : float

The duty cycle of the square wave jitter (between 0 and 1)

Returns:
np.real(st) : ndarray

The response of a TES to a square wave jitter in time domain with the given fit parameters. The real part is taken in order to ensure that the trace is real

qetpy.detcal.didv.squarewaveguessparams(trace, sgamp, rshunt)

Function to guess the fit parameters for the 1-pole fit.

Parameters:
trace : array_like

The trace in time domain (in Amps).

sgamp : float

The peak-to-peak size of the square wave jitter (in Amps)

rshunt : float

Shunt resistance of the TES electronics (in Ohms)

Returns:
A0 : float

Guess of the fit parameter A (in Ohms)

tau20 : float

Guess of the fit parameter tau2 (in s)

qetpy.detcal.didv.guessdidvparams(trace, flatpts, sgamp, rshunt, L0=1e-07)

Function to find the fit parameters for either the 1-pole (A, tau2, dt), 2-pole (A, B, tau1, tau2, dt), or 3-pole (A, B, C, tau1, tau2, tau3, dt) fit.

Parameters:
trace : array_like

The trace in time domain (in Amps)

flatpts : array_like

The flat parts of the trace (in Amps)

sgamp : float

The peak-to-peak size of the square wave jitter (in Amps)

rshunt : float

Shunt resistance of the TES electronics (in Ohms)

L0 : float, optional

The guess of the inductance (in Henries)

Returns:
A0 : float

Guess of the fit parameter A (in Ohms)

B0 : float

Guess of the fit parameter B (in Ohms)

tau10 : float

Guess of the fit parameter tau1 (in s)

tau20 : float

Guess of the fit parameter tau2 (in s)

isloopgainsub1 : boolean

Boolean flag that gives whether the loop gain is greater than one (False) or less than one (True)

qetpy.detcal.didv.fitdidv(freq, didv, yerr=None, A0=0.25, B0=-0.6, C0=-0.6, tau10=-0.0003183098861837907, tau20=1.5915494309189535e-06, tau30=0.0, dt=-1e-05, poles=2, isloopgainsub1=None)

Function to find the fit parameters for either the 1-pole (A, tau2, dt), 2-pole (A, B, tau1, tau2, dt), or 3-pole (A, B, C, tau1, tau2, tau3, dt) fit.

Parameters:
freq : ndarray

Frequencies corresponding to the didv

didv : ndarray

Complex impedance extracted from the trace in frequency space

yerr : ndarray, NoneType, optional

Error at each frequency of the didv. Should be a complex number, e.g. yerr = yerr_real + 1.0j * yerr_imag, where yerr_real is the standard deviation of the real part of the didv, and yerr_imag is the standard deviation of the imaginary part of the didv. If left as None, then each frequency will be assumed to be equally weighted.

A0 : float, optional

Guess of the fit parameter A (in Ohms). Default is 0.25.

B0 : float, optional

Guess of the fit parameter B (in Ohms). Default is -0.6.

C0 : float, optional

Guess of the fit parameter C (unitless). Default is -0.6.

tau10 : float, optional

Guess of the fit parameter tau1 (in s). Default is -1.0/(2*pi*5e2).

tau20 : float, optional

Guess of the fit parameter tau2 (in s). Default is 1.0/(2*pi*1e5).

tau30 : float, optional

Guess of the fit parameter tau3 (in s). Default is 0.0.

dt : float, optional

Guess of the time shift (in s). Default is -10.0e-6.

poles : int, optional

The number of poles to use in the fit (should be 1, 2, or 3). Default is 2.

isloopgainsub1 : boolean, NoneType, optional

If set, should be used to specify if the fit should be done assuming that the Irwin loop gain is less than 1 (True) or greater than 1 (False). Default is None, in which case loop gain less than 1 and greater than 1 fits will be done, returning the one with a lower Chi^2.

Returns:
popt : ndarray

The fitted parameters for the specificed number of poles

pcov : ndarray

The corresponding covariance matrix for the fitted parameters

cost : float

The cost of the the fit

qetpy.detcal.didv.converttotesvalues(popt, pcov, r0, rload, r0_err=0.001, rload_err=0.001)

Function to convert the fit parameters for either 1-pole (A, tau2, dt), 2-pole (A, B, tau1, tau2, dt), or 3-pole (A, B, C, tau1, tau2, tau3, dt) fit to the corresponding TES parameters: 1-pole (rtot, L, r0, rload, dt), 2-pole (rload, r0, beta, l, L, tau0, dt), and 3-pole (no conversion done).

Parameters:
popt : ndarray

The fit parameters for either the 1-pole, 2-pole, or 3-pole fit

pcov : ndarray

The corresponding covariance matrix for the fit parameters

r0 : float

The resistance of the TES (in Ohms)

rload : float

The load resistance of the TES circuit (in Ohms)

r0_err : float, optional

The error in the r0 value (in Ohms). Default is 0.001.

rload_err : float, optional

The error in the rload value (in Ohms). Default is 0.001.

Returns:
popt_out : ndarray

The TES parameters for the specified fit

pcov_out : ndarray

The corresponding covariance matrix for the TES parameters

qetpy.detcal.didv.fitdidvpriors(freq, didv, priors, invpriorscov, yerr=None, rload=0.35, r0=0.13, beta=0.5, l=10.0, L=5e-07, tau0=0.0005, dt=-1e-05)

Function to directly fit Irwin’s TES parameters (rload, r0, beta, l, L, tau0, dt) with the knowledge of prior known values any number of the parameters. In order for the degeneracy of the parameters to be broken, at least 2 fit parameters should have priors knowledge. This is usually rload and r0, as these can be known from IV data.

Parameters:
freq : ndarray

Frequencies corresponding to the didv

didv : ndarray

Complex impedance extracted from the trace in frequency space

priors : ndarray

Prior known values of Irwin’s TES parameters for the trace. Should be in the order of (rload,r0,beta,l,L,tau0,dt)

invpriorscov : ndarray

Inverse of the covariance matrix of the prior known values of Irwin’s TES parameters for the trace (any values that are set to zero mean that we have no knowledge of that parameter)

yerr : ndarray, optional

Error at each frequency of the didv. Should be a complex number, e.g. yerr = yerr_real + 1.0j * yerr_imag, where yerr_real is the standard deviation of the real part of the didv, and yerr_imag is the standard deviation of the imaginary part of the didv. If left as None, then each frequency will be assumed to be equally weighted.

rload : float, optional

Guess of the load resistance of the TES circuit (in Ohms). Default is 0.35.

r0 : float, optional

Guess of the resistance of the TES (in Ohms). Default is 0.130.

beta : float, optional

Guess of the current sensitivity beta (unitless). Default is 0.5.

l : float, optional

Guess of Irwin’s loop gain (unitless). Default is 10.0.

L : float, optional

Guess of the inductance (in Henrys). Default is 500.0e-9.

tau0 : float, optional

Guess of the thermal time constant (in s). Default is 500.0e-6.

dt : float, optional

Guess of the time shift (in s). Default is -10.0e-6.

Returns:
popt : ndarray

The fitted parameters in the order of (rload, r0, beta, l, L, tau0, dt)

pcov : ndarray

The corresponding covariance matrix for the fitted parameters

cost : float

The cost of the the fit

qetpy.detcal.didv.convertfromtesvalues(popt, pcov)

Function to convert from Irwin’s TES parameters (rload, r0, beta, l, L, tau0, dt) to the fit parameters (A, B, tau1, tau2, dt)

Parameters:
popt : ndarray

Irwin’s TES parameters in the order of (rload, r0, beta, l, L, tau0, dt), should be a 1-dimensional np.array of length 7

pcov : ndarray

The corresponding covariance matrix for Irwin’s TES parameters. Should be a 2-dimensional, 7-by-7 np.array

Returns:
popt_out : ndarray

The fit parameters in the order of (A, B, tau1, tau2, dt)

pcov_out : ndarray

The corresponding covariance matrix for the fit parameters

qetpy.detcal.didv.findpolefalltimes(params)

Function for taking TES params from a 1-pole, 2-pole, or 3-pole didv and calculating the falltimes (i.e. the values of the poles in the complex plane)

Parameters:
params : ndarray

TES parameters for either 1-pole, 2-pole, or 3-pole didv. This will be a 1-dimensional np.array of varying length, depending on the fit. 1-pole fit has 3 parameters (A,tau2,dt), 2-pole fit has 5 parameters (A,B,tau1,tau2,dt), and 3-pole fit has 7 parameters (A,B,C,tau1,tau2,tau3,dt). The parameters should be in that order, and any other number of parameters will print a warning and return zero.

Returns:
np.sort(falltimes) : ndarray

The falltimes for the didv fit, sorted from fastest to slowest.

qetpy.detcal.didv.deconvolvedidv(x, trace, rshunt, sgamp, sgfreq, dutycycle)

Function for taking a trace with a known square wave jitter and extracting the complex impedance via deconvolution of the square wave and the TES response in frequency space.

Parameters:
x : ndarray

Time values for the trace

trace : ndarray

The trace in time domain (in Amps)

rshunt : float

Shunt resistance for electronics (in Ohms)

sgamp : float

Peak to peak value of square wave jitter (in Amps, jitter in QET bias)

sgfreq : float

Frequency of square wave jitter

dutycycle : float

duty cycle of square wave jitter

Returns:
freq : ndarray

The frequencies that each point of the trace corresponds to

didv : ndarray

Complex impedance of the trace in frequency space

zeroinds : ndarray

Indices of the frequencies where the trace’s Fourier Transform is zero. Since we divide by the FT of the trace, we need to know which values should be zero, so that we can ignore these points in the complex impedance.

class qetpy.detcal.didv.DIDV(rawtraces, fs, sgfreq, sgamp, rshunt, tracegain=1.0, r0=0.3, r0_err=0.001, rload=0.01, rload_err=0.001, dutycycle=0.5, add180phase=False, priors=None, invpriorscov=None, dt0=1e-05)

Bases: object

Class for fitting a didv curve for different types of models of the didv. Also gives various other useful values pertaining to the didv. This class supports doing 1, 2, and 3 pole fits, as well as a 2 pole priors fit. This is supported in a way that does one dataset at a time.

Attributes:
rawtraces : ndarray

The array of rawtraces to use when fitting the didv. Should be of shape (number of traces, length of trace in bins). This can be any units, as long as tracegain will convert this to Amps.

fs : float

Sample rate of the data taken, in Hz

sgfreq : float

Frequency of the signal generator, in Hz

sgamp : float

Amplitude of the signal generator, in Amps (equivalent to jitter in the QET bias)

r0 : float

Resistance of the TES in Ohms

r0_err : float

Error in the resistance of the TES (Ohms)

rload : float

Load resistance of the circuit (rload = rshunt + rparasitic), Ohms

rload_err : float

Error in the load resistance, Ohms

rshunt : float

Shunt resistance in the circuit, Ohms

tracegain : float

The factor that the rawtraces should be divided by to convert the units to Amps. If rawtraces already has units of Amps, then this should be set to 1.0

dutycycle : float

The duty cycle of the signal generator, should be a float between 0 and 1. Set to 0.5 by default

add180phase : boolean

If the signal generator is out of phase (i.e. if it looks like –__ instead of __–), then this should be set to True. Adds half a period of the signal generator to the dt0 attribute

priors : ndarray

Prior known values of Irwin’s TES parameters for the trace. Should be in the order of (rload,r0,beta,l,L,tau0,dt)

invpriorscov : ndarray

Inverse of the covariance matrix of the prior known values of Irwin’s TES parameters for the trace (any values that are set to zero mean that we have no knowledge of that parameter)

dt0 : float

The value of the starting guess for the time offset of the didv when fitting. The best way to use this value if it isn’t converging well is to run the fit multiple times, setting dt0 equal to the fit’s next value, and seeing where the dt0 value converges. The fit can have a difficult time finding the value on the first run if it the initial value is far from the actual value, so a solution is to do this iteratively.

freq : ndarray

The frequencies of the didv fit

time : ndarray

The times the didv trace

ntraces : float

The number of traces in the data

traces : ndarray

The traces being used in units of Amps and also truncated so as to include only an integer number of signal generator periods

flatinds : ndarray

The indices where the traces are flat

tmean : ndarray

The average trace in time domain, units of Amps

zeroinds : ndarray

The indices of the didv fit in frequency space where the values should be zero

didvstd : ndarray

The complex standard deviation of the didv in frequency space for each frequency

didvmean : ndarray

The average trace converted to didv

offset : float

The offset (i.e. baseline value) of the didv trace, in Amps

offset_err : float

The error in the offset of the didv trace, in Amps

fitparams1 : ndarray

The fit parameters of the 1-pole fit, in order of (A, tau2, dt)

fitcov1 : ndarray

The corresponding covariance for the 1-pole fit parameters

fitcost1 : float

The cost of the 1-pole fit

irwinparams1 : ndarray

The Irwin parameters of the 1-pole fit, in order of (rtot, L , r0, rload, dt)

irwincov1 : ndarray

The corresponding covariance for the Irwin parameters for the 1-pole fit

falltimes1 : ndarray

The fall times of the 1-pole fit, same as tau2, in s

didvfit1_timedomain : ndarray

The 1-pole fit in time domain

didvfit1_freqdomain : ndarray

The 1-pole fit in frequency domain

fitparams2 : ndarray

The fit parameters of the 2-pole fit, in order of (A, B, tau1, tau2, dt)

fitcov2 : ndarray

The corresponding covariance for the 2-pole fit parameters

fitcost2 : float

The cost of the 2-pole fit

irwinparams2 : ndarray

The Irwin parameters of the 2-pole fit, in order of (rload, r0, beta, l, L, tau0, dt)

irwincov2 : ndarray

The corresponding covariance for the Irwin parameters for the 2-pole fit

falltimes2 : ndarray

The fall times of the 2-pole fit, tau_plus and tau_minus, in s

didvfit2_timedomain : ndarray

The 2-pole fit in time domain

didvfit2_freqdomain : ndarray

The 2-pole fit in frequency domain

fitparams3 : ndarray

The fit parameters of the 3-pole fit, in order of (A, B, C, tau1, tau2, tau3, dt)

fitcov3 : ndarray

The corresponding covariance for the 3-pole fit parameters

fitcost3 : float

The cost of the 3-pole fit

irwinparams3 : NoneType

The Irwin parameters of the 3-pole fit, this returns None now, as there is no model that we convert to

irwincov3 : NoneType

The corresponding covariance for the Irwin parameters for the 3-pole fit, also returns None

falltimes3 : ndarray

The fall times of the 3-pole fit in s

didvfit3_timedomain : ndarray

The 3-pole fit in time domain

didvfit3_freqdomain : ndarray

The 3-pole fit in frequency domain

fitparams2priors : ndarray

The fit parameters of the 2-pole priors fit, in order of (A, B, tau1, tau2, dt), converted from the Irwin parameters

fitcov2priors : ndarray

The corresponding covariance for the 2-pole priors fit parameters

fitcost2priors : float

The cost of the 2-pole priors fit

irwinparams2priors : ndarray

The Irwin parameters of the 2-pole priors fit, in order of (rload, r0, beta, l, L, tau0, dt)

irwincov2priors : ndarray

The corresponding covariance for the Irwin parameters for the 2-pole priors fit

falltimes2priors : ndarray

The fall times of the 2-pole priors fit, tau_plus and tau_minus, in s

didvfit2priors_timedomain : ndarray

The 2-pole priors fit in time domain

didvfit2priors_freqdomain : ndarray

The 2-pole priors fit in frequency domain

Methods

doallfits() This module does all of the fits consecutively.
dofit(poles) This method does the fit that is specified by the variable poles.
dopriorsfit() This module runs the priorsfit, assuming that the priors and invpriorscov attributes have been set to the proper values.
get_irwinparams_dict(poles[, lgcpriors]) Returns a dictionary with the irwin fit parameters for a given number of poles
plot_didv_flipped([poles, plotpriors, …]) Module to plot the flipped trace in time domain.
plot_full_trace([poles, plotpriors, …]) Module to plot the entire trace in time domain
plot_re_im_didv([poles, plotpriors, …]) Module to plot the real and imaginary parts of the didv in frequency space.
plot_single_period_of_trace([poles, …]) Module to plot a single period of the trace in time domain
plot_zoomed_in_trace([poles, zoomfactor, …]) Module to plot a zoomed in portion of the trace in time domain.
processtraces() This method processes the traces loaded to the DIDV class object.
doallfits()

This module does all of the fits consecutively. The priors fit is not done if the attributes priors and invpriorscov have not yet been set.

dofit(poles)

This method does the fit that is specified by the variable poles. If the processtraces module has not been run yet, then this module will run that first. This module does not do the priors fit.

Parameters:
poles : int

The fit that should be run. Should be 1, 2, or 3.

dopriorsfit()

This module runs the priorsfit, assuming that the priors and invpriorscov attributes have been set to the proper values.

get_irwinparams_dict(poles, lgcpriors=False)

Returns a dictionary with the irwin fit parameters for a given number of poles

Parameters:
poles: int

The number of poles used for the fit

lgcpriors: bool, optional

If true, the values from the priors fit are returned

Returns:
return_dict: dictionary

The irwim parameters stored in a dictionary

plot_didv_flipped(poles='all', plotpriors=True, lgcsave=False, savepath='', savename='')

Module to plot the flipped trace in time domain. This function should be used to test if there are nonlinearities in the didv

Parameters:
poles : int, string, array_like, optional

The pole fits that we want to plot. If set to “all”, then plots all of the fits. Can also be set to just one of the fits. Can be set as an array of different fits, e.g. [1, 2]

plotpriors : boolean, optional

Boolean value on whether or not the priors fit should be plotted.

lgcsave : boolean, optional

Boolean value on whether or not the figure should be saved

savepath : string, optional

Where the figure should be saved. Saved in the current directory by default.

savename : string, optional

A string to append to the end of the file name if saving. Empty string by default.

plot_full_trace(poles='all', plotpriors=True, lgcsave=False, savepath='', savename='')

Module to plot the entire trace in time domain

Parameters:
poles : int, string, array_like, optional

The pole fits that we want to plot. If set to “all”, then plots all of the fits. Can also be set to just one of the fits. Can be set as an array of different fits, e.g. [1, 2]

plotpriors : boolean, optional

Boolean value on whether or not the priors fit should be plotted.

lgcsave : boolean, optional

Boolean value on whether or not the figure should be saved

savepath : string, optional

Where the figure should be saved. Saved in the current directory by default.

savename : string, optional

A string to append to the end of the file name if saving. Empty string by default.

plot_re_im_didv(poles='all', plotpriors=True, lgcsave=False, savepath='', savename='')

Module to plot the real and imaginary parts of the didv in frequency space. Currently creates two different plots.

Parameters:
poles : int, string, array_like, optional

The pole fits that we want to plot. If set to “all”, then plots all of the fits. Can also be set to just one of the fits. Can be set as an array of different fits, e.g. [1, 2]

plotpriors : boolean, optional

Boolean value on whether or not the priors fit should be plotted.

lgcsave : boolean, optional

Boolean value on whether or not the figure should be saved

savepath : string, optional

Where the figure should be saved. Saved in the current directory by default.

savename : string, optional

A string to append to the end of the file name if saving. Empty string by default.

plot_single_period_of_trace(poles='all', plotpriors=True, lgcsave=False, savepath='', savename='')

Module to plot a single period of the trace in time domain

Parameters:
poles : int, string, array_like, optional

The pole fits that we want to plot. If set to “all”, then plots all of the fits. Can also be set to just one of the fits. Can be set as an array of different fits, e.g. [1, 2]

plotpriors : boolean, optional

Boolean value on whether or not the priors fit should be plotted.

lgcsave : boolean, optional

Boolean value on whether or not the figure should be saved

savepath : string, optional

Where the figure should be saved. Saved in the current directory by default.

savename : string, optional

A string to append to the end of the file name if saving. Empty string by default.

plot_zoomed_in_trace(poles='all', zoomfactor=0.1, plotpriors=True, lgcsave=False, savepath='', savename='')

Module to plot a zoomed in portion of the trace in time domain. This plot zooms in on the overshoot of the didv.

Parameters:
poles : int, string, array_like, optional

The pole fits that we want to plot. If set to “all”, then plots all of the fits. Can also be set to just one of the fits. Can be set as an array of different fits, e.g. [1, 2]

zoomfactor : float, optional, optional

Number between zero and 1 to show different amounts of the zoomed in trace.

plotpriors : boolean, optional

Boolean value on whether or not the priors fit should be plotted.

lgcsave : boolean, optional

Boolean value on whether or not the figure should be saved

savepath : string, optional

Where the figure should be saved. Saved in the current directory by default.

savename : string, optional

A string to append to the end of the file name if saving. Empty string by default.

processtraces()

This method processes the traces loaded to the DIDV class object. This sets up the object for fitting.

qetpy.detcal.iv module

qetpy.detcal.iv.fitfunc(x, b, m)

Function to use for fitting to a straight line

Parameters:
x : array_like

x-values of the data

b : float

y-intercept of line

m : float

slope of line

Returns:
linfunc : array_like

Outputted line for x with slope m and intercept b

qetpy.detcal.iv.findnormalinds(vb, dites, dites_err, tol=10)

Function to determine the indices of the normal data in an IV curve

Parameters:
vb : array_like

Bias voltage, should be a 1d array or list

dites : array_like

The current read out by the electronics with some offset from the true current

dites_err : array_like

The error in the current

tol : float, optional

The tolerance in the reduced chi-squared for the cutoff on which points to use

Returns:
normalinds : iterable

The iterable which stores the range of the normal indices

class qetpy.detcal.iv.IV(dites, dites_err, vb, vb_err, rload, rload_err, chan_names, normalinds=None)

Bases: object

Class for creating the IV curve and calculating various values, such as the normal resistance, the resistance of the TES, the power, etc., as well as the corresponding errors. This class supports data for multple bath temperatures, multiple channels, and multiple bias points.

Note: If different bath temperatures have different numbers of bias points (iters), then the user should pad the end of the arrays with NaN so that the data can be put into an ndarray and loaded into this class.

Attributes:
dites : ndarray

The current read out by the electronics

dites_err : ndarray

The error in the current read out by the electronics

vb : ndarray

The bias voltage (vb = qet bias * rshunt)

vb_err : ndarray

The corresponding error in the bias voltage

rload : scalar, ndarray

The load resistance, this can be scalar if using the same rload for all values. If 1-dimensional, then this should be the load resistance for each channel. If 2-dimensional, this should be the load resistance for each bath temperature and each bias point, where the shape is (ntemps, nch). If 3-dimensional, then this should be with shape (ntemps, nch, niters)

rload_err : scalar, ndarray

The corresponding error in the load resistance, should be the same type as rload

chan_names : array_like

Array of strings corresponding to the names of each channel in the data. Should have the same length as the nch axis in dites

ioff : ndarray

The current offset calculated from the fit, shape (ntemps, nch)

ioff_err : ndarray

The corresponding error in the current offset

rfit : ndarray

The total resistance (rnorm + rload) from the fit, shape (ntemps, nch)

rfit_err : ndarray

The corresponding error in the fit resistance

rnorm : ndarray

The normal resistance of the TES, using the fit, shape (ntemps, nch)

rnorm_err : ndarray

The corresponding error in the normal resistance

ites : ndarray

The calculated current through the TES, shape (ntemps, nch, niters)

ites_err : ndarray

The corresponding error in the current through the TES

r0 : ndarray

The calculated resistance of the TES, shape (ntemps, nch, niters)

r0_err : ndarray

The corresponding error in the resistance of the TES

ptes : ndarray

The calculated power of the TES, shape (ntemps, nch, niters)

ptes_err : ndarray

The corresponding error in the power of the TES

Methods

calc_iv() Method to calculate the IV curve for the intialized object.
plot_all_curves([temps, chans, showfit, …]) Function to plot the IV, resistance, and power curves for the data in an IV object.
plot_iv([temps, chans, showfit, lgcsave, …]) Function to plot the IV curves for the data in an IV object.
plot_pv([temps, chans, lgcsave, savepath, …]) Function to plot the power curves for the data in an IV object.
plot_rv([temps, chans, lgcsave, savepath, …]) Function to plot the resistance curves for the data in an IV object.
calc_iv()

Method to calculate the IV curve for the intialized object. Calculates the power and resistance of each bias point, as well as saving the fit parameters from the fit to the normal points and the calculated normal reistance from these points.

plot_all_curves(temps='all', chans='all', showfit=True, lgcsave=False, savepath='', savename='')

Function to plot the IV, resistance, and power curves for the data in an IV object.

Parameters:
temps : string, array_like, int, optional

Which bath temperatures to plot. Setting to “all” plots all of them. Can also set to a subset of bath temperatures, or just one

chans : string, array_like, int, optional

Which bath temperatures to plot. Setting to “all” plots all of them. Can also set to a subset of bath temperatures, or just one

showfit : boolean, optional

Boolean flag to also plot the linear fit to the normal data

lgcsave : boolean, optional

Boolean flag to save the plot

savepath : string, optional

Path to save the plot to, saves it to the current directory by default

savename : string, optional

Name to append to the plot file name, if saving

plot_iv(temps='all', chans='all', showfit=True, lgcsave=False, savepath='', savename='')

Function to plot the IV curves for the data in an IV object.

Parameters:
temps : string, array_like, int, optional

Which bath temperatures to plot. Setting to “all” plots all of them. Can also set to a subset of bath temperatures, or just one

chans : string, array_like, int, optional

Which bath temperatures to plot. Setting to “all” plots all of them. Can also set to a subset of bath temperatures, or just one

showfit : boolean, optional

Boolean flag to also plot the linear fit to the normal data

lgcsave : boolean, optional

Boolean flag to save the plot

savepath : string, optional

Path to save the plot to, saves it to the current directory by default

savename : string, optional

Name to append to the plot file name, if saving

plot_pv(temps='all', chans='all', lgcsave=False, savepath='', savename='')

Function to plot the power curves for the data in an IV object.

Parameters:
temps : string, array_like, int, optional

Which bath temperatures to plot. Setting to “all” plots all of them. Can also set to a subset of bath temperatures, or just one

chans : string, array_like, int, optional

Which bath temperatures to plot. Setting to “all” plots all of them. Can also set to a subset of bath temperatures, or just one

lgcsave : boolean, optional

Boolean flag to save the plot

savepath : string, optional

Path to save the plot to, saves it to the current directory by default

savename : string, optional

Name to append to the plot file name, if saving

plot_rv(temps='all', chans='all', lgcsave=False, savepath='', savename='')

Function to plot the resistance curves for the data in an IV object.

Parameters:
temps : string, array_like, int, optional

Which bath temperatures to plot. Setting to “all” plots all of them. Can also set to a subset of bath temperatures, or just one

chans : string, array_like, int, optional

Which bath temperatures to plot. Setting to “all” plots all of them. Can also set to a subset of bath temperatures, or just one

lgcsave : boolean, optional

Boolean flag to save the plot

savepath : string, optional

Path to save the plot to, saves it to the current directory by default

savename : string, optional

Name to append to the plot file name, if saving

qetpy.detcal.noise module

qetpy.detcal.noise.slope(x, y, removemeans=True)

Computes the maximum likelihood slope of a set of x and y points.

Parameters:
x : array_like

Array of real-valued independent variables.

y : array_like

Array of real-valued dependent variables.

removemeans : boolean

Boolean flag for if the mean of x should be subtracted. This should be set to True if x has not already had its mean subtracted. Set to False if the mean has been subtracted. Default is True.

Returns:
slope : float

Maximum likelihood slope estimate, calculated as sum((x-<x>)(y-<y>))/sum((x-<x>)**2)

qetpy.detcal.noise.fill_negatives(arr)

Simple helper function to remove negative and zero values from PSD’s and replace them with interpolated values.

Parameters:
arr: ndarray

Array of values to replace neagive values on

Returns:
arr: ndarray

Modified input array with the negative and zero values replace by interpelate values

class qetpy.detcal.noise.Noise(traces, fs, channames, tracegain=1.0, fname=None, time=None)

Bases: object

This class allows the user to calculate the power spectral densities of signals from detectors, study correlations, and decouple the intrinsic noise from cross channel correlated noise.

Attributes:
traces : ndarray

Array of the traces to use in the noise analysis. Should be shape (# of traces, # of channels, # of bins)

fs : float

The digitization rate of the data in Hz.

channames : list

A list of strings that name each of the channels.

time : ndarray

The time values for each bin in each trace.

fname : str

The file name of the data, this will be used when saving the file.

tracegain : float

The factor that traces should be divided by to convert the units to Amps. If rawtraces already has units of Amps, then this should be set to 1.0

freqs : ndarray

The frequencies that correspond to each value in the spectral densities

psd : ndarray

The power spectral density of the data in A^2/Hz

realpsd : ndarray

The real power spectral density of the data in A^2/Hz

imagpsd : ndarray

The imaginary power spectral density of the data in A^2/Hz

corrcoeff : ndarray

The array of the correlation coefficients between each of the channels

uncorrnoise : ndarray

The uncorrelated noise psd in A^2/Hz

corrnoise : ndarray

The correlated noise psd in A^2/Hz

realcsd : ndarray

The real part of the cross spectral density in A^2/Hz

imagcsd : ndarray

The imaginary part of the cross spectral density in A^2/Hz

realcsdstd : ndarray

The standard deviation of the real part of the cross spectral density at each frequency

imagcsdstd : ndarray

The standard deviation of the imaginary part of the cross spectral density at each frequency

csd : ndarray

The cross spectral density of the traces

chandict : dict

A dictionary that stores the channel number for each channel name.

Methods

calculate_corrcoeff() Calculates the correlations between channels as a function of frequency.
calculate_csd() Calculates the csd for each channel in traces.
calculate_psd() Calculates the psd for each channel in traces.
calculate_uncorr_noise() Calculates the uncorrelated and correlated total noise.
plot_corrcoeff([lgcsmooth, nwindow, …]) Function to plot the cross channel correlation coefficients.
plot_csd([whichcsd, lgcreal, lgcsave, savepath]) Function to plot the cross channel noise spectrum referenced to the TES line in units of Amperes^2/Hz
plot_decorrelatednoise([lgcoverlay, …]) Function to plot the de-correlated noise spectrum referenced to the TES line in units of Amperes/sqrt(Hz) from fitted parameters calculated calculate_deCorrelated_noise
plot_psd([lgcoverlay, lgcsave, savepath]) Function to plot the noise spectrum referenced to the TES line in units of Amperes/sqrt(Hz).
plot_reim_psd([lgcsave, savepath]) Function to plot the real vs imaginary noise spectrum referenced to the TES line in units of Amperes/sqrt(Hz).
remove_trace_slope() Function to remove the slope from each trace.
save(path) Saves the noise object as a pickle file
calculate_corrcoeff()

Calculates the correlations between channels as a function of frequency. Stores results in self.corrcoeff

calculate_csd()

Calculates the csd for each channel in traces. Stores csd in self.csd

calculate_psd()

Calculates the psd for each channel in traces. Stores psd in self.psd

calculate_uncorr_noise()

Calculates the uncorrelated and correlated total noise.

plot_corrcoeff(lgcsmooth=True, nwindow=7, lgcsave=False, savepath=None)

Function to plot the cross channel correlation coefficients. Since there are typically few traces, the correlations are often noisy. a savgol_filter is used to smooth out some of the noise

Parameters:
lgcsmooth : boolean, optional

If True, a savgol_filter will be used when plotting.

nwindow : int, optional

the number of bins used for the window in the savgol_filter

lgcsave : boolean, optional

If True, the figure is saved in the user provided directory

savepath : str, optional

Absolute path for the figure to be saved

plot_csd(whichcsd=['01'], lgcreal=True, lgcsave=False, savepath=None)

Function to plot the cross channel noise spectrum referenced to the TES line in units of Amperes^2/Hz

Parameters:
whichcsd : list, optional

a list of strings, where each element of the list refers to the pair of indices of the desired csd plot

lgcreal : boolean, optional

If True, the Re(csd) is plotted. If False, the Im(csd) is plotted

lgcsave : boolean, optional

If True, the figure is saved in the user provided directory

savepath : str, optional

Absolute path for the figure to be saved

plot_decorrelatednoise(lgcoverlay=False, lgcdata=True, lgcuncorrnoise=True, lgccorrelated=False, lgcsum=False, lgcsave=False, savepath=None)

Function to plot the de-correlated noise spectrum referenced to the TES line in units of Amperes/sqrt(Hz) from fitted parameters calculated calculate_deCorrelated_noise

Parameters:
lgcoverlay : boolean, optional

If True, de-correlated for all channels are overlayed in a single plot, If False, the noise for each channel is plotted in a seperate subplot

lgcdata : boolean, optional

Only applies when lgcoverlay = False. If True, the csd data is plotted

lgcuncorrnoise : boolean, optional

Only applies when lgcoverlay = False. If True, the de-correlated noise is plotted

lgccorrelated : boolean, optional

Only applies when lgcoverlay = False. If True, the correlated component of the fitted noise is plotted

lgcsum : boolean, optional

Only applies when lgcoverlay = False. If True, the sum of the fitted de-correlated noise and and correlated noise is plotted

lgcsave : boolean, optional

If True, the figure is saved in the user provided directory

savepath : str, optional

Absolute path for the figure to be saved

plot_psd(lgcoverlay=True, lgcsave=False, savepath=None)

Function to plot the noise spectrum referenced to the TES line in units of Amperes/sqrt(Hz).

Parameters:
lgcoverlay : boolean, optional

If True, psd’s for all channels are overlayed in a single plot, If False, each psd for each channel is plotted in a seperate subplot

lgcsave : boolean, optional

If True, the figure is saved in the user provided directory

savepath : str, optional

Absolute path for the figure to be saved

plot_reim_psd(lgcsave=False, savepath=None)

Function to plot the real vs imaginary noise spectrum referenced to the TES line in units of Amperes/sqrt(Hz). This is done to check for thermal muon tails making it passed the quality cuts

Parameters:
lgcsave : boolean, optional

If True, the figure is saved in the user provided directory

savepath : str, optional

Absolute path for the figure to be saved

remove_trace_slope()

Function to remove the slope from each trace. self.traces is changed to be the slope subtracted traces.

save(path)

Saves the noise object as a pickle file

Parameters:
path : str

Path where the noise object should be saved.

qetpy.detcal.sim module

qetpy.detcal.sim.loadfromdidv(DIDVobj, G=5e-10, qetbias=0.00016, tc=0.04, tload=0.9, tbath=0.02, squiddc=2.5e-12, squidpole=0.0, squidn=1.0, noisetype='transition', lgcpriors=False)

Function for loading the parameters from a DIDV class object.

Parameters:
DIDVobj : Object

A DIDV class object after a fit has been run, such that there are Irwin parameters that can be used to model the noise.

G : float, optional

The thermal conductance of the TES in W/K

qetbias : float, optional

The QET bias in Amps

tc : float

The critical temperature of the TES in K

tload : float

The effective temperature of the load resistor in K

tbath : float

The bath temperature in K

squiddc : float, optional

The DC value of the SQUID and downstream electronics noise, in Amps/rtHz. The SQUID/electronics noise should have been fit beforehand, using the following model:

(squiddc*(1.0+(squidpole/f)**squidn))**2.0

squidpole : float, optional

The frequency pole for the SQUID and downstream electronics noise, in Hz. The SQUID/electronics noise should have been fit beforehand, using the following model:

(squiddc*(1.0+(squidpole/f)**squidn))**2.0

squidn : float, optional

The power of the SQUID and downstream electronics noise, in Hz. The SQUID/electronics noise should have been fit beforehand, using the following model:

(squiddc*(1.0+(squidpole/f)**squidn))**2.0

noisetype : str, optional

The type of the noise that is to be loaded. The options are transition : Use the Irwin parameters from the two pole fit as the transition noise model superconducting : Use the Irwin parameters from the one pole fit as the superconducting noise model normal : Use the Irwin parameters from the one pole fit as the normal noise model

lgcpriors : bool, optional

If True, the priors fit values are loaded from the didv object, if False, the regular fit values are loaded

Returns:
TESobj : Object

A TESnoise class object with all of the fit parameters loaded.

class qetpy.detcal.sim.TESnoise(freqs=None, rload=0.012, r0=0.15, rshunt=0.005, beta=1.0, loopgain=10.0, inductance=4e-07, tau0=0.0005, G=5e-10, qetbias=0.00016, tc=0.04, tload=0.9, tbath=0.02, n=5.0, lgcb=True, squiddc=2.5e-12, squidpole=0.0, squidn=1.0)

Bases: object

Class for the simulation of the TES noise using the simple Irwin theory. Supports noise simulation for in transition, superconducting, and normal.

Attributes:
freqs : float, array_like

The frequencies for which we will calculate the noise simulation

rload : float

The load resistance of the TES (sum of shunt and parasitic resistances) in Ohms

r0 : float

The bias resistance of the TES in Ohms

rshunt : float

The shunt resistance of the TES circuit in Ohms

beta : float

The current sensitivity of the TES (dlogR/dlogI), unitless

loopgain : float

The Irwin loop gain of the TES, unitless

inductance : float

The inductance of the TES circuit in Henries

tau0 : float

The thermal time constant (equals C/G) in s

G : float

The thermal conductance of the TES in W/K

qetbias : float

The QET bias in Amps

tc : float

The critical temperature of the TES in K

tload : float

The effective temperature of the load resistor in K

tbath : float

The bath temperature in K

n : float

The power-law dependence of the power flow to the heat bath

lgcb : boolean

Boolean flag that determines whether we use the ballistic (True) or diffusive limit when calculating TFN power noise

squiddc : float

The frequency pole for the SQUID and downstream electronics noise, in Hz. The SQUID/electronics noise should have been fit beforehand, using the following model:

(squiddc*(1.0+(squidpole/f)**squidn))**2.0

squidpole : float

The frequency pole for the SQUID and downstream electronics noise, in Hz. The SQUID/electronics noise should have been fit beforehand, using the following model:

(squiddc*(1.0+(squidpole/f)**squidn))**2.0

squidn : float

The power of the SQUID and downstream electronics noise, in Hz. The SQUID/electronics noise should have been fit beforehand, using the following model:

(squiddc*(1.0+(squidpole/f)**squidn))**2.0

f_tfn : float

Function that estimates the noise suppression of the thermal fluctuation noise due to the difference in temperature between the bath and the TES. Supports the ballistic and diffusive limits, which is chosen via lgcb

Methods

dIdP([freqs]) The two-pole dIdP function determined from the TES parameters.
dIdV([freqs]) The two-pole dIdV function determined from the TES parameters.
dIdVnormal([freqs]) The one-pole dIdV function determined from the TES parameters for when the TES is normal.
dIdVsc([freqs]) The one-pole dIdV function determined from the TES parameters for when the TES is superconducting.
s_iload([freqs]) The Johnson load current noise determined from the TES parameters for in transition.
s_iloadnormal([freqs]) The Johnson load current noise determined from the TES parameters for normal.
s_iloadsc([freqs]) The Johnson load current noise determined from the TES parameters for superconducting.
s_isquid([freqs]) The SQUID and downstream electronics current noise, currently is using a 1/f model that must be specified when initializing the class.
s_ites([freqs]) The Johnson TES current noise determined from the TES parameters for in transition.
s_itesnormal([freqs]) The Johnson TES current noise determined from the TES parameters for normal.
s_itfn([freqs]) The thermal fluctuation noise in current determined from the TES parameters for in transition.
s_itot([freqs]) The total current noise for the TES in transition.
s_itotnormal([freqs]) The total current noise for the TES when normal.
s_itotsc([freqs]) The total current noise for the TES when superconducting.
s_pload([freqs]) The Johnson load power noise determined from the TES parameters for in transition.
s_psquid([freqs]) The SQUID and downstream electronics power noise, currently is using a 1/f model that must be specified when initializing the class.
s_ptes([freqs]) The Johnson TES power noise determined from the TES parameters for in transition.
s_ptfn([freqs]) The thermal fluctuation noise in power determined from the TES parameters for in transition.
s_ptot([freqs]) The total power noise for the TES in transition.
s_vload([freqs]) The Johnson load voltage noise determined from the TES parameters.
s_vtes([freqs]) The Johnson TES voltage noise determined from the TES parameters for in transition.
s_vtesnormal([freqs]) The Johnson TES voltage noise determined from the TES parameters for normal.
dIdP(freqs=None)

The two-pole dIdP function determined from the TES parameters.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
dIdP : float, ndarray

The two-pole dIdP function

dIdV(freqs=None)

The two-pole dIdV function determined from the TES parameters.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
dIdV : float, ndarray

The two-pole dIdV function

dIdVnormal(freqs=None)

The one-pole dIdV function determined from the TES parameters for when the TES is normal.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
dIdVnormal : float, ndarray

The one-pole dIdV function for when the TES is normal.

dIdVsc(freqs=None)

The one-pole dIdV function determined from the TES parameters for when the TES is superconducting.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
dIdVsc : float, ndarray

The one-pole dIdV function for when the TES is superconducting.

s_iload(freqs=None)

The Johnson load current noise determined from the TES parameters for in transition.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_iload : float, ndarray

The Johnson load current noise at the specified frequencies

s_iloadnormal(freqs=None)

The Johnson load current noise determined from the TES parameters for normal.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_iloadnormal : float, ndarray

The Johnson load current noise at the specified frequencies

s_iloadsc(freqs=None)

The Johnson load current noise determined from the TES parameters for superconducting.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_iloadsc : float, ndarray

The Johnson load current noise at the specified frequencies

s_isquid(freqs=None)

The SQUID and downstream electronics current noise, currently is using a 1/f model that must be specified when initializing the class.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_isquid : float, ndarray

The SQUID and downstream electronics current noise at the specified frequencies

s_ites(freqs=None)

The Johnson TES current noise determined from the TES parameters for in transition. This noise has both an electronic and thermal component.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_ites : float, ndarray

The Johnson TES current noise at the specified frequencies

s_itesnormal(freqs=None)

The Johnson TES current noise determined from the TES parameters for normal.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_itesnormal : float, ndarray

The Johnson TES current noise at the specified frequencies

s_itfn(freqs=None)

The thermal fluctuation noise in current determined from the TES parameters for in transition.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_itfn : float, ndarray

The thermal fluctuation noise in current at the specified frequencies

s_itot(freqs=None)

The total current noise for the TES in transition. This is calculated by summing each of the current noise sources together. Units are [A^2/Hz].

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_itot : float, ndarray

The total current noise at the specified frequencies

s_itotnormal(freqs=None)

The total current noise for the TES when normal. This is calculated by summing each of the current noise sources together. Units are [A^2/Hz].

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_itotnormal : float, ndarray

The total current noise at the specified frequencies

s_itotsc(freqs=None)

The total current noise for the TES when superconducting. This is calculated by summing each of the current noise sources together. Units are [A^2/Hz].

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_itotsc : float, ndarray

The total current noise at the specified frequencies

s_pload(freqs=None)

The Johnson load power noise determined from the TES parameters for in transition.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_pload : float, ndarray

The Johnson load power noise at the specified frequencies

s_psquid(freqs=None)

The SQUID and downstream electronics power noise, currently is using a 1/f model that must be specified when initializing the class. This is only used for when the TES is in transition.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_psquid : float, ndarray

The SQUID and downstream electronics power noise at the specified frequencies

s_ptes(freqs=None)

The Johnson TES power noise determined from the TES parameters for in transition.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_ptes : float, ndarray

The Johnson TES power noise at the specified frequencies

s_ptfn(freqs=None)

The thermal fluctuation noise in power determined from the TES parameters for in transition.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_ptfn : float, ndarray

The thermal fluctuation noise in power at the specified frequencies

s_ptot(freqs=None)

The total power noise for the TES in transition. This is calculated by summing each of the current noise sources together and using dIdP to convert to power noise. Units are [W^2/Hz].

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_ptot : float, ndarray

The total power noise at the specified frequencies

s_vload(freqs=None)

The Johnson load voltage noise determined from the TES parameters. This formula holds no matter where we are in transition.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_vload : float, ndarray

The Johnson load voltage noise at the specified frequencies

s_vtes(freqs=None)

The Johnson TES voltage noise determined from the TES parameters for in transition.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_vtes : float, ndarray

The Johnson TES voltage noise at the specified frequencies

s_vtesnormal(freqs=None)

The Johnson TES voltage noise determined from the TES parameters for normal.

Parameters:
freqs : float, ndarray, optional

The frequencies for which we will calculate the noise simulation. If left as None, the function will use the values from the initialization.

Returns:
s_vtesnormal : float, ndarray

The Johnson TES voltage noise at the specified frequencies

Module contents