Belle II Software  release-05-02-19
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 <tracking/dataobjects/RecoTrack.h>
38 #include <ecl/digitization/EclConfiguration.h>
39 
40 #include <analysis/utility/PCmsLabTransform.h>
41 #include <analysis/ClusterUtility/ClusterUtils.h>
42 #include <boost/optional.hpp>
43 
44 #include <TH2F.h>
45 #include <TTree.h>
46 #include <TFile.h>
47 
48 using namespace Belle2;
49 using namespace ECL;
50 using namespace std;
51 
52 //-----------------------------------------------------------------
53 // Register the Module
54 //-----------------------------------------------------------------
55 REG_MODULE(eclBhabhaTimeCalibrationValidationCollector)
56 
57 //-----------------------------------------------------------------
58 // Implementation
59 //-----------------------------------------------------------------
60 
63  m_dbg_tree_electronClusters(0),
64  m_dbg_tree_event(0),
65  m_dbg_tree_run(0),
66  m_CrateTimeDB("ECLCrateTimeOffset")//,
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 
84 
85  // specify this flag if you need parallel processing
86  setPropertyFlags(c_ParallelProcessingCertified);
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  StoreObjPtr<EventMetaData> evtMetaData ;
142  B2INFO("eclBhabhaTimeCalibrationValidationCollector: Experiment = " << evtMetaData->getExperiment() <<
143  " run = " << evtMetaData->getRun()) ;
144 
145  //=== Create histograms and register them in the data store
146 
147  // Define the bin size, which is equivalent to the
148  double binSize = 2000.0 / pow(2, 12);
149  double halfBinSize = binSize / 2.0;
150 
151  /* Determine the number of bins required to go from the edge of the bin centred
152  on zero to a value just larger than the negative cut off */
153  double nBinsNeg = floor((m_timeAbsMax - halfBinSize) / binSize);
154  double min_t = -nBinsNeg * binSize - halfBinSize; // lower edge value of left most bin
155  int nbins = nBinsNeg * 2 + 1; // number of negative time bins + t=0 bin + number of positive time bins
156  double max_t = min_t + nbins * binSize; // upper edge value of right most bin
157 
158  /* Variable bin width information for the time information vs energy since
159  the time width should vary as a function of 1/E */
160  const Int_t N_E_BIN_EDGES = 64;
161  const Int_t N_E_BINS = N_E_BIN_EDGES - 1;
162  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};
163 
164 
165  auto cutflow = new TH1F("cutflow", " ;Cut label number ;Number of events passing cut", 10, 0, 10) ;
166  registerObject<TH1F>("cutflow", cutflow) ;
167 
168  auto clusterTime = new TH1F("clusterTime", ";Electron ECL cluster time [ns]; number of ECL clusters", nbins, min_t, max_t) ;
169  registerObject<TH1F>("clusterTime", clusterTime) ;
170 
171  auto clusterTime_cid = new TH2F("clusterTime_cid",
172  ";crystal Cell ID ;Electron ECL cluster time [ns]", 8736, 1, 8736 + 1, nbins, 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 
206 {
207  int cutIndexPassed = 0;
208  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
209  B2DEBUG(22, "Cutflow: no cuts: index = " << cutIndexPassed);
210 
211 
212  /* Use ECLChannelMapper to get other detector indices for the crystals */
213  /* For conversion from CellID to crate, shaper, and channel ids. */
214 
215  // Use smart pointer to avoid memory leak when the ECLChannelMapper object needs destroying at the end of the event.
216  shared_ptr< ECL::ECLChannelMapper > crystalMapper(new ECL::ECLChannelMapper());
217  crystalMapper->initFromDB();
218 
219  B2DEBUG(29, "Finished checking if previous crystal time payload has changed");
220 
221  if (m_CrateTimeDB.hasChanged()) {
222  m_CrateTime = m_CrateTimeDB->getCalibVector();
223  m_CrateTimeUnc = m_CrateTimeDB->getCalibUncVector();
224  }
225 
226  B2DEBUG(25, "eclBhabhaTimeCalibrationValidationCollector:: loaded ECLCrateTimeOffset from the database"
227  << LogVar("IoV", m_CrateTimeDB.getIoV())
228  << LogVar("Checksum", m_CrateTimeDB.getChecksum()));
229 
230  // Conversion coefficient from ADC ticks to nanoseconds
231  // TICKS_TO_NS ~ 0.4931 ns/clock tick
232  // 1/(4fRF) = 0.4913 ns/clock tick, where fRF is the accelerator RF frequency, fRF=508.889 MHz.
233  const double TICKS_TO_NS = 1.0 / (4.0 * EclConfiguration::m_rf) * 1e3;
234 
235 
236  vector<float> Crate_time_ns(52, 0.0);
237  vector<float> Crate_time_unc_ns(52, 0.0);
239  // Make a crate time offset vector with an entry per crate (instead of per crystal) and convert from ADC counts to ns.
240  for (int crysID = 1; crysID <= 8736; crysID++) {
241  int crateID_temp = crystalMapper->getCrateID(crysID);
242  Crate_time_ns[crateID_temp - 1] = m_CrateTime[crysID] * TICKS_TO_NS;
243  Crate_time_unc_ns[crateID_temp - 1] = m_CrateTimeUnc[crysID] * TICKS_TO_NS;
244  }
245 
246 
247 
248  // Storage crystal energies
249  m_EperCrys.resize(8736);
250  for (auto& eclCalDigit : m_eclCalDigitArray) {
251  int tempCrysID = eclCalDigit.getCellId() - 1;
252  m_EperCrys[tempCrysID] = eclCalDigit.getEnergy();
253  }
254 
255 
256  // Getting the event t0 using the full event t0 rather than from the CDC specifically
257 
258  double evt_t0 = -1000 ;
259  double evt_t0_unc = -1000 ;
260 
261  // Determine if there is an event t0 to use and then extract the information about it
262  if (m_eventT0.isOptional()) {
263  if (!m_eventT0.isValid()) {
264  return;
265  }
266  if (!m_eventT0->hasEventT0()) {
267  return;
268  } else {
269  // Overall event t0 (combination of multiple event t0s from different detectors)
270  evt_t0 = m_eventT0->getEventT0() ;
271  evt_t0_unc = m_eventT0->getEventT0Uncertainty() ;
272  }
273  B2DEBUG(26, "Found event t0") ;
274  }
275 
276 
277  //---------------------------------------------------------------------
278  //..Some utilities
279  ClusterUtils cUtil;
280  const TVector3 clustervertex = cUtil.GetIPPosition();
281  PCmsLabTransform boostrotate;
282 
283  //---------------------------------------------------------------------
284  //..Track properties. Use pion (211) mass hypothesis,
285  // which is the only particle hypothesis currently available???
286  double maxp[2] = {0., 0.};
287  int maxiTrk[2] = { -1, -1};
288  int nTrkAll = tracks.getEntries() ;
289 
290  int nTrkLoose = 0 ;
291  int nTrkTight = 0 ;
294  /* Loop over all the tracks to define the tight and loose selection tracks
295  We will select events with only 2 tight tracks and no additional loose tracks.
296  Tight tracks are a subset of looses tracks. */
297  for (int iTrk = 0 ; iTrk < nTrkAll ; iTrk++) {
298  // Get track biasing towards the particle being a pion based on what particle types
299  // are used for reconstruction at this stage.
300  const TrackFitResult* tempTrackFit = tracks[iTrk]->getTrackFitResultWithClosestMass(Const::pion);
301  if (not tempTrackFit) {continue ;}
302 
303  // Collect track info to be used for categorizing
304  short charge = tempTrackFit->getChargeSign() ;
305  double z0 = tempTrackFit->getZ0() ;
306  double d0 = tempTrackFit->getD0() ;
307  int nCDChits = tempTrackFit->getHitPatternCDC().getNHits() ;
308  double p = tempTrackFit->getMomentum().Mag() ;
309 
310  /* Test if loose track */
311 
312  // d0 and z0 cuts
313  if (fabs(d0) > m_looseTrkD0) {
314  continue;
315  }
316  if (fabs(z0) > m_looseTrkZ0) {
317  continue;
318  }
319  // Number of hits in the CDC
320  if (nCDChits < 1) {
321  continue;
322  }
323  nTrkLoose++;
324 
325 
326 
327  /* Test if the loose track is also a tight track */
328 
329  // Number of hits in the CDC
330  if (nCDChits < 20) {
331  continue;
332  }
333  // d0 and z0 cuts
334  if (fabs(d0) > m_tightTrkD0) {
335  continue;
336  }
337  if (fabs(z0) > m_tightTrkZ0) {
338  continue;
339  }
340  nTrkTight++;
341 
342  // Sorting of tight tracks. Not really required as we only want two tight tracks (at the moment) but okay.
343  //..Find the maximum p negative [0] and positive [1] tracks
344  int icharge = 0;
345  if (charge > 0) {icharge = 1;}
346  if (p > maxp[icharge]) {
347  maxp[icharge] = p;
348  maxiTrk[icharge] = iTrk;
349  }
350 
351  }
352  /* After that last section the numbers of loose and tight tracks are known as well as the
353  index of the loose tracks that have the highest p negatively charged and highest p positively
354  charged tracks as measured in the centre of mass frame */
355 
356 
357  if (nTrkTight != 2) {
358  return;
359  }
360  // There are exactly two tight tracks
361  cutIndexPassed++;
362  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
363  B2DEBUG(22, "Cutflow: Two tight tracks: index = " << cutIndexPassed);
364 
365 
366  if (nTrkLoose != 2) {
367  return;
368  }
369  // There are exactly two loose tracks as well, i.e. no additional loose tracks
370  cutIndexPassed++ ;
371  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed) ;
372  B2DEBUG(22, "Cutflow: No additional loose tracks: index = " << cutIndexPassed) ;
373  /* Determine if the two tracks have the opposite electric charge.
374  We know this because the track indices stores the max pt track in [0] for negatively charged track
375  and [1] fo the positively charged track. If both are filled then both a negatively charged
376  and positively charged track were found. */
377  bool oppositelyChargedTracksPassed = maxiTrk[0] != -1 && maxiTrk[1] != -1;
378  if (!oppositelyChargedTracksPassed) {
379  return;
380  }
381  // The two tracks have the opposite electric charges.
382  cutIndexPassed++;
383  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
384  B2DEBUG(22, "Cutflow: Oppositely charged tracks: index = " << cutIndexPassed);
385 
386 
387 
388  //---------------------------------------------------------------------
389  /* Determine associated energy clusters to each of the two tracks. Sum the energies of the
390  multiple clusters to each track and find the crystal with the maximum energy within all
391  the sets of clusters associated to the tracks. Extract the good cluster times.*/
392  double trkEClustLab[2] = {0., 0.};
393  double trkEClustCOM[2] = {0., 0.};
394  double trkpLab[2];
395  double trkpCOM[2];
396  TLorentzVector trkp4Lab[2];
397  TLorentzVector trkp4COM[2];
398 
399  // Index of the cluster and the crystal that has the highest energy crystal for the two tracks
400  int numClustersPerTrack[2] = { 0, 0 };
401  double E_DIV_p[2];
402 
403  vector<double> goodClustTimes ;
404  vector<double> goodClust_dt99 ;
405  vector<double> goodClustE ;
406  vector<int> goodClustMaxEcrys_cid ;
407 
408  for (int icharge = 0; icharge < 2; icharge++) {
409  if (maxiTrk[icharge] > -1) {
410  B2DEBUG(22, "looping over the 2 max pt tracks");
411 
412  const TrackFitResult* tempTrackFit = tracks[maxiTrk[icharge]]->getTrackFitResultWithClosestMass(Const::pion);
413  if (not tempTrackFit) {continue ;}
414 
415  trkp4Lab[icharge] = tempTrackFit->get4Momentum();
416  trkp4COM[icharge] = boostrotate.rotateLabToCms() * trkp4Lab[icharge];
417  trkpLab[icharge] = trkp4Lab[icharge].Rho();
418  trkpCOM[icharge] = trkp4COM[icharge].Rho();
419 
420 
421  /* For each cluster associated to the current track, sum up the energies to get the total
422  energy of all clusters associated to the track and find which crystal has the highest
423  energy from all those clusters*/
424  auto eclClusterRelationsFromTracks = tracks[maxiTrk[icharge]]->getRelationsTo<ECLCluster>();
425  for (unsigned int clusterIdx = 0; clusterIdx < eclClusterRelationsFromTracks.size(); clusterIdx++) {
426 
427  B2DEBUG(22, "Looking at clusters. index = " << clusterIdx);
428  auto cluster = eclClusterRelationsFromTracks[clusterIdx];
429 
430  if (cluster->hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons)) {
431  numClustersPerTrack[icharge]++;
432  double eClust = cluster->getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons);
433  double electronTime = cluster->getTime();
434  bool badElectronTime = cluster->hasFailedFitTime();
435  bool badElectronTimeResolution = cluster->hasFailedTimeResolution();
436  if ((fabs(electronTime) < m_timeAbsMax) &&
437  (!badElectronTime) &&
438  (!badElectronTimeResolution)) {
439  trkEClustLab[icharge] += eClust ;
440  short cid = cluster->getMaxECellId() ;
441  goodClustMaxEcrys_cid.push_back(cid) ;
442  goodClustTimes.push_back(electronTime) ;
443  goodClust_dt99.push_back(cluster->getDeltaTime99()) ;
444  goodClustE.push_back(eClust);
445  }
446  }
447  }
448  trkEClustCOM[icharge] = trkEClustLab[icharge] * trkpCOM[icharge] / trkpLab[icharge];
449 
450  // Check both electrons to see if their cluster energy / track momentum is good.
451  // The Belle II physics book shows that this is the main way of separating electrons from other particles
452  // Done in the centre of mass reference frame although I believe E/p is invariant under a boost.
453  E_DIV_p[icharge] = trkEClustCOM[icharge] / trkpCOM[icharge];
454 
455  }
456  }
457  /* At the end of this section the 3-momenta magnitudes and the cluster energies are known
458  for the two saved track indices for both the lab and COM frames.
459  The crystal with the maximum energy, one associated to each track, is recorded*/
460  B2DEBUG(26, "Extracted time information and E/p for the tracks") ;
461 
462 
463 
464  /* Cut on the number of ECL cluster connected to tracks
465 
466  THIS IS DIFFERENT FROM THE CODE THAT PERFORMS THE CALIBRATIONS. THIS VALIDATIONS REQUIRES
467  THAT THERE ARE EXACTLY TWO CLUSTERS ASSOCIATED TO THE TRACKS WHILE THE CALIBRATION
468  CODE ALLOWS FOR MORE THAN ONE CLUSTER PER TRACK. THIS VALIDATION ALSO DOES NOT CUT ON THE
469  NUMBER OF EXTRA CLUSTERS NOT ASSOCIATED TO THE 2 TRACKS, WHICH IS A LOOSER CUT THAN USED
470  TO PERFORM THE CALIBRATION. */
471  long unsigned int numGoodElectronClusters_cut = 2 ;
472  if (goodClustTimes.size() != numGoodElectronClusters_cut) {
473  return ;
474  }
475  // There is exactly two ECL clusters connected to tracks in the event
476  cutIndexPassed++ ;
477  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed) ;
478  B2DEBUG(22, "Cutflow: Exactly " << numGoodElectronClusters_cut
479  << " good clusters connected to tracks: index = " << cutIndexPassed);
480 
481 
482  // Check both electrons to see if their cluster energy / track momentum is good.
483  // The Belle II physics book shows that this is the main way of separating electrons from other particles
484  // Done in the centre of mass reference frame although I believe E/p is invariant under a boost.
485  bool E_DIV_p_instance_passed[2] = {false, false};
486  double E_DIV_p_CUT = 0.7;
487  for (int icharge = 0; icharge < 2; icharge++) {
488  E_DIV_p_instance_passed[icharge] = E_DIV_p[icharge] > E_DIV_p_CUT;
489  }
490  if (!E_DIV_p_instance_passed[0] || !E_DIV_p_instance_passed[1]) {
491  return;
492  }
493  // E/p sufficiently large
494  cutIndexPassed++;
495  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
496  B2DEBUG(22, "Cutflow: E_i/p_i > " << E_DIV_p_CUT << ": index = " << cutIndexPassed);
497 
498 
499 
500  // Cut on the invariant mass of the tracks in the event
501  double invMassTrk = (trkp4Lab[0] + trkp4Lab[1]).M();
502  double invMass_CUT = 0.9;
503 
504  bool invMassCutsPassed = invMassTrk > (invMass_CUT * boostrotate.getCMSEnergy());
505  if (!invMassCutsPassed) {
506  return;
507  }
508  // Invariable mass of the two tracks are above the minimum
509  cutIndexPassed++;
510  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
511  B2DEBUG(22, "Cutflow: m(track 1+2) > " << invMass_CUT << "*E_COM = " << invMass_CUT << " * " << boostrotate.getCMSEnergy() <<
512  " : index = " << cutIndexPassed);
513 
514 
515  B2DEBUG(22, "Event passed all cuts");
516 
517 
518  // Fill the histogram for the event level variables
519  getObjectPtr<TH1F>("eventT0")->Fill(evt_t0) ;
520 
521  bool isCDCt0 = (static_cast<EventT0::EventT0Component>(*m_eventT0->getEventT0Component())).detectorSet.contains(Const::CDC);
522  bool isECLt0 = (static_cast<EventT0::EventT0Component>(*m_eventT0->getEventT0Component())).detectorSet.contains(Const::ECL);
523  string t0Detector = "UNKNOWN... WHY?";
524  if (isCDCt0) {
525  t0Detector = "CDC" ;
526  } else if (isECLt0) {
527  t0Detector = "ECL" ;
528  }
529 
530  B2DEBUG(26, "t0 = " << evt_t0 << " ns. t0 is from CDC?=" << isCDCt0 << ", t0 is from ECL?=" << isECLt0 << " t0 from " <<
531  t0Detector);
532 
533 
534  //=== For each good electron cluster in the processed event and fill histogram.
535 
536  StoreObjPtr<EventMetaData> evtMetaData ;
537  for (long unsigned int i = 0 ; i < goodClustTimes.size() ; i++) {
538  getObjectPtr<TH1F>("clusterTime")->Fill(goodClustTimes[i]) ;
539  getObjectPtr<TH2F>("clusterTime_cid")->Fill(goodClustMaxEcrys_cid[i] + 0.001, goodClustTimes[i] , 1) ;
540  getObjectPtr<TH2F>("clusterTime_run")->Fill(evtMetaData->getRun() + 0.001, goodClustTimes[i] , 1) ;
541  getObjectPtr<TH2F>("clusterTimeClusterE")->Fill(goodClustE[i], goodClustTimes[i], 1) ;
542  getObjectPtr<TH2F>("dt99_clusterE")->Fill(goodClustE[i], goodClust_dt99[i], 1) ;
543 
544 
545  //== Save debug TTree with detailed information if necessary.
546  if (m_saveTree) {
547 
548  m_tree_time = goodClustTimes[i] ;
549  m_tree_t0 = evt_t0 ;
550  m_tree_t0_unc = evt_t0_unc ;
551  m_E_electron_clust = goodClustE[i] ;
552  m_NtightTracks = nTrkTight ;
553  m_tree_evt_num = evtMetaData->getEvent() ;
554  m_tree_run = evtMetaData->getRun() ;
555  m_tree_cid = goodClustMaxEcrys_cid[i] ;
556  m_tree_dt99 = goodClust_dt99[i] ;
557 
558  m_dbg_tree_electronClusters->Fill() ;
559 
560  }
561  }
562  B2DEBUG(26, "Filled cluster tree") ;
563 
564  //=== Fill histogram for cluster time difference of the two electrons
565  double tDiff;
566  if (goodClustE[0] > goodClustE[1]) {
567  tDiff = goodClustTimes[0] - goodClustTimes[1];
568  } else {
569  tDiff = goodClustTimes[1] - goodClustTimes[0];
570  }
571 
572  getObjectPtr<TH1F>("clusterTimeE0E1diff")->Fill(tDiff) ;
573 
574 
575 
576  if (m_saveTree) {
577  m_tree_t0 = evt_t0 ;
578  m_tree_t0_unc = evt_t0_unc ;
579  m_tree_evt_num = evtMetaData->getEvent() ;
580  m_tree_run = evtMetaData->getRun() ;
581  m_NtightTracks = nTrkTight ;
582  m_tree_E0 = goodClustE[0] ;
583  m_tree_E1 = goodClustE[1] ;
584  m_tree_time_fromE0 = goodClustTimes[0] ;
585  m_tree_time_fromE1 = goodClustTimes[1] ;
586 
587  m_dbg_tree_event->Fill() ;
588 
589 
590  int runNum = evtMetaData->getRun();
591  if (m_tree_PreviousRun != runNum) {
592  for (int icrate = 1; icrate <= 52; icrate++) {
593  m_tree_run = runNum ;
594  m_tree_crateid = icrate ;
595  m_tree_tcrate = Crate_time_ns[icrate] ;
596  m_tree_tcrate_unc = Crate_time_unc_ns[icrate] ;
597 
598  m_dbg_tree_run->Fill() ;
599  }
600  m_tree_PreviousRun = m_tree_run ;
601  }
602  }
603 
604  B2DEBUG(26, "Filled event tree") ;
605 
606 }
Belle2::CalibrationCollectorModule
Calibration collector module base class.
Definition: CalibrationCollectorModule.h:44
Belle2::eclBhabhaTimeCalibrationValidationCollectorModule::prepare
void prepare() override
Define histograms and read payloads from DB.
Definition: eclBhabhaTimeCalibrationValidationCollectorModule.cc:98
Belle2::ECLCluster
ECL cluster data.
Definition: ECLCluster.h:39
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::ECL::ECLChannelMapper::getCrateID
int getCrateID(int iCOPPERNode, int iFINESSE)
get crate number by given COPPER node number and FINESSE number
Definition: ECLChannelMapper.cc:211
Belle2::PCmsLabTransform::getCMSEnergy
double getCMSEnergy() const
Returns CMS energy of e+e- (aka.
Definition: PCmsLabTransform.h:57
Belle2::ECLCluster::EHypothesisBit::c_nPhotons
@ c_nPhotons
CR is split into n photons (N1)
Belle2::TrackFitResult
Values of the result of a track fit with a given particle hypothesis.
Definition: TrackFitResult.h:59
Belle2::EventT0::EventT0Component
Structure for storing the extracted event t0s together with its detector and its uncertainty.
Definition: EventT0.h:44
Belle2::ClusterUtils::GetIPPosition
const TVector3 GetIPPosition()
Returns default IP position from beam parameters.
Definition: ClusterUtils.cc:182
Belle2::eclBhabhaTimeCalibrationValidationCollectorModule
This module generates 'TimevsCrys' histogram to later (in eclBhabhaTAlgorithm) find time offset from ...
Definition: eclBhabhaTimeCalibrationValidationCollectorModule.h:38
Belle2::Const::pion
static const ChargedStable pion
charged pion particle
Definition: Const.h:535
Belle2::ECL::ECLChannelMapper
This class provides access to ECL channel map that is either a) Loaded from the database (see ecl/dbo...
Definition: ECLChannelMapper.h:36
Belle2::ClusterUtils
Class to provide momentum-related information from ECLClusters.
Definition: ClusterUtils.h:44
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
Belle2::ECL::EclConfiguration::m_rf
static constexpr double m_rf
accelerating RF, http://ptep.oxfordjournals.org/content/2013/3/03A006.full.pdf
Definition: EclConfiguration.h:45
Belle2::eclBhabhaTimeCalibrationValidationCollectorModule::~eclBhabhaTimeCalibrationValidationCollectorModule
virtual ~eclBhabhaTimeCalibrationValidationCollectorModule()
Module destructor.
Definition: eclBhabhaTimeCalibrationValidationCollectorModule.cc:89
Belle2::eclBhabhaTimeCalibrationValidationCollectorModule::collect
void collect() override
Select events and crystals and accumulate histograms.
Definition: eclBhabhaTimeCalibrationValidationCollectorModule.cc:205
Belle2::PCmsLabTransform
Class to hold Lorentz transformations from/to CMS and boost vector.
Definition: PCmsLabTransform.h:37
Belle2::PCmsLabTransform::rotateLabToCms
const TLorentzRotation rotateLabToCms() const
Returns Lorentz transformation from Lab to CMS.
Definition: PCmsLabTransform.h:74
Belle2::ECL::ECLChannelMapper::initFromDB
bool initFromDB()
Initialize channel mapper from the conditions database.
Definition: ECLChannelMapper.cc:105
Belle2::eclBhabhaTimeCalibrationValidationCollectorModule::inDefineHisto
void inDefineHisto() override
Replacement for defineHisto() in CalibrationCollector modules.
Definition: eclBhabhaTimeCalibrationValidationCollectorModule.cc:93