Belle II Software prerelease-11-00-00a
TOPLocalCalFitter.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
9#pragma once
10#include <array>
11#include <string>
12#include <calibration/CalibrationAlgorithm.h>
13
14class TH1;
15
16namespace Belle2 {
21 namespace TOP {
22
31 public:
32
35
37 ~TOPLocalCalFitter() override;
38
41
44
46 void setMinEntries(int minEntries)
47 {
48 m_minEntries = minEntries;
49 }
50
52 void setOutputFileName(const std::string& output)
53 {
54 m_output = output;
55 }
56
60 void setFitConstraintsFileName(const std::string& fitConstraints)
61 {
62 m_fitConstraints = fitConstraints;
63 }
64
66 void setTTSFileName(const std::string& TTSData)
67 {
68 m_TTSData = TTSData;
69 }
70
79 void setFitMode(const std::string& fitterMode)
80 {
81 if (fitterMode == "calibration")
82 B2INFO("Fitter set to calibration mode");
83 else if (fitterMode == "monitoring")
84 B2INFO("Fitter set to monitoring mode");
85 else if (fitterMode == "MC")
86 B2INFO("Fitter set to MC mode");
87 else
88 B2ERROR("Unknown fitter type " << fitterMode << ". The valid options are calibration, monitoring or MC");
89
90 m_fitterMode = fitterMode;
91 }
92
94 void fitInAmpliduteBins(bool isFitInAmplitudeBins)
95 {
96 m_isFitInAmplitudeBins = isFitInAmplitudeBins;
97 }
98
99 protected:
100
103
105 void loadMCInfoTrees();
106
108 void determineFitStatus();
109
111 void fitChannel(short slot, short channel, TH1* h);
112
120 void fitChannel(short slot, short channel, TH1* h, bool inBins, double frac);
121
123 void fitPulser(TH1*, TH1*);
124
128 void calculateChannelT0();
129
132 EResult calibrate() override;
133
134 private:
135
136 int m_minEntries = 50;
137 std::string m_output = "laserFitResult.root";
138 std::string m_fitConstraints =
139 "/group/belle2/group/detector/TOP/calibration/MCreferences/LaserMCParameters.root";
140 std::string m_TTSData =
141 "/group/belle2/group/detector/TOP/calibration/MCreferences/TTSParametrization.root";
142 std::string m_fitterMode = "calibration";
144 std::vector<float> m_binEdges = {50, 100, 150, 200, 250, 300, 350, 400, 500, 600, 800, 1000, 1500, 2000};
145 TFile* m_inputTTS = nullptr;
146 TFile* m_inputConstraints = nullptr;
147 TTree* m_treeTTS = nullptr;
149 nullptr;
150 TFile* m_histFile = nullptr;
151 TTree* m_fitTree = nullptr;
153 nullptr;
155
156 bool m_detectCrosstalk = false;
157
158 TTree* m_crosstalkTree = nullptr;
159 TTree* m_fitTree_noXtalk = nullptr;
160
161 // Per-candidate (crosstalk) variables
162 short m_sl0 = -1;
163 short m_sl1 = -1;
164 short m_ch0 = -1;
165 short m_ch1 = -1;
166 float m_ht0 = NAN;
167 float m_ht1 = NAN;
168 float m_a0 = NAN;
169 float m_a1 = NAN;
170 float m_w0 = NAN;
171 float m_w1 = NAN;
172 float m_q0 = NAN;
173 float m_q1 = NAN;
174 float m_f_q0 = NAN;
175
176 // Variables for the TTS parametrization tree
177 float m_mean2 = 0;
178 float m_sigma1 = 0;
179 float m_sigma2 = 0;
180 float m_f1 = 0;
181 float m_f2 = 0;
182 short m_pixelRow = 0;
183 short m_pixelCol = 0;
184
185 // Variables for the MC truth infos
195
196
197 // Variables for the output tree
198 float m_binLowerEdge = 0;
199 float m_binUpperEdge = 0;
200 short m_channel = 0;
201 short m_slot = 0;
202 short m_row = 0;
203 short m_col = 0;
204 float m_peakTime = 0;
205 float m_deltaT =
206 0;
207 float m_sigma = 0.;
208 float m_fraction = 0.;
209 float m_yieldLaser = 0.;
210 float m_histoIntegral = 0.;
211
212 float m_peakTimeErr = 0;
213 float m_deltaTErr = 0;
214 float m_sigmaErr = 0.;
215 float m_fractionErr = 0.;
216 float m_yieldLaserErr = 0.;
217
218 float m_timeExtra = 0.;
219 float m_sigmaExtra = 0.;
220 float m_yieldLaserExtra = 0.;
221 float m_alphaExtra = 0.;
222 float m_nExtra = 0.;
223
224 float m_timeBackground = 0.;
225 float m_sigmaBackground = 0.;
227
228 float m_fractionMC = 0.;
229 float m_deltaTMC = 0.;
231 0.;
233
234 float m_chi2 = 0;
235 float m_rms = 0;
236
238 0.;
242 float m_channelT0Err = 0.;
243
244 float m_firstPulserTime = 0.;
246
248 0.;
250
251 short m_fitStatus = 1;
252
253 double m_width = 0;
254 double m_amplitude = 0;
255
256 // Row/Col per (slot, channel). Slots are 0..15, channels 0..511 (0-based inside code).
257 std::array<std::array<short, 512>, 16> m_rowOf{};
258 std::array<std::array<short, 512>, 16> m_colOf{};
259 bool m_hasChannelMaps{false};
260
262 inline short rowOf(short slot, short ch) const noexcept
263 {
264 return (slot >= 0 && slot < 16 && ch >= 0 && ch < 512) ? m_rowOf[slot][ch] : short(-1);
265 }
266
267 inline short colOf(short slot, short ch) const noexcept
268 {
269 return (slot >= 0 && slot < 16 && ch >= 0 && ch < 512) ? m_colOf[slot][ch] : short(-1);
270 }
271
276 inline bool areNeighbors(short slot, short a, short b, int drMax = 1, int dcMax = 1) const noexcept
277 {
278 const int dr = std::abs(rowOf(slot, a) - rowOf(slot, b));
279 const int dc = std::abs(colOf(slot, a) - colOf(slot, b));
280 return (dr + dc > 0) && (dr <= drMax) && (dc <= dcMax);
281 }
282
284 void buildChannelMaps();
285
286 };
287 } // namespace TOP
289} // namespace Belle2
EResult
The result of calibration.
CalibrationAlgorithm(const std::string &collectorModuleName)
Constructor - sets the prefix for collected objects (won't be accesses until execute(....
short m_fitStatus
Fit quality flag, propagated to the constants.
short colOf(short slot, short ch) const noexcept
Column index for (slot,channel), or -1 if out of bounds.
TTree * m_crosstalkTree
Output tree for crosstalk candidates.
float m_yieldLaserErr
Statistical error on yield.
float m_binLowerEdge
Lower edge of the amplitude bin in which this fit is performed.
float m_nExtraConstraints
parameter n of the tail of the extra peak
float m_chi2
Reduced chi2 of the fit.
float m_timeExtraConstraints
Position of the gaussian used to describe the extra peak on the timing distribution tail.
void fitChannel(short slot, short channel, TH1 *h)
Fits the laser light on one channel.
int m_minEntries
Minimum number of entries to perform the fit.
float m_f1
Fraction of the first gaussian on the TTS parametrization.
float m_fraction
Fraction of events in the secondary peak.
float m_sigmaExtra
Gaussian sigma of the extra peak in the timing tail.
float m_secondPulserSigma
Time resolution from the fit of the first electronic pulse, from a Gaussian fit.
float m_yieldLaserBackground
Integral of the background gaussian.
float m_sigma
Gaussian time resolution, fitted.
void setMinEntries(int minEntries)
Sets the minimum number of entries to perform the calibration in one channel.
float m_peakTimeMC
Time of the main peak in the MC simulation, i.e.
std::string m_output
Name of the output file.
float m_q0
Integrated charge for channel 0 in pair.
float m_timeBackground
Position of the gaussian used to describe the background, w/ respect to peakTime.
void loadMCInfoTrees()
loads the TTS parameters and the MC truth info
float m_deltaTMC
Time difference between the main peak and the secondary peak in the MC simulation.
float m_peakTimeErr
Statistical error on peakTime.
short m_channel
Channel number (0-511)
float m_channelT0Err
Statistical error on channelT0.
bool m_detectCrosstalk
Enables the crosstalk detection algorithm.
std::string m_TTSData
File with the TTS parametrization.
std::vector< float > m_binEdges
Amplitude bins.
float m_sigmaBackground
Sigma of the gaussian used to describe the background.
float m_peakTime
Fitted time of the main (i.e.
float m_fractionErr
Statistical error on fraction.
bool m_isFitInAmplitudeBins
Enables the fit in amplitude bins.
TOPLocalCalFitter & operator=(const TOPLocalCalFitter &)=delete
Assignment operator (disabled)
TOPLocalCalFitter(const TOPLocalCalFitter &)=delete
Copy constructor (disabled)
float m_alphaExtra
alpha parameter of the tail of the extra peak.
float m_a0
Amplitude for channel 0 in pair.
float m_rms
RMS of the histogram used for the fit.
float m_mean2
Position of the second gaussian of the TTS parametrization with respect to the first one.
std::string m_fitConstraints
File with the Fit constraints.
bool m_hasChannelMaps
Flag indicating if channel->(row,col) maps have been built.
TTree * m_timewalkTree
Output of the fitter.
void setFitConstraintsFileName(const std::string &fitConstraints)
Sets the name of the root file containing the laser MC time corrections and the fit constraints.
float m_w0
Width for channel 0 in pair.
float m_yieldLaser
Total number of laser hits from the fitting function integral.
bool areNeighbors(short slot, short a, short b, int drMax=1, int dcMax=1) const noexcept
Return true if channels a and b are neighbors on the same slot in row/col space.
TFile * m_inputTTS
File containing m_treeTTS.
float m_nExtra
parameter n of the tail of the extra peak
float m_ht1
Hit time for channel 1 in pair.
void fitInAmpliduteBins(bool isFitInAmplitudeBins)
Enables the fit amplitude bins.
float m_alphaExtraConstraints
alpha parameter of the tail of the extra peak.
float m_channelT0
Raw, channelT0 calibration, defined as peakTime-peakTimeMC.
short m_ch1
Channel number (0-511)
float m_sigmaErr
Statistical error on sigma.
float m_binUpperEdge
Upper edge of the amplitude bin in which this fit is performed.
float m_f_q0
Fraction of charge on channel 0 in pair.
TFile * m_histFile
Output of the fitter.
TTree * m_treeConstraints
Input to the fitter.
float m_sigmaExtraConstraints
Width of the gaussian used to describe the extra peak on the timing distribution tail.
float m_w1
Width for channel 1 in pair.
float m_firstPulserSigma
Time resolution from the fit of the first electronic pulse, from a Gaussian fit.
float m_ht0
Hit time for channel 0 in pair.
std::array< std::array< short, 512 >, 16 > m_colOf
Column index for (slot,channel), or -1 if out of bounds.
short rowOf(short slot, short ch) const noexcept
Row index for (slot,channel), or -1 if out of bounds.
float m_q1
Integrated charge for channel 1 in pair.
float m_deltaT
Time difference between the main peak and the secondary peak.
TTree * m_fitTree
Output of the fitter.
float m_f2
Fraction of the second gaussian on the TTS parametrization.
void setOutputFileName(const std::string &output)
Sets the name of the output root file.
void determineFitStatus()
determines if the constant obtained by the fit are good or not
void buildChannelMaps()
Build (row,col) lookup tables from the TTS tree; call after opening m_treeTTS.
void setFitMode(const std::string &fitterMode)
Sets the fitter mode.
float m_secondPulserTime
Average time of the second electronic pulse respect to the reference pulse, from a gaussian fit.
float m_peakTimeConstraints
Time of the main laser peak in the MC simulation (aka MC correction)
float m_sigma1
Width of the first gaussian on the TTS parametrization.
float m_fractionMC
Fraction of events in the secondary peak form the MC simulation.
float m_sigmaBackgroundConstraints
Sigma of the gaussian used to describe the background.
TTree * m_treeTTS
Input to the fitter.
void setupOutputTreeAndFile()
prepares the output tree
~TOPLocalCalFitter() override
Destructor.
TTree * m_fitTree_noXtalk
Output tree for non-crosstalk candidates.
float m_timeBackgroundConstraints
Position of the gaussian used to describe the background, w/ respect to peakTime.
short m_ch0
Channel number (0-511)
TFile * m_inputConstraints
File containing m_treeConstraints.
float m_deltaTErr
Statistical error on deltaT.
void calculateChannelT0()
Calculates the commonT0 calibration after the fits have been done.
std::array< std::array< short, 512 >, 16 > m_rowOf
Row index for (slot,channel), or -1 if out of bounds.
float m_firstPulserTime
Average time of the first electronic pulse respect to the reference pulse, from a Gaussian fit.
float m_histoIntegral
Integral of the fitted histogram.
float m_yieldLaserExtra
Integral of the extra peak.
float m_sigma2
Width of the second gaussian on the TTS parametrization.
float m_a1
Amplitude for channel 1 in pair.
void fitPulser(TH1 *, TH1 *)
Fits the two pulsers.
EResult calibrate() override
Runs the algorithm on events.
float m_timeExtra
Position of the extra peak seen in the timing tail, w/ respect to peakTime.
float m_fractionConstraints
Fraction of the main peak.
float m_deltaTConstraints
Distance between the main and the secondary laser peak.
void setTTSFileName(const std::string &TTSData)
Sets the name of the root file containing the TTS parameters.
Abstract base class for different kinds of events.