|
DDASToys for NSCLDAQ
6.2-000
|
Plugin library for analytic fitting. More...

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... | |
Plugin library for analytic fitting.
| 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.
| A1 | Amplitude of pulse |
| k1 | Steepness of pulse rise. |
| k2 | Decay time of pulse fall. |
| x1 | Position of the pulse. |
| C | Constant offset of the trace. |
| points | Set of (x, y) data points which contribute to the total chi-square value |
Neyman's chi-square. The chi-square value is calculated from a passed set of (x, y) points.
| 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.
| A1 | Amplitude of pulse |
| k1 | Steepness of pulse rise. |
| k2 | Decay time of pulse fall. |
| x1 | Position of the pulse. |
| C | Constant offset of the trace. |
| trace | Trace to compute the chi-square value with respect to. |
| low,high | Region of interest over which to compute the chi-square value. |
Neyman's chi-square. The chi-square value is calculated based on the passed limits low, high.
| 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.
| A1 | Amplitude of the first pulse. |
| k1 | Steepness of first pulse rise. |
| k2 | Decay time of the first pulse. |
| x1 | Position of the first pulse. |
| A2 | Amplitude of the second pulse. |
| k3 | Steepness of second pulse rise. |
| k4 | Decay time of second pulse. |
| x2 | Position of second pulse. |
| C | Constant offset the pulses sit on. |
| points | Set of x, y data points which contribute to the total chi-square value |
Neyman's chi-square. The chi-square value is calculated from a passed set of (x, y) points.
| 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.
| A1 | Amplitude of the first pulse. |
| k1 | Steepness of first pulse rise. |
| k2 | Decay time of the first pulse. |
| x1 | Position of the first pulse. |
| A2 | Amplitude of the second pulse. |
| k3 | Steepness of second pulse rise. |
| k4 | Decay time of second pulse. |
| x2 | Position of second pulse. |
| C | Constant offset the pulses sit on. |
| trace | Trace to compute the chisquare with respect to. |
| low,high | Region of interest over which to compute the chi-square value. |
Neyman's chi-square. The chi-square values is calculated based on the passed limits low, high.
| double ddastoys::analyticfit::decay | ( | double | A, |
| double | k, | ||
| double | x1, | ||
| double | x | ||
| ) |
Evaluate an exponential decay for the specific parameters and point.
| A | Amplitude of the signal |
| k | Decay constant of the signal. |
| x1 | Position of the pulse. |
| x | Where to evaluate the signal. |
Signals from detectors usually have a falling shape that approximates an exponential. This function evaluates this decay at some point.
| 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.
| A1 | Amplitude of the first pulse. |
| k1 | Steepness of first pulse rise. |
| k2 | Decay time of the first pulse. |
| x1 | Position of the first pulse. |
| A2 | Amplitude of the second pulse. |
| k3 | Steepness of second pulse rise. |
| k4 | Decay time of second pulse. |
| x2 | Position of second pulse. |
| C | Constant offset the pulses sit on. |
| x | Position at which to evaluate the pulse. |
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.
| 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.
| pResult | Struct that will get the results of the fit. |
| trace | Vector of trace points. |
| limits | Limits of the trace over which to conduct the fit. |
| saturation | Value 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.
| 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.
| pResult | Results will be stored here. |
| trace | References the trace to fit. |
| limits | The limits of the trace that can be fit. |
| pSinglePulseFit | The fit for a single pulse, used to seed initial guesses if present. Otherwise a single pulse fit is done for that. |
| saturation | ADC 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.
| 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.
| pResult | Results will be stored here. |
| trace | References the trace to fit. |
| limits | The limits of the trace that can be fit. |
| pSinglePulseFit | The fit for a single pulse, used to seed initial uesses if present. Otherwise a single pulse fit is done for that. |
| saturation | ADC 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.
| double ddastoys::analyticfit::logistic | ( | double | A, |
| double | k, | ||
| double | x1, | ||
| double | x | ||
| ) |
Evaluate a logistic function for the specified parameters and point.
| A | Amplitude of the signal. |
| k | Steepness of the signal (related to the rise time). |
| x1 | Mid-point of the rise of the logistic. |
| x | Location at which to evaluate the function. |
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.
| 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.
| A | The scaling term of the pulse. |
| k1 | The steepness term of the logistic. |
| k2 | The fall time term of the decay. |
| x0 | The position of the pulse. |
| -1 | If |
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:
![]()
We plug that position back into the pulse to get the amplitude.
| 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.
| A1 | Pulse amplitiude. |
| k1 | Logistic rise steepness. |
| k2 | Exponential decay time constant. |
| x1 | Logistic position. |
| C | Constant offset. |
| x | Position at which to evaluate the function. |
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.
| double ddastoys::analyticfit::switchOn | ( | double | x1, |
| double | x | ||
| ) |
Logistic switch to turn on a function evaluation at a given point.
| x1 | Switch on point. The switch is "on" for values greater than x1. |
| x | Point in space at which to evaluate the switch. |
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.
| void ddastoys::analyticfit::writeTrace | ( | const char * | filename, |
| const char * | title, | ||
| const std::vector< uint16_t > & | trace | ||
| ) |
Write a single trace to a file.
| filename | Where to write. |
| title | Title string. |
| trace | The trace data. |
| 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.
| filename | Where to write. |
| title | Title string. |
| t1 | The first trace. |
| t2 | The second trace. |