Belle II Software  release-06-02-00
KLMTimeAlgorithm.cc
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 /* Own header. */
10 #include <klm/calibration/KLMTimeAlgorithm.h>
11 
12 /* KLM headers. */
13 #include <klm/dataobjects/bklm/BKLMElementNumbers.h>
14 
15 /* Belle 2 headers. */
16 #include <framework/database/Database.h>
17 #include <framework/database/DBObjPtr.h>
18 #include <framework/database/DBStore.h>
19 #include <framework/dataobjects/EventMetaData.h>
20 #include <framework/gearbox/Const.h>
21 #include <framework/logging/Logger.h>
22 
23 /* ROOT headers. */
24 #include <Math/MinimizerOptions.h>
25 #include <TFile.h>
26 #include <TFitResult.h>
27 #include <TMinuit.h>
28 #include <TROOT.h>
29 #include <TString.h>
30 #include <TTree.h>
31 #include <TVector3.h>
32 
33 using namespace Belle2;
34 using namespace ROOT::Math;
35 
37 const int c_NBinsTime = 100;
38 
40 const int c_NBinsDistance = 100;
41 
43 static double s_BinnedData[c_NBinsTime][c_NBinsDistance];
44 
46 static double s_LowerTimeBoundary = 0;
47 
49 static double s_UpperTimeBoundary = 0;
50 
52 static double s_StripLength = 0;
53 
54 static bool compareEventNumber(const std::pair<KLMChannelNumber, unsigned int>& pair1,
55  const std::pair<KLMChannelNumber, unsigned int>& pair2)
56 {
57  return pair1.second < pair2.second;
58 }
59 
60 static double timeDensity(double x[2], double* par)
61 {
62  double polynomial, t0, gauss;
63  polynomial = par[0];
64  t0 = par[2] + par[4] * x[1];
65  gauss = par[1] / (sqrt(2.0 * M_PI) * par[3]) *
66  exp(-0.5 * pow((x[0] - t0) / par[3], 2));
67  return fabs(polynomial + gauss);
68 }
69 
70 static void fcn(int& npar, double* grad, double& fval, double* par, int iflag)
71 {
72  (void)npar;
73  (void)grad;
74  (void)iflag;
75  double x[2];
76  fval = 0;
77  for (int i = 0; i < c_NBinsTime; ++i) {
78  x[0] = s_LowerTimeBoundary +
79  (s_UpperTimeBoundary - s_LowerTimeBoundary) *
80  (double(i) + 0.5) / c_NBinsTime;
81  for (int j = 0; j < c_NBinsDistance; ++j) {
82  x[1] = s_StripLength * (double(j) + 0.5) / c_NBinsDistance;
83  double f = timeDensity(x, par);
84  if (s_BinnedData[i][j] == 0)
85  fval = fval + 2.0 * f;
86  else
87  fval = fval + 2.0 * (f - s_BinnedData[i][j] *
88  (1.0 - log(s_BinnedData[i][j] / f)));
89  }
90  }
91 }
92 
94  CalibrationAlgorithm("KLMTimeCollector")
95 {
97  m_minimizerOptions = ROOT::Math::MinimizerOptions();
98 }
99 
101 {
102 }
103 
105 {
106  const std::vector<Calibration::ExpRun>& runs = getRunList();
107  int firstExperiment = runs[0].first;
108  int lastExperiment = runs[runs.size() - 1].first;
109  if (firstExperiment != lastExperiment) {
110  B2FATAL("Runs from different experiments are used "
111  "for KLM time calibration (single algorithm run).");
112  }
113  /* DataStore. */
115  StoreObjPtr<EventMetaData> eventMetaData;
116  eventMetaData.registerInDataStore();
118  /* Database. */
119  if (eventMetaData.isValid()) {
120  if (eventMetaData->getExperiment() != firstExperiment) {
121  B2FATAL("Runs from different experiments are used "
122  "for KLM time calibration (consecutive algorithm runs).");
123  }
124  eventMetaData->setExperiment(firstExperiment);
125  eventMetaData->setRun(runs[0].second);
126  } else {
127  eventMetaData.construct(1, runs[0].second, firstExperiment);
128  }
129  DBStore& dbStore = DBStore::Instance();
130  dbStore.update();
131  dbStore.updateEvent();
132  /*
133  * For calibration algorithms, the database is not initialized on class
134  * creation. Do not move the database object to class members.
135  */
136  DBObjPtr<BKLMGeometryPar> bklmGeometry;
139 }
140 
142 {
143  B2INFO("Read tree entries and seprate events by module id.");
144  Event event;
145  std::shared_ptr<TTree> timeCalibrationData;
146  timeCalibrationData = getObjectPtr<TTree>("time_calibration_data");
147  timeCalibrationData->SetBranchAddress("t0", &event.t0);
148  timeCalibrationData->SetBranchAddress("flyTime", &event.flyTime);
149  timeCalibrationData->SetBranchAddress("recTime", &event.recTime);
150  timeCalibrationData->SetBranchAddress("dist", &event.dist);
151  timeCalibrationData->SetBranchAddress("diffDistX", &event.diffDistX);
152  timeCalibrationData->SetBranchAddress("diffDistY", &event.diffDistY);
153  timeCalibrationData->SetBranchAddress("diffDistZ", &event.diffDistZ);
154  timeCalibrationData->SetBranchAddress("eDep", &event.eDep);
155  timeCalibrationData->SetBranchAddress("nPE", &event.nPE);
156  timeCalibrationData->SetBranchAddress("channelId", &event.channelId);
157  timeCalibrationData->SetBranchAddress("inRPC", &event.inRPC);
158  timeCalibrationData->SetBranchAddress("isFlipped", &event.isFlipped);
159 
160  B2INFO(LogVar("Total number of digits:", timeCalibrationData->GetEntries()));
161  m_evts.clear();
162 
163  int n = timeCalibrationData->GetEntries();
164  if (n < m_MinimalDigitNumber)
166  for (int i = 0; i < n; ++i) {
167  timeCalibrationData->GetEntry(i);
168  m_evts[event.channelId].push_back(event);
169  }
170  B2INFO("Events packing finish.");
172 }
173 
175 {
176  if (m_mc) {
177  m_LowerTimeBoundaryRPC = -10.0;
178  m_UpperTimeBoundaryRPC = 10.0;
183  } else {
184  m_LowerTimeBoundaryRPC = -800.0;
185  m_UpperTimeBoundaryRPC = -600.0;
190  }
191  int nBin = 200;
192  int nBin_scint = 400;
193 
194  TString iFstring[2] = {"Backward", "Forward"};
195  TString iPstring[2] = {"ZReadout", "PhiReadout"};
196  TString hn, ht;
197 
198  h_diff = new TH1F("h_diff", "Position difference between bklmHit2d and extHit;position difference", 100, 0, 10);
199  h_calibrated = new TH1I("h_calibrated_summary", "h_calibrated_summary;calibrated or not", 3, 0, 3);
200 
201  gre_time_channel_scint = new TGraphErrors();
202  gre_time_channel_rpc = new TGraphErrors();
203  gre_time_channel_scint_end = new TGraphErrors();
204 
205  gr_timeShift_channel_scint = new TGraph();
206  gr_timeShift_channel_rpc = new TGraph();
207  gr_timeShift_channel_scint_end = new TGraph();
208 
209  double maximalPhiStripLengthBKLM =
211  double maximalZStripLengthBKLM =
213  double maximalStripLengthEKLM =
215 
216  m_ProfileRpcPhi = new TProfile("hprf_rpc_phi_effC",
217  "Time over propagation length for RPCs (Phi_Readout); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]", 400, 0.0,
218  400.0);
219  m_ProfileRpcZ = new TProfile("hprf_rpc_z_effC",
220  "Time over propagation length for RPCs (Z_Readout); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]", 400, 0.0,
221  400.0);
222  m_ProfileBKLMScintillatorPhi = new TProfile("hprf_scint_phi_effC",
223  "Time over propagation length for scintillators (Phi_Readout); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]",
224  200, 0.0, maximalPhiStripLengthBKLM);
225  m_ProfileBKLMScintillatorZ = new TProfile("hprf_scint_z_effC",
226  "Time over propagation length for scintillators (Z_Readout); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]",
227  200, 0.0, maximalZStripLengthBKLM);
228  m_ProfileEKLMScintillatorPlane1 = new TProfile("hprf_scint_plane1_effC_end",
229  "Time over propagation length for scintillators (plane1, Endcap); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]",
230  200, 0.0, maximalStripLengthEKLM);
231  m_ProfileEKLMScintillatorPlane2 = new TProfile("hprf_scint_plane2_effC_end",
232  "Time over propagation length for scintillators (plane2, Endcap); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]",
233  200, 0.0, maximalStripLengthEKLM);
234 
235  m_Profile2RpcPhi = new TProfile("hprf2_rpc_phi_effC",
236  "Time over propagation length for RPCs (Phi_Readout); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]", 400, 0.0,
237  400.0);
238  m_Profile2RpcZ = new TProfile("hprf2_rpc_z_effC",
239  "Time over propagation length for RPCs (Z_Readout); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]", 400, 0.0,
240  400.0);
241  m_Profile2BKLMScintillatorPhi = new TProfile("hprf2_scint_phi_effC",
242  "Time over propagation length for scintillators (Phi_Readout); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]",
243  200, 0.0, maximalPhiStripLengthBKLM);
244  m_Profile2BKLMScintillatorZ = new TProfile("hprf2_scint_z_effC",
245  "Time over propagation length for scintillators (Z_Readout); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]",
246  200, 0.0, maximalZStripLengthBKLM);
247  m_Profile2EKLMScintillatorPlane1 = new TProfile("hprf2_scint_plane1_effC_end",
248  "Time over propagation length for scintillators (plane1, Endcap); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]",
249  200, 0.0, maximalStripLengthEKLM);
250  m_Profile2EKLMScintillatorPlane2 = new TProfile("hprf2_scint_plane2_effC_end",
251  "Time over propagation length for scintillators (plane2, Endcap); propagation distance[cm]; T_rec-T_0-T_fly-'T_calibration'[ns]",
252  200, 0.0, maximalStripLengthEKLM);
253 
254  h_time_rpc_tc = new TH1F("h_time_rpc_tc", "time distribtution for RPC", nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
255  h_time_scint_tc = new TH1F("h_time_scint_tc", "time distribtution for Scintillator", nBin_scint,
257  h_time_scint_tc_end = new TH1F("h_time_scint_tc_end", "time distribtution for Scintillator (Endcap)", nBin_scint,
260 
262  h_time_rpc = new TH1F("h_time_rpc", "time distribtution for RPC; T_rec-T_0-T_fly-T_propagation[ns]", nBin, m_LowerTimeBoundaryRPC,
264  h_time_scint = new TH1F("h_time_scint", "time distribtution for Scintillator; T_rec-T_0-T_fly-T_propagation[ns]", nBin_scint,
266  h_time_scint_end = new TH1F("h_time_scint_end", "time distribtution for Scintillator (Endcap); T_rec-T_0-T_fly-T_propagation[ns]",
268 
269  hc_time_rpc = new TH1F("hc_time_rpc", "Calibrated time distribtution for RPC; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
271  hc_time_scint = new TH1F("hc_time_scint",
272  "Calibrated time distribtution for Scintillator; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]", nBin_scint,
275  hc_time_scint_end = new TH1F("hc_time_scint_end",
276  "Calibrated time distribtution for Scintillator (Endcap); T_rec-T_0-T_fly-T_propagation-T_calibration[ns]", nBin_scint,
278 
279  for (int iF = 0; iF < 2; ++iF) {
280  hn = Form("h_timeF%d_rpc", iF);
281  ht = Form("Time distribtution for RPC of %s; T_rec-T_0-T_fly-T_propagation[ns]", iFstring[iF].Data());
282  h_timeF_rpc[iF] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
283  hn = Form("h_timeF%d_scint", iF);
284  ht = Form("Time distribtution for Scintillator of %s; T_rec-T_0-T_fly-T_propagation[ns]", iFstring[iF].Data());
285  h_timeF_scint[iF] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsBKLM,
287  hn = Form("h_timeF%d_scint_end", iF);
288  ht = Form("Time distribtution for Scintillator of %s (Endcap); T_rec-T_0-T_fly-T_propagation[ns]", iFstring[iF].Data());
289  h_timeF_scint_end[iF] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
291 
292  hn = Form("h2_timeF%d_rpc", iF);
293  ht = Form("Time distribtution for RPC of %s; Sector Index; T_rec-T_0-T_fly-T_propagation[ns]", iFstring[iF].Data());
294  h2_timeF_rpc[iF] = new TH2F(hn.Data(), ht.Data(), 8, 0, 8, nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
295  hn = Form("h2_timeF%d_scint", iF);
296  ht = Form("Time distribtution for Scintillator of %s; Sector Index; T_rec-T_0-T_fly-T_propagation[ns]", iFstring[iF].Data());
297  h2_timeF_scint[iF] = new TH2F(hn.Data(), ht.Data(), 8, 0, 8, nBin_scint, m_LowerTimeBoundaryScintilltorsBKLM,
299  hn = Form("h2_timeF%d_scint_end", iF);
300  ht = Form("Time distribtution for Scintillator of %s (Endcap); Sector Index; T_rec-T_0-T_fly-T_propagation[ns]",
301  iFstring[iF].Data());
302  h2_timeF_scint_end[iF] = new TH2F(hn.Data(), ht.Data(), 4, 0, 4, nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
304 
305  hn = Form("hc_timeF%d_rpc", iF);
306  ht = Form("Calibrated time distribtution for RPC of %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]", iFstring[iF].Data());
307  hc_timeF_rpc[iF] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryCalibratedRPC, m_UpperTimeBoundaryCalibratedRPC);
308  hn = Form("hc_timeF%d_scint", iF);
309  ht = Form("Calibrated time distribtution for Scintillator of %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
310  iFstring[iF].Data());
311  hc_timeF_scint[iF] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsBKLM,
313  hn = Form("hc_timeF%d_scint_end", iF);
314  ht = Form("Calibrated time distribtution for Scintillator of %s (Endcap); T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
315  iFstring[iF].Data());
316  hc_timeF_scint_end[iF] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsEKLM,
318 
319  hn = Form("h2c_timeF%d_rpc", iF);
320  ht = Form("Calibrated time distribtution for RPC of %s; Sector Index; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
321  iFstring[iF].Data());
322  h2c_timeF_rpc[iF] = new TH2F(hn.Data(), ht.Data(), 8, 0, 8, nBin, m_LowerTimeBoundaryCalibratedRPC,
324  hn = Form("h2c_timeF%d_scint", iF);
325  ht = Form("Calibrated time distribtution for Scintillator of %s; Sector Index; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
326  iFstring[iF].Data());
327  h2c_timeF_scint[iF] = new TH2F(hn.Data(), ht.Data(), 8, 0, 8, nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsBKLM,
329  hn = Form("h2c_timeF%d_scint_end", iF);
330  ht = Form("Calibrated time distribtution for Scintillator of %s (Endcap) ; Sector Index; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
331  iFstring[iF].Data());
332  h2c_timeF_scint_end[iF] = new TH2F(hn.Data(), ht.Data(), 4, 0, 4, nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsEKLM,
334 
335  for (int iS = 0; iS < 8; ++iS) {
336  // Barrel parts
337  hn = Form("h_timeF%d_S%d_scint", iF, iS);
338  ht = Form("Time distribtution for Scintillator of Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iS, iFstring[iF].Data());
339  h_timeFS_scint[iF][iS] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsBKLM,
341  hn = Form("h_timeF%d_S%d_rpc", iF, iS);
342  ht = Form("Time distribtution for RPC of Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iS, iFstring[iF].Data());
343  h_timeFS_rpc[iF][iS] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
344  hn = Form("h2_timeF%d_S%d", iF, iS);
345  ht = Form("Time distribtution of Sector%d, %s; Layer Index; T_rec-T_0-T_fly-T_propagation[ns]", iS, iFstring[iF].Data());
346  h2_timeFS[iF][iS] = new TH2F(hn.Data(), ht.Data(), 15, 0, 15, nBin_scint, m_LowerTimeBoundaryRPC,
348 
349  hn = Form("hc_timeF%d_S%d_scint", iF, iS);
350  ht = Form("Calibrated time distribtution for Scintillator of Sector%d, %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]", iS,
351  iFstring[iF].Data());
352  hc_timeFS_scint[iF][iS] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsBKLM,
354  hn = Form("hc_timeF%d_S%d_rpc", iF, iS);
355  ht = Form("Calibrated time distribtution for RPC of Sector%d, %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]", iS,
356  iFstring[iF].Data());
357  hc_timeFS_rpc[iF][iS] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedRPC,
359  hn = Form("h2c_timeF%d_S%d", iF, iS);
360  ht = Form("Calibrated time distribtution of Sector%d, %s; Layer Index; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]", iS,
361  iFstring[iF].Data());
362  h2c_timeFS[iF][iS] = new TH2F(hn.Data(), ht.Data(), 15, 0, 15, nBin_scint, m_LowerTimeBoundaryCalibratedRPC,
364 
365  // Inner 2 layers --> Scintillators
366  for (int iL = 0; iL < 2; ++iL) {
367  hn = Form("h_timeF%d_S%d_L%d", iF, iS, iL);
368  ht = Form("Time distribtution for Scintillator of Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iL, iS,
369  iFstring[iF].Data());
370  h_timeFSL[iF][iS][iL] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsBKLM,
372  hn = Form("hc_timeF%d_S%d_L%d", iF, iS, iL);
373  ht = Form("Calibrated time distribtution for Scintillator of Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
374  iL, iS, iFstring[iF].Data());
375  hc_timeFSL[iF][iS][iL] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsBKLM,
377 
378  for (int iP = 0; iP < 2; ++iP) {
379  hn = Form("h_timeF%d_S%d_L%d_P%d", iF, iS, iL, iP);
380  ht = Form("Time distribtution for Scintillator of %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]",
381  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
382  h_timeFSLP[iF][iS][iL][iP] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsBKLM,
384  hn = Form("h2_timeF%d_S%d_L%d_P%d", iF, iS, iL, iP);
385  ht = Form("Time distribtution for Scintillator of %s, Layer%d, Sector%d, %s; Channel Index; T_rec-T_0-T_fly-T_propagation[ns]",
386  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
387  h2_timeFSLP[iF][iS][iL][iP] = new TH2F(hn.Data(), ht.Data(), 54, 0, 54, nBin_scint, m_LowerTimeBoundaryScintilltorsBKLM,
389 
390  hn = Form("hc_timeF%d_S%d_L%d_P%d", iF, iS, iL, iP);
391  ht = Form("Calibrated time distribtution for Scintillator of %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
392  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
393  hc_timeFSLP[iF][iS][iL][iP] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsBKLM,
395  hn = Form("h2c_timeF%d_S%d_L%d_P%d", iF, iS, iL, iP);
396  ht = Form("Calibrated time distribtution for Scintillator of %s, Layer%d, Sector%d, %s; Channel Index; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
397  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
398  h2c_timeFSLP[iF][iS][iL][iP] = new TH2F(hn.Data(), ht.Data(), 54, 0, 54, nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsBKLM,
400 
401  int nchannel_max = BKLMElementNumbers::getNStrips(iF, iS + 1, iL + 1, iP);
402  for (int iC = 0; iC < nchannel_max; ++iC) {
403  hn = Form("h_timeF%d_S%d_L%d_P%d_C%d_tc", iF, iS, iL, iP, iC);
404  ht = Form("time distribtution for Scintillator of Channel%d, %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iC,
405  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
406  h_timeFSLPC_tc[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsBKLM,
408 
409  hn = Form("h_timeF%d_S%d_L%d_P%d_C%d", iF, iS, iL, iP, iC);
410  ht = Form("time distribtution for Scintillator of Channel%d, %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iC,
411  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
412  h_timeFSLPC[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsBKLM,
414 
415  hn = Form("hc_timeF%d_S%d_L%d_P%d_C%d", iF, iS, iL, iP, iC);
416  ht = Form("Calibrated time distribtution for Scintillator of Channel%d, %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
417  iC, iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
418  hc_timeFSLPC[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsBKLM,
420  hn = Form("time_length_bklm_F%d_S%d_L%d_P%d_C%d", iF, iS, iL, iP, iC);
421  double stripLength = 200;
422  m_HistTimeLengthBKLM[iF][iS][iL][iP][iC] =
423  new TH2F(hn.Data(),
424  "Time versus propagation length; "
425  "propagation distance[cm]; "
426  "T_rec-T_0-T_fly-'T_calibration'[ns]",
427  200, 0.0, stripLength,
430  }
431  }
432  }
433 
434  for (int iL = 2; iL < 15; ++iL) {
435  hn = Form("h_timeF%d_S%d_L%d", iF, iS, iL);
436  ht = Form("time distribtution for RPC of Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iL, iS, iFstring[iF].Data());
437  h_timeFSL[iF][iS][iL] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
438 
439  hn = Form("hc_timeF%d_S%d_L%d", iF, iS, iL);
440  ht = Form("Calibrated time distribtution for RPC of Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]", iL, iS,
441  iFstring[iF].Data());
442  hc_timeFSL[iF][iS][iL] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryCalibratedRPC, m_UpperTimeBoundaryCalibratedRPC);
443 
444  for (int iP = 0; iP < 2; ++iP) {
445  hn = Form("h_timeF%d_S%d_L%d_P%d", iF, iS, iL, iP);
446  ht = Form("time distribtution for RPC of %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iPstring[iP].Data(), iL, iS,
447  iFstring[iF].Data());
448  h_timeFSLP[iF][iS][iL][iP] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
449 
450  hn = Form("h2_timeF%d_S%d_L%d_P%d", iF, iS, iL, iP);
451  ht = Form("time distribtution for RPC of %s, Layer%d, Sector%d, %s; Channel Index; T_rec-T_0-T_fly-T_propagation[ns]",
452  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
453  h2_timeFSLP[iF][iS][iL][iP] = new TH2F(hn.Data(), ht.Data(), 48, 0, 48, nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
454 
455  hn = Form("hc_timeF%d_S%d_L%d_P%d", iF, iS, iL, iP);
456  ht = Form("Calibrated time distribtution for RPC of %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
457  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
458  hc_timeFSLP[iF][iS][iL][iP] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryCalibratedRPC,
460 
461  hn = Form("h2c_timeF%d_S%d_L%d_P%d", iF, iS, iL, iP);
462  ht = Form("Calibrated time distribtution for RPC of %s, Layer%d, Sector%d, %s; Channel Index; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
463  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
464  h2c_timeFSLP[iF][iS][iL][iP] = new TH2F(hn.Data(), ht.Data(), 48, 0, 48, nBin, m_LowerTimeBoundaryCalibratedRPC,
466 
467  int nchannel_max = BKLMElementNumbers::getNStrips(iF, iS + 1, iL + 1, iP);
468  for (int iC = 0; iC < nchannel_max; ++iC) {
469  hn = Form("h_timeF%d_S%d_L%d_P%d_C%d_tc", iF, iS, iL, iP, iC);
470  ht = Form("Time distribtution for RPC of Channel%d, %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iC,
471  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
472  h_timeFSLPC_tc[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
473 
474  hn = Form("h_timeF%d_S%d_L%d_P%d_C%d", iF, iS, iL, iP, iC);
475  ht = Form("Time distribtution for RPC of Channel%d, %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation[ns]", iC,
476  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
477  h_timeFSLPC[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryRPC, m_UpperTimeBoundaryRPC);
478 
479  hn = Form("hc_timeF%d_S%d_L%d_P%d_C%d", iF, iS, iL, iP, iC);
480  ht = Form("Calibrated time distribtution for RPC of Channel%d, %s, Layer%d, Sector%d, %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
481  iC, iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
482  hc_timeFSLPC[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin, m_LowerTimeBoundaryCalibratedRPC,
484  }
485  }
486  }
487  }
488  // Endcap part
489  int maxLay = 12 + 2 * iF;
490  for (int iS = 0; iS < 4; ++iS) {
491  hn = Form("h_timeF%d_S%d_scint_end", iF, iS);
492  ht = Form("Time distribtution for Scintillator of Sector%d, %s (Endcap); T_rec-T_0-T_fly-T_propagation[ns]", iS,
493  iFstring[iF].Data());
494  h_timeFS_scint_end[iF][iS] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
496  hn = Form("h2_timeF%d_S%d_end", iF, iS);
497  ht = Form("Time distribtution of Sector%d, %s (Endcap); Layer Index; T_rec-T_0-T_fly-T_propagation[ns]", iS, iFstring[iF].Data());
498  h2_timeFS_end[iF][iS] = new TH2F(hn.Data(), ht.Data(), maxLay, 0, maxLay, nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
500  hn = Form("hc_timeF%d_S%d_scint_end", iF, iS);
501  ht = Form("Calibrated time distribtution for Scintillator of Sector%d (Endcap), %s; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
502  iS, iFstring[iF].Data());
503  hc_timeFS_scint_end[iF][iS] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsEKLM,
505  hn = Form("h2c_timeF%d_S%d_end", iF, iS);
506  ht = Form("Calibrated time distribtution of Sector%d, %s (Endcap); Layer Index; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
507  iS, iFstring[iF].Data());
508  h2c_timeFS_end[iF][iS] = new TH2F(hn.Data(), ht.Data(), maxLay, 0, maxLay, nBin_scint,
511 
512  for (int iL = 0; iL < maxLay; ++iL) {
513  hn = Form("h_timeF%d_S%d_L%d_end", iF, iS, iL);
514  ht = Form("Time distribtution for Scintillator of Layer%d, Sector%d, %s (Endcap); T_rec-T_0-T_fly-T_propagation[ns]", iL, iS,
515  iFstring[iF].Data());
516  h_timeFSL_end[iF][iS][iL] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
518  hn = Form("hc_timeF%d_S%d_L%d_end", iF, iS, iL);
519  ht = Form("Calibrated time distribtution for Scintillator of Layer%d, Sector%d, %s (Endcap); T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
520  iL, iS, iFstring[iF].Data());
521  hc_timeFSL_end[iF][iS][iL] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsEKLM,
523 
524  for (int iP = 0; iP < 2; ++iP) {
525  hn = Form("h_timeF%d_S%d_L%d_P%d_end", iF, iS, iL, iP);
526  ht = Form("Time distribtution for Scintillator of %s, Layer%d, Sector%d, %s (Endcap); T_rec-T_0-T_fly-T_propagation[ns]",
527  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
528  h_timeFSLP_end[iF][iS][iL][iP] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
530 
531  hn = Form("h2_timeF%d_S%d_L%d_P%d_end", iF, iS, iL, iP);
532  ht = Form("Time distribtution for Scintillator of %s, Layer%d, Sector%d, %s (Endcap); Channel Index; T_rec-T_0-T_fly-T_propagation[ns]",
533  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
534  h2_timeFSLP_end[iF][iS][iL][iP] = new TH2F(hn.Data(), ht.Data(), 75, 0, 75, nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
536 
537  hn = Form("hc_timeF%d_S%d_L%d_P%d_end", iF, iS, iL, iP);
538  ht = Form("Calibrated time distribtution for Scintillator of %s, Layer%d, Sector%d, %s (Endcap); T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
539  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
540  hc_timeFSLP_end[iF][iS][iL][iP] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsEKLM,
542 
543  hn = Form("h2c_timeF%d_S%d_L%d_P%d_end", iF, iS, iL, iP);
544  ht = Form("Calibrated time distribtution for Scintillator of %s, Layer%d, Sector%d, %s (Endcap); Channel Index; T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
545  iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
546  h2c_timeFSLP_end[iF][iS][iL][iP] = new TH2F(hn.Data(), ht.Data(), 75, 0, 75, nBin_scint,
549 
550  for (int iC = 0; iC < 75; ++iC) {
551  hn = Form("h_timeF%d_S%d_L%d_P%d_C%d_tc_end", iF, iS, iL, iP, iC);
552  ht = Form("Time distribtution for Scintillator of Channel%d, %s, Layer%d, Sector%d, %s (Endcap); T_rec-T_0-T_fly-T_propagation[ns]",
553  iC, iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
554  h_timeFSLPC_tc_end[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
556 
557  hn = Form("h_timeF%d_S%d_L%d_P%d_C%d_end", iF, iS, iL, iP, iC);
558  ht = Form("Time distribtution for Scintillator of Channel%d, %s, Layer%d, Sector%d, %s (Endcap); T_rec-T_0-T_fly-T_propagation[ns]",
559  iC, iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
560  h_timeFSLPC_end[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryScintilltorsEKLM,
562 
563  hn = Form("hc_timeF%d_S%d_L%d_P%d_C%d_end", iF, iS, iL, iP, iC);
564  ht = Form("Calibrated time distribtution for Scintillator of Channel%d, %s, Layer%d, Sector%d, %s (Endcap); T_rec-T_0-T_fly-T_propagation-T_calibration[ns]",
565  iC, iPstring[iP].Data(), iL, iS, iFstring[iF].Data());
566  hc_timeFSLPC_end[iF][iS][iL][iP][iC] = new TH1F(hn.Data(), ht.Data(), nBin_scint, m_LowerTimeBoundaryCalibratedScintilltorsEKLM,
568  hn = Form("time_length_eklm_F%d_S%d_L%d_P%d_C%d", iF, iS, iL, iP, iC);
569  double stripLength = m_EKLMGeometry->getStripLength(iC + 1) /
570  CLHEP::cm * Unit::cm;
571  m_HistTimeLengthEKLM[iF][iS][iL][iP][iC] =
572  new TH2F(hn.Data(),
573  "Time versus propagation length; "
574  "propagation distance[cm]; "
575  "T_rec-T_0-T_fly-'T_calibration'[ns]",
576  200, 0.0, stripLength,
579  }
580  }
581  }
582  }
583  }
584 }
585 
587  TProfile* profileRpcPhi, TProfile* profileRpcZ,
588  TProfile* profileBKLMScintillatorPhi, TProfile* profileBKLMScintillatorZ,
589  TProfile* profileEKLMScintillatorPlane1,
590  TProfile* profileEKLMScintillatorPlane2, bool fill2dHistograms)
591 {
592  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
593  KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
594  if (m_cFlag[channel] == ChannelCalibrationStatus::c_NotEnoughData)
595  continue;
596 
597  std::vector<struct Event>::iterator it;
598  std::vector<struct Event> eventsChannel;
599  eventsChannel = m_evts[channel];
600  int iSub = klmChannel.getSubdetector();
601 
602  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
603  double timeHit = it->time() - m_timeShift[channel];
604  if (m_useEventT0)
605  timeHit = timeHit - it->t0;
606  double distHit = it->dist;
607 
608  if (iSub == KLMElementNumbers::c_BKLM) {
609  int iF = klmChannel.getSection();
610  int iS = klmChannel.getSector() - 1;
611  int iL = klmChannel.getLayer() - 1;
612  int iP = klmChannel.getPlane();
613  int iC = klmChannel.getStrip() - 1;
614  if (iL > 1) {
615  if (iP) {
616  profileRpcPhi->Fill(distHit, timeHit);
617  } else {
618  profileRpcZ->Fill(distHit, timeHit);
619  }
620  } else {
621  if (fill2dHistograms)
622  m_HistTimeLengthBKLM[iF][iS][iL][iP][iC]->Fill(distHit, timeHit);
623  if (iP) {
624  profileBKLMScintillatorPhi->Fill(distHit, timeHit);
625  } else {
626  profileBKLMScintillatorZ->Fill(distHit, timeHit);
627  }
628  }
629  } else {
630  int iF = klmChannel.getSection() - 1;
631  int iS = klmChannel.getSector() - 1;
632  int iL = klmChannel.getLayer() - 1;
633  int iP = klmChannel.getPlane() - 1;
634  int iC = klmChannel.getStrip() - 1;
635  if (fill2dHistograms)
636  m_HistTimeLengthEKLM[iF][iS][iL][iP][iC]->Fill(distHit, timeHit);
637  if (iP) {
638  profileEKLMScintillatorPlane1->Fill(distHit, timeHit);
639  } else {
640  profileEKLMScintillatorPlane2->Fill(distHit, timeHit);
641  }
642  }
643  }
644  }
645 }
646 
648  const std::vector< std::pair<KLMChannelNumber, unsigned int> >& channels,
649  double& delay, double& delayError)
650 {
651  std::vector<struct Event>::iterator it;
652  std::vector<struct Event> eventsChannel;
653  int nFits = 1000;
654  int nConvergedFits = 0;
655  delay = 0;
656  delayError = 0;
657  if (nFits > (int)channels.size())
658  nFits = channels.size();
659  for (int i = 0; i < nFits; ++i) {
660  int subdetector, section, sector, layer, plane, strip;
662  channels[i].first, &subdetector, &section, &sector, &layer, &plane,
663  &strip);
664  if (subdetector == KLMElementNumbers::c_BKLM) {
665  s_LowerTimeBoundary = m_LowerTimeBoundaryScintilltorsBKLM;
666  s_UpperTimeBoundary = m_UpperTimeBoundaryScintilltorsBKLM;
667  const bklm::Module* module =
668  m_BKLMGeometry->findModule(section, sector, layer);
669  s_StripLength = module->getStripLength(plane, strip);
670  } else {
671  s_LowerTimeBoundary = m_LowerTimeBoundaryScintilltorsEKLM;
672  s_UpperTimeBoundary = m_UpperTimeBoundaryScintilltorsEKLM;
673  s_StripLength = m_EKLMGeometry->getStripLength(strip) /
674  CLHEP::cm * Unit::cm;
675  }
676  for (int j = 0; j < c_NBinsTime; ++j) {
677  for (int k = 0; k < c_NBinsDistance; ++k)
678  s_BinnedData[j][k] = 0;
679  }
680  eventsChannel = m_evts[channels[i].first];
681  double averageTime = 0;
682  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
683  double timeHit = it->time();
684  if (m_useEventT0)
685  timeHit = timeHit - it->t0;
686  averageTime = averageTime + timeHit;
687  int timeBin = std::floor((timeHit - s_LowerTimeBoundary) * c_NBinsTime /
688  (s_UpperTimeBoundary - s_LowerTimeBoundary));
689  if (timeBin < 0 || timeBin >= c_NBinsTime)
690  continue;
691  int distanceBin = std::floor(it->dist * c_NBinsDistance / s_StripLength);
692  if (distanceBin < 0 || distanceBin >= c_NBinsDistance) {
693  B2ERROR("The distance to SiPM is greater than the strip length.");
694  continue;
695  }
696  s_BinnedData[timeBin][distanceBin] += 1;
697  }
698  averageTime = averageTime / eventsChannel.size();
699  TMinuit minuit(5);
700  minuit.SetPrintLevel(-1);
701  int minuitResult;
702  minuit.mnparm(0, "P0", 1, 0.001, 0, 0, minuitResult);
703  minuit.mnparm(1, "N", 10, 0.001, 0, 0, minuitResult);
704  minuit.mnparm(2, "T0", averageTime, 0.001, 0, 0, minuitResult);
705  minuit.mnparm(3, "SIGMA", 10, 0.001, 0, 0, minuitResult);
706  minuit.mnparm(4, "DELAY", 0.0, 0.001, 0, 0, minuitResult);
707  minuit.SetFCN(fcn);
708  minuit.mncomd("FIX 2 3 4 5", minuitResult);
709  minuit.mncomd("MIGRAD 10000", minuitResult);
710  minuit.mncomd("RELEASE 2", minuitResult);
711  minuit.mncomd("MIGRAD 10000", minuitResult);
712  minuit.mncomd("RELEASE 3", minuitResult);
713  minuit.mncomd("MIGRAD 10000", minuitResult);
714  minuit.mncomd("RELEASE 4", minuitResult);
715  minuit.mncomd("MIGRAD 10000", minuitResult);
716  minuit.mncomd("RELEASE 5", minuitResult);
717  minuit.mncomd("MIGRAD 10000", minuitResult);
718  /* Require converged fit with accurate error matrix. */
719  if (minuit.fISW[1] != 3)
720  continue;
721  nConvergedFits++;
722  double channelDelay, channelDelayError;
723  minuit.GetParameter(4, channelDelay, channelDelayError);
724  delay = delay + channelDelay;
725  delayError = delayError + channelDelayError * channelDelayError;
726  }
727  delay = delay / nConvergedFits;
728  delayError = sqrt(delayError) / (nConvergedFits - 1);
729 }
730 
732 {
733  int channelId;
734  gROOT->SetBatch(kTRUE);
735  setupDatabase();
738 
739  fcn_gaus = new TF1("fcn_gaus", "gaus");
740  fcn_land = new TF1("fcn_land", "landau");
741  fcn_pol1 = new TF1("fcn_pol1", "pol1");
742  fcn_const = new TF1("fcn_const", "pol0");
743 
745  if (result != CalibrationAlgorithm::c_OK)
746  return result;
747 
748  /* Choose non-existing file name. */
749  std::string name = "time_calibration.root";
750  int i = 1;
751  while (1) {
752  struct stat buffer;
753  if (stat(name.c_str(), &buffer) != 0)
754  break;
755  name = "time_calibration_" + std::to_string(i) + ".root";
756  i = i + 1;
757  /* Overflow. */
758  if (i < 0)
759  break;
760  }
761  m_outFile = new TFile(name.c_str(), "recreate");
763 
764  std::vector<struct Event>::iterator it;
765  std::vector<struct Event> eventsChannel;
766 
767  eventsChannel.clear();
768  m_cFlag.clear();
769  m_minimizerOptions.SetDefaultStrategy(2);
770 
771  /* Sort channels by number of events. */
772  std::vector< std::pair<KLMChannelNumber, unsigned int> > channelsBKLM;
773  std::vector< std::pair<KLMChannelNumber, unsigned int> > channelsEKLM;
774  KLMChannelIndex klmChannels;
775  for (KLMChannelIndex& klmChannel : klmChannels) {
776  KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
777  m_cFlag[channel] = ChannelCalibrationStatus::c_NotEnoughData;
778  if (m_evts.find(channel) == m_evts.end())
779  continue;
780  int nEvents = m_evts[channel].size();
781  if (nEvents < m_lower_limit_counts) {
782  B2WARNING("Not enough calibration data collected."
783  << LogVar("channel", channel)
784  << LogVar("number of digit", nEvents));
785  continue;
786  }
787  m_cFlag[channel] = ChannelCalibrationStatus::c_FailedFit;
788  if (klmChannel.getSubdetector() == KLMElementNumbers::c_BKLM &&
789  klmChannel.getLayer() < BKLMElementNumbers::c_FirstRPCLayer) {
790  channelsBKLM.push_back(
791  std::pair<KLMChannelNumber, unsigned int>(channel, nEvents));
792  }
793  if (klmChannel.getSubdetector() == KLMElementNumbers::c_EKLM) {
794  channelsEKLM.push_back(
795  std::pair<KLMChannelNumber, unsigned int>(channel, nEvents));
796  }
797  }
798  std::sort(channelsBKLM.begin(), channelsBKLM.end(), compareEventNumber);
799  std::sort(channelsEKLM.begin(), channelsEKLM.end(), compareEventNumber);
800 
801  /* Two-dimensional fit for the channel with the maximal number of events. */
802  double delayBKLM, delayBKLMError;
803  double delayEKLM, delayEKLMError;
804  timeDistance2dFit(channelsBKLM, delayBKLM, delayBKLMError);
805  timeDistance2dFit(channelsEKLM, delayEKLM, delayEKLMError);
806 
807  /**********************************
808  * First loop
809  * Estimation of effective light speed for Scintillators and RPCs, separately.
810  **********************************/
811  B2INFO("Effective light speed Estimation.");
812  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
813  channelId = klmChannel.getKLMChannelNumber();
814  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
815  continue;
816 
817  eventsChannel = m_evts[channelId];
818  int iSub = klmChannel.getSubdetector();
819 
820  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
821  TVector3 diffD = TVector3(it->diffDistX, it->diffDistY, it->diffDistZ);
822  h_diff->Fill(diffD.Mag());
823  double timeHit = it->time();
824  if (m_useEventT0)
825  timeHit = timeHit - it->t0;
826  if (iSub == KLMElementNumbers::c_BKLM) {
827  int iF = klmChannel.getSection();
828  int iS = klmChannel.getSector() - 1;
829  int iL = klmChannel.getLayer() - 1;
830  int iP = klmChannel.getPlane();
831  int iC = klmChannel.getStrip() - 1;
832  h_timeFSLPC_tc[iF][iS][iL][iP][iC]->Fill(timeHit);
833  if (iL > 1) {
834  h_time_rpc_tc->Fill(timeHit);
835  } else {
836  h_time_scint_tc->Fill(timeHit);
837  }
838  } else {
839  int iF = klmChannel.getSection() - 1;
840  int iS = klmChannel.getSector() - 1;
841  int iL = klmChannel.getLayer() - 1;
842  int iP = klmChannel.getPlane() - 1;
843  int iC = klmChannel.getStrip() - 1;
844  h_timeFSLPC_tc_end[iF][iS][iL][iP][iC]->Fill(timeHit);
845  h_time_scint_tc_end->Fill(timeHit);
846  }
847  }
848  }
849  B2INFO("Effective light speed Estimation! Hists and Graph filling done.");
850 
851  m_timeShift.clear();
852 
853  double tmpMean_rpc_global = h_time_rpc_tc->GetMean();
854  double tmpMean_scint_global = h_time_scint_tc->GetMean();
855  double tmpMean_scint_global_end = h_time_scint_tc_end->GetMean();
856 
857  B2INFO("Global Mean for Raw." << LogVar("RPC", tmpMean_rpc_global) << LogVar("Scint BKLM",
858  tmpMean_scint_global) << LogVar("Scint EKLM", tmpMean_scint_global_end));
859 
860  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
861  channelId = klmChannel.getKLMChannelNumber();
862  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
863  continue;
864 
865  int iSub = klmChannel.getSubdetector();
866  if (iSub == KLMElementNumbers::c_BKLM) {
867  int iF = klmChannel.getSection();
868  int iS = klmChannel.getSector() - 1;
869  int iL = klmChannel.getLayer() - 1;
870  int iP = klmChannel.getPlane();
871  int iC = klmChannel.getStrip() - 1;
872  h_timeFSLPC_tc[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
873  double tmpMean_channel = fcn_gaus->GetParameter(1);
874  if (iL > 1) {
875  m_timeShift[channelId] = tmpMean_channel - tmpMean_rpc_global;
876  } else {
877  m_timeShift[channelId] = tmpMean_channel - tmpMean_scint_global;
878  }
879  } else {
880  int iF = klmChannel.getSection() - 1;
881  int iS = klmChannel.getSector() - 1;
882  int iL = klmChannel.getLayer() - 1;
883  int iP = klmChannel.getPlane() - 1;
884  int iC = klmChannel.getStrip() - 1;
885  h_timeFSLPC_tc_end[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
886  double tmpMean_channel = fcn_gaus->GetParameter(1);
887  m_timeShift[channelId] = tmpMean_channel - tmpMean_scint_global_end;
888  }
889  }
890 
891  delete h_time_scint_tc;
892  delete h_time_scint_tc_end;
893  delete h_time_rpc_tc;
894  B2INFO("Effective Light m_timeShift obtained. done.");
895 
900 
901  B2INFO("Effective light speed fitting.");
902  m_ProfileRpcPhi->Fit("fcn_pol1", "EMQ");
903  double delayRPCPhi = fcn_pol1->GetParameter(1);
904  double e_slope_rpc_phi = fcn_pol1->GetParError(1);
905 
906  m_ProfileRpcZ->Fit("fcn_pol1", "EMQ");
907  double delayRPCZ = fcn_pol1->GetParameter(1);
908  double e_slope_rpc_z = fcn_pol1->GetParError(1);
909 
910  m_ProfileBKLMScintillatorPhi->Fit("fcn_pol1", "EMQ");
911  double slope_scint_phi = fcn_pol1->GetParameter(1);
912  double e_slope_scint_phi = fcn_pol1->GetParError(1);
913 
914  m_ProfileBKLMScintillatorZ->Fit("fcn_pol1", "EMQ");
915  double slope_scint_z = fcn_pol1->GetParameter(1);
916  double e_slope_scint_z = fcn_pol1->GetParError(1);
917 
918  m_ProfileEKLMScintillatorPlane1->Fit("fcn_pol1", "EMQ");
919  double slope_scint_plane1_end = fcn_pol1->GetParameter(1);
920  double e_slope_scint_plane1_end = fcn_pol1->GetParError(1);
921 
922  m_ProfileEKLMScintillatorPlane2->Fit("fcn_pol1", "EMQ");
923  double slope_scint_plane2_end = fcn_pol1->GetParameter(1);
924  double e_slope_scint_plane2_end = fcn_pol1->GetParError(1);
925 
926  TString logStr_phi, logStr_z;
927  logStr_phi = Form("%.4f +/- %.4f ns/cm", delayRPCPhi, e_slope_rpc_phi);
928  logStr_z = Form("%.4f +/- %.4f ns/cm", delayRPCZ, e_slope_rpc_z);
929  B2INFO("Delay in RPCs:"
930  << LogVar("Fitted Value (phi readout) ", logStr_phi.Data())
931  << LogVar("Fitted Value (z readout) ", logStr_z.Data()));
932  logStr_phi = Form("%.4f +/- %.4f ns/cm", slope_scint_phi, e_slope_scint_phi);
933  logStr_z = Form("%.4f +/- %.4f ns/cm", slope_scint_z, e_slope_scint_z);
934  B2INFO("Delay in BKLM scintillators:"
935  << LogVar("Fitted Value (phi readout) ", logStr_phi.Data())
936  << LogVar("Fitted Value (z readout) ", logStr_z.Data()));
937  logStr_phi = Form("%.4f +/- %.4f ns/cm", slope_scint_plane1_end,
938  e_slope_scint_plane1_end);
939  logStr_z = Form("%.4f +/- %.4f ns/cm", slope_scint_plane2_end,
940  e_slope_scint_plane2_end);
941  B2INFO("Delay in EKLM scintillators:"
942  << LogVar("Fitted Value (plane1 readout) ", logStr_phi.Data())
943  << LogVar("Fitted Value (plane2 readout) ", logStr_z.Data()));
944 
945  logStr_z = Form("%.4f +/- %.4f ns/cm", delayBKLM, delayBKLMError);
946  B2INFO("Delay in BKLM scintillators:"
947  << LogVar("Fitted Value (2d fit) ", logStr_z.Data()));
948  logStr_z = Form("%.4f +/- %.4f ns/cm", delayEKLM, delayEKLMError);
949  B2INFO("Delay in EKLM scintillators:"
950  << LogVar("Fitted Value (2d fit) ", logStr_z.Data()));
951 
952  // Default Effective light speed in current Database
953  //delayEKLM = 0.5 * (slope_scint_plane1_end + slope_scint_plane2_end);
954  //delayBKLM = 0.5 * (slope_scint_phi + slope_scint_z);
955 
960 
962  B2INFO("Time distribution filling begins.");
963  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
964  channelId = klmChannel.getKLMChannelNumber();
965  int iSub = klmChannel.getSubdetector();
966 
967  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
968  continue;
969  eventsChannel = m_evts[channelId];
970 
971  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
972  double timeHit = it->time();
973  if (m_useEventT0)
974  timeHit = timeHit - it->t0;
975  if (iSub == KLMElementNumbers::c_BKLM) {
976  int iF = klmChannel.getSection();
977  int iS = klmChannel.getSector() - 1;
978  int iL = klmChannel.getLayer() - 1;
979  int iP = klmChannel.getPlane();
980  int iC = klmChannel.getStrip() - 1;
981  if (iL > 1) {
982  double propgationT;
984  propgationT = it->dist * delayRPCZ;
985  else
986  propgationT = it->dist * delayRPCPhi;
987  double time = timeHit - propgationT;
988  h_time_rpc->Fill(time);
989  h_timeF_rpc[iF]->Fill(time);
990  h_timeFS_rpc[iF][iS]->Fill(time);
991  h_timeFSL[iF][iS][iL]->Fill(time);
992  h_timeFSLP[iF][iS][iL][iP]->Fill(time);
993  h_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
994  h2_timeF_rpc[iF]->Fill(iS, time);
995  h2_timeFS[iF][iS]->Fill(iL, time);
996  h2_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
997  } else {
998  double propgationT = it->dist * delayBKLM;
999  double time = timeHit - propgationT;
1000  h_time_scint->Fill(time);
1001  h_timeF_scint[iF]->Fill(time);
1002  h_timeFS_scint[iF][iS]->Fill(time);
1003  h_timeFSL[iF][iS][iL]->Fill(time);
1004  h_timeFSLP[iF][iS][iL][iP]->Fill(time);
1005  h_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1006  h2_timeF_scint[iF]->Fill(iS, time);
1007  h2_timeFS[iF][iS]->Fill(iL, time);
1008  h2_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1009  }
1010  } else {
1011  int iF = klmChannel.getSection() - 1;
1012  int iS = klmChannel.getSector() - 1;
1013  int iL = klmChannel.getLayer() - 1;
1014  int iP = klmChannel.getPlane() - 1;
1015  int iC = klmChannel.getStrip() - 1;
1016  double propgationT = it->dist * delayEKLM;
1017  double time = timeHit - propgationT;
1018  h_time_scint_end->Fill(time);
1019  h_timeF_scint_end[iF]->Fill(time);
1020  h_timeFS_scint_end[iF][iS]->Fill(time);
1021  h_timeFSL_end[iF][iS][iL]->Fill(time);
1022  h_timeFSLP_end[iF][iS][iL][iP]->Fill(time);
1023  h_timeFSLPC_end[iF][iS][iL][iP][iC]->Fill(time);
1024  h2_timeF_scint_end[iF]->Fill(iS, time);
1025  h2_timeFS_end[iF][iS]->Fill(iL, time);
1026  h2_timeFSLP_end[iF][iS][iL][iP]->Fill(iC, time);
1027  }
1028  }
1029  }
1030 
1031  B2INFO("Orignal filling done.");
1032 
1033  int iChannel_rpc = 0;
1034  int iChannel = 0;
1035  int iChannel_end = 0;
1036  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1037  channelId = klmChannel.getKLMChannelNumber();
1038  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
1039  continue;
1040  int iSub = klmChannel.getSubdetector();
1041 
1042  if (iSub == KLMElementNumbers::c_BKLM) {
1043  int iF = klmChannel.getSection();
1044  int iS = klmChannel.getSector() - 1;
1045  int iL = klmChannel.getLayer() - 1;
1046  int iP = klmChannel.getPlane();
1047  int iC = klmChannel.getStrip() - 1;
1048 
1049  TFitResultPtr r = h_timeFSLPC[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
1050  if (int(r) != 0)
1051  continue;
1052  if (int(r) == 0)
1053  m_cFlag[channelId] = ChannelCalibrationStatus::c_SuccessfulCalibration;
1054  m_time_channel[channelId] = fcn_gaus->GetParameter(1);
1055  m_etime_channel[channelId] = fcn_gaus->GetParError(1);
1056  if (iL > 1) {
1057  gre_time_channel_rpc->SetPoint(iChannel, channelId, m_time_channel[channelId]);
1058  gre_time_channel_rpc->SetPointError(iChannel, 0., m_etime_channel[channelId]);
1059  iChannel++;
1060  } else {
1061  gre_time_channel_scint->SetPoint(iChannel_rpc, channelId, m_time_channel[channelId]);
1062  gre_time_channel_scint->SetPointError(iChannel_rpc, 0., m_etime_channel[channelId]);
1063  iChannel_rpc++;
1064  }
1065  } else {
1066  int iF = klmChannel.getSection() - 1;
1067  int iS = klmChannel.getSector() - 1;
1068  int iL = klmChannel.getLayer() - 1;
1069  int iP = klmChannel.getPlane() - 1;
1070  int iC = klmChannel.getStrip() - 1;
1071 
1072  TFitResultPtr r = h_timeFSLPC_end[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
1073  if (int(r) != 0)
1074  continue;
1075  if (int(r) == 0)
1076  m_cFlag[channelId] = ChannelCalibrationStatus::c_SuccessfulCalibration;
1077  m_time_channel[channelId] = fcn_gaus->GetParameter(1);
1078  m_etime_channel[channelId] = fcn_gaus->GetParError(1);
1079  gre_time_channel_scint_end->SetPoint(iChannel_end, channelId, m_time_channel[channelId]);
1080  gre_time_channel_scint_end->SetPointError(iChannel_end, 0., m_etime_channel[channelId]);
1081  iChannel_end++;
1082  }
1083  }
1084 
1085  gre_time_channel_scint->Fit("fcn_const", "EMQ");
1086  m_time_channelAvg_scint = fcn_const->GetParameter(0);
1087  m_etime_channelAvg_scint = fcn_const->GetParError(0);
1088 
1089  gre_time_channel_scint_end->Fit("fcn_const", "EMQ");
1090  m_time_channelAvg_scint_end = fcn_const->GetParameter(0);
1091  m_etime_channelAvg_scint_end = fcn_const->GetParError(0);
1092 
1093  gre_time_channel_rpc->Fit("fcn_const", "EMQ");
1094  m_time_channelAvg_rpc = fcn_const->GetParameter(0);
1095  m_etime_channelAvg_rpc = fcn_const->GetParError(0);
1096 
1097  B2INFO("Channel's time distribution fitting done.");
1098  B2DEBUG(20, LogVar("Average time (RPC)", m_time_channelAvg_rpc)
1099  << LogVar("Average time (BKLM scintillators)", m_time_channelAvg_scint)
1100  << LogVar("Average time (EKLM scintillators)", m_time_channelAvg_scint_end));
1101 
1102  B2INFO("Calibrated channel's time distribution filling begins.");
1103 
1104  m_timeShift.clear();
1105  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1106  channelId = klmChannel.getKLMChannelNumber();
1107  h_calibrated->Fill(m_cFlag[channelId]);
1108  if (m_time_channel.find(channelId) == m_time_channel.end())
1109  continue;
1110  double timeShift = m_time_channel[channelId];
1111  int iSub = klmChannel.getSubdetector();
1112  if (iSub == KLMElementNumbers::c_BKLM) {
1113  int iL = klmChannel.getLayer() - 1;
1114  if (iL > 1)
1115  m_timeShift[channelId] = timeShift;
1116  else
1117  m_timeShift[channelId] = timeShift;
1118  } else {
1119  m_timeShift[channelId] = timeShift;
1120  }
1121  m_timeCableDelay->setTimeDelay(channelId, m_timeShift[channelId]);
1122  }
1123 
1124  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1125  channelId = klmChannel.getKLMChannelNumber();
1126  if (m_timeShift.find(channelId) != m_timeShift.end())
1127  continue;
1128  m_timeShift[channelId] = esti_timeShift(klmChannel);
1129  m_timeCableDelay->setTimeDelay(channelId, m_timeShift[channelId]);
1130  B2DEBUG(20, "Uncalibrated Estimation " << LogVar("Channel", channelId) << LogVar("Estimated value", m_timeShift[channelId]));
1131  }
1132 
1133  iChannel_rpc = 0;
1134  iChannel = 0;
1135  iChannel_end = 0;
1136  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1137  channelId = klmChannel.getKLMChannelNumber();
1138  if (m_timeShift.find(channelId) == m_timeShift.end()) {
1139  B2ERROR("!!! Not All Channels Calibration Constant Set. Error Happended on " << LogVar("Channel", channelId));
1140  continue;
1141  }
1142  int iSub = klmChannel.getSubdetector();
1143  if (iSub == KLMElementNumbers::c_BKLM) {
1144  int iL = klmChannel.getLayer();
1145  if (iL > 2) {
1146  gr_timeShift_channel_rpc->SetPoint(iChannel_rpc, channelId, m_timeShift[channelId]);
1147  iChannel_rpc++;
1148  } else {
1149  gr_timeShift_channel_scint->SetPoint(iChannel, channelId, m_timeShift[channelId]);
1150  iChannel++;
1151  }
1152  } else {
1153  gr_timeShift_channel_scint_end->SetPoint(iChannel_end, channelId, m_timeShift[channelId]);
1154  iChannel_end++;
1155  }
1156  }
1157 
1162  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1163  channelId = klmChannel.getKLMChannelNumber();
1164  int iSub = klmChannel.getSubdetector();
1165  eventsChannel = m_evts[channelId];
1166  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
1167  double timeHit = it->time();
1168  if (m_useEventT0)
1169  timeHit = timeHit - it->t0;
1170  if (iSub == KLMElementNumbers::c_BKLM) {
1171  int iF = klmChannel.getSection();
1172  int iS = klmChannel.getSector() - 1;
1173  int iL = klmChannel.getLayer() - 1;
1174  int iP = klmChannel.getPlane();
1175  int iC = klmChannel.getStrip() - 1;
1176  if (iL > 1) {
1177  double propgationT;
1178  if (iP == BKLMElementNumbers::c_ZPlane)
1179  propgationT = it->dist * delayRPCZ;
1180  else
1181  propgationT = it->dist * delayRPCPhi;
1182  double time = timeHit - propgationT - m_timeShift[channelId];
1183  hc_time_rpc->Fill(time);
1184  hc_timeF_rpc[iF]->Fill(time);
1185  hc_timeFS_rpc[iF][iS]->Fill(time);
1186  hc_timeFSL[iF][iS][iL]->Fill(time);
1187  hc_timeFSLP[iF][iS][iL][iP]->Fill(time);
1188  hc_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1189  h2c_timeF_rpc[iF]->Fill(iS, time);
1190  h2c_timeFS[iF][iS]->Fill(iL, time);
1191  h2c_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1192  } else {
1193  double propgationT = it->dist * delayBKLM;
1194  double time = timeHit - propgationT - m_timeShift[channelId];
1195  hc_time_scint->Fill(time);
1196  hc_timeF_scint[iF]->Fill(time);
1197  hc_timeFS_scint[iF][iS]->Fill(time);
1198  hc_timeFSL[iF][iS][iL]->Fill(time);
1199  hc_timeFSLP[iF][iS][iL][iP]->Fill(time);
1200  hc_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1201  h2c_timeF_scint[iF]->Fill(iS, time);
1202  h2c_timeFS[iF][iS]->Fill(iL, time);
1203  h2c_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1204  }
1205  } else {
1206  int iF = klmChannel.getSection() - 1;
1207  int iS = klmChannel.getSector() - 1;
1208  int iL = klmChannel.getLayer() - 1;
1209  int iP = klmChannel.getPlane() - 1;
1210  int iC = klmChannel.getStrip() - 1;
1211  double propgationT = it->dist * delayEKLM;
1212  double time = timeHit - propgationT - m_timeShift[channelId];
1213  hc_time_scint_end->Fill(time);
1214  hc_timeF_scint_end[iF]->Fill(time);
1215  hc_timeFS_scint_end[iF][iS]->Fill(time);
1216  hc_timeFSL_end[iF][iS][iL]->Fill(time);
1217  hc_timeFSLP_end[iF][iS][iL][iP]->Fill(time);
1218  hc_timeFSLPC_end[iF][iS][iL][iP][iC]->Fill(time);
1219  h2c_timeF_scint_end[iF]->Fill(iS, time);
1220  h2c_timeFS_end[iF][iS]->Fill(iL, time);
1221  h2c_timeFSLP_end[iF][iS][iL][iP]->Fill(iC, time);
1222  }
1223  }
1224  }
1225 
1226  saveHist();
1227 
1228  delete fcn_const;
1229  m_evts.clear();
1230  m_timeShift.clear();
1231  m_cFlag.clear();
1232 
1233  saveCalibration(m_timeCableDelay, "KLMTimeCableDelay");
1234  saveCalibration(m_timeConstants, "KLMTimeConstants");
1236 }
1237 
1238 
1240 {
1241  m_outFile->cd();
1242  B2INFO("Save Histograms into Files.");
1243  TDirectory* dir_monitor = m_outFile->mkdir("monitor_Hists");
1244  dir_monitor->cd();
1245  h_calibrated->SetDirectory(dir_monitor);
1246  h_diff->SetDirectory(dir_monitor);
1247 
1248  m_outFile->cd();
1249  TDirectory* dir_effC = m_outFile->mkdir("effC_Hists");
1250  dir_effC->cd();
1251  m_ProfileRpcPhi->SetDirectory(dir_effC);
1252  m_ProfileRpcZ->SetDirectory(dir_effC);
1253  m_ProfileBKLMScintillatorPhi->SetDirectory(dir_effC);
1254  m_ProfileBKLMScintillatorZ->SetDirectory(dir_effC);
1255  m_ProfileEKLMScintillatorPlane1->SetDirectory(dir_effC);
1256  m_ProfileEKLMScintillatorPlane2->SetDirectory(dir_effC);
1257  m_Profile2RpcPhi->SetDirectory(dir_effC);
1258  m_Profile2RpcZ->SetDirectory(dir_effC);
1259  m_Profile2BKLMScintillatorPhi->SetDirectory(dir_effC);
1260  m_Profile2BKLMScintillatorZ->SetDirectory(dir_effC);
1261  m_Profile2EKLMScintillatorPlane1->SetDirectory(dir_effC);
1262  m_Profile2EKLMScintillatorPlane2->SetDirectory(dir_effC);
1263 
1264  m_outFile->cd();
1265  TDirectory* dir_time = m_outFile->mkdir("time");
1266  dir_time->cd();
1267 
1268  h_time_scint->SetDirectory(dir_time);
1269  hc_time_scint->SetDirectory(dir_time);
1270 
1271  h_time_scint_end->SetDirectory(dir_time);
1272  hc_time_scint_end->SetDirectory(dir_time);
1273 
1274  h_time_rpc->SetDirectory(dir_time);
1275  hc_time_rpc->SetDirectory(dir_time);
1276 
1277  gre_time_channel_rpc->Write("gre_time_channel_rpc");
1278  gre_time_channel_scint->Write("gre_time_channel_scint");
1279  gre_time_channel_scint_end->Write("gre_time_channel_scint_end");
1280  gr_timeShift_channel_rpc->Write("gr_timeShift_channel_rpc");
1281  gr_timeShift_channel_scint->Write("gr_timeShift_channel_scint");
1282  gr_timeShift_channel_scint_end->Write("gr_timeShift_channel_scint_end");
1283 
1284  B2INFO("Top file setup Done.");
1285 
1286  TDirectory* dir_time_F[2];
1287  TDirectory* dir_time_FS[2][8];
1288  TDirectory* dir_time_FSL[2][8][15];
1289  TDirectory* dir_time_FSLP[2][8][15][2];
1290  TDirectory* dir_time_F_end[2];
1291  TDirectory* dir_time_FS_end[2][4];
1292  TDirectory* dir_time_FSL_end[2][4][14];
1293  TDirectory* dir_time_FSLP_end[2][4][14][2];
1294  char dirname[50];
1295  B2INFO("Sub files declare Done.");
1296  for (int iF = 0; iF < 2; ++iF) {
1297  h_timeF_rpc[iF]->SetDirectory(dir_time);
1298  hc_timeF_rpc[iF]->SetDirectory(dir_time);
1299 
1300  h2_timeF_rpc[iF]->SetDirectory(dir_time);
1301  h2c_timeF_rpc[iF]->SetDirectory(dir_time);
1302 
1303  h_timeF_scint[iF]->SetDirectory(dir_time);
1304  hc_timeF_scint[iF]->SetDirectory(dir_time);
1305 
1306  h2_timeF_scint[iF]->SetDirectory(dir_time);
1307  h2c_timeF_scint[iF]->SetDirectory(dir_time);
1308 
1309  h_timeF_scint_end[iF]->SetDirectory(dir_time);
1310  hc_timeF_scint_end[iF]->SetDirectory(dir_time);
1311 
1312  h2_timeF_scint_end[iF]->SetDirectory(dir_time);
1313  h2c_timeF_scint_end[iF]->SetDirectory(dir_time);
1314 
1315  sprintf(dirname, "isForward_%d", iF);
1316  dir_time_F[iF] = dir_time->mkdir(dirname);
1317  dir_time_F[iF]->cd();
1318 
1319  for (int iS = 0; iS < 8; ++iS) {
1320  h_timeFS_rpc[iF][iS]->SetDirectory(dir_time_F[iF]);
1321  hc_timeFS_rpc[iF][iS]->SetDirectory(dir_time_F[iF]);
1322 
1323  h_timeFS_scint[iF][iS]->SetDirectory(dir_time_F[iF]);
1324  hc_timeFS_scint[iF][iS]->SetDirectory(dir_time_F[iF]);
1325 
1326  h2_timeFS[iF][iS]->SetDirectory(dir_time_F[iF]);
1327  h2c_timeFS[iF][iS]->SetDirectory(dir_time_F[iF]);
1328 
1329  sprintf(dirname, "Sector_%d", iS + 1);
1330  dir_time_FS[iF][iS] = dir_time_F[iF]->mkdir(dirname);
1331  dir_time_FS[iF][iS]->cd();
1332 
1333  for (int iL = 0; iL < 15; ++iL) {
1334  h_timeFSL[iF][iS][iL]->SetDirectory(dir_time_FS[iF][iS]);
1335  hc_timeFSL[iF][iS][iL]->SetDirectory(dir_time_FS[iF][iS]);
1336 
1337  sprintf(dirname, "Layer_%d", iL + 1);
1338  dir_time_FSL[iF][iS][iL] = dir_time_FS[iF][iS]->mkdir(dirname);
1339  dir_time_FSL[iF][iS][iL]->cd();
1340  for (int iP = 0; iP < 2; ++iP) {
1341  h_timeFSLP[iF][iS][iL][iP]->SetDirectory(dir_time_FSL[iF][iS][iL]);
1342  hc_timeFSLP[iF][iS][iL][iP]->SetDirectory(dir_time_FSL[iF][iS][iL]);
1343  h2_timeFSLP[iF][iS][iL][iP]->SetDirectory(dir_time_FSL[iF][iS][iL]);
1344  h2c_timeFSLP[iF][iS][iL][iP]->SetDirectory(dir_time_FSL[iF][iS][iL]);
1345 
1346  sprintf(dirname, "Plane_%d", iP);
1347  dir_time_FSLP[iF][iS][iL][iP] = dir_time_FSL[iF][iS][iL]->mkdir(dirname);
1348  dir_time_FSLP[iF][iS][iL][iP]->cd();
1349 
1350  int nchannel_max = BKLMElementNumbers::getNStrips(iF, iS + 1, iL + 1, iP);
1351  for (int iC = 0; iC < nchannel_max; ++iC) {
1352  if (iL < 2)
1353  m_HistTimeLengthBKLM[iF][iS][iL][iP][iC]->SetDirectory(dir_time_FSLP[iF][iS][iL][iP]);
1354  h_timeFSLPC[iF][iS][iL][iP][iC]->SetDirectory(dir_time_FSLP[iF][iS][iL][iP]);
1355  hc_timeFSLPC[iF][iS][iL][iP][iC]->SetDirectory(dir_time_FSLP[iF][iS][iL][iP]);
1356  delete h_timeFSLPC_tc[iF][iS][iL][iP][iC];
1357  }
1358  }
1359  }
1360  }
1361 
1362  sprintf(dirname, "isForward_%d_end", iF + 1);
1363  dir_time_F_end[iF] = dir_time->mkdir(dirname);
1364  dir_time_F_end[iF]->cd();
1365  int maxLayer = 12 + 2 * iF;
1366  for (int iS = 0; iS < 4; ++iS) {
1367  h_timeFS_scint_end[iF][iS]->SetDirectory(dir_time_F_end[iF]);
1368  hc_timeFS_scint_end[iF][iS]->SetDirectory(dir_time_F_end[iF]);
1369 
1370  h2_timeFS_end[iF][iS]->SetDirectory(dir_time_F_end[iF]);
1371  h2c_timeFS_end[iF][iS]->SetDirectory(dir_time_F_end[iF]);
1372 
1373  sprintf(dirname, "Sector_%d_end", iS + 1);
1374  dir_time_FS_end[iF][iS] = dir_time_F_end[iF]->mkdir(dirname);
1375  dir_time_FS_end[iF][iS]->cd();
1376  for (int iL = 0; iL < maxLayer; ++iL) {
1377  h_timeFSL_end[iF][iS][iL]->SetDirectory(dir_time_FS_end[iF][iS]);
1378  hc_timeFSL_end[iF][iS][iL]->SetDirectory(dir_time_FS_end[iF][iS]);
1379 
1380  sprintf(dirname, "Layer_%d_end", iL + 1);
1381  dir_time_FSL_end[iF][iS][iL] = dir_time_FS_end[iF][iS]->mkdir(dirname);
1382  dir_time_FSL_end[iF][iS][iL]->cd();
1383  for (int iP = 0; iP < 2; ++iP) {
1384  h_timeFSLP_end[iF][iS][iL][iP]->SetDirectory(dir_time_FSL_end[iF][iS][iL]);
1385  hc_timeFSLP_end[iF][iS][iL][iP]->SetDirectory(dir_time_FSL_end[iF][iS][iL]);
1386  h2_timeFSLP_end[iF][iS][iL][iP]->SetDirectory(dir_time_FSL_end[iF][iS][iL]);
1387  h2c_timeFSLP_end[iF][iS][iL][iP]->SetDirectory(dir_time_FSL_end[iF][iS][iL]);
1388 
1389  sprintf(dirname, "plane_%d_end", iP);
1390  dir_time_FSLP_end[iF][iS][iL][iP] = dir_time_FSL_end[iF][iS][iL]->mkdir(dirname);
1391  dir_time_FSLP_end[iF][iS][iL][iP]->cd();
1392 
1393  for (int iC = 0; iC < 75; ++iC) {
1394  h_timeFSLPC_end[iF][iS][iL][iP][iC]->SetDirectory(dir_time_FSLP_end[iF][iS][iL][iP]);
1395  hc_timeFSLPC_end[iF][iS][iL][iP][iC]->SetDirectory(dir_time_FSLP_end[iF][iS][iL][iP]);
1396  m_HistTimeLengthEKLM[iF][iS][iL][iP][iC]->SetDirectory(dir_time_FSLP_end[iF][iS][iL][iP]);
1397  delete h_timeFSLPC_tc_end[iF][iS][iL][iP][iC];
1398  }
1399  }
1400  }
1401  }
1402  }
1403  m_outFile->cd();
1404  m_outFile->Write();
1405  m_outFile->Close();
1406  B2INFO("File Write and Close. Done.");
1407 }
1408 
1410 {
1411  double tS = 0.0;
1412  int iSub = klmChannel.getSubdetector();
1413  int iF = klmChannel.getSection();
1414  int iS = klmChannel.getSector();
1415  int iL = klmChannel.getLayer();
1416  int iP = klmChannel.getPlane();
1417  int iC = klmChannel.getStrip();
1418  int totNStrips = EKLMElementNumbers::getMaximalStripNumber();
1419  if (iSub == KLMElementNumbers::c_BKLM)
1420  totNStrips = BKLMElementNumbers::getNStrips(iF, iS, iL, iP);
1421  if (iC == 1) {
1422  KLMChannelIndex kCIndex_upper(iSub, iF, iS, iL, iP, iC + 1);
1423  tS = tS_upperStrip(kCIndex_upper).second;
1424  } else if (iC == totNStrips) {
1425  KLMChannelIndex kCIndex_lower(iSub, iF, iS, iL, iP, iC - 1);
1426  tS = tS_lowerStrip(kCIndex_lower).second;
1427  } else {
1428  KLMChannelIndex kCIndex_upper(iSub, iF, iS, iL, iP, iC + 1);
1429  KLMChannelIndex kCIndex_lower(iSub, iF, iS, iL, iP, iC - 1);
1430  std::pair<int, double> tS_upper = tS_upperStrip(kCIndex_upper);
1431  std::pair<int, double> tS_lower = tS_lowerStrip(kCIndex_lower);
1432  unsigned int td_upper = tS_upper.first - iC;
1433  unsigned int td_lower = iC - tS_lower.first;
1434  unsigned int td = tS_upper.first - tS_lower.first;
1435  tS = (double(td_upper) * tS_lower.second + double(td_lower) * tS_upper.second) / double(td);
1436  }
1437  return tS;
1438 }
1439 
1440 std::pair<int, double> KLMTimeAlgorithm::tS_upperStrip(const KLMChannelIndex& klmChannel)
1441 {
1442  std::pair<int, double> tS;
1443  int cId = klmChannel.getKLMChannelNumber();
1444  int iSub = klmChannel.getSubdetector();
1445  int iF = klmChannel.getSection();
1446  int iS = klmChannel.getSector();
1447  int iL = klmChannel.getLayer();
1448  int iP = klmChannel.getPlane();
1449  int iC = klmChannel.getStrip();
1450  int totNStrips = EKLMElementNumbers::getMaximalStripNumber();
1451  if (iSub == KLMElementNumbers::c_BKLM)
1452  totNStrips = BKLMElementNumbers::getNStrips(iF, iS, iL, iP);
1453  if (m_timeShift.find(cId) != m_timeShift.end()) {
1454  tS.first = iC;
1455  tS.second = m_timeShift[cId];
1456  } else if (iC == totNStrips) {
1457  tS.first = iC;
1458  tS.second = 0.0;
1459  } else {
1460  KLMChannelIndex kCIndex(iSub, iF, iS, iL, iP, iC + 1);
1461  tS = tS_upperStrip(kCIndex);
1462  }
1463  return tS;
1464 }
1465 
1466 std::pair<int, double> KLMTimeAlgorithm::tS_lowerStrip(const KLMChannelIndex& klmChannel)
1467 {
1468  std::pair<int, double> tS;
1469  int cId = klmChannel.getKLMChannelNumber();
1470  int iSub = klmChannel.getSubdetector();
1471  int iF = klmChannel.getSection();
1472  int iS = klmChannel.getSector();
1473  int iL = klmChannel.getLayer();
1474  int iP = klmChannel.getPlane();
1475  int iC = klmChannel.getStrip();
1476  if (m_timeShift.find(cId) != m_timeShift.end()) {
1477  tS.first = iC;
1478  tS.second = m_timeShift[cId];
1479  } else if (iC == 1) {
1480  tS.first = iC;
1481  tS.second = 0.0;
1482  } else {
1483  KLMChannelIndex kCIndex(iSub, iF, iS, iL, iP, iC - 1);
1484  tS = tS_lowerStrip(kCIndex);
1485  }
1486  return tS;
1487 }
@ c_FirstRPCLayer
First RPC layer.
static int getNStrips(int section, int sector, int layer, int plane)
Get number of strips.
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Singleton class to cache database objects.
Definition: DBStore.h:32
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:52
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:92
static constexpr int getMaximalStripNumber()
Get maximal strip number.
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:33
double getStripLength(int strip) const
Get strip length.
Definition: GeometryData.h:71
double getMaximalStripLength() const
Get maximal strip length.
Definition: GeometryData.h:79
KLM channel index.
int getSubdetector() const
Get subdetector.
int getLayer() const
Get layer.
KLMChannelIndex begin()
First channel.
KLMChannelIndex & end()
Last channel.
int getSection() const
Get section.
int getPlane() const
Get plane.
int getStrip() const
Get strip.
int getSector() const
Get sector.
KLMChannelNumber getKLMChannelNumber() const
Get KLM channel number.
static const KLMElementNumbers & Instance()
Instantiation.
void channelNumberToElementNumbers(KLMChannelNumber channel, int *subdetector, int *section, int *sector, int *layer, int *plane, int *strip) const
Get element numbers by channel number.
TProfile * m_Profile2EKLMScintillatorPlane2
For EKLM scintillator plane2.
TH2F * m_HistTimeLengthEKLM[2][4][14][2][75]
Two-dimensional distributions of time versus propagation length.
TH1F * h_timeFSLPC_tc[2][8][15][2][54]
BKLM part, used for effective light speed estimation.
TH2F * h2c_timeF_scint_end[2]
EKLM part.
TH1F * h_timeFSLPC_tc_end[2][4][14][2][75]
EKLM part, used for effective light speed estimation.
TH1F * h_time_scint_tc_end
EKLM part.
void createHistograms()
Create histograms.
TGraphErrors * gre_time_channel_scint
BKLM Scintillator.
double m_LowerTimeBoundaryCalibratedScintilltorsEKLM
Lower time boundary for EKLM scintillators (calibrated data).
double m_LowerTimeBoundaryScintilltorsEKLM
Lower time boundary for EKLM scintillators.
TH1F * h_timeFSL[2][8][15]
BKLM part.
TH1F * hc_timeFSLPC_end[2][4][14][2][75]
EKLM part.
TH1F * hc_timeFSL_end[2][4][14]
EKLM part.
TH1F * h_timeFSLP_end[2][4][14][2]
EKLM part.
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 ecah 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 m_UpperTimeBoundaryCalibratedRPC
Upper time boundary for RPC (calibrated data).
KLMTimeConstants * m_timeConstants
DBObject of time cost on some parts of the detector.
void setupDatabase()
Setup the database.
double m_LowerTimeBoundaryCalibratedScintilltorsBKLM
Lower time boundary for BKLM scintillators (calibrated data).
double m_UpperTimeBoundaryCalibratedScintilltorsBKLM
Upper time boundary for BKLM scintillators (calibrated data).
TH1F * hc_timeFS_scint[2][8]
BKLM scintillator part.
std::map< KLMChannelNumber, double > m_time_channel
Time distribution central value of each channel.
CalibrationAlgorithm::EResult readCalibrationData()
Read calibration data.
TGraph * gr_timeShift_channel_scint_end
EKLM.
TH1F * hc_timeF_scint[2]
BKLM scintillator part.
TH1F * h_timeFS_scint[2][8]
BKLM scintillator part.
const KLMElementNumbers * m_ElementNumbers
Element numbers.
TH1F * hc_timeF_rpc[2]
BKLM RPC part.
TH2F * h2c_timeFS_end[2][4]
EKLM part.
TH1F * h_timeFSLPC_end[2][4][14][2][75]
EKLM part.
const EKLM::GeometryData * m_EKLMGeometry
EKLM geometry data.
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.
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.
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 avaiable channels with increasing strip number.
TH1F * hc_timeFSLP[2][8][15][2]
BKLM part.
void saveHist()
Save histograms to file.
const bklm::GeometryPar * m_BKLMGeometry
BKLM geometry data.
TH2F * h2c_timeFSLP[2][8][15][2]
BKLM part.
TF1 * fcn_const
Const function.
TProfile * m_Profile2RpcZ
For BKLM RPC z plane.
TH1F * h_timeFSLPC[2][8][15][2][54]
BKLM part.
TH2F * m_HistTimeLengthBKLM[2][8][15][2][54]
Two-dimensional distributions of time versus propagation length.
TH2F * h2_timeFSLP_end[2][4][14][2]
EKLM part.
void timeDistance2dFit(const std::vector< std::pair< KLMChannelNumber, unsigned int > > &channels, double &delay, double &delayError)
Two-dimensional fit for individual channels.
TH1I * h_calibrated
Calibration statistics for each channel.
TProfile * m_ProfileBKLMScintillatorZ
For BKLM scintillator z plane.
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.
TH2F * h2c_timeFSLP_end[2][4][14][2]
EKLM part.
TH2F * h2_timeF_scint_end[2]
EKLM part.
double m_UpperTimeBoundaryScintilltorsBKLM
Upper time boundary for BKLM scintillators.
TF1 * fcn_gaus
Gaussian function.
TH1F * h_timeF_rpc[2]
BKLM RPC part.
TProfile * m_ProfileBKLMScintillatorPhi
For BKLM scintillator phi plane.
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).
double m_UpperTimeBoundaryCalibratedScintilltorsEKLM
Upper time boundary for BKLM scintillators (calibrated data).
TH1F * h_timeF_scint_end[2]
EKLM part.
TProfile * m_ProfileRpcPhi
For BKLM RPC phi plane.
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.
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.
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).
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 avaiable channels with decreasing strip number.
TProfile * m_ProfileRpcZ
For BKLM RPC z plane.
double m_LowerTimeBoundaryScintilltorsBKLM
Lower time boundary for BKLM scintillators.
TH1F * hc_timeFSLPC[2][8][15][2][54]
BKLM part.
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).
std::map< KLMChannelNumber, double > m_etime_channel
Time distribution central value Error of each channel.
double m_UpperTimeBoundaryScintilltorsEKLM
Upper time boundary for BKLM scintillators.
int m_lower_limit_counts
Lower limit of hits collected for on single channel.
Class to store BKLM delay time coused by cable in the database.
void setTimeDelay(KLMChannelNumber channel, float timeDelay)
Set the time calibration constant.
Class to store KLM constants related to time.
@ c_BKLM
BKLM scintillator.
@ c_EKLM
EKLM scintillator.
void setDelay(float delay, int cType)
Set effective light speed of scintillators.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:110
bool construct(Args &&... params)
Construct an object of type T in this StoreObjPtr, using the provided constructor arguments.
Definition: StoreObjPtr.h:118
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
Definition: GeometryPar.cc:727
double getMaximalZStripLength() const
Get maximal Z strip length (for scintillators).
Definition: GeometryPar.h:296
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Definition: GeometryPar.cc:28
double getMaximalPhiStripLength() const
Get maximal phi strip length (for scintillators).
Definition: GeometryPar.h:290
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition: Module.h:76
Class to store variables with their name which were sent to the logging service.
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:26
void updateEvent()
Updates all intra-run dependent objects.
Definition: DBStore.cc:140
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:77
uint16_t KLMChannelNumber
Channel number.
Abstract base class for different kinds of events.