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