Belle II Software  release-08-00-10
ECLBhabhaTCollectorModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 /* Own header. */
10 #include <ecl/modules/eclBhabhaTCollector/ECLBhabhaTCollectorModule.h>
11 
12 /* ECL headers. */
13 #include <ecl/dataobjects/ECLCalDigit.h>
14 #include <ecl/dataobjects/ECLDigit.h>
15 #include <ecl/dataobjects/ECLElementNumbers.h>
16 #include <ecl/dbobjects/ECLCrystalCalib.h>
17 #include <ecl/dbobjects/ECLReferenceCrystalPerCrateCalib.h>
18 #include <ecl/digitization/EclConfiguration.h>
19 #include <ecl/geometry/ECLGeometryPar.h>
20 
21 /* Basf2 headers. */
22 #include <analysis/ClusterUtility/ClusterUtils.h>
23 #include <analysis/utility/PCmsLabTransform.h>
24 #include <framework/dataobjects/EventMetaData.h>
25 #include <framework/gearbox/Const.h>
26 #include <mdst/dataobjects/ECLCluster.h>
27 #include <mdst/dataobjects/HitPatternCDC.h>
28 #include <mdst/dataobjects/Track.h>
29 #include <tracking/dataobjects/RecoTrack.h>
30 
31 /* ROOT headers. */
32 #include <TFile.h>
33 #include <TH2F.h>
34 #include <TTree.h>
35 
36 using namespace Belle2;
37 using namespace ECL;
38 using namespace std;
39 
40 //-----------------------------------------------------------------
41 // Register the Module
42 //-----------------------------------------------------------------
43 REG_MODULE(ECLBhabhaTCollector);
44 
45 //-----------------------------------------------------------------
46 // Implementation
47 //-----------------------------------------------------------------
48 
50  m_ElectronicsDB("ECLCrystalElectronics"),
51  m_ElectronicsTimeDB("ECLCrystalElectronicsTime"),
52  m_FlightTimeDB("ECLCrystalFlightTime"),
53  m_PreviousCrystalTimeDB("ECLCrystalTimeOffset"),
54  m_CrateTimeDB("ECLCrateTimeOffset"),
55  m_RefCrystalsCalibDB("ECLReferenceCrystalPerCrateCalib"),
56  m_channelMapDB("ECLChannelMap")//,
57 {
58  setDescription("This module generates sum of all event times per crystal");
59 
60  addParam("timeAbsMax", m_timeAbsMax, // (Time in ns)
61  "Events with fabs(getTimeFit) > m_timeAbsMax "
62  "are excluded", (short)80);
63 
64  addParam("minCrystal", m_minCrystal,
65  "First CellId to handle.", 1);
66  addParam("maxCrystal", m_maxCrystal,
67  "Last CellId to handle.", ECLElementNumbers::c_NCrystals);
68 
69  addParam("saveTree", m_saveTree,
70  "If true, TTree 'tree' with more detailed event info is saved in "
71  "the output file specified by HistoManager",
72  false);
73 
74  addParam("looseTrkZ0", m_looseTrkZ0, "max Z0 for loose tracks (cm)", 10.);
75  addParam("tightTrkZ0", m_tightTrkZ0, "max Z0 for tight tracks (cm)", 2.);
76  addParam("looseTrkD0", m_looseTrkD0, "max D0 for loose tracks (cm)", 2.);
77  addParam("tightTrkD0", m_tightTrkD0, "max D0 for tight tracks (cm)", 0.5); // beam pipe radius = 1cm in 2019
78  addParam("skipTrgSel", skipTrgSel, "boolean to skip the trigger skim selection", false);
79 
80  addParam("hadronEventT0_TO_bhabhaEventT0_correction", m_hadronEventT0_TO_bhabhaEventT0_correction,
81  "CDC bhabha t0 bias correction (ns)", 0.);
82 
83  // specify this flag if you need parallel processing
85 }
86 
88 {
89 }
90 
92 {
93  //=== Prepare TTree for debug output
94  if (m_saveTree) { // /* >>>>>>>>>>>>>>>>>>>>>>>>>>>> if boolean true for saving debug trees >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
95  // Per electron
96  m_dbgTree_electrons = new TTree("tree_electrons", "Debug data for bhabha time calibration - one entry per electron");
97  m_dbgTree_electrons->Branch("EventNum", &m_tree_evtNum)->SetTitle("Event number");
98  m_dbgTree_electrons->Branch("CrystalCellID", &m_tree_cid)->SetTitle("Cell ID, 1..8736");
99  m_dbgTree_electrons->Branch("ADCamplitude", &m_tree_amp)->SetTitle("Amplitude, ADC units");
100  m_dbgTree_electrons->Branch("maxEcrystalEnergy", &m_tree_en)->SetTitle("Max Energy Crystal Energy, GeV");
101  m_dbgTree_electrons->Branch("maxEcrystalEnergyDIVclusterE",
102  &m_tree_E1Etot)->SetTitle("Max Energy Crystal Fraction Energy of Cluster");
103  m_dbgTree_electrons->Branch("E1divE2crystal",
104  &m_tree_E1E2)->SetTitle("Max Energy Crystal DIV second max energy crystal in cluster");
105  m_dbgTree_electrons->Branch("E1crystal_DIV_p", &m_tree_E1p)->SetTitle("Max Energy Crystal in cluster DIV track p");
106  m_dbgTree_electrons->Branch("timetsPreviousTimeCalibs",
107  &m_tree_timetsPreviousTimeCalibs)->SetTitle("Time t_psi after application of previous Ts, ns");
108  m_dbgTree_electrons->Branch("E_DIV_p", &m_E_DIV_p)->SetTitle("E DIV p");
109  m_dbgTree_electrons->Branch("timeF", &m_tree_timeF)->SetTitle("Time F, ns");
110  m_dbgTree_electrons->Branch("time_t_psi", &m_tree_time)->SetTitle("Time t_psi for Ts, ns");
111  m_dbgTree_electrons->Branch("quality", &m_tree_quality)->SetTitle("ECL FPGA fit quality, see Confluence article");
112  m_dbgTree_electrons->Branch("t0", &m_tree_t0)->SetTitle("T0, ns");
113  m_dbgTree_electrons->Branch("t0_unc", &m_tree_t0_unc)->SetTitle("T0 uncertainty, ns");
114  m_dbgTree_electrons->Branch("CrateID", &m_crystalCrate)->SetTitle("Crate id for crystal");
115  m_dbgTree_electrons->Branch("runNum", &m_runNum)->SetTitle("Run number");
116 
117  m_dbgTree_electrons->SetAutoSave(10);
118 
119 
120  // Per track
121  m_dbgTree_tracks = new TTree("tree_tracks", "Debug data for bhabha time calibration - one entry per track");
122  m_dbgTree_tracks->Branch("d0", &m_tree_d0)->SetTitle("d0, cm");
123  m_dbgTree_tracks->Branch("z0", &m_tree_z0)->SetTitle("z0, cm");
124  m_dbgTree_tracks->Branch("p", &m_tree_p)->SetTitle("track momentum, GeV");
125  m_dbgTree_tracks->Branch("charge", &m_charge)->SetTitle("track electric charge");
126  m_dbgTree_tracks->Branch("Num_CDC_hits", &m_tree_nCDChits)->SetTitle("Num CDC hits");
127 
128  m_dbgTree_tracks->SetAutoSave(10);
129 
130  // Per crystal
131  m_dbgTree_crystals = new TTree("tree_crystals",
132  "Debug data for bhabha time calibration - one entry per electron - one entry per crystal");
133  m_dbgTree_crystals->Branch("clustCrysE_DIV_maxEcrys",
134  &m_tree_clustCrysE_DIV_maxEcrys)->SetTitle("E of crystal i from cluster / E of max E crystal");
135  m_dbgTree_crystals->Branch("Crystal_E", &m_tree_clustCrysE) ->SetTitle("E of crystal i from cluster");
136  m_dbgTree_crystals->Branch("time_t_psi", &m_tree_time)->SetTitle("Time for Ts, ns");
137  m_dbgTree_crystals->Branch("Crystal_cell_ID", &m_tree_cid)->SetTitle("Cell ID, 1..8736");
138  m_dbgTree_crystals->Branch("quality", &m_tree_quality)->SetTitle("ECL FPGA fit quality, see Confluence article");
139 
140  m_dbgTree_crystals->SetAutoSave(10);
141 
142 
143  // Per event
144  m_dbgTree_event = new TTree("tree_event", "Debug data for bhabha time calibration - one entry per event");
145  m_dbgTree_event->Branch("massInvTracks", &m_massInvTracks)->SetTitle("Invariant mass of the two tracks");
146 
147  m_dbgTree_event->SetAutoSave(10);
148 
149 
150  m_dbgTree_evt_allCuts = new TTree("tree_evt_allCuts",
151  "Debug data for bhabha time calibration - one entry per event after all the cuts");
152  m_dbgTree_evt_allCuts->Branch("EclustPlus", &m_tree_enPlus)->SetTitle("Energy of cluster with +ve charge, GeV");
153  m_dbgTree_evt_allCuts->Branch("EclustNeg", &m_tree_enNeg)->SetTitle("Energy of cluster with -ve charge, GeV");
154  m_dbgTree_evt_allCuts->Branch("clustTimePos", &m_tree_tClustPos)->SetTitle("Cluster time of cluster with +ve charge, GeV");
155  m_dbgTree_evt_allCuts->Branch("clustTimeNeg", &m_tree_tClustNeg)->SetTitle("Cluster time of cluster with -ve charge, GeV");
156  m_dbgTree_evt_allCuts->Branch("maxEcrysTimePosClust",
157  &m_tree_maxEcrystPosClust)->SetTitle("Time of maximum energy crystal in cluster with +ve charge, GeV");
158  m_dbgTree_evt_allCuts->Branch("maxEcrysTimeNegClust",
159  &m_tree_maxEcrystNegClust)->SetTitle("Time of maximum energy crystal in cluster with -ve charge, GeV");
160  m_dbgTree_evt_allCuts->Branch("t0", &m_tree_t0)->SetTitle("T0, ns");
161  m_dbgTree_evt_allCuts->Branch("t0_ECL_closestCDC", &m_tree_t0_ECLclosestCDC)->SetTitle("T0 ECL closest to CDC t0, ns");
162  m_dbgTree_evt_allCuts->Branch("t0_ECL_minChi2", &m_tree_t0_ECL_minChi2)->SetTitle("T0 ECL with smallest chi squared, ns");
163 
164  m_dbgTree_evt_allCuts->SetAutoSave(10);
165 
166 
167  // Per crystal within each cluster, entry after all the cuts
168  m_dbgTree_crys_allCuts = new TTree("m_dbgTree_crys_allCuts",
169  "Debug data for bhabha time calibration - one entry per crystal per cluster entry after all cuts");
170 
171  m_dbgTree_crys_allCuts->Branch("runNum", &m_runNum)->SetTitle("Run number");
172  m_dbgTree_crys_allCuts->Branch("EventNum", &m_tree_evtNum)->SetTitle("Event number");
173  m_dbgTree_crys_allCuts->Branch("m_tree_ECLCalDigitTime",
175  ->SetTitle("Time of a crystal within the cluster after application of previous calibrations except t0, ns");
176  m_dbgTree_crys_allCuts->Branch("m_tree_ECLCalDigitE", &m_tree_ECLCalDigitE)->SetTitle("Energy of crystal, GeV");
177  m_dbgTree_crys_allCuts->Branch("m_tree_ECLDigitAmplitude",
178  &m_tree_ECLDigitAmplitude)->SetTitle("Amplitude of crystal signal pulse");
179  m_dbgTree_crys_allCuts->Branch("timetsPreviousTimeCalibs",
180  &m_tree_timetsPreviousTimeCalibs)->SetTitle("Time t_psi after application of previous Ts, ns");
181  m_dbgTree_crys_allCuts->Branch("t0", &m_tree_t0)->SetTitle("T0, ns");
182  m_dbgTree_crys_allCuts->Branch("t0_ECL_closestCDC", &m_tree_t0_ECLclosestCDC)->SetTitle("T0 ECL closest to CDC t0, ns");
183  m_dbgTree_crys_allCuts->Branch("t0_ECL_minChi2", &m_tree_t0_ECL_minChi2)->SetTitle("T0 ECL with smallest chi squared, ns");
184  m_dbgTree_crys_allCuts->Branch("CrystalCellID", &m_tree_cid)->SetTitle("Cell ID, 1..8736");
185 
186  m_dbgTree_crys_allCuts->SetAutoSave(10);
187 
188 
189  } // <<<<<<<<<<<<<<<<<<<<< if boolean true for saving debug trees <<<<<<<<<<<<<<<<<<<<<<<<<< */
190 
191 
192  // Per max E crystal entry after all the cuts
193  // this tree is always saved
194  m_dbgTree_allCuts = new TTree("tree_allCuts",
195  "Debug data for bhabha time calibration - one entry per max E crystal entry after cuts");
196 
197  m_dbgTree_allCuts->Branch("time_t_psi", &m_tree_time)->SetTitle("Time t_psi for Ts, ns");
198  m_dbgTree_allCuts->Branch("crateID", &m_crystalCrate)->SetTitle("Crate id for crystal");
199  m_dbgTree_allCuts->Branch("EventNum", &m_tree_evtNum)->SetTitle("Event number");
200  m_dbgTree_allCuts->Branch("runNum", &m_runNum)->SetTitle("Run number");
201  m_dbgTree_allCuts->Branch("CrystalCellID", &m_tree_cid)->SetTitle("Cell ID, 1..8736");
202  m_dbgTree_allCuts->Branch("maxEcrystalEnergy", &m_tree_en)->SetTitle("Max Energy Crystal Energy, GeV");
203  m_dbgTree_allCuts->Branch("maxEcrystalEnergyDIVclusterE",
204  &m_tree_E1Etot)->SetTitle("Max Energy Crystal Fraction Energy of Cluster");
205  m_dbgTree_allCuts->Branch("E1divE2crystal",
206  &m_tree_E1E2)->SetTitle("Max Energy Crystal DIV second max energy crystal in cluster");
207  m_dbgTree_allCuts->Branch("E1crystalDIVp", &m_tree_E1p)->SetTitle("Max Energy Crystal in cluster DIV track p");
208  m_dbgTree_allCuts->Branch("timetsPreviousTimeCalibs",
209  &m_tree_timetsPreviousTimeCalibs)->SetTitle("Time t_psi after application of previous Ts, ns");
210  m_dbgTree_allCuts->Branch("massInvTracks", &m_massInvTracks)->SetTitle("Invariant mass of the two tracks");
211  m_dbgTree_allCuts->Branch("t0", &m_tree_t0)->SetTitle("T0, ns");
212  m_dbgTree_allCuts->Branch("t0_ECL_closestCDC", &m_tree_t0_ECLclosestCDC)->SetTitle("T0 ECL closest to CDC t0, ns");
213  m_dbgTree_allCuts->Branch("t0_ECL_minChi2", &m_tree_t0_ECL_minChi2)->SetTitle("T0 ECL with smallest chi squared, ns");
214 
215  m_dbgTree_allCuts->Branch("clusterTime", &m_tree_tClust)->SetTitle("Cluster time of cluster with +ve charge, GeV");
216 
217  m_dbgTree_allCuts->SetAutoSave(10);
218 
219 }
220 
222 {
223  //=== MetaData
224  B2INFO("ECLBhabhaTCollector: Experiment = " << m_EventMetaData->getExperiment() <<
225  " run = " << m_EventMetaData->getRun());
226 
227 
228  //=== Create histograms and register them in the data store
229  int nbins = m_timeAbsMax * 8;
230  int max_t = m_timeAbsMax;
231  int min_t = -m_timeAbsMax;
232 
233 
234  auto TimevsCrysPrevCrateCalibPrevCrystCalib = new TH2F("TimevsCrysPrevCrateCalibPrevCrystCalib",
235  "Time t psi - ts - tcrate (previous calibs) vs crystal cell ID;crystal cell ID;Time t_psi with previous calib (ns)",
237  registerObject<TH2F>("TimevsCrysPrevCrateCalibPrevCrystCalib", TimevsCrysPrevCrateCalibPrevCrystCalib);
238 
239  auto TimevsCratePrevCrateCalibPrevCrystCalib = new TH2F("TimevsCratePrevCrateCalibPrevCrystCalib",
240  "Time t psi - ts - tcrate (previous calibs) vs crate ID;crate ID;Time t_psi previous calib (ns)",
241  52, 1, 52 + 1, nbins, min_t, max_t);
242  registerObject<TH2F>("TimevsCratePrevCrateCalibPrevCrystCalib", TimevsCratePrevCrateCalibPrevCrystCalib);
243 
244  auto TimevsCrysNoCalibrations = new TH2F("TimevsCrysNoCalibrations",
245  "Time tpsi vs crystal cell ID;crystal cell ID;Time t_psi (ns)", ECLElementNumbers::c_NCrystals, 1,
246  ECLElementNumbers::c_NCrystals + 1, nbins, min_t, max_t);
247  registerObject<TH2F>("TimevsCrysNoCalibrations", TimevsCrysNoCalibrations);
248 
249  auto TimevsCrateNoCalibrations = new TH2F("TimevsCrateNoCalibrations",
250  "Time tpsi vs crate ID;crate ID;Time t_psi (ns)", 52, 1, 52 + 1, nbins, min_t, max_t);
251  registerObject<TH2F>("TimevsCrateNoCalibrations", TimevsCrateNoCalibrations);
252 
253  auto TimevsCrysPrevCrateCalibNoCrystCalib = new TH2F("TimevsCrysPrevCrateCalibNoCrystCalib",
254  "Time tpsi - tcrate (previous calib) vs crystal cell ID;crystal cell ID;Time t_psi including previous crate calib (ns)",
256  registerObject<TH2F>("TimevsCrysPrevCrateCalibNoCrystCalib", TimevsCrysPrevCrateCalibNoCrystCalib);
257 
258  auto TimevsCrateNoCrateCalibPrevCrystCalib = new TH2F("TimevsCrateNoCrateCalibPrevCrystCalib",
259  "Time tpsi - ts (previous calib) vs crate ID;crate ID;Time t_psi including previous crystal calib (ns)",
260  52, 1, 52 + 1, nbins, min_t, max_t);
261  registerObject<TH2F>("TimevsCrateNoCrateCalibPrevCrystCalib", TimevsCrateNoCrateCalibPrevCrystCalib);
262 
263 
264  auto TsDatabase = new TH1F("TsDatabase", ";cell id;Ts from database", ECLElementNumbers::c_NCrystals, 1,
266  registerObject<TH1F>("TsDatabase", TsDatabase);
267 
268  auto TsDatabaseUnc = new TH1F("TsDatabaseUnc", ";cell id;Ts uncertainty from database", ECLElementNumbers::c_NCrystals, 1,
270  registerObject<TH1F>("TsDatabaseUnc", TsDatabaseUnc);
271 
272  auto TcrateDatabase = new TH1F("TcrateDatabase", ";cell id;Tcrate from database", ECLElementNumbers::c_NCrystals, 1,
274  registerObject<TH1F>("TcrateDatabase", TcrateDatabase);
275 
276  auto TcrateUncDatabase = new TH1F("TcrateUncDatabase", ";cell id;Tcrate uncertainty from database", ECLElementNumbers::c_NCrystals,
278  registerObject<TH1F>("TcrateUncDatabase", TcrateUncDatabase);
279 
280 
281  auto tcrateDatabase_ns = new TH1F("tcrateDatabase_ns", ";crate id;tcrate derived from database", 52, 1, 52 + 1);
282  registerObject<TH1F>("tcrateDatabase_ns", tcrateDatabase_ns);
283 
284 
285  auto databaseCounter = new TH1I("databaseCounter",
286  ";A database was read in;Number of times database was saved to histogram", 1, 1, 2);
287  registerObject<TH1I>("databaseCounter", databaseCounter);
288 
289 
290  auto numCrystalEntriesPerEvent = new TH1F("numCrystalEntriesPerEvent",
291  ";Number crystal entries;Number of events", 15, 0, 15);
292  registerObject<TH1F>("numCrystalEntriesPerEvent", numCrystalEntriesPerEvent);
293 
294  auto cutflow = new TH1F("cutflow", ";Cut label number;Number of events passing cut", 20, 0, 20);
295  registerObject<TH1F>("cutflow", cutflow);
296 
297  auto maxEcrsytalEnergyFraction = new TH1F("maxEcrsytalEnergyFraction",
298  ";Maximum energy crystal energy / (sum) cluster energy;Number", 22, 0, 1.1);
299  registerObject<TH1F>("maxEcrsytalEnergyFraction", maxEcrsytalEnergyFraction);
300 
301  auto refCrysIDzeroingCrate = new TH1F("refCrysIDzeroingCrate", ";cell id;Boolean - is reference crystal",
303  registerObject<TH1F>("refCrysIDzeroingCrate", refCrysIDzeroingCrate);
304 
305  auto CDCEventT0Correction = new TH1F("CDCEventT0Correction", ";;CDC event t0 offset correction [ns]", 1, 1, 2);
306  registerObject<TH1F>("CDCEventT0Correction", CDCEventT0Correction);
307 
308 
309  //=== Required data objects
310  m_eventT0.isRequired();
311  tracks.isRequired();
312  m_eclClusterArray.isRequired();
313  m_eclCalDigitArray.isRequired();
314  m_eclDigitArray.isRequired();
315 
316  B2INFO("hadronEventT0_TO_bhabhaEventT0_correction = " << m_hadronEventT0_TO_bhabhaEventT0_correction <<
317  " ns correction to CDC event t0 will be applied");
318 
319  B2INFO("skipTrgSel = " << skipTrgSel);
320 
321 }
322 
324 {
325  int cutIndexPassed = 0;
326  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
327  B2DEBUG(22, "Cutflow: no cuts: index = " << cutIndexPassed);
328  B2DEBUG(22, "Event number = " << m_EventMetaData->getEvent());
329 
330 
331  // --- Check the trigger skim is the type that has two tracks
332 
333  /* If we skip the trigger skim selection then still fill the cutflow histogram
334  just so that the positions don't change. */
335  if (!skipTrgSel) {
336  if (!m_TrgResult.isValid()) {
337  B2WARNING("SoftwareTriggerResult required to select bhabha event is not found");
338  return;
339  }
340 
341  /* Release05: bhabha_all is grand skim = bhabha+bhabhaecl+radee. We only want
342  to look at the 2 track bhabha events. */
343  const std::map<std::string, int>& fresults = m_TrgResult->getResults();
344  if (fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end()) {
345  B2WARNING("Can't find required bhabha trigger identifier");
346  return;
347  }
348 
349  const bool eBhabha = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
351  B2DEBUG(22, "eBhabha (trigger passed) = " << eBhabha);
352 
353  if (!eBhabha) {
354  return;
355  }
356  }
357 
358  /* Fill the histgram showing that the trigger skim cut passed OR that we
359  are skipping this selection. */
360  cutIndexPassed++;
361  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
362  B2DEBUG(22, "Cutflow: Trigger cut passed: index = " << cutIndexPassed);
363 
364 
365 
366  /* Use ECLChannelMapper to get other detector indices for the crystals
367  For conversion from CellID to crate, shaper, and channel ids.
368  The initialization function automatically checks to see if the
369  object has been initialized and ifthe payload has changed and
370  thus needs updating. */
371  bool ECLchannelMapHasChanged = m_channelMapDB.hasChanged();
372  if (ECLchannelMapHasChanged) {
373  B2INFO("ECLBhabhaTCollectorModule::collect() " << LogVar("ECLchannelMapHasChanged", ECLchannelMapHasChanged));
374  if (!m_crystalMapper->initFromDB()) {
375  B2FATAL("ECLBhabhaTCollectorModule::collect() : Can't initialize eclChannelMapper!");
376  }
377  }
378 
379 
380  //== Get expected energies and calibration constants from DB. Need to call
381  // hasChanged() for later comparison
382  if (m_ElectronicsDB.hasChanged()) {
383  m_Electronics = m_ElectronicsDB->getCalibVector();
384  }
385  if (m_ElectronicsTimeDB.hasChanged()) {
386  m_ElectronicsTime = m_ElectronicsTimeDB->getCalibVector();
387  }
388  if (m_FlightTimeDB.hasChanged()) {
389  m_FlightTime = m_FlightTimeDB->getCalibVector();
390  }
391 
392  // Get the previous crystal time offset (the same thing that this calibration is meant to calculate).
393  // This can be used for testing purposes, and for the crate time offset.
394  if (m_PreviousCrystalTimeDB.hasChanged()) {
396  m_PreviousCrystalTimeUnc = m_PreviousCrystalTimeDB->getCalibUncVector();
397  }
398 
399  B2DEBUG(29, "Finished checking if previous crystal time payload has changed");
400  if (m_CrateTimeDB.hasChanged()) {
401  m_CrateTime = m_CrateTimeDB->getCalibVector();
402  m_CrateTimeUnc = m_CrateTimeDB->getCalibUncVector();
403  }
404  B2DEBUG(29, "Finished checking if previous crate time payload has changed");
405  B2DEBUG(29, "m_CrateTime size = " << m_CrateTime.size());
406  B2DEBUG(25, "Crate time +- uncertainty [0]= " << m_CrateTime[0] << " +- " << m_CrateTimeUnc[0]);
407  B2DEBUG(25, "Crate time +- uncertainty [8735]= " << m_CrateTime[8735] << " +- " << m_CrateTimeUnc[8735]);
408 
409  B2DEBUG(29, "Finished checking if previous crate time payload has changed");
410  if (m_RefCrystalsCalibDB.hasChanged()) {
411  m_RefCrystalsCalib = m_RefCrystalsCalibDB->getReferenceCrystals();
412  }
413  B2DEBUG(29, "Finished checking if reference crystal cell ids payload has changed");
414 
415 
416  B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrystalTimeOffset from the database"
417  << LogVar("IoV", m_PreviousCrystalTimeDB.getIoV())
418  << LogVar("Checksum", m_PreviousCrystalTimeDB.getChecksum()));
419  B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrateTimeOffset from the database"
420  << LogVar("IoV", m_CrateTimeDB.getIoV())
421  << LogVar("Checksum", m_CrateTimeDB.getChecksum()));
422  B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrystalElectronics from the database"
423  << LogVar("IoV", m_ElectronicsDB.getIoV())
424  << LogVar("Checksum", m_ElectronicsDB.getChecksum()));
425  B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrystalElectronicsTime from the database"
426  << LogVar("IoV", m_ElectronicsTimeDB.getIoV())
427  << LogVar("Checksum", m_ElectronicsTimeDB.getChecksum()));
428  B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrystalFlightTime from the database"
429  << LogVar("IoV", m_FlightTimeDB.getIoV())
430  << LogVar("Checksum", m_FlightTimeDB.getChecksum()));
431  B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLReferenceCrystalPerCrateCalib from the database"
432  << LogVar("IoV", m_RefCrystalsCalibDB.getIoV())
433  << LogVar("Checksum", m_RefCrystalsCalibDB.getChecksum()));
434 
435 
436 
437  // Conversion coefficient from ADC ticks to nanoseconds
438  // TICKS_TO_NS ~ 0.4913 ns/clock tick
439  // 1/(4fRF) = 0.4913 ns/clock tick, where fRF is the accelerator RF frequency
440  const double TICKS_TO_NS = 1.0 / (4.0 * EclConfiguration::getRF()) * 1e3;
441 
442 
443  vector<float> Crate_time_ns(52, 0.0);
445  // Make a crate time offset vector with an entry per crate (instead of per crystal) and convert from ADC counts to ns.
446  for (int crysID = 1; crysID <= ECLElementNumbers::c_NCrystals; crysID++) {
447  int crateID_temp = m_crystalMapper->getCrateID(crysID);
448  Crate_time_ns[crateID_temp - 1] = m_CrateTime[crysID] * TICKS_TO_NS;
449  }
450 
451 
452 
458  for (int crateID_temp = 1; crateID_temp <= 52; crateID_temp++) {
459  getObjectPtr<TH1F>("refCrysIDzeroingCrate")->Fill(m_RefCrystalsCalib[crateID_temp - 1] + 0.001);
460  }
461 
462 
464  for (int crysID = 1; crysID <= ECLElementNumbers::c_NCrystals; crysID++) {
465  getObjectPtr<TH1F>("TsDatabase")->SetBinContent(crysID + 0.001, m_PreviousCrystalTime[crysID - 1]);
466  getObjectPtr<TH1F>("TsDatabaseUnc")->SetBinContent(crysID + 0.001, m_PreviousCrystalTimeUnc[crysID - 1]);
467  getObjectPtr<TH1F>("TcrateDatabase")->SetBinContent(crysID + 0.001, m_CrateTime[crysID - 1]);
468  getObjectPtr<TH1F>("TcrateUncDatabase")->SetBinContent(crysID + 0.001, m_CrateTimeUnc[crysID - 1]);
469  }
470  if (m_storeCalib) {
471  B2INFO("ECLBhabhaTCollector:: ECLCrystalTimeOffset from the database information:"
472  << LogVar("IoV", m_PreviousCrystalTimeDB.getIoV())
473  << LogVar("Checksum", m_PreviousCrystalTimeDB.getChecksum()));
474  B2INFO("First event so print out previous ts values");
475  for (int crysID = 1; crysID <= ECLElementNumbers::c_NCrystals; crysID++) {
476  B2INFO("cid = " << crysID << ", Ts previous = " << m_PreviousCrystalTime[crysID - 1]);
477  }
478  m_storeCalib = false;
479  }
480 
481 
482 
483 
484  for (int crateID_temp = 1; crateID_temp <= 52; crateID_temp++) {
485  getObjectPtr<TH1F>("tcrateDatabase_ns")->SetBinContent(crateID_temp + 0.001, Crate_time_ns[crateID_temp - 1]);
486  }
487 
488  // Use a histogram with only one bin as a counter to know the number of times the database histograms were filled.
489  // This is mostly useful for the talg when running over multiple runs and trying to read ts values.
490  getObjectPtr<TH1I>("databaseCounter")->SetBinContent(1, 1);
491 
492 
493 
494  // Save what CDC event t0 correction was applied
495  getObjectPtr<TH1F>("CDCEventT0Correction")->SetBinContent(1, m_hadronEventT0_TO_bhabhaEventT0_correction);
496 
497 
498 
499 
500  /* Getting the event t0 using the full event t0 rather than from the CDC specifically */
501  double evt_t0 = -1;
502  double evt_t0_unc = -1;
503  double evt_t0_ECL_closestCDC = -1;
504  double evt_t0_ECL_minChi2 = -1;
505 
506  // Determine if there is an event t0 to use and then extract the information about it
507  if (!m_eventT0.isValid()) {
508  //cout << "event t0 not valid\n";
509  return;
510  } else if (!m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
511  //cout << "no event t0\n";
512  return;
513  } else {
514  // Event has a t0 from CDC
515  cutIndexPassed++;
516  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
517  B2DEBUG(22, "Cutflow: Event t0 exists: index = " << cutIndexPassed);
518 
519 
520  // Get event t0 from CDC. We don't want event t0 from ECL as we are calibrating the ECL wrt the more accurately measured time measurements of the time. Start with the CDC since it has an event t0 but in the future we may switch to the TOP detector.
521  // Based on the information from Thomas Hauth <Thomas.Hauth@kit.edu> (leaving physics) we should take the last event t0 in the list of event t0's from the CDC as the later event t0 measurements are calculated in slower but more accurate ways.
522  vector<EventT0::EventT0Component> evt_t0_list = m_eventT0->getTemporaryEventT0s(Const::EDetector::CDC);
523  evt_t0 = evt_t0_list.back().eventT0; // time value
524  evt_t0_unc = evt_t0_list.back().eventT0Uncertainty; // uncertainty on event t0
525 
526 
527  // Correct the CDC event t0 value for the bhabha bias
528  evt_t0 = evt_t0 + m_hadronEventT0_TO_bhabhaEventT0_correction; // Bias not yet fixed in CDC t0 reco.
529 
530 
531  // Get the ECL event t0 for comparison - validations
532  if (m_eventT0->hasTemporaryEventT0(Const::EDetector::ECL)) {
533  vector<EventT0::EventT0Component> evt_t0_list_ECL = m_eventT0->getTemporaryEventT0s(Const::EDetector::ECL);
534 
535 
536  double smallest_CDC_ECL_t0_diff = fabs(evt_t0_list_ECL[0].eventT0 - evt_t0);
537  int smallest_CDC_ECL_t0_diff_idx = 0;
538  for (long unsigned int ECLi = 0; ECLi < evt_t0_list_ECL.size(); ECLi++) {
539  double tempt_ECL_t0 = evt_t0_list_ECL[ECLi].eventT0;
540  if (fabs(tempt_ECL_t0 - evt_t0) < smallest_CDC_ECL_t0_diff) {
541  smallest_CDC_ECL_t0_diff = fabs(tempt_ECL_t0 - evt_t0);
542  smallest_CDC_ECL_t0_diff_idx = ECLi;
543  }
544  }
545 
546  evt_t0_ECL_closestCDC = evt_t0_list_ECL[smallest_CDC_ECL_t0_diff_idx].eventT0; // time value
547  B2DEBUG(26, "evt_t0_ECL_closestCDC = " << evt_t0_ECL_closestCDC);
548 
549 
550 
551  double smallest_ECL_t0_minChi2 = evt_t0_list_ECL[0].quality;
552  int smallest_ECL_t0_minChi2_idx = 0;
553 
554  B2DEBUG(26, "evt_t0_list_ECL[0].quality = " << evt_t0_list_ECL[0].quality
555  << ", with ECL event t0 = " << evt_t0_list_ECL[0].eventT0);
556 
557  for (long unsigned int ECLi = 0; ECLi < evt_t0_list_ECL.size(); ECLi++) {
558  B2DEBUG(26, "evt_t0_list_ECL[" << ECLi << "].quality = " << evt_t0_list_ECL[ECLi].quality
559  << ", with ECL event t0 = " <<
560  evt_t0_list_ECL[ECLi].eventT0);
561  if (evt_t0_list_ECL[ECLi].quality < smallest_ECL_t0_minChi2) {
562  smallest_ECL_t0_minChi2 = evt_t0_list_ECL[ECLi].quality;
563  smallest_ECL_t0_minChi2_idx = ECLi;
564  }
565  }
566 
567  evt_t0_ECL_minChi2 = evt_t0_list_ECL[smallest_ECL_t0_minChi2_idx].eventT0; // time value
568 
569  B2DEBUG(26, "evt_t0_ECL_minChi2 = " << evt_t0_ECL_minChi2);
570  B2DEBUG(26, "smallest_ECL_t0_minChi2_idx = " << smallest_ECL_t0_minChi2_idx);
571  }
572  }
573 
574 
575 
576  /* Determine the energies for each of the crystals since this isn't naturally connected to the cluster.
577  Also determine the indexing of the ecl cal digits and the ecl digits
578  Taken from Chris's ec/modules/eclGammaGammaECollector */
579 
580  // Resize vectors
584 
585 
586  int idx = 0;
587  for (auto& eclCalDigit : m_eclCalDigitArray) {
588  int tempCrysID = eclCalDigit.getCellId() - 1;
589  m_EperCrys[tempCrysID] = eclCalDigit.getEnergy();
590  m_eclCalDigitID[tempCrysID] = idx;
591  idx++;
592  }
593 
594  idx = 0;
595  for (auto& eclDigit : m_eclDigitArray) {
596  int tempCrysID = eclDigit.getCellId() - 1;
597  m_eclDigitID[tempCrysID] = idx;
598  idx++;
599  }
600 
601 
602 
603 
604  //---------------------------------------------------------------------
605  //..Some utilities
606  PCmsLabTransform boostrotate;
607 
608  //---------------------------------------------------------------------
609  //..Track properties, including 2 maxp tracks. Use pion (211) mass hypothesis,
610  // which is the only particle hypothesis currently available???
611  double maxp[2] = {0., 0.};
612  int maxiTrk[2] = { -1, -1};
613  int nTrkAll = tracks.getEntries();
614 
615  int nTrkLoose = 0;
616  int nTrkTight = 0;
618  /* Loop over all the tracks to define the tight and loose selection tracks.
619  We will select events with only 2 tight tracks and no additional loose tracks.
620  Tight tracks are a subset of looses tracks. */
621  for (int iTrk = 0; iTrk < nTrkAll; iTrk++) {
622  // Get track biasing towards the particle being a pion based on what particle types
623  // are used for reconstruction at this stage.
624  const TrackFitResult* tempTrackFit = tracks[iTrk]->getTrackFitResultWithClosestMass(Const::pion);
625  if (not tempTrackFit) {continue;}
626 
627  // Collect track info to be used for categorizing
628  short charge = tempTrackFit->getChargeSign();
629  double z0 = tempTrackFit->getZ0();
630  double d0 = tempTrackFit->getD0();
631  int nCDChits = tempTrackFit->getHitPatternCDC().getNHits();
632  double p = tempTrackFit->getMomentum().R();
633 
634  // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
635  //== Save debug TTree with detailed information if necessary.
636  m_tree_d0 = d0;
637  m_tree_z0 = z0;
638  m_tree_p = p;
639  m_charge = charge;
640  m_tree_nCDChits = nCDChits;
641 
642  if (m_saveTree) {
643  m_dbgTree_tracks->Fill();
644  }
645  //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
646 
647 
648  /* Test if loose track */
649 
650  // d0 and z0 cuts
651  if (fabs(d0) > m_looseTrkD0) {
652  continue;
653  }
654  if (fabs(z0) > m_looseTrkZ0) {
655  continue;
656  }
657  // Number of hits in the CDC
658  if (nCDChits < 1) {
659  continue;
660  }
661  nTrkLoose++;
662 
663 
664 
665  /* Test if the loose track is also a tight track */
666 
667  // Number of hits in the CDC
668  if (nCDChits < 20) {
669  continue;
670  }
671 
672 
673  // d0 and z0 cuts
674  if (fabs(d0) > m_tightTrkD0) {
675  continue;
676  }
677  if (fabs(z0) > m_tightTrkZ0) {
678  continue;
679  }
680  nTrkTight++;
681 
682  // Sorting of tight tracks. Not really required as we only want two tight tracks (at the moment) but okay.
683  //..Find the maximum p negative [0] and positive [1] tracks
684  int icharge = 0;
685  if (charge > 0) {icharge = 1;}
686  if (p > maxp[icharge]) {
687  maxp[icharge] = p;
688  maxiTrk[icharge] = iTrk;
689  }
690 
691  }
692  /* After that last section the numbers of loose and tight tracks are known as well as the
693  index of the loose tracks that have the highest p negatively charged and highest p positively
694  charged tracks as measured in the centre of mass frame */
695 
696 
697  if (nTrkTight != 2) {
698  return;
699  }
700  // There are exactly two tight tracks
701  cutIndexPassed++;
702  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
703  B2DEBUG(22, "Cutflow: Two tight tracks: index = " << cutIndexPassed);
704 
705 
706  if (nTrkLoose != 2) {
707  return;
708  }
709  // There are exactly two loose tracks as well, i.e. no additional loose tracks
710  cutIndexPassed++;
711  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
712  B2DEBUG(22, "Cutflow: No additional loose tracks: index = " << cutIndexPassed);
713 
714  /* Determine if the two tracks have the opposite electric charge.
715  We know this because the track indices stores the max pt track in [0] for negatively charged track
716  and [1] fo the positively charged track. If both are filled then both a negatively charged
717  and positively charged track were found. */
718  bool oppositelyChargedTracksPassed = maxiTrk[0] != -1 && maxiTrk[1] != -1;
719  if (!oppositelyChargedTracksPassed) {
720  return;
721  }
722  // The two tracks have the opposite electric charges.
723  cutIndexPassed++;
724  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
725  B2DEBUG(22, "Cutflow: Oppositely charged tracks: index = " << cutIndexPassed);
726 
727 
728 
729 
730  //---------------------------------------------------------------------
731  /* Determine associated energy clusters to each of the two tracks. Sum the energies of the
732  multiple clusters to each track and find the crystal with the maximum energy within all
733  the sets of clusters associated to the tracks*/
734  double trkEClustLab[2] = {0., 0.};
735  double trkEClustCOM[2] = {0., 0.};
736  double trkpLab[2];
737  double trkpCOM[2];
738  ROOT::Math::PxPyPzEVector trkp4Lab[2];
739  ROOT::Math::PxPyPzEVector trkp4COM[2];
740 
741  // Index of the cluster and the crystal that has the highest energy crystal for the two tracks
742  int crysIDMax[2] = { -1, -1 };
743  double crysEMax[2] = { -1, -1 };
744  double crysE2Max[2] = { -1, -1 };
745  int numClustersPerTrack[2] = { 0, 0 };
746 
747  double clusterTime[2] = {0, 0};
748 
749  double E_DIV_p[2];
750 
751  vector<double> time_ECLCaldigits_bothClusters;
752  vector<int> cid_ECLCaldigits_bothClusters;
753  vector<double> E_ECLCaldigits_bothClusters;
754  vector<double> amp_ECLDigits_bothClusters;
755  vector<int> chargeID_ECLCaldigits_bothClusters;
756 
757  for (int icharge = 0; icharge < 2; icharge++) {
758  if (maxiTrk[icharge] > -1) {
759  B2DEBUG(22, "looping over the 2 max pt tracks");
760 
761  const TrackFitResult* tempTrackFit = tracks[maxiTrk[icharge]]->getTrackFitResultWithClosestMass(Const::pion);
762  if (not tempTrackFit) {continue;}
763  trkp4Lab[icharge] = tempTrackFit->get4Momentum();
764  trkp4COM[icharge] = boostrotate.rotateLabToCms() * trkp4Lab[icharge];
765  trkpLab[icharge] = trkp4Lab[icharge].P();
766  trkpCOM[icharge] = trkp4COM[icharge].P();
767 
768 
769  /* For each cluster associated to the current track, sum up the energies to get the total
770  energy of all clusters associated to the track and find which crystal has the highest
771  energy from all those clusters*/
772  auto eclClusterRelationsFromTracks = tracks[maxiTrk[icharge]]->getRelationsTo<ECLCluster>();
773  for (unsigned int clusterIdx = 0; clusterIdx < eclClusterRelationsFromTracks.size(); clusterIdx++) {
774 
775  B2DEBUG(22, "Looking at clusters. index = " << clusterIdx);
776  auto cluster = eclClusterRelationsFromTracks[clusterIdx];
777  bool goodClusterType = false;
778 
779  if (cluster->hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons)) {
780  trkEClustLab[icharge] += cluster->getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons);
781  goodClusterType = true;
782  numClustersPerTrack[icharge]++;
783  }
784 
785  if (goodClusterType) {
786 
787  clusterTime[icharge] = cluster->getTime();
788 
789  auto eclClusterRelations = cluster->getRelationsTo<ECLCalDigit>("ECLCalDigits");
790 
791  // Find the crystal that has the largest energy
792  for (unsigned int ir = 0; ir < eclClusterRelations.size(); ir++) {
793  const auto calDigit = eclClusterRelations.object(ir);
794  int tempCrysID = calDigit->getCellId() - 1;
795  double tempE = m_EperCrys[tempCrysID];
796 
797  int eclDigitIndex = m_eclDigitID[tempCrysID];
798  ECLDigit* ecl_dig = m_eclDigitArray[eclDigitIndex];
799 
800  // for the max E crystal
801  if (tempE > crysEMax[icharge]) {
802  // Set 2nd highest E crystal to the info from the highest E crystal
803  crysE2Max[icharge] = crysEMax[icharge];
804  // Set the highest E crystal to the current crystal
805  crysEMax[icharge] = tempE;
806  crysIDMax[icharge] = tempCrysID;
807  }
808  // for the 2nd highest E crystal
809  if (tempE > crysE2Max[icharge] && tempCrysID != crysIDMax[icharge]) {
810  crysE2Max[icharge] = tempE;
811  }
812 
813  // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
814  // If we drop the information about the second highest energy crystal we could use
815  // m_eclClusterArray[ic]->getMaxECellId()
816 
817  B2DEBUG(26, "calDigit(ir" << ir << ") time = " << calDigit->getTime() << "ns , with E = " << tempE << " GeV");
818  time_ECLCaldigits_bothClusters.push_back(calDigit->getTime());
819  cid_ECLCaldigits_bothClusters.push_back(tempCrysID);
820  E_ECLCaldigits_bothClusters.push_back(tempE);
821  amp_ECLDigits_bothClusters.push_back(ecl_dig->getAmp());
822  chargeID_ECLCaldigits_bothClusters.push_back(icharge);
823 
824  }
825  }
826  }
827  trkEClustCOM[icharge] = trkEClustLab[icharge] * trkpCOM[icharge] / trkpLab[icharge];
828 
829  // Check both electrons to see if their cluster energy / track momentum is good.
830  // The Belle II physics book shows that this is the main way of separating electrons from other particles
831  // Done in the centre of mass reference frame although I believe E/p is invariant under a boost.
832  E_DIV_p[icharge] = trkEClustCOM[icharge] / trkpCOM[icharge];
833 
834  }
835  }
836  /* At the end of this section the 3-momenta magnitudes and the cluster energies are known
837  for the two saved track indices for both the lab and COM frames.
838  The crystal with the maximum energy, one associated to each track, is recorded*/
839 
840 
841 
842  //=== Check each crystal in the processed event and fill histogram.
843 
844  int numCrystalsPassingCuts = 0;
845 
846  int crystalIDs[2] = { -1, -1};
847  int crateIDs[2] = { -1, -1};
848  double ts_prevCalib[2] = { -1, -1};
849  double tcrate_prevCalib[2] = { -1, -1};
850  double times_noPreviousCalibrations[2] = { -1, -1};
851  bool crystalCutsPassed[2] = {false, false};
852  double crystalEnergies[2] = { -1, -1};
853  double crystalEnergies2[2] = { -1, -1};
854 
855  for (int iCharge = 0; iCharge < 2; iCharge++) {
856  int crystal_idx = crysIDMax[iCharge];
857  int eclCalDigitIndex = m_eclCalDigitID[crystal_idx];
858  int eclDigitIndex = m_eclDigitID[crystal_idx];
859 
860  ECLDigit* ecl_dig = m_eclDigitArray[eclDigitIndex];
861  ECLCalDigit* ecl_cal = m_eclCalDigitArray[eclCalDigitIndex];
862 
863  //== Check whether specific ECLDigits should be excluded.
864 
865  auto en = ecl_cal->getEnergy();
866  auto amplitude = ecl_dig->getAmp();
867  crystalEnergies[iCharge] = en;
868 
869  int cid = ecl_dig->getCellId();
870  double time = ecl_dig->getTimeFit() * TICKS_TO_NS - evt_t0;
871 
872  // Offset time by electronics calibration and flight time calibration.
873  time -= m_ElectronicsTime[cid - 1] * TICKS_TO_NS;
874  time -= m_FlightTime[cid - 1];
875 
876 
877  // Apply the time walk correction: time shift as a function of the amplitude corrected by the electronics calibration.
878  // The electronics calibration also accounts for crystals that have a dead pre-amp and thus half the normal amplitude.
879  double energyTimeShift = m_ECLTimeUtil->energyDependentTimeOffsetElectronic(amplitude * m_Electronics[cid - 1]) * TICKS_TO_NS;
880 
881  B2DEBUG(29, "cellid = " << cid << ", amplitude = " << amplitude << ", time before t(E) shift = " << time <<
882  ", t(E) shift = " << energyTimeShift << " ns");
883  time -= energyTimeShift;
884 
885 
886  // Cell ID should be within specified range.
887  if (cid < m_minCrystal || cid > m_maxCrystal) continue;
888 
889  // Absolute time should be in specified range condition.
890  if (fabs(time) > m_timeAbsMax) continue;
891 
892  // Fit quality flag -- choose only events with best fit quality
893  if (ecl_dig->getQuality() != 0) continue;
894 
895  //== Save time and crystal information. Fill plot after both electrons are tested
896  crystalIDs[iCharge] = cid;
897  crateIDs[iCharge] = m_crystalMapper->getCrateID(ecl_cal->getCellId());
898 
899 
900  ts_prevCalib[iCharge] = m_PreviousCrystalTime[cid - 1] * TICKS_TO_NS;
901  tcrate_prevCalib[iCharge] = m_CrateTime[cid - 1] * TICKS_TO_NS;
902  times_noPreviousCalibrations[iCharge] = time;
903 
904 
905  B2DEBUG(26, "iCharge = " << iCharge);
906  B2DEBUG(26, "crateIDs[iCharge] = " << crateIDs[iCharge]);
907  B2DEBUG(26, "times_noPreviousCalibrations[iCharge] = " << times_noPreviousCalibrations[iCharge]);
908  B2DEBUG(26, "tcrate_prevCalib[iCharge] = " << tcrate_prevCalib[iCharge]);
909  B2DEBUG(26, "ts_prevCalib[iCharge] = " << ts_prevCalib[iCharge]);
910 
911 
912  crystalCutsPassed[iCharge] = true;
913 
914 
915  // For second most energetic energy crystal
916  crystalEnergies2[iCharge] = crysE2Max[iCharge];
917 
918 
919 // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
920  //== Save debug TTree with detailed information if necessary.
921  m_tree_cid = ecl_dig->getCellId();
922  m_tree_amp = ecl_dig->getAmp();
923  m_tree_en = en;
924  m_tree_E1Etot = en / trkEClustLab[iCharge];
925  m_tree_E1E2 = en / crystalEnergies2[iCharge];
926  m_tree_E1p = en / trkpLab[iCharge];
927  m_tree_timetsPreviousTimeCalibs = time - ts_prevCalib[iCharge] - tcrate_prevCalib[iCharge];
928  m_tree_timeF = ecl_dig->getTimeFit() * TICKS_TO_NS;
929  m_tree_time = time;
930  m_tree_quality = ecl_dig->getQuality();
931  m_tree_t0 = evt_t0;
932  m_tree_t0_unc = evt_t0_unc;
933  m_E_DIV_p = E_DIV_p[iCharge];
934  m_tree_evtNum = m_EventMetaData->getEvent();
935  m_crystalCrate = m_crystalMapper->getCrateID(ecl_cal->getCellId());
936  m_runNum = m_EventMetaData->getRun();
937 
938  if (m_saveTree) {
939  m_dbgTree_electrons->Fill();
940  }
941 // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
942 
943  // Fill histogram with information about maximum energy crystal energy fraction
944  getObjectPtr<TH1F>("maxEcrsytalEnergyFraction")->Fill(en / trkEClustLab[iCharge]);
945 
946 
947  }
948 
949 
950 
951  // Check both electrons to see if their cluster energy / track momentum is good.
952  // The Belle II physics book shows that this is the main way of separating electrons from other particles
953  // Done in the centre of mass reference frame although I believe E/p is invariant under a boost.
954  bool E_DIV_p_instance_passed[2] = {false, false};
955  double E_DIV_p_CUT = 0.7;
956  for (int icharge = 0; icharge < 2; icharge++) {
957  E_DIV_p_instance_passed[icharge] = E_DIV_p[icharge] > E_DIV_p_CUT;
958  }
959  if (!E_DIV_p_instance_passed[0] || !E_DIV_p_instance_passed[1]) {
960  return;
961  }
962  // E/p sufficiently large
963  cutIndexPassed++;
964  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
965  B2DEBUG(22, "Cutflow: E_i/p_i > " << E_DIV_p_CUT << ": index = " << cutIndexPassed);
966 
967 
968 
969  // Start of cuts on both the combined system of tracks and energy clusters
970 
971  double invMassTrk = (trkp4Lab[0] + trkp4Lab[1]).M();
972  double invMass_CUT = 0.9;
973  m_massInvTracks = invMassTrk; // invariant mass of the two tracks
974 
975 // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
976  if (m_saveTree) {
977  m_dbgTree_event->Fill();
978  }
979 // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
980 
981  bool invMassCutsPassed = invMassTrk > (invMass_CUT * boostrotate.getCMSEnergy());
982  if (!invMassCutsPassed) {
983  return;
984  }
985  // Invariable mass of the two tracks are above the minimum
986  cutIndexPassed++;
987  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
988  B2DEBUG(22, "Cutflow: m(track 1+2) > " << invMass_CUT << "*E_COM = " << invMass_CUT << " * " << boostrotate.getCMSEnergy() <<
989  " : index = " << cutIndexPassed);
990 
991 
992 
993  //== Fill output histogram.
994  for (int iCharge = 0; iCharge < 2; iCharge++) {
995  if (crystalCutsPassed[iCharge]) {
996  getObjectPtr<TH2F>("TimevsCrysPrevCrateCalibPrevCrystCalib")->Fill((crystalIDs[iCharge]) + 0.001,
997  times_noPreviousCalibrations[iCharge] - ts_prevCalib[iCharge] - tcrate_prevCalib[iCharge], 1);
998  getObjectPtr<TH2F>("TimevsCratePrevCrateCalibPrevCrystCalib")->Fill((crateIDs[iCharge]) + 0.001,
999  times_noPreviousCalibrations[iCharge] - ts_prevCalib[iCharge] - tcrate_prevCalib[iCharge], 1);
1000  getObjectPtr<TH2F>("TimevsCrysNoCalibrations")->Fill((crystalIDs[iCharge]) + 0.001, times_noPreviousCalibrations[iCharge], 1);
1001  getObjectPtr<TH2F>("TimevsCrateNoCalibrations")->Fill((crateIDs[iCharge]) + 0.001, times_noPreviousCalibrations[iCharge], 1);
1002  getObjectPtr<TH2F>("TimevsCrysPrevCrateCalibNoCrystCalib")->Fill((crystalIDs[iCharge]) + 0.001,
1003  times_noPreviousCalibrations[iCharge] - tcrate_prevCalib[iCharge], 1);
1004  getObjectPtr<TH2F>("TimevsCrateNoCrateCalibPrevCrystCalib")->Fill((crateIDs[iCharge]) + 0.001,
1005  times_noPreviousCalibrations[iCharge] - ts_prevCalib[iCharge], 1);
1006 
1007  // Record number of crystals used from the event. Should be exactly two.
1008  numCrystalsPassingCuts++;
1009 
1010  }
1011  }
1012 
1013 
1014  // Change cutflow method for this bit ... don't call return because we used to call the hadron cluster stuff afterwards
1015  //
1016  if (crystalCutsPassed[0] || crystalCutsPassed[1]) {
1017  // At least one ECL crystal time and quality cuts passed
1018  cutIndexPassed++;
1019  getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
1020  B2DEBUG(22, "Cutflow: At least one crystal time and quality cuts passed: index = " << cutIndexPassed);
1021 
1022  getObjectPtr<TH1F>("numCrystalEntriesPerEvent")->Fill(numCrystalsPassingCuts);
1023  }
1024 
1025 
1026  // Save final information to the tree after all cuts are applied
1027  for (int iCharge = 0; iCharge < 2; iCharge++) {
1028  if (crystalCutsPassed[iCharge]) {
1029  // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
1030  m_tree_evtNum = m_EventMetaData->getEvent();
1031  m_tree_cid = crystalIDs[iCharge];
1032  //m_tree_time = times[iCharge];
1033  m_tree_time = times_noPreviousCalibrations[iCharge];
1034  m_crystalCrate = crateIDs[iCharge];
1035  m_runNum = m_EventMetaData->getRun();
1036  m_tree_en = crystalEnergies[iCharge]; // for studies of ts as a function of energy
1037  m_tree_E1Etot = crystalEnergies[iCharge] / trkEClustLab[iCharge];
1038  m_tree_E1E2 = crystalEnergies[iCharge] / crystalEnergies2[iCharge];
1039  m_tree_E1p = crystalEnergies[iCharge] / trkpLab[iCharge];
1040  m_tree_timetsPreviousTimeCalibs = times_noPreviousCalibrations[iCharge] - ts_prevCalib[iCharge] - tcrate_prevCalib[iCharge];
1041  m_tree_t0 = evt_t0;
1042  m_tree_t0_ECLclosestCDC = evt_t0_ECL_closestCDC;
1043  m_tree_t0_ECL_minChi2 = evt_t0_ECL_minChi2;
1044  m_tree_tClust = clusterTime[iCharge];
1045 
1046  m_massInvTracks = invMassTrk; // This is probably already set but I'll set it again anyways just so that it is clear
1047 
1048  if (m_saveTree) {
1049  m_dbgTree_allCuts->Fill();
1050  }
1051  // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
1052  }
1053  }
1054 
1055 
1056 
1057 
1058  if (crystalCutsPassed[0] && crystalCutsPassed[1] &&
1059  numClustersPerTrack[0] == 1 && numClustersPerTrack[1] == 1) {
1060  m_tree_enNeg = trkEClustLab[0];
1061  m_tree_enPlus = trkEClustLab[1];
1062  m_tree_tClustNeg = clusterTime[0];
1063  m_tree_tClustPos = clusterTime[1];
1064  m_tree_maxEcrystPosClust = times_noPreviousCalibrations[0] - ts_prevCalib[0] - tcrate_prevCalib[0];
1065  m_tree_maxEcrystNegClust = times_noPreviousCalibrations[1] - ts_prevCalib[1] - tcrate_prevCalib[1];
1066  m_tree_t0 = evt_t0;
1067  m_tree_t0_ECLclosestCDC = evt_t0_ECL_closestCDC;
1068  m_tree_t0_ECL_minChi2 = evt_t0_ECL_minChi2;
1069 
1070  // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
1071  if (m_saveTree) {
1072  m_dbgTree_evt_allCuts->Fill();
1073  }
1074  // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
1075  }
1076 
1077 
1078  B2DEBUG(26, "m_tree_maxEcrystPosClust + evt_t0 = " << m_tree_maxEcrystPosClust + evt_t0);
1079  B2DEBUG(26, "m_tree_maxEcrystNegClust + evt_t0 = " << m_tree_maxEcrystNegClust + evt_t0);
1080  B2DEBUG(26, "CDC evt_t0 = " << evt_t0);
1081  B2DEBUG(26, "ECL min chi2 even t0, m_tree_t0_ECL_minChi2 = " << m_tree_t0_ECL_minChi2);
1082 
1083 
1084 
1085  for (long unsigned int digit_i = 0; digit_i < time_ECLCaldigits_bothClusters.size(); digit_i++) {
1086  m_runNum = m_EventMetaData->getRun();
1087  m_tree_evtNum = m_EventMetaData->getEvent();
1088  m_tree_ECLCalDigitTime = time_ECLCaldigits_bothClusters[digit_i];
1089  m_tree_ECLCalDigitE = E_ECLCaldigits_bothClusters[digit_i];
1090  m_tree_ECLDigitAmplitude = amp_ECLDigits_bothClusters[digit_i];
1091  m_tree_t0 = evt_t0;
1092  m_tree_t0_ECLclosestCDC = evt_t0_ECL_closestCDC;
1093  m_tree_t0_ECL_minChi2 = evt_t0_ECL_minChi2;
1094  m_tree_timetsPreviousTimeCalibs = times_noPreviousCalibrations[chargeID_ECLCaldigits_bothClusters[digit_i]] -
1095  ts_prevCalib[chargeID_ECLCaldigits_bothClusters[digit_i]] -
1096  tcrate_prevCalib[chargeID_ECLCaldigits_bothClusters[digit_i]];
1097  m_tree_cid = cid_ECLCaldigits_bothClusters[digit_i];
1098 
1099  // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
1100  if (m_saveTree) {
1101  m_dbgTree_crys_allCuts->Fill();
1102  }
1103  // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
1104 
1105  }
1106 
1107 
1108  B2DEBUG(26, "This was for event number = " << m_tree_evtNum);
1109 
1110 }
1111 
1112 
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).
double m_tree_d0
Track d0 for debug TTree output.
double m_tree_maxEcrystNegClust
Time of the highest energy crystal in the cluster associated to negatively charged track,...
std::vector< float > m_CrateTime
vector obtained from DB object
std::vector< int > m_eclDigitID
ECL digit id sorter.
double m_tree_timeF
ECL fitting time for debug TTree output.
double m_tree_E1Etot
Energy of crystal with maximum energy within ECL cluster divided by total cluster energy,...
double m_tree_tClustNeg
Cluster time of cluster associated to negatively charged track, ns for debug TTree output.
DBObjPtr< ECLCrystalCalib > m_CrateTimeDB
Time offset from crate time calibration (also this calibration) from database.
TTree * m_dbgTree_crystals
Debug TTree output per crystal.
StoreObjPtr< EventT0 > m_eventT0
StoreObjPtr for T0.
StoreArray< ECLDigit > m_eclDigitArray
Required input array of ECLDigits.
double m_tree_clustCrysE
crystal energy, only for the crystals that meet all the selection criteria for debug TTree output
std::vector< float > m_FlightTime
vector obtained from DB object
double m_tree_timetsPreviousTimeCalibs
Time for Ts distribution after application of previous time calibrations for debug TTree output.
double m_looseTrkD0
Loose track d0 minimum cut.
double m_tree_E1E2
Energy of crystal with maximum energy within ECL cluster divided by second most energetic crystal in ...
double m_tree_ECLDigitAmplitude
Amplitude (used to calculate energy) of an ECLDigit within a cluster, for debug TTree output.
double m_tree_tClustPos
Cluster time of cluster associated to positively charged track, ns for debug TTree output.
int m_tree_quality
ECL fit quality for debug TTree output.
bool m_saveTree
If true, save TTree with more detailed event info.
double m_tree_enNeg
Energy of cluster associated to negatively charged track, GeV for debug TTree output.
DBObjPtr< ECLCrystalCalib > m_ElectronicsTimeDB
Time offset from electronics calibration from database.
double m_tree_tClust
Cluster time of a cluster, ns for debug TTree output.
double m_E_DIV_p
Energy divided by momentum, for debug TTree output.
double m_tree_en
Energy of crystal with maximum energy within ECL cluster, GeV for debug TTree output.
std::vector< int > m_eclCalDigitID
ECL cal digit id sorter.
bool m_storeCalib
Boolean for whether or not to store the previous calibration calibration constants.
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.
DBObjPtr< ECLReferenceCrystalPerCrateCalib > m_RefCrystalsCalibDB
Crystal IDs of the one reference crystal per crate from database.
int m_crystalCrate
Crate id for the crystal.
int m_tree_cid
ECL Cell ID (1..ECLElementNumbers::c_NCrystals) for debug TTree output.
double m_tree_time
Time for Ts distribution for debug TTree output.
double m_massInvTracks
invariant mass of the two tracks, for debug TTree output
double m_tree_nCDChits
Number of CDC hits along the track for debug TTree output.
double m_tree_enPlus
Energy of cluster associated to positively charged track, GeV for debug TTree output.
int m_tree_evtNum
Event number for debug TTree output.
double m_tightTrkZ0
Tight track z0 minimum cut.
DBObjPtr< Belle2::ECLChannelMap > m_channelMapDB
Mapper of ecl channels to various other objects, like crates.
std::vector< short > m_RefCrystalsCalib
vector obtained from DB object
double m_tree_t0_ECL_minChi2
EventT0 (from ECL) min chi2 for debug TTree output.
double m_tree_maxEcrystPosClust
Time of the highest energy crystal in the cluster associated to positively charged track,...
void collect() override
Select events and crystals and accumulate histograms.
std::vector< float > m_CrateTimeUnc
uncertainty vector obtained from DB object
TTree * m_dbgTree_allCuts
Debug TTree output after all cuts.
double m_tree_ECLCalDigitTime
Time of an ECLCalDigit within a cluster, ns for debug TTree output.
StoreArray< ECLCalDigit > m_eclCalDigitArray
Required input array of ECLCalDigits.
TTree * m_dbgTree_event
Debug TTree output per event.
double m_tree_z0
Track z0 for debug TTree output.
std::vector< float > m_EperCrys
ECL cal digit energy for each crystal.
bool skipTrgSel
flag to skip the trigger skim selection in the module
double m_tree_p
Track momentum for debug TTree output.
void prepare() override
Define histograms and read payloads from DB.
std::vector< float > m_PreviousCrystalTime
vector obtained from DB object
DBObjPtr< ECLCrystalCalib > m_PreviousCrystalTimeDB
Time offset from previous crystal time calibration (this calibration) from database.
double m_tightTrkD0
Tight track d0 minimum cut.
double m_tree_t0
EventT0 (not from ECL) for debug TTree output.
int m_tree_amp
Fitting amplitude from ECL for debug TTree output.
double m_tree_E1p
Energy of crystal with maximum energy within ECL cluster divided by total cluster energy divided by t...
TTree * m_dbgTree_evt_allCuts
Debug TTree output per event after all cuts.
std::vector< float > m_Electronics
vector obtained from DB object
void inDefineHisto() override
Replacement for defineHisto() in CalibrationCollector modules.
int m_charge
particle charge, for debug TTree output
TTree * m_dbgTree_electrons
Output tree with detailed event data.
std::vector< float > m_PreviousCrystalTimeUnc
vector obtained from DB object
std::unique_ptr< Belle2::ECL::ECLTimingUtilities > m_ECLTimeUtil
ECL timing tools.
StoreArray< Track > tracks
StoreArray for tracks.
double m_tree_t0_ECLclosestCDC
EventT0 (from ECL) closest to CDC for debug TTree output.
TTree * m_dbgTree_crys_allCuts
Debug TTree output per crystal after all cuts.
double m_tree_clustCrysE_DIV_maxEcrys
ratio of crystal energy to energy of the crystal that has the maximum energy, only for the crystals t...
DBObjPtr< ECLCrystalCalib > m_FlightTimeDB
Time offset from flight time b/w IP and crystal from database.
DBObjPtr< ECLCrystalCalib > m_ElectronicsDB
electronics amplitude calibration from database Scale amplitudefor each crystal and for dead pre-amps
double m_tree_ECLCalDigitE
Energy of an ECLCalDigit within a cluster, GeV for debug TTree output.
StoreObjPtr< EventMetaData > m_EventMetaData
Event metadata.
double m_tree_t0_unc
EventT0 uncertainty for debug TTree output.
std::unique_ptr< Belle2::ECL::ECLChannelMapper > m_crystalMapper
ECL object for keeping track of mapping between crystals and crates etc.
virtual ~ECLBhabhaTCollectorModule()
Module destructor.
double m_hadronEventT0_TO_bhabhaEventT0_correction
correction to apply to CDC event t0 values in bhabha events to correct for CDC event t0 bias compared...
std::vector< float > m_ElectronicsTime
vector obtained from DB object
StoreObjPtr< SoftwareTriggerResult > m_TrgResult
Store array for Trigger selection.
TTree * m_dbgTree_tracks
Debug TTree output per track.
double m_looseTrkZ0
Loose track z0 minimum cut.
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
int getCellId() const
Get Cell ID.
Definition: ECLCalDigit.h:114
double getEnergy() const
Get Calibrated Energy.
Definition: ECLCalDigit.h:119
ECL cluster data.
Definition: ECLCluster.h:27
@ c_nPhotons
CR is split into n photons (N1)
Class to store ECL digitized hits (output of ECLDigi) relation to ECLHit filled in ecl/modules/eclDig...
Definition: ECLDigit.h:24
int getAmp() const
Get Fitting Amplitude.
Definition: ECLDigit.h:70
int getQuality() const
Get Fitting Quality.
Definition: ECLDigit.h:80
int getCellId() const
Get Cell ID.
Definition: ECLDigit.h:64
int getTimeFit() const
Get Fitting Time.
Definition: ECLDigit.h:75
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;.
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.