DDASToys for NSCLDAQ  6.2-000
Classes
libFitEditorAnalytic.so

Plugin library for analytic fitting. More...

Collaboration diagram for libFitEditorAnalytic.so:

Classes

class  ddastoys::Configuration
 Manage fit configuration information. More...
 
struct  ddastoys::analyticfit::GslFitParameters
 Data passed around the fitting subsystem to Jacobian and function evaluators. More...
 
double ddastoys::analyticfit::logistic (double A, double k, double x1, double x)
 Evaluate a logistic function for the specified parameters and point. More...
 
double ddastoys::analyticfit::decay (double A, double k, double x1, double x)
 Evaluate an exponential decay for the specific parameters and point. More...
 
double ddastoys::analyticfit::switchOn (double x1, double x)
 Logistic switch to turn on a function evaluation at a given point. More...
 
double ddastoys::analyticfit::singlePulse (double A1, double k1, double k2, double x1, double C, double x)
 Evaluate the single pulse fit function at a given point. More...
 
double ddastoys::analyticfit::doublePulse (double A1, double k1, double k2, double x1, double A2, double k3, double k4, double x2, double C, double x)
 Evaluate the double pulse fit function at a given point. More...
 
double ddastoys::analyticfit::pulseAmplitude (double A, double k1, double k2, double x0)
 Calculate the pulse amplitude corrected for the ballistic deficit imposed by the exponential decay. More...
 
double ddastoys::analyticfit::chiSquare1 (double A1, double k1, double k2, double x1, double C, const std::vector< uint16_t > &trace, int low=0, int high=-1)
 Computes the chi-square goodness-of-fit for a specific parameterization of a single pulse canonical form with respect to a trace. More...
 
double ddastoys::analyticfit::chiSquare1 (double A1, double k1, double k2, double x1, double C, const std::vector< std::pair< uint16_t, uint16_t > > &points)
 Computes the chi-square goodness-of-fit for a specific parameterization of a single pulse canonical form with respect to a trace. More...
 
double ddastoys::analyticfit::chiSquare2 (double A1, double k1, double k2, double x1, double A2, double k3, double k4, double x2, double C, const std::vector< uint16_t > &trace, int low=0, int high=-1)
 Computes the chi-square goodness of a specific parameterization of a double pulse canonical form with respect to a trace. More...
 
double ddastoys::analyticfit::chiSquare2 (double A1, double k1, double k2, double x1, double A2, double k3, double k4, double x2, double C, const std::vector< std::pair< uint16_t, uint16_t > > &points)
 Computes the chi-square goodness of a specific parameterization of a double pulse canonical form with respect to a trace. More...
 
void ddastoys::analyticfit::writeTrace (const char *filename, const char *title, const std::vector< uint16_t > &trace)
 Write a single trace to a file. More...
 
void ddastoys::analyticfit::writeTrace2 (const char *filename, const char *title, const std::vector< uint16_t > &trace1, const std::vector< uint16_t > &trace2)
 Write two traces to a file. More...
 
void ddastoys::analyticfit::lmfit1 (fit1Info *pResult, std::vector< uint16_t > &trace, const std::pair< unsigned, unsigned > &limits, uint16_t saturation=0xffff)
 Driver for the GSL LM fitter for single pulses. More...
 
void ddastoys::analyticfit::lmfit2 (fit2Info *pResult, std::vector< uint16_t > &trace, const std::pair< unsigned, unsigned > &limits, fit1Info *pSinglePulseFit=nullptr, uint16_t saturation=0xffff)
 Driver for the GSL LM fitter for double pulses. More...
 
void ddastoys::analyticfit::lmfit2fixedT (fit2Info *pResult, std::vector< uint16_t > &trace, const std::pair< unsigned, unsigned > &limits, fit1Info *pSinglePulseFit=nullptr, uint16_t saturation=0xffff)
 Driver for the GSL LM fitter for double pulses, constraining the two timing parameters (rise time and fall time) to be the same for both pulses. More...
 

Detailed Description

Plugin library for analytic fitting.

Function Documentation

◆ chiSquare1() [1/2]

double ddastoys::analyticfit::chiSquare1 ( double  A1,
double  k1,
double  k2,
double  x1,
double  C,
const std::vector< std::pair< uint16_t, uint16_t > > &  points 
)

Computes the chi-square goodness-of-fit for a specific parameterization of a single pulse canonical form with respect to a trace.

Parameters
A1Amplitude of pulse
k1Steepness of pulse rise.
k2Decay time of pulse fall.
x1Position of the pulse.
CConstant offset of the trace.
pointsSet of (x, y) data points which contribute to the total chi-square value
Returns
The chi-square goodness-of-fit statistic.

Neyman's chi-square. The chi-square value is calculated from a passed set of (x, y) points.

◆ chiSquare1() [2/2]

double ddastoys::analyticfit::chiSquare1 ( double  A1,
double  k1,
double  k2,
double  x1,
double  C,
const std::vector< uint16_t > &  trace,
int  low = 0,
int  high = -1 
)

Computes the chi-square goodness-of-fit for a specific parameterization of a single pulse canonical form with respect to a trace.

Parameters
A1Amplitude of pulse
k1Steepness of pulse rise.
k2Decay time of pulse fall.
x1Position of the pulse.
CConstant offset of the trace.
traceTrace to compute the chi-square value with respect to.
low,highRegion of interest over which to compute the chi-square value.
Note
high = -1 will set the high limit to the last sample in the trace.
Returns
The chi-square goodness-of-fit statistic.

Neyman's chi-square. The chi-square value is calculated based on the passed limits low, high.

◆ chiSquare2() [1/2]

double ddastoys::analyticfit::chiSquare2 ( double  A1,
double  k1,
double  k2,
double  x1,
double  A2,
double  k3,
double  k4,
double  x2,
double  C,
const std::vector< std::pair< uint16_t, uint16_t > > &  points 
)

Computes the chi-square goodness of a specific parameterization of a double pulse canonical form with respect to a trace.

Parameters
A1Amplitude of the first pulse.
k1Steepness of first pulse rise.
k2Decay time of the first pulse.
x1Position of the first pulse.
A2Amplitude of the second pulse.
k3Steepness of second pulse rise.
k4Decay time of second pulse.
x2Position of second pulse.
CConstant offset the pulses sit on.
pointsSet of x, y data points which contribute to the total chi-square value
Returns
The chi-square goodness-of-fit statistic.

Neyman's chi-square. The chi-square value is calculated from a passed set of (x, y) points.

◆ chiSquare2() [2/2]

double ddastoys::analyticfit::chiSquare2 ( double  A1,
double  k1,
double  k2,
double  x1,
double  A2,
double  k3,
double  k4,
double  x2,
double  C,
const std::vector< uint16_t > &  trace,
int  low = 0,
int  high = -1 
)

Computes the chi-square goodness of a specific parameterization of a double pulse canonical form with respect to a trace.

Parameters
A1Amplitude of the first pulse.
k1Steepness of first pulse rise.
k2Decay time of the first pulse.
x1Position of the first pulse.
A2Amplitude of the second pulse.
k3Steepness of second pulse rise.
k4Decay time of second pulse.
x2Position of second pulse.
CConstant offset the pulses sit on.
traceTrace to compute the chisquare with respect to.
low,highRegion of interest over which to compute the chi-square value.
Note
high = -1 will set the high limit to the last sample in the trace.
Returns
The chi-square goodness-of-fit statistic.

Neyman's chi-square. The chi-square values is calculated based on the passed limits low, high.

◆ decay()

double ddastoys::analyticfit::decay ( double  A,
double  k,
double  x1,
double  x 
)

Evaluate an exponential decay for the specific parameters and point.

Parameters
AAmplitude of the signal
kDecay constant of the signal.
x1Position of the pulse.
xWhere to evaluate the signal.
Returns
Function value at x.

Signals from detectors usually have a falling shape that approximates an exponential. This function evaluates this decay at some point.

◆ doublePulse()

double ddastoys::analyticfit::doublePulse ( double  A1,
double  k1,
double  k2,
double  x1,
double  A2,
double  k3,
double  k4,
double  x2,
double  C,
double  x 
)

Evaluate the double pulse fit function at a given point.

Parameters
A1Amplitude of the first pulse.
k1Steepness of first pulse rise.
k2Decay time of the first pulse.
x1Position of the first pulse.
A2Amplitude of the second pulse.
k3Steepness of second pulse rise.
k4Decay time of second pulse.
x2Position of second pulse.
CConstant offset the pulses sit on.
xPosition at which to evaluate the pulse.
Returns
Double pulse function evaluated at x.

Evaluate the canonical form of a double pulse. This is done by summing two single pulses. The constant term is thrown into the first pulse. The second pulse gets a constant term of 0.

◆ lmfit1()

void ddastoys::analyticfit::lmfit1 ( fit1Info pResult,
std::vector< uint16_t > &  trace,
const std::pair< unsigned, unsigned > &  limits,
uint16_t  saturation = 0xffff 
)

Driver for the GSL LM fitter for single pulses.

Parameters
pResultStruct that will get the results of the fit.
traceVector of trace points.
limitsLimits of the trace over which to conduct the fit.
saturationValue at which the ADC is saturated (points at or above this value are removed from the fit.)

Implements GSL's Levenburg-Marquardt fitter to determine best-fit parameters by chi-square minimization.

◆ lmfit2()

void ddastoys::analyticfit::lmfit2 ( fit2Info pResult,
std::vector< uint16_t > &  trace,
const std::pair< unsigned, unsigned > &  limits,
fit1Info pSinglePulseFit = nullptr,
uint16_t  saturation = 0xffff 
)

Driver for the GSL LM fitter for double pulses.

Parameters
pResultResults will be stored here.
traceReferences the trace to fit.
limitsThe limits of the trace that can be fit.
pSinglePulseFitThe fit for a single pulse, used to seed initial guesses if present. Otherwise a single pulse fit is done for that.
saturationADC saturation value. Points with values at or above this value are removed from the fit.

Implements GSL's Levenburg-Marquardt fitter to determine best-fit parameters by chi-square minimization.

◆ lmfit2fixedT()

void ddastoys::analyticfit::lmfit2fixedT ( fit2Info pResult,
std::vector< uint16_t > &  trace,
const std::pair< unsigned, unsigned > &  limits,
fit1Info pSinglePulseFit = nullptr,
uint16_t  saturation = 0xffff 
)

Driver for the GSL LM fitter for double pulses, constraining the two timing parameters (rise time and fall time) to be the same for both pulses.

Parameters
pResultResults will be stored here.
traceReferences the trace to fit.
limitsThe limits of the trace that can be fit.
pSinglePulseFitThe fit for a single pulse, used to seed initial uesses if present. Otherwise a single pulse fit is done for that.
saturationADC saturation value. Points with values at or above this value are removed from the fit.

Implements GSL's Levenburg-Marquardt fitter to determine best-fit parameters by chi-square minimization.

◆ logistic()

double ddastoys::analyticfit::logistic ( double  A,
double  k,
double  x1,
double  x 
)

Evaluate a logistic function for the specified parameters and point.

Parameters
AAmplitude of the signal.
kSteepness of the signal (related to the rise time).
x1Mid-point of the rise of the logistic.
xLocation at which to evaluate the function.
Returns
Function value at x.

A logistic function is a function with a sigmoidal shape. We use it to fit the rising edge of signals from DDAS digitizes. See e.g. https://en.wikipedia.org/wiki/Logistic_function for a discussion of this function.

◆ pulseAmplitude()

double ddastoys::analyticfit::pulseAmplitude ( double  A,
double  k1,
double  k2,
double  x0 
)

Calculate the pulse amplitude corrected for the ballistic deficit imposed by the exponential decay.

Parameters
AThe scaling term of the pulse.
k1The steepness term of the logistic.
k2The fall time term of the decay.
x0The position of the pulse.
Returns
The corrected amplitude of the fitted pulse.
Return values
-1If $k_1/k_2 \leq 1$.

The A1 term in a pulse fit is not actually the "true" pulse amplitude. The effect of the exponential decay is already important causing A1 to over-estimate the amplitude.

This function computes the "true" amplitude of a pulse given its parameters. This is done by noting that the deriviative of the pulse has a zero at:

$x = x_0 + \log(k_1/k_2 - 1)/k_1$

We plug that position back into the pulse to get the amplitude.

◆ singlePulse()

double ddastoys::analyticfit::singlePulse ( double  A1,
double  k1,
double  k2,
double  x1,
double  C,
double  x 
)

Evaluate the single pulse fit function at a given point.

Parameters
A1Pulse amplitiude.
k1Logistic rise steepness.
k2Exponential decay time constant.
x1Logistic position.
CConstant offset.
xPosition at which to evaluate the function.
Returns
double Single pulse function evaluated at x.

Evaluate the value of a single pulse in accordance with our canonical functional form. The form is a logistic rise with an exponential decay that sits on top of a constant offset.

◆ switchOn()

double ddastoys::analyticfit::switchOn ( double  x1,
double  x 
)

Logistic switch to turn on a function evaluation at a given point.

Parameters
x1Switch on point. The switch is "on" for values greater than x1.
xPoint in space at which to evaluate the switch.
Returns
The switch value: very nearly 0.0 left of x1 and very nearly 1.0 to the right of x1.

Provides a multipier that can be used to turn on a function at a point in space. This is provided for GPUs that may want to define functions with conditional chunks without the associated flow control divergence conditionals would require. This is implemented as a very steep logistic with rise centered at the point in question.

◆ writeTrace()

void ddastoys::analyticfit::writeTrace ( const char *  filename,
const char *  title,
const std::vector< uint16_t > &  trace 
)

Write a single trace to a file.

Parameters
filenameWhere to write.
titleTitle string.
traceThe trace data.

◆ writeTrace2()

void ddastoys::analyticfit::writeTrace2 ( const char *  filename,
const char *  title,
const std::vector< uint16_t > &  t1,
const std::vector< uint16_t > &  t2 
)

Write two traces to a file.

Parameters
filenameWhere to write.
titleTitle string.
t1The first trace.
t2The second trace.
Note
The traces must be the same length.