Belle II Software development
WaveFitter.h
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#ifndef WAVEFITTER_H
9#define WAVEFITTER_H
10
11#include <array>
12#include <vector>
13#include <algorithm>
14#include <functional>
15#include <cmath>
16
17namespace Belle2 {
22 namespace SVD {
23
26 public:
27
31 DefaultWave(double tau = s_default_tau) : m_tau(tau) {}
32
33 static double s_default_tau;
40 inline double getValue(double t) const
41 {
42 if (t < 0)
43 return 0.0;
44 else {
45 double z = t / m_tau;
46 return z * std::exp(1.0 - z);
47 }
48 }
49
51 double getTau() const { return m_tau; }
52
53 private:
54 double m_tau;
55 };
56
57
59 class WaveFitter {
60 public:
61
63 typedef std::function<double(double)> wave_function_type;
65 typedef std::vector< std::array<double, 6> > strip_data_type;
66
71 WaveFitter(wave_function_type wave, std::array<double, 6> times):
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 }
76
81 void addData(const std::array<double, 6>& inpData, double rmsNoise)
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 }
90
92 void reset()
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 }
102
107 {
108 if (!m_hasFit) doFit();
109 return m_fittedTime;
110 }
111
116 {
117 if (!m_hasFit) doFit();
118 return m_fittedTimeError;
119 }
120
124 const std::vector<double>& getFittedAmplitudes()
125 {
126 if (!m_hasFit) doFit();
127 return m_fittedAmplitudes;
128 }
129
133 const std::vector<double>& getFittedAmplitudeErrors()
134 {
135 if (!m_hasFit) doFit();
137 }
138
143 {
144 if (!m_hasFit) doFit();
145 return m_fittedLik;
146 }
147
151 const strip_data_type& getData() const { return m_data; }
152
157 {
158 if (!m_hasFit) doFit();
159 return m_fittedData;
160 }
161
163 double getDt() const { return m_dt; }
165 const std::array<double, 6>& getTimes() const { return m_times; }
166
171 double Chi(double t);
176 double negLogLikelihood(double t);
179 double pSignal();
185 double lrSignal(double a, double b);
186
188 double wave(double t) const { return m_wave(t); }
190 static void setPrecision(int ndigits);
191
192 private:
193
196
198 void doFit();
200 void calculateAmplitudes();
204 void calculateFitErrors();
206 void calculateFittedData();
208 // cppcheck-suppress unusedPrivateFunction
209 double integral12(double lower, double upper, std::function<double(double)> f);
211 double integral20(double lower, double upper, std::function<double(double)> f);
212 bool m_hasFit;
213 double m_dt;
215 std::vector<double> m_noises;
216 std::array<double, 6> m_times;
220 std::vector<double> m_fittedAmplitudes;
221 std::vector<double> m_fittedAmplitudeErrors;
222 double m_fittedLik;
223 int m_ndf;
225 };
226
227 } // namespace SVD
229} // namespace Belle2
230#endif // WAVEFITTER_H
A functor to provide a simple model of APV25 strip response.
Definition: WaveFitter.h:25
static double s_default_tau
Default waveform decay time.
Definition: WaveFitter.h:33
double getTau() const
return the decay time
Definition: WaveFitter.h:51
DefaultWave(double tau=s_default_tau)
Constructor takes waveform decay time as parameter.
Definition: WaveFitter.h:31
double getValue(double t) const
getValue() returns the value at desired time and time shift.
Definition: WaveFitter.h:40
double m_tau
Waveform decay time.
Definition: WaveFitter.h:54
Waveform fitter class.
Definition: WaveFitter.h:59
double m_fittedTimeError
Error estimate of the fitted time.
Definition: WaveFitter.h:219
const std::vector< double > & getFittedAmplitudes()
Retrieve strip amplitudes after a fit.
Definition: WaveFitter.h:124
std::vector< std::array< double, 6 > > strip_data_type
Function type for cluster data.
Definition: WaveFitter.h:65
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
double getFittedTime()
Retrieve fitted time shift value after a fit.
Definition: WaveFitter.h:106
std::vector< double > m_noises
RMS noise for strips.
Definition: WaveFitter.h:215
double getDt() const
Get sampling interval.
Definition: WaveFitter.h:163
const strip_data_type & getFitData()
Get fitted values for all samples.
Definition: WaveFitter.h:156
double getFitLikelihood()
Return negative log-likelihood for the fit.
Definition: WaveFitter.h:142
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
const strip_data_type & getData() const
Get original fitted data.
Definition: WaveFitter.h:151
std::vector< double > m_fittedAmplitudeErrors
Fitted amplitude errors.
Definition: WaveFitter.h:221
const std::vector< double > & getFittedAmplitudeErrors()
Retrieve errors of strip amplitudes after a fit.
Definition: WaveFitter.h:133
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
const std::array< double, 6 > & getTimes() const
Get sample times.
Definition: WaveFitter.h:165
void calculateFitErrors()
Calculate fit errors.
Definition: WaveFitter.cc:134
double getFittedTimeError()
Retrieve error of time shift value after a fit.
Definition: WaveFitter.h:115
void addData(const std::array< double, 6 > &inpData, double rmsNoise)
Add strip data to the fitter.
Definition: WaveFitter.h:81
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
WaveFitter(wave_function_type wave, std::array< double, 6 > times)
Constructor creates an empty fitter with default settings.
Definition: WaveFitter.h:71
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
void reset()
Reset fitter data.
Definition: WaveFitter.h:92
std::function< double(double)> wave_function_type
Function type for waveform.
Definition: WaveFitter.h:63
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
Abstract base class for different kinds of events.