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