Belle II Software development
ECLVariables.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 <analysis/variables/ECLVariables.h>
11
12//framework
13#include <framework/logging/Logger.h>
14
15//analysis
16#include <analysis/dataobjects/Particle.h>
17#include <analysis/dataobjects/ParticleList.h>
18#include <analysis/dataobjects/ECLEnergyCloseToTrack.h>
19#include <analysis/utility/ReferenceFrame.h>
20#include <analysis/ClusterUtility/ClusterUtils.h>
21#include <analysis/VariableManager/Utility.h>
22
23//MDST
24#include <mdst/dataobjects/KlId.h>
25#include <mdst/dataobjects/ECLCluster.h>
26#include <mdst/dataobjects/Track.h>
27#include <mdst/dataobjects/EventLevelClusteringInfo.h>
28
29#include <Math/Vector4D.h>
30
31#include <cmath>
32#include <stack>
33
34
35
36namespace Belle2 {
41 namespace Variable {
42 double distanceToMcKl(const Particle* particle)
43 {
44 if (particle->hasExtraInfo("mcdistanceKL")) {
45 return particle->getExtraInfo("mcdistanceKL");
46 } else {
47 B2WARNING("The extraInfo mcdistanceKL is not registered! \n"
48 "This variable is only available for ECL based lists, and you have to run the function getNeutralHadronGeomMatches to use it");
49 return Const::doubleNaN;
50 }
51 }
52
53 double distanceToMcNeutron(const Particle* particle)
54 {
55 if (particle->hasExtraInfo("mcdistanceNeutron")) {
56 return particle->getExtraInfo("mcdistanceNeutron");
57 } else {
58 B2WARNING("The extraInfo mcdistanceNeutron is not registered! \n"
59 "This variable is only available for ECL based lists, and you have to run the function getNeutralHadronGeomMatches to use it");
60 return Const::doubleNaN;
61 }
62 }
63
64 int mdstIndexMcKl(const Particle* particle)
65 {
66 if (particle->hasExtraInfo("mdstIndexTruthKL")) {
67 return int(particle->getExtraInfo("mdstIndexTruthKL") + 0.1);
68 } else {
69 B2WARNING("The extraInfo mdstIndexTruthKL is not registered! \n"
70 "This variable is only available for ECL based lists, and you have to run the function getNeutralHadronGeomMatches to use it");
71 return -1;
72 }
73 }
74
75 int mdstIndexMcNeutron(const Particle* particle)
76 {
77 if (particle->hasExtraInfo("mdstIndexTruthNeutron")) {
78 return int(particle->getExtraInfo("mdstIndexTruthNeutron") + 0.1);
79 } else {
80 B2WARNING("The extraInfo mdstIndexTruthNeutron is not registered! \n"
81 "This variable is only available for ECL based lists, and you have to run the function getNeutralHadronGeomMatches to use it");
82 return -1;
83 }
84 }
85
86
87 double beamBackgroundSuppression(const Particle* particle)
88 {
89 if (particle->hasExtraInfo("beamBackgroundSuppression")) {
90 return particle->getExtraInfo("beamBackgroundSuppression");
91 } else {
92 B2WARNING("The extraInfo beamBackgroundSuppression is not registered! \n"
93 "This variable is only available for photons, and you either have to use the standard particle lists (stdPhotons or stdPi0s) or run getBeamBackgroundProbability on a photon list.");
94 return Const::doubleNaN;
95 }
96 }
97
98 double fakePhotonSuppression(const Particle* particle)
99 {
100 if (particle->hasExtraInfo("fakePhotonSuppression")) {
101 return particle->getExtraInfo("fakePhotonSuppression");
102 } else {
103 B2WARNING("The extraInfo fakePhotonSuppression is not registered! \n"
104 "This variable is only available for photons, and you either have to use the standard particle lists (stdPhotons or stdPi0s) or run getFakePhotonProbability on a photon list.");
105 return Const::doubleNaN;
106 }
107 }
108
109 double hadronicSplitOffSuppression(const Particle* particle)
110 {
111 B2WARNING("This variable has been deprecated since light-2302-genetta and is no longer maintained with up to date weights. Please use the variable fakePhotonSuppression instead.");
112 return fakePhotonSuppression(particle);
113 }
114
115 double eclClusterKlId(const Particle* particle)
116 {
117 const ECLCluster* cluster = particle->getECLCluster();
118 if (!cluster) {
119 return Const::doubleNaN;
120 }
121 const KlId* klid = cluster->getRelatedTo<KlId>();
122 if (!klid) {
123 return Const::doubleNaN;
124 }
125 return klid->getKlId();
126 }
127
128
129 double eclPulseShapeDiscriminationMVA(const Particle* particle)
130 {
131 const ECLCluster* cluster = particle->getECLCluster();
132 if (cluster) {
133 if (eclClusterHasPulseShapeDiscrimination(particle)) {
134 return cluster->getPulseShapeDiscriminationMVA();
135 } else {
136 return Const::doubleNaN;
137 }
138 }
139 return Const::doubleNaN;
140 }
141
142 double eclClusterNumberOfHadronDigits(const Particle* particle)
143 {
144
145 const ECLCluster* cluster = particle->getECLCluster();
146 if (cluster) {
147 if (eclClusterHasPulseShapeDiscrimination(particle)) {
148 return cluster->getNumberOfHadronDigits();
149 } else
150 return Const::doubleNaN;
151 }
152 return Const::doubleNaN;
153 }
154
155 double eclClusterDetectionRegion(const Particle* particle)
156 {
157
158 const ECLCluster* cluster = particle->getECLCluster();
159 if (cluster)
160 return cluster->getDetectorRegion();
161
162 return Const::doubleNaN;
163 }
164
165 double eclClusterIsolation(const Particle* particle)
166 {
167
168 const ECLCluster* cluster = particle->getECLCluster();
169 if (cluster) {
170 auto minDist = cluster->getMinTrkDistance();
171 if (minDist > 0)
172 return minDist;
173 }
174 return Const::doubleNaN;
175 }
176
177 double eclClusterIsolationID(const Particle* particle)
178 {
179
180 const ECLCluster* cluster = particle->getECLCluster();
181 if (cluster)
182 return cluster->getMinTrkDistanceID();
183
184 return Const::doubleNaN;
185 }
186
187 Manager::FunctionPtr eclClusterIsolationVar(const std::vector<std::string>& arguments)
188 {
189 if (arguments.size() > 2 or arguments.size() == 0)
190 B2FATAL("Wrong number of arguments (2 required) for meta variable minC2TDistVar");
191 std::string listName = "pi-:all";
192 std::string variableName = arguments[0];
193 if (arguments.size() == 2)
194 listName = arguments[1];
195
196
197 auto func = [listName, variableName](const Particle * particle) -> double {
198 StoreObjPtr<ParticleList> particleList(listName);
199 if (!(particleList.isValid()))
200 {
201 B2FATAL("Invalid Listname " << listName << " given to minC2TDistVar!");
202 }
203 const Variable::Manager::Var* var = Manager::Instance().getVariable(variableName);
204 const ECLCluster* cluster = particle->getECLCluster();
205 if (!cluster)
206 return Const::doubleNaN;
207 auto trackID = cluster->getMinTrkDistanceID();
208 double result = Const::doubleNaN;
209 // Find particle with that track ID:
210 for (unsigned int i = 0; i < particleList->getListSize(); i++)
211 {
212 const Particle* listParticle = particleList->getParticle(i);
213 if (listParticle and listParticle->getTrack() and listParticle->getTrack()->getArrayIndex() == trackID) {
214 result = std::get<double>(var->function(listParticle));
215 break;
216 }
217 }
218 return result;
219 };
220 return func;
221 }
222
223 double eclClusterConnectedRegionID(const Particle* particle)
224 {
225
226 const ECLCluster* cluster = particle->getECLCluster();
227 if (cluster)
228 return cluster->getConnectedRegionId();
229
230 return Const::doubleNaN;
231 }
232
233 double eclClusterDeltaL(const Particle* particle)
234 {
235
236 const ECLCluster* cluster = particle->getECLCluster();
237 if (cluster)
238 return cluster->getDeltaL();
239
240 return Const::doubleNaN;
241 }
242
243 double eclClusterErrorE(const Particle* particle)
244 {
245
246 const ECLCluster* cluster = particle->getECLCluster();
247 if (cluster) {
248 ClusterUtils clutls;
249 const auto EPhiThetaCov = clutls.GetCovarianceMatrix3x3FromCluster(cluster);
250 return sqrt(fabs(EPhiThetaCov[0][0]));
251 }
252 return Const::doubleNaN;
253 }
254
255 double eclClusterUncorrectedE(const Particle* particle)
256 {
257
258 const ECLCluster* cluster = particle->getECLCluster();
259 if (cluster) {
260 return cluster->getEnergyRaw();
261 }
262 return Const::doubleNaN;
263 }
264
265 double eclClusterE(const Particle* particle)
266 {
267 const auto& frame = ReferenceFrame::GetCurrent();
268
269 const ECLCluster* cluster = particle->getECLCluster();
270 if (cluster) {
271 ClusterUtils clutls;
272 ROOT::Math::PxPyPzEVector p4Cluster = clutls.GetCluster4MomentumFromCluster(cluster, particle->getECLClusterEHypothesisBit());
273
274 return frame.getMomentum(p4Cluster).E();
275 }
276 return Const::doubleNaN;
277 }
278
279 double eclClusterHighestE(const Particle* particle)
280 {
281
282 const ECLCluster* cluster = particle->getECLCluster();
283 if (cluster) {
284 return cluster->getEnergyHighestCrystal();
285 }
286 return Const::doubleNaN;
287 }
288
289 double eclClusterCellId(const Particle* particle)
290 {
291
292 const ECLCluster* cluster = particle->getECLCluster();
293 if (cluster) {
294 return cluster->getMaxECellId();
295 }
296 return Const::doubleNaN;
297 }
298
299 // An array with each number representing the last number of the cellID per thetaID. There are 69 thetaIDs in total.
300 const std::array<int, 69> lastCellIDperThetaID{48, 96, 160, 224, 288, 384, 480, 576, 672, 768, 864,
301 1008, 1152, 1296, 1440, 1584, 1728, 1872, 2016, 2160, 2304, 2448,
302 2592, 2736, 2880, 3024, 3168, 3312, 3456, 3600, 3744, 3888, 4032,
303 4176, 4320, 4464, 4608, 4752, 4896, 5040, 5184, 5328, 5472, 5616,
304 5760, 5904, 6048, 6192, 6336, 6480, 6624, 6768, 6912, 7056, 7200,
305 7344, 7488, 7632, 7776, 7920, 8064, 8160, 8256, 8352, 8448, 8544,
306 8608, 8672, 8736};
307
308 double eclClusterThetaId(const Particle* particle)
309 {
310
311 const ECLCluster* cluster = particle->getECLCluster();
312 if (cluster) {
313 int cellID = cluster->getMaxECellId();
314 return std::distance(lastCellIDperThetaID.begin(), std::lower_bound(lastCellIDperThetaID.begin(), lastCellIDperThetaID.end(),
315 cellID));
316 }
317 return Const::doubleNaN;
318 }
319
320 double eclClusterPhiId(const Particle* particle)
321 {
322
323 const ECLCluster* cluster = particle->getECLCluster();
324 if (cluster) {
325 int cellID = cluster->getMaxECellId();
326 if (cellID <= 48) {
327 return cellID - 1;
328 } else {
329 int closestinlist = lastCellIDperThetaID[std::distance(lastCellIDperThetaID.begin(), std::lower_bound(lastCellIDperThetaID.begin(),
330 lastCellIDperThetaID.end(), cellID)) - 1];
331 return cellID - closestinlist - 1;
332 }
333 }
334 return Const::doubleNaN;
335 }
336
337 double eclClusterTiming(const Particle* particle)
338 {
339
340 const ECLCluster* cluster = particle->getECLCluster();
341 if (cluster) {
342 return cluster->getTime();
343 }
344 return Const::doubleNaN;
345 }
346
347 double eclClusterHasFailedTiming(const Particle* particle)
348 {
349 const ECLCluster* cluster = particle->getECLCluster();
350 if (cluster) {
351 return cluster->hasFailedFitTime();
352 }
353 return Const::doubleNaN;
354 }
355
356 double eclClusterErrorTiming(const Particle* particle)
357 {
358
359 const ECLCluster* cluster = particle->getECLCluster();
360 if (cluster) {
361 return cluster->getDeltaTime99();
362 }
363 return Const::doubleNaN;
364 }
365
366 double eclClusterHasFailedErrorTiming(const Particle* particle)
367 {
368 const ECLCluster* cluster = particle->getECLCluster();
369 if (cluster) {
370 return cluster->hasFailedTimeResolution();
371 }
372 return Const::doubleNaN;
373 }
374
375 double eclClusterTheta(const Particle* particle)
376 {
377 const auto& frame = ReferenceFrame::GetCurrent();
378
379 const ECLCluster* cluster = particle->getECLCluster();
380 if (cluster) {
381 ClusterUtils clutls;
382 ROOT::Math::PxPyPzEVector p4Cluster = clutls.Get4MomentumFromCluster(cluster, particle->getECLClusterEHypothesisBit());
383
384 return frame.getMomentum(p4Cluster).Theta();
385 }
386 return Const::doubleNaN;
387 }
388
389 double eclClusterErrorTheta(const Particle* particle)
390 {
391
392 const ECLCluster* cluster = particle->getECLCluster();
393 if (cluster) {
394 ClusterUtils clutls;
395 const auto EPhiThetaCov = clutls.GetCovarianceMatrix3x3FromCluster(cluster);
396 return sqrt(fabs(EPhiThetaCov[2][2]));
397 }
398 return Const::doubleNaN;
399 }
400
401 double eclClusterErrorPhi(const Particle* particle)
402 {
403
404 const ECLCluster* cluster = particle->getECLCluster();
405 if (cluster) {
406 ClusterUtils clutls;
407 const auto EPhiThetaCov = clutls.GetCovarianceMatrix3x3FromCluster(cluster);
408 return sqrt(fabs(EPhiThetaCov[1][1]));
409 }
410 return Const::doubleNaN;
411 }
412
413 double eclClusterPhi(const Particle* particle)
414 {
415 const auto& frame = ReferenceFrame::GetCurrent();
416
417 const ECLCluster* cluster = particle->getECLCluster();
418 if (cluster) {
419 ClusterUtils clutls;
420 ROOT::Math::PxPyPzEVector p4Cluster = clutls.Get4MomentumFromCluster(cluster, particle->getECLClusterEHypothesisBit());
421
422 return frame.getMomentum(p4Cluster).Phi();
423 }
424 return Const::doubleNaN;
425 }
426
427 double eclClusterR(const Particle* particle)
428 {
429
430 const ECLCluster* cluster = particle->getECLCluster();
431 if (cluster) {
432 return cluster->getR();
433 }
434 return Const::doubleNaN;
435 }
436
437 double eclClusterE1E9(const Particle* particle)
438 {
439
440 const ECLCluster* cluster = particle->getECLCluster();
441 if (cluster) {
442 return cluster->getE1oE9();
443 }
444 return Const::doubleNaN;
445 }
446
447 double eclClusterE9E21(const Particle* particle)
448 {
449
450 const ECLCluster* cluster = particle->getECLCluster();
451 if (cluster) {
452 return cluster->getE9oE21();
453 }
454 return Const::doubleNaN;
455 }
456
457 double eclClusterAbsZernikeMoment40(const Particle* particle)
458 {
459
460 const ECLCluster* cluster = particle->getECLCluster();
461 if (cluster) {
462 return cluster->getAbsZernike40();
463 }
464 return Const::doubleNaN;
465 }
466
467 double eclClusterAbsZernikeMoment51(const Particle* particle)
468 {
469
470 const ECLCluster* cluster = particle->getECLCluster();
471 if (cluster) {
472 return cluster->getAbsZernike51();
473 }
474 return Const::doubleNaN;
475 }
476
477 double eclClusterZernikeMVA(const Particle* particle)
478 {
479
480 const ECLCluster* cluster = particle->getECLCluster();
481 if (cluster) {
482 return cluster->getZernikeMVA();
483 }
484 return Const::doubleNaN;
485 }
486
487 double eclClusterSecondMoment(const Particle* particle)
488 {
489
490 const ECLCluster* cluster = particle->getECLCluster();
491 if (cluster) {
492 return cluster->getSecondMoment();
493 }
494 return Const::doubleNaN;
495 }
496
497 double eclClusterLAT(const Particle* particle)
498 {
499
500 const ECLCluster* cluster = particle->getECLCluster();
501 if (cluster) {
502 return cluster->getLAT();
503 }
504 return Const::doubleNaN;
505 }
506
507 double eclClusterNHits(const Particle* particle)
508 {
509
510 const ECLCluster* cluster = particle->getECLCluster();
511 if (cluster) {
512 return cluster->getNumberOfCrystals();
513 }
514 return Const::doubleNaN;
515 }
516
517 double eclClusterTrackMatched(const Particle* particle)
518 {
519
520 const ECLCluster* cluster = particle->getECLCluster();
521 if (cluster) {
522 const Track* track = cluster->getRelatedFrom<Track>();
523
524 if (track)
525 return 1.0;
526 else
527 return 0.0;
528 }
529 return Const::doubleNaN;
530 }
531
532 double nECLClusterTrackMatches(const Particle* particle)
533 {
534 // if no ECL cluster then nan
535 const ECLCluster* cluster = particle->getECLCluster();
536 if (!cluster)
537 return Const::doubleNaN;
538
539 // one or more tracks may be matched to charged particles
540 size_t out = cluster->getRelationsFrom<Track>().size();
541 return double(out);
542 }
543
544 double eclClusterId(const Particle* particle)
545 {
546 const ECLCluster* cluster = particle->getECLCluster();
547 if (cluster) {
548 return cluster->getClusterId();
549 }
550 return Const::doubleNaN;
551 }
552
553 double eclClusterHasNPhotonsHypothesis(const Particle* particle)
554 {
555 const ECLCluster* cluster = particle->getECLCluster();
556 if (cluster) {
557 return cluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons);
558 }
559 return Const::doubleNaN;
560 }
561
562 double eclClusterHasNeutralHadronHypothesis(const Particle* particle)
563 {
564 const ECLCluster* cluster = particle->getECLCluster();
565 if (cluster) {
566 return cluster->hasHypothesis(ECLCluster::EHypothesisBit::c_neutralHadron);
567 }
568 return Const::doubleNaN;
569 }
570
571 double eclClusterHasPulseShapeDiscrimination(const Particle* particle)
572 {
573 const ECLCluster* cluster = particle->getECLCluster();
574 if (cluster) {
575 return cluster->hasPulseShapeDiscrimination();
576 }
577 return Const::doubleNaN;
578 }
579
580 double eclExtTheta(const Particle* particle)
581 {
582
583 const Track* track = particle->getTrack();
584 if (track) {
585
586 auto* eclinfo = track->getRelatedTo<ECLEnergyCloseToTrack>();
587
588 if (eclinfo) {
589 return eclinfo->getExtTheta();
590 } else {
591 B2WARNING("Relation to ECLEnergyCloseToTrack not found, did you forget to run ECLTrackCalDigitMatchModule?");
592 return Const::doubleNaN;
593 }
594 }
595
596 return Const::doubleNaN;
597 }
598
599 double eclExtPhi(const Particle* particle)
600 {
601
602 const Track* track = particle->getTrack();
603 if (track) {
604
605 auto* eclinfo = track->getRelatedTo<ECLEnergyCloseToTrack>();
606
607 if (eclinfo) {
608 return eclinfo->getExtPhi();
609 } else {
610 B2WARNING("Relation to ECLEnergyCloseToTrack not found, did you forget to run ECLTrackCalDigitMatchModule?");
611 return Const::doubleNaN;
612 }
613 }
614
615 return Const::doubleNaN;
616 }
617
618 double eclExtPhiId(const Particle* particle)
619 {
620 const Track* track = particle->getTrack();
621 if (track) {
622
623 auto* eclinfo = track->getRelatedTo<ECLEnergyCloseToTrack>();
624
625 if (eclinfo) {
626 return eclinfo->getExtPhiId();
627 } else {
628 B2WARNING("Relation to ECLEnergyCloseToTrack not found, did you forget to run ECLTrackCalDigitMatchModule?");
629 return Const::doubleNaN;
630 }
631 }
632
633 return Const::doubleNaN;
634 }
635
636 double weightedAverageECLTime(const Particle* particle)
637 {
638 int nDaughters = particle->getNDaughters();
639 if (nDaughters < 1) {
640 B2WARNING("The provided particle has no daughters!");
641 return Const::doubleNaN;
642 }
643
644 double numer = 0, denom = 0;
645 int numberOfClusterDaughters = 0;
646
647 auto weightedECLTimeAverage = [&numer, &denom, &numberOfClusterDaughters](const Particle * p) {
648 const ECLCluster* cluster = p->getECLCluster();
649 if (cluster and not cluster->hasFailedFitTime()) {
650 numberOfClusterDaughters ++;
651
652 double time = cluster->getTime();
653 B2DEBUG(10, "time[" << numberOfClusterDaughters << "] = " << time);
654 double deltatime = cluster->getDeltaTime99();
655 B2DEBUG(10, "deltatime[" << numberOfClusterDaughters << "] = " << deltatime);
656 numer += time / pow(deltatime, 2);
657 B2DEBUG(11, "numer[" << numberOfClusterDaughters << "] = " << numer);
658 denom += 1 / pow(deltatime, 2);
659 B2DEBUG(11, "denom[" << numberOfClusterDaughters << "] = " << denom);
660 }
661 return false;
662 };
663
664 particle->forEachDaughter(weightedECLTimeAverage, true, true);
665
666 if (numberOfClusterDaughters < 1) {
667 B2WARNING("There are no clusters or cluster matches amongst the daughters of the provided particle!");
668 return Const::doubleNaN;
669 }
670
671 if (denom == 0) {
672 B2WARNING("The denominator of the weighted mean is zero!");
673 return Const::doubleNaN;
674 } else {
675 B2DEBUG(10, "numer/denom = " << numer / denom);
676 return numer / denom;
677 }
678 }
679
680 double maxWeightedDistanceFromAverageECLTime(const Particle* particle)
681 {
682 int nDaughters = particle->getNDaughters();
683 if (nDaughters < 1) {
684 B2WARNING("The provided particle has no daughters!");
685 return Const::doubleNaN;
686 }
687
688 double maxTimeDiff = -DBL_MAX;
689 int numberOfClusterDaughters = 0;
690
691 double averageECLTime = weightedAverageECLTime(particle);
692
693 auto maxTimeDifference = [&maxTimeDiff, &numberOfClusterDaughters, &averageECLTime](const Particle * p) {
694
695 const ECLCluster* cluster = p->getECLCluster();
696 if (cluster) {
697 numberOfClusterDaughters ++;
698
699 double time = cluster->getTime();
700 B2DEBUG(10, "time[" << numberOfClusterDaughters << "] = " << time);
701 double deltatime = cluster->getDeltaTime99();
702 B2DEBUG(10, "deltatime[" << numberOfClusterDaughters << "] = " << deltatime);
703 double maxTimeDiff_temp = fabs((time - averageECLTime) / deltatime);
704 B2DEBUG(11, "maxTimeDiff_temp[" << numberOfClusterDaughters << "] = " << maxTimeDiff_temp);
705 if (maxTimeDiff_temp > maxTimeDiff)
706 maxTimeDiff = maxTimeDiff_temp;
707 B2DEBUG(11, "maxTimeDiff[" << numberOfClusterDaughters << "] = " << maxTimeDiff);
708 }
709 return false;
710 };
711
712 particle->forEachDaughter(maxTimeDifference, true, true);
713
714 if (numberOfClusterDaughters < 1) {
715 B2WARNING("There are no clusters or cluster matches amongst the daughters of the provided particle!");
716 return Const::doubleNaN;
717 }
718
719 if (maxTimeDiff < 0) {
720 B2WARNING("The max time difference is negative!");
721 return Const::doubleNaN;
722 } else {
723 B2DEBUG(10, "maxTimeDiff = " << maxTimeDiff);
724 return maxTimeDiff;
725 }
726 }
727
728 double eclClusterMdstIndex(const Particle* particle)
729 {
730 const ECLCluster* cluster = particle->getECLCluster();
731 if (cluster) {
732 return cluster->getArrayIndex();
733 } else return Const::doubleNaN;
734
735 return Const::doubleNaN;
736 }
737
738
739 /*************************************************************
740 * Event-based ECL clustering information
741 */
742 double nECLOutOfTimeCrystalsFWDEndcap(const Particle*)
743 {
744 StoreObjPtr<EventLevelClusteringInfo> elci;
745 if (!elci) return Const::doubleNaN;
746 return (double)elci->getNECLCalDigitsOutOfTimeFWD();
747 }
748
749 double nECLOutOfTimeCrystalsBarrel(const Particle*)
750 {
751 StoreObjPtr<EventLevelClusteringInfo> elci;
752 if (!elci) return Const::doubleNaN;
753 return (double)elci->getNECLCalDigitsOutOfTimeBarrel();
754 }
755
756 double nECLOutOfTimeCrystalsBWDEndcap(const Particle*)
757 {
758 StoreObjPtr<EventLevelClusteringInfo> elci;
759 if (!elci) return Const::doubleNaN;
760 return (double)elci->getNECLCalDigitsOutOfTimeBWD();
761 }
762
763 double nECLOutOfTimeCrystals(const Particle*)
764 {
765 StoreObjPtr<EventLevelClusteringInfo> elci;
766 if (!elci) return Const::doubleNaN;
767 return (double)elci->getNECLCalDigitsOutOfTime();
768 }
769
770 double nRejectedECLShowersFWDEndcap(const Particle*)
771 {
772 StoreObjPtr<EventLevelClusteringInfo> elci;
773 if (!elci) return Const::doubleNaN;
774 return (double)elci->getNECLShowersRejectedFWD();
775 }
776
777 double nRejectedECLShowersBarrel(const Particle*)
778 {
779 StoreObjPtr<EventLevelClusteringInfo> elci;
780 if (!elci) return Const::doubleNaN;
781 return (double)elci->getNECLShowersRejectedBarrel();
782 }
783
784 double nRejectedECLShowersBWDEndcap(const Particle*)
785 {
786 StoreObjPtr<EventLevelClusteringInfo> elci;
787 if (!elci) return Const::doubleNaN;
788 return (double)elci->getNECLShowersRejectedBWD();
789 }
790
791 double nRejectedECLShowers(const Particle*)
792 {
793 StoreObjPtr<EventLevelClusteringInfo> elci;
794 if (!elci) return Const::doubleNaN;
795 return (double) elci->getNECLShowersRejected();
796 }
797
798 double nKLMMultistripHitsFWDEndcap(const Particle*)
799 {
800 StoreObjPtr<EventLevelClusteringInfo> elci;
801 if (!elci) return Const::doubleNaN;
802 return (double) elci->getNKLMDigitsMultiStripFWD();
803 }
804
805 double nKLMMultistripHitsBarrel(const Particle*)
806 {
807 StoreObjPtr<EventLevelClusteringInfo> elci;
808 if (!elci) return Const::doubleNaN;
809 return (double) elci->getNKLMDigitsMultiStripBarrel();
810 }
811
812 double nKLMMultistripHitsBWDEndcap(const Particle*)
813 {
814 StoreObjPtr<EventLevelClusteringInfo> elci;
815 if (!elci) return Const::doubleNaN;
816 return (double) elci->getNKLMDigitsMultiStripBWD();
817 }
818
819 double nKLMMultistripHits(const Particle*)
820 {
821 StoreObjPtr<EventLevelClusteringInfo> elci;
822 if (!elci) return Const::doubleNaN;
823 return (double) elci->getNKLMDigitsMultiStrip();
824 }
825
826 double nECLShowersFWDEndcap(const Particle*)
827 {
828 StoreObjPtr<EventLevelClusteringInfo> elci;
829 if (!elci) return Const::doubleNaN;
830 return (double) elci->getNECLShowersFWD();
831 }
832
833 double nECLShowersBarrel(const Particle*)
834 {
835 StoreObjPtr<EventLevelClusteringInfo> elci;
836 if (!elci) return Const::doubleNaN;
837 return (double) elci->getNECLShowersBarrel();
838 }
839
840 double nECLShowersBWDEndcap(const Particle*)
841 {
842 StoreObjPtr<EventLevelClusteringInfo> elci;
843 if (!elci) return Const::doubleNaN;
844 return (double) elci->getNECLShowersBWD();
845 }
846
847 double nECLShowers(const Particle*)
848 {
849 StoreObjPtr<EventLevelClusteringInfo> elci;
850 if (!elci) return Const::doubleNaN;
851 return (double) elci->getNECLShowers();
852 }
853
854 double nECLLocalMaximumsFWDEndcap(const Particle*)
855 {
856 StoreObjPtr<EventLevelClusteringInfo> elci;
857 if (!elci) return Const::doubleNaN;
858 return (double) elci->getNECLLocalMaximumsFWD();
859 }
860
861 double nECLLocalMaximumsBarrel(const Particle*)
862 {
863 StoreObjPtr<EventLevelClusteringInfo> elci;
864 if (!elci) return Const::doubleNaN;
865 return (double) elci->getNECLLocalMaximumsBarrel();
866 }
867
868 double nECLLocalMaximumsBWDEndcap(const Particle*)
869 {
870 StoreObjPtr<EventLevelClusteringInfo> elci;
871 if (!elci) return Const::doubleNaN;
872 return (double) elci->getNECLLocalMaximumsBWD();
873 }
874
875 double nECLLocalMaximums(const Particle*)
876 {
877 StoreObjPtr<EventLevelClusteringInfo> elci;
878 if (!elci) return Const::doubleNaN;
879 return (double) elci->getNECLLocalMaximums();
880 }
881
882 double nECLTriggerCellsFWDEndcap(const Particle*)
883 {
884 StoreObjPtr<EventLevelClusteringInfo> elci;
885 if (!elci) return Const::doubleNaN;
886 return (double) elci->getNECLTriggerCellsFWD();
887 }
888
889 double nECLTriggerCellsBarrel(const Particle*)
890 {
891 StoreObjPtr<EventLevelClusteringInfo> elci;
892 if (!elci) return Const::doubleNaN;
893 return (double) elci->getNECLTriggerCellsBarrel();
894 }
895
896 double nECLTriggerCellsBWDEndcap(const Particle*)
897 {
898 StoreObjPtr<EventLevelClusteringInfo> elci;
899 if (!elci) return Const::doubleNaN;
900 return (double) elci->getNECLTriggerCellsBWD();
901 }
902
903 double nECLTriggerCells(const Particle*)
904 {
905 StoreObjPtr<EventLevelClusteringInfo> elci;
906 if (!elci) return Const::doubleNaN;
907 return (double) elci->getNECLTriggerCells();
908 }
909
910 double eclClusterEoP(const Particle* part)
911 {
912 double E = eclClusterE(part);
913 if (part->hasExtraInfo("bremsCorrectedPhotonEnergy")) {
914 E += part->getExtraInfo("bremsCorrectedPhotonEnergy");
915 }
916 const double p = part->getMomentumMagnitude();
917 if (0 == p) { return Const::doubleNaN;}
918 return E / p;
919 }
920
921 double eclClusterOnlyInvariantMass(const Particle* part)
922 {
923 int nDaughters = part->getNDaughters();
924 ROOT::Math::PxPyPzEVector sum;
925
926 if (nDaughters < 1) {
927 return part->getMass();
928 } else {
929 int nClusterDaughters = 0;
930 std::stack<const Particle*> stacked;
931 stacked.push(part);
932 while (!stacked.empty()) {
933 const Particle* current = stacked.top();
934 stacked.pop();
935
936 const ECLCluster* cluster = current->getECLCluster();
937 if (cluster) {
938 const ECLCluster::EHypothesisBit clusterBit = current->getECLClusterEHypothesisBit();
939 nClusterDaughters ++;
940 ClusterUtils clutls;
941 ROOT::Math::PxPyPzEVector p4Cluster = clutls.Get4MomentumFromCluster(cluster, clusterBit);
942 sum += p4Cluster;
943 } else {
944 const std::vector<Particle*> daughters = current->getDaughters();
945 nDaughters = current->getNDaughters();
946 for (int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
947 stacked.push(daughters[iDaughter]);
948 }
949 }
950 }
951
952 if (nClusterDaughters < 1) {
953 B2WARNING("There are no clusters amongst the daughters of the provided particle!");
954 return Const::doubleNaN;
955 }
956 B2DEBUG(10, "Number of daughters with cluster associated = " << nClusterDaughters);
957 return sum.M();
958 }
959 }
960
961 Manager::FunctionPtr photonHasOverlap(const std::vector<std::string>& arguments)
962 {
963 std::string cutString = "";
964 if (arguments.size() > 0) {
965 cutString = arguments[0];
966 }
967 std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
968
969 std::string photonlistname = "gamma:all";
970 if (arguments.size() > 1) {
971 photonlistname = arguments[1];
972 }
973
974 std::string tracklistname = "e-:all";
975 if (arguments.size() > 2) {
976 tracklistname = arguments[2];
977 }
978
979 auto func = [cut, photonlistname, tracklistname](const Particle * particle) -> double {
980
981 if (particle->getPDGCode() != Const::photon.getPDGCode())
982 {
983 B2WARNING("The variable photonHasOverlap is supposed to be calculated for photons. Returning NaN.");
984 return Const::doubleNaN;
985 }
986
987 StoreObjPtr<ParticleList> photonlist(photonlistname);
988 if (!(photonlist.isValid()))
989 {
990 B2WARNING("The provided particle list " << photonlistname << " does not exist."
991 " Therefore, the variable photonHasOverlap can not be calculated. Returning NaN.");
992 return Const::doubleNaN;
993 }
994 if (photonlist->getPDGCode() != Const::photon.getPDGCode())
995 {
996 B2WARNING("The list " << photonlistname << " does not contain photons."
997 " Therefore, the variable photonHasOverlap can not be calculated reliably. Returning NaN.");
998 return Const::doubleNaN;
999 }
1000
1001 StoreObjPtr<ParticleList> tracklist(tracklistname);
1002 if (!(tracklist.isValid()))
1003 {
1004 B2WARNING("The provided particle list " << tracklistname << " does not exist."
1005 " Therefore, the variable photonHasOverlap can not be calculated. Returning NaN.");
1006 return Const::doubleNaN;
1007 }
1008 if (!Const::chargedStableSet.contains(Const::ParticleType(abs(tracklist->getPDGCode()))))
1009 {
1010 B2WARNING("The list " << tracklistname << " does not contain charged final state particles."
1011 " Therefore, the variable photonHasOverlap can not be calculated reliably. Returning NaN.");
1012 return Const::doubleNaN;
1013 }
1014
1015 double connectedRegionID = eclClusterConnectedRegionID(particle);
1016 unsigned mdstArrayIndex = particle->getMdstArrayIndex();
1017
1018 for (unsigned int i = 0; i < photonlist->getListSize(); i++)
1019 {
1020 const Particle* part = photonlist->getParticle(i);
1021
1022 // skip the particle itself
1023 if (part->getMdstArrayIndex() == mdstArrayIndex) {
1024 continue;
1025 }
1026
1027 // skip photons that do not fulfill the provided criteria
1028 if (!cut->check(part)) {
1029 continue;
1030 }
1031
1032 if (connectedRegionID == eclClusterConnectedRegionID(part)) {
1033 return 1;
1034 }
1035 }
1036
1037 for (unsigned int i = 0; i < tracklist->getListSize(); i++)
1038 {
1039 const Particle* part = tracklist->getParticle(i);
1040
1041 // skip tracks that do not fulfill the provided criteria
1042 if (!cut->check(part)) {
1043 continue;
1044 }
1045
1046 if (connectedRegionID == eclClusterConnectedRegionID(part)) {
1047 return 1;
1048 }
1049 }
1050
1051 return 0;
1052 };
1053 return func;
1054 }
1055
1056 VARIABLE_GROUP("ECL Cluster related");
1057 REGISTER_VARIABLE("clusterEoP", eclClusterEoP, R"DOC(
1058Returns ratio of uncorrelated energy E over momentum p, a convenience
1059alias for (clusterE / p).
1060)DOC");
1061 REGISTER_VARIABLE("clusterReg", eclClusterDetectionRegion, R"DOC(
1062Returns an integer code for the ECL region of a cluster.
1063
1064- 1: forward, 2: barrel, 3: backward,
1065- 11: between FWD and barrel, 13: between BWD and barrel,
1066- 0: otherwise
1067)DOC");
1068 REGISTER_VARIABLE("clusterDeltaLTemp", eclClusterDeltaL, R"DOC(
1069| Returns DeltaL for the shower shape.
1070| A cluster comprises the energy depositions of several crystals. All these crystals have slightly
1071 different orientations in space. A shower direction can be constructed by calculating the weighted
1072 average of these orientations using the corresponding energy depositions as weights. The intersection
1073 (more precisely the point of closest approach) of the vector with this direction originating from the
1074 cluster center and an extrapolated track can be used as reference for the calculation of the shower
1075 depth. It is defined as the distance between this intersection and the cluster center.
1076
1077.. warning::
1078 This distance is calculated on the reconstructed level and is temporarily
1079 included to the ECL cluster MDST data format for studying purposes. If it is found
1080 that it is not crucial for physics analysis then this variable will be removed
1081 in future releases.
1082 Therefore, keep in mind that this variable might be removed in the future!
1083
1084.. note::
1085 | Please read `this <importantNoteECL>` first.
1086 | Lower limit: :math:`-250.0`
1087 | Upper limit: :math:`250.0`
1088 | Precision: :math:`10` bit
1089..
1090
1091)DOC","cm");
1092
1093 REGISTER_VARIABLE("minC2TDist", eclClusterIsolation, R"DOC(
1094Returns the distance between the ECL cluster and its nearest track.
1095
1096For all tracks in the event, the distance between each of their extrapolated hits in the ECL and the ECL shower
1097position is calculated, and the overall smallest distance is returned. The track array index of the track that is
1098closest to the ECL cluster can be retrieved using `minC2TDistID`.
1099
1100If the calculated distance is greater than :math:`250.0`, the returned distance will be capped at :math:`250.0`.
1101If there are no extrapolated hits found in the ECL for the event, NaN will be returned.
1102
1103.. note::
1104 This distance is calculated on the reconstructed level.
1105
1106.. note::
1107 | Please read `this <importantNoteECL>` first.
1108 | Lower limit: :math:`0.0`
1109 | Upper limit: :math:`250.0`
1110 | Precision: :math:`10` bit
1111..
1112
1113)DOC","cm");
1114 REGISTER_VARIABLE("minC2TDistID", eclClusterIsolationID, R"DOC(
1115Returns the track array index of the nearest track to the ECL cluster. The nearest track is calculated
1116using the `minC2TDist` variable.
1117)DOC");
1118 REGISTER_METAVARIABLE("minC2TDistVar(variable,particleList=pi-:all)", eclClusterIsolationVar, R"DOC(
1119Returns the variable value of the nearest track to the given ECL cluster as calculated by `minC2TDist`. The
1120first argument is the variable name, e.g. `nCDCHits`, while the second (optional) argument is the particle list name which
1121will be used to pick up the nearest track in the calculation of `minC2TDist`. The default particle list used
1122is ``pi-:all``.
1123)DOC", Manager::VariableDataType::c_double);
1124 REGISTER_VARIABLE("clusterE", eclClusterE, R"DOC(
1125Returns ECL cluster's energy corrected for leakage and background.
1126
1127The raw photon energy is given by the weighted sum of all ECL crystal energies within the ECL cluster.
1128The weights per crystals are :math:`\leq 1` after cluster energy splitting in the case of overlapping
1129clusters. The number of crystals that are included in the sum depends on a initial energy estimation
1130and local beam background levels at the highest energy crystal position. It is optimized to minimize
1131the core width (resolution) of true photons. Photon energy distributions always show a low energy tail
1132due to unavoidable longitudinal and transverse leakage that can be further modified by the clustering
1133algorithm and beam backgrounds.The peak position of the photon energy distributions are corrected to
1134match the true photon energy in MC:
1135
1136- Leakage correction: Using large MC samples of mono-energetic single photons, a correction factor
1137 :math:`f` as function of reconstructed detector position, reconstructed photon energy and beam backgrounds
1138 is determined via :math:`f = \frac{\text{peak_reconstructed}}{\text{energy_true}}`.
1139
1140- Cluster energy calibration (data only): To reach the target precision of :math:`< 1.8\%` energy
1141 resolution for high energetic photons, the remaining difference between MC and data must be calibrated
1142 using kinematically fit muon pairs. This calibration is only applied to data and not to MC and will
1143 take time to develop.
1144
1145- Energy Bias Correction module, sub-percent correction, is NOT applied on clusterE, but on photon energy
1146 and momentum. Only applied to data.
1147
1148It is important to note that after perfect leakage correction and cluster energy calibration,
1149the :math:`\pi^{0}` mass peak will be shifted slightly to smaller values than the PDG average
1150due to the low energy tails of photons. The :math:`\pi^{0}` mass peak must not be corrected
1151to the PDG value by adjusting the reconstructed photon energies. Selection criteria based on
1152the mass for :math:`\pi^{0}` candidates must be based on the biased value. Most analysis
1153will used mass constrained :math:`\pi^{0}` s anyhow.
1154
1155.. warning::
1156 We only store clusters with :math:`E > 20\,` MeV.
1157
1158.. note::
1159 | Please read `this <importantNoteECL>` first.
1160 | Lower limit: :math:`-5` (:math:`e^{-5} = 0.00674\,` GeV)
1161 | Upper limit: :math:`3.0` (:math:`e^3 = 20.08553\,` GeV)
1162 | Precision: :math:`18` bit
1163 | This value can be changed to a different reference frame with :b2:var:`useCMSFrame`.
1164..
1165
1166)DOC","GeV");
1167 REGISTER_VARIABLE("clusterErrorE", eclClusterErrorE, R"DOC(
1168Returns ECL cluster's uncertainty on energy
1169(from background level and energy dependent tabulation).
1170
1171)DOC","GeV");
1172 REGISTER_VARIABLE("clusterErrorPhi", eclClusterErrorPhi, R"DOC(
1173Returns ECL cluster's uncertainty on :math:`\phi`
1174(from background level and energy dependent tabulation).
1175
1176)DOC","rad");
1177 REGISTER_VARIABLE("clusterErrorTheta", eclClusterErrorTheta, R"DOC(
1178Returns ECL cluster's uncertainty on :math:`\theta`
1179(from background level and energy dependent tabulation).
1180
1181)DOC","rad");
1182
1183 REGISTER_VARIABLE("clusterR", eclClusterR, R"DOC(
1184Returns ECL cluster's centroid distance from :math:`(0,0,0)`.
1185
1186)DOC","cm");
1187 REGISTER_VARIABLE("clusterPhi", eclClusterPhi, R"DOC(
1188Returns ECL cluster's azimuthal angle :math:`\phi`
1189(this is not generally equal to a photon azimuthal angle).
1190
1191| The direction of a cluster is given by the connecting line of :math:`\,(0,0,0)\,` and
1192 cluster centroid position in the ECL.
1193| The cluster centroid position is calculated using up to 21 crystals (5x5 excluding corners)
1194 after cluster energy splitting in the case of overlapping clusters.
1195| The centroid position is the logarithmically weighted average of all crystals evaluated at
1196 the crystal centers. Cluster centroids are generally biased towards the centers of the
1197 highest energetic crystal. This effect is larger for low energetic photons.
1198| Beam backgrounds slightly decrease the position resolution, mainly for low energetic photons.
1199
1200.. note::
1201 Radius of a cluster is almost constant in the barrel and should not be used directly in any selection.
1202
1203Unlike for charged tracks, the uncertainty (covariance) of the photon directions is not determined
1204based on individual cluster properties but taken from on MC-based parametrizations of the resolution
1205as function of true photon energy, true photon direction and beam background level.
1206
1207.. warning::
1208 Users must use the actual particle direction (done automatically in the modularAnalysis using the average
1209 IP position (can be changed if needed)) and not the ECL Cluster direction (position in the ECL measured
1210 from :math:`(0,0,0)`) for particle kinematics.
1211
1212.. note::
1213 | Please read `this <importantNoteECL>` first.
1214 | Lower limit: :math:`-\pi`
1215 | Upper limit: :math:`\pi`
1216 | Precision: :math:`16` bit
1217..
1218
1219)DOC","rad");
1220 REGISTER_VARIABLE("clusterConnectedRegionID", eclClusterConnectedRegionID, R"DOC(
1221Returns ECL cluster's connected region ID.
1222)DOC");
1223 REGISTER_VARIABLE("clusterTheta", eclClusterTheta, R"DOC(
1224Returns ECL cluster's polar angle :math:`\theta`
1225(this is not generally equal to a photon polar angle).
1226
1227| The direction of a cluster is given by the connecting line of :math:`\,(0,0,0)\,` and
1228 cluster centroid position in the ECL.
1229| The cluster centroid position is calculated using up to 21 crystals (5x5 excluding corners)
1230 after cluster energy splitting in the case of overlapping clusters.
1231| The centroid position is the logarithmically weighted average of all crystals evaluated at
1232 the crystal centers. Cluster centroids are generally biased towards the centers of the
1233 highest energetic crystal. This effect is larger for low energetic photons.
1234| Beam backgrounds slightly decrease the position resolution, mainly for low energetic photons.
1235
1236.. note::
1237 Radius of a cluster is almost constant in the barrel and should not be used directly in any selection.
1238
1239Unlike for charged tracks, the uncertainty (covariance) of the photon directions is not determined
1240based on individual cluster properties but taken from on MC-based parametrizations of the resolution
1241as function of true photon energy, true photon direction and beam background level.
1242
1243.. warning::
1244 Users must use the actual particle direction (done automatically in the modularAnalysis using the average
1245 IP position (can be changed if needed)) and not the ECL Cluster direction (position in the ECL measured
1246 from :math:`(0,0,0)`) for particle kinematics.
1247
1248.. note::
1249 | Please read `this <importantNoteECL>` first.
1250 | Lower limit: :math:`0.0`
1251 | Upper limit: :math:`\pi`
1252 | Precision: :math:`16` bit
1253..
1254
1255)DOC","rad");
1256 REGISTER_VARIABLE("clusterTiming", eclClusterTiming, R"DOC(
1257**In Belle II:**
1258Returns the time of the ECL cluster. It is calculated as the Photon timing minus the Event t0.
1259Photon timing is given by the fitted time of the recorded waveform of the highest energy crystal in the
1260cluster. After all calibrations and corrections (including Time-Of-Flight), photons from the interaction
1261point (IP) should have a Photon timing that corresponds to the Event t0, :math:`t_{0}`. The Event t0 is the
1262time of the event and may be measured by a different sub-detector (see Event t0 documentation). For an ECL
1263cluster produced at the interaction point in time with the event, the cluster time should be consistent with zero
1264within the uncertainties. Special values are returned if the fit for the Photon timing fails (see
1265documentation for `clusterHasFailedTiming`). (For MC, the calibrations and corrections are not fully simulated).
1266
1267.. note::
1268 | Please read `this <importantNoteECL>` first.
1269 | Lower limit: :math:`-1000.0`
1270 | Upper limit: :math:`1000.0`
1271 | Precision: :math:`12` bit
1272..
1273
1274**In Belle:**
1275Returns the trigger cell (TC) time of the ECL cluster (photon).
1276This information is available only in Belle data since experiment 31, and not available in Belle MC.
1277Clusters produced at the interaction point in time with the event, have TC time in the range of 9000-11000
1278Calculated based on the Appendix of Belle note 831.
1279
1280.. note::
1281 | In case this variable is obtained from Belle data that is stored in Belle II mdst/udst format, it will be truncated to:
1282 | Lower limit: :math:`-1000.0`
1283 | Upper limit: :math:`1000.0`
1284 | Precision: :math:`12` bit
1285..
1286
1287)DOC","ns");
1288 REGISTER_VARIABLE("clusterHasFailedTiming", eclClusterHasFailedTiming, R"DOC(
1289Status bit for if the ECL cluster's timing fit failed. Photon timing is given by the fitted time
1290of the recorded waveform of the highest energetic crystal in a cluster; however, that fit can fail and so
1291this variable tells the user if that has happened.
1292)DOC");
1293 REGISTER_VARIABLE("clusterErrorTiming", eclClusterErrorTiming, R"DOC(
1294Returns ECL cluster's timing uncertainty that contains :math:`99\%` of true photons (dt99).
1295
1296The photon timing uncertainty is currently determined using MC. The resulting parametrization depends on
1297the true energy deposition in the highest energetic crystal and the local beam background level in that crystal.
1298The resulting timing distribution is non-Gaussian and for each photon the value dt99 is stored,
1299where :math:`|\text{timing}| / \text{dt99} < 1` is designed to give a :math:`99\%`
1300timing efficiency for true photons from the IP.
1301The resulting efficiency is approximately flat in energy and independent of beam background levels.
1302
1303Very large values of dt99 are an indication of failed waveform fits in the ECL.
1304We remove such clusters in most physics photon lists.
1305
1306.. note::
1307 | Please read `this <importantNoteECL>` first.
1308 | Lower limit: :math:`0.0`
1309 | Upper limit: :math:`1000.0`
1310 | Precision: :math:`12` bit
1311
1312.. warning::
1313 In real data there will be a sizeable number of high energetic Bhabha events
1314 (from previous or later bunch collisions) that can easily be rejected by timing cuts.
1315 However, these events create large ECL clusters that can overlap with other ECL clusters
1316 and it is not clear that a simple rejection is the correction strategy.
1317..
1318
1319)DOC","ns");
1320 REGISTER_VARIABLE("clusterHasFailedErrorTiming", eclClusterHasFailedErrorTiming, R"DOC(
1321Status bit for if the ECL cluster's timing uncertainty calculation failed. Photon timing is given by the fitted time
1322of the recorded waveform of the highest energetic crystal in a cluster; however, that fit can fail and so
1323this variable tells the user if that has happened.
1324)DOC");
1325 REGISTER_VARIABLE("clusterHighestE", eclClusterHighestE, R"DOC(
1326Returns energy of the highest energetic crystal in the ECL cluster after reweighting.
1327
1328.. warning::
1329 This variable must be used carefully since it can bias shower selection
1330 towards photons that hit crystals in the center and hence have a large energy
1331 deposition in the highest energy crystal.
1332
1333.. note::
1334 | Please read `this <importantNoteECL>` first.
1335 | Lower limit: :math:`-5` (:math:`e^{-5} = 0.00674\,` GeV)
1336 | Upper limit: :math:`3.0` (:math:`e^3 = 20.08553\,` GeV)
1337 | Precision: :math:`18` bit
1338..
1339
1340)DOC","GeV");
1341 REGISTER_VARIABLE("clusterCellID", eclClusterCellId,
1342 "Returns cellId of the crystal with highest energy in the ECLCluster.");
1343 REGISTER_VARIABLE("clusterThetaID", eclClusterThetaId,
1344 "Returns thetaId of the crystal with highest energy in the ECLCluster.");
1345 REGISTER_VARIABLE("clusterPhiID", eclClusterPhiId,
1346 "Returns phiId of the crystal with highest energy in the ECLCluster.");
1347 REGISTER_VARIABLE("clusterE1E9", eclClusterE1E9, R"DOC(
1348Returns ratio of energies of the central crystal, E1, and 3x3 crystals, E9, around the central crystal.
1349Since :math:`E1 \leq E9`, this ratio is :math:`\leq 1` and tends towards larger values for photons
1350and smaller values for hadrons.
1351
1352.. note::
1353 | Please read `this <importantNoteECL>` first.
1354 | Lower limit: :math:`0.0`
1355 | Upper limit: :math:`1.0`
1356 | Precision: :math:`10` bit
1357)DOC");
1358 REGISTER_VARIABLE("clusterE9E25", eclClusterE9E25, R"DOC(
1359Deprecated - kept for backwards compatibility - returns clusterE9E21.
1360)DOC");
1361 REGISTER_VARIABLE("clusterE9E21", eclClusterE9E21, R"DOC(
1362Returns ratio of energies in inner 3x3 crystals, E9, and 5x5 crystals around the central crystal without corners.
1363Since :math:`E9 \leq E21`, this ratio is :math:`\leq 1` and tends towards larger values for photons
1364and smaller values for hadrons.
1365
1366.. note::
1367 | Please read `this <importantNoteECL>` first.
1368 | Lower limit: :math:`0.0`
1369 | Upper limit: :math:`1.0`
1370 | Precision: :math:`10` bit
1371)DOC");
1372 REGISTER_VARIABLE("clusterAbsZernikeMoment40", eclClusterAbsZernikeMoment40, R"DOC(
1373Returns absolute value of Zernike moment 40 (:math:`|Z_{40}|`). (shower shape variable).
1374
1375.. note::
1376 | Please read `this <importantNoteECL>` first.
1377 | Lower limit: :math:`0.0`
1378 | Upper limit: :math:`1.7`
1379 | Precision: :math:`10` bit
1380)DOC");
1381 REGISTER_VARIABLE("clusterAbsZernikeMoment51", eclClusterAbsZernikeMoment51, R"DOC(
1382Returns absolute value of Zernike moment 51 (:math:`|Z_{51}|`). (shower shape variable).
1383
1384.. note::
1385 | Please read `this <importantNoteECL>` first.
1386 | Lower limit: :math:`0.0`
1387 | Upper limit: :math:`1.2`
1388 | Precision: :math:`10` bit
1389)DOC");
1390 REGISTER_VARIABLE("clusterZernikeMVA", eclClusterZernikeMVA, R"DOC(
1391Returns output of a MVA using eleven Zernike moments of the cluster. Zernike moments are calculated per
1392shower in a plane perpendicular to the shower direction via
1393
1394.. math::
1395 |Z_{nm}| = \frac{n+1}{\pi} \frac{1}{\sum_{i} w_{i} E_{i}} \left|\sum_{i} R_{nm}(\rho_{i}) e^{-im\alpha_{i}} w_{i} E_{i} \right|
1396
1397where n, m are the integers, :math:`i` runs over the crystals in the shower,
1398:math:`E_{i}` is the energy of the i-th crystal in the shower,
1399:math:`R_{nm}` is a polynomial of degree :math:`n`,
1400:math:`\rho_{i}` is the radial distance of the :math:`i`-th crystal in the perpendicular plane,
1401and :math:`\alpha_{i}` is the polar angle of the :math:`i`-th crystal in the perpendicular plane.
1402As a crystal can be related to more than one shower, :math:`w_{i}` is the fraction of the
1403energy of the :math:`i`-th crystal associated with the shower.
1404
1405More details about the implementation can be found in `BELLE2-NOTE-TE-2017-001 <https://docs.belle2.org/record/454?ln=en>`_ .
1406
1407More details about Zernike polynomials can be found in `Wikipedia <https://en.wikipedia.org/wiki/Zernike_polynomials>`_ .
1408
1409| For cluster with hypothesisId==N1: raw MVA output.
1410| For cluster with hypothesisId==N2: 1 - prod{clusterZernikeMVA}, where the product is on all N1 showers
1411 belonging to the same connected region (shower shape variable).
1412
1413.. note::
1414 | Please read `this <importantNoteECL>` first.
1415 | Lower limit: :math:`0.0`
1416 | Upper limit: :math:`1.0`
1417 | Precision: :math:`10` bit
1418)DOC");
1419 REGISTER_VARIABLE("clusterSecondMoment", eclClusterSecondMoment, R"DOC(
1420Returns second moment :math:`S`. It is defined as:
1421
1422.. math::
1423 S = \frac{1}{S_{0}(\theta)}\frac{\sum_{i=0}^{n} w_{i} E_{i} r^2_{i}}{\sum_{i=0}^{n} w_{i} E_{i}}
1424
1425where :math:`E_{i} = (E_0, E_1, ...)` are the single crystal energies sorted by energy, :math:`w_{i}` is
1426the crystal weight, and :math:`r_{i}` is the distance of the :math:`i`-th digit to the shower center projected
1427to a plane perpendicular to the shower axis. :math:`S_{0}(\theta)` normalizes :math:`S` to 1 for isolated photons.
1428
1429.. note::
1430 | Please read `this <importantNoteECL>` first.
1431 | Lower limit: :math:`0.0`
1432 | Upper limit: :math:`40.0`
1433 | Precision: :math:`10` bit
1434..
1435
1436)DOC","dimensionless");
1437 REGISTER_VARIABLE("clusterLAT", eclClusterLAT, R"DOC(
1438Returns lateral energy distribution (shower variable). It is defined as following:
1439
1440.. math::
1441 S = \frac{\sum_{i=2}^{n} w_{i} E_{i} r^2_{i}}{(w_{0} E_{0} + w_{1} E_{1}) r^2_{0} + \sum_{i=2}^{n} w_{i} E_{i} r^2_{i}}
1442
1443where :math:`E_{i} = (E_{0}, E_{1}, ...)` are the single crystal energies sorted by energy
1444(:math:`E_{0}` is the highest energy and :math:`E_{1}` the second highest), :math:`w_{i}`
1445is the crystal weight, :math:`r_{i}` is the distance of the :math:`i`-th digit to the
1446shower center projected to a plane perpendicular to the shower axis,
1447and :math:`r_{0} \approx 5\,cm` is the distance between two crystals.
1448
1449clusterLAT peaks around 0.3 for radially symmetrical electromagnetic showers and is larger
1450for hadronic events, and electrons with a close-by radiative or Bremsstrahlung photon.
1451
1452.. note::
1453 | Please read `this <importantNoteECL>` first.
1454 | Lower limit: :math:`0.0`
1455 | Upper limit: :math:`1.0`
1456 | Precision: :math:`10` bit
1457)DOC");
1458 REGISTER_VARIABLE("clusterNHits", eclClusterNHits, R"DOC(
1459Returns sum of weights :math:`w_{i}` (:math:`w_{i} \leq 1`) of all crystals in an ECL cluster.
1460For non-overlapping clusters this is equal to the number of crystals in the cluster.
1461In case of energy splitting among nearby clusters, this can be a non-integer value.
1462
1463.. note::
1464 | Please read `this <importantNoteECL>` first.
1465 | Lower limit: :math:`0.0`
1466 | Upper limit: :math:`200.0`
1467 | Precision: :math:`10` bit
1468 | If fractional weights are not of interest, this value should be cast to the nearest integer.
1469)DOC");
1470 REGISTER_VARIABLE("clusterTrackMatch", eclClusterTrackMatched, R"DOC(
1471Returns 1.0 if at least one reconstructed charged track is matched to the ECL cluster.
1472
1473Every reconstructed charged track is extrapolated into the ECL.
1474Every ECL crystal that is crossed by the track extrapolation is marked.
1475Each ECL cluster that contains any marked crystal is matched to the track.
1476Multiple tracks can be matched to one cluster and multiple clusters can be matched to one track.
1477It is conceptually correct to have two tracks matched to the same cluster.
1478)DOC");
1479 REGISTER_VARIABLE("nECLClusterTrackMatches", nECLClusterTrackMatches, R"DOC(
1480Returns number of charged tracks matched to this cluster.
1481
1482.. note::
1483 Sometimes (perfectly correctly) two tracks are extrapolated into the same cluster.
1484
1485 - For charged particles, this should return at least 1 (but sometimes 2 or more).
1486 - For neutrals, this should always return 0.
1487 - Returns NaN if there is no cluster.
1488)DOC");
1489 REGISTER_VARIABLE("clusterHasPulseShapeDiscrimination", eclClusterHasPulseShapeDiscrimination, R"DOC(
1490Status bit to indicate if cluster has digits with waveforms that passed energy and :math:`\chi^2`
1491thresholds for computing PSD variables.
1492)DOC");
1493 REGISTER_VARIABLE("beamBackgroundSuppression", beamBackgroundSuppression, R"DOC(
1494Returns the output of an MVA classifier that uses shower-related variables to distinguish true photon clusters from beam background clusters.
1495Class 1 is for true photon clusters while class 0 is for beam background clusters.
1496
1497The MVA has been trained using MC and the features used are:
1498
1499- `clusterTiming`
1500- `clusterPulseShapeDiscriminationMVA`
1501- `clusterE`
1502- `clusterTheta`
1503- `clusterZernikeMVA`
1504
1505Both run-dependent and run-independent weights are available. For more information on this, and for usage recommendations, please see
1506the `Neutrals Performance XWiki Page <https://xwiki.desy.de/xwiki/rest/p/e23c8>`_.
1507)DOC");
1508 REGISTER_VARIABLE("fakePhotonSuppression", fakePhotonSuppression, R"DOC(
1509Returns the output of an MVA classifier that uses shower-related variables to distinguish true photon clusters from fake photon clusters (e.g. split-offs,
1510track-cluster matching failures etc.). Class 1 is for true photon clusters while class 0 is for fake photon clusters.
1511
1512The MVA has been trained using MC and the features are:
1513
1514- `clusterPulseShapeDiscriminationMVA`
1515- `minC2TDist`
1516- `clusterZernikeMVA`
1517- `clusterE`
1518- `clusterTiming`
1519- `clusterTheta`
1520
1521This MVA is the same as the one used for `hadronicSplitOffSuppression` but that variable should not be used as it is deprecated and does not use the new weights.
1522
1523Both run-dependent and run-independent weights are available. For more information on this, and for usage recommendations, please see
1524the `Neutrals Performance XWiki Page <https://xwiki.desy.de/xwiki/rest/p/e23c8>`_.
1525)DOC");
1526 REGISTER_VARIABLE("hadronicSplitOffSuppression", hadronicSplitOffSuppression, R"DOC(
1527Returns the output of an MVA classifier that uses shower-related variables to distinguish true photon clusters from hadronic splitoff clusters.
1528The classes are:
1529
1530- 1 for true photon clusters
1531- 0 for hadronic splitoff clusters
1532
1533The MVA has been trained using samples of signal photons and hadronic splitoff photons coming from MC. The features used are (in decreasing order of significance):
1534
1535- `clusterPulseShapeDiscriminationMVA`
1536- `minC2TDist`
1537- `clusterZernikeMVA`
1538- `clusterE`
1539- `clusterLAT`
1540- `clusterE1E9`
1541- `clusterSecondMoment`
1542)DOC");
1543 MAKE_DEPRECATED("hadronicSplitOffSuppression", false, "light-2302-genetta", R"DOC(
1544 Use the variable `fakePhotonSuppression` instead which is maintained and uses the latest weight files.
1545)DOC");
1546 REGISTER_VARIABLE("clusterKlId", eclClusterKlId, R"DOC(
1547Returns MVA classifier that uses ECL clusters variables to discriminate Klong clusters from em background.
1548
1549- 1 for Kl
1550- 0 for background
1551)DOC");
1552 REGISTER_VARIABLE("clusterPulseShapeDiscriminationMVA", eclPulseShapeDiscriminationMVA, R"DOC(
1553Returns MVA classifier that uses pulse shape discrimination to identify electromagnetic vs hadronic showers.
1554
1555- 1 for electromagnetic showers
1556- 0 for hadronic showers
1557)DOC");
1558 REGISTER_VARIABLE("clusterNumberOfHadronDigits", eclClusterNumberOfHadronDigits, R"DOC(
1559Returns ECL cluster's number of hadron digits in cluster (pulse shape discrimination variable).
1560Weighted sum of digits in cluster with significant scintillation emission (:math:`> 3\,` MeV)
1561in the hadronic scintillation component.
1562Computed only using cluster digits with energy :math:`> 50\,` MeV and good offline waveform fit :math:`\chi^2`.
1563
1564.. note::
1565 | Please read `this <importantNoteECL>` first.
1566 | Lower limit: :math:`0.0`
1567 | Upper limit: :math:`255.0`
1568 | Precision: :math:`18` bit
1569)DOC");
1570 REGISTER_VARIABLE("clusterClusterID", eclClusterId, R"DOC(
1571Returns ECL cluster ID of this ECL cluster within the connected region (CR) to which it belongs to.
1572)DOC");
1573 REGISTER_VARIABLE("clusterHasNPhotons", eclClusterHasNPhotonsHypothesis, R"DOC(
1574Returns 1.0 if cluster has the 'N photons' hypothesis (historically called 'N1'),
15750.0 if not, and NaN if no cluster is associated to the particle.
1576)DOC");
1577 REGISTER_VARIABLE("clusterHasNeutralHadron", eclClusterHasNeutralHadronHypothesis, R"DOC(
1578Returns 1.0 if the cluster has the 'neutral hadrons' hypothesis (historically called 'N2'),
15790.0 if not, and NaN if no cluster is associated to the particle.
1580)DOC");
1581 REGISTER_VARIABLE("eclExtTheta", eclExtTheta, R"DOC(
1582Returns extrapolated :math:`\theta` of particle track associated to the cluster (if any). Requires module ECLTrackCalDigitMatch to be executed.
1583
1584)DOC","rad");
1585 REGISTER_VARIABLE("eclExtPhi", eclExtPhi, R"DOC(
1586Returns extrapolated :math:`\phi` of particle track associated to the cluster (if any). Requires module ECLTrackCalDigitMatch to be executed..
1587
1588)DOC","rad");
1589 REGISTER_VARIABLE("eclExtPhiId", eclExtPhiId, R"DOC(
1590Returns extrapolated :math:`\phi` ID of particle track associated to the cluster (if any). Requires module ECLTrackCalDigitMatch to be executed..
1591)DOC");
1592 REGISTER_VARIABLE("weightedAverageECLTime", weightedAverageECLTime, R"DOC(
1593Returns ECL weighted average time of all clusters (neutrals) and matched clusters (charged) of daughters
1594(of any generation) of the provided particle.
1595
1596)DOC", "ns");
1597 REGISTER_VARIABLE("maxWeightedDistanceFromAverageECLTime", maxWeightedDistanceFromAverageECLTime, R"DOC(
1598Returns maximum weighted distance between time of the cluster of a photon and the ECL average time, amongst
1599the clusters (neutrals) and matched clusters (charged) of daughters (of all generations) of the provided particle.
1600
1601)DOC", "ns");
1602 REGISTER_VARIABLE("clusterMdstIndex", eclClusterMdstIndex, R"DOC(
1603StoreArray index(0 - based) of the MDST ECLCluster (useful for track-based particles matched to a cluster).
1604)DOC");
1605
1606 REGISTER_VARIABLE("nECLOutOfTimeCrystals", nECLOutOfTimeCrystals, R"DOC(
1607[Eventbased] Returns the number of crystals (ECLCalDigits) that are out of time.
1608)DOC");
1609
1610 REGISTER_VARIABLE("nECLOutOfTimeCrystalsFWDEndcap", nECLOutOfTimeCrystalsFWDEndcap, R"DOC(
1611[Eventbased] Returns the number of crystals (ECLCalDigits) that are out of time in the forward endcap.
1612)DOC");
1613
1614 REGISTER_VARIABLE("nECLOutOfTimeCrystalsBarrel", nECLOutOfTimeCrystalsBarrel, R"DOC(
1615[Eventbased] Returns the number of crystals (ECLCalDigits) that are out of time in the barrel.
1616)DOC");
1617
1618 REGISTER_VARIABLE("nECLOutOfTimeCrystalsBWDEndcap", nECLOutOfTimeCrystalsBWDEndcap, R"DOC(
1619[Eventbased] Returns the number of crystals (ECLCalDigits) that are out of time in the backward endcap.
1620)DOC");
1621
1622 REGISTER_VARIABLE("nRejectedECLShowers", nRejectedECLShowers, R"DOC(
1623[Eventbased] Returns the number of showers in the ECL that do not become clusters.
1624)DOC");
1625
1626 REGISTER_VARIABLE("nRejectedECLShowersFWDEndcap", nRejectedECLShowersFWDEndcap, R"DOC(
1627[Eventbased] Returns the number of showers in the ECL that do not become clusters, from the forward endcap.
1628If the number exceeds 255 (uint8_t maximum value) the variable is set to 255.
1629)DOC");
1630
1631 REGISTER_VARIABLE("nRejectedECLShowersBarrel", nRejectedECLShowersBarrel, R"DOC(
1632[Eventbased] Returns the number of showers in the ECL that do not become clusters, from the barrel.
1633If the number exceeds 255 (uint8_t maximum value) the variable is set to 255.
1634)DOC");
1635
1636 REGISTER_VARIABLE("nRejectedECLShowersBWDEndcap", nRejectedECLShowersBWDEndcap, R"DOC(
1637[Eventbased] Returns the number of showers in the ECL that do not become clusters, from the backward endcap.
1638If the number exceeds 255 (uint8_t maximum value) the variable is set to 255.
1639)DOC");
1640
1641 REGISTER_VARIABLE("nKLMMultistripHitsFWDEndcap", nKLMMultistripHitsFWDEndcap, R"DOC(
1642[Eventbased] Returns the number of multi-strip hits in the KLM forward endcap.
1643)DOC");
1644
1645 REGISTER_VARIABLE("nKLMMultistripHitsBarrel", nKLMMultistripHitsBarrel, R"DOC(
1646[Eventbased] Returns the number of multi-strip hits in the KLM barrel.
1647)DOC");
1648
1649 REGISTER_VARIABLE("nKLMMultistripHitsBWDEndcap", nKLMMultistripHitsBWDEndcap, R"DOC(
1650[Eventbased] Returns the number of multi-strip hits in the KLM backward endcap.
1651)DOC");
1652
1653 REGISTER_VARIABLE("nKLMMultistripHits", nKLMMultistripHits, R"DOC(
1654[Eventbased] Returns the number of multi-strip hits in the KLM.
1655)DOC");
1656
1657 REGISTER_VARIABLE("nECLShowersFWDEndcap", nECLShowersFWDEndcap, R"DOC(
1658[Eventbased] Returns the number of ECLShowers in the forward endcap.
1659)DOC");
1660
1661 REGISTER_VARIABLE("nECLShowersBarrel", nECLShowersBarrel, R"DOC(
1662[Eventbased] Returns the number of ECLShowers in the barrel.
1663)DOC");
1664
1665 REGISTER_VARIABLE("nECLShowersBWDEndcap", nECLShowersBWDEndcap, R"DOC(
1666[Eventbased] Returns the number of ECLShowers in the backward endcap.
1667)DOC");
1668
1669 REGISTER_VARIABLE("nECLShowers", nECLShowers, R"DOC(
1670[Eventbased] Returns the number of ECLShowers.
1671)DOC");
1672
1673 REGISTER_VARIABLE("nECLLocalMaximumsFWDEndcap", nECLLocalMaximumsFWDEndcap, R"DOC(
1674[Eventbased] Returns the number of LocalMaximums in the ECL forward endcap.
1675)DOC");
1676
1677 REGISTER_VARIABLE("nECLLocalMaximumsBarrel", nECLLocalMaximumsBarrel, R"DOC(
1678[Eventbased] Returns the number of LocalMaximums in the ECL barrel.
1679)DOC");
1680
1681 REGISTER_VARIABLE("nECLLocalMaximumsBWDEndcap", nECLLocalMaximumsBWDEndcap, R"DOC(
1682[Eventbased] Returns the number of LocalMaximums in the ECL backward endcap.
1683)DOC");
1684
1685 REGISTER_VARIABLE("nECLLocalMaximums", nECLLocalMaximums, R"DOC(
1686[Eventbased] Returns the number of LocalMaximums in the ECL.
1687)DOC");
1688
1689 REGISTER_VARIABLE("nECLTriggerCellsFWDEndcap", nECLTriggerCellsFWDEndcap, R"DOC(
1690[Eventbased] Returns the number of ECL trigger cells above 100 MeV in the forward endcap.
1691)DOC");
1692
1693 REGISTER_VARIABLE("nECLTriggerCellsBarrel", nECLTriggerCellsBarrel, R"DOC(
1694[Eventbased] Returns the number of ECL trigger cells above 100 MeV in the barrel.
1695)DOC");
1696
1697 REGISTER_VARIABLE("nECLTriggerCellsBWDEndcap", nECLTriggerCellsBWDEndcap, R"DOC(
1698[Eventbased] Returns the number of ECL trigger cells above 100 MeV in the backward endcap.
1699)DOC");
1700
1701 REGISTER_VARIABLE("nECLTriggerCells", nECLTriggerCells, R"DOC(
1702[Eventbased] Returns the number of ECL trigger cells above 100 MeV.
1703)DOC");
1704
1705 REGISTER_VARIABLE("eclClusterOnlyInvariantMass", eclClusterOnlyInvariantMass, R"DOC(
1706[Expert] The invariant mass calculated from all ECLCluster daughters (i.e. photons) and
1707cluster-matched tracks using the cluster 4-momenta.
1708
1709Used for ECL-based dark sector physics and debugging track-cluster matching.
1710
1711)DOC","GeV/:math:`\\text{c}^2`");
1712
1713 REGISTER_METAVARIABLE("photonHasOverlap(cutString, photonlistname, tracklistname)", photonHasOverlap, R"DOC(
1714 Returns true if the connected ECL region of the particle's cluster is shared by another particle's cluster.
1715 Neutral and charged cluster are considered.
1716 A cut string can be provided to ignore cluster that do not satisfy the given criteria.
1717 By default, the ParticleList ``gamma:all`` is used for the check of neutral ECL cluster,
1718 and the ParticleList ``e-:all`` for the check of charged ECL cluster.
1719 However, one can customize the name of the ParticleLists via additional arguments.
1720 If no argument or only a cut string is provided and ``gamma:all`` or ``e-:all`` does not exist
1721 or if the variable is requested for a particle that is not a photon, NaN is returned.
1722 )DOC", Manager::VariableDataType::c_double);
1723
1724 REGISTER_VARIABLE("clusterUncorrE", eclClusterUncorrectedE, R"DOC(
1725[Expert] [Calibration] Returns ECL cluster's uncorrected energy. That is, before leakage corrections.
1726This variable should only be used for study of the ECL. Please see :b2:var:`clusterE`.
1727
1728)DOC","GeV");
1729
1730 REGISTER_VARIABLE("distanceToMcKl",distanceToMcKl,R"DOC(
1731Returns the distance to the nearest truth KL particle, extrapolated to the cluster radius. To use
1732this variable, it is required to run getNeutralHadronGeomMatches function. Optionally, it can return
1733negative values to indicate that the ECL cluster should be removed from the analysis to correct for data
1734to MC difference in KL efficiency.
1735
1736)DOC", "cm");
1737
1738 REGISTER_VARIABLE("distanceToMcNeutron",distanceToMcNeutron,R"DOC(
1739Returns the distance to the nearest truth (anti)neutron, extrapolated to the cluster radius. To use
1740this variable, it is required to run getNeutralHadronGeomMatches function. Optionally, it can return
1741negative values to indicate that the ECL cluster should be removed from the analysis to correct for data
1742to MC difference in KL efficiency.
1743
1744)DOC", "cm");
1745
1746 REGISTER_VARIABLE("mdstIndexMcKl",mdstIndexMcKl,R"DOC(
1747 Returns the mdst index of the nearest truth KL, extrapolated to the cluster radius, if it is
1748 within the matching cone. To use this variable, it is required to run getNeutralHadronGeomMatches function.
1749)DOC");
1750
1751 REGISTER_VARIABLE("mdstIndexMcNeutron",mdstIndexMcNeutron,R"DOC(
1752 Returns the mdst index of the nearest truth (anti)neutron, extrapolated to the cluster radius, if it is
1753 within the matching cone. To use this variable, it is required to run getNeutralHadronGeomMatches function.
1754)DOC");
1755
1756 }
1758}
R E
internal precision of FFTW codelets
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:618
static const double doubleNaN
quiet_NaN
Definition: Const.h:703
static const ParticleType photon
photon particle
Definition: Const.h:673
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
Definition: ECLCluster.h:31
@ c_nPhotons
CR is split into n photons (N1)
@ c_neutralHadron
CR is reconstructed as a neutral hadron (N2)
static std::unique_ptr< GeneralCut > compile(const std::string &cut)
Creates an instance of a cut and returns a unique_ptr to it, if you need a copy-able object instead y...
Definition: GeneralCut.h:84
static const ReferenceFrame & GetCurrent()
Get current rest frame.
std::function< VarVariant(const Particle *)> FunctionPtr
functions stored take a const Particle* and return VarVariant.
Definition: Manager.h:113
const Var * getVariable(std::string name)
Get the variable belonging to the given key.
Definition: Manager.cc:57
static Manager & Instance()
get singleton instance.
Definition: Manager.cc:25
#define MAKE_DEPRECATED(name, make_fatal, version, description)
Registers a variable as deprecated.
Definition: Manager.h:443
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.