Belle II Software  release-08-01-10
WaveFitter Class Reference

Waveform fitter class. More...

#include <WaveFitter.h>

Collaboration diagram for WaveFitter:

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. More...
 
void addData (const std::array< double, 6 > &inpData, double rmsNoise)
 Add strip data to the fitter. More...
 
void reset ()
 Reset fitter data.
 
double getFittedTime ()
 Retrieve fitted time shift value after a fit. More...
 
double getFittedTimeError ()
 Retrieve error of time shift value after a fit. More...
 
const std::vector< double > & getFittedAmplitudes ()
 Retrieve strip amplitudes after a fit. More...
 
const std::vector< double > & getFittedAmplitudeErrors ()
 Retrieve errors of strip amplitudes after a fit. More...
 
double getFitLikelihood ()
 Return negative log-likelihood for the fit. More...
 
const strip_data_typegetData () const
 Get original fitted data. More...
 
const strip_data_typegetFitData ()
 Get fitted values for all samples. More...
 
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. More...
 
double negLogLikelihood (double t)
 Get negative profile log-likelihood for a given time shift. More...
 
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. More...
 
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. More...
 
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 60 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 72 of file WaveFitter.h.

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 82 of file WaveFitter.h.

◆ 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 }
double m_fittedTimeError
Error estimate of the fitted time.
Definition: WaveFitter.h:220
std::vector< double > m_fittedAmplitudes
Fitted amplitudes of the sample.
Definition: WaveFitter.h:221
std::vector< double > m_noises
RMS noise for strips.
Definition: WaveFitter.h:216
strip_data_type m_data
Vector of sextets of APV samples.
Definition: WaveFitter.h:215
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
double m_fittedTime
Time from interval (-dt,+dt).
Definition: WaveFitter.h:219
wave_function_type m_wave
Waveform function.
Definition: WaveFitter.h:218
std::array< double, 6 > m_times
Vector of sample times.
Definition: WaveFitter.h:217
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ 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.

◆ getData()

const strip_data_type& getData ( void  ) const
inline

Get original fitted data.

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

Definition at line 152 of file WaveFitter.h.

◆ 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 157 of file WaveFitter.h.

◆ getFitLikelihood()

double getFitLikelihood ( )
inline

Return negative log-likelihood for the fit.

Returns
negative log-likelihood for the fit

Definition at line 143 of file WaveFitter.h.

◆ 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 134 of file WaveFitter.h.

◆ getFittedAmplitudes()

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

Retrieve strip amplitudes after a fit.

Returns
vector<double> of strip amplitude estimates

Definition at line 125 of file WaveFitter.h.

◆ getFittedTime()

double getFittedTime ( )
inline

Retrieve fitted time shift value after a fit.

Returns
time shift estimate

Definition at line 107 of file WaveFitter.h.

◆ getFittedTimeError()

double getFittedTimeError ( )
inline

Retrieve error of time shift value after a fit.

Returns
time shift estimate error

Definition at line 116 of file WaveFitter.h.

◆ 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.

◆ 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.


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