Belle II Software  release-08-01-10
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 
17 namespace Belle2 {
22  namespace SVD {
23 
25  class DefaultWave {
26  public:
27 
31  DefaultWave(double tau = s_default_tau) : m_tau(tau) {}
32 
33  static double s_default_tau;
41  inline double getValue(double t) const
42  {
43  if (t < 0)
44  return 0.0;
45  else {
46  double z = t / m_tau;
47  return z * std::exp(1.0 - z);
48  }
49  }
50 
52  double getTau() const { return m_tau; }
53 
54  private:
55  double m_tau;
56  };
57 
58 
60  class WaveFitter {
61  public:
62 
64  typedef std::function<double(double)> wave_function_type;
66  typedef std::vector< std::array<double, 6> > strip_data_type;
67 
72  WaveFitter(wave_function_type wave, std::array<double, 6> times):
73  m_hasFit(false), m_times(times), m_wave(wave), m_fittedTime(0), m_fittedTimeError(), m_fittedLik(1.0e10), m_ndf(0)
74  {
75  m_dt = m_times[1] - m_times[0];
76  }
77 
82  void addData(const std::array<double, 6>& inpData, double rmsNoise)
83  {
84  std::array<double, 6> v;
85  std::copy(inpData.begin(), inpData.end(), v.begin());
86  m_data.push_back(v);
87  m_noises.push_back(rmsNoise);
88  m_hasFit = false;
89  m_ndf = 6 * m_data.size() - m_data.size() - 2;
90  }
91 
93  void reset()
94  {
95  m_data.clear();
96  m_noises.clear();
97  m_fittedData.clear();
98  m_fittedAmplitudes.clear();
100  m_ndf = 0;
101  m_hasFit = false;
102  }
103 
107  double getFittedTime()
108  {
109  if (!m_hasFit) doFit();
110  return m_fittedTime;
111  }
112 
117  {
118  if (!m_hasFit) doFit();
119  return m_fittedTimeError;
120  }
121 
125  const std::vector<double>& getFittedAmplitudes()
126  {
127  if (!m_hasFit) doFit();
128  return m_fittedAmplitudes;
129  }
130 
134  const std::vector<double>& getFittedAmplitudeErrors()
135  {
136  if (!m_hasFit) doFit();
138  }
139 
144  {
145  if (!m_hasFit) doFit();
146  return m_fittedLik;
147  }
148 
152  const strip_data_type& getData() const { return m_data; }
153 
158  {
159  if (!m_hasFit) doFit();
160  return m_fittedData;
161  }
162 
164  double getDt() const { return m_dt; }
166  const std::array<double, 6>& getTimes() const { return m_times; }
167 
172  double Chi(double t);
177  double negLogLikelihood(double t);
180  double pSignal();
186  double lrSignal(double a, double b);
187 
189  double wave(double t) const { return m_wave(t); }
191  static void setPrecision(int ndigits);
192 
193  private:
194 
197 
199  void doFit();
201  void calculateAmplitudes();
205  void calculateFitErrors();
207  void calculateFittedData();
209  // cppcheck-suppress unusedPrivateFunction
210  double integral12(double lower, double upper, std::function<double(double)> f);
212  double integral20(double lower, double upper, std::function<double(double)> f);
213  bool m_hasFit;
214  double m_dt;
216  std::vector<double> m_noises;
217  std::array<double, 6> m_times;
219  double m_fittedTime;
221  std::vector<double> m_fittedAmplitudes;
222  std::vector<double> m_fittedAmplitudeErrors;
223  double m_fittedLik;
224  int m_ndf;
226  };
227 
228  } // namespace SVD
230 } // namespace Belle2
231 #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:52
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:41
double m_tau
Waveform decay time.
Definition: WaveFitter.h:55
Waveform fitter class.
Definition: WaveFitter.h:60
double m_fittedTimeError
Error estimate of the fitted time.
Definition: WaveFitter.h:220
const std::array< double, 6 > & getTimes() const
Get sample times.
Definition: WaveFitter.h:166
std::vector< std::array< double, 6 > > strip_data_type
Function type for cluster data.
Definition: WaveFitter.h:66
std::vector< double > m_fittedAmplitudes
Fitted amplitudes of the sample.
Definition: WaveFitter.h:221
strip_data_type m_fittedData
Fitted values from current fit.
Definition: WaveFitter.h:225
static int s_minimizer_precision
minimizer precision
Definition: WaveFitter.h:196
double getFittedTime()
Retrieve fitted time shift value after a fit.
Definition: WaveFitter.h:107
std::vector< double > m_noises
RMS noise for strips.
Definition: WaveFitter.h:216
const std::vector< double > & getFittedAmplitudes()
Retrieve strip amplitudes after a fit.
Definition: WaveFitter.h:125
double getDt() const
Get sampling interval.
Definition: WaveFitter.h:164
const strip_data_type & getFitData()
Get fitted values for all samples.
Definition: WaveFitter.h:157
double getFitLikelihood()
Return negative log-likelihood for the fit.
Definition: WaveFitter.h:143
strip_data_type m_data
Vector of sextets of APV samples.
Definition: WaveFitter.h:215
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:222
double m_dt
Time interval between samples.
Definition: WaveFitter.h:214
double wave(double t) const
Get unit waveform value at a given time.
Definition: WaveFitter.h:189
void calculateFitErrors()
Calculate fit errors.
Definition: WaveFitter.cc:134
double getFittedTimeError()
Retrieve error of time shift value after a fit.
Definition: WaveFitter.h:116
void addData(const std::array< double, 6 > &inpData, double rmsNoise)
Add strip data to the fitter.
Definition: WaveFitter.h:82
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:224
static void setPrecision(int ndigits)
Set precision of the minimzer.
Definition: WaveFitter.cc:22
const std::vector< double > & getFittedAmplitudeErrors()
Retrieve errors of strip amplitudes after a fit.
Definition: WaveFitter.h:134
double m_fittedTime
Time from interval (-dt,+dt).
Definition: WaveFitter.h:219
wave_function_type m_wave
Waveform function.
Definition: WaveFitter.h:218
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:72
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:217
const strip_data_type & getData() const
Get original fitted data.
Definition: WaveFitter.h:152
double m_fittedLik
Chi2 from the current fit.
Definition: WaveFitter.h:223
void reset()
Reset fitter data.
Definition: WaveFitter.h:93
std::function< double(double)> wave_function_type
Function type for waveform.
Definition: WaveFitter.h:64
bool m_hasFit
Are fit results available?
Definition: WaveFitter.h:213
double Chi(double t)
Get standardized chi2 for a given time shift.
Definition: WaveFitter.cc:24
Abstract base class for different kinds of events.