10#include <analysis/variables/Variables.h>
13#include <analysis/VariableManager/Manager.h>
15#include <analysis/utility/PCmsLabTransform.h>
16#include <analysis/utility/ReferenceFrame.h>
17#include <analysis/utility/MCMatching.h>
18#include <analysis/ClusterUtility/ClusterUtils.h>
21#include <framework/datastore/StoreArray.h>
22#include <framework/datastore/StoreObjPtr.h>
25#include <analysis/dataobjects/Particle.h>
26#include <analysis/dataobjects/ParticleList.h>
27#include <framework/dataobjects/EventExtraInfo.h>
28#include <analysis/dataobjects/EventShapeContainer.h>
30#include <mdst/dataobjects/MCParticle.h>
31#include <mdst/dataobjects/Track.h>
32#include <mdst/dataobjects/ECLCluster.h>
33#include <mdst/dataobjects/V0.h>
35#include <mdst/dbobjects/BeamSpot.h>
38#include <framework/logging/Logger.h>
39#include <framework/geometry/B2Vector3.h>
40#include <framework/geometry/BFieldManager.h>
41#include <framework/gearbox/Const.h>
42#include <framework/utilities/Conversion.h>
44#include <Math/Vector4D.h>
50#include <boost/algorithm/string.hpp>
64 double particleP(
const Particle* part)
67 return frame.getMomentum(part).P();
70 double particleE(
const Particle* part)
73 return frame.getMomentum(part).E();
76 double particleClusterEUncertainty(
const Particle* part)
78 const ECLCluster* cluster = part->getECLCluster();
81 const auto EPhiThetaCov = clutls.GetCovarianceMatrix3x3FromCluster(cluster);
82 return sqrt(fabs(EPhiThetaCov[0][0]));
87 double particlePx(
const Particle* part)
90 return frame.getMomentum(part).Px();
93 double particlePy(
const Particle* part)
96 return frame.getMomentum(part).Py();
99 double particlePz(
const Particle* part)
102 return frame.getMomentum(part).Pz();
105 double particlePt(
const Particle* part)
108 return frame.getMomentum(part).Pt();
111 double covMatrixElement(
const Particle* part,
const std::vector<double>& element)
113 int elementI = int(std::lround(element[0]));
114 int elementJ = int(std::lround(element[1]));
116 if (elementI < 0 || elementI > 6) {
117 B2WARNING(
"Requested particle's momentumVertex covariance matrix element is out of boundaries [0 - 6]:" <<
LogVar(
"i", elementI));
120 if (elementJ < 0 || elementJ > 6) {
121 B2WARNING(
"Requested particle's momentumVertex covariance matrix element is out of boundaries [0 - 6]:" <<
LogVar(
"j", elementJ));
125 return part->getMomentumVertexErrorMatrix()(elementI, elementJ);
128 double particleEUncertainty(
const Particle* part)
132 double errorSquared = frame.getMomentumErrorMatrix(part)(3, 3);
134 if (errorSquared >= 0.0)
135 return sqrt(errorSquared);
140 double particlePErr(
const Particle* part)
142 TMatrixD jacobianRot(3, 3);
145 double cosPhi = cos(particlePhi(part));
146 double sinPhi = sin(particlePhi(part));
147 double cosTheta = particleCosTheta(part);
148 double sinTheta = sin(acos(cosTheta));
149 double p = particleP(part);
151 jacobianRot(0, 0) = sinTheta * cosPhi;
152 jacobianRot(0, 1) = sinTheta * sinPhi;
153 jacobianRot(1, 0) = cosTheta * cosPhi / p;
154 jacobianRot(1, 1) = cosTheta * sinPhi / p;
155 jacobianRot(0, 2) = cosTheta;
156 jacobianRot(2, 0) = -sinPhi / sinTheta / p;
157 jacobianRot(1, 2) = -sinTheta / p;
158 jacobianRot(2, 1) = cosPhi / sinTheta / p;
162 double errorSquared = frame.getMomentumErrorMatrix(part).GetSub(0, 2, 0, 2,
" ").Similarity(jacobianRot)(0, 0);
164 if (errorSquared >= 0.0)
165 return sqrt(errorSquared);
170 double particlePxErr(
const Particle* part)
174 double errorSquared = frame.getMomentumErrorMatrix(part)(0, 0);
176 if (errorSquared >= 0.0)
177 return sqrt(errorSquared);
182 double particlePyErr(
const Particle* part)
185 double errorSquared = frame.getMomentumErrorMatrix(part)(1, 1);
187 if (errorSquared >= 0.0)
188 return sqrt(errorSquared);
193 double particlePzErr(
const Particle* part)
196 double errorSquared = frame.getMomentumErrorMatrix(part)(2, 2);
198 if (errorSquared >= 0.0)
199 return sqrt(errorSquared);
204 double particlePtErr(
const Particle* part)
206 TMatrixD jacobianRot(3, 3);
209 double px = particlePx(part);
210 double py = particlePy(part);
211 double pt = particlePt(part);
213 jacobianRot(0, 0) = px / pt;
214 jacobianRot(0, 1) = py / pt;
215 jacobianRot(1, 0) = -py / (pt * pt);
216 jacobianRot(1, 1) = px / (pt * pt);
217 jacobianRot(2, 2) = 1;
220 double errorSquared = frame.getMomentumErrorMatrix(part).GetSub(0, 2, 0, 2,
" ").Similarity(jacobianRot)(0, 0);
222 if (errorSquared >= 0.0)
223 return sqrt(errorSquared);
228 double momentumDeviationChi2(
const Particle* part)
233 if (part->getPValue() < 0.0)
242 result += TMath::Power(part->getPx() - mcp->getMomentum().X(), 2.0) / part->getMomentumVertexErrorMatrix()(0, 0);
243 result += TMath::Power(part->getPy() - mcp->getMomentum().Y(), 2.0) / part->getMomentumVertexErrorMatrix()(1, 1);
244 result += TMath::Power(part->getPz() - mcp->getMomentum().Z(), 2.0) / part->getMomentumVertexErrorMatrix()(2, 2);
249 double particleTheta(
const Particle* part)
252 return frame.getMomentum(part).Theta();
255 double particleThetaErr(
const Particle* part)
257 TMatrixD jacobianRot(3, 3);
260 double cosPhi = cos(particlePhi(part));
261 double sinPhi = sin(particlePhi(part));
262 double cosTheta = particleCosTheta(part);
263 double sinTheta = sin(acos(cosTheta));
264 double p = particleP(part);
266 jacobianRot(0, 0) = sinTheta * cosPhi;
267 jacobianRot(0, 1) = sinTheta * sinPhi;
268 jacobianRot(1, 0) = cosTheta * cosPhi / p;
269 jacobianRot(1, 1) = cosTheta * sinPhi / p;
270 jacobianRot(0, 2) = cosTheta;
271 jacobianRot(2, 0) = -sinPhi / sinTheta / p;
272 jacobianRot(1, 2) = -sinTheta / p;
273 jacobianRot(2, 1) = cosPhi / sinTheta / p;
276 double errorSquared = frame.getMomentumErrorMatrix(part).GetSub(0, 2, 0, 2,
" ").Similarity(jacobianRot)(1, 1);
278 if (errorSquared >= 0.0)
279 return sqrt(errorSquared);
284 double particleCosTheta(
const Particle* part)
287 return cos(frame.getMomentum(part).Theta());
290 double particleCosThetaErr(
const Particle* part)
292 return fabs(particleThetaErr(part) * sin(particleTheta(part)));
295 double particlePhi(
const Particle* part)
298 return frame.getMomentum(part).Phi();
301 double particlePhiErr(
const Particle* part)
303 TMatrixD jacobianRot(3, 3);
306 double px = particlePx(part);
307 double py = particlePy(part);
308 double pt = particlePt(part);
310 jacobianRot(0, 0) = px / pt;
311 jacobianRot(0, 1) = py / pt;
312 jacobianRot(1, 0) = -py / (pt * pt);
313 jacobianRot(1, 1) = px / (pt * pt);
314 jacobianRot(2, 2) = 1;
317 double errorSquared = frame.getMomentumErrorMatrix(part).GetSub(0, 2, 0, 2,
" ").Similarity(jacobianRot)(1, 1);
319 if (errorSquared >= 0.0)
320 return sqrt(errorSquared);
325 double particleXp(
const Particle* part)
328 ROOT::Math::PxPyPzEVector p4 = part -> get4Vector();
329 ROOT::Math::PxPyPzEVector p4CMS = T.rotateLabToCms() * p4;
330 float s = T.getCMSEnergy();
331 float M = part->getMass();
332 return p4CMS.P() / TMath::Sqrt(s * s / 4 - M * M);
335 int particlePDGCode(
const Particle* part)
337 return part->getPDGCode();
340 double cosAngleBetweenMomentumAndVertexVectorInXYPlane(
const Particle* part)
342 static DBObjPtr<BeamSpot> beamSpotDB;
343 double px = part->getPx();
344 double py = part->getPy();
346 double xV = part->getX();
347 double yV = part->getY();
349 double xIP = (beamSpotDB->getIPPosition()).X();
350 double yIP = (beamSpotDB->getIPPosition()).Y();
355 double cosangle = (px * x + py * y) / (sqrt(px * px + py * py) * sqrt(x * x + y * y));
359 double cosAngleBetweenMomentumAndVertexVector(
const Particle* part)
361 static DBObjPtr<BeamSpot> beamSpotDB;
362 return cos((
B2Vector3D(part->getVertex()) - beamSpotDB->getIPPosition()).Angle(
B2Vector3D(part->getMomentum())));
365 double cosThetaBetweenParticleAndNominalB(
const Particle* part)
368 int particlePDG = abs(part->getPDGCode());
369 if (particlePDG != 511 and particlePDG != 521)
370 B2FATAL(
"The Variables cosThetaBetweenParticleAndNominalB is only meant to be used on B mesons!");
373 double e_Beam = T.getCMSEnergy() / 2.0;
374 double m_B = part->getPDGMass();
376 if (e_Beam * e_Beam - m_B * m_B < 0) {
377 e_Beam = 1.0579400E+1 / 2.0;
379 double p_B = std::sqrt(e_Beam * e_Beam - m_B * m_B);
381 ROOT::Math::PxPyPzEVector p = T.rotateLabToCms() * part->get4Vector();
386 double theta_BY = (2 * e_Beam * e_d - m_B * m_B - m_d * m_d)
391 double cosToThrustOfEvent(
const Particle* part)
393 StoreObjPtr<EventShapeContainer> evtShape;
395 B2WARNING(
"Cannot find thrust of event information, did you forget to load the event shape calculation?");
400 B2Vector3D particleMomentum = (T.rotateLabToCms() * part -> get4Vector()).Vect();
401 return std::cos(th.Angle(particleMomentum));
404 double ImpactXY(
const Particle* particle)
406 static DBObjPtr<BeamSpot> beamSpotDB;
408 ROOT::Math::XYZVector mom = particle->getMomentum();
410 ROOT::Math::XYZVector r = particle->getVertex() - ROOT::Math::XYZVector(beamSpotDB->getIPPosition());
415 double T = TMath::Sqrt(mom.Perp2() - 2.0 * curvature.Dot(r.Cross(mom)) + curvature.Mag2() * r.Perp2());
417 return TMath::Abs((-2 * r.Cross(mom).Z() + curvature.R() * r.Perp2()) / (T + mom.Rho()));
420 double ArmenterosLongitudinalMomentumAsymmetry(
const Particle* part)
423 int nDaughters = part -> getNDaughters();
425 B2FATAL(
"You are trying to use an Armenteros variable. The mother particle is required to have exactly two daughters");
427 const auto& daughters = part -> getDaughters();
428 B2Vector3D motherMomentum = frame.getMomentum(part).Vect();
429 B2Vector3D daughter1Momentum = frame.getMomentum(daughters[0]).Vect();
430 B2Vector3D daughter2Momentum = frame.getMomentum(daughters[1]).Vect();
432 int daughter1Charge = daughters[0] -> getCharge();
433 int daughter2Charge = daughters[1] -> getCharge();
434 double daughter1Ql = daughter1Momentum.Dot(motherMomentum) / motherMomentum.Mag();
435 double daughter2Ql = daughter2Momentum.Dot(motherMomentum) / motherMomentum.Mag();
438 if (daughter2Charge > daughter1Charge)
439 Arm_alpha = (daughter2Ql - daughter1Ql) / (daughter2Ql + daughter1Ql);
441 Arm_alpha = (daughter1Ql - daughter2Ql) / (daughter1Ql + daughter2Ql);
446 double ArmenterosDaughter1Qt(
const Particle* part)
449 int nDaughters = part -> getNDaughters();
451 B2FATAL(
"You are trying to use an Armenteros variable. The mother particle is required to have exactly two daughters.");
453 const auto& daughters = part -> getDaughters();
454 B2Vector3D motherMomentum = frame.getMomentum(part).Vect();
455 B2Vector3D daughter1Momentum = frame.getMomentum(daughters[0]).Vect();
456 double qt = daughter1Momentum.
Perp(motherMomentum);
461 double ArmenterosDaughter2Qt(
const Particle* part)
464 int nDaughters = part -> getNDaughters();
466 B2FATAL(
"You are trying to use an Armenteros variable. The mother particle is required to have exactly two daughters.");
468 const auto& daughters = part -> getDaughters();
469 B2Vector3D motherMomentum = frame.getMomentum(part).Vect();
470 B2Vector3D daughter2Momentum = frame.getMomentum(daughters[1]).Vect();
471 double qt = daughter2Momentum.
Perp(motherMomentum);
479 double particleMass(
const Particle* part)
481 return part->getMass();
484 double particleDMass(
const Particle* part)
486 return part->getMass() - part->getPDGMass();
489 double particleInvariantMassFromDaughters(
const Particle* part)
491 const std::vector<Particle*> daughters = part->getDaughters();
492 if (daughters.size() > 0) {
493 ROOT::Math::PxPyPzEVector sum;
494 for (
auto daughter : daughters)
495 sum += daughter->get4Vector();
499 return part->getMass();
503 double particleInvariantMassFromDaughtersDisplaced(
const Particle* part)
505 ROOT::Math::XYZVector vertex = part->getVertex();
506 if (part->getParticleSource() != Particle::EParticleSourceObject::c_V0
507 && vertex.Rho() < 0.5)
return particleInvariantMassFromDaughters(part);
509 const std::vector<Particle*> daughters = part->getDaughters();
510 if (daughters.size() == 0)
return particleInvariantMassFromDaughters(part);
513 ROOT::Math::PxPyPzMVector sum;
514 for (
auto daughter : daughters) {
515 const TrackFitResult* tfr = daughter->getTrackFitResult();
517 sum += daughter->get4Vector();
520 Helix helix = tfr->getHelix();
521 helix.passiveMoveBy(vertex);
522 double scalingFactor = daughter->getEffectiveMomentumScale();
523 double momX = scalingFactor * helix.getMomentumX(bField);
524 double momY = scalingFactor * helix.getMomentumY(bField);
525 double momZ = scalingFactor * helix.getMomentumZ(bField);
526 float mPDG = daughter->getPDGMass();
527 sum += ROOT::Math::PxPyPzMVector(momX, momY, momZ, mPDG);
532 double particleInvariantMassLambda(
const Particle* part)
534 const std::vector<Particle*> daughters = part->getDaughters();
535 if (daughters.size() == 2) {
536 ROOT::Math::PxPyPzEVector dt1;
537 ROOT::Math::PxPyPzEVector dt2;
538 ROOT::Math::PxPyPzEVector dtsum;
541 dt1 = daughters[0]->get4Vector();
542 dt2 = daughters[1]->get4Vector();
543 double E1 = hypot(mpi, dt1.P());
544 double E2 = hypot(mpr, dt2.P());
546 return sqrt((E1 + E2) * (E1 + E2) - dtsum.P() * dtsum.P());
549 return part->getMass();
553 double particleInvariantMassError(
const Particle* part)
555 float invMass = part->getMass();
556 TMatrixFSym covarianceMatrix = part->getMomentumErrorMatrix();
558 TVectorF jacobian(Particle::c_DimMomentum);
559 jacobian[0] = -1.0 * part->getPx() / invMass;
560 jacobian[1] = -1.0 * part->getPy() / invMass;
561 jacobian[2] = -1.0 * part->getPz() / invMass;
562 jacobian[3] = 1.0 * part->getEnergy() / invMass;
564 double result = jacobian * (covarianceMatrix * jacobian);
569 return TMath::Sqrt(result);
572 double particleInvariantMassSignificance(
const Particle* part)
574 return particleDMass(part) / particleInvariantMassError(part);
577 double particleMassSquared(
const Particle* part)
579 ROOT::Math::PxPyPzEVector p4 = part->get4Vector();
583 double b2bTheta(
const Particle* part)
586 ROOT::Math::PxPyPzEVector pcms = T.rotateLabToCms() * part->get4Vector();
587 ROOT::Math::PxPyPzEVector b2bcms(-pcms.Px(), -pcms.Py(), -pcms.Pz(), pcms.E());
588 ROOT::Math::PxPyPzEVector b2blab = T.rotateCmsToLab() * b2bcms;
589 return b2blab.Theta();
592 double b2bPhi(
const Particle* part)
595 ROOT::Math::PxPyPzEVector pcms = T.rotateLabToCms() * part->get4Vector();
596 ROOT::Math::PxPyPzEVector b2bcms(-pcms.Px(), -pcms.Py(), -pcms.Pz(), pcms.E());
597 ROOT::Math::PxPyPzEVector b2blab = T.rotateCmsToLab() * b2bcms;
601 double b2bClusterTheta(
const Particle* part)
604 const ECLCluster* cluster = part->getECLCluster();
610 ROOT::Math::PxPyPzEVector p4Cluster = clutls.Get4MomentumFromCluster(cluster, clusterHypothesis);
614 ROOT::Math::PxPyPzEVector pcms = T.rotateLabToCms() * p4Cluster;
615 ROOT::Math::PxPyPzEVector b2bcms(-pcms.Px(), -pcms.Py(), -pcms.Pz(), pcms.E());
616 ROOT::Math::PxPyPzEVector b2blab = T.rotateCmsToLab() * b2bcms;
617 return b2blab.Theta();
620 double b2bClusterPhi(
const Particle* part)
623 const ECLCluster* cluster = part->getECLCluster();
629 ROOT::Math::PxPyPzEVector p4Cluster = clutls.Get4MomentumFromCluster(cluster, clusterHypothesis);
633 ROOT::Math::PxPyPzEVector pcms = T.rotateLabToCms() * p4Cluster;
634 ROOT::Math::PxPyPzEVector b2bcms(-pcms.Px(), -pcms.Py(), -pcms.Pz(), pcms.E());
635 ROOT::Math::PxPyPzEVector b2blab = T.rotateCmsToLab() * b2bcms;
642 double particleQ(
const Particle* part)
644 float m = part->getMass();
645 for (
unsigned i = 0; i < part->getNDaughters(); i++) {
646 const Particle* child = part->getDaughter(i);
648 m -= child->getMass();
653 double particleDQ(
const Particle* part)
655 float m = part->getMass() - part->getPDGMass();
656 for (
unsigned i = 0; i < part->getNDaughters(); i++) {
657 const Particle* child = part->getDaughter(i);
659 m -= (child->getMass() - child->getPDGMass());
666 double particleMbc(
const Particle* part)
669 ROOT::Math::PxPyPzEVector vec = T.rotateLabToCms() * part->get4Vector();
670 double E = T.getCMSEnergy() / 2;
671 double m2 = E * E - vec.P2();
672 double mbc = m2 >= 0 ? sqrt(m2) : Const::doubleNaN;
676 double particleDeltaE(
const Particle* part)
679 ROOT::Math::PxPyPzEVector vec = T.rotateLabToCms() * part->get4Vector();
680 return vec.E() - T.getCMSEnergy() / 2;
686 void printParticleInternal(
const Particle* p,
int depth)
689 for (
int i = 0; i < depth; i++) {
692 s << p->getPDGCode();
693 const MCParticle* mcp = p->getMCParticle();
696 s <<
" -> MC: " << mcp->getPDG() <<
", mcErrors: " << flags <<
" ("
698 s <<
", mc-index " << mcp->getIndex();
700 s <<
" (no MC match)";
702 s <<
", mdst-source " << p->getMdstSource();
704 for (
const auto* daughter : p->getDaughters()) {
705 printParticleInternal(daughter, depth + 1);
709 bool printParticle(
const Particle* p)
711 printParticleInternal(p, 0);
716 double particleMCMomentumTransfer2(
const Particle* part)
724 ROOT::Math::PxPyPzEVector pB = mcB->get4Vector();
726 std::vector<MCParticle*> mcDaug = mcB->getDaughters();
733 ROOT::Math::PxPyPzEVector pX;
735 for (
auto mcTemp : mcDaug) {
736 if (abs(mcTemp->getPDG()) <= 16)
739 pX += mcTemp->get4Vector();
742 ROOT::Math::PxPyPzEVector q = pB - pX;
748 double recoilPx(
const Particle* particle)
753 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
757 return frame.getMomentum(pIN - particle->get4Vector()).Px();
760 double recoilPy(
const Particle* particle)
765 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
769 return frame.getMomentum(pIN - particle->get4Vector()).Py();
772 double recoilPz(
const Particle* particle)
777 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
781 return frame.getMomentum(pIN - particle->get4Vector()).Pz();
784 double recoilMomentum(
const Particle* particle)
789 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
793 return frame.getMomentum(pIN - particle->get4Vector()).P();
796 double recoilMomentumTheta(
const Particle* particle)
801 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
805 return frame.getMomentum(pIN - particle->get4Vector()).Theta();
808 double recoilMomentumPhi(
const Particle* particle)
813 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
817 return frame.getMomentum(pIN - particle->get4Vector()).Phi();
820 double recoilEnergy(
const Particle* particle)
825 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
829 return frame.getMomentum(pIN - particle->get4Vector()).E();
832 double recoilMass(
const Particle* particle)
837 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
841 return frame.getMomentum(pIN - particle->get4Vector()).M();
844 double recoilMassSquared(
const Particle* particle)
849 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
853 return frame.getMomentum(pIN - particle->get4Vector()).M2();
856 double m2RecoilSignalSide(
const Particle* part)
859 double beamEnergy = T.getCMSEnergy() / 2.;
861 ROOT::Math::PxPyPzEVector tagVec = T.rotateLabToCms() * part->getDaughter(0)->get4Vector();
862 ROOT::Math::PxPyPzEVector sigVec = T.rotateLabToCms() * part->getDaughter(1)->get4Vector();
863 tagVec.SetE(-beamEnergy);
864 return (-tagVec - sigVec).M2();
867 double recoilMCDecayType(
const Particle* particle)
869 auto* mcp = particle->getMCParticle();
874 MCParticle* mcMother = mcp->getMother();
879 std::vector<MCParticle*> daughters = mcMother->getDaughters();
881 if (daughters.size() != 2)
884 MCParticle* recoilMC =
nullptr;
885 if (daughters[0]->getArrayIndex() == mcp->getArrayIndex())
886 recoilMC = daughters[1];
888 recoilMC = daughters[0];
894 checkMCParticleDecay(recoilMC, decayType,
false);
897 checkMCParticleDecay(recoilMC, decayType,
true);
902 void checkMCParticleDecay(MCParticle* mcp,
int& decayType,
bool recursive)
904 int nHadronicParticles = 0;
905 int nPrimaryParticleDaughters = 0;
906 std::vector<MCParticle*> daughters = mcp->getDaughters();
909 for (
auto& daughter : daughters) {
913 nPrimaryParticleDaughters++;
915 nHadronicParticles++;
918 if (nPrimaryParticleDaughters > 1) {
919 for (
auto& daughter : daughters) {
923 if (abs(daughter->getPDG()) == 12 or abs(daughter->getPDG()) == 14 or abs(daughter->getPDG()) == 16) {
925 if (nHadronicParticles == 0) {
939 checkMCParticleDecay(daughter, decayType, recursive);
944 int nRemainingTracksInEvent(
const Particle* particle)
947 StoreArray<Track> tracks;
948 int event_tracks = tracks.getEntries();
951 const auto& daughters = particle->getFinalStateDaughters();
952 for (
const auto& daughter : daughters) {
953 int pdg = abs(daughter->getPDGCode());
958 return event_tracks - par_tracks;
961 double trackMatchType(
const Particle* particle)
966 const ECLCluster* cluster = particle->getECLCluster();
970 if (cluster->isTrack()) {
979 bool False(
const Particle*)
984 bool True(
const Particle*)
989 double infinity(
const Particle*)
991 double inf = std::numeric_limits<double>::infinity();
995 double random(
const Particle*)
997 return gRandom->Uniform(0, 1);
1000 double eventRandom(
const Particle*)
1002 std::string key =
"__eventRandom";
1003 StoreObjPtr<EventExtraInfo> eventExtraInfo;
1004 if (not eventExtraInfo.isValid())
1005 eventExtraInfo.create();
1006 if (eventExtraInfo->hasExtraInfo(key)) {
1007 return eventExtraInfo->getExtraInfo(key);
1009 double value = gRandom->Uniform(0, 1);
1010 eventExtraInfo->addExtraInfo(key, value);
1017 if (arguments.size() != 3 && arguments.size() != 4) {
1018 B2ERROR(
"Wrong number of arguments (3 or 4 required) for meta variable minET2ETDist");
1021 bool useHighestProbMassForExt(
true);
1022 if (arguments.size() == 4) {
1024 useHighestProbMassForExt =
static_cast<bool>(Belle2::convertString<int>(arguments[3]));
1025 }
catch (std::invalid_argument& e) {
1026 B2ERROR(
"Fourth (optional) argument of minET2ETDist must be an integer flag.");
1031 std::string detName = arguments[0];
1032 std::string detLayer = arguments[1];
1033 std::string referenceListName = arguments[2];
1034 std::string extraSuffix = (useHighestProbMassForExt) ?
"__useHighestProbMassForExt" :
"";
1036 std::string extraInfo =
"distToClosestTrkAt" + detName + detLayer +
"_VS_" + referenceListName + extraSuffix;
1038 auto func = [ = ](
const Particle * part) ->
double {
1039 auto dist = (part->hasExtraInfo(extraInfo)) ? part->getExtraInfo(extraInfo) :
Const::doubleNaN;
1050 if (arguments.size() != 4) {
1051 B2ERROR(
"Wrong number of arguments (4 required) for meta variable minET2ETDistVar");
1055 std::string detName = arguments[0];
1056 std::string detLayer = arguments[1];
1057 std::string referenceListName = arguments[2];
1058 std::string variableName = arguments[3];
1060 std::string extraInfo =
"idxOfClosestPartAt" + detName + detLayer +
"In_" + referenceListName;
1062 auto func = [ = ](
const Particle * part) ->
double {
1064 StoreObjPtr<ParticleList> refPartList(referenceListName);
1065 if (!refPartList.isValid())
1067 B2FATAL(
"Invalid Listname " << referenceListName <<
" given to minET2ETDistVar!");
1070 if (!part->hasExtraInfo(extraInfo))
1076 auto refPart = refPartList->getParticleWithMdstIdx(part->getExtraInfo(extraInfo));
1078 return std::get<double>(var->function(refPart));
1088 if (arguments.size() < 3) {
1089 B2ERROR(
"Wrong number of arguments (at least 3 required) for meta variable minET2ETIsoScore");
1093 std::string referenceListName = arguments[0];
1094 bool useHighestProbMassForExt;
1096 useHighestProbMassForExt =
static_cast<bool>(Belle2::convertString<int>(arguments[1]));
1097 }
catch (std::invalid_argument& e) {
1098 B2ERROR(
"Second argument must be an integer flag.");
1101 std::string extraSuffix = (useHighestProbMassForExt) ?
"__useHighestProbMassForExt" :
"";
1103 std::vector<std::string> detectorNames(arguments.begin() + 2, arguments.end());
1105 std::string detNamesConcat(
"");
1106 for (
auto& detName : detectorNames) {
1107 boost::to_upper(detName);
1108 detNamesConcat +=
"_" + detName;
1111 std::string extraInfo =
"trkIsoScore" + detNamesConcat +
"_VS_" + referenceListName + extraSuffix;
1113 auto func = [ = ](
const Particle * part) ->
double {
1115 StoreObjPtr<ParticleList> refPartList(referenceListName);
1116 if (!refPartList.isValid())
1118 B2FATAL(
"Invalid Listname " << referenceListName <<
" given to minET2ETIsoScore!");
1121 if (!part->hasExtraInfo(extraInfo))
1125 auto scoreDet = part->getExtraInfo(extraInfo);
1126 if (std::isnan(scoreDet))
1139 Manager::FunctionPtr particleExtTrkIsoScoreVarAsWeightedAvg(
const std::vector<std::string>& arguments)
1142 if (arguments.size() < 3) {
1143 B2ERROR(
"Wrong number of arguments (at least 3 required) for meta variable minET2ETIsoScoreAsWeightedAvg");
1147 std::string referenceListName = arguments[0];
1148 bool useHighestProbMassForExt;
1150 useHighestProbMassForExt =
static_cast<bool>(Belle2::convertString<int>(arguments[1]));
1151 }
catch (std::invalid_argument& e) {
1152 B2ERROR(
"Second argument must be an integer flag.");
1155 std::string extraSuffix = (useHighestProbMassForExt) ?
"__useHighestProbMassForExt" :
"";
1157 std::vector<std::string> detectorNames(arguments.begin() + 2, arguments.end());
1159 std::string detNamesConcat(
"");
1160 for (
auto& detName : detectorNames) {
1161 boost::to_upper(detName);
1162 detNamesConcat +=
"_" + detName;
1165 std::string extraInfo =
"trkIsoScoreAsWeightedAvg" + detNamesConcat +
"_VS_" + referenceListName + extraSuffix;
1167 auto func = [ = ](
const Particle * part) ->
double {
1169 StoreObjPtr<ParticleList> refPartList(referenceListName);
1170 if (!refPartList.isValid())
1172 B2FATAL(
"Invalid Listname " << referenceListName <<
" given to minET2ETIsoScoreAsWeightedAvg!");
1175 if (!part->hasExtraInfo(extraInfo))
1179 auto scoreDet = part->getExtraInfo(extraInfo);
1180 if (std::isnan(scoreDet))
1192 VARIABLE_GROUP(
"Kinematics");
1193 REGISTER_VARIABLE(
"p", particleP,
"momentum magnitude\n\n",
"GeV/c");
1194 REGISTER_VARIABLE(
"E", particleE,
"energy\n\n",
"GeV");
1196 REGISTER_VARIABLE(
"E_uncertainty", particleEUncertainty, R
"DOC(
1197 energy uncertainty (:math:`\sqrt{\sigma^2}`)
1200 REGISTER_VARIABLE(
"ECLClusterE_uncertainty", particleClusterEUncertainty,
1201 "energy uncertainty as given by the underlying ECL cluster\n\n",
"GeV");
1202 REGISTER_VARIABLE(
"px", particlePx,
"momentum component x\n\n",
"GeV/c");
1203 REGISTER_VARIABLE(
"py", particlePy,
"momentum component y\n\n",
"GeV/c");
1204 REGISTER_VARIABLE(
"pz", particlePz,
"momentum component z\n\n",
"GeV/c");
1205 REGISTER_VARIABLE(
"pt", particlePt,
"transverse momentum\n\n",
"GeV/c");
1206 REGISTER_VARIABLE(
"xp", particleXp,
1207 "scaled momentum: the momentum of the particle in the CMS as a fraction of its maximum available momentum in the collision");
1208 REGISTER_VARIABLE(
"pErr", particlePErr,
"error of momentum magnitude\n\n",
"GeV/c");
1209 REGISTER_VARIABLE(
"pxErr", particlePxErr,
"error of momentum component x\n\n",
"GeV/c");
1210 REGISTER_VARIABLE(
"pyErr", particlePyErr,
"error of momentum component y\n\n",
"GeV/c");
1211 REGISTER_VARIABLE(
"pzErr", particlePzErr,
"error of momentum component z\n\n",
"GeV/c");
1212 REGISTER_VARIABLE(
"ptErr", particlePtErr,
"error of transverse momentum\n\n",
"GeV/c");
1213 REGISTER_VARIABLE(
"momVertCovM(i,j)", covMatrixElement,
1214 "returns the (i,j)-th element of the MomentumVertex Covariance Matrix (7x7).\n"
1215 "Order of elements in the covariance matrix is: px, py, pz, E, x, y, z.\n\n",
"GeV/c, GeV/c, GeV/c, GeV, cm, cm, cm");
1216 REGISTER_VARIABLE(
"momDevChi2", momentumDeviationChi2, R
"DOC(
1217momentum deviation :math:`\chi^2` value calculated as :math:`\chi^2 = \sum_i (p_i - mc(p_i))^2/\sigma(p_i)^2`,
1218where :math:`\sum` runs over i = px, py, pz and :math:`mc(p_i)` is the mc truth value and :math:`\sigma(p_i)` is the estimated error of i-th component of momentum vector
1220 REGISTER_VARIABLE("theta", particleTheta,
"polar angle\n\n",
"rad");
1221 REGISTER_VARIABLE(
"thetaErr", particleThetaErr,
"error of polar angle\n\n",
"rad");
1222 REGISTER_VARIABLE(
"cosTheta", particleCosTheta,
"momentum cosine of polar angle");
1223 REGISTER_VARIABLE(
"cosThetaErr", particleCosThetaErr,
"error of momentum cosine of polar angle");
1224 REGISTER_VARIABLE(
"phi", particlePhi,
"momentum azimuthal angle\n\n",
"rad");
1225 REGISTER_VARIABLE(
"phiErr", particlePhiErr,
"error of momentum azimuthal angle\n\n",
"rad");
1226 REGISTER_VARIABLE(
"PDG", particlePDGCode,
"PDG code");
1227 REGISTER_VARIABLE(
"cosAngleBetweenMomentumAndVertexVectorInXYPlane",
1228 cosAngleBetweenMomentumAndVertexVectorInXYPlane,
1229 "cosine of the angle between momentum and vertex vector (vector connecting ip and fitted vertex) of this particle in xy-plane");
1230 REGISTER_VARIABLE(
"cosAngleBetweenMomentumAndVertexVector",
1231 cosAngleBetweenMomentumAndVertexVector,
1232 "cosine of the angle between momentum and vertex vector (vector connecting ip and fitted vertex) of this particle");
1233 REGISTER_VARIABLE(
"cosThetaBetweenParticleAndNominalB",
1234 cosThetaBetweenParticleAndNominalB,
1235 "cosine of the angle in CMS between momentum the particle and a nominal B particle. It is somewhere between -1 and 1 if only a massless particle like a neutrino is missing in the reconstruction.");
1236 REGISTER_VARIABLE(
"cosToThrustOfEvent", cosToThrustOfEvent,
1237 "Returns the cosine of the angle between the particle and the thrust axis of the event, as calculate by the EventShapeCalculator module. buildEventShape() must be run before calling this variable");
1239 REGISTER_VARIABLE(
"ImpactXY" , ImpactXY ,
"The impact parameter of the given particle in the xy plane\n\n",
"cm");
1241 REGISTER_VARIABLE(
"M", particleMass, R
"DOC(
1244Note that this is context-dependent variable and can take different values depending on the situation. This should be the "best"
1245value possible with the information provided.
1247- If this particle is track- or cluster- based, then this is the value of the mass hypothesis.
1248- If this particle is an MC particle then this is the mass of that particle.
1249- If this particle is composite, then *initially* this takes the value of the invariant mass of the daughters.
1250- If this particle is composite and a *mass or vertex fit* has been performed then this may be updated by the fit.
1252* You will see a difference between this mass and the :b2:var:`InvM`.
1254)DOC", "GeV/:math:`\\text{c}^2`");
1255 REGISTER_VARIABLE(
"dM", particleDMass,
"mass minus nominal mass\n\n",
"GeV/:math:`\\text{c}^2`");
1256 REGISTER_VARIABLE(
"Q", particleQ,
"energy released in decay\n\n",
"GeV");
1257 REGISTER_VARIABLE(
"dQ", particleDQ,
":b2:var:`Q` minus nominal energy released in decay\n\n",
"GeV");
1258 REGISTER_VARIABLE(
"Mbc", particleMbc,
"beam constrained mass\n\n",
"GeV/:math:`\\text{c}^2`");
1259 REGISTER_VARIABLE(
"deltaE", particleDeltaE,
"difference between :b2:var:`E` and half the center of mass energy\n\n",
"GeV");
1260 REGISTER_VARIABLE(
"M2", particleMassSquared,
"The particle's mass squared.\n\n",
":math:`[\\text{GeV}/\\text{c}^2]^2`");
1262 REGISTER_VARIABLE(
"InvM", particleInvariantMassFromDaughtersDisplaced,
1263 "invariant mass (determined from particle's daughter 4-momentum vectors). If this particle is V0 or decays at rho > 5 mm, its daughter 4-momentum vectors at fitted vertex are taken.\n"
1264 "If this particle has no daughters, defaults to :b2:var:`M`.\n\n",
"GeV/:math:`\\text{c}^2`");
1265 REGISTER_VARIABLE(
"InvMLambda", particleInvariantMassLambda,
1266 "Invariant mass (determined from particle's daughter 4-momentum vectors), assuming the first daughter is a pion and the second daughter is a proton.\n"
1267 "If the particle has not 2 daughters, it returns just the mass value.\n\n",
"GeV/:math:`\\text{c}^2`");
1269 REGISTER_VARIABLE(
"ErrM", particleInvariantMassError,
1270 "uncertainty of invariant mass\n\n",
"GeV/:math:`\\text{c}^2`");
1271 REGISTER_VARIABLE(
"SigM", particleInvariantMassSignificance,
1272 "signed deviation of particle's invariant mass from its nominal mass in units of the uncertainty on the invariant mass (:b2:var:`dM`/:b2:var:`ErrM`)");
1274 REGISTER_VARIABLE(
"pxRecoil", recoilPx,
1275 "component x of 3-momentum recoiling against given Particle\n\n",
"GeV/c");
1276 REGISTER_VARIABLE(
"pyRecoil", recoilPy,
1277 "component y of 3-momentum recoiling against given Particle\n\n",
"GeV/c");
1278 REGISTER_VARIABLE(
"pzRecoil", recoilPz,
1279 "component z of 3-momentum recoiling against given Particle\n\n",
"GeV/c");
1281 REGISTER_VARIABLE(
"pRecoil", recoilMomentum,
1282 "magnitude of 3 - momentum recoiling against given Particle\n\n",
"GeV/c");
1283 REGISTER_VARIABLE(
"pRecoilTheta", recoilMomentumTheta,
1284 "Polar angle of a particle's missing momentum\n\n",
"rad");
1285 REGISTER_VARIABLE(
"pRecoilPhi", recoilMomentumPhi,
1286 "Azimuthal angle of a particle's missing momentum\n\n",
"rad");
1287 REGISTER_VARIABLE(
"eRecoil", recoilEnergy,
1288 "energy recoiling against given Particle\n\n",
"GeV");
1289 REGISTER_VARIABLE(
"mRecoil", recoilMass,
1290 "Invariant mass of the system recoiling against given Particle\n\n",
"GeV/:math:`\\text{c}^2`");
1291 REGISTER_VARIABLE(
"m2Recoil", recoilMassSquared,
1292 "invariant mass squared of the system recoiling against given Particle\n\n",
":math:`[\\text{GeV}/\\text{c}^2]^2`");
1293 REGISTER_VARIABLE(
"m2RecoilSignalSide", m2RecoilSignalSide, R
"DOC(
1294 Squared recoil mass of the signal side which is calculated in the CMS frame under the assumption that the
1295 signal and tag side are produced back to back and the tag side energy equals the beam energy. The variable
1296 must be applied to the Upsilon and the tag side must be the first, the signal side the second daughter
1298 )DOC", ":math:`[\\text{GeV}/\\text{c}^2]^2`");
1300 REGISTER_VARIABLE(
"b2bTheta", b2bTheta,
1301 "Polar angle in the lab system that is back-to-back to the particle in the CMS. Useful for low multiplicity studies.\n\n",
"rad");
1302 REGISTER_VARIABLE(
"b2bPhi", b2bPhi,
1303 "Azimuthal angle in the lab system that is back-to-back to the particle in the CMS. Useful for low multiplicity studies.\n\n",
1305 REGISTER_VARIABLE(
"b2bClusterTheta", b2bClusterTheta,
1306 "Polar angle in the lab system that is back-to-back to the particle's associated ECLCluster in the CMS. Returns NAN if no cluster is found. Useful for low multiplicity studies.\n\n",
1308 REGISTER_VARIABLE(
"b2bClusterPhi", b2bClusterPhi,
1309 "Azimuthal angle in the lab system that is back-to-back to the particle's associated ECLCluster in the CMS. Returns NAN if no cluster is found. Useful for low multiplicity studies.\n\n",
1311 REGISTER_VARIABLE(
"ArmenterosLongitudinalMomentumAsymmetry", ArmenterosLongitudinalMomentumAsymmetry,
1312 "Longitudinal momentum asymmetry of V0's daughters.\n"
1313 "The mother (V0) is required to have exactly two daughters");
1314 REGISTER_VARIABLE(
"ArmenterosDaughter1Qt", ArmenterosDaughter1Qt, R
"DOC(
1315 Transverse momentum of the first daughter with respect to the V0 mother.
1316 The mother is required to have exactly two daughters
1319 REGISTER_VARIABLE(
"ArmenterosDaughter2Qt", ArmenterosDaughter2Qt, R
"DOC(
1320 Transverse momentum of the second daughter with respect to the V0 mother.
1321 The mother is required to have exactly two daughters
1325 VARIABLE_GROUP(
"Miscellaneous");
1326 REGISTER_VARIABLE(
"nRemainingTracksInEvent", nRemainingTracksInEvent,
1327 "Number of tracks in the event - Number of tracks( = charged FSPs) of particle.");
1328 REGISTER_VARIABLE(
"trackMatchType", trackMatchType, R
"DOC(
1330* -1 particle has no ECL cluster
1331* 0 particle has no associated track
1332* 1 there is a matched track called connected - region(CR) track match
1336 Use better variables like `trackNECLClusters`, `clusterTrackMatch`, and `nECLClusterTrackMatches`.)DOC");
1338 REGISTER_VARIABLE("decayTypeRecoil", recoilMCDecayType,
1339 "type of the particle decay(no related mcparticle = -1, hadronic = 0, direct leptonic = 1, direct semileptonic = 2,"
1340 "lower level leptonic = 3.");
1342 REGISTER_VARIABLE(
"printParticle", printParticle,
1343 "For debugging, print Particle and daughter PDG codes, plus MC match. Returns 0.");
1344 REGISTER_VARIABLE(
"mcMomTransfer2", particleMCMomentumTransfer2,
1345 "Return the true momentum transfer to lepton pair in a B(semi -) leptonic B meson decay.\n\n",
"GeV/c");
1346 REGISTER_VARIABLE(
"False", False,
1347 "returns always 0, used for testing and debugging.");
1348 REGISTER_VARIABLE(
"True", True,
1349 "returns always 1, used for testing and debugging.");
1350 REGISTER_VARIABLE(
"infinity", infinity,
1351 "returns std::numeric_limits<double>::infinity()");
1352 REGISTER_VARIABLE(
"random", random,
1353 "return a random number between 0 and 1 for each candidate. Can be used, e.g. for picking a random"
1354 "candidate in the best candidate selection.");
1355 REGISTER_VARIABLE(
"eventRandom", eventRandom,
1356 "[Eventbased] Returns a random number between 0 and 1 for this event. Can be used, e.g. for applying an event prescale.");
1357 REGISTER_METAVARIABLE(
"minET2ETDist(detName, detLayer, referenceListName, useHighestProbMassForExt=1)", particleDistToClosestExtTrk,
1358 R
"DOC(Returns the distance :math:`d_{\mathrm{i}}` in [cm] between the particle and the nearest particle in the reference list at the given detector :math:`i`-th layer surface.
1359The definition is based on the track helices extrapolation.
1361* The first argument is the name of the detector to consider.
1362* The second argument is the detector layer on whose surface we search for the nearest neighbour.
1363* The third argument is the reference particle list name used to search for the nearest neighbour.
1364* The fourth (optional) argument is an integer ("boolean") flag: if 1 (the default, if nothing is set), it is assumed the extrapolation was done with the most probable mass hypothesis for the track fit;
1365 if 0, it is assumed the mass hypothesis matching the particle lists' PDG was used.
1368 This variable requires to run the ``TrackIsoCalculator`` module first.
1369 Note that the choice of input parameters of this metafunction must correspond to the settings used to configure the module!
1371 Manager::VariableDataType::c_double);
1373 REGISTER_METAVARIABLE("minET2ETDistVar(detName, detLayer, referenceListName, variableName)", particleDistToClosestExtTrkVar,
1374 R
"DOC(Returns the value of the variable for the nearest neighbour to this particle as taken from the reference list at the given detector :math:`i`-th layer surface
1375, according to the distance definition of `minET2ETDist`.
1377* The first argument is the name of the detector to consider.
1378* The second argument is the detector layer on whose surface we search for the nearest neighbour.
1379* The third argument is the reference particle list name used to search for the nearest neighbour.
1380* The fourth argument is a variable name, e.g. `nCDCHits`.
1383 This variable requires to run the ``TrackIsoCalculator`` module first.
1384 Note that the choice of input parameters of this metafunction must correspond to the settings used to configure the module!
1386 Manager::VariableDataType::c_double);
1388 REGISTER_METAVARIABLE("minET2ETIsoScore(referenceListName, useHighestProbMassForExt, detectorList)", particleExtTrkIsoScoreVar,
1389 R
"DOC(Returns a particle's isolation score :math:`s` defined as:
1396 s &= \sum_{\mathrm{det}} 1 - \left(-w_{\mathrm{det}} \cdot \frac{\sum_{i}^{N_{\mathrm{det}}^{\mathrm{layers}}} H(i)}{N_{\mathrm{det}}^{\mathrm{layers}}}\right), \\
1400 0 & d_{\mathrm{i}} > D_{\mathrm{det}}^{\mathrm{thresh}} \\
1401 1 & d_{\mathrm{i}} <= D_{\mathrm{det}}^{\mathrm{thresh}}, \\
1406where :math:`d_{\mathrm{i}}` is the distance to the closest neighbour at the :math:`i`-th layer of the given detector (c.f., `minET2ETDist`), :math:`N_{\mathrm{det}}^{\mathrm{layers}}` is the
1407number of layers of the detector, :math:`D_{\mathrm{det}}^{\mathrm{thresh}}` is a threshold length related to the detector's granularity defined in the ``TrackIsoCalculator`` module,
1408and :math:`w_{\mathrm{det}}` are (negative) weights associated to the detector's impact on PID for this particle type, read from a CDB payload.
1410The score is normalised in [0, 1], where values closer to 1 indicates a well-isolated particle.
1412* The first argument is the reference particle list name used to search for the nearest neighbour.
1413* The second argument is an integer ("boolean") flag: if 1, it is assumed the extrapolation was done with the most probable mass hypothesis for the track fit;
1414 if 0, it is assumed the mass hypothesis matching the particle lists' PDG was used.
1415* The remaining arguments are a comma-separated list of detector names, which must correspond to the one given to the `TrackIsoCalculator` module.
1418 The PID detector weights :math:`w_{\mathrm{det}}` are non-trivial only if ``excludePIDDetWeights=false`` in the ``TrackIsoCalculator`` module configuration.
1419 Otherwise :math:`w_{\mathrm{det}} = -1`.
1422 This variable requires to run the `TrackIsoCalculator` module first.
1423 Note that the choice of input parameters of this metafunction must correspond to the settings used to configure the module!
1426 Manager::VariableDataType::c_double);
1428 REGISTER_METAVARIABLE("minET2ETIsoScoreAsWeightedAvg(referenceListName, useHighestProbMassForExt, detectorList)", particleExtTrkIsoScoreVarAsWeightedAvg,
1429 R
"DOC(Returns a particle's isolation score :math:`s` based on the weighted average:
1433 s = 1 - \frac{\sum_{\mathrm{det}} \sum_{i}^{N_{\mathrm{det}}^{\mathrm{layers}}} w_{\mathrm{det}} \cdot \frac{D_{\mathrm{det}}^{\mathrm{thresh}}}{d_{\mathrm{i}}} }{ \sum_{\mathrm{det}} w_{\mathrm{det}}},
1435where :math:`d_{\mathrm{i}}` is the distance to the closest neighbour at the :math:`i`-th layer of the given detector (c.f., `minET2ETDist`), :math:`N_{\mathrm{det}}^{\mathrm{layers}}` is the
1436number of layers of the detector, :math:`D_{\mathrm{det}}^{\mathrm{thresh}}` is a threshold length related to the detector's granularity defined in the ``TrackIsoCalculator`` module,
1437and :math:`w_{\mathrm{det}}` are (negative) weights associated to the detector's impact on PID for this particle type, read from a CDB payload.
1439The score is normalised in [0, 1], where values closer to 1 indicates a well-isolated particle.
1441* The first argument is the reference particle list name used to search for the nearest neighbour.
1442* The second argument is an integer ("boolean") flag: if 1, it is assumed the extrapolation was done with the most probable mass hypothesis for the track fit;
1443 if 0, it is assumed the mass hypothesis matching the particle lists' PDG was used.
1444* The remaining arguments are a comma-separated list of detector names, which must correspond to the one given to the `TrackIsoCalculator` module.
1447 The PID detector weights :math:`w_{\mathrm{det}}` are non-trivial only if ``excludePIDDetWeights=false`` in the ``TrackIsoCalculator`` module configuration.
1448 Otherwise :math:`\lvert w_{\mathrm{det}} \rvert = 1`.
1451 This variable requires to run the `TrackIsoCalculator` module first.
1452 Note that the choice of input parameters of this metafunction must correspond to the settings used to configure the module!
1455 Manager::VariableDataType::c_double);
DataType Perp() const
The transverse component (R in cylindrical coordinate system).
static BFieldManager & getInstance()
Return the instance of the magnetic field manager.
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
int getPDGCode() const
PDG code.
static const ChargedStable muon
muon particle
static const ChargedStable pion
charged pion particle
static const double pionMass
charged pion mass
static const double speedOfLight
[cm/ns]
static const ChargedStable proton
proton particle
static const double protonMass
proton mass
static const double doubleNaN
quiet_NaN
static const ChargedStable kaon
charged kaon particle
static const ParticleType photon
photon particle
static const ChargedStable electron
electron particle
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
@ c_PrimaryParticle
bit 0: Particle is primary particle.
const MCParticle * getMCParticle() const
Returns the pointer to the MCParticle object that was used to create this Particle (ParticleType == c...
static const ReferenceFrame & GetCurrent()
Get current rest frame.
std::function< VarVariant(const Particle *)> FunctionPtr
functions stored take a const Particle* and return VarVariant.
const Var * getVariable(std::string name)
Get the variable belonging to the given key.
static Manager & Instance()
get singleton instance.
Class to store variables with their name which were sent to the logging service.
#define MAKE_DEPRECATED(name, make_fatal, version, description)
Registers a variable as deprecated.
B2Vector3< double > B2Vector3D
typedef for common usage with double
Abstract base class for different kinds of events.
static std::string explainFlags(unsigned int flags)
Return string with all human-readable flags, e.g.
static int getMCErrors(const Belle2::Particle *particle, const Belle2::MCParticle *mcParticle=nullptr)
Returns quality indicator of the match as a bit pattern where the individual bits indicate the the ty...