Belle II Software  release-06-01-15
TOPDigitVariables.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 #include <top/variables/TOPDigitVariables.h>
9 
10 // // framework
11 #include <framework/gearbox/Const.h>
12 #include <framework/datastore/StoreObjPtr.h>
13 #include <framework/datastore/StoreArray.h>
14 
15 // // dataobjects
16 #include <tracking/dataobjects/ExtHit.h>
17 #include <mdst/dataobjects/Track.h>
18 #include <analysis/dataobjects/Particle.h>
19 #include <top/dataobjects/TOPDigit.h>
20 #include <top/geometry/TOPGeometryPar.h>
21 #include <top/dataobjects/TOPLikelihood.h>
22 #include <top/dataobjects/TOPRecBunch.h>
23 #include <top/dataobjects/TOPLikelihoodScanResult.h>
24 
25 #include <top/dataobjects/TOPBarHit.h>
26 #include <mdst/dataobjects/MCParticle.h>
27 
28 #include <top/reconstruction_cpp/TOPRecoManager.h>
29 #include <top/reconstruction_cpp/TOPTrack.h>
30 #include <top/reconstruction_cpp/PDFConstructor.h>
31 
32 
33 #include <algorithm> // for sort
34 using namespace std;
35 
36 namespace Belle2 {
41  // contains a couple of helper functions that are related to TOP variables
42  namespace Variable {
43  namespace TOPVariable {
44  // returns the TOP likelihood that is associated with a given particle
45  const TOPLikelihood* getTOPLikelihood(const Particle* particle)
46  {
47  if (not particle) return nullptr;
48  const auto* track = particle->getTrack();
49  return track ? track->getRelated<TOPLikelihood>() : nullptr;
50  }
51 
52  // returns the ExtHit that is associated with a given particle
53  const ExtHit* getExtHit(const Particle* particle)
54  {
55  const auto* topLikelihood = getTOPLikelihood(particle);
56  return topLikelihood ? topLikelihood->getRelated<ExtHit>() : nullptr;
57  }
58 
59  // returns the TOP slot ID of the particle
60  double getSlotID(const Particle* particle)
61  {
62  const auto* extHit = getExtHit(particle);
63  // zero is invalid slot ID, e.g. particle didn't hit the TOP
64  return extHit ? extHit->getCopyID() : 0;
65  }
66 
67  // returns the local coordinate of the particle's entry point to the TOP
68  TVector3 getLocalPosition(const Particle* particle)
69  {
70  const auto* extHit = getExtHit(particle);
71  if (not extHit) return TVector3(0, 0, 0);
72  int slotID = extHit->getCopyID();
73  const auto& position = extHit->getPosition(); // TVector3
74  const auto* geo = TOP::TOPGeometryPar::Instance()->getGeometry();
75  if (not geo or not geo->isModuleIDValid(slotID)) return TVector3(0, 0, 0);
76  const auto& module = geo->getModule(slotID);
77  return module.pointToLocal(position); // TVector3
78  }
79 
80  // returns the local coordinate of the MC particle's entry point to the TOP
81  TVector3 getLocalPositionMCMatch(const Particle* particle)
82  {
83  const MCParticle* mcparticle = particle->getRelatedTo<MCParticle>();
84  if (mcparticle == nullptr) {
85  return TVector3(0, 0, 0);
86  }
87  const auto* barHit = mcparticle->getRelated<TOPBarHit>();
88  if (!barHit) {
89  return TVector3(0, 0, 0);
90  }
91  int slotID = barHit->getModuleID();
92  const auto& position = barHit->getPosition(); // TVector3
93  const auto* geo = TOP::TOPGeometryPar::Instance()->getGeometry();
94  if (not geo or not geo->isModuleIDValid(slotID)) return TVector3(0, 0, 0);
95  const auto& module = geo->getModule(slotID);
96  return module.pointToLocal(position); // TVector3
97  }
98 
99  // returns the local coordinates of the particles momentum in the TOP
100  TVector3 getLocalMomentum(const Particle* particle)
101  {
102  const auto* extHit = getExtHit(particle);
103  if (not extHit) return TVector3(0, 0, 0);
104  int slotID = extHit->getCopyID();
105  const auto& momentum = extHit->getMomentum(); // TVector3
106  const auto* geo = TOP::TOPGeometryPar::Instance()->getGeometry();
107  if ((not geo) or (not geo->isModuleIDValid(slotID))) return TVector3(0, 0, 0);
108  const auto& module = geo->getModule(slotID);
109  return module.momentumToLocal(momentum); // TVector3
110  }
111 
112  // helper function to compute the TOF for an arbitrary hypothesis
113  double computeTOF(const Particle* particle, int pdg)
114  {
115  const auto* extHit = getExtHit(particle);
116  if (not extHit) return 0;
117  auto extPDGCode = abs(extHit->getPdgCode());
118  double pmom = particle->getMomentumMagnitude();
119  double massExtHit = Const::ChargedStable(extPDGCode).getMass();
120  double betaExtHit = pmom / sqrt(pmom * pmom + massExtHit * massExtHit);
121  double mass = pdg == 0 ? particle->getMass() : Const::ChargedStable(abs(pdg)).getMass();
122  double beta = pmom / sqrt(pmom * pmom + mass * mass);
123  return extHit->getTOF() * betaExtHit / beta;
124  }
125 
126  // returns the time of flight from the origin to the TOP
127  double getTOF(const Particle* particle)
128  {
129  return computeTOF(particle, 0);
130  }
131 
132  // returns the time of flight from the origin to the TOP under a given hypothesis
133  double getTOFExpert(const Particle* particle, const vector<double>& vars)
134  {
135  if (vars.size() != 1) {
136  B2FATAL("Need exactly one parameter (pdg id).");
137  }
138  int pdg = static_cast<int>(vars[0]);
139  return computeTOF(particle, pdg);
140  }
141 
142  // returns the average time of the first 5 (good) digits
143  double getAverageTimeOfFirst5(const Particle* particle)
144  {
145  int slotID = static_cast<int>(getSlotID(particle));
146  StoreArray<TOPDigit> digits;
147  vector<double> digitTimes;
148  for (const auto& digit : digits) {
149  if (digit.getModuleID() != slotID) continue;
150  // skip bad digits only when we want to clean
151  if (digit.getHitQuality() != TOPDigit::c_Good) continue;
152  digitTimes.push_back(digit.getTime());
153  }
154  if (digitTimes.empty()) return 0;
155  sort(digitTimes.begin(), digitTimes.end());
156  double T0 = 0;
157  size_t count = 0;
158  for (auto t : digitTimes) {
159  T0 += t;
160  count += 1;
161  if (count == 5) break;
162  }
163  return T0 / count;
164  }
165 
166  // counts the number of photons in the TOP in a given time frame
167  // if tmin < 0, count from the time of the first photon
168  int countHits(const Particle* particle, double tmin, double tmax, bool clean)
169  {
170  int slotID = static_cast<int>(getSlotID(particle));
171  StoreArray<TOPDigit> digits;
172  vector<double> digitTimes;
173  for (const auto& digit : digits) {
174  if (digit.getModuleID() != slotID) continue;
175  // skip bad digits only when we want to clean
176  if (clean && digit.getHitQuality() != TOPDigit::c_Good) continue;
177  digitTimes.push_back(digit.getTime());
178  }
179  if (digitTimes.empty()) return 0;
180  sort(digitTimes.begin(), digitTimes.end());
181  int count = 0;
182  if (tmin < 0) tmin = digitTimes[0];
183  for (auto t : digitTimes) {
184  if (t > tmax) break;
185  if (t >= tmin) ++count;
186  }
187  return count;
188  }
189 
190  // counts the number of photons regardless of hit quality
191  int countRawHits(const Particle* particle, double tmin, double tmax)
192  {
193  return countHits(particle, tmin, tmax, false);
194  }
195 
196  // returns the expected number of photons for a given hypothesis
197  double getExpectedPhotonCount(const Particle* particle, int pdg)
198  {
199  const auto* topLikelihood = getTOPLikelihood(particle);
200  if (not topLikelihood) return 0;
201  // if the user does select a hypothesis, use the particle's pdg code
202  pdg = pdg != 0 ? pdg : particle->getPDGCode();
203  const auto& chargedStable = Const::chargedStableSet.find(abs(pdg));
204  if (chargedStable == Const::invalidParticle) return 0; // PDG code not one of e, mu, pi, K, p, d
205  return topLikelihood->getEstPhot(chargedStable);
206  }
208  double topDigitCount(const Particle* particle)
209  {
210  auto trk = particle->getTrack();
211  if (not trk) {
212  return -1.0;
213  }
214  auto extHits = trk->getRelationsWith<ExtHit>();
215  int thisModuleID = static_cast<int>(getSlotID(particle));
216  if (thisModuleID == 0) return 0;
217  StoreArray<TOPDigit> topDigits;
218  int count = 0;
219  for (const auto& t : topDigits) {
220  if (t.getModuleID() != thisModuleID) continue; // catch the case where one of the module IDs is negative
221  if (t.getHitQuality() != TOPDigit::c_Good) continue;
222  count += 1;
223  }
224  return count;
225  }
226 
228  double topBackgroundDigitCount(const Particle* particle)
229  {
230  auto trk = particle->getTrack();
231  if (not trk) {
232  return -1.0;
233  }
234  auto extHits = trk->getRelationsWith<ExtHit>();
235  int thisModuleID = static_cast<int>(getSlotID(particle));
236  if (thisModuleID == 0) return 0;
237  StoreArray<TOPDigit> topDigits;
238  int count = 0;
239  for (const auto& t : topDigits) {
240  if (abs(t.getModuleID()) == abs(thisModuleID)) continue; // catch the case where one of the module IDs is negative
241  if (t.getHitQuality() != TOPDigit::c_Good) continue;
242  count += 1;
243  }
244  return count;
245  }
246 
248  double topBackgroundDigitCountRaw(const Particle* particle)
249  {
250  auto trk = particle->getTrack();
251  if (not trk) {
252  return -1.0;
253  }
254  auto extHits = trk->getRelationsWith<ExtHit>();
255  int thisModuleID = static_cast<int>(getSlotID(particle));
256  if (thisModuleID == 0) return 0;
257  StoreArray<TOPDigit> topDigits;
258  int count = 0;
259  for (const auto& t : topDigits) {
260  if (abs(t.getModuleID()) == abs(thisModuleID)) continue; // catch the case where one of the module IDs is negative
261  count += 1;
262  }
263  return count;
264  }
265 
267  double topRawDigitCount(const Particle* particle)
268  {
269  auto trk = particle->getTrack();
270  if (not trk) {
271  return -1.0;
272  }
273  auto extHits = trk->getRelationsWith<ExtHit>();
274  int thisModuleID = static_cast<int>(getSlotID(particle));
275  if (thisModuleID == 0) return 0;
276  StoreArray<TOPDigit> topDigits;
277  int count = 0;
278  for (const auto& t : topDigits) {
279  if (abs(t.getModuleID()) != abs(thisModuleID)) continue; // catch the case where one of the module IDs is negative
280  count += 1;
281  }
282  return count;
283  }
284 
286  double topDigitGapSize(const Particle* particle)
287  {
288  auto trk = particle->getTrack();
289  if (not trk) {
290  return -1.0;
291  }
292  auto extHits = trk->getRelationsWith<ExtHit>();
293  int thisModuleID = static_cast<int>(getSlotID(particle));
294  if (thisModuleID == 0) return 0;
295  StoreArray<TOPDigit> topDigits;
296  double maxGap = 0; // the largest time difference between two consecutive hits
297  vector<double> digitTimes; // all digits in the module that the track entered
298  for (const auto& t : topDigits) {
299  if (abs(t.getModuleID()) != abs(thisModuleID)) continue;
300  if (t.getHitQuality() != TOPDigit::c_Good) continue;
301  digitTimes.push_back(t.getTime());
302  }
303  if (digitTimes.empty()) {
304  return -1.0;
305  }
306  sort(digitTimes.begin(), digitTimes.end());
307  for (size_t i = 0; i < digitTimes.size() - 1; ++i) {
308  double gap = digitTimes[i + 1] - digitTimes[i];
309  if (gap > maxGap) {
310  maxGap = gap;
311  }
312  }
313  return maxGap;
314  }
315 
316  // The number of reflected digits is defined as the number of digits after the gap
317  // This method is a helper function to count the top digits after the largest gap
318  // between subsequent hits, under the constraints gap > minGap and gap < maxGap
319  double topCountPhotonsAfterLargesGapWithin(const Particle* particle, double minGap, double maxGap)
320  {
321  auto trk = particle->getTrack();
322  if (not trk) {
323  return -1.0;
324  }
325  auto extHits = trk->getRelationsWith<ExtHit>();
326  int thisModuleID = static_cast<int>(getSlotID(particle));
327  if (thisModuleID == 0) return 0;
328  StoreArray<TOPDigit> topDigits;
329  vector<double> digitTimes; // the times for all digits in the module that the track entered
330  for (const auto& t : topDigits) {
331  if (abs(t.getModuleID()) != abs(thisModuleID)) continue;
332  if (t.getHitQuality() != TOPDigit::c_Good) continue;
333  digitTimes.push_back(t.getTime());
334  }
335  if (digitTimes.empty()) {
336  return -1.0;
337  }
338 
339  double currentMaxGap = -1; // the largest time difference between two consecutive hits
340  size_t maxGapIndex = 0; // the index of the first hit *after* the gap
341  sort(digitTimes.begin(), digitTimes.end());
342  for (size_t i = 0; i < digitTimes.size() - 1; ++i) {
343  double gap = digitTimes[i + 1] - digitTimes[i];
344  if ((gap > minGap) and (gap < maxGap) and (gap > currentMaxGap)) {
345  currentMaxGap = gap;
346  maxGapIndex = i + 1;
347  }
348  }
349  return digitTimes.size() - maxGapIndex;
350  }
351 
353  double extrapTrackToTOPz(const Particle* particle)
354  {
355  auto trk = particle->getTrack();
356  if (not trk) {
357  return std::numeric_limits<double>::quiet_NaN();
358  }
359  auto trkfit = trk->getTrackFitResultWithClosestMass(Belle2::Const::ChargedStable(std::abs(particle->getPDGCode())));
360  auto top = trkfit->getHelix();
361  double arcLength = top.getArcLength2DAtCylindricalR(120);
362  const auto& result = top.getPositionAtArcLength2D(arcLength);
363  return result.z();
364  }
365 
367  double extrapTrackToTOPtheta(const Particle* particle)
368  {
369  auto trk = particle->getTrack();
370  if (not trk) {
371  return std::numeric_limits<double>::quiet_NaN();
372  }
373  auto trkfit = trk->getTrackFitResultWithClosestMass(Belle2::Const::ChargedStable(std::abs(particle->getPDGCode())));
374  auto top = trkfit->getHelix();
375  double arcLength = top.getArcLength2DAtCylindricalR(120);
376  const auto& result = top.getPositionAtArcLength2D(arcLength);
377  return result.Theta();
378  }
379 
381  double extrapTrackToTOPphi(const Particle* particle)
382  {
383  auto trk = particle->getTrack();
384  if (not trk) {
385  return std::numeric_limits<double>::quiet_NaN();
386  }
387  auto trkfit = trk->getTrackFitResultWithClosestMass(Belle2::Const::ChargedStable(std::abs(particle->getPDGCode())));
388  auto top = trkfit->getHelix();
389  double arcLength = top.getArcLength2DAtCylindricalR(120);
390  const auto& result = top.getPositionAtArcLength2D(arcLength);
391  return result.Phi();
392  }
393 
395  double topReflectedDigitCount(const Particle* particle)
396  {
397  return topCountPhotonsAfterLargesGapWithin(particle, 0, 10000);
398  }
399 
400  double topReflectedDigitCountExpert(const Particle* particle, const vector<double>& vars)
401  {
402  if (vars.size() != 2) {
403  B2FATAL("Need exactly two parameters (min, max)");
404  }
405  return topCountPhotonsAfterLargesGapWithin(particle, vars[0], vars[1]);
406  }
407 
409  double getTOPLocalX(const Particle* particle)
410  {
411  return TOPVariable::getLocalPosition(particle).X();
412  }
413 
415  double getTOPLocalY(const Particle* particle)
416  {
417  return TOPVariable::getLocalPosition(particle).Y();
418  }
419 
421  double getTOPLocalZ(const Particle* particle)
422  {
423  return TOPVariable::getLocalPosition(particle).Z();
424  }
425 
427  double getTOPLocalXMCMatch(const Particle* particle)
428  {
429  return TOPVariable::getLocalPositionMCMatch(particle).X();
430  }
431 
433  double getTOPLocalYMCMatch(const Particle* particle)
434  {
435  return TOPVariable::getLocalPositionMCMatch(particle).Y();
436  }
437 
439  double getTOPLocalZMCMatch(const Particle* particle)
440  {
441  return TOPVariable::getLocalPositionMCMatch(particle).Z();
442  }
443 
445  double getTOPLocalPhi(const Particle* particle)
446  {
447  return TOPVariable::getLocalMomentum(particle).Phi();
448  }
449 
451  double getTOPLocalTheta(const Particle* particle)
452  {
453  return TOPVariable::getLocalMomentum(particle).Theta();
454  }
455 
457  double getTOPPhotonCount(const Particle* particle)
458  {
459  const auto* topLikelihood = TOPVariable::getTOPLikelihood(particle);
460  return topLikelihood ? topLikelihood->getNphot() : 0;
461  }
462 
464  double getExpectedTOPPhotonCount(const Particle* particle, const vector<double>& vars)
465  {
466  if (vars.size() != 1) {
467  B2FATAL("Need exactly one parameter (pdg id).");
468  }
469  return TOPVariable::getExpectedPhotonCount(particle, static_cast<int>(vars[0]));
470  }
471 
473  double countTOPHitsInInterval(const Particle* particle, const vector<double>& vars)
474  {
475  if (vars.size() != 2) {
476  B2FATAL("Need exactly two parameters (tmin, tmax)");
477  }
478  return TOPVariable::countHits(particle, vars[0], vars[1]);
479  }
480 
482  double countTOPHitsInFirst20ns(const Particle* particle)
483  {
484  return TOPVariable::countHits(particle, -1.0, 20.0);
485  }
486 
488  double countRawTOPHitsInInterval([[maybe_unused]] const Particle* particle, const vector<double>& vars)
489  {
490  if (vars.size() != 2) {
491  B2FATAL("Need exactly two parameters (tmin, tmax)");
492  }
493  return TOPVariable::countRawHits(particle, vars[0], vars[1]);
494  }
495 
496  double getFlag(const Particle* particle)
497  {
498  const auto* topLikelihood = getTOPLikelihood(particle);
499  if (not topLikelihood) return 0;
500  return topLikelihood->getFlag();
501  }
502 
503  double getElectronLogL(const Particle* particle)
504  {
505  const auto* topLikelihood = getTOPLikelihood(particle);
506  if (not topLikelihood) return 0;
507  return topLikelihood->getLogL_e();
508  }
509 
510  double getMuonLogL(const Particle* particle)
511  {
512  const auto* topLikelihood = getTOPLikelihood(particle);
513  if (not topLikelihood) return 0;
514  return topLikelihood->getLogL_mu();
515  }
516 
517  double getPionLogL(const Particle* particle)
518  {
519  const auto* topLikelihood = getTOPLikelihood(particle);
520  if (not topLikelihood) return 0;
521  return topLikelihood->getLogL_pi();
522  }
523 
524  double getKaonLogL(const Particle* particle)
525  {
526  const auto* topLikelihood = getTOPLikelihood(particle);
527  if (not topLikelihood) return 0;
528  return topLikelihood->getLogL_K();
529  }
530 
531  double getProtonLogL(const Particle* particle)
532  {
533  const auto* topLikelihood = getTOPLikelihood(particle);
534  if (not topLikelihood) return 0;
535  return topLikelihood->getLogL_p();
536  }
537 
538 
539  double getLogLScanMass(const Particle* particle)
540  {
541  const auto* track = particle->getTrack();
542  if (!track) return -1;
543  auto scanRes = track->getRelated<TOPLikelihoodScanResult>();
544  if (!scanRes) {
545  B2WARNING("No TOPLikelihoodScanResult objcte found. Are you sure you added TOPLLScanner to the path?");
546  return -1;
547  }
548  return scanRes->getMostLikelyMass();
549  }
550 
551  double getLogLScanMassUpperInterval(const Particle* particle)
552  {
553  const auto* track = particle->getTrack();
554  if (!track) return -1;
555  auto scanRes = track->getRelated<TOPLikelihoodScanResult>();
556  if (!scanRes) {
557  B2WARNING("No TOPLikelihoodScanResult object found. Are you sure you added TOPLLScanner to the path?");
558  return -1;
559  }
560  return scanRes->getMostLikelyMassIntervalUp();
561  }
562 
563 
564  double getLogLScanMassLowerInterval(const Particle* particle)
565  {
566  const auto* track = particle->getTrack();
567  if (!track) return -1;
568  auto scanRes = track->getRelated<TOPLikelihoodScanResult>();
569  if (!scanRes) {
570  B2WARNING("No TOPLikelihoodScanResult object found. Are you sure you added TOPLLScanner to the path?");
571  return -1;
572  }
573  return scanRes->getMostLikelyMassIntervalLow();
574  }
575 
576  double getLogLScanThreshold(const Particle* particle)
577  {
578  const auto* track = particle->getTrack();
579  if (!track) return -1;
580  auto scanRes = track->getRelated<TOPLikelihoodScanResult>();
581  if (!scanRes) {
582  B2WARNING("No TOPLikelihoodScanResult objcte found. Are you sure you added TOPLLScanner to the path?");
583  return -1;
584  }
585  return scanRes->getThreshold();
586  }
587 
588 
589  double getLogLScanExpectedSignalPhotons(const Particle* particle)
590  {
591  const auto* track = particle->getTrack();
592  if (!track) return -1;
593  auto scanRes = track->getRelated<TOPLikelihoodScanResult>();
594  if (!scanRes) {
595  B2WARNING("No TOPLikelihoodScanResult objcte found. Are you sure you added TOPLLScanner to the path?");
596  return -1;
597  }
598  return scanRes->getMostLikelySignalPhotonCount();
599  }
600 
601 
602  //---------------- TOPRecBunch related --------------------
603 
605  double isTOPRecBunchReconstructed([[maybe_unused]] const Particle* particle)
606  {
607  StoreObjPtr<TOPRecBunch> recBunch;
608  // Attention! 0.0 is false, everything else is true
609  // returning -1, like for the others will most likely lead to bugs
610  // if the caller is not careful about return values.
611  if (not recBunch.isValid()) return 0.0;
612  return recBunch->isReconstructed();
613  }
614 
616  double TOPRecBunchNumber([[maybe_unused]] const Particle* particle)
617  {
618  StoreObjPtr<TOPRecBunch> recBunch;
619  if (not recBunch.isValid()) return -9999.0;
620  return recBunch->getBunchNo();
621  }
622 
624  double TOPRecBunchCurrentOffset([[maybe_unused]] const Particle* particle)
625  {
626  StoreObjPtr<TOPRecBunch> recBunch;
627  if (not recBunch.isValid()) return -9999.0;
628  return recBunch->getCurrentOffset();
629  }
630 
632  double TOPRecBunchTrackCount([[maybe_unused]] const Particle* particle)
633  {
634  StoreObjPtr<TOPRecBunch> recBunch;
635  if (not recBunch.isValid()) return -9999.0;
636  return recBunch->getNumTracks();
637  }
638 
640  double TOPRecBunchUsedTrackCount([[maybe_unused]] const Particle* particle)
641  {
642  StoreObjPtr<TOPRecBunch> recBunch;
643  if (not recBunch.isValid()) return -9999.0;
644  return recBunch->getUsedTracks();
645  }
646  //-------------- Event based -----------------------------------
647 
649  double TOPRawPhotonsInSlot([[maybe_unused]] const Particle* particle, const vector<double>& vars)
650  {
651  if (vars.size() != 1) { B2FATAL("Need exactly one parameter (slot id).");}
652  StoreArray<TOPDigit> topDigits;
653  int thisModuleID = static_cast<int>(vars[0]);
654  size_t count = 0;
655  for (const auto& t : topDigits) {
656  if (abs(t.getModuleID()) != abs(thisModuleID)) continue;
657  count += 1;
658  }
659  return count;
660  }
661 
663  double TOPGoodPhotonsInSlot([[maybe_unused]] const Particle* particle, const vector<double>& vars)
664  {
665  if (vars.size() != 1) { B2FATAL("Need exactly one parameter (slot id).");}
666  StoreArray<TOPDigit> topDigits;
667  int thisModuleID = static_cast<int>(vars[0]);
668  size_t count = 0;
669  for (const auto& t : topDigits) {
670  if (abs(t.getModuleID()) != abs(thisModuleID)) continue;
671  if (t.getHitQuality() != TOPDigit::c_Good) continue;
672  count += 1;
673  }
674  return count;
675  }
676 
678  double TOPTracksInSlot([[maybe_unused]] const Particle* particle)
679  {
680  const auto* trk = particle->getTrack();
681  if (not trk) {
682  return -1.0;
683  }
684  int thisModuleID = static_cast<int>(getSlotID(particle));
685  if (thisModuleID == 0) return 0;
686  StoreArray<Track> tracks;
687  int nTracks = 0;
688  for (const auto& t : tracks) {
689  const auto* tl = t.getRelated<TOPLikelihood>();
690  if (not tl) continue;
691  const auto* te = tl->getRelated<ExtHit>();
692  if (not te) continue;
693  if (te->getCopyID() != thisModuleID) continue;
694  nTracks += 1;
695  }
696  return nTracks;
697  }
698 
699  } // TOPVariable
700 
701  VARIABLE_GROUP("TOP Calibration");
702  REGISTER_VARIABLE("extrapTrackToTOPimpactZ", TOPVariable::extrapTrackToTOPz,
703  "[calibration] z coordinate of the impact point of the track extrapolated to TOP using helix data from TrackFitResult");
704  REGISTER_VARIABLE("extrapTrackToTOPimpactTheta", TOPVariable::extrapTrackToTOPtheta,
705  "[calibration] theta coordinate of the impact point of the track extrapolated to TOP using helix data from TrackFitResult");
706  REGISTER_VARIABLE("extrapTrackToTOPimpactPhi", TOPVariable::extrapTrackToTOPphi,
707  "[calibration] phi coordinate of the impact point of the track extrapolated to TOP using helix data from TrackFitResult");
708  REGISTER_VARIABLE("topDigitCount", TOPVariable::topDigitCount,
709  "[calibration] The number of TOPDigits in the module to which the track was extrapolated");
710  REGISTER_VARIABLE("topBackgroundDigitCount", TOPVariable::topBackgroundDigitCount,
711  "[calibration] The number of TOPDigits in all modules except the one to which the track was extrapolated");
712  REGISTER_VARIABLE("topBackgroundDigitCountRaw", TOPVariable::topBackgroundDigitCountRaw,
713  "[calibration] The number of TOPDigits in all modules except the one to which the track was extrapolated");
714  REGISTER_VARIABLE("topDigitCountRaw", TOPVariable::topDigitCount,
715  "[calibration] The number of TOPDigits in the module to which the track was extrapolated, regardless of hit quality");
716  REGISTER_VARIABLE("topReflectedDigitCount", TOPVariable::topReflectedDigitCount,
717  "[calibration] The number of reflected photons in the same module");
718  REGISTER_VARIABLE("topReflectedDigitCountExpert(minGap, maxGap)", TOPVariable::topReflectedDigitCountExpert,
719  "[calibration] The number of photons after the largest gap between minGap and maxGap");
720  REGISTER_VARIABLE("topDigitGapSize", TOPVariable::topDigitGapSize,
721  "[calibration] The largest time difference between two consecutive hits in the same module");
722  REGISTER_VARIABLE("topLocalX", TOPVariable::getTOPLocalX,
723  "[calibration] The local x coordinate of the particle's entry point to the TOP module");
724  REGISTER_VARIABLE("topLocalY", TOPVariable::getTOPLocalY,
725  "[calibration] The local y coordinate of the particle's entry point to the TOP module");
726  REGISTER_VARIABLE("topLocalZ", TOPVariable::getTOPLocalZ,
727  "[calibration] The local z coordinate of the particle's entry point to the TOP module");
728  REGISTER_VARIABLE("topLocalXMCMatch", TOPVariable::getTOPLocalXMCMatch,
729  "[calibration] The local x coordinate of the MC particle's entry point to the TOP module");
730  REGISTER_VARIABLE("topLocalYMCMatch", TOPVariable::getTOPLocalYMCMatch,
731  "[calibration] The local y coordinate of the MC particle's entry point to the TOP module");
732  REGISTER_VARIABLE("topLocalZMCMatch", TOPVariable::getTOPLocalZMCMatch,
733  "[calibration] The local z coordinate of the MC particle's entry point to the TOP module");
734  REGISTER_VARIABLE("topLocalPhi", TOPVariable::getTOPLocalPhi,
735  "[calibration] The local phi coordinate of the particle's momentum in the TOP module");
736  REGISTER_VARIABLE("topLocalTheta", TOPVariable::getTOPLocalTheta,
737  "[calibration] The local phi coordinate of the particle's momentum in the TOP module");
738  REGISTER_VARIABLE("topTOF", TOPVariable::getTOF,
739  "[calibration] The time of flight from the origin to the TOP");
740  REGISTER_VARIABLE("topTOFExpert(pdg)", TOPVariable::getTOFExpert,
741  "[calibration] The time of flight from the origin to the TOP under the given hypothesis");
742  REGISTER_VARIABLE("topAverageTimeOfFirst5", TOPVariable::getAverageTimeOfFirst5,
743  "[calibration] The average time of the first (up to) 5 hits in the module with the track");
744  REGISTER_VARIABLE("topSlotID", TOPVariable::getSlotID,
745  "[calibration] The ID of the TOP slot that was hit by the particle");
746  REGISTER_VARIABLE("topExpectedPhotonCount(pdg)", TOPVariable::getExpectedTOPPhotonCount,
747  "[calibration] The expected number of photons in the TOP for the particle under the given hypothesis");
748  REGISTER_VARIABLE("topPhotonCount", TOPVariable::getTOPPhotonCount,
749  "[calibration] The number of (bg-subtracted) TOP photons in for the given particle");
750  REGISTER_VARIABLE("countTOPHitsInInterval(tmin, tmax)", TOPVariable::countTOPHitsInInterval,
751  "[calibration] The number of photons in the given interval");
752  REGISTER_VARIABLE("countTOPHitsInFirst20ns", TOPVariable::countTOPHitsInFirst20ns,
753  "[calibration] The number of photons in the first 20 ns after the first photon");
754  REGISTER_VARIABLE("countRawTOPHitsInInterval(tmin, tmax)", TOPVariable::countRawTOPHitsInInterval,
755  "[calibration] The number of photons in the given interval");
756  REGISTER_VARIABLE("topFlag", TOPVariable::getFlag,
757  "[calibration] reconstruction flag, log likelihoods are valid if flag==1");
758  REGISTER_VARIABLE("topElectronLogL", TOPVariable::getElectronLogL,
759  "[calibration] electron log likelihood");
760  REGISTER_VARIABLE("topMuonLogL", TOPVariable::getMuonLogL,
761  "[calibration] muon log likelihood");
762  REGISTER_VARIABLE("topPionLogL", TOPVariable::getPionLogL,
763  "[calibration] pion log likelihood");
764  REGISTER_VARIABLE("topKaonLogL", TOPVariable::getKaonLogL,
765  "[calibration] kaon log likelihood");
766  REGISTER_VARIABLE("topProtonLogL", TOPVariable::getProtonLogL,
767  "[calibration] proton log likelihood");
768  REGISTER_VARIABLE("logLScanMass", TOPVariable::getLogLScanMass,
769  "[calibration] mass at the logL maximum from the LL scan");
770  REGISTER_VARIABLE("logLScanMassUpperInterval", TOPVariable::getLogLScanMassUpperInterval,
771  "[calibration] Upper edge of the mass interval determined by the LL scan");
772  REGISTER_VARIABLE("logLScanMassLowerInterval", TOPVariable::getLogLScanMassLowerInterval,
773  "[calibration] Lower edge of the mass interval determined by the LL scan");
774  REGISTER_VARIABLE("logLScanThreshold", TOPVariable::getLogLScanThreshold,
775  "[calibration] Cherenkov threshold determind by the LL scan");
776  REGISTER_VARIABLE("logLScanExpectedSignalPhotons", TOPVariable::getLogLScanExpectedSignalPhotons,
777  "[calibration] Expected signal photon yeild at the LL maximum");
778  REGISTER_VARIABLE("topRecBunchUsedTrackCount", TOPVariable::TOPRecBunchUsedTrackCount,
779  "[calibration] The number of tracks used in the bunch reconstruction");
780  REGISTER_VARIABLE("topRecBunchTrackCount", TOPVariable::TOPRecBunchTrackCount,
781  "[calibration] The number of tracks in the TOP acceptance");
782  REGISTER_VARIABLE("topRecBunchCurrentOffset", TOPVariable::TOPRecBunchCurrentOffset,
783  "[calibration] The current offset");
784  REGISTER_VARIABLE("topRecBunchNumber", TOPVariable::TOPRecBunchNumber,
785  "[calibration] The number of the bunch relative to the interaction");
786  REGISTER_VARIABLE("isTopRecBunchReconstructed", TOPVariable::isTOPRecBunchReconstructed,
787  "[calibration] Flag to indicate whether the bunch was reconstructed");
788  REGISTER_VARIABLE("topRawPhotonsInSlot(id)", TOPVariable::TOPRawPhotonsInSlot,
789  "[calibration] The number of all photons in the given slot");
790  REGISTER_VARIABLE("topGoodPhotonsInSlot(id)", TOPVariable::TOPGoodPhotonsInSlot,
791  "[calibration] The number of good photons in the given slot");
792  REGISTER_VARIABLE("topTracksInSlot", TOPVariable::TOPTracksInSlot,
793  "[calibration] The number of tracks in the same slot as the particle");
794  } // Variable
796 } // Belle2
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:470
Abstract base class for different kinds of events.