Belle II Software development
WaveFitter Class Reference

Waveform fitter class. More...

#include <WaveFitter.h>

Public Types

typedef std::function< double(double)> wave_function_type
 Function type for waveform.
 
typedef std::vector< std::array< double, 6 > > strip_data_type
 Function type for cluster data.
 

Public Member Functions

 WaveFitter (wave_function_type wave, std::array< double, 6 > times)
 Constructor creates an empty fitter with default settings.
 
void addData (const std::array< double, 6 > &inpData, double rmsNoise)
 Add strip data to the fitter.
 
void reset ()
 Reset fitter data.
 
double getFittedTime ()
 Retrieve fitted time shift value after a fit.
 
double getFittedTimeError ()
 Retrieve error of time shift value after a fit.
 
const std::vector< double > & getFittedAmplitudes ()
 Retrieve strip amplitudes after a fit.
 
const std::vector< double > & getFittedAmplitudeErrors ()
 Retrieve errors of strip amplitudes after a fit.
 
double getFitLikelihood ()
 Return negative log-likelihood for the fit.
 
const strip_data_typegetData () const
 Get original fitted data.
 
const strip_data_typegetFitData ()
 Get fitted values for all samples.
 
double getDt () const
 Get sampling interval.
 
const std::array< double, 6 > & getTimes () const
 Get sample times.
 
double Chi (double t)
 Get standardized chi2 for a given time shift.
 
double negLogLikelihood (double t)
 Get negative profile log-likelihood for a given time shift.
 
double pSignal ()
 Calculate probability of the hit being a signal, with given priors for signal and background.
 
double lrSignal (double a, double b)
 Calculate likelihood ratio Lsignal/Lbackground for acceptance (signal) window given as arguments.
 
double wave (double t) const
 Get unit waveform value at a given time.
 

Static Public Member Functions

static void setPrecision (int ndigits)
 Set precision of the minimzer.
 

Private Member Functions

void doFit ()
 Perform fit on the data.
 
void calculateAmplitudes ()
 Calculate amplitudes.
 
void calculateFitErrors ()
 Calculate fit errors.
 
void calculateFittedData ()
 Calculate fitted data.
 
double integral12 (double lower, double upper, std::function< double(double)> f)
 High-order Gauss-Legendre quadrature for likelihood integrals.
 
double integral20 (double lower, double upper, std::function< double(double)> f)
 High-order Gauss-Legendre quadrature for likelihood integrals.
 

Private Attributes

bool m_hasFit
 Are fit results available?
 
double m_dt
 Time interval between samples.
 
strip_data_type m_data
 Vector of sextets of APV samples.
 
std::vector< double > m_noises
 RMS noise for strips.
 
std::array< double, 6 > m_times
 Vector of sample times.
 
wave_function_type m_wave
 Waveform function.
 
double m_fittedTime
 Time from interval (-dt,+dt).
 
double m_fittedTimeError
 Error estimate of the fitted time.
 
std::vector< double > m_fittedAmplitudes
 Fitted amplitudes of the sample.
 
std::vector< double > m_fittedAmplitudeErrors
 Fitted amplitude errors.
 
double m_fittedLik
 Chi2 from the current fit.
 
int m_ndf
 Degrees of freedom of the fit.
 
strip_data_type m_fittedData
 Fitted values from current fit.
 

Static Private Attributes

static int s_minimizer_precision = 6
 minimizer precision
 

Detailed Description

Waveform fitter class.

Definition at line 59 of file WaveFitter.h.

Member Typedef Documentation

◆ strip_data_type

typedef std::vector< std::array<double, 6> > strip_data_type

Function type for cluster data.

Definition at line 65 of file WaveFitter.h.

◆ wave_function_type

typedef std::function<double(double)> wave_function_type

Function type for waveform.

Definition at line 63 of file WaveFitter.h.

Constructor & Destructor Documentation

◆ WaveFitter()

WaveFitter ( wave_function_type  wave,
std::array< double, 6 >  times 
)
inline

Constructor creates an empty fitter with default settings.

Parameters
waveWave function to use for fitting data
timesarray of points in time at which data are taken

Definition at line 71 of file WaveFitter.h.

71 :
72 m_hasFit(false), m_times(times), m_wave(wave), m_fittedTime(0), m_fittedTimeError(), m_fittedLik(1.0e10), m_ndf(0)
73 {
74 m_dt = m_times[1] - m_times[0];
75 }
double m_fittedTimeError
Error estimate of the fitted time.
Definition: WaveFitter.h:219
double m_dt
Time interval between samples.
Definition: WaveFitter.h:213
double wave(double t) const
Get unit waveform value at a given time.
Definition: WaveFitter.h:188
int m_ndf
Degrees of freedom of the fit.
Definition: WaveFitter.h:223
double m_fittedTime
Time from interval (-dt,+dt).
Definition: WaveFitter.h:218
wave_function_type m_wave
Waveform function.
Definition: WaveFitter.h:217
std::array< double, 6 > m_times
Vector of sample times.
Definition: WaveFitter.h:216
double m_fittedLik
Chi2 from the current fit.
Definition: WaveFitter.h:222
bool m_hasFit
Are fit results available?
Definition: WaveFitter.h:212

Member Function Documentation

◆ addData()

void addData ( const std::array< double, 6 > &  inpData,
double  rmsNoise 
)
inline

Add strip data to the fitter.

Parameters
inpDataarray of 6 APV signals.
rmsNoiseRMS noise for the strip.

Definition at line 81 of file WaveFitter.h.

82 {
83 std::array<double, 6> v;
84 std::copy(inpData.begin(), inpData.end(), v.begin());
85 m_data.push_back(v);
86 m_noises.push_back(rmsNoise);
87 m_hasFit = false;
88 m_ndf = 6 * m_data.size() - m_data.size() - 2;
89 }
std::vector< double > m_noises
RMS noise for strips.
Definition: WaveFitter.h:215
strip_data_type m_data
Vector of sextets of APV samples.
Definition: WaveFitter.h:214

◆ calculateAmplitudes()

void calculateAmplitudes ( )
private

Calculate amplitudes.

Definition at line 115 of file WaveFitter.cc.

116{
117 // This we only need to calculate once
118 array<double, 6> waves;
119 double sumTemplateByTemplate = 0;
120 for (int i = 0; i < 6; ++i) {
121 waves[i] = m_wave(m_times[i] - m_fittedTime);
122 sumTemplateByTemplate += waves[i] * waves[i];
123 }
124 // Calculate amplitudes
125 for (unsigned int row = 0; row < m_data.size(); ++row) {
126 double sumDataByTemplate = 0;
127 for (int i = 0; i < 6; i++) {
128 sumDataByTemplate += m_data[row][i] * waves[i];
129 }
130 m_fittedAmplitudes.push_back(sumDataByTemplate / sumTemplateByTemplate);
131 }
132}
std::vector< double > m_fittedAmplitudes
Fitted amplitudes of the sample.
Definition: WaveFitter.h:220

◆ calculateFitErrors()

void calculateFitErrors ( )
private

Calculate fit errors.

This only calculates local (linearization) errors.

Definition at line 134 of file WaveFitter.cc.

135{
136 // Errors on amplitudes - simple from the linear model
137 double sum_w2 = 0;
138 for (int i = 0; i < 6; ++i) {
139 double wave = m_wave(m_times[i] - m_fittedTime);
140 sum_w2 += wave * wave;
141 }
142 double one_by_sqrtSumw2 = 1.0 / sqrt(sum_w2);
143 for (unsigned int row = 0; row < m_data.size(); ++row)
144 m_fittedAmplitudeErrors.push_back(one_by_sqrtSumw2 * m_noises[row]);
145 // For time: ignore correlation with amplitudes
146 // We will need derivatives, calculate them numerically
147 double timestep = 0.1 * m_dt;
148 double dwave_sq = 0;
149 for (int i = 0; i < 6; ++i) {
150 double dwave = (wave(m_times[i] - m_fittedTime + 0.5 * timestep) - wave(m_times[i] - m_fittedTime - 0.5 * timestep)) / timestep;
151 dwave_sq += dwave * dwave;
152 }
153 double sn_sq = 0;
154 for (unsigned int row = 0; row < m_data.size(); ++row) {
155 double sn = m_fittedAmplitudes[row] / m_noises[row];
156 sn_sq += sn * sn;
157 }
158 m_fittedTimeError = 1.0 / sqrt(sn_sq * dwave_sq);
159}
std::vector< double > m_fittedAmplitudeErrors
Fitted amplitude errors.
Definition: WaveFitter.h:221
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ calculateFittedData()

void calculateFittedData ( )
private

Calculate fitted data.

Definition at line 161 of file WaveFitter.cc.

162{
163 m_fittedData.resize(m_data.size());
164 for (int i = 0; i < 6; ++i) {
165 double wave = m_wave(m_times[i] - m_fittedTime);
166 for (unsigned int row = 0; row < m_data.size(); ++row)
167 m_fittedData[row][i] = m_fittedAmplitudes[row] * wave;
168 }
169}
strip_data_type m_fittedData
Fitted values from current fit.
Definition: WaveFitter.h:224

◆ Chi()

double Chi ( double  t)

Get standardized chi2 for a given time shift.

Parameters
ttime shift
Returns
chi2/ndf at a given time shift

Definition at line 24 of file WaveFitter.cc.

25{
26 // This we only need to calculate once
27 array<double, 6> waves;
28 double sumTemplateByTemplate = 0;
29 for (int i = 0; i < 6; ++i) {
30 waves[i] = m_wave(m_times[i] - t0);
31 sumTemplateByTemplate += waves[i] * waves[i];
32 }
33 // Calculate amplitudes
34 vector<double> amplitudes(m_data.size());
35 for (unsigned int row = 0; row < m_data.size(); ++row) {
36 double sumDataByTemplate = 0;
37 for (int i = 0; i < 6; i++) {
38 sumDataByTemplate += m_data[row][i] * waves[i];
39 }
40 amplitudes[row] = sumDataByTemplate / sumTemplateByTemplate;
41 }
42 double result = 0;
43 for (unsigned int row = 0; row < m_data.size(); ++row)
44 for (int i = 0; i < 6; i++) {
45 double diff =
46 (m_data[row][i] - amplitudes[row] * waves[i]) / m_noises[row];
47 result += diff * diff;
48 }
49 return result / m_ndf;
50}

◆ doFit()

void doFit ( )
private

Perform fit on the data.

Definition at line 73 of file WaveFitter.cc.

74{
75 // Degrees of freedom:
76 // # samples - # strips (amplitudes) - 1 (timeshift) - 1 (sum = 0)
77 m_ndf = 6 * m_data.size() - m_data.size() - 2;
78 if (m_ndf < 1) { // We have no data, so there'll be no fit and an error
79 B2ERROR("SVD::WaveFitterter : Fit required with no data provided.");
80 m_fittedTime = 0;
82 m_fittedAmplitudes.clear();
84 m_fittedLik = 1.0e10;
85 m_fittedData.clear();
86 // We don't need to repeat this
87 m_hasFit = true;
88 return;
89 }
90 // Normal operation: minimize wrt. t0 in (-dt,0) and (0,dt),
91 // take the solution giving a better fit
92 auto result_low = boost::math::tools::brent_find_minima(
93 [this](double t)->double { return negLogLikelihood(t); }, -m_dt, 0.0, s_minimizer_precision
94 );
95 double t0_low = result_low.first;
96 double lik_low = result_low.second;
97 auto result_hi = boost::math::tools::brent_find_minima(
98 [this](double t)->double { return negLogLikelihood(t); }, 0.0, m_dt, s_minimizer_precision
99 );
100 double t0_hi = result_hi.first;
101 double lik_hi = result_hi.second;
102 if (lik_low < lik_hi) {
103 m_fittedTime = t0_low;
104 m_fittedLik = lik_low;
105 } else {
106 m_fittedTime = t0_hi;
107 m_fittedLik = lik_hi;
108 }
112 m_hasFit = true;
113}
static int s_minimizer_precision
minimizer precision
Definition: WaveFitter.h:195
void calculateFitErrors()
Calculate fit errors.
Definition: WaveFitter.cc:134
double negLogLikelihood(double t)
Get negative profile log-likelihood for a given time shift.
Definition: WaveFitter.cc:52
void calculateFittedData()
Calculate fitted data.
Definition: WaveFitter.cc:161
void calculateAmplitudes()
Calculate amplitudes.
Definition: WaveFitter.cc:115

◆ getData()

const strip_data_type & getData ( ) const
inline

Get original fitted data.

Returns
vector of arrays of strip data used in the fit.

Definition at line 151 of file WaveFitter.h.

151{ return m_data; }

◆ getDt()

double getDt ( ) const
inline

Get sampling interval.

Definition at line 163 of file WaveFitter.h.

163{ return m_dt; }

◆ getFitData()

const strip_data_type & getFitData ( )
inline

Get fitted values for all samples.

Returns
vector of arrays of fit estimates of fitted data.

Definition at line 156 of file WaveFitter.h.

157 {
158 if (!m_hasFit) doFit();
159 return m_fittedData;
160 }
void doFit()
Perform fit on the data.
Definition: WaveFitter.cc:73

◆ getFitLikelihood()

double getFitLikelihood ( )
inline

Return negative log-likelihood for the fit.

Returns
negative log-likelihood for the fit

Definition at line 142 of file WaveFitter.h.

143 {
144 if (!m_hasFit) doFit();
145 return m_fittedLik;
146 }

◆ getFittedAmplitudeErrors()

const std::vector< double > & getFittedAmplitudeErrors ( )
inline

Retrieve errors of strip amplitudes after a fit.

Returns
vector<double> of errors of strip amplitude estimates

Definition at line 133 of file WaveFitter.h.

134 {
135 if (!m_hasFit) doFit();
137 }

◆ getFittedAmplitudes()

const std::vector< double > & getFittedAmplitudes ( )
inline

Retrieve strip amplitudes after a fit.

Returns
vector<double> of strip amplitude estimates

Definition at line 124 of file WaveFitter.h.

125 {
126 if (!m_hasFit) doFit();
127 return m_fittedAmplitudes;
128 }

◆ getFittedTime()

double getFittedTime ( )
inline

Retrieve fitted time shift value after a fit.

Returns
time shift estimate

Definition at line 106 of file WaveFitter.h.

107 {
108 if (!m_hasFit) doFit();
109 return m_fittedTime;
110 }

◆ getFittedTimeError()

double getFittedTimeError ( )
inline

Retrieve error of time shift value after a fit.

Returns
time shift estimate error

Definition at line 115 of file WaveFitter.h.

116 {
117 if (!m_hasFit) doFit();
118 return m_fittedTimeError;
119 }

◆ getTimes()

const std::array< double, 6 > & getTimes ( ) const
inline

Get sample times.

Definition at line 165 of file WaveFitter.h.

165{ return m_times; }

◆ integral12()

double integral12 ( double  lower,
double  upper,
std::function< double(double)>  f 
)
private

High-order Gauss-Legendre quadrature for likelihood integrals.

Definition at line 171 of file WaveFitter.cc.

172{
173 const unsigned int half_order = 6;
174 const double knots[half_order] = {
175 0.1252334085114689154724414,
176 0.3678314989981801937526915,
177 0.5873179542866174472967024,
178 0.7699026741943046870368938,
179 0.9041172563704748566784659,
180 0.9815606342467192506905491
181 };
182 const double weights[half_order] = {
183 0.2491470458134027850005624,
184 0.2334925365383548087608499,
185 0.2031674267230659217490645,
186 0.1600783285433462263346525,
187 0.1069393259953184309602547,
188 0.0471753363865118271946160
189 };
190 double span = 0.5 * (upper - lower);
191 double center = 0.5 * (upper + lower);
192 double result = 0;
193 for (unsigned int iknot = 0; iknot < half_order; ++iknot) {
194 result += weights[iknot] * f(center + span * knots[iknot]);
195 result += weights[iknot] * f(center - span * knots[iknot]);
196 }
197 result *= span;
198 return result;
199}

◆ integral20()

double integral20 ( double  lower,
double  upper,
std::function< double(double)>  f 
)
private

High-order Gauss-Legendre quadrature for likelihood integrals.

Definition at line 201 of file WaveFitter.cc.

202{
203 const unsigned int half_order = 10;
204 const double knots[half_order] = {
205 0.0765265211,
206 0.2277858511,
207 0.3737060887,
208 0.510867002,
209 0.6360536807,
210 0.7463319065,
211 0.8391169718,
212 0.9122344283,
213 0.9639719273,
214 0.9931285992
215 };
216 const double weights[half_order] = {
217 0.1527533871,
218 0.1491729865,
219 0.1420961093,
220 0.1316886384,
221 0.118194532,
222 0.1019301198,
223 0.0832767416,
224 0.0626720483,
225 0.0406014298,
226 0.0176140071
227 };
228 double span = 0.5 * (upper - lower);
229 double center = 0.5 * (upper + lower);
230 double result = 0;
231 for (unsigned int iknot = 0; iknot < half_order; ++iknot) {
232 result += weights[iknot] * f(center + span * knots[iknot]);
233 result += weights[iknot] * f(center - span * knots[iknot]);
234 }
235 result *= span;
236 return result;
237}

◆ lrSignal()

double lrSignal ( double  a,
double  b 
)

Calculate likelihood ratio Lsignal/Lbackground for acceptance (signal) window given as arguments.

Parameters
alower bound of the acceptance window
bupper bound of the acceptance window
Returns
L(in acc. w.)/L(out of acc. w.)

Definition at line 257 of file WaveFitter.cc.

258{
259 const double outer = m_dt;
260 if (!m_hasFit) doFit();
261 auto integrand = [this](double t)->double { return exp(-negLogLikelihood(t));};
262 // Calculate total integral and integral between a and b:
263 double int_total = integral20(-outer, 0, integrand);
264 int_total += integral20(0, outer, integrand);
265 double int_inner = 0;
266 if (low * high < 0) {
267 int_inner += integral20(low, 0, integrand);
268 int_inner += integral20(0, high, integrand);
269 } else {
270 int_inner += integral20(low, high, integrand);
271 }
272 double int_outer = int_total - int_inner;
273 double LR = int_inner / int_outer;
274 //cout << L_background << " <><> " << L_signal << " <> " << prob << endl;
275 return LR;
276}
double integral20(double lower, double upper, std::function< double(double)> f)
High-order Gauss-Legendre quadrature for likelihood integrals.
Definition: WaveFitter.cc:201

◆ negLogLikelihood()

double negLogLikelihood ( double  t)

Get negative profile log-likelihood for a given time shift.

Parameters
tthe time shift at which to calculate the likelihood
Returns
log-likelihood integrated over amplitudes at the given t

Definition at line 52 of file WaveFitter.cc.

53{
54 // wave values
55 array<double, 6> w;
56 for (int i = 0; i < 6; ++i) w[i] = m_wave(m_times[i] - t);
57 double sum_w2 = 0;
58 for (int i = 0; i < 6; ++i) sum_w2 += w[i] * w[i];
59 double lognorm = 0.5 * m_data.size() * log(sum_w2);
60 double sum_ksi = 0.0;
61 for (unsigned int row = 0; row < m_data.size(); ++row) {
62 double ksi = 0;
63 for (int i = 1; i < 6; ++i)
64 for (int j = 0; j < i; ++j) {
65 double cross = w[i] * m_data[row][j] - w[j] * m_data[row][i];
66 ksi += cross * cross;
67 }
68 sum_ksi += ksi;
69 }
70 return lognorm + 0.5 * sum_ksi / sum_w2;
71}

◆ pSignal()

double pSignal ( )

Calculate probability of the hit being a signal, with given priors for signal and background.

Definition at line 239 of file WaveFitter.cc.

240{
241 const double outer = m_dt;
242 const double inner = 5; // ns
243 if (!m_hasFit) doFit();
244 auto integrand = [this](double t)->double { return exp(-negLogLikelihood(t));};
245
246 double i1 = integral20(-outer, -inner, integrand);
247 double i2 = integral20(-inner, 0, integrand);
248 double i3 = integral20(0, inner, integrand);
249 double i4 = integral20(inner, outer, integrand);
250 double L_background = 1.0 / (2 * outer) * (i1 + i2 + i3 + i4);
251 double L_signal = 1.0 / (2 * inner) * (i2 + i3);
252 double prob = (L_signal / (L_background + L_signal));
253 //cout << L_background << " <><> " << L_signal << " <> " << prob << endl;
254 return prob;
255}

◆ reset()

void reset ( )
inline

Reset fitter data.

Definition at line 92 of file WaveFitter.h.

93 {
94 m_data.clear();
95 m_noises.clear();
96 m_fittedData.clear();
97 m_fittedAmplitudes.clear();
99 m_ndf = 0;
100 m_hasFit = false;
101 }

◆ setPrecision()

void setPrecision ( int  ndigits)
static

Set precision of the minimzer.

Definition at line 22 of file WaveFitter.cc.

22{ s_minimizer_precision = ndigits; }

◆ wave()

double wave ( double  t) const
inline

Get unit waveform value at a given time.

Definition at line 188 of file WaveFitter.h.

188{ return m_wave(t); }

Member Data Documentation

◆ m_data

strip_data_type m_data
private

Vector of sextets of APV samples.

Definition at line 214 of file WaveFitter.h.

◆ m_dt

double m_dt
private

Time interval between samples.

Definition at line 213 of file WaveFitter.h.

◆ m_fittedAmplitudeErrors

std::vector<double> m_fittedAmplitudeErrors
private

Fitted amplitude errors.

Definition at line 221 of file WaveFitter.h.

◆ m_fittedAmplitudes

std::vector<double> m_fittedAmplitudes
private

Fitted amplitudes of the sample.

Definition at line 220 of file WaveFitter.h.

◆ m_fittedData

strip_data_type m_fittedData
private

Fitted values from current fit.

Definition at line 224 of file WaveFitter.h.

◆ m_fittedLik

double m_fittedLik
private

Chi2 from the current fit.

Definition at line 222 of file WaveFitter.h.

◆ m_fittedTime

double m_fittedTime
private

Time from interval (-dt,+dt).

Definition at line 218 of file WaveFitter.h.

◆ m_fittedTimeError

double m_fittedTimeError
private

Error estimate of the fitted time.

Definition at line 219 of file WaveFitter.h.

◆ m_hasFit

bool m_hasFit
private

Are fit results available?

Definition at line 212 of file WaveFitter.h.

◆ m_ndf

int m_ndf
private

Degrees of freedom of the fit.

Definition at line 223 of file WaveFitter.h.

◆ m_noises

std::vector<double> m_noises
private

RMS noise for strips.

Definition at line 215 of file WaveFitter.h.

◆ m_times

std::array<double, 6> m_times
private

Vector of sample times.

Definition at line 216 of file WaveFitter.h.

◆ m_wave

wave_function_type m_wave
private

Waveform function.

Definition at line 217 of file WaveFitter.h.

◆ s_minimizer_precision

int s_minimizer_precision = 6
staticprivate

minimizer precision

Definition at line 195 of file WaveFitter.h.


The documentation for this class was generated from the following files: