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