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