Belle II Software development
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
32using namespace Belle2;
33using namespace ECL;
34using namespace std;
35
36//-----------------------------------------------------------------
37// Register the Module
38//-----------------------------------------------------------------
39REG_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 "SVD 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 XWiki 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 XWiki 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_closestSVD", &m_tree_t0_ECLclosestSVD)->SetTitle("T0 ECL closest to SVD 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_closestSVD", &m_tree_t0_ECLclosestSVD)->SetTitle("T0 ECL closest to SVD 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_closestSVD", &m_tree_t0_ECLclosestSVD)->SetTitle("T0 ECL closest to SVD 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 svdEventT0 = new TH1D("svdEventT0", ";EventT0 from SVD", 200, -50, 50);
294 registerObject<TH1D>("svdEventT0", svdEventT0);
295
296 auto maxEcrsytalEnergyFraction = new TH1F("maxEcrsytalEnergyFraction",
297 ";Maximum energy crystal energy / (sum) cluster energy;Number", 22, 0, 1.1);
298 registerObject<TH1F>("maxEcrsytalEnergyFraction", maxEcrsytalEnergyFraction);
299
300 auto refCrysIDzeroingCrate = new TH1F("refCrysIDzeroingCrate", ";cell id;Boolean - is reference crystal",
302 registerObject<TH1F>("refCrysIDzeroingCrate", refCrysIDzeroingCrate);
303
304 auto SVDEventT0Correction = new TH1F("SVDEventT0Correction", ";;SVD event t0 offset correction [ns]", 1, 1, 2);
305 registerObject<TH1F>("SVDEventT0Correction", SVDEventT0Correction);
306
307
308 //=== Required data objects
309 m_eventT0.isRequired();
310 tracks.isRequired();
311 m_eclClusterArray.isRequired();
312 m_eclCalDigitArray.isRequired();
313 m_eclDigitArray.isRequired();
314
315 B2INFO("hadronEventT0_TO_bhabhaEventT0_correction = " << m_hadronEventT0_TO_bhabhaEventT0_correction <<
316 " ns correction to SVD event t0 will be applied");
317
318 B2INFO("skipTrgSel = " << skipTrgSel);
319
320}
321
323{
324 int cutIndexPassed = 0;
325 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
326 B2DEBUG(22, "Cutflow: no cuts: index = " << cutIndexPassed);
327 B2DEBUG(22, "Event number = " << m_EventMetaData->getEvent());
328
329
330 // --- Check the trigger skim is the type that has two tracks
331
332 /* If we skip the trigger skim selection then still fill the cutflow histogram
333 just so that the positions don't change. */
334 if (!skipTrgSel) {
335 if (!m_TrgResult.isValid()) {
336 B2WARNING("SoftwareTriggerResult required to select bhabha event is not found");
337 return;
338 }
339
340 /* Release05: bhabha_all is grand skim = bhabha+bhabhaecl+radee. We only want
341 to look at the 2 track bhabha events. */
342 const std::map<std::string, int>& fresults = m_TrgResult->getResults();
343 if (fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end()) {
344 B2WARNING("Can't find required bhabha trigger identifier");
345 return;
346 }
347
348 const bool eBhabha = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
350 B2DEBUG(22, "eBhabha (trigger passed) = " << eBhabha);
351
352 if (!eBhabha) {
353 return;
354 }
355 }
356
357 /* Fill the histgram showing that the trigger skim cut passed OR that we
358 are skipping this selection. */
359 cutIndexPassed++;
360 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
361 B2DEBUG(22, "Cutflow: Trigger cut passed: index = " << cutIndexPassed);
362
363
364
365 /* Use ECLChannelMapper to get other detector indices for the crystals
366 For conversion from CellID to crate, shaper, and channel ids.
367 The initialization function automatically checks to see if the
368 object has been initialized and ifthe payload has changed and
369 thus needs updating. */
370 bool ECLchannelMapHasChanged = m_channelMapDB.hasChanged();
371 if (ECLchannelMapHasChanged) {
372 B2INFO("ECLBhabhaTCollectorModule::collect() " << LogVar("ECLchannelMapHasChanged", ECLchannelMapHasChanged));
373 if (!m_crystalMapper->initFromDB()) {
374 B2FATAL("ECLBhabhaTCollectorModule::collect() : Can't initialize eclChannelMapper!");
375 }
376 }
377
378
379 //== Get expected energies and calibration constants from DB. Need to call
380 // hasChanged() for later comparison
381 if (m_ElectronicsDB.hasChanged()) {
382 m_Electronics = m_ElectronicsDB->getCalibVector();
383 }
384 if (m_ElectronicsTimeDB.hasChanged()) {
385 m_ElectronicsTime = m_ElectronicsTimeDB->getCalibVector();
386 }
387 if (m_FlightTimeDB.hasChanged()) {
388 m_FlightTime = m_FlightTimeDB->getCalibVector();
389 }
390
391 // Get the previous crystal time offset (the same thing that this calibration is meant to calculate).
392 // This can be used for testing purposes, and for the crate time offset.
393 if (m_PreviousCrystalTimeDB.hasChanged()) {
396 }
397
398 B2DEBUG(29, "Finished checking if previous crystal time payload has changed");
399 if (m_CrateTimeDB.hasChanged()) {
400 m_CrateTime = m_CrateTimeDB->getCalibVector();
401 m_CrateTimeUnc = m_CrateTimeDB->getCalibUncVector();
402 }
403 B2DEBUG(29, "Finished checking if previous crate time payload has changed");
404 B2DEBUG(29, "m_CrateTime size = " << m_CrateTime.size());
405 B2DEBUG(25, "Crate time +- uncertainty [0]= " << m_CrateTime[0] << " +- " << m_CrateTimeUnc[0]);
406 B2DEBUG(25, "Crate time +- uncertainty [8735]= " << m_CrateTime[8735] << " +- " << m_CrateTimeUnc[8735]);
407
408 B2DEBUG(29, "Finished checking if previous crate time payload has changed");
409 if (m_RefCrystalsCalibDB.hasChanged()) {
410 m_RefCrystalsCalib = m_RefCrystalsCalibDB->getReferenceCrystals();
411 }
412 B2DEBUG(29, "Finished checking if reference crystal cell ids payload has changed");
413
414
415 B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrystalTimeOffset from the database"
416 << LogVar("IoV", m_PreviousCrystalTimeDB.getIoV())
417 << LogVar("Checksum", m_PreviousCrystalTimeDB.getChecksum()));
418 B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrateTimeOffset from the database"
419 << LogVar("IoV", m_CrateTimeDB.getIoV())
420 << LogVar("Checksum", m_CrateTimeDB.getChecksum()));
421 B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrystalElectronics from the database"
422 << LogVar("IoV", m_ElectronicsDB.getIoV())
423 << LogVar("Checksum", m_ElectronicsDB.getChecksum()));
424 B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrystalElectronicsTime from the database"
425 << LogVar("IoV", m_ElectronicsTimeDB.getIoV())
426 << LogVar("Checksum", m_ElectronicsTimeDB.getChecksum()));
427 B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLCrystalFlightTime from the database"
428 << LogVar("IoV", m_FlightTimeDB.getIoV())
429 << LogVar("Checksum", m_FlightTimeDB.getChecksum()));
430 B2DEBUG(25, "ECLBhabhaTCollector:: loaded ECLReferenceCrystalPerCrateCalib from the database"
431 << LogVar("IoV", m_RefCrystalsCalibDB.getIoV())
432 << LogVar("Checksum", m_RefCrystalsCalibDB.getChecksum()));
433
434
435
436 // Conversion coefficient from ADC ticks to nanoseconds
437 // TICKS_TO_NS ~ 0.4913 ns/clock tick
438 // 1/(4fRF) = 0.4913 ns/clock tick, where fRF is the accelerator RF frequency
439 const double TICKS_TO_NS = 1.0 / (4.0 * EclConfiguration::getRF()) * 1e3;
440
441
442 vector<float> Crate_time_ns(52, 0.0);
444 // Make a crate time offset vector with an entry per crate (instead of per crystal) and convert from ADC counts to ns.
445 for (int crysID = 1; crysID <= ECLElementNumbers::c_NCrystals; crysID++) {
446 int crateID_temp = m_crystalMapper->getCrateID(crysID);
447 Crate_time_ns[crateID_temp - 1] = m_CrateTime[crysID] * TICKS_TO_NS;
448 }
449
450
451
457 for (int crateID_temp = 1; crateID_temp <= 52; crateID_temp++) {
458 getObjectPtr<TH1F>("refCrysIDzeroingCrate")->Fill(m_RefCrystalsCalib[crateID_temp - 1] + 0.001);
459 }
460
461
463 for (int crysID = 1; crysID <= ECLElementNumbers::c_NCrystals; crysID++) {
464 getObjectPtr<TH1F>("TsDatabase")->SetBinContent(crysID + 0.001, m_PreviousCrystalTime[crysID - 1]);
465 getObjectPtr<TH1F>("TsDatabaseUnc")->SetBinContent(crysID + 0.001, m_PreviousCrystalTimeUnc[crysID - 1]);
466 getObjectPtr<TH1F>("TcrateDatabase")->SetBinContent(crysID + 0.001, m_CrateTime[crysID - 1]);
467 getObjectPtr<TH1F>("TcrateUncDatabase")->SetBinContent(crysID + 0.001, m_CrateTimeUnc[crysID - 1]);
468 }
469 if (m_storeCalib) {
470 B2INFO("ECLBhabhaTCollector:: ECLCrystalTimeOffset from the database information:"
471 << LogVar("IoV", m_PreviousCrystalTimeDB.getIoV())
472 << LogVar("Checksum", m_PreviousCrystalTimeDB.getChecksum()));
473 B2INFO("First event so print out previous ts values");
474 for (int crysID = 1; crysID <= ECLElementNumbers::c_NCrystals; crysID++) {
475 B2INFO("cid = " << crysID << ", Ts previous = " << m_PreviousCrystalTime[crysID - 1]);
476 }
477 m_storeCalib = false;
478 }
479
480
481
482
483 for (int crateID_temp = 1; crateID_temp <= 52; crateID_temp++) {
484 getObjectPtr<TH1F>("tcrateDatabase_ns")->SetBinContent(crateID_temp + 0.001, Crate_time_ns[crateID_temp - 1]);
485 }
486
487 // Use a histogram with only one bin as a counter to know the number of times the database histograms were filled.
488 // This is mostly useful for the talg when running over multiple runs and trying to read ts values.
489 getObjectPtr<TH1I>("databaseCounter")->SetBinContent(1, 1);
490
491
492
493 // Save what SVD event t0 correction was applied
494 getObjectPtr<TH1F>("SVDEventT0Correction")->SetBinContent(1, m_hadronEventT0_TO_bhabhaEventT0_correction);
495
496
497
498
499 /* Getting the event t0 using the full event t0 rather than from the SVD specifically */
500 double evt_t0 = -1;
501 double evt_t0_unc = -1;
502 double evt_t0_ECL_closestSVD = -1;
503 double evt_t0_ECL_minChi2 = -1;
504
505 // Determine if the event t0 comes from the SVD then extract the information about it
506 if (!m_eventT0.isValid()) {
507 //cout << "event t0 not valid\n";
508 return;
509 } else if (!m_eventT0->hasTemporaryEventT0(Const::EDetector::SVD)) {
510 //cout << "no event t0\n";
511 return;
512 } else {
513 // SVD eventT0 is available
514 cutIndexPassed++;
515 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
516 B2DEBUG(22, "Cutflow: Event t0 exists: index = " << cutIndexPassed);
517
518 // Get eventT0 from SVD
519 const auto bestSVDEventT0Candidate = m_eventT0->getBestSVDTemporaryEventT0();
520 evt_t0 = bestSVDEventT0Candidate->eventT0; // time value
521 evt_t0_unc = bestSVDEventT0Candidate->eventT0Uncertainty; // uncertainty on event t0
522 getObjectPtr<TH1D>("svdEventT0")->Fill(evt_t0);
523
524 // Correct the SVD event t0 value for the bhabha bias, if any
526
527
528 // Get the ECL event t0 for comparison - validations
529 if (m_eventT0->hasTemporaryEventT0(Const::EDetector::ECL)) {
530 vector<EventT0::EventT0Component> evt_t0_list_ECL = m_eventT0->getTemporaryEventT0s(Const::EDetector::ECL);
531
532
533 double smallest_SVD_ECL_t0_diff = fabs(evt_t0_list_ECL[0].eventT0 - evt_t0);
534 int smallest_SVD_ECL_t0_diff_idx = 0;
535 for (long unsigned int ECLi = 0; ECLi < evt_t0_list_ECL.size(); ECLi++) {
536 double tempt_ECL_t0 = evt_t0_list_ECL[ECLi].eventT0;
537 if (fabs(tempt_ECL_t0 - evt_t0) < smallest_SVD_ECL_t0_diff) {
538 smallest_SVD_ECL_t0_diff = fabs(tempt_ECL_t0 - evt_t0);
539 smallest_SVD_ECL_t0_diff_idx = ECLi;
540 }
541 }
542
543 evt_t0_ECL_closestSVD = evt_t0_list_ECL[smallest_SVD_ECL_t0_diff_idx].eventT0; // time value
544 B2DEBUG(26, "evt_t0_ECL_closestSVD = " << evt_t0_ECL_closestSVD);
545
546
547
548 double smallest_ECL_t0_minChi2 = evt_t0_list_ECL[0].quality;
549 int smallest_ECL_t0_minChi2_idx = 0;
550
551 B2DEBUG(26, "evt_t0_list_ECL[0].quality = " << evt_t0_list_ECL[0].quality
552 << ", with ECL event t0 = " << evt_t0_list_ECL[0].eventT0);
553
554 for (long unsigned int ECLi = 0; ECLi < evt_t0_list_ECL.size(); ECLi++) {
555 B2DEBUG(26, "evt_t0_list_ECL[" << ECLi << "].quality = " << evt_t0_list_ECL[ECLi].quality
556 << ", with ECL event t0 = " <<
557 evt_t0_list_ECL[ECLi].eventT0);
558 if (evt_t0_list_ECL[ECLi].quality < smallest_ECL_t0_minChi2) {
559 smallest_ECL_t0_minChi2 = evt_t0_list_ECL[ECLi].quality;
560 smallest_ECL_t0_minChi2_idx = ECLi;
561 }
562 }
563
564 evt_t0_ECL_minChi2 = evt_t0_list_ECL[smallest_ECL_t0_minChi2_idx].eventT0; // time value
565
566 B2DEBUG(26, "evt_t0_ECL_minChi2 = " << evt_t0_ECL_minChi2);
567 B2DEBUG(26, "smallest_ECL_t0_minChi2_idx = " << smallest_ECL_t0_minChi2_idx);
568 }
569 }
570
571
572
573 /* Determine the energies for each of the crystals since this isn't naturally connected to the cluster.
574 Also determine the indexing of the ecl cal digits and the ecl digits
575 Taken from Chris's ec/modules/eclGammaGammaECollector */
576
577 // Resize vectors
581
582
583 int idx = 0;
584 for (auto& eclCalDigit : m_eclCalDigitArray) {
585 int tempCrysID = eclCalDigit.getCellId() - 1;
586 m_EperCrys[tempCrysID] = eclCalDigit.getEnergy();
587 m_eclCalDigitID[tempCrysID] = idx;
588 idx++;
589 }
590
591 idx = 0;
592 for (auto& eclDigit : m_eclDigitArray) {
593 int tempCrysID = eclDigit.getCellId() - 1;
594 m_eclDigitID[tempCrysID] = idx;
595 idx++;
596 }
597
598
599
600
601 //---------------------------------------------------------------------
602 //..Some utilities
603 PCmsLabTransform boostrotate;
604
605 //---------------------------------------------------------------------
606 //..Track properties, including 2 maxp tracks. Use pion (211) mass hypothesis,
607 // which is the only particle hypothesis currently available???
608 double maxp[2] = {0., 0.};
609 int maxiTrk[2] = { -1, -1};
610 int nTrkAll = tracks.getEntries();
611
612 int nTrkLoose = 0;
613 int nTrkTight = 0;
615 /* Loop over all the tracks to define the tight and loose selection tracks.
616 We will select events with only 2 tight tracks and no additional loose tracks.
617 Tight tracks are a subset of looses tracks. */
618 for (int iTrk = 0; iTrk < nTrkAll; iTrk++) {
619 // Get track biasing towards the particle being a pion based on what particle types
620 // are used for reconstruction at this stage.
621 const TrackFitResult* tempTrackFit = tracks[iTrk]->getTrackFitResultWithClosestMass(Const::pion);
622 if (not tempTrackFit) {continue;}
623
624 // Collect track info to be used for categorizing
625 short charge = tempTrackFit->getChargeSign();
626 double z0 = tempTrackFit->getZ0();
627 double d0 = tempTrackFit->getD0();
628 int nCDChits = tempTrackFit->getHitPatternCDC().getNHits();
629 double p = tempTrackFit->getMomentum().R();
630
631 // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
632 //== Save debug TTree with detailed information if necessary.
633 m_tree_d0 = d0;
634 m_tree_z0 = z0;
635 m_tree_p = p;
636 m_charge = charge;
637 m_tree_nCDChits = nCDChits;
638
639 if (m_saveTree) {
640 m_dbgTree_tracks->Fill();
641 }
642 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
643
644
645 /* Test if loose track */
646
647 // d0 and z0 cuts
648 if (fabs(d0) > m_looseTrkD0) {
649 continue;
650 }
651 if (fabs(z0) > m_looseTrkZ0) {
652 continue;
653 }
654 // Number of hits in the CDC
655 if (nCDChits < 1) {
656 continue;
657 }
658 nTrkLoose++;
659
660
661
662 /* Test if the loose track is also a tight track */
663
664 // Number of hits in the CDC
665 if (nCDChits < 20) {
666 continue;
667 }
668
669
670 // d0 and z0 cuts
671 if (fabs(d0) > m_tightTrkD0) {
672 continue;
673 }
674 if (fabs(z0) > m_tightTrkZ0) {
675 continue;
676 }
677 nTrkTight++;
678
679 // Sorting of tight tracks. Not really required as we only want two tight tracks (at the moment) but okay.
680 //..Find the maximum p negative [0] and positive [1] tracks
681 int icharge = 0;
682 if (charge > 0) {icharge = 1;}
683 if (p > maxp[icharge]) {
684 maxp[icharge] = p;
685 maxiTrk[icharge] = iTrk;
686 }
687
688 }
689 /* After that last section the numbers of loose and tight tracks are known as well as the
690 index of the loose tracks that have the highest p negatively charged and highest p positively
691 charged tracks as measured in the centre of mass frame */
692
693
694 if (nTrkTight != 2) {
695 return;
696 }
697 // There are exactly two tight tracks
698 cutIndexPassed++;
699 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
700 B2DEBUG(22, "Cutflow: Two tight tracks: index = " << cutIndexPassed);
701
702
703 if (nTrkLoose != 2) {
704 return;
705 }
706 // There are exactly two loose tracks as well, i.e. no additional loose tracks
707 cutIndexPassed++;
708 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
709 B2DEBUG(22, "Cutflow: No additional loose tracks: index = " << cutIndexPassed);
710
711 /* Determine if the two tracks have the opposite electric charge.
712 We know this because the track indices stores the max pt track in [0] for negatively charged track
713 and [1] fo the positively charged track. If both are filled then both a negatively charged
714 and positively charged track were found. */
715 bool oppositelyChargedTracksPassed = maxiTrk[0] != -1 && maxiTrk[1] != -1;
716 if (!oppositelyChargedTracksPassed) {
717 return;
718 }
719 // The two tracks have the opposite electric charges.
720 cutIndexPassed++;
721 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
722 B2DEBUG(22, "Cutflow: Oppositely charged tracks: index = " << cutIndexPassed);
723
724
725
726
727 //---------------------------------------------------------------------
728 /* Determine associated energy clusters to each of the two tracks. Sum the energies of the
729 multiple clusters to each track and find the crystal with the maximum energy within all
730 the sets of clusters associated to the tracks*/
731 double trkEClustLab[2] = {0., 0.};
732 double trkEClustCOM[2] = {0., 0.};
733 double trkpLab[2];
734 double trkpCOM[2];
735 ROOT::Math::PxPyPzEVector trkp4Lab[2];
736 ROOT::Math::PxPyPzEVector trkp4COM[2];
737
738 // Index of the cluster and the crystal that has the highest energy crystal for the two tracks
739 int crysIDMax[2] = { -1, -1 };
740 double crysEMax[2] = { -1, -1 };
741 double crysE2Max[2] = { -1, -1 };
742 int numClustersPerTrack[2] = { 0, 0 };
743
744 double clusterTime[2] = {0, 0};
745
746 double E_DIV_p[2];
747
748 vector<double> time_ECLCaldigits_bothClusters;
749 vector<int> cid_ECLCaldigits_bothClusters;
750 vector<double> E_ECLCaldigits_bothClusters;
751 vector<double> amp_ECLDigits_bothClusters;
752 vector<int> chargeID_ECLCaldigits_bothClusters;
753
754 for (int icharge = 0; icharge < 2; icharge++) {
755 if (maxiTrk[icharge] > -1) {
756 B2DEBUG(22, "looping over the 2 max pt tracks");
757
758 const TrackFitResult* tempTrackFit = tracks[maxiTrk[icharge]]->getTrackFitResultWithClosestMass(Const::pion);
759 if (not tempTrackFit) {continue;}
760 trkp4Lab[icharge] = tempTrackFit->get4Momentum();
761 trkp4COM[icharge] = boostrotate.rotateLabToCms() * trkp4Lab[icharge];
762 trkpLab[icharge] = trkp4Lab[icharge].P();
763 trkpCOM[icharge] = trkp4COM[icharge].P();
764
765
766 /* For each cluster associated to the current track, sum up the energies to get the total
767 energy of all clusters associated to the track and find which crystal has the highest
768 energy from all those clusters*/
769 auto eclClusterRelationsFromTracks = tracks[maxiTrk[icharge]]->getRelationsTo<ECLCluster>();
770 for (unsigned int clusterIdx = 0; clusterIdx < eclClusterRelationsFromTracks.size(); clusterIdx++) {
771
772 B2DEBUG(22, "Looking at clusters. index = " << clusterIdx);
773 auto cluster = eclClusterRelationsFromTracks[clusterIdx];
774 bool goodClusterType = false;
775
776 if (cluster->hasHypothesis(Belle2::ECLCluster::EHypothesisBit::c_nPhotons)) {
777 trkEClustLab[icharge] += cluster->getEnergy(Belle2::ECLCluster::EHypothesisBit::c_nPhotons);
778 goodClusterType = true;
779 numClustersPerTrack[icharge]++;
780 }
781
782 if (goodClusterType) {
783
784 clusterTime[icharge] = cluster->getTime();
785
786 auto eclClusterRelations = cluster->getRelationsTo<ECLCalDigit>("ECLCalDigits");
787
788 // Find the crystal that has the largest energy
789 for (unsigned int ir = 0; ir < eclClusterRelations.size(); ir++) {
790 const auto calDigit = eclClusterRelations.object(ir);
791 int tempCrysID = calDigit->getCellId() - 1;
792 double tempE = m_EperCrys[tempCrysID];
793
794 int eclDigitIndex = m_eclDigitID[tempCrysID];
795 ECLDigit* ecl_dig = m_eclDigitArray[eclDigitIndex];
796
797 // for the max E crystal
798 if (tempE > crysEMax[icharge]) {
799 // Set 2nd highest E crystal to the info from the highest E crystal
800 crysE2Max[icharge] = crysEMax[icharge];
801 // Set the highest E crystal to the current crystal
802 crysEMax[icharge] = tempE;
803 crysIDMax[icharge] = tempCrysID;
804 }
805 // for the 2nd highest E crystal
806 if (tempE > crysE2Max[icharge] && tempCrysID != crysIDMax[icharge]) {
807 crysE2Max[icharge] = tempE;
808 }
809
810 // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
811 // If we drop the information about the second highest energy crystal we could use
812 // m_eclClusterArray[ic]->getMaxECellId()
813
814 B2DEBUG(26, "calDigit(ir" << ir << ") time = " << calDigit->getTime() << "ns , with E = " << tempE << " GeV");
815 time_ECLCaldigits_bothClusters.push_back(calDigit->getTime());
816 cid_ECLCaldigits_bothClusters.push_back(tempCrysID);
817 E_ECLCaldigits_bothClusters.push_back(tempE);
818 amp_ECLDigits_bothClusters.push_back(ecl_dig->getAmp());
819 chargeID_ECLCaldigits_bothClusters.push_back(icharge);
820
821 }
822 }
823 }
824 trkEClustCOM[icharge] = trkEClustLab[icharge] * trkpCOM[icharge] / trkpLab[icharge];
825
826 // Check both electrons to see if their cluster energy / track momentum is good.
827 // The Belle II physics book shows that this is the main way of separating electrons from other particles
828 // Done in the centre of mass reference frame although I believe E/p is invariant under a boost.
829 E_DIV_p[icharge] = trkEClustCOM[icharge] / trkpCOM[icharge];
830
831 }
832 }
833 /* At the end of this section the 3-momenta magnitudes and the cluster energies are known
834 for the two saved track indices for both the lab and COM frames.
835 The crystal with the maximum energy, one associated to each track, is recorded*/
836
837
838
839 //=== Check each crystal in the processed event and fill histogram.
840
841 int numCrystalsPassingCuts = 0;
842
843 int crystalIDs[2] = { -1, -1};
844 int crateIDs[2] = { -1, -1};
845 double ts_prevCalib[2] = { -1, -1};
846 double tcrate_prevCalib[2] = { -1, -1};
847 double times_noPreviousCalibrations[2] = { -1, -1};
848 bool crystalCutsPassed[2] = {false, false};
849 double crystalEnergies[2] = { -1, -1};
850 double crystalEnergies2[2] = { -1, -1};
851
852 for (int iCharge = 0; iCharge < 2; iCharge++) {
853 int crystal_idx = crysIDMax[iCharge];
854 int eclCalDigitIndex = m_eclCalDigitID[crystal_idx];
855 int eclDigitIndex = m_eclDigitID[crystal_idx];
856
857 ECLDigit* ecl_dig = m_eclDigitArray[eclDigitIndex];
858 ECLCalDigit* ecl_cal = m_eclCalDigitArray[eclCalDigitIndex];
859
860 //== Check whether specific ECLDigits should be excluded.
861
862 auto en = ecl_cal->getEnergy();
863 auto amplitude = ecl_dig->getAmp();
864 crystalEnergies[iCharge] = en;
865
866 int cid = ecl_dig->getCellId();
867 double time = ecl_dig->getTimeFit() * TICKS_TO_NS - evt_t0;
868
869 // Offset time by electronics calibration and flight time calibration.
870 time -= m_ElectronicsTime[cid - 1] * TICKS_TO_NS;
871 time -= m_FlightTime[cid - 1];
872
873
874 // Apply the time walk correction: time shift as a function of the amplitude corrected by the electronics calibration.
875 // The electronics calibration also accounts for crystals that have a dead pre-amp and thus half the normal amplitude.
876 double energyTimeShift = m_ECLTimeUtil->energyDependentTimeOffsetElectronic(amplitude * m_Electronics[cid - 1]) * TICKS_TO_NS;
877
878 B2DEBUG(29, "cellid = " << cid << ", amplitude = " << amplitude << ", time before t(E) shift = " << time <<
879 ", t(E) shift = " << energyTimeShift << " ns");
880 time -= energyTimeShift;
881
882
883 // Cell ID should be within specified range.
884 if (cid < m_minCrystal || cid > m_maxCrystal) continue;
885
886 // Absolute time should be in specified range condition.
887 if (fabs(time) > m_timeAbsMax) continue;
888
889 // Fit quality flag -- choose only events with best fit quality
890 if (ecl_dig->getQuality() != 0) continue;
891
892 //== Save time and crystal information. Fill plot after both electrons are tested
893 crystalIDs[iCharge] = cid;
894 crateIDs[iCharge] = m_crystalMapper->getCrateID(ecl_cal->getCellId());
895
896
897 ts_prevCalib[iCharge] = m_PreviousCrystalTime[cid - 1] * TICKS_TO_NS;
898 tcrate_prevCalib[iCharge] = m_CrateTime[cid - 1] * TICKS_TO_NS;
899 times_noPreviousCalibrations[iCharge] = time;
900
901
902 B2DEBUG(26, "iCharge = " << iCharge);
903 B2DEBUG(26, "crateIDs[iCharge] = " << crateIDs[iCharge]);
904 B2DEBUG(26, "times_noPreviousCalibrations[iCharge] = " << times_noPreviousCalibrations[iCharge]);
905 B2DEBUG(26, "tcrate_prevCalib[iCharge] = " << tcrate_prevCalib[iCharge]);
906 B2DEBUG(26, "ts_prevCalib[iCharge] = " << ts_prevCalib[iCharge]);
907
908
909 crystalCutsPassed[iCharge] = true;
910
911
912 // For second most energetic energy crystal
913 crystalEnergies2[iCharge] = crysE2Max[iCharge];
914
915
916// >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
917 //== Save debug TTree with detailed information if necessary.
918 m_tree_cid = ecl_dig->getCellId();
919 m_tree_amp = ecl_dig->getAmp();
920 m_tree_en = en;
921 m_tree_E1Etot = en / trkEClustLab[iCharge];
922 m_tree_E1E2 = en / crystalEnergies2[iCharge];
923 m_tree_E1p = en / trkpLab[iCharge];
924 m_tree_timetsPreviousTimeCalibs = time - ts_prevCalib[iCharge] - tcrate_prevCalib[iCharge];
925 m_tree_timeF = ecl_dig->getTimeFit() * TICKS_TO_NS;
926 m_tree_time = time;
927 m_tree_quality = ecl_dig->getQuality();
928 m_tree_t0 = evt_t0;
929 m_tree_t0_unc = evt_t0_unc;
930 m_E_DIV_p = E_DIV_p[iCharge];
931 m_tree_evtNum = m_EventMetaData->getEvent();
932 m_crystalCrate = m_crystalMapper->getCrateID(ecl_cal->getCellId());
933 m_runNum = m_EventMetaData->getRun();
934
935 if (m_saveTree) {
936 m_dbgTree_electrons->Fill();
937 }
938// <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
939
940 // Fill histogram with information about maximum energy crystal energy fraction
941 getObjectPtr<TH1F>("maxEcrsytalEnergyFraction")->Fill(en / trkEClustLab[iCharge]);
942
943
944 }
945
946
947
948 // Check both electrons to see if their cluster energy / track momentum is good.
949 // The Belle II physics book shows that this is the main way of separating electrons from other particles
950 // Done in the centre of mass reference frame although I believe E/p is invariant under a boost.
951 bool E_DIV_p_instance_passed[2] = {false, false};
952 double E_DIV_p_CUT = 0.7;
953 for (int icharge = 0; icharge < 2; icharge++) {
954 E_DIV_p_instance_passed[icharge] = E_DIV_p[icharge] > E_DIV_p_CUT;
955 }
956 if (!E_DIV_p_instance_passed[0] || !E_DIV_p_instance_passed[1]) {
957 return;
958 }
959 // E/p sufficiently large
960 cutIndexPassed++;
961 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
962 B2DEBUG(22, "Cutflow: E_i/p_i > " << E_DIV_p_CUT << ": index = " << cutIndexPassed);
963
964
965
966 // Start of cuts on both the combined system of tracks and energy clusters
967
968 double invMassTrk = (trkp4Lab[0] + trkp4Lab[1]).M();
969 double invMass_CUT = 0.9;
970 m_massInvTracks = invMassTrk; // invariant mass of the two tracks
971
972// >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
973 if (m_saveTree) {
974 m_dbgTree_event->Fill();
975 }
976// <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
977
978 bool invMassCutsPassed = invMassTrk > (invMass_CUT * boostrotate.getCMSEnergy());
979 if (!invMassCutsPassed) {
980 return;
981 }
982 // Invariable mass of the two tracks are above the minimum
983 cutIndexPassed++;
984 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
985 B2DEBUG(22, "Cutflow: m(track 1+2) > " << invMass_CUT << "*E_COM = " << invMass_CUT << " * " << boostrotate.getCMSEnergy() <<
986 " : index = " << cutIndexPassed);
987
988
989
990 //== Fill output histogram.
991 for (int iCharge = 0; iCharge < 2; iCharge++) {
992 if (crystalCutsPassed[iCharge]) {
993 getObjectPtr<TH2F>("TimevsCrysPrevCrateCalibPrevCrystCalib")->Fill((crystalIDs[iCharge]) + 0.001,
994 times_noPreviousCalibrations[iCharge] - ts_prevCalib[iCharge] - tcrate_prevCalib[iCharge], 1);
995 getObjectPtr<TH2F>("TimevsCratePrevCrateCalibPrevCrystCalib")->Fill((crateIDs[iCharge]) + 0.001,
996 times_noPreviousCalibrations[iCharge] - ts_prevCalib[iCharge] - tcrate_prevCalib[iCharge], 1);
997 getObjectPtr<TH2F>("TimevsCrysNoCalibrations")->Fill((crystalIDs[iCharge]) + 0.001, times_noPreviousCalibrations[iCharge], 1);
998 getObjectPtr<TH2F>("TimevsCrateNoCalibrations")->Fill((crateIDs[iCharge]) + 0.001, times_noPreviousCalibrations[iCharge], 1);
999 getObjectPtr<TH2F>("TimevsCrysPrevCrateCalibNoCrystCalib")->Fill((crystalIDs[iCharge]) + 0.001,
1000 times_noPreviousCalibrations[iCharge] - tcrate_prevCalib[iCharge], 1);
1001 getObjectPtr<TH2F>("TimevsCrateNoCrateCalibPrevCrystCalib")->Fill((crateIDs[iCharge]) + 0.001,
1002 times_noPreviousCalibrations[iCharge] - ts_prevCalib[iCharge], 1);
1003
1004 // Record number of crystals used from the event. Should be exactly two.
1005 numCrystalsPassingCuts++;
1006
1007 }
1008 }
1009
1010
1011 // Change cutflow method for this bit ... don't call return because we used to call the hadron cluster stuff afterwards
1012 //
1013 if (crystalCutsPassed[0] || crystalCutsPassed[1]) {
1014 // At least one ECL crystal time and quality cuts passed
1015 cutIndexPassed++;
1016 getObjectPtr<TH1F>("cutflow")->Fill(cutIndexPassed);
1017 B2DEBUG(22, "Cutflow: At least one crystal time and quality cuts passed: index = " << cutIndexPassed);
1018
1019 getObjectPtr<TH1F>("numCrystalEntriesPerEvent")->Fill(numCrystalsPassingCuts);
1020 }
1021
1022
1023 // Save final information to the tree after all cuts are applied
1024 for (int iCharge = 0; iCharge < 2; iCharge++) {
1025 if (crystalCutsPassed[iCharge]) {
1026 // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
1027 m_tree_evtNum = m_EventMetaData->getEvent();
1028 m_tree_cid = crystalIDs[iCharge];
1029 //m_tree_time = times[iCharge];
1030 m_tree_time = times_noPreviousCalibrations[iCharge];
1031 m_crystalCrate = crateIDs[iCharge];
1032 m_runNum = m_EventMetaData->getRun();
1033 m_tree_en = crystalEnergies[iCharge]; // for studies of ts as a function of energy
1034 m_tree_E1Etot = crystalEnergies[iCharge] / trkEClustLab[iCharge];
1035 m_tree_E1E2 = crystalEnergies[iCharge] / crystalEnergies2[iCharge];
1036 m_tree_E1p = crystalEnergies[iCharge] / trkpLab[iCharge];
1037 m_tree_timetsPreviousTimeCalibs = times_noPreviousCalibrations[iCharge] - ts_prevCalib[iCharge] - tcrate_prevCalib[iCharge];
1038 m_tree_t0 = evt_t0;
1039 m_tree_t0_ECLclosestSVD = evt_t0_ECL_closestSVD;
1040 m_tree_t0_ECL_minChi2 = evt_t0_ECL_minChi2;
1041 m_tree_tClust = clusterTime[iCharge];
1042
1043 m_massInvTracks = invMassTrk; // This is probably already set but I'll set it again anyways just so that it is clear
1044
1045 if (m_saveTree) {
1046 m_dbgTree_allCuts->Fill();
1047 }
1048 // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
1049 }
1050 }
1051
1052
1053
1054
1055 if (crystalCutsPassed[0] && crystalCutsPassed[1] &&
1056 numClustersPerTrack[0] == 1 && numClustersPerTrack[1] == 1) {
1057 m_tree_enNeg = trkEClustLab[0];
1058 m_tree_enPlus = trkEClustLab[1];
1059 m_tree_tClustNeg = clusterTime[0];
1060 m_tree_tClustPos = clusterTime[1];
1061 m_tree_maxEcrystPosClust = times_noPreviousCalibrations[0] - ts_prevCalib[0] - tcrate_prevCalib[0];
1062 m_tree_maxEcrystNegClust = times_noPreviousCalibrations[1] - ts_prevCalib[1] - tcrate_prevCalib[1];
1063 m_tree_t0 = evt_t0;
1064 m_tree_t0_ECLclosestSVD = evt_t0_ECL_closestSVD;
1065 m_tree_t0_ECL_minChi2 = evt_t0_ECL_minChi2;
1066
1067 // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
1068 if (m_saveTree) {
1069 m_dbgTree_evt_allCuts->Fill();
1070 }
1071 // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
1072 }
1073
1074
1075 B2DEBUG(26, "m_tree_maxEcrystPosClust + evt_t0 = " << m_tree_maxEcrystPosClust + evt_t0);
1076 B2DEBUG(26, "m_tree_maxEcrystNegClust + evt_t0 = " << m_tree_maxEcrystNegClust + evt_t0);
1077 B2DEBUG(26, "SVD evt_t0 = " << evt_t0);
1078 B2DEBUG(26, "ECL min chi2 even t0, m_tree_t0_ECL_minChi2 = " << m_tree_t0_ECL_minChi2);
1079
1080
1081
1082 for (long unsigned int digit_i = 0; digit_i < time_ECLCaldigits_bothClusters.size(); digit_i++) {
1083 m_runNum = m_EventMetaData->getRun();
1084 m_tree_evtNum = m_EventMetaData->getEvent();
1085 m_tree_ECLCalDigitTime = time_ECLCaldigits_bothClusters[digit_i];
1086 m_tree_ECLCalDigitE = E_ECLCaldigits_bothClusters[digit_i];
1087 m_tree_ECLDigitAmplitude = amp_ECLDigits_bothClusters[digit_i];
1088 m_tree_t0 = evt_t0;
1089 m_tree_t0_ECLclosestSVD = evt_t0_ECL_closestSVD;
1090 m_tree_t0_ECL_minChi2 = evt_t0_ECL_minChi2;
1091 m_tree_timetsPreviousTimeCalibs = times_noPreviousCalibrations[chargeID_ECLCaldigits_bothClusters[digit_i]] -
1092 ts_prevCalib[chargeID_ECLCaldigits_bothClusters[digit_i]] -
1093 tcrate_prevCalib[chargeID_ECLCaldigits_bothClusters[digit_i]];
1094 m_tree_cid = cid_ECLCaldigits_bothClusters[digit_i];
1095
1096 // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> Tree saving
1097 if (m_saveTree) {
1098 m_dbgTree_crys_allCuts->Fill();
1099 }
1100 // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Tree saving
1101
1102 }
1103
1104
1105 B2DEBUG(26, "This was for event number = " << m_tree_evtNum);
1106
1107}
1108
1109
Calibration collector module base class.
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
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.
double m_tree_t0_ECLclosestSVD
EventT0 (from ECL) closest to SVD for debug TTree output.
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.
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 SVD event t0 values in bhabha events to correct for SVD 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:118
double getEnergy() const
Get Calibrated Energy.
Definition: ECLCalDigit.h:123
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.
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
@ c_accept
Accept this event.
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.
STL namespace.