Belle II Software development
KLMTimeAlgorithm.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
11/* KLM headers. */
12#include <klm/bklm/geometry/GeometryPar.h>
13#include <klm/dataobjects/KLMChannelIndex.h>
14#include <klm/dataobjects/KLMElementNumbers.h>
15#include <klm/dbobjects/KLMTimeCableDelay.h>
16#include <klm/dbobjects/KLMTimeConstants.h>
17#include <klm/dbobjects/KLMTimeResolution.h>
18#include <klm/eklm/geometry/GeometryData.h>
19
20/* Basf2 headers. */
21#include <calibration/CalibrationAlgorithm.h>
22
23/* ROOT headers. */
24#include <Math/MinimizerOptions.h>
25#include <TF1.h>
26#include <TGraphErrors.h>
27#include <TH1F.h>
28#include <TH1I.h>
29#include <TH2F.h>
30#include <TProfile.h>
31
32/* C++ STL headers. */
33#include <functional>
34#include <map> // For std::map in readCalibrationDataCounts
35#include <utility> // For std::pair in readCalibrationDataFor2DFit
36#include <vector> // For std::vector in readCalibrationDataFor2DFit
37#include <functional> // ADD THIS - needed for std::function in readCalibrationDataBatch
38
39class TFile; // forward declaration
40
41namespace Belle2 {
46
51
52 public:
53
57 struct Event {
58
60 double t0 = 0;
61
63 double flyTime = 0;
64
66 double recTime = 0;
67
69 double dist = 0;
70
72 double diffDistX = 0;
73
75 double diffDistY = 0;
76
78 double diffDistZ = 0;
79
81 double eDep = 0;
82
84 double nPE = 0;
85
87 int channelId = 0;
88
90 bool inRPC = 0;
91
93 bool isFlipped = 0;
94
98 double time() const
99 {
100 return recTime - flyTime;
101 }
102
103 };
104
109
110 /* Not enough data. */
111 c_NotEnoughData = 0,
112
113 /* Failed fit. */
114 c_FailedFit = 1,
115
116 /* Successful calibration. */
117 c_SuccessfulCalibration = 2,
118
119 };
120
125
130
134 void setDebug()
135 {
136 m_debug = true;
137 }
138
144 void setMC(bool mc)
145 {
146 m_mc = mc;
147 }
148
152 void useEvtT0()
153 {
154 m_useEventT0 = true;
155 }
156
160 void setMinimalDigitNumber(int minimalDigitNumber)
161 {
162 m_MinimalDigitNumber = minimalDigitNumber;
163 }
164
170 void setLowerLimit(int counts)
171 {
172 m_lower_limit_counts = counts;
173 }
174
179 void setSaveAllPlots(bool on)
180 {
181 m_saveAllPlots = on;
182 }
183
190 {
192 }
193
197 void saveHist();
198
203 double esti_timeShift(const KLMChannelIndex& klmChannel);
204
209 std::pair<int, double> tS_upperStrip(const KLMChannelIndex& klmChannel);
210
215 std::pair<int, double> tS_lowerStrip(const KLMChannelIndex& klmChannel);
216
221 double esti_timeRes(const KLMChannelIndex& klmChannel);
222
227 std::pair<int, double> tR_upperStrip(const KLMChannelIndex& klmChannel);
228
233 std::pair<int, double> tR_lowerStrip(const KLMChannelIndex& klmChannel);
234
235
236 protected:
237
241 virtual EResult calibrate() override;
242
243 private:
244
248 void setupDatabase();
249
257
262 void readCalibrationDataCounts(std::map<KLMChannelNumber, unsigned int>& eventCounts);
263
270 const std::vector<std::pair<KLMChannelNumber, unsigned int>>& channelsBKLM,
271 const std::vector<std::pair<KLMChannelNumber, unsigned int>>& channelsEKLM);
272
277 void readCalibrationDataBatch(std::function<bool(const KLMChannelIndex&)> channelFilter);
278
282 void createHistograms();
283
288 void writeThenDelete_(TH1* h, bool write);
289
291 void writeThenDelete_(TH2* h, bool write);
292
318 TProfile* profileRpcPhi, TProfile* profileRpcZ,
319 TProfile* profileBKLMScintillatorPhi, TProfile* profileBKLMScintillatorZ,
320 TProfile* profileEKLMScintillatorPlane1,
321 TProfile* profileEKLMScintillatorPlane2, bool fill2dHistograms);
322
330 const std::vector< std::pair<KLMChannelNumber, unsigned int> >& channels,
331 double& delay, double& delayError);
332
337 std::map<KLMChannelNumber, std::vector<struct Event> > m_evts;
338
343 std::map<KLMChannelNumber, int> m_cFlag;
344
346 std::map<KLMChannelNumber, double> m_timeShift;
347
349 std::map<KLMChannelNumber, double> m_timeRes;
350
352 std::map<KLMChannelNumber, double> m_time_channel;
353
355 std::map<KLMChannelNumber, double> m_etime_channel;
356
358 std::map<KLMChannelNumber, double> m_ctime_channel;
359
361 std::map<KLMChannelNumber, double> mc_etime_channel;
362
365
368
371
374
377
380
383
386
389
392
395
398
401
404
408
412
416
420
423
426
430
434
438
442
444 int m_MinimalDigitNumber = 100000000;
445
448
451
454
457
460
462 ROOT::Math::MinimizerOptions m_minimizerOptions;
463
466
472
475
477 bool m_debug = false;
478
480 bool m_mc = false;
481
483 bool m_useEventT0 = true;
484
486 TH1I* h_calibrated = nullptr;
487
489 TH1I* hc_calibrated = nullptr;
490
492 TH1F* h_diff = nullptr;
493
494 /* Monitor graphs of peak value of time distribution for each channel. */
495
497 TGraphErrors* gre_time_channel_rpc = nullptr;
498
500 TGraphErrors* gre_time_channel_scint = nullptr;
501
503 TGraphErrors* gre_time_channel_scint_end = nullptr;
504
505 /* Monitor graphs of peak value of calibrated time distribution for each channel. */
506
508 TGraphErrors* gre_ctime_channel_rpc = nullptr;
509
511 TGraphErrors* gre_ctime_channel_scint = nullptr;
512
514 TGraphErrors* gre_ctime_channel_scint_end = nullptr;
515
516 /* Monitor graphs of calibration constant value of each channel. */
517
519 TGraph* gr_timeShift_channel_rpc = nullptr;
520
522 TGraph* gr_timeShift_channel_scint = nullptr;
523
526
527 /* Monitor graphs of calibrated time resolution value of each channel. */
528
530 TGraph* gr_timeRes_channel_rpc = nullptr;
531
533 TGraph* gr_timeRes_channel_scint = nullptr;
534
537
538 /* Profiles used for effective light speed estimation. */
539
541 TProfile* m_ProfileRpcPhi = nullptr;
542
544 TProfile* m_ProfileRpcZ = nullptr;
545
547 TProfile* m_ProfileBKLMScintillatorPhi = nullptr;
548
550 TProfile* m_ProfileBKLMScintillatorZ = nullptr;
551
554
557
558 /* Profiles of time versus distance (after fit). */
559
561 TProfile* m_Profile2RpcPhi = nullptr;
562
564 TProfile* m_Profile2RpcZ = nullptr;
565
567 TProfile* m_Profile2BKLMScintillatorPhi = nullptr;
568
570 TProfile* m_Profile2BKLMScintillatorZ = nullptr;
571
574
577
578 /*
579 * Histograms of global time distribution used for effective light speed
580 * estimation.
581 */
582
584 TH1F* h_time_rpc_tc = nullptr;
585
587 TH1F* h_time_scint_tc = nullptr;
588
590 TH1F* h_time_scint_tc_end = nullptr;
591
592 /* Histograms of global time distribution before calibration. */
593
595 TH1F* h_time_rpc = nullptr;
596
598 TH1F* h_time_scint = nullptr;
599
601 TH1F* h_time_scint_end = nullptr;
602
603 /* Histograms of global time distribution after calibration. */
604
606 TH1F* hc_time_rpc = nullptr;
607
609 TH1F* hc_time_scint = nullptr;
610
612 TH1F* hc_time_scint_end = nullptr;
613
614 /*
615 * Histograms of time distribution for forward (backward)
616 * before calibration.
617 */
618
620 TH1F* h_timeF_rpc[2] = {nullptr};
621
623 TH1F* h_timeF_scint[2] = {nullptr};
624
626 TH1F* h_timeF_scint_end[2] = {nullptr};
627
628 /*
629 * Histograms of time distribution for forward (backward)
630 * after calibration.
631 */
632
634 TH1F* hc_timeF_rpc[2] = {nullptr};
635
637 TH1F* hc_timeF_scint[2] = {nullptr};
638
640 TH1F* hc_timeF_scint_end[2] = {nullptr};
641
642 /*
643 * Histograms of time dependent on sector for forward (backward)
644 * before calibration.
645 */
646
648 TH2F* h2_timeF_rpc[2] = {nullptr};
649
651 TH2F* h2_timeF_scint[2] = {nullptr};
652
654 TH2F* h2_timeF_scint_end[2] = {nullptr};
655
656 /*
657 * Histograms of time dependent on sector for forward (backward)
658 * after calibration.
659 */
660
662 TH2F* h2c_timeF_rpc[2] = {nullptr};
663
665 TH2F* h2c_timeF_scint[2] = {nullptr};
666
668 TH2F* h2c_timeF_scint_end[2] = {nullptr};
669
670 /* Histograms of time distribution for sectors before calibration. */
671
673 TH1F* h_timeFS_rpc[2][8] = {{nullptr}};
674
676 TH1F* h_timeFS_scint[2][8] = {{nullptr}};
677
679 TH1F* h_timeFS_scint_end[2][4] = {{nullptr}};
680
681 /* Histograms of time distribution for sectors after calibration. */
682
684 TH1F* hc_timeFS_rpc[2][8] = {{nullptr}};
685
687 TH1F* hc_timeFS_scint[2][8] = {{nullptr}};
688
690 TH1F* hc_timeFS_scint_end[2][4] = {{nullptr}};
691
692 /*
693 * Histograms of time distribution dependent on layer of sectors
694 * before calibration.
695 */
696
698 TH2F* h2_timeFS[2][8] = {{nullptr}};
699
701 TH2F* h2_timeFS_end[2][4] = {{nullptr}};
702
703 /*
704 * Histograms of time distribution dependent on layer of sectors
705 * after calibration.
706 */
707
709 TH2F* h2c_timeFS[2][8] = {{nullptr}};
710
712 TH2F* h2c_timeFS_end[2][4] = {{nullptr}};
713
714 /* Histograms of time distribution of one layer before calibration. */
715
717 TH1F* h_timeFSL[2][8][15] = {{{nullptr}}};
718
720 TH1F* h_timeFSL_end[2][4][14] = {{{nullptr}}};
721
722 /* Histograms of time distribution of one layer after calibration. */
723
725 TH1F* hc_timeFSL[2][8][15] = {{{nullptr}}};
726
728 TH1F* hc_timeFSL_end[2][4][14] = {{{nullptr}}};
729
730 /* Histograms of time distribution of one plane before calibration. */
731
733 TH1F* h_timeFSLP[2][8][15][2] = {{{{nullptr}}}};
734
736 TH1F* h_timeFSLP_end[2][4][14][2] = {{{{nullptr}}}};
737
738 /* Histograms of time distribution of one plane after calibration. */
739
741 TH1F* hc_timeFSLP[2][8][15][2] = {{{{nullptr}}}};
742
744 TH1F* hc_timeFSLP_end[2][4][14][2] = {{{{nullptr}}}};
745
746 /*
747 * Histograms of time distribution dependent on channels
748 * before calibration.
749 */
750
752 TH2F* h2_timeFSLP[2][8][15][2] = {{{{nullptr}}}};
753
755 TH2F* h2_timeFSLP_end[2][4][14][2] = {{{{nullptr}}}};
756
757 /*
758 * Histograms of time distribution dependent on channels
759 * after calibration.
760 */
761
763 TH2F* h2c_timeFSLP[2][8][15][2] = {{{{nullptr}}}};
764
766 TH2F* h2c_timeFSLP_end[2][4][14][2] = {{{{nullptr}}}};
767
768 /* Formulas used for fitting. */
769
771 TF1* fcn_pol1 = nullptr;
772
774 TF1* fcn_const = nullptr;
775
777 TF1* fcn_gaus = nullptr;
778
780 TF1* fcn_land = nullptr;
781
783 TFile* m_outFile = nullptr;
784
786 bool m_saveAllPlots = false;
787
789 bool m_saveChannelHists = false;
790
791 };
792
794}
EResult
The result of calibration.
CalibrationAlgorithm(const std::string &collectorModuleName)
Constructor - sets the prefix for collected objects (won't be accesses until execute(....
EKLM geometry data.
KLM channel index.
TProfile * m_Profile2EKLMScintillatorPlane2
For EKLM scintillator plane2.
double mc_etime_channelAvg_rpc
Calibrated central value error of the global time distribution (BKLM RPC part).
TH2F * h2c_timeF_scint_end[2]
EKLM part.
KLMTimeResolution * m_timeResolution
DBObject of time resolution.
TH1F * h_time_scint_tc_end
EKLM part.
void createHistograms()
Create histograms.
TGraphErrors * gre_time_channel_scint
BKLM Scintillator.
TH1F * h_timeFSL[2][8][15]
BKLM part.
TH1F * hc_timeFSL_end[2][4][14]
EKLM part.
TH1F * h_timeFSLP_end[2][4][14][2]
EKLM part.
TGraph * gr_timeRes_channel_rpc
BKLM RPC.
TH1F * hc_timeFSLP_end[2][4][14][2]
EKLM part.
TH1F * h_timeFSLP[2][8][15][2]
BKLM part.
TH1F * hc_timeF_scint_end[2]
EKLM part.
std::map< KLMChannelNumber, double > m_timeShift
Shift values of each channel.
TH1F * h_time_scint
BKLM scintillator part.
double m_time_channelAvg_scint
Central value of the global time distribution (BKLM scintillator part).
TH1F * hc_timeFS_scint_end[2][4]
EKLM part.
double esti_timeRes(const KLMChannelIndex &klmChannel)
Estimate value of calibration constant for calibrated channels.
double m_UpperTimeBoundaryCalibratedRPC
Upper time boundary for RPC (calibrated data).
double m_ctime_channelAvg_rpc
Calibrated central value of the global time distribution (BKLM RPC part).
KLMTimeConstants * m_timeConstants
DBObject of time cost on some parts of the detector.
void setupDatabase()
Setup the database.
void setMinimalDigitNumber(int minimalDigitNumber)
Set minimal digit number (total).
TH1F * hc_timeFS_scint[2][8]
BKLM scintillator part.
std::map< KLMChannelNumber, double > m_time_channel
Time distribution central value of each channel.
double m_ctime_channelAvg_scint_end
Calibrated central value of the global time distribution (EKLM scintillator part).
CalibrationAlgorithm::EResult readCalibrationData()
Read calibration data.
TGraph * gr_timeShift_channel_scint_end
EKLM.
TGraph * gr_timeRes_channel_scint
BKLM scintillator.
TH1F * hc_timeF_scint[2]
BKLM scintillator part.
TH1F * h_timeFS_scint[2][8]
BKLM scintillator part.
bool m_saveChannelHists
Write per-channel temporary histograms (tc/raw/hc) in minimal mode.
double m_UpperTimeBoundaryScintillatorsBKLM
Upper time boundary for BKLM scintillators.
const KLMElementNumbers * m_ElementNumbers
Element numbers.
std::pair< int, double > tR_upperStrip(const KLMChannelIndex &klmChannel)
Tracing available channels with increasing strip number.
TH1F * hc_timeF_rpc[2]
BKLM RPC part.
TH2F * h2c_timeFS_end[2][4]
EKLM part.
const EKLM::GeometryData * m_EKLMGeometry
EKLM geometry data.
TGraphErrors * gre_ctime_channel_scint_end
EKLM.
TProfile * m_Profile2BKLMScintillatorPhi
For BKLM scintillator phi plane.
TH1F * hc_time_scint_end
EKLM part.
TGraphErrors * gre_time_channel_scint_end
EKLM.
TH2F * h2_timeFSLP[2][8][15][2]
BKLM part.
double m_UpperTimeBoundaryCalibratedScintillatorsEKLM
Upper time boundary for BKLM scintillators (calibrated data).
TGraph * gr_timeShift_channel_scint
BKLM scintillator.
double m_time_channelAvg_scint_end
Central value of the global time distribution (EKLM scintillator part).
TProfile * m_Profile2EKLMScintillatorPlane1
For EKLM scintillator plane1.
TH1F * hc_timeFS_rpc[2][8]
BKLM RPC part.
double m_UpperTimeBoundaryScintillatorsEKLM
Upper time boundary for BKLM scintillators.
void fillTimeDistanceProfiles(TProfile *profileRpcPhi, TProfile *profileRpcZ, TProfile *profileBKLMScintillatorPhi, TProfile *profileBKLMScintillatorZ, TProfile *profileEKLMScintillatorPlane1, TProfile *profileEKLMScintillatorPlane2, bool fill2dHistograms)
Fill profiles of time versus distance.
TFile * m_outFile
Output file.
double esti_timeShift(const KLMChannelIndex &klmChannel)
Estimate value of calibration constant for uncalibrated channels.
std::pair< int, double > tS_upperStrip(const KLMChannelIndex &klmChannel)
Tracing available channels with increasing strip number.
double m_LowerTimeBoundaryCalibratedScintillatorsEKLM
Lower time boundary for EKLM scintillators (calibrated data).
TH1F * hc_timeFSLP[2][8][15][2]
BKLM part.
TGraphErrors * gre_ctime_channel_rpc
BKLM RPC.
void saveHist()
Save histograms to file.
const bklm::GeometryPar * m_BKLMGeometry
BKLM geometry data.
bool m_saveAllPlots
Default minimal unless you set true in your header script.
TH2F * h2c_timeFSLP[2][8][15][2]
BKLM part.
double m_ctime_channelAvg_scint
Calibrated central value of the global time distribution (BKLM scintillator part).
TF1 * fcn_const
Const function.
double m_UpperTimeBoundaryCalibratedScintillatorsBKLM
Upper time boundary for BKLM scintillators (calibrated data).
TProfile * m_Profile2RpcZ
For BKLM RPC z plane.
void setDebug()
Turn on debug mode (prints histograms and output running log).
TH2F * h2_timeFSLP_end[2][4][14][2]
EKLM part.
TH1I * hc_calibrated
Calibration statistics for each channel.
void timeDistance2dFit(const std::vector< std::pair< KLMChannelNumber, unsigned int > > &channels, double &delay, double &delayError)
Two-dimensional fit for individual channels.
TGraph * gr_timeRes_channel_scint_end
EKLM.
TH1I * h_calibrated
Calibration statistics for each channel.
TProfile * m_ProfileBKLMScintillatorZ
For BKLM scintillator z plane.
double m_LowerTimeBoundaryScintillatorsBKLM
Lower time boundary for BKLM scintillators.
TH1F * h_time_rpc_tc
BKLM RPC part.
TH1F * h_time_scint_end
EKLM part.
TH2F * h2c_timeF_scint[2]
BKLM scintillator part.
TF1 * fcn_pol1
Pol1 function.
double m_etime_channelAvg_scint_end
Central value error of the global time distribution (EKLM scintillator part).
TH1F * hc_time_rpc
BKLM RPC part.
double mc_etime_channelAvg_scint
Calibrated central value error of the global time distribution (BKLM scintillator part).
TH2F * h2c_timeFSLP_end[2][4][14][2]
EKLM part.
double mc_etime_channelAvg_scint_end
Calibrated central value error of the global time distribution (EKLM scintillator part).
void writeThenDelete_(TH1 *h, bool write)
Optionally write a histogram, then delete it to free memory.
TH2F * h2_timeF_scint_end[2]
EKLM part.
void useEvtT0()
Use event T0 as the initial time point or not.
TF1 * fcn_gaus
Gaussian function.
double m_LowerTimeBoundaryCalibratedScintillatorsBKLM
Lower time boundary for BKLM scintillators (calibrated data).
TH1F * h_timeF_rpc[2]
BKLM RPC part.
TProfile * m_ProfileBKLMScintillatorPhi
For BKLM scintillator phi plane.
void setSaveAllPlots(bool on)
Save every preallocated debug histogram family (sectors/layers/planes/2D maps).
TH1F * hc_time_scint
BKLM scintillator part.
TH2F * h2_timeFS[2][8]
BKLM part.
double m_etime_channelAvg_scint
Central value error of the global time distribution (BKLM scintillator part).
TH1F * h_timeF_scint_end[2]
EKLM part.
TProfile * m_ProfileRpcPhi
For BKLM RPC phi plane.
TGraphErrors * gre_ctime_channel_scint
BKLM Scintillator.
TProfile * m_ProfileEKLMScintillatorPlane2
For EKLM scintillator plane2.
TH1F * h_time_rpc
BKLM RPC part.
TH1F * h_timeFSL_end[2][4][14]
EKLM part.
TProfile * m_Profile2RpcPhi
For BKLM RPC phi plane.
TH1F * h_timeF_scint[2]
BKLM scintillator part.
TH1F * hc_timeFSL[2][8][15]
BKLM part.
TProfile * m_Profile2BKLMScintillatorZ
For BKLM scintillator z plane.
TH1F * h_timeFS_rpc[2][8]
BKLM RPC part.
KLMChannelIndex m_klmChannels
KLM ChannelIndex object.
TGraph * gr_timeShift_channel_rpc
BKLM RPC.
std::map< KLMChannelNumber, double > m_timeRes
Resolution values of each channel.
TH1F * h_diff
Distance between global and local position.
TH2F * h2_timeF_scint[2]
BKLM scintillator part.
TH1F * h_time_scint_tc
BKLM scintillator part.
double m_LowerTimeBoundaryRPC
Lower time boundary for RPC.
virtual EResult calibrate() override
Run algorithm on data.
std::map< KLMChannelNumber, double > m_ctime_channel
Calibrated time distribution central value of each channel.
double m_LowerTimeBoundaryCalibratedRPC
Lower time boundary for RPC (calibrated data).
bool m_useEventT0
Whether to use event T0 from CDC.
int m_MinimalDigitNumber
Minimal digit number (total).
void setLowerLimit(int counts)
Set the lower number of hits collected on one single strip.
ChannelCalibrationStatus
Channel calibration status.
void setSaveChannelHists(bool on)
When running in minimal mode, also write the per-channel temporary histograms (the vital tc,...
TProfile * m_ProfileEKLMScintillatorPlane1
For EKLM scintillator plane1.
double m_UpperTimeBoundaryRPC
Upper time boundary for RPC.
TH2F * h2c_timeF_rpc[2]
BKLM RPC part.
TF1 * fcn_land
Landau function.
KLMTimeCableDelay * m_timeCableDelay
DBObject of the calibration constant of each channel due to cable decay.
std::pair< int, double > tS_lowerStrip(const KLMChannelIndex &klmChannel)
Tracing available channels with decreasing strip number.
TProfile * m_ProfileRpcZ
For BKLM RPC z plane.
std::map< KLMChannelNumber, double > mc_etime_channel
Calibrated time distribution central value Error of each channel.
std::pair< int, double > tR_lowerStrip(const KLMChannelIndex &klmChannel)
Tracing available channels with decreasing strip number.
void setMC(bool mc)
Set flag indicating whether the input is MC sample.
void readCalibrationDataBatch(std::function< bool(const KLMChannelIndex &)> channelFilter)
Load calibration data for a specific batch of channels.
ROOT::Math::MinimizerOptions m_minimizerOptions
Minimization options.
std::map< KLMChannelNumber, int > m_cFlag
Calibration flag if the channel has enough hits collected and fitted OK.
TGraphErrors * gre_time_channel_rpc
BKLM RPC.
TH2F * h2_timeFS_end[2][4]
EKLM part.
TH2F * h2_timeF_rpc[2]
BKLM RPC part.
std::map< KLMChannelNumber, std::vector< struct Event > > m_evts
Container of hit information.
TH1F * h_timeFS_scint_end[2][4]
EKLM part.
double m_time_channelAvg_rpc
Central value of the global time distribution (BKLM RPC part).
TH2F * h2c_timeFS[2][8]
BKLM part.
double m_etime_channelAvg_rpc
Central value error of the global time distribution (BKLM RPC part).
double m_LowerTimeBoundaryScintillatorsEKLM
Lower time boundary for EKLM scintillators.
void readCalibrationDataCounts(std::map< KLMChannelNumber, unsigned int > &eventCounts)
Count events per channel (lightweight scan without loading full data).
std::map< KLMChannelNumber, double > m_etime_channel
Time distribution central value Error of each channel.
int m_lower_limit_counts
Lower limit of hits collected for on single channel.
void readCalibrationDataFor2DFit(const std::vector< std::pair< KLMChannelNumber, unsigned int > > &channelsBKLM, const std::vector< std::pair< KLMChannelNumber, unsigned int > > &channelsEKLM)
Load calibration data only for channels needed for 2D fit.
Class to store BKLM delay time coused by cable in the database.
Class to store KLM constants related to time.
Class to store KLM time resolution in the database.
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition GeometryPar.h:37
Abstract base class for different kinds of events.
double flyTime
Particle flying time.
double diffDistX
Global position difference between klmHit2d and ExtHit (X).
double t0
EventT0 for the digit.
int channelId
Unique channel id Barral and endcap merged.
bool inRPC
BKLM RPC flag, used for testing and not necessary.
double eDep
Collect energy eV.
double diffDistY
Global position difference between klmHit2d and ExtHit (Y).
double diffDistZ
Global position difference between klmHit2d and ExtHit (Z).
double time() const
Get propagation time + cableDelay time.
double dist
Propagation distance from hit to FEE.
double nPE
Number of photon electron.
double recTime
Recosntruction time respected to the trigger time.
bool isFlipped
If phi and z plane flipped, used for testing and not necessary.