DDASToys for NSCLDAQ  6.2-000
Classes
libFitEditorTemplate.so

Plugin library for template fitting. More...

Collaboration diagram for libFitEditorTemplate.so:

Classes

class  ddastoys::Configuration
 Manage fit configuration information. More...
 
class  ddastoys::FitEditorTemplate
 Fit trace data with the template fitting functions. More...
 
struct  ddastoys::templatefit::GslFitParameters
 Data passed around the fitting subsystem to Jacobian and function evaluators. More...
 
double ddastoys::templatefit::singlePulse (double S1, double x1, double C, double x, const std::vector< double > &trace_template)
 Evaluate the single pulse fit function at a given point. More...
 
double ddastoys::templatefit::doublePulse (double S1, double x1, double S2, double x2, double C, double x, const std::vector< double > &trace_template)
 Evaluate the double pulse fit function at a given point. More...
 
double ddastoys::templatefit::chiSquare1 (double S1, double x1, double C, const std::vector< std::pair< uint16_t, uint16_t > > &points, const std::vector< double > &trace_template)
 Computes the chi-square goodness of a specific parameterization of a single pulse canonical form with respect to a trace. More...
 
double ddastoys::templatefit::chiSquare2 (double S1, double x1, double S2, double x2, double C, const std::vector< std::pair< uint16_t, uint16_t > > &points, const std::vector< double > &trace_template)
 Computes the chi-square goodness of a specific parameterization of a double pulse canonical form with respect to a trace. More...
 
void ddastoys::templatefit::lmfit1 (fit1Info *pResult, std::vector< uint16_t > &trace, std::vector< double > &traceTemplate, unsigned alignPoint, const std::pair< unsigned, unsigned > &limits, uint16_t saturation=0xffff)
 Driver for the GSL LM fitter for single pulses. More...
 
void ddastoys::templatefit::lmfit2 (fit2Info *pResult, std::vector< uint16_t > &trace, std::vector< double > &traceTemplate, unsigned alignPoint, const std::pair< unsigned, unsigned > &limits, fit1Info *pSinglePulseFit=nullptr, uint16_t saturation=0xffff)
 Driver for the GSL LM fitter for double pulses. More...
 

Detailed Description

Plugin library for template fitting.

Function Documentation

◆ chiSquare1()

double ddastoys::templatefit::chiSquare1 ( double  A1,
double  x1,
double  C,
const std::vector< std::pair< uint16_t, uint16_t > > &  points,
const std::vector< double > &  trace_template 
)

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

Parameters
A1Scaling factor for the template fit.
x1Offset with respect to the template.
CConstant baseline.
pointsSet of (x, y) data points which contribute to the total chi-square value.
trace_templateTemplate trace used to fit the data.
Returns
The chi-square goodness-of-fit statistic.

Neyman's chi-square value is computed from a passed set of (x, y) data.

◆ chiSquare2()

double ddastoys::templatefit::chiSquare2 ( double  A1,
double  x1,
double  A2,
double  x2,
double  C,
const std::vector< std::pair< uint16_t, uint16_t > > &  points,
const std::vector< double > &  trace_template 
)

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

Parameters
A1Scaling factor for the template fit
x1Offset with respect to the template
A2Scaling factor for the template fit
x2Offset with respect to the template
CConstant baseline
pointsSet of x, y data points which contribute to the total chi-square value
trace_templateTemplate trace used to fit the data.
Returns
The chi-square goodness-of-fit statistic.

Neyman's chi-square value is computed from a passed set of (x, y) data.

◆ doublePulse()

double ddastoys::templatefit::doublePulse ( double  A1,
double  x1,
double  A2,
double  x2,
double  C,
double  x,
const std::vector< double > &  trace_template 
)

Evaluate the double pulse fit function at a given point.

Parameters
A1Scaling factor for the template fit.
x1Offset with respect to the template.
A2Scaling factor for the template fit.
x2Offset with respect to the template.
CConstant baseline.
xPosition at which to evaluate this function.
trace_templateTemplate trace used to fit the data.
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::templatefit::lmfit1 ( fit1Info pResult,
std::vector< uint16_t > &  trace,
std::vector< double > &  traceTemplate,
unsigned  alignPoint,
const std::pair< unsigned, unsigned > &  limits,
uint16_t  saturation = 0xffff 
)

Driver for the GSL LM fitter for single pulses.

Parameters
[in,out]pResultStruct that will get the results of the fit.
[in]traceReferences the trace to fit.
[in]traceTemplateReferences the template trace for the fit. The template is assumed to have a baseline of 0 and a maximum value of 1.
[in]alignPointThe internal alignment point of the template.
[in]limitsLimits of the trace over which to conduct the fit.
[in]saturationValue at which the ADC is saturated (points at or above this value are removed from the fit.)

GSL non-linear least-squares fitter with Levenburg-Marquardt method used to solve the "trust-region problem" of approximating the objective function around its minimum. See https://www.gnu.org/software/gsl/doc/html/index.html.

◆ lmfit2()

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

Driver for the GSL LM fitter for double pulses.

Parameters
[in,out]pResultStruct that will get the results of the fit.
[in]traceReferences the trace to fit.
[in]traceTemplateReferences the template trace for the fit. The template is assumed to have a baseline of 0 and a maximum value of 1.
[in]alignPointThe internal alignment point of the template.
[in]limitsThe limits of the trace that can be fit.
[in]pSinglePulseFitPointer to the fit for a single pulse, used to seed initial guesses if present. Otherwise a single pulse fit is done for that.
[in]saturationADC saturation value. Points with values at or above this value are removed from the fit.

GSL non-linear least-squares fitter with Levenburg-Marquardt method used to solve the "trust-region problem" of approximating the objective function around its minimum. See https://www.gnu.org/software/gsl/doc/html/index.html.

◆ singlePulse()

double ddastoys::templatefit::singlePulse ( double  A1,
double  x1,
double  C,
double  x,
const std::vector< double > &  trace_template 
)

Evaluate the single pulse fit function at a given point.

Parameters
A1Scaling factor for the template fit.
x1Offset with respect to the template.
CConstant baseline.
xPosition at which to evaluate this function.
trace_templateTemplate trace used to fit the data.
Returns
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.