Belle II Software  release-05-02-19
eclHadronTimeCalibrationValidationCollectorModule.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 
10 #include <ecl/modules/eclHadronTimeCalibrationValidationCollector/eclHadronTimeCalibrationValidationCollectorModule.h>
11 #include <framework/dataobjects/EventMetaData.h>
12 #include <framework/gearbox/Const.h>
13 #include <ecl/dbobjects/ECLCrystalCalib.h>
14 #include <ecl/dataobjects/ECLDigit.h>
15 #include <ecl/dataobjects/ECLCalDigit.h>
16 #include <ecl/dataobjects/ECLTrig.h>
17 #include <mdst/dataobjects/ECLCluster.h>
18 #include <mdst/dataobjects/Track.h>
19 #include <mdst/dataobjects/HitPatternCDC.h>
20 #include <tracking/dataobjects/RecoTrack.h>
21 #include <ecl/digitization/EclConfiguration.h>
22 
23 #include <analysis/utility/PCmsLabTransform.h>
24 #include <analysis/ClusterUtility/ClusterUtils.h>
25 #include <boost/optional.hpp>
26 
27 #include <TH2F.h>
28 #include <TTree.h>
29 #include <TFile.h>
30 
31 using namespace Belle2;
32 using namespace ECL;
33 using namespace std;
34 
35 //-----------------------------------------------------------------
36 // Register the Module
37 //-----------------------------------------------------------------
38 REG_MODULE(eclHadronTimeCalibrationValidationCollector)
39 
40 //-----------------------------------------------------------------
41 // Implementation
42 //-----------------------------------------------------------------
43 
46  m_dbg_tree_photonClusters(0),
47  m_dbg_tree_event(0),
48  m_tree_evt_num(0)//,
49  //m_GammaGammaECalib("ECLCrystalEnergyGammaGamma")
50 {
51  setDescription("This module validates the ECL cluster times");
52 
53  addParam("timeAbsMax", m_timeAbsMax, // (Time in ns)
54  "Events with fabs(getTimeFit) > m_timeAbsMax "
55  "are excluded", (short)80);
56 
57  addParam("saveTree", m_saveTree,
58  "If true, TTree 'tree' with more detailed event info is saved in "
59  "the output file specified by HistoManager",
60  false);
61 
62  addParam("looseTrkZ0", m_looseTrkZ0, "max Z0 for loose tracks (cm)", 10.);
63  addParam("tightTrkZ0", m_tightTrkZ0, "max Z0 for tight tracks (cm)", 2.);
64  addParam("looseTrkD0", m_looseTrkD0, "max D0 for loose tracks (cm)", 2.);
65  addParam("tightTrkD0", m_tightTrkD0, "max D0 for tight tracks (cm)", 0.5); // beam pipe radius = 1cm in 2019
66 
67 
68  // specify this flag if you need parallel processing
69  setPropertyFlags(c_ParallelProcessingCertified);
70 }
71 
73 {
74 }
75 
77 {
78 
79  //=== Prepare TTree for debug output
80  if (m_saveTree) {
81  // Per photon cluster
82  m_dbg_tree_photonClusters = new TTree("tree_photonClusters",
83  "Validating crystal and crate time calibrations using photon clusters in events with lots of tracks and clusters") ;
84  m_dbg_tree_photonClusters->Branch("EventNum" , &m_tree_evt_num) ->SetTitle("Event number") ;
85  m_dbg_tree_photonClusters->Branch("cluster_time" , &m_tree_time) ->SetTitle("Cluster time t (calibrated), ns") ;
86  m_dbg_tree_photonClusters->Branch("clust_E" , &m_E_photon_clust) ->SetTitle("Photon type cluster energy, GeV") ;
87  m_dbg_tree_photonClusters->Branch("Ntracks" , &m_NtightTracks) ->SetTitle("Number of tracks") ;
88  m_dbg_tree_photonClusters->Branch("NphotonClusters" , &m_NphotonClusters) ->SetTitle("Number of photons") ;
89  m_dbg_tree_photonClusters->Branch("NGoodClusters" , &m_NGoodClusters) ->SetTitle("Number of good ECL clusters") ;
90  m_dbg_tree_photonClusters->Branch("t0" , &m_tree_t0) ->SetTitle("T0, ns") ;
91  m_dbg_tree_photonClusters->Branch("t0_unc" , &m_tree_t0_unc) ->SetTitle("T0 uncertainty, ns") ;
92  m_dbg_tree_photonClusters->Branch("runNum" , &m_tree_run) ->SetTitle("Run number") ;
93  m_dbg_tree_photonClusters->Branch("CrystalCellID" , &m_tree_cid) ->SetTitle("Cell ID, 1..8736") ;
94  m_dbg_tree_photonClusters->Branch("dt99" , &m_tree_dt99) ->SetTitle("Cluster dt99, ns") ;
95  m_dbg_tree_photonClusters->SetAutoSave(10) ;
96 
97  // Per event
98  m_dbg_tree_event = new TTree("tree_event",
99  "Validating crystal and crate time calibrations using photon clusters in events with lots of tracks and clusters") ;
100  m_dbg_tree_event->Branch("EventNum" , &m_tree_evt_num) ->SetTitle("Event number") ;
101  m_dbg_tree_event->Branch("t0" , &m_tree_t0) ->SetTitle("T0, ns") ;
102  m_dbg_tree_event->Branch("t0_unc" , &m_tree_t0_unc) ->SetTitle("T0 uncertainty, ns") ;
103  m_dbg_tree_event->Branch("runNum" , &m_tree_run) ->SetTitle("Run number") ;
104  m_dbg_tree_event->Branch("Ntracks" , &m_NtightTracks) ->SetTitle("Number of tracks") ;
105  m_dbg_tree_event->Branch("NphotonClusters" , &m_NphotonClusters) ->SetTitle("Number of photons") ;
106  m_dbg_tree_event->Branch("NGoodClusters" , &m_NGoodClusters) ->SetTitle("Number of good ECL clusters") ;
107  m_dbg_tree_event->Branch("E0" , &m_tree_E0) ->SetTitle("Highest E cluster E") ;
108  m_dbg_tree_event->Branch("time_E0" , &m_tree_time_fromE0) ->SetTitle("Cluster time of highest E cluster") ;
109  m_dbg_tree_event->SetAutoSave(10) ;
110  }
111 }
112 
114 {
115  //=== MetaData
116  StoreObjPtr<EventMetaData> evtMetaData ;
117  B2INFO("eclHadronTimeCalibrationValidationCollector: Experiment = " << evtMetaData->getExperiment() <<
118  " run = " << evtMetaData->getRun()) ;
119 
120  //=== Create histograms and register them in the data store
121 
122  // Define the bin size, which is equivalent to the
123  double binSize = 2000.0 / pow(2, 12);
124  double halfBinSize = binSize / 2.0;
125 
126  /* Determine the number of bins required to go from the edge of the bin centred
127  on zero to a value just larger than the negative cut off */
128  double nBinsNeg = floor((m_timeAbsMax - halfBinSize) / binSize);
129  double min_t = -nBinsNeg * binSize - halfBinSize; // lower edge value of left most bin
130  int nbins = nBinsNeg * 2 + 1; // number of negative time bins + t=0 bin + number of positive time bins
131  double max_t = min_t + nbins * binSize; // upper edge value of right most bin
132 
133  /* Variable bin width information for the time information vs energy since
134  the time width should vary as a function of 1/E */
135  const Int_t N_E_BIN_EDGES = 64;
136  const Int_t N_E_BINS = N_E_BIN_EDGES - 1;
137  Double_t energyBinEdges[N_E_BIN_EDGES] = {0, 0.05, 0.051, 0.052, 0.053, 0.054, 0.055, 0.056, 0.057, 0.058, 0.059, 0.06, 0.062, 0.064, 0.066, 0.068, 0.07, 0.075, 0.08, 0.085, 0.09, 0.095, 0.1, 0.11, 0.12, 0.13, 0.14, 0.15, 0.16, 0.17, 0.18, 0.19, 0.2, 0.25, 0.3, 0.35, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1, 1.2, 1.4, 1.6, 1.8, 2, 2.25, 2.5, 2.8, 3.2, 3.6, 4, 4.4, 4.8, 5.2, 5.6, 6, 6.4, 6.8, 7.2, 7.6, 8};
138 
139 
140  auto cutflow = new TH1F("cutflow", " ;Cut label number ;Number of events passing cut", 10, 0, 10) ;
141  registerObject<TH1F>("cutflow", cutflow) ;
142 
143  auto clusterTime = new TH1F("clusterTime", ";Photon ECL cluster time [ns]; number of photon ECL clusters", nbins, min_t, max_t) ;
144  registerObject<TH1F>("clusterTime", clusterTime) ;
145 
146  auto clusterTime_cid = new TH2F("clusterTime_cid",
147  ";crystal Cell ID ;Photon ECL cluster time [ns]", 8736, 1, 8736 + 1, nbins, min_t, max_t) ;
148  registerObject<TH2F>("clusterTime_cid", clusterTime_cid) ;
149 
150  auto clusterTime_run = new TH2F("clusterTime_run",
151  ";Run number ;Photon ECL cluster time [ns]", 7000, 0, 7000, nbins, min_t, max_t) ;
152  registerObject<TH2F>("clusterTime_run", clusterTime_run) ;
153 
154 
155  auto clusterTimeClusterE = new TH2F("clusterTimeClusterE",
156  ";Photon cluster energy [GeV];Photon cluster time [ns]", N_E_BINS, energyBinEdges, nbins, min_t, max_t) ;
157  registerObject<TH2F>("clusterTimeClusterE", clusterTimeClusterE) ;
158 
159 
160  auto dt99_clusterE = new TH2F("dt99_clusterE",
161  ";Photon cluster energy [GeV];dt99 [ns]", N_E_BINS, energyBinEdges, nbins, 0, max_t) ;
162  registerObject<TH2F>("dt99_clusterE", dt99_clusterE) ;
163 
164 
165  auto eventT0 = new TH1F("eventT0", ";event t0 [ns]; number of events", nbins, min_t, max_t) ;
166  registerObject<TH1F>("eventT0", eventT0) ;
167 
168 
169  auto clusterTimeE0E1diff = new TH1F("clusterTimeE0E1diff",
170  ";ECL cluster time of max E photon - ECL cluster time of 2nd max E photon [ns]; number of photon ECL cluster time differences",
171  nbins, min_t, max_t) ;
172  registerObject<TH1F>("clusterTimeE0E1diff", clusterTimeE0E1diff) ;
173 
174 
175  //=== Required data objects
176  tracks.isRequired() ;
177  m_eclClusterArray.isRequired() ;
178  m_eclCalDigitArray.isRequired() ;
179 
180 }
181 
183 {
184  int cutIndexPassed = 0;
185  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
186  B2DEBUG(22, "Cutflow: no cuts: index = " << cutIndexPassed);
187 
188 
189  /* Use ECLChannelMapper to get other detector indices for the crystals */
190  /* For conversion from CellID to crate, shaper, and channel ids. */
191 
192  // Use smart pointer to avoid memory leak when the ECLChannelMapper object needs destroying at the end of the event.
193  shared_ptr< ECL::ECLChannelMapper > crystalMapper(new ECL::ECLChannelMapper());
194  crystalMapper->initFromDB();
195 
196 
197  // Storage crystal energies
198  m_EperCrys.resize(8736);
199  for (auto& eclCalDigit : m_eclCalDigitArray) {
200  int tempCrysID = eclCalDigit.getCellId() - 1;
201  m_EperCrys[tempCrysID] = eclCalDigit.getEnergy();
202  }
203 
204  // Getting the event t0 using the full event t0 rather than from the CDC specifically
205 
206  double evt_t0 = -1000 ;
207  double evt_t0_unc = -1000 ;
208 
209  // Determine if there is an event t0 to use and then extract the information about it
210  if (m_eventT0.isOptional()) {
211  if (!m_eventT0.isValid()) {
212  return;
213  }
214  if (!m_eventT0->hasEventT0()) {
215  return;
216  } else {
217  // Overall event t0 (combination of multiple event t0s from different detectors)
218  evt_t0 = m_eventT0->getEventT0() ;
219  evt_t0_unc = m_eventT0->getEventT0Uncertainty() ;
220  }
221  B2DEBUG(26, "Found event t0") ;
222  }
223 
224  //---------------------------------------------------------------------
225  //..Track properties. Use pion (211) mass hypothesis,
226  // which is the only particle hypothesis currently available???
227  int nTrkAll = tracks.getEntries() ;
228 
229  int nTrkLoose = 0 ;
230  int nTrkTight = 0 ;
233  /* Loop over all the tracks to define the tight and loose selection tracks.
234  We will select events with only a few tight tracks and no additional loose tracks.
235  Tight tracks are a subset of looses tracks. */
236  for (int iTrk = 0 ; iTrk < nTrkAll ; iTrk++) {
237 
238  // Get track biasing towards the particle being a pion
239  const TrackFitResult* tempTrackFit = tracks[iTrk]->getTrackFitResultWithClosestMass(Const::pion) ;
240  if (not tempTrackFit) {continue ;}
241 
242  // Collect track info to be used for categorizing
243  //short charge = tempTrackFit->getChargeSign() ;
244  double z0 = tempTrackFit->getZ0() ;
245  double d0 = tempTrackFit->getD0() ;
246  int nCDChits = tempTrackFit->getHitPatternCDC().getNHits() ;
247  //double pt = tempTrackFit->getTransverseMomentum() ;
248  //double p = tempTrackFit->getMomentum().Mag() ;
249 
250  /* Test if loose track */
251 
252  // d0 and z0 cuts
253  if (fabs(d0) > m_looseTrkD0) {
254  continue;
255  }
256  if (fabs(z0) > m_looseTrkZ0) {
257  continue;
258  }
259  // Number of hits in the CDC
260  if (nCDChits < 1) {
261  continue;
262  }
263  nTrkLoose++;
264 
265 
266 
267  /* Test if the loose track is also a tight track */
268 
269  // Number of hits in the CDC
270  if (nCDChits < 20) {
271  continue;
272  }
273  // d0 and z0 cuts
274  if (fabs(d0) > m_tightTrkD0) {
275  continue;
276  }
277  if (fabs(z0) > m_tightTrkZ0) {
278  continue;
279  }
280  nTrkTight++;
281 
282  }
283  // After that last section the numbers of loose and tight tracks are known
284  B2DEBUG(26, "Found loose and tight tracks") ;
285 
286 
287  int numGoodTightTracks_minCut = 4 ;
288  if (nTrkTight < numGoodTightTracks_minCut) {
289  return ;
290  }
291  // There are at least X tight tracks
292  cutIndexPassed++ ;
293  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed) ;
294  B2DEBUG(22, "Cutflow: At least " << numGoodTightTracks_minCut << " tight tracks: index = " << cutIndexPassed) ;
295 
296 
297  int numGoodLooseTracks_minCut = numGoodTightTracks_minCut ;
298  if (nTrkLoose < numGoodLooseTracks_minCut) {
299  return ;
300  }
301  // There are more loose tracks than tight tracks then veto the event. If there are fewer loose tracks than tight tracks then veto the event, although this should be impossible
302  cutIndexPassed++ ;
303  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed) ;
304  B2DEBUG(22, "Cutflow: No additional loose tracks: index = " << cutIndexPassed) ;
305 
306 
307  //------------------------------------------------------------------------
308  // Find the good ECL clusters
309  double clusterE_minCut = 0.1 ; // GeV
310  int nclust = m_eclClusterArray.getEntries();
311  int nGoodClusts = 0 ;
312  vector<int> goodClusterIdxs ;
313  for (int ic = 0; ic < nclust; ic++) {
314  double eClust = m_eclClusterArray[ic]->getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons) ;
315  if (eClust > clusterE_minCut) {
316  goodClusterIdxs.push_back(ic) ;
317  nGoodClusts++ ;
318  }
319  }
320 
321 
322  // Cut on the minimum number of good clusters
323  int numGoodEMclusters_minCut = 5 ;
324  if (nGoodClusts < numGoodEMclusters_minCut) {
325  return ;
326  }
327  // There are at least 5 good EM clusters (photon = basically all clusters)
328  cutIndexPassed++ ;
329  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed) ;
330  B2DEBUG(22, "Cutflow: At least " << numGoodEMclusters_minCut << " ECL clusters: index = " << cutIndexPassed) ;
331 
332 
333  //------------------------------------------------------------------------
334  // Find the good photons first before doing anything with them
335 
336  double photonE_minCut = 0.05 ; // GeV
337  double zernikeMVA_minCut = 0.2 ;
338  int nPhotons = 0 ;
339 
340  vector<int> goodPhotonClusterIdxs ;
341  for (int ic = 0; ic < nclust; ic++) {
342  if (m_eclClusterArray[ic]->hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons)) {
343  double eClust = m_eclClusterArray[ic]->getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons);
344  double photonTime = m_eclClusterArray[ic]->getTime();
345  double zernikeMVA = m_eclClusterArray[ic]->getZernikeMVA();
346  bool badPhotonTime = m_eclClusterArray[ic]->hasFailedFitTime();
347  bool badPhotonTimeResolution = m_eclClusterArray[ic]->hasFailedTimeResolution();
348  bool hasTrack = m_eclClusterArray[ic]->isTrack();
349  if ((eClust > photonE_minCut) &&
350  (fabs(photonTime) < m_timeAbsMax) &&
351  (!badPhotonTime) &&
352  (!badPhotonTimeResolution) &&
353  (zernikeMVA > zernikeMVA_minCut) &&
354  (!hasTrack)) {
355  goodPhotonClusterIdxs.push_back(ic) ;
356  nPhotons++;
357  }
358  }
359  }
360 
361 
362  // Cut on the minimum number of good photon clusters
363  int numGoodPhotonclusters_minCut = 1 ;
364  if (nPhotons < numGoodPhotonclusters_minCut) {
365  return ;
366  }
367  // There is at least one good photon in the event
368  cutIndexPassed++ ;
369  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed) ;
370  B2DEBUG(22, "Cutflow: At least " << numGoodPhotonclusters_minCut << " good photon: index = " << cutIndexPassed) ;
371 
372 
373 
374  //------------------------------------------------------------------------
375  /* Extract the times of the good clusters and
376  save the maximum energy crystal information (cid) */
377  vector<double> goodClustTimes ;
378  vector<double> goodClust_dt99 ;
379  vector<double> goodClusters_crysE ;
380  vector<double> goodClustE ;
381  vector<int> goodClustMaxEcrys_cid ;
382  for (long unsigned int i = 0; i < goodPhotonClusterIdxs.size(); i++) {
383  int ic = goodPhotonClusterIdxs[i] ;
384 
385  if (m_eclClusterArray[ic]->hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons)) {
386  double eClust = m_eclClusterArray[ic]->getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons);
387  short cid = m_eclClusterArray[ic]->getMaxECellId() ;
388 
389  goodClustMaxEcrys_cid.push_back(cid) ;
390  goodClustTimes.push_back(m_eclClusterArray[ic]->getTime()) ;
391  goodClust_dt99.push_back(m_eclClusterArray[ic]->getDeltaTime99()) ;
392  goodClusters_crysE.push_back(m_EperCrys[cid - 1]) ;
393  goodClustE.push_back(eClust);
394  }
395  }
396 
397 
398  // Define a pair (energy,time) so that we can quickly and easily sort the cluster information
399  // based on the energy of the clusters
400  vector< pair<double, double> > pair_energy_time ;
401  for (long unsigned int ic = 0; ic < goodClusters_crysE.size(); ic++) {
402  pair_energy_time.push_back(make_pair(goodClusters_crysE[ic], goodClustTimes[ic])) ;
403  }
404 
405  // sorts pairs in decreasing order of their first value (energy)
406  // i.e. highest energy first
407  sort(pair_energy_time.begin(), pair_energy_time.end(), greater<>()) ;
408 
409 
410 
411  B2DEBUG(22, "Event passed all cuts");
412 
413 
414  // Fill the histogram for the event level variables
415  getObjectPtr<TH1F>("eventT0")->Fill(evt_t0) ;
416 
417  bool isCDCt0 = (static_cast<EventT0::EventT0Component>(*m_eventT0->getEventT0Component())).detectorSet.contains(Const::CDC);
418  bool isECLt0 = (static_cast<EventT0::EventT0Component>(*m_eventT0->getEventT0Component())).detectorSet.contains(Const::ECL);
419  string t0Detector = "UNKNOWN... WHY?";
420  if (isCDCt0) {
421  t0Detector = "CDC" ;
422  } else if (isECLt0) {
423  t0Detector = "ECL" ;
424  }
425 
426  B2DEBUG(26, "t0 = " << evt_t0 << " ns. t0 is from CDC?=" << isCDCt0 << ", t0 is from ECL?=" << isECLt0 << " t0 from " <<
427  t0Detector);
428 
429 
430 
431  //=== For each good photon cluster in the processed event and fill histogram.
432 
433  StoreObjPtr<EventMetaData> evtMetaData ;
434  for (long unsigned int i = 0 ; i < goodPhotonClusterIdxs.size() ; i++) {
435  getObjectPtr<TH1F>("clusterTime")->Fill(goodClustTimes[i]) ;
436  getObjectPtr<TH2F>("clusterTime_cid")->Fill(goodClustMaxEcrys_cid[i] + 0.001, goodClustTimes[i] , 1) ;
437  getObjectPtr<TH2F>("clusterTime_run")->Fill(evtMetaData->getRun() + 0.001, goodClustTimes[i] , 1) ;
438  getObjectPtr<TH2F>("clusterTimeClusterE")->Fill(goodClustE[i], goodClustTimes[i], 1) ;
439  getObjectPtr<TH2F>("dt99_clusterE")->Fill(goodClustE[i], goodClust_dt99[i], 1) ;
440 
441  //== Save debug TTree with detailed information if necessary.
442  if (m_saveTree) {
443 
444  m_tree_time = goodClustTimes[i] ;
445  m_E_photon_clust = goodClusters_crysE[i] ;
446  m_tree_t0 = evt_t0 ;
447  m_tree_t0_unc = evt_t0_unc ;
448  m_NtightTracks = nTrkTight ;
449  m_NphotonClusters = nPhotons ;
450  m_NGoodClusters = nGoodClusts ;
451  m_tree_evt_num = evtMetaData->getEvent() ;
452  m_tree_run = evtMetaData->getRun() ;
453  m_tree_cid = goodClustMaxEcrys_cid[i] ;
454  m_tree_dt99 = goodClust_dt99[i] ;
455 
456  m_dbg_tree_photonClusters->Fill() ;
457 
458  }
459  }
460  B2DEBUG(26, "Filled cluster tree") ;
461 
462  //=== Fill histogram for cluster time difference of the two max E photons
463  if (pair_energy_time.size() >= 2) {
464  getObjectPtr<TH1F>("clusterTimeE0E1diff")->Fill(pair_energy_time[0].second - pair_energy_time[1].second) ;
465  }
466 
467 
468 
469  if (m_saveTree) {
470  m_tree_t0 = evt_t0 ;
471  m_tree_t0_unc = evt_t0_unc ;
472  m_tree_evt_num = evtMetaData->getEvent() ;
473  m_tree_run = evtMetaData->getRun() ;
474  m_NtightTracks = nTrkTight ;
475  m_NphotonClusters = nPhotons ;
476  m_NGoodClusters = nGoodClusts ;
477 
478  m_tree_E0 = pair_energy_time[0].first ;
479  m_tree_time_fromE0 = pair_energy_time[0].second ;
480  m_dbg_tree_event->Fill() ;
481  }
482 
483  B2DEBUG(26, "Filled event tree") ;
484 
485 }
Belle2::CalibrationCollectorModule
Calibration collector module base class.
Definition: CalibrationCollectorModule.h:44
Belle2::eclHadronTimeCalibrationValidationCollectorModule
This module generates 'TimevsCrys' histogram to later (in eclBhabhaTAlgorithm) find time offset from ...
Definition: eclHadronTimeCalibrationValidationCollectorModule.h:44
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::ECLCluster::EHypothesisBit::c_nPhotons
@ c_nPhotons
CR is split into n photons (N1)
Belle2::eclHadronTimeCalibrationValidationCollectorModule::prepare
void prepare() override
Define histograms and read payloads from DB.
Definition: eclHadronTimeCalibrationValidationCollectorModule.cc:113
Belle2::eclHadronTimeCalibrationValidationCollectorModule::collect
void collect() override
Select events and crystals and accumulate histograms.
Definition: eclHadronTimeCalibrationValidationCollectorModule.cc:182
Belle2::TrackFitResult
Values of the result of a track fit with a given particle hypothesis.
Definition: TrackFitResult.h:59
Belle2::EventT0::EventT0Component
Structure for storing the extracted event t0s together with its detector and its uncertainty.
Definition: EventT0.h:44
Belle2::Const::pion
static const ChargedStable pion
charged pion particle
Definition: Const.h:535
Belle2::ECL::ECLChannelMapper
This class provides access to ECL channel map that is either a) Loaded from the database (see ecl/dbo...
Definition: ECLChannelMapper.h:36
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::eclHadronTimeCalibrationValidationCollectorModule::inDefineHisto
void inDefineHisto() override
Replacement for defineHisto() in CalibrationCollector modules.
Definition: eclHadronTimeCalibrationValidationCollectorModule.cc:76
Belle2::eclHadronTimeCalibrationValidationCollectorModule::~eclHadronTimeCalibrationValidationCollectorModule
virtual ~eclHadronTimeCalibrationValidationCollectorModule()
Module destructor.
Definition: eclHadronTimeCalibrationValidationCollectorModule.cc:72
Belle2::ECL::ECLChannelMapper::initFromDB
bool initFromDB()
Initialize channel mapper from the conditions database.
Definition: ECLChannelMapper.cc:105