8 #include <framework/logging/Logger.h>
9 #include <svd/reconstruction/WaveFitter.h>
13 #include <boost/math/tools/minima.hpp>
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];
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];
40 amplitudes[row] = sumDataByTemplate / sumTemplateByTemplate;
43 for (
unsigned int row = 0; row < m_data.size(); ++row)
44 for (
int i = 0; i < 6; i++) {
46 (m_data[row][i] - amplitudes[row] * waves[i]) / m_noises[row];
47 result += diff * diff;
49 return result / m_ndf;
56 for (
int i = 0; i < 6; ++i) w[i] = m_wave(m_times[i] - t);
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);
61 for (
unsigned int row = 0; row < m_data.size(); ++row) {
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];
70 return lognorm + 0.5 * sum_ksi / sum_w2;
77 m_ndf = 6 * m_data.size() - m_data.size() - 2;
79 B2ERROR(
"SVD::WaveFitterter : Fit required with no data provided.");
81 m_fittedTimeError = 0;
82 m_fittedAmplitudes.clear();
83 m_fittedAmplitudeErrors.clear();
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
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
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;
106 m_fittedTime = t0_hi;
107 m_fittedLik = lik_hi;
109 calculateAmplitudes();
110 calculateFitErrors();
111 calculateFittedData();
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];
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];
130 m_fittedAmplitudes.push_back(sumDataByTemplate / sumTemplateByTemplate);
138 for (
int i = 0; i < 6; ++i) {
139 double wave = m_wave(m_times[i] - m_fittedTime);
140 sum_w2 += wave * wave;
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]);
147 double timestep = 0.1 * m_dt;
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;
154 for (
unsigned int row = 0; row < m_data.size(); ++row) {
155 double sn = m_fittedAmplitudes[row] / m_noises[row];
158 m_fittedTimeError = 1.0 /
sqrt(sn_sq * dwave_sq);
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;
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
182 const double weights[half_order] = {
183 0.2491470458134027850005624,
184 0.2334925365383548087608499,
185 0.2031674267230659217490645,
186 0.1600783285433462263346525,
187 0.1069393259953184309602547,
188 0.0471753363865118271946160
190 double span = 0.5 * (upper - lower);
191 double center = 0.5 * (upper + lower);
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]);
203 const unsigned int half_order = 10;
204 const double knots[half_order] = {
216 const double weights[half_order] = {
228 double span = 0.5 * (upper - lower);
229 double center = 0.5 * (upper + lower);
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]);
241 const double outer = m_dt;
242 const double inner = 5;
243 if (!m_hasFit) doFit();
244 auto integrand = [
this](
double t)->
double {
return exp(-negLogLikelihood(t));};
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));
259 const double outer = m_dt;
260 if (!m_hasFit) doFit();
261 auto integrand = [
this](
double t)->
double {
return exp(-negLogLikelihood(t));};
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);
270 int_inner += integral20(low, high, integrand);
272 double int_outer = int_total - int_inner;
273 double LR = int_inner / int_outer;
static double s_default_tau
Default waveform decay time.
static int s_minimizer_precision
minimizer precision
double pSignal()
Calculate probability of the hit being a signal, with given priors for signal and background.
void doFit()
Perform fit on the data.
void calculateFitErrors()
Calculate fit errors.
double integral12(double lower, double upper, std::function< double(double)> f)
High-order Gauss-Legendre quadrature for likelihood integrals.
double negLogLikelihood(double t)
Get negative profile log-likelihood for a given time shift.
double integral20(double lower, double upper, std::function< double(double)> f)
High-order Gauss-Legendre quadrature for likelihood integrals.
static void setPrecision(int ndigits)
Set precision of the minimzer.
void calculateFittedData()
Calculate fitted data.
void calculateAmplitudes()
Calculate amplitudes.
double lrSignal(double a, double b)
Calculate likelihood ratio Lsignal/Lbackground for acceptance (signal) window given as arguments.
double Chi(double t)
Get standardized chi2 for a given time shift.
double sqrt(double a)
sqrt for double
Namespace to encapsulate code needed for simulation and reconstrucion of the SVD.