Belle II Software  release-08-00-10
eclBhabhaTimeCalibrationValidationCollectorModule.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  * Description: *
11  * Tool to perform validations using bhabha events. *
12  * The selections here are very similar to the selections used to *
13  * perform the bhabha based timing calibrations; however, *
14  * there are some small differences so that the cuts are slightly *
15  * different. Further differences could appear if the two codes *
16  * are not updated simultaneously if/when there are changes. *
17  * The intention of this code is just to make sure the calibrations *
18  * results are self consistent and to validate the calibrations *
19  * as a whole for a collection of bhabha events rather than to *
20  * validate the cluster times of each individual cluster from the *
21  * original calibration code. Hopefully the selection is looser *
22  * so that fewer events can be used to validate than are required *
23  * to perform the calibration itself. *
24  **************************************************************************/
25 
26 /* Own header. */
27 #include <ecl/modules/eclBhabhaTimeCalibrationValidationCollector/eclBhabhaTimeCalibrationValidationCollectorModule.h>
28 
29 /* ECL headers. */
30 #include <ecl/dataobjects/ECLCalDigit.h>
31 #include <ecl/dataobjects/ECLDigit.h>
32 #include <ecl/dataobjects/ECLElementNumbers.h>
33 #include <ecl/dataobjects/ECLTrig.h>
34 #include <ecl/dbobjects/ECLCrystalCalib.h>
35 #include <ecl/digitization/EclConfiguration.h>
36 
37 /* Basf2 headers. */
38 #include <analysis/ClusterUtility/ClusterUtils.h>
39 #include <analysis/utility/PCmsLabTransform.h>
40 #include <framework/gearbox/Const.h>
41 #include <mdst/dataobjects/ECLCluster.h>
42 #include <mdst/dataobjects/HitPatternCDC.h>
43 #include <mdst/dataobjects/Track.h>
44 
45 /* ROOT headers. */
46 #include <TH2F.h>
47 #include <TTree.h>
48 #include <TFile.h>
49 
50 using namespace Belle2;
51 using namespace ECL;
52 using namespace std;
53 
54 //-----------------------------------------------------------------
55 // Register the Module
56 //-----------------------------------------------------------------
57 REG_MODULE(eclBhabhaTimeCalibrationValidationCollector);
58 
59 //-----------------------------------------------------------------
60 // Implementation
61 //-----------------------------------------------------------------
62 
65  m_dbg_tree_electronClusters(0),
66  m_dbg_tree_event(0),
67  m_dbg_tree_run(0),
68  m_CrateTimeDB("ECLCrateTimeOffset"),
69  m_channelMapDB("ECLChannelMap")//,
70 {
71  setDescription("This module validates the ECL cluster times");
72 
73  addParam("timeAbsMax", m_timeAbsMax, // (Time in ns)
74  "Events with fabs(getTimeFit) > m_timeAbsMax "
75  "are excluded", (short)80);
76 
77  addParam("saveTree", m_saveTree,
78  "If true, TTree 'tree' with more detailed event info is saved in "
79  "the output file specified by HistoManager",
80  false);
81 
82  addParam("looseTrkZ0", m_looseTrkZ0, "max Z0 for loose tracks (cm)", 10.);
83  addParam("tightTrkZ0", m_tightTrkZ0, "max Z0 for tight tracks (cm)", 2.);
84  addParam("looseTrkD0", m_looseTrkD0, "max D0 for loose tracks (cm)", 2.);
85  addParam("tightTrkD0", m_tightTrkD0, "max D0 for tight tracks (cm)", 0.5); // beam pipe radius = 1cm in 2019
86  addParam("skipTrgSel", skipTrgSel, "boolean to skip the trigger skim selection", false);
87 
88 
89  // specify this flag if you need parallel processing
91 }
92 
94 {
95 }
96 
98 {
99 
100 }
101 
103 {
104  //=== Prepare TTree for debug output
105  if (m_saveTree) {
106  // Per electron cluster
107  m_dbg_tree_electronClusters = new TTree("tree_electronClusters",
108  "Validating crystal and crate time calibrations using electron clusters in events with lots of tracks and clusters") ;
109  m_dbg_tree_electronClusters->Branch("EventNum", &m_tree_evt_num) ->SetTitle("Event number") ;
110  m_dbg_tree_electronClusters->Branch("cluster_time", &m_tree_time) ->SetTitle("Cluster time t (calibrated), ns") ;
111  m_dbg_tree_electronClusters->Branch("clust_E", &m_E_electron_clust) ->SetTitle("Electron type cluster energy, GeV") ;
112  m_dbg_tree_electronClusters->Branch("t0", &m_tree_t0) ->SetTitle("T0, ns") ;
113  m_dbg_tree_electronClusters->Branch("t0_unc", &m_tree_t0_unc) ->SetTitle("T0 uncertainty, ns") ;
114  m_dbg_tree_electronClusters->Branch("runNum", &m_tree_run) ->SetTitle("Run number") ;
115  m_dbg_tree_electronClusters->Branch("CrystalCellID", &m_tree_cid) ->SetTitle("Cell ID, 1..8736") ;
116  m_dbg_tree_electronClusters->Branch("dt99", &m_tree_dt99) ->SetTitle("Cluster dt99, ns") ;
117  m_dbg_tree_electronClusters->SetAutoSave(10) ;
118 
119  // Per event
120  m_dbg_tree_event = new TTree("tree_event",
121  "Validating crystal and crate time calibrations using electron clusters in events with lots of tracks and clusters") ;
122  m_dbg_tree_event->Branch("EventNum", &m_tree_evt_num) ->SetTitle("Event number") ;
123  m_dbg_tree_event->Branch("t0", &m_tree_t0) ->SetTitle("T0, ns") ;
124  m_dbg_tree_event->Branch("t0_unc", &m_tree_t0_unc) ->SetTitle("T0 uncertainty, ns") ;
125  m_dbg_tree_event->Branch("runNum", &m_tree_run) ->SetTitle("Run number") ;
126  m_dbg_tree_event->Branch("E0", &m_tree_E0) ->SetTitle("Highest E cluster E") ;
127  m_dbg_tree_event->Branch("E1", &m_tree_E1) ->SetTitle("2nd highest E cluster E") ;
128  m_dbg_tree_event->Branch("time_E0", &m_tree_time_fromE0) ->SetTitle("Cluster time of highest E cluster") ;
129  m_dbg_tree_event->Branch("time_E1", &m_tree_time_fromE1) ->SetTitle("Cluster time of 2nd highest E cluster") ;
130  m_dbg_tree_event->SetAutoSave(10) ;
131 
132 
133  // Per run
134  m_dbg_tree_run = new TTree("tree_run", "Storing crate time constants") ;
135  m_dbg_tree_run->Branch("runNum", &m_tree_run) ->SetTitle("Run number") ;
136  m_dbg_tree_run->Branch("crateid", &m_tree_crateid) ->SetTitle("Crate ID") ;
137  m_dbg_tree_run->Branch("tcrate", &m_tree_tcrate) ->SetTitle("Crate time") ;
138  m_dbg_tree_run->Branch("tcrate_unc", &m_tree_tcrate_unc) ->SetTitle("Crate time uncertainty") ;
139  m_dbg_tree_run->SetAutoSave(10) ;
140 
141  }
142 
143 
144  //=== MetaData
145  B2INFO("eclBhabhaTimeCalibrationValidationCollector: Experiment = " << m_EventMetaData->getExperiment() <<
146  " run = " << m_EventMetaData->getRun()) ;
147 
148  //=== Create histograms and register them in the data store
149 
150  // Define the bin size, which is equivalent to the
151  double binSize = 2000.0 / pow(2, 12);
152  double halfBinSize = binSize / 2.0;
153 
154  /* Determine the number of bins required to go from the edge of the bin centred
155  on zero to a value just larger than the negative cut off */
156  double nBinsNeg = floor((m_timeAbsMax - halfBinSize) / binSize);
157  double min_t = -nBinsNeg * binSize - halfBinSize; // lower edge value of left most bin
158  int nbins = nBinsNeg * 2 + 1; // number of negative time bins + t=0 bin + number of positive time bins
159  double max_t = min_t + nbins * binSize; // upper edge value of right most bin
160 
161  /* Variable bin width information for the time information vs energy since
162  the time width should vary as a function of 1/E */
163  const Int_t N_E_BIN_EDGES = 64;
164  const Int_t N_E_BINS = N_E_BIN_EDGES - 1;
165  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};
166 
167 
168  auto cutflow = new TH1F("cutflow", " ;Cut label number ;Number of events passing cut", 10, 0, 10) ;
169  registerObject<TH1F>("cutflow", cutflow) ;
170 
171  auto clusterTime = new TH1F("clusterTime", ";Electron ECL cluster time [ns]; number of ECL clusters", nbins, min_t, max_t) ;
172  registerObject<TH1F>("clusterTime", clusterTime) ;
173 
174  auto clusterTime_cid = new TH2F("clusterTime_cid",
175  ";crystal Cell ID ;Electron ECL cluster time [ns]", ECLElementNumbers::c_NCrystals, 1, ECLElementNumbers::c_NCrystals + 1, nbins,
176  min_t, max_t) ;
177  registerObject<TH2F>("clusterTime_cid", clusterTime_cid) ;
178 
179  auto clusterTime_run = new TH2F("clusterTime_run",
180  ";Run number ;Electron ECL cluster time [ns]", 7000, 0, 7000, nbins, min_t, max_t) ;
181  registerObject<TH2F>("clusterTime_run", clusterTime_run) ;
182 
183 
184  auto clusterTimeClusterE = new TH2F("clusterTimeClusterE",
185  ";Electron cluster energy [GeV];Electron cluster time [ns]", N_E_BINS, energyBinEdges, nbins, min_t, max_t) ;
186  registerObject<TH2F>("clusterTimeClusterE", clusterTimeClusterE) ;
187 
188  auto dt99_clusterE = new TH2F("dt99_clusterE",
189  ";Electron cluster energy [GeV];dt99 [ns]", N_E_BINS, energyBinEdges, nbins, 0, max_t) ;
190  registerObject<TH2F>("dt99_clusterE", dt99_clusterE) ;
191 
192 
193  auto eventT0 = new TH1F("eventT0", ";event t0 [ns]; number of events", nbins, min_t, max_t) ;
194  registerObject<TH1F>("eventT0", eventT0) ;
195 
196  auto clusterTimeE0E1diff = new TH1F("clusterTimeE0E1diff",
197  ";ECL cluster time of max E electron - ECL cluster time of 2nd max E electron [ns]; number of electron ECL cluster time differences",
198  nbins, min_t, max_t) ;
199  registerObject<TH1F>("clusterTimeE0E1diff", clusterTimeE0E1diff) ;
200 
201 
202  //=== Required data objects
203  tracks.isRequired() ;
204  m_eclClusterArray.isRequired() ;
205  m_eclCalDigitArray.isRequired() ;
206 
207 
208  B2INFO("skipTrgSel = " << skipTrgSel);
209 
210 }
211 
213 {
214  int cutIndexPassed = 0;
215  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
216  B2DEBUG(22, "Cutflow: no cuts: index = " << cutIndexPassed);
217 
218 
219  // --- Check the trigger skim is the type that has two tracks
220 
221  /* If we skip the trigger skim selection then still fill the cutflow histogram
222  just so that the positions don't change. */
223  if (!skipTrgSel) {
224  if (!m_TrgResult.isValid()) {
225  B2WARNING("SoftwareTriggerResult required to select bhabha event is not found");
226  return;
227  }
228 
229  /* Release05: bhabha_all is grand skim = bhabha+bhabhaecl+radee. We only want
230  to look at the 2 track bhabha events. */
231  const std::map<std::string, int>& fresults = m_TrgResult->getResults();
232  if (fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end()) {
233  B2WARNING("Can't find required bhabha trigger identifier");
234  return;
235  }
236 
237  const bool eBhabha = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
239  B2DEBUG(22, "eBhabha (trigger passed) = " << eBhabha);
240 
241  if (!eBhabha) {
242  return;
243  }
244  }
245 
246  /* Fill the histgram showing that the trigger skim cut passed OR that we
247  are skipping this selection. */
248  cutIndexPassed++;
249  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
250  B2DEBUG(22, "Cutflow: Trigger cut passed: index = " << cutIndexPassed);
251 
252 
253 
254 
255 
256  /* Use ECLChannelMapper to get other detector indices for the crystals
257  For conversion from CellID to crate, shaper, and channel ids.
258  The initialization function automatically checks to see if the
259  object has been initialized and ifthe payload has changed and
260  thus needs updating. */
261  bool ECLchannelMapHasChanged = m_channelMapDB.hasChanged();
262  if (ECLchannelMapHasChanged) {
263  B2INFO("eclBhabhaTimeCalibrationValidationCollectorModule::collect() " << LogVar("ECLchannelMapHasChanged",
264  ECLchannelMapHasChanged));
265  if (!m_crystalMapper->initFromDB()) {
266  B2FATAL("eclBhabhaTimeCalibrationValidationCollectorModule::collect() : Can't initialize eclChannelMapper!");
267  }
268  }
269 
270 
271 
272 
273 
274  B2DEBUG(29, "Finished checking if previous crystal time payload has changed");
275 
276  if (m_CrateTimeDB.hasChanged()) {
277  m_CrateTime = m_CrateTimeDB->getCalibVector();
278  m_CrateTimeUnc = m_CrateTimeDB->getCalibUncVector();
279  }
280 
281  B2DEBUG(25, "eclBhabhaTimeCalibrationValidationCollector:: loaded ECLCrateTimeOffset from the database"
282  << LogVar("IoV", m_CrateTimeDB.getIoV())
283  << LogVar("Checksum", m_CrateTimeDB.getChecksum()));
284 
285  // Conversion coefficient from ADC ticks to nanoseconds
286  // TICKS_TO_NS ~ 0.4931 ns/clock tick
287  // 1/(4fRF) = 0.4913 ns/clock tick, where fRF is the accelerator RF frequency
288  const double TICKS_TO_NS = 1.0 / (4.0 * EclConfiguration::getRF()) * 1e3;
289 
290 
291  vector<float> Crate_time_ns(52, 0.0);
292  vector<float> Crate_time_unc_ns(52, 0.0);
294  // Make a crate time offset vector with an entry per crate (instead of per crystal) and convert from ADC counts to ns.
295  for (int crysID = 1; crysID <= ECLElementNumbers::c_NCrystals; crysID++) {
296  int crateID_temp = m_crystalMapper->getCrateID(crysID);
297  Crate_time_ns[crateID_temp - 1] = m_CrateTime[crysID] * TICKS_TO_NS;
298  Crate_time_unc_ns[crateID_temp - 1] = m_CrateTimeUnc[crysID] * TICKS_TO_NS;
299  }
300 
301 
302 
303  // Storage crystal energies
305  for (auto& eclCalDigit : m_eclCalDigitArray) {
306  int tempCrysID = eclCalDigit.getCellId() - 1;
307  m_EperCrys[tempCrysID] = eclCalDigit.getEnergy();
308  }
309 
310 
311  // Getting the event t0 using the full event t0 rather than from the CDC specifically
312 
313  double evt_t0 = -1000 ;
314  double evt_t0_unc = -1000 ;
315 
316  // Determine if there is an event t0 to use and then extract the information about it
317  if (m_eventT0.isOptional()) {
318  if (!m_eventT0.isValid()) {
319  return;
320  }
321  if (!m_eventT0->hasEventT0()) {
322  return;
323  } else {
324  // Overall event t0 (combination of multiple event t0s from different detectors)
325  evt_t0 = m_eventT0->getEventT0() ;
326  evt_t0_unc = m_eventT0->getEventT0Uncertainty() ;
327  }
328  B2DEBUG(26, "Found event t0") ;
329  }
330 
331 
332  //---------------------------------------------------------------------
333  //..Some utilities
334  PCmsLabTransform boostrotate;
335 
336  //---------------------------------------------------------------------
337  //..Track properties. Use pion (211) mass hypothesis,
338  // which is the only particle hypothesis currently available???
339  double maxp[2] = {0., 0.};
340  int maxiTrk[2] = { -1, -1};
341  int nTrkAll = tracks.getEntries() ;
342 
343  int nTrkLoose = 0 ;
344  int nTrkTight = 0 ;
347  /* Loop over all the tracks to define the tight and loose selection tracks
348  We will select events with only 2 tight tracks and no additional loose tracks.
349  Tight tracks are a subset of looses tracks. */
350  for (int iTrk = 0 ; iTrk < nTrkAll ; iTrk++) {
351  // Get track biasing towards the particle being a pion based on what particle types
352  // are used for reconstruction at this stage.
353  const TrackFitResult* tempTrackFit = tracks[iTrk]->getTrackFitResultWithClosestMass(Const::pion);
354  if (not tempTrackFit) {continue ;}
355 
356  // Collect track info to be used for categorizing
357  short charge = tempTrackFit->getChargeSign() ;
358  double z0 = tempTrackFit->getZ0() ;
359  double d0 = tempTrackFit->getD0() ;
360  int nCDChits = tempTrackFit->getHitPatternCDC().getNHits() ;
361  double p = tempTrackFit->getMomentum().R() ;
362 
363  /* Test if loose track */
364 
365  // d0 and z0 cuts
366  if (fabs(d0) > m_looseTrkD0) {
367  continue;
368  }
369  if (fabs(z0) > m_looseTrkZ0) {
370  continue;
371  }
372  // Number of hits in the CDC
373  if (nCDChits < 1) {
374  continue;
375  }
376  nTrkLoose++;
377 
378 
379 
380  /* Test if the loose track is also a tight track */
381 
382  // Number of hits in the CDC
383  if (nCDChits < 20) {
384  continue;
385  }
386  // d0 and z0 cuts
387  if (fabs(d0) > m_tightTrkD0) {
388  continue;
389  }
390  if (fabs(z0) > m_tightTrkZ0) {
391  continue;
392  }
393  nTrkTight++;
394 
395  // Sorting of tight tracks. Not really required as we only want two tight tracks (at the moment) but okay.
396  //..Find the maximum p negative [0] and positive [1] tracks
397  int icharge = 0;
398  if (charge > 0) {icharge = 1;}
399  if (p > maxp[icharge]) {
400  maxp[icharge] = p;
401  maxiTrk[icharge] = iTrk;
402  }
403 
404  }
405  /* After that last section the numbers of loose and tight tracks are known as well as the
406  index of the loose tracks that have the highest p negatively charged and highest p positively
407  charged tracks as measured in the centre of mass frame */
408 
409 
410  if (nTrkTight != 2) {
411  return;
412  }
413  // There are exactly two tight tracks
414  cutIndexPassed++;
415  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
416  B2DEBUG(22, "Cutflow: Two tight tracks: index = " << cutIndexPassed);
417 
418 
419  if (nTrkLoose != 2) {
420  return;
421  }
422  // There are exactly two loose tracks as well, i.e. no additional loose tracks
423  cutIndexPassed++ ;
424  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed) ;
425  B2DEBUG(22, "Cutflow: No additional loose tracks: index = " << cutIndexPassed) ;
426  /* Determine if the two tracks have the opposite electric charge.
427  We know this because the track indices stores the max pt track in [0] for negatively charged track
428  and [1] fo the positively charged track. If both are filled then both a negatively charged
429  and positively charged track were found. */
430  bool oppositelyChargedTracksPassed = maxiTrk[0] != -1 && maxiTrk[1] != -1;
431  if (!oppositelyChargedTracksPassed) {
432  return;
433  }
434  // The two tracks have the opposite electric charges.
435  cutIndexPassed++;
436  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
437  B2DEBUG(22, "Cutflow: Oppositely charged tracks: index = " << cutIndexPassed);
438 
439 
440 
441  //---------------------------------------------------------------------
442  /* Determine associated energy clusters to each of the two tracks. Sum the energies of the
443  multiple clusters to each track and find the crystal with the maximum energy within all
444  the sets of clusters associated to the tracks. Extract the good cluster times.*/
445  double trkEClustLab[2] = {0., 0.};
446  double trkEClustCOM[2] = {0., 0.};
447  double trkpLab[2];
448  double trkpCOM[2];
449  ROOT::Math::PxPyPzEVector trkp4Lab[2];
450  ROOT::Math::PxPyPzEVector trkp4COM[2];
451 
452  // Index of the cluster and the crystal that has the highest energy crystal for the two tracks
453  int numClustersPerTrack[2] = { 0, 0 };
454  double E_DIV_p[2];
455 
456  vector<double> goodClustTimes ;
457  vector<double> goodClust_dt99 ;
458  vector<double> goodClustE ;
459  vector<int> goodClustMaxEcrys_cid ;
460 
461  for (int icharge = 0; icharge < 2; icharge++) {
462  if (maxiTrk[icharge] > -1) {
463  B2DEBUG(22, "looping over the 2 max pt tracks");
464 
465  const TrackFitResult* tempTrackFit = tracks[maxiTrk[icharge]]->getTrackFitResultWithClosestMass(Const::pion);
466  if (not tempTrackFit) {continue ;}
467 
468  trkp4Lab[icharge] = tempTrackFit->get4Momentum();
469  trkp4COM[icharge] = boostrotate.rotateLabToCms() * trkp4Lab[icharge];
470  trkpLab[icharge] = trkp4Lab[icharge].P();
471  trkpCOM[icharge] = trkp4COM[icharge].P();
472 
473 
474  /* For each cluster associated to the current track, sum up the energies to get the total
475  energy of all clusters associated to the track and find which crystal has the highest
476  energy from all those clusters*/
477  auto eclClusterRelationsFromTracks = tracks[maxiTrk[icharge]]->getRelationsTo<ECLCluster>();
478  for (unsigned int clusterIdx = 0; clusterIdx < eclClusterRelationsFromTracks.size(); clusterIdx++) {
479 
480  B2DEBUG(22, "Looking at clusters. index = " << clusterIdx);
481  auto cluster = eclClusterRelationsFromTracks[clusterIdx];
482 
483  if (cluster->hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons)) {
484  numClustersPerTrack[icharge]++;
485  double eClust = cluster->getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons);
486  double electronTime = cluster->getTime();
487  bool badElectronTime = cluster->hasFailedFitTime();
488  bool badElectronTimeResolution = cluster->hasFailedTimeResolution();
489  if ((fabs(electronTime) < m_timeAbsMax) &&
490  (!badElectronTime) &&
491  (!badElectronTimeResolution)) {
492  trkEClustLab[icharge] += eClust ;
493  short cid = cluster->getMaxECellId() ;
494  goodClustMaxEcrys_cid.push_back(cid) ;
495  goodClustTimes.push_back(electronTime) ;
496  goodClust_dt99.push_back(cluster->getDeltaTime99()) ;
497  goodClustE.push_back(eClust);
498  }
499  }
500  }
501  trkEClustCOM[icharge] = trkEClustLab[icharge] * trkpCOM[icharge] / trkpLab[icharge];
502 
503  // Check both electrons to see if their cluster energy / track momentum is good.
504  // The Belle II physics book shows that this is the main way of separating electrons from other particles
505  // Done in the centre of mass reference frame although I believe E/p is invariant under a boost.
506  E_DIV_p[icharge] = trkEClustCOM[icharge] / trkpCOM[icharge];
507 
508  }
509  }
510  /* At the end of this section the 3-momenta magnitudes and the cluster energies are known
511  for the two saved track indices for both the lab and COM frames.
512  The crystal with the maximum energy, one associated to each track, is recorded*/
513  B2DEBUG(26, "Extracted time information and E/p for the tracks") ;
514 
515 
516 
517  /* Cut on the number of ECL cluster connected to tracks
518 
519  THIS IS DIFFERENT FROM THE CODE THAT PERFORMS THE CALIBRATIONS. THIS VALIDATIONS REQUIRES
520  THAT THERE ARE EXACTLY TWO CLUSTERS ASSOCIATED TO THE TRACKS WHILE THE CALIBRATION
521  CODE ALLOWS FOR MORE THAN ONE CLUSTER PER TRACK. THIS VALIDATION ALSO DOES NOT CUT ON THE
522  NUMBER OF EXTRA CLUSTERS NOT ASSOCIATED TO THE 2 TRACKS, WHICH IS A LOOSER CUT THAN USED
523  TO PERFORM THE CALIBRATION. */
524  long unsigned int numGoodElectronClusters_cut = 2 ;
525  if (goodClustTimes.size() != numGoodElectronClusters_cut) {
526  return ;
527  }
528  // There is exactly two ECL clusters connected to tracks in the event
529  cutIndexPassed++ ;
530  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed) ;
531  B2DEBUG(22, "Cutflow: Exactly " << numGoodElectronClusters_cut
532  << " good clusters connected to tracks: index = " << cutIndexPassed);
533 
534 
535  // Check both electrons to see if their cluster energy / track momentum is good.
536  // The Belle II physics book shows that this is the main way of separating electrons from other particles
537  // Done in the centre of mass reference frame although I believe E/p is invariant under a boost.
538  bool E_DIV_p_instance_passed[2] = {false, false};
539  double E_DIV_p_CUT = 0.7;
540  for (int icharge = 0; icharge < 2; icharge++) {
541  E_DIV_p_instance_passed[icharge] = E_DIV_p[icharge] > E_DIV_p_CUT;
542  }
543  if (!E_DIV_p_instance_passed[0] || !E_DIV_p_instance_passed[1]) {
544  return;
545  }
546  // E/p sufficiently large
547  cutIndexPassed++;
548  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
549  B2DEBUG(22, "Cutflow: E_i/p_i > " << E_DIV_p_CUT << ": index = " << cutIndexPassed);
550 
551 
552 
553  // Cut on the invariant mass of the tracks in the event
554  double invMassTrk = (trkp4Lab[0] + trkp4Lab[1]).M();
555  double invMass_CUT = 0.9;
556 
557  bool invMassCutsPassed = invMassTrk > (invMass_CUT * boostrotate.getCMSEnergy());
558  if (!invMassCutsPassed) {
559  return;
560  }
561  // Invariable mass of the two tracks are above the minimum
562  cutIndexPassed++;
563  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
564  B2DEBUG(22, "Cutflow: m(track 1+2) > " << invMass_CUT << "*E_COM = " << invMass_CUT << " * " << boostrotate.getCMSEnergy() <<
565  " : index = " << cutIndexPassed);
566 
567 
568  B2DEBUG(22, "Event passed all cuts");
569 
570 
571  // Fill the histogram for the event level variables
572  getObjectPtr<TH1F>("eventT0")->Fill(evt_t0) ;
573 
574  bool isCDCt0 = (static_cast<EventT0::EventT0Component>(*m_eventT0->getEventT0Component())).detectorSet.contains(Const::CDC);
575  bool isECLt0 = (static_cast<EventT0::EventT0Component>(*m_eventT0->getEventT0Component())).detectorSet.contains(Const::ECL);
576  string t0Detector = "UNKNOWN... WHY?";
577  if (isCDCt0) {
578  t0Detector = "CDC" ;
579  } else if (isECLt0) {
580  t0Detector = "ECL" ;
581  }
582 
583  B2DEBUG(26, "t0 = " << evt_t0 << " ns. t0 is from CDC?=" << isCDCt0 << ", t0 is from ECL?=" << isECLt0 << " t0 from " <<
584  t0Detector);
585 
586 
587  //=== For each good electron cluster in the processed event and fill histogram.
588  for (long unsigned int i = 0 ; i < goodClustTimes.size() ; i++) {
589  getObjectPtr<TH1F>("clusterTime")->Fill(goodClustTimes[i]) ;
590  getObjectPtr<TH2F>("clusterTime_cid")->Fill(goodClustMaxEcrys_cid[i] + 0.001, goodClustTimes[i], 1) ;
591  getObjectPtr<TH2F>("clusterTime_run")->Fill(m_EventMetaData->getRun() + 0.001, goodClustTimes[i], 1) ;
592  getObjectPtr<TH2F>("clusterTimeClusterE")->Fill(goodClustE[i], goodClustTimes[i], 1) ;
593  getObjectPtr<TH2F>("dt99_clusterE")->Fill(goodClustE[i], goodClust_dt99[i], 1) ;
594 
595 
596  //== Save debug TTree with detailed information if necessary.
597  if (m_saveTree) {
598 
599  m_tree_time = goodClustTimes[i] ;
600  m_tree_t0 = evt_t0 ;
601  m_tree_t0_unc = evt_t0_unc ;
602  m_E_electron_clust = goodClustE[i] ;
603  m_NtightTracks = nTrkTight ;
604  m_tree_evt_num = m_EventMetaData->getEvent() ;
605  m_tree_run = m_EventMetaData->getRun() ;
606  m_tree_cid = goodClustMaxEcrys_cid[i] ;
607  m_tree_dt99 = goodClust_dt99[i] ;
608 
610 
611  }
612  }
613  B2DEBUG(26, "Filled cluster tree") ;
614 
615  //=== Fill histogram for cluster time difference of the two electrons
616  double tDiff;
617  if (goodClustE[0] > goodClustE[1]) {
618  tDiff = goodClustTimes[0] - goodClustTimes[1];
619  } else {
620  tDiff = goodClustTimes[1] - goodClustTimes[0];
621  }
622 
623  getObjectPtr<TH1F>("clusterTimeE0E1diff")->Fill(tDiff) ;
624 
625 
626 
627  if (m_saveTree) {
628  m_tree_t0 = evt_t0 ;
629  m_tree_t0_unc = evt_t0_unc ;
630  m_tree_evt_num = m_EventMetaData->getEvent() ;
631  m_tree_run = m_EventMetaData->getRun() ;
632  m_NtightTracks = nTrkTight ;
633  m_tree_E0 = goodClustE[0] ;
634  m_tree_E1 = goodClustE[1] ;
635  m_tree_time_fromE0 = goodClustTimes[0] ;
636  m_tree_time_fromE1 = goodClustTimes[1] ;
637 
638  m_dbg_tree_event->Fill() ;
639 
640 
641  int runNum = m_EventMetaData->getRun();
642  if (m_tree_PreviousRun != runNum) {
643  for (int icrate = 1; icrate <= 52; icrate++) {
644  m_tree_run = runNum ;
645  m_tree_crateid = icrate ;
646  m_tree_tcrate = Crate_time_ns[icrate] ;
647  m_tree_tcrate_unc = Crate_time_unc_ns[icrate] ;
648 
649  m_dbg_tree_run->Fill() ;
650  }
652  }
653  }
654 
655  B2DEBUG(26, "Filled event tree") ;
656 
657 }
Calibration collector module base class.
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
bool hasChanged()
Check whether the object has changed since the last call to hasChanged of the accessor).
ECL cluster data.
Definition: ECLCluster.h:27
@ c_nPhotons
CR is split into n photons (N1)
static double getRF()
See m_rf.
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
Class to hold Lorentz transformations from/to CMS and boost vector.
double getCMSEnergy() const
Returns CMS energy of e+e- (aka.
const ROOT::Math::LorentzRotation rotateLabToCms() const
Returns Lorentz transformation from Lab to CMS.
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
ROOT::Math::PxPyPzEVector get4Momentum() const
Getter for the 4Momentum at the closest approach of the track in the r/phi projection.
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
int m_tree_PreviousRun
Run number for the previous run for debug TTree output.
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.
DBObjPtr< Belle2::ECLChannelMap > m_channelMapDB
Mapper of ecl channels to various other objects, like crates.
void collect() override
Select events and crystals and accumulate histograms.
std::vector< float > m_CrateTimeUnc
uncertainty vector obtained from DB object
StoreArray< ECLCalDigit > m_eclCalDigitArray
Required input array of ECLCalDigits.
void inDefineHisto() override
Replacement for defineHisto() in CalibrationCollector modules.
std::unique_ptr< Belle2::ECL::ECLChannelMapper > m_crystalMapper
ECL object for keeping track of mapping between crystals and crates etc.
StoreObjPtr< SoftwareTriggerResult > m_TrgResult
Store array for Trigger selection.
Class to store variables with their name which were sent to the logging service.
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
@ c_accept
Accept this event.
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