Belle II Software development
WaveFitter.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8#include <framework/logging/Logger.h>
9#include <svd/reconstruction/WaveFitter.h>
10#include <iostream>
11#include <math.h>
12#include <utility>
13#include <boost/math/tools/minima.hpp>
14
15using namespace Belle2::SVD;
16using namespace std;
17
19
21
22void WaveFitter::setPrecision(int ndigits) { s_minimizer_precision = ndigits; }
23
24double WaveFitter::Chi(double t0)
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}
51
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}
72
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}
114
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}
133
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}
160
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}
170
171double WaveFitter::integral12(double lower, double upper, std::function<double(double)> f)
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}
200
201double WaveFitter::integral20(double lower, double upper, std::function<double(double)> f)
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}
238
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}
256
257double WaveFitter::lrSignal(double low, double high)
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}
277
static double s_default_tau
Default waveform decay time.
Definition: WaveFitter.h:33
double m_fittedTimeError
Error estimate of the fitted time.
Definition: WaveFitter.h:219
std::vector< double > m_fittedAmplitudes
Fitted amplitudes of the sample.
Definition: WaveFitter.h:220
strip_data_type m_fittedData
Fitted values from current fit.
Definition: WaveFitter.h:224
static int s_minimizer_precision
minimizer precision
Definition: WaveFitter.h:195
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
double pSignal()
Calculate probability of the hit being a signal, with given priors for signal and background.
Definition: WaveFitter.cc:239
void doFit()
Perform fit on the data.
Definition: WaveFitter.cc:73
std::vector< double > m_fittedAmplitudeErrors
Fitted amplitude errors.
Definition: WaveFitter.h:221
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
void calculateFitErrors()
Calculate fit errors.
Definition: WaveFitter.cc:134
double integral12(double lower, double upper, std::function< double(double)> f)
High-order Gauss-Legendre quadrature for likelihood integrals.
Definition: WaveFitter.cc:171
double negLogLikelihood(double t)
Get negative profile log-likelihood for a given time shift.
Definition: WaveFitter.cc:52
double integral20(double lower, double upper, std::function< double(double)> f)
High-order Gauss-Legendre quadrature for likelihood integrals.
Definition: WaveFitter.cc:201
int m_ndf
Degrees of freedom of the fit.
Definition: WaveFitter.h:223
static void setPrecision(int ndigits)
Set precision of the minimzer.
Definition: WaveFitter.cc:22
double m_fittedTime
Time from interval (-dt,+dt).
Definition: WaveFitter.h:218
wave_function_type m_wave
Waveform function.
Definition: WaveFitter.h:217
void calculateFittedData()
Calculate fitted data.
Definition: WaveFitter.cc:161
void calculateAmplitudes()
Calculate amplitudes.
Definition: WaveFitter.cc:115
double lrSignal(double a, double b)
Calculate likelihood ratio Lsignal/Lbackground for acceptance (signal) window given as arguments.
Definition: WaveFitter.cc:257
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
double Chi(double t)
Get standardized chi2 for a given time shift.
Definition: WaveFitter.cc:24
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Namespace to encapsulate code needed for simulation and reconstrucion of the SVD.
Definition: GeoSVDCreator.h:23
STL namespace.