Belle II Software development
KLMEventT0EstimatorModule.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 <klm/modules/KLMEventT0Estimator/KLMEventT0EstimatorModule.h>
11
12/* KLM headers — only those not already pulled in via the own header. */
13#include <klm/bklm/geometry/Module.h>
14#include <klm/dataobjects/bklm/BKLMElementNumbers.h>
15#include <klm/dataobjects/eklm/EKLMElementNumbers.h>
16
17/* Basf2 headers — only those not already pulled in via the own header. */
18#include <framework/gearbox/Const.h>
19#include <framework/logging/Logger.h>
20#include <framework/datastore/DataStore.h>
21
22using namespace Belle2;
23using namespace Belle2::bklm;
24using namespace Belle2::EKLM;
25
26REG_MODULE(KLMEventT0Estimator);
27
28/* Constructor and destructor. */
29
30KLMEventT0EstimatorModule::KLMEventT0EstimatorModule() :
32 m_geoParB(nullptr),
33 m_geoParE(nullptr),
34 m_transformE(nullptr)
35{
36 setDescription("Estimate per-event T0 using KLM digits matched to extrapolated tracks (BKLM/EKLM scintillators and RPC) with per-event track averages, uncertainties, and final combined KLM value.");
37 setPropertyFlags(c_ParallelProcessingCertified);
38
39 addParam("MuonListName", m_MuonListName,
40 "Muon (or generic) ParticleList name used to access tracks and KLM relations (e.g. 'mu+:forT0').",
41 std::string("mu+:forT0"));
42 addParam("useCDCTemporaryT0", m_useCDCTemporaryT0,
43 "Read CDC temporary EventT0 as a seed/diagnostic (not applied to the mean).",
44 true);
45 addParam("IgnoreBackwardPropagation", m_ignoreBackward,
46 "Ignore backward-propagated ExtHits when forming entry/exit pairs.",
47 false);
48 addParam("histogramDirectoryName", m_histDirName,
49 "Top directory for KLMEventT0Estimator histograms inside the ROOT file.",
50 std::string("KLMEventT0Estimator"));
51 addParam("histogramSubdirUncorrected", m_histSubdirUncorr,
52 "Subdirectory name for uncorrected timing histograms.",
53 std::string("uncorrected"));
54
55 // ADC cut parameters
56 addParam("ADCCut_BKLM_Scint_Min", m_ADCCut_BKLM_Scint_Min,
57 "Minimum ADC charge cut for BKLM scintillator. Set to 0 to disable lower cut.",
58 30.0);
59 addParam("ADCCut_BKLM_Scint_Max", m_ADCCut_BKLM_Scint_Max,
60 "Maximum ADC charge cut for BKLM scintillator. Set to large value to disable upper cut.",
61 320.0);
62 addParam("ADCCut_EKLM_Scint_Min", m_ADCCut_EKLM_Scint_Min,
63 "Minimum ADC charge cut for EKLM scintillator. Set to 0 to disable lower cut.",
64 40.0);
65 addParam("ADCCut_EKLM_Scint_Max", m_ADCCut_EKLM_Scint_Max,
66 "Maximum ADC charge cut for EKLM scintillator. Set to large value to disable upper cut.",
67 350.0);
68
69}
70
71KLMEventT0EstimatorModule::~KLMEventT0EstimatorModule()
72{
73 delete m_transformE;
74}
75
76/* Histogram definition. */
77
79{
80 // Parent directory for this module
81 TDirectory* topdir = gDirectory->mkdir(m_histDirName.c_str());
82 TDirectory::TContext ctxTop{gDirectory, topdir};
83
84 // Subdirectory for uncorrected timing histograms (always created)
85 TDirectory* d_unc = topdir->mkdir(m_histSubdirUncorr.c_str());
86
87 auto H1 = [](const char* n, const char* t, int nb, double lo, double hi) {
88 return new TH1D(n, t, nb, lo, hi);
89 };
90
91 /* Uncorrected timing histograms. */
92 {
93 TDirectory::TContext ctxUnc{gDirectory, d_unc};
94
95 // --- per_track/ subdirectory ---
96 TDirectory* d_per_track = d_unc->mkdir("per_track");
97 {
98 TDirectory::TContext ctxTrk{gDirectory, d_per_track};
99 m_hT0Trk_BKLM_Scint = H1("h_t0trk_bklm_scint", "Per-track T0 (BKLM Scint);T0 [ns]", 800, -100, 100);
100 m_hT0Trk_BKLM_RPC = H1("h_t0trk_bklm_rpc", "Per-track T0 (BKLM RPC);T0 [ns]", 800, -100, 100);
101 m_hT0Trk_EKLM_Scint = H1("h_t0trk_eklm_scint", "Per-track T0 (EKLM Scint);T0 [ns]", 800, -100, 100);
102 }
103
104 // --- per_event/ subdirectory ---
105 TDirectory* d_per_event = d_unc->mkdir("per_event");
106 {
107 TDirectory::TContext ctxEvt{gDirectory, d_per_event};
108
109 // Track-average
110 m_hT0Evt_TrkAvg_BKLM_Scint = H1("h_t0evt_trkavg_bklm_scint", "Per-event T0 (track-avg, BKLM Scint);T0 [ns]", 800, -100,
111 100);
112 m_hT0Evt_TrkAvg_BKLM_RPC = H1("h_t0evt_trkavg_bklm_rpc", "Per-event T0 (track-avg, BKLM RPC);T0 [ns]", 800, -100,
113 100);
114 m_hT0Evt_TrkAvg_EKLM_Scint = H1("h_t0evt_trkavg_eklm_scint", "Per-event T0 (track-avg, EKLM Scint);T0 [ns]", 800, -100,
115 100);
116 m_hT0Evt_TrkAvg_All = H1("h_t0evt_trkavg_all", "Per-event T0 (track-avg, all categories);T0 [ns]", 800, -100,
117 100);
118 m_hT0Evt_TrkAvg_BKLM_Scint_SEM = H1("h_t0evt_trkavg_bklm_scint_sem", "SEM (track-avg, BKLM Scint);SEM [ns]", 800, 0.0, 20.0);
119 m_hT0Evt_TrkAvg_BKLM_RPC_SEM = H1("h_t0evt_trkavg_bklm_rpc_sem", "SEM (track-avg, BKLM RPC);SEM [ns]", 800, 0.0, 20.0);
120 m_hT0Evt_TrkAvg_EKLM_Scint_SEM = H1("h_t0evt_trkavg_eklm_scint_sem", "SEM (track-avg, EKLM Scint);SEM [ns]", 800, 0.0, 20.0);
121 m_hT0Evt_TrkAvg_All_SEM = H1("h_t0evt_trkavg_all_sem", "SEM (track-avg, all categories);SEM [ns]", 800, 0.0, 20.0);
122
123 // Final-source audit: all possible combinations of B(KLM Scint), E(KLM Scint), R(PC)
124 m_hFinalSource = new TH1I("h_final_source", "Final KLM source;;events", 7, 0.5, 7.5);
125 m_hFinalSource->GetXaxis()->SetBinLabel(1, "B only");
126 m_hFinalSource->GetXaxis()->SetBinLabel(2, "E only");
127 m_hFinalSource->GetXaxis()->SetBinLabel(3, "R only");
128 m_hFinalSource->GetXaxis()->SetBinLabel(4, "B+E");
129 m_hFinalSource->GetXaxis()->SetBinLabel(5, "B+R");
130 m_hFinalSource->GetXaxis()->SetBinLabel(6, "E+R");
131 m_hFinalSource->GetXaxis()->SetBinLabel(7, "B+E+R");
132 }
133 } // end uncorrected/ directory
134
135}
136
137/* Lifecycle. */
138
140{
141 // Register that we define histograms
142 REG_HISTOGRAM;
143
144 // Inputs
145 m_MuonList.isRequired(m_MuonListName);
146 m_tracks.isRequired();
147
148 // Geometry
152
153 // Log ADC cut settings
154 B2DEBUG(20, "KLMEventT0Estimator: ADC cuts configured:"
155 << LogVar("BKLM Scint min", m_ADCCut_BKLM_Scint_Min)
156 << LogVar("BKLM Scint max", m_ADCCut_BKLM_Scint_Max)
157 << LogVar("EKLM Scint min", m_ADCCut_EKLM_Scint_Min)
158 << LogVar("EKLM Scint max", m_ADCCut_EKLM_Scint_Max));
159}
160
162{
163 if (!m_eventT0HitResolution.isValid())
164 B2FATAL("KLMEventT0Estimator: KLM EventT0 hit resolution data are not available.");
165
166 B2DEBUG(20, "KLMEventT0Estimator: Using calibrated per-hit resolution."
167 << LogVar("sigma_RPC (ns)", m_eventT0HitResolution->getSigmaRPC())
168 << LogVar("sigma_BKLM_Scint (ns)", m_eventT0HitResolution->getSigmaBKLMScint())
169 << LogVar("sigma_EKLM_Scint (ns)", m_eventT0HitResolution->getSigmaEKLMScint()));
170}
171
175
177
178/* Helper methods. */
179
180bool KLMEventT0EstimatorModule::passesADCCut(double charge, int subdetector, int layer, bool inRPC) const
181{
182 // RPC: No ADC cut applied
183 if (subdetector == KLMElementNumbers::c_BKLM && (inRPC || layer >= BKLMElementNumbers::c_FirstRPCLayer)) {
184 return true;
185 }
186
187 // BKLM Scintillator
188 if (subdetector == KLMElementNumbers::c_BKLM) {
189 return (charge >= m_ADCCut_BKLM_Scint_Min && charge <= m_ADCCut_BKLM_Scint_Max);
190 }
191
192 // EKLM Scintillator
193 if (subdetector == KLMElementNumbers::c_EKLM) {
194 return (charge >= m_ADCCut_EKLM_Scint_Min && charge <= m_ADCCut_EKLM_Scint_Max);
195 }
196
197 return true;
198}
199
200double KLMEventT0EstimatorModule::getHitSigma(int subdetector, int layer, bool inRPC, int plane) const
201{
202 if (!m_eventT0HitResolution.isValid()) {
203 B2ERROR("KLMEventT0Estimator: Calibrated hit resolution payload not available!");
204 return 1.0; // Fallback
205 }
206
207 if (subdetector == KLMElementNumbers::c_BKLM) {
208 if (inRPC || layer >= BKLMElementNumbers::c_FirstRPCLayer) {
209 // Use direction-specific RPC resolution if available (version 2+)
210 // Fall back to combined RPC resolution for backward compatibility
211 if (plane == BKLMElementNumbers::c_ZPlane) {
212 float sigmaZ = m_eventT0HitResolution->getSigmaRPCZ();
213 return (sigmaZ > 0.0) ? sigmaZ : m_eventT0HitResolution->getSigmaRPC();
214 } else {
215 float sigmaPhi = m_eventT0HitResolution->getSigmaRPCPhi();
216 return (sigmaPhi > 0.0) ? sigmaPhi : m_eventT0HitResolution->getSigmaRPC();
217 }
218 } else {
219 return m_eventT0HitResolution->getSigmaBKLMScint();
220 }
221 } else { // EKLM
222 return m_eventT0HitResolution->getSigmaEKLMScint();
223 }
224}
225
227KLMEventT0EstimatorModule::matchExt(unsigned int key, ExtMap& v_ExtHits)
228{
229 ExtHit* entryHit = nullptr;
230 ExtHit* exitHit = nullptr;
231 auto itlow = v_ExtHits.lower_bound(key);
232 auto itup = v_ExtHits.upper_bound(key);
233 for (auto it = itlow; it != itup; ++it) {
234 if (!entryHit || it->second.getTOF() < entryHit->getTOF()) entryHit = &(it->second);
235 if (!exitHit || it->second.getTOF() > exitHit->getTOF()) exitHit = &(it->second);
236 }
237 return std::make_pair(entryHit, exitHit);
238}
239
241 ExtMap& scintMap,
242 ExtMap& rpcMap)
243{
244 scintMap.clear();
245 rpcMap.clear();
246
247 RelationVector<ExtHit> extHits = track->getRelationsTo<ExtHit>();
248 KLMMuidLikelihood* muidLikelihood = track->getRelatedTo<KLMMuidLikelihood>();
249
250 for (const ExtHit& eHit : extHits) {
251 if (eHit.getStatus() != EXT_EXIT) continue;
252 if (m_ignoreBackward && eHit.isBackwardPropagated()) continue;
253
254 const bool isB = (eHit.getDetectorID() == Const::EDetector::BKLM);
255 const bool isE = (eHit.getDetectorID() == Const::EDetector::EKLM);
256 if (!isB && !isE) continue;
257
258 int copyId = eHit.getCopyID();
259 int tFor, tSec, tLay, tPla, tStr;
260 int tSub = -1;
261
262 if (isE) {
265 &tFor, &tLay, &tSec, &tPla, &tStr);
266 }
267 if (isB) {
270 &tFor, &tSec, &tLay, &tPla, &tStr);
271 }
272 if (tSub < 0) continue;
273
274 bool crossed = false;
275
276 if (isB) {
277 crossed = muidLikelihood
278 ? muidLikelihood->isExtrapolatedBarrelLayerCrossed(tLay - 1)
279 : true;
280 if (!crossed) continue;
281
282 const bool isRPC = (tLay >= BKLMElementNumbers::c_FirstRPCLayer);
283
284 if (isRPC) {
285 // RPC: match by module
286 unsigned int moduleKey =
287 m_elementNum->moduleNumber(tSub, tFor, tSec, tLay);
288 rpcMap.insert(std::make_pair(moduleKey, eHit));
289 } else {
290 // BKLM scintillator: match by channel
291 unsigned int channelKey =
292 m_elementNum->channelNumber(tSub, tFor, tSec, tLay, tPla, tStr);
293 if (m_channelStatus.isValid() &&
294 m_channelStatus->getChannelStatus(channelKey) != KLMChannelStatus::c_Normal)
295 continue;
296 scintMap.insert(std::make_pair(channelKey, eHit));
297 }
298 }
299
300 if (isE) {
301 crossed = muidLikelihood
302 ? muidLikelihood->isExtrapolatedEndcapLayerCrossed(tLay - 1)
303 : true;
304 if (!crossed) continue;
305
306 unsigned int channelKey =
307 m_elementNum->channelNumber(tSub, tFor, tSec, tLay, tPla, tStr);
308 if (m_channelStatus.isValid() &&
309 m_channelStatus->getChannelStatus(channelKey) != KLMChannelStatus::c_Normal)
310 continue;
311 scintMap.insert(std::make_pair(channelKey, eHit));
312 }
313 }
314}
315
316/* Per-subdetector accumulators. */
317
318namespace {
319 // Weighted accumulation using per-hit sigma
320 inline void acc_stat_weighted(double t, double sigma, double& sumW, double& sumWT)
321 {
322 if (sigma <= 0.0 || !std::isfinite(sigma)) return;
323 const double w = 1.0 / (sigma * sigma);
324 sumW += w;
325 sumWT += w * t;
326 }
327}
328
329/* EKLM scintillator */
331 const ExtMap& scintMap,
332 double& sumW, double& sumWT)
333{
334 DBObjPtr<KLMTimeConstants> timeConstants;
335 DBObjPtr<KLMTimeCableDelay> timeCableDelay;
336
337 const double delayScint = timeConstants.isValid()
338 ? timeConstants->getDelay(KLMTimeConstants::c_EKLM)
339 : 0.0;
340
341 HepGeom::Point3D<double> hitGlobal_ext, hitLocal_ext;
342
343 for (const KLMHit2d& hit2d : klmHit2ds) {
344 if (hit2d.getSubdetector() != KLMElementNumbers::c_EKLM) continue;
345
346 RelationVector<KLMDigit> digits = hit2d.getRelationsTo<KLMDigit>();
347 if (digits.size() == 0) continue;
348
349 for (const KLMDigit& d : digits) {
350 if (!d.isGood()) continue;
351
352 unsigned int cid = d.getUniqueChannelID();
353 if (m_channelStatus.isValid() &&
354 m_channelStatus->getChannelStatus(cid) != KLMChannelStatus::c_Normal) continue;
355
356 // Apply ADC cut
357 if (!passesADCCut(d.getCharge(), KLMElementNumbers::c_EKLM, d.getLayer(), false)) {
358 continue;
359 }
360
361 // Match using channel ID
362 ExtPair ex = const_cast<KLMEventT0EstimatorModule*>(this)->matchExt(cid, const_cast<ExtMap&>(scintMap));
363 if (!ex.first || !ex.second) continue;
364 const double flyTime = 0.5 * (ex.first->getTOF() + ex.second->getTOF());
365
366 // Distance along strip to readout
367 const ROOT::Math::XYZVector posGlobExt = 0.5 * (ex.first->getPosition() + ex.second->getPosition());
368 hitGlobal_ext.setX(posGlobExt.X() / Unit::mm * CLHEP::mm);
369 hitGlobal_ext.setY(posGlobExt.Y() / Unit::mm * CLHEP::mm);
370 hitGlobal_ext.setZ(posGlobExt.Z() / Unit::mm * CLHEP::mm);
371
372 const double Lmm = m_geoParE->getStripLength(d.getStrip()) / CLHEP::mm * Unit::mm;
373 const HepGeom::Transform3D* tr = m_transformE->getStripGlobalToLocal(const_cast<KLMDigit*>(&d));
374 hitLocal_ext = (*tr) * hitGlobal_ext;
375 const double dist_mm = 0.5 * Lmm - hitLocal_ext.x() / CLHEP::mm * Unit::mm;
376
377 // Components
378 const double Trec = d.getTime();
379 const double Tcable = timeCableDelay.isValid() ? timeCableDelay->getTimeDelay(cid) : 0.0;
380 const double Tprop = dist_mm * delayScint;
381 const double Tfly = flyTime;
382
383 // Correct digit time
384 double t = Trec;
385 if (timeCableDelay.isValid()) t -= Tcable;
386 t -= Tprop;
387
388 const double t0_est = t - Tfly;
389 if (!std::isfinite(t0_est)) continue;
390
391 // Weighted accumulation using calibrated sigma
392 const double sigma = getHitSigma(KLMElementNumbers::c_EKLM, d.getLayer(), false);
393 acc_stat_weighted(t0_est, sigma, sumW, sumWT);
394 }
395 }
396}
397
398/* BKLM scintillator */
400 const ExtMap& scintMap,
401 double& sumW, double& sumWT)
402{
403 DBObjPtr<KLMTimeConstants> timeConstants;
404 DBObjPtr<KLMTimeCableDelay> timeCableDelay;
405
406 const double delayScint = timeConstants.isValid()
407 ? timeConstants->getDelay(KLMTimeConstants::c_BKLM)
408 : 0.0;
409
410 for (KLMHit2d& hit2d : klmHit2ds) {
411 if (hit2d.getSubdetector() != KLMElementNumbers::c_BKLM) continue;
412 if (hit2d.inRPC()) continue;
413 if (hit2d.isOutOfTime()) continue;
414
415 RelationVector<BKLMHit1d> b1ds = hit2d.getRelationsTo<BKLMHit1d>();
416 if (b1ds.size() == 0) continue;
417
418 const bklm::Module* mod = m_geoParB->findModule(hit2d.getSection(), hit2d.getSector(), hit2d.getLayer());
419 const ROOT::Math::XYZVector posG2d = hit2d.getPosition();
420
421 for (const BKLMHit1d& h1d : b1ds) {
422 RelationVector<KLMDigit> digits = h1d.getRelationsTo<KLMDigit>();
423
424 for (const KLMDigit& d : digits) {
425 if (d.inRPC() || !d.isGood()) continue;
426
427 unsigned int cid = d.getUniqueChannelID();
428 if (m_channelStatus.isValid() &&
429 m_channelStatus->getChannelStatus(cid) != KLMChannelStatus::c_Normal) continue;
430
431 // Apply ADC cut
432 if (!passesADCCut(d.getCharge(), KLMElementNumbers::c_BKLM, d.getLayer(), false)) {
433 continue;
434 }
435
436 // Match using channel ID
437 ExtPair p = const_cast<KLMEventT0EstimatorModule*>(this)->matchExt(cid, const_cast<ExtMap&>(scintMap));
438 if (!p.first || !p.second) continue;
439
440 const double flyTime = 0.5 * (p.first->getTOF() + p.second->getTOF());
441 const ROOT::Math::XYZVector posGext = 0.5 * (p.first->getPosition() + p.second->getPosition());
442
443 // Gate in local
444 const CLHEP::Hep3Vector locExt = mod->globalToLocal(CLHEP::Hep3Vector(posGext.X(), posGext.Y(), posGext.Z()), true);
445 const CLHEP::Hep3Vector locHit2 = mod->globalToLocal(CLHEP::Hep3Vector(posG2d.X(), posG2d.Y(), posG2d.Z()), true);
446 const CLHEP::Hep3Vector diff = locExt - locHit2;
447 if (std::fabs(diff.z()) > mod->getZStripWidth() || std::fabs(diff.y()) > mod->getPhiStripWidth()) continue;
448
449 // Prop distance
450 const bool isPhiReadout = h1d.isPhiReadout();
451 double propaLen = mod->getPropagationDistance(locExt, d.getStrip(), isPhiReadout);
452
453 // Components
454 const double Trec = d.getTime();
455 const double Tcable = timeCableDelay.isValid() ? timeCableDelay->getTimeDelay(cid) : 0.0;
456 const double Tprop = propaLen * delayScint;
457 const double Tfly = flyTime;
458
459 double t = Trec;
460 if (timeCableDelay.isValid()) t -= Tcable;
461 t -= Tprop;
462
463 const double t0_est = t - Tfly;
464 if (!std::isfinite(t0_est)) continue;
465
466 // Weighted accumulation using calibrated sigma
467 const double sigma = getHitSigma(KLMElementNumbers::c_BKLM, d.getLayer(), false);
468 acc_stat_weighted(t0_est, sigma, sumW, sumWT);
469 }
470 }
471 }
472}
473
474/* BKLM RPC */
476 const ExtMap& rpcMap,
477 double& sumW, double& sumWT)
478{
479 DBObjPtr<KLMTimeConstants> timeConstants;
480 DBObjPtr<KLMTimeCableDelay> timeCableDelay;
481
482 const double delayPhi = timeConstants.isValid()
483 ? timeConstants->getDelay(KLMTimeConstants::c_RPCPhi)
484 : 0.0;
485 const double delayZ = timeConstants.isValid()
486 ? timeConstants->getDelay(KLMTimeConstants::c_RPCZ)
487 : 0.0;
488
489 for (KLMHit2d& hit2d : klmHit2ds) {
490 if (hit2d.getSubdetector() != KLMElementNumbers::c_BKLM) continue;
491 if (!hit2d.inRPC()) continue;
492 if (hit2d.isOutOfTime()) continue;
493
494 RelationVector<BKLMHit1d> b1ds = hit2d.getRelationsTo<BKLMHit1d>();
495 if (b1ds.size() == 0) continue;
496
497 const bklm::Module* mod = m_geoParB->findModule(hit2d.getSection(), hit2d.getSector(), hit2d.getLayer());
498 const ROOT::Math::XYZVector posG2d = hit2d.getPosition();
499
500 for (const BKLMHit1d& h1d : b1ds) {
501 const bool isPhi = h1d.isPhiReadout();
502 RelationVector<KLMDigit> digits = h1d.getRelationsTo<KLMDigit>();
503
504 for (const KLMDigit& d : digits) {
505 if (!d.inRPC()) continue;
506
507 unsigned int cid = d.getUniqueChannelID();
508 if (m_channelStatus.isValid() &&
509 m_channelStatus->getChannelStatus(cid) != KLMChannelStatus::c_Normal) continue;
510
511 if (!d.isGood()) continue;
512
513 // RPC matched by module key
514 unsigned int moduleKey = m_elementNum->moduleNumber(d.getSubdetector(), d.getSection(), d.getSector(), d.getLayer());
515 ExtPair p = const_cast<KLMEventT0EstimatorModule*>(this)->matchExt(moduleKey, const_cast<ExtMap&>(rpcMap));
516 if (!p.first || !p.second) continue;
517
518 const double flyTime = 0.5 * (p.first->getTOF() + p.second->getTOF());
519 const ROOT::Math::XYZVector posGext = 0.5 * (p.first->getPosition() + p.second->getPosition());
520
521 const CLHEP::Hep3Vector locExt = mod->globalToLocal(CLHEP::Hep3Vector(posGext.X(), posGext.Y(), posGext.Z()), true);
522 const CLHEP::Hep3Vector locHit2 = mod->globalToLocal(CLHEP::Hep3Vector(posG2d.X(), posG2d.Y(), posG2d.Z()), true);
523
524 const CLHEP::Hep3Vector diff = locExt - locHit2;
525 if (std::fabs(diff.z()) > mod->getZStripWidth() || std::fabs(diff.y()) > mod->getPhiStripWidth()) continue;
526
527 const CLHEP::Hep3Vector propaV = mod->getPropagationDistance(locExt);
528 const double propaDist = isPhi ? propaV.y() : propaV.z();
529
530 // Components
531 const double Trec = d.getTime();
532 const double Tcable = timeCableDelay.isValid() ? timeCableDelay->getTimeDelay(cid) : 0.0;
533 const double Tprop = propaDist * (isPhi ? delayPhi : delayZ);
534 const double Tfly = flyTime;
535
536 double t = Trec;
537 if (timeCableDelay.isValid()) t -= Tcable;
538 t -= Tprop;
539
540 const double t0_est = t - Tfly;
541 if (!std::isfinite(t0_est)) continue;
542
543 // Weighted accumulation using calibrated sigma with direction-specific resolution
545 const double sigma = getHitSigma(KLMElementNumbers::c_BKLM, d.getLayer(), true, plane);
546 acc_stat_weighted(t0_est, sigma, sumW, sumWT);
547 }
548 }
549 }
550}
551
552/* Event. */
553
555{
556 // CDC seed for logging (not used to compute the means)
557 m_seedT0 = 0.0;
560 if (evtT0.isValid() && evtT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
561 const auto best = evtT0->getBestCDCTemporaryEventT0();
562 if (best) m_seedT0 = best->eventT0;
563 }
564 }
565
566 if (!m_MuonList.isValid()) { B2WARNING("KLMEventT0Estimator: ParticleList '" << m_MuonListName << "' not found."); return; }
567 const unsigned nTracks = m_MuonList->getListSize();
568 if (nTracks == 0u) return;
569
570 // Weighted mean and uncertainty from running sums
571 auto weighted_result = [](double wsum, double wtsum) -> std::pair<double, double> {
572 if (wsum <= 0.0) return {NAN, NAN};
573 return {wtsum / wsum, std::sqrt(1.0 / wsum)};
574 };
575
576 // Weighted track averaging using inverse-variance (1/SEM²) weighting
577 auto mean_sem_tracks = [](const std::vector<std::pair<double, double>>& v) -> std::pair<double, double> {
578 if (v.empty()) return {NAN, NAN};
579 if (v.size() == 1)
580 {
581 return {v[0].first, std::isfinite(v[0].second) ? v[0].second : 0.0};
582 }
583
584 bool allValid = true;
585 for (const auto& [t0, sem] : v)
586 {
587 if (!std::isfinite(sem) || sem <= 0.0) { allValid = false; break; }
588 }
589
590 if (allValid)
591 {
592 double wsum = 0.0, wtsum = 0.0;
593 for (const auto& [t0, sem] : v) {
594 const double w = 1.0 / (sem * sem);
595 wsum += w;
596 wtsum += w * t0;
597 }
598 if (wsum > 0.0) {
599 return {wtsum / wsum, std::sqrt(1.0 / wsum)};
600 }
601 }
602
603 // Fallback to simple average if weights not valid
604 double s = 0.0;
605 for (const auto& [t0, sem] : v) s += t0;
606 const double mu = s / v.size();
607 double ss = 0.0;
608 for (const auto& [t0, sem] : v) { const double d = t0 - mu; ss += d * d; }
609 const double var = (v.size() > 1) ? ss / (v.size() - 1) : 0.0;
610 return {mu, std::sqrt(var / v.size())};
611 };
612
613 // For per-event track-averages: pairs of (T0, SEM) for weighted averaging
614 std::vector<std::pair<double, double>> vTrk_B, vTrk_R, vTrk_E, vTrk_All;
615
616 for (unsigned i = 0; i < nTracks; ++i) {
617 const Particle* particle = m_MuonList->getParticle(i);
618 if (!particle) continue;
619 const Track* track = particle->getTrack();
620 if (!track) continue;
621
622 RelationVector<KLMHit2d> hit2ds = track->getRelationsTo<KLMHit2d>();
623 if (hit2ds.size() == 0) continue;
624
625 // Build ExtHit maps for this track
626 m_extScint.clear();
627 m_extRPC.clear();
629
630 // Per-track digit sums per category (weighted)
631 double wE = 0, wTE = 0;
632 accumulateEKLM(hit2ds, m_extScint, wE, wTE);
633
634 double wB = 0, wTB = 0;
635 accumulateBKLMScint(hit2ds, m_extScint, wB, wTB);
636
637 double wR = 0, wTR = 0;
638 accumulateBKLMRPC(hit2ds, m_extRPC, wR, wTR);
639
640 // Per-track means and SEMs by category
641 if (wB > 0.0) {
642 auto [muB, seB] = weighted_result(wB, wTB);
643 if (m_hT0Trk_BKLM_Scint && std::isfinite(muB)) m_hT0Trk_BKLM_Scint->Fill(muB);
644 if (std::isfinite(muB)) vTrk_B.push_back({muB, seB});
645 }
646
647 if (wR > 0.0) {
648 auto [muR, seR] = weighted_result(wR, wTR);
649 if (m_hT0Trk_BKLM_RPC && std::isfinite(muR)) m_hT0Trk_BKLM_RPC->Fill(muR);
650 if (std::isfinite(muR)) vTrk_R.push_back({muR, seR});
651 }
652
653 if (wE > 0.0) {
654 auto [muE, seE] = weighted_result(wE, wTE);
655 if (m_hT0Trk_EKLM_Scint && std::isfinite(muE)) m_hT0Trk_EKLM_Scint->Fill(muE);
656 if (std::isfinite(muE)) vTrk_E.push_back({muE, seE});
657 }
658
659 // Per-track overall (if any category present)
660 {
661 const double wAll = wB + wE + wR;
662 const double wtAll = wTB + wTE + wTR;
663 if (wAll > 0.0) {
664 const double t0 = wtAll / wAll;
665 const double se = std::sqrt(1.0 / wAll);
666 vTrk_All.push_back({t0, se});
667 }
668 }
669 }
670
671 if (vTrk_All.empty()) {
672 B2DEBUG(20, "KLMEventT0Estimator: no usable KLM timing residuals for this event.");
673 return;
674 }
675
676 // Per-event track-averages using inverse-variance (1/SEM²) weighting
677 const auto [muB_trk, seB_trk] = mean_sem_tracks(vTrk_B);
678 const auto [muR_trk, seR_trk] = mean_sem_tracks(vTrk_R);
679 const auto [muE_trk, seE_trk] = mean_sem_tracks(vTrk_E);
680 const auto [muAll_trk, seAll_trk] = mean_sem_tracks(vTrk_All);
681
682 if (m_hT0Evt_TrkAvg_BKLM_Scint && std::isfinite(muB_trk)) m_hT0Evt_TrkAvg_BKLM_Scint->Fill(muB_trk);
683 if (m_hT0Evt_TrkAvg_BKLM_RPC && std::isfinite(muR_trk)) m_hT0Evt_TrkAvg_BKLM_RPC->Fill(muR_trk);
684 if (m_hT0Evt_TrkAvg_EKLM_Scint && std::isfinite(muE_trk)) m_hT0Evt_TrkAvg_EKLM_Scint->Fill(muE_trk);
685 if (m_hT0Evt_TrkAvg_All && std::isfinite(muAll_trk)) m_hT0Evt_TrkAvg_All->Fill(muAll_trk);
686
687 if (m_hT0Evt_TrkAvg_BKLM_Scint_SEM && std::isfinite(seB_trk)) m_hT0Evt_TrkAvg_BKLM_Scint_SEM->Fill(seB_trk);
688 if (m_hT0Evt_TrkAvg_BKLM_RPC_SEM && std::isfinite(seR_trk)) m_hT0Evt_TrkAvg_BKLM_RPC_SEM->Fill(seR_trk);
689 if (m_hT0Evt_TrkAvg_EKLM_Scint_SEM && std::isfinite(seE_trk)) m_hT0Evt_TrkAvg_EKLM_Scint_SEM->Fill(seE_trk);
690 if (m_hT0Evt_TrkAvg_All_SEM && std::isfinite(seAll_trk)) m_hT0Evt_TrkAvg_All_SEM->Fill(seAll_trk);
691
692 // ---------------- Final KLM combination (single saved component) ----------------
693 const bool useB = std::isfinite(muB_trk);
694 const bool useE = std::isfinite(muE_trk);
695 const bool useR = std::isfinite(muR_trk);
696
697 double finalT0 = NAN, finalSE = NAN;
698 int sourceBin = -1;
699
700 {
701 std::vector<std::pair<double, double>> parts;
702 if (useB) parts.emplace_back(muB_trk, seB_trk);
703 if (useE) parts.emplace_back(muE_trk, seE_trk);
704 if (useR) parts.emplace_back(muR_trk, seR_trk);
705
706 auto [t0, se] = mean_sem_tracks(parts);
707 finalT0 = t0;
708 finalSE = se;
709
710 if (useB && useE && useR) sourceBin = 7;
711 else if (useB && useE) sourceBin = 4;
712 else if (useB && useR) sourceBin = 5;
713 else if (useE && useR) sourceBin = 6;
714 else if (useB) sourceBin = 1;
715 else if (useE) sourceBin = 2;
716 else if (useR) sourceBin = 3;
717 }
718
719 B2DEBUG(20, "KLMEventT0Estimator: "
720 << "T0_trkavg_all=" << muAll_trk << " ns (seed CDC=" << m_seedT0 << " ns)"
721 << " | E=" << muE_trk << " | Bsc=" << muB_trk << " | Brpc=" << muR_trk
722 << (std::isfinite(finalT0) ? (std::string(" | FINAL KLM=") + std::to_string(finalT0) + " ns") : std::string("")));
723
725 if (!outT0.isValid()) outT0.construct();
726
727 if (std::isfinite(finalT0)) {
728 if (m_hFinalSource && sourceBin > 0) m_hFinalSource->Fill(sourceBin);
729 const double quality = static_cast<double>((useB ? 1 : 0) + (useE ? 1 : 0) + (useR ? 1 : 0));
730 EventT0::EventT0Component klmT0Component(finalT0, std::isfinite(finalSE) ? finalSE : 0.0,
731 Const::KLM, "KLM", quality);
732 outT0->addTemporaryEventT0(klmT0Component);
733 outT0->setEventT0(klmT0Component);
734 }
735}
static void channelNumberToElementNumbers(KLMChannelNumber channel, int *section, int *sector, int *layer, int *plane, int *strip)
Get element numbers by channel number.
@ c_FirstRPCLayer
First RPC layer.
Store one reconstructed BKLM 1D hit as a ROOT object.
Definition BKLMHit1d.h:30
Class for accessing objects in the database.
Definition DBObjPtr.h:21
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition DataStore.h:59
static const EKLMElementNumbers & Instance()
Instantiation.
void stripNumberToElementNumbers(int stripGlobal, int *section, int *layer, int *sector, int *plane, int *strip) const
Get element numbers by strip global number.
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Transformation data.
@ c_None
Displacement is not used.
Store one Ext hit as a ROOT object.
Definition ExtHit.h:31
double getTOF() const
Get time of flight from the point of closest approach near the origin to this hit.
Definition ExtHit.h:136
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition HistoModule.h:29
@ c_Normal
Normally operating channel.
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition KLMDigit.h:29
Belle2::EKLM::TransformData * m_transformE
EKLM strip transformation data.
TH1D * m_hT0Evt_TrkAvg_BKLM_RPC_SEM
Per-event T0 track-average for BKLM RPC (SEM) [ns].
TH1D * m_hT0Evt_TrkAvg_EKLM_Scint_SEM
Per-event T0 track-average for EKLM scintillator (SEM) [ns].
TH1D * m_hT0Trk_BKLM_RPC
Per-track T0 for BKLM RPC [ns].
double m_ADCCut_BKLM_Scint_Min
Minimum ADC cut for BKLM scintillator.
TH1D * m_hT0Trk_EKLM_Scint
Per-track T0 for EKLM scintillator [ns].
TH1I * m_hFinalSource
Final EventT0 source selection (7 bins).
std::string m_histDirName
Parent directory inside the ROOT file (HistoManager) for this module.
void accumulateEKLM(const RelationVector< KLMHit2d > &, const ExtMap &, double &sumW, double &sumWT)
Accumulate EKLM scintillator per-digit T0 estimates (weighted).
ExtPair matchExt(unsigned int key, ExtMap &v_ExtHits)
Find earliest (entry) and latest (exit) ExtHits matching a key (channel or module).
const Belle2::EKLM::GeometryData * m_geoParE
EKLM geometry data.
TH1D * m_hT0Trk_BKLM_Scint
Per-track T0 for BKLM scintillator [ns].
bool m_useCDCTemporaryT0
Use CDC temporary EventT0 as a diagnostic seed (not applied to averaging).
void initialize() override
Register inputs/params; get geometry; call REG_HISTOGRAM.
void event() override
Per-event algorithm: collect hits, compute residuals, fill outputs.
std::string m_histSubdirUncorr
Subdirectory name for uncorrected timing histograms.
DBObjPtr< KLMEventT0HitResolution > m_eventT0HitResolution
Per-hit time resolution for EventT0 estimation.
std::multimap< unsigned int, Belle2::ExtHit > ExtMap
Multimap of ExtHit objects keyed by channel or module number.
void endRun() override
Called when the current run ends.
ExtMap m_extRPC
Extrapolated hits keyed by module number (RPC).
void terminate() override
Called at the end of processing.
TH1D * m_hT0Evt_TrkAvg_BKLM_Scint_SEM
Per-event T0 track-average for BKLM scintillator (SEM) [ns].
void collectExtrapolatedHits(const Track *track, ExtMap &scintMap, ExtMap &rpcMap)
Build maps of extrapolated hits for a track (scint: channel key; RPC: module key).
double m_ADCCut_BKLM_Scint_Max
Maximum ADC cut for BKLM scintillator.
DBObjPtr< KLMChannelStatus > m_channelStatus
Channel status (Normal/Dead/etc.).
StoreObjPtr< ParticleList > m_MuonList
Selected muon particle list.
void beginRun() override
Per-run resets if desired (histos remain booked).
StoreArray< Track > m_tracks
Reconstructed tracks.
double m_ADCCut_EKLM_Scint_Max
Maximum ADC cut for EKLM scintillator.
bool passesADCCut(double charge, int subdetector, int layer, bool inRPC) const
Check if a digit passes the ADC charge cut.
double m_ADCCut_EKLM_Scint_Min
Minimum ADC cut for EKLM scintillator.
Belle2::bklm::GeometryPar * m_geoParB
BKLM geometry.
bool m_ignoreBackward
Ignore backward-propagated ExtHits when forming entry/exit pairs.
TH1D * m_hT0Evt_TrkAvg_BKLM_Scint
Per-event T0 track-average for BKLM scintillator (mean) [ns].
std::pair< Belle2::ExtHit *, Belle2::ExtHit * > ExtPair
Pair of entry and exit ExtHit pointers.
double m_seedT0
Optional seed from CDC (for logging only).
std::string m_MuonListName
Input ParticleList (e.g.
double getHitSigma(int subdetector, int layer, bool inRPC, int plane=0) const
Get per-hit sigma for a digit based on detector category.
TH1D * m_hT0Evt_TrkAvg_All
Per-event T0 track-average combined (mean) [ns].
const KLMElementNumbers * m_elementNum
Element numbering helpers.
void accumulateBKLMRPC(RelationVector< KLMHit2d > &klmHit2ds, const ExtMap &rpcMap, double &sumW, double &sumWT)
Accumulate BKLM RPC per-digit T0 estimates (weighted, both readout directions).
TH1D * m_hT0Evt_TrkAvg_All_SEM
Per-event T0 track-average combined (SEM) [ns].
TH1D * m_hT0Evt_TrkAvg_BKLM_RPC
Per-event T0 track-average for BKLM RPC (mean) [ns].
ExtMap m_extScint
Extrapolated hits keyed by channel number (scintillator).
TH1D * m_hT0Evt_TrkAvg_EKLM_Scint
Per-event T0 track-average for EKLM scintillator (mean) [ns].
void accumulateBKLMScint(RelationVector< KLMHit2d > &, const ExtMap &, double &sumW, double &sumWT)
Accumulate BKLM scintillator per-digit T0 estimates (weighted).
void defineHisto() override
Definition of histograms (called once by HistoManager).
KLM 2d hit.
Definition KLMHit2d.h:33
Class to store the likelihoods from KLM with additional information related to the extrapolation.
@ c_BKLM
BKLM scintillator.
@ c_EKLM
EKLM scintillator.
Class to store reconstructed particles.
Definition Particle.h:76
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
Type-safe access to single objects in the data store.
Definition StoreObjPtr.h:96
bool isValid() const
Check whether the object was created.
bool construct(Args &&... params)
Construct an object of type T in this StoreObjPtr, using the provided constructor arguments.
Class that bundles various TrackFitResults.
Definition Track.h:25
static const double mm
[millimeters]
Definition Unit.h:70
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition Module.h:76
Class to store variables with their name which were sent to the logging service.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
Abstract base class for different kinds of events.
Structure for storing the extracted event t0s together with its detector and its uncertainty.
Definition EventT0.h:33