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