10 #include <analysis/variables/ECLVariables.h> 
   13 #include <framework/logging/Logger.h> 
   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> 
   24 #include <mdst/dataobjects/KlId.h> 
   25 #include <mdst/dataobjects/ECLCluster.h> 
   26 #include <mdst/dataobjects/Track.h> 
   27 #include <mdst/dataobjects/EventLevelClusteringInfo.h> 
   29 #include <Math/Vector4D.h> 
   42     double beamBackgroundSuppression(
const Particle* particle)
 
   44       if (particle->hasExtraInfo(
"beamBackgroundSuppression")) {
 
   45         return particle->getExtraInfo(
"beamBackgroundSuppression");
 
   47         B2WARNING(
"The extraInfo beamBackgroundSuppression is not registered! \n" 
   48                   "This variable is only available for photons, and you either have to run the function getBeamBackgroundProbability or turn the argument loadPhotonBeamBackgroundMVA to True when using fillParticleList.");
 
   49         return std::numeric_limits<float>::quiet_NaN();
 
   53     double fakePhotonSuppression(
const Particle* particle)
 
   55       if (particle->hasExtraInfo(
"fakePhotonSuppression")) {
 
   56         return particle->getExtraInfo(
"fakePhotonSuppression");
 
   58         B2WARNING(
"The extraInfo fakePhotonSuppression is not registered! \n" 
   59                   "This variable is only available for photons, and you either have to run the function getFakePhotonProbability or turn the argument loadFakePhotonMVA to True when using fillParticleList.");
 
   60         return std::numeric_limits<float>::quiet_NaN();
 
   64     double hadronicSplitOffSuppression(
const Particle* particle)
 
   66       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.");
 
   67       return fakePhotonSuppression(particle);
 
   70     double eclClusterKlId(
const Particle* particle)
 
   72       const ECLCluster* cluster = particle->getECLCluster();
 
   74         return std::numeric_limits<double>::quiet_NaN();
 
   76       const KlId* klid = cluster->getRelatedTo<KlId>();
 
   78         return std::numeric_limits<double>::quiet_NaN();
 
   80       return klid->getKlId();
 
   84     double eclPulseShapeDiscriminationMVA(
const Particle* particle)
 
   86       const ECLCluster* cluster = particle->getECLCluster();
 
   88         if (eclClusterHasPulseShapeDiscrimination(particle)) {
 
   89           return cluster->getPulseShapeDiscriminationMVA();
 
   91           return std::numeric_limits<float>::quiet_NaN();
 
   94       return std::numeric_limits<float>::quiet_NaN();
 
   97     double eclClusterNumberOfHadronDigits(
const Particle* particle)
 
  100       const ECLCluster* cluster = particle->getECLCluster();
 
  102         if (eclClusterHasPulseShapeDiscrimination(particle)) {
 
  103           return cluster->getNumberOfHadronDigits();
 
  105           return std::numeric_limits<float>::quiet_NaN();
 
  107       return std::numeric_limits<float>::quiet_NaN();
 
  110     double eclClusterDetectionRegion(
const Particle* particle)
 
  113       const ECLCluster* cluster = particle->getECLCluster();
 
  115         return cluster->getDetectorRegion();
 
  117       return std::numeric_limits<float>::quiet_NaN();
 
  120     double eclClusterIsolation(
const Particle* particle)
 
  123       const ECLCluster* cluster = particle->getECLCluster();
 
  125         auto minDist = cluster->getMinTrkDistance();
 
  129       return std::numeric_limits<float>::quiet_NaN();
 
  132     double eclClusterIsolationID(
const Particle* particle)
 
  135       const ECLCluster* cluster = particle->getECLCluster();
 
  137         return cluster->getMinTrkDistanceID();
 
  139       return std::numeric_limits<float>::quiet_NaN();
 
  144       if (arguments.size() > 2 or arguments.size() == 0)
 
  145         B2FATAL(
"Wrong number of arguments (2 required) for meta variable minC2TDistVar");
 
  146       std::string listName = 
"pi-:all";
 
  147       std::string variableName = arguments[0];
 
  148       if (arguments.size() == 2)
 
  149         listName = arguments[1];
 
  152       auto func = [listName, variableName](
const Particle * particle) -> 
double {
 
  153         StoreObjPtr<ParticleList> particleList(listName);
 
  154         if (!(particleList.isValid()))
 
  156           B2FATAL(
"Invalid Listname " << listName << 
" given to minC2TDistVar!");
 
  159         const ECLCluster* cluster = particle->getECLCluster();
 
  161           return std::numeric_limits<float>::quiet_NaN();
 
  162         auto trackID = cluster->getMinTrkDistanceID();
 
  163         double result = std::numeric_limits<float>::quiet_NaN();
 
  165         for (
unsigned int i = 0; i < particleList->getListSize(); i++)
 
  167           const Particle* listParticle = particleList->getParticle(i);
 
  168           if (listParticle and listParticle->getTrack() and listParticle->getTrack()->getArrayIndex() == trackID) {
 
  169             result = std::get<double>(var->function(listParticle));
 
  178     double eclClusterConnectedRegionID(
const Particle* particle)
 
  181       const ECLCluster* cluster = particle->getECLCluster();
 
  183         return  cluster->getConnectedRegionId();
 
  185       return std::numeric_limits<float>::quiet_NaN();
 
  188     double eclClusterDeltaL(
const Particle* particle)
 
  191       const ECLCluster* cluster = particle->getECLCluster();
 
  193         return cluster->getDeltaL();
 
  195       return std::numeric_limits<float>::quiet_NaN();
 
  198     double eclClusterErrorE(
const Particle* particle)
 
  201       const ECLCluster* cluster = particle->getECLCluster();
 
  203         return cluster->getUncertaintyEnergy();
 
  205       return std::numeric_limits<float>::quiet_NaN();
 
  208     double eclClusterUncorrectedE(
const Particle* particle)
 
  211       const ECLCluster* cluster = particle->getECLCluster();
 
  213         return cluster->getEnergyRaw();
 
  215       return std::numeric_limits<float>::quiet_NaN();
 
  218     double eclClusterE(
const Particle* particle)
 
  222       const ECLCluster* cluster = particle->getECLCluster();
 
  225         ROOT::Math::PxPyPzEVector p4Cluster = clutls.GetCluster4MomentumFromCluster(cluster, particle->getECLClusterEHypothesisBit());
 
  227         return frame.getMomentum(p4Cluster).E();
 
  229       return std::numeric_limits<float>::quiet_NaN();
 
  232     double eclClusterHighestE(
const Particle* particle)
 
  235       const ECLCluster* cluster = particle->getECLCluster();
 
  237         return cluster->getEnergyHighestCrystal();
 
  239       return std::numeric_limits<float>::quiet_NaN();
 
  242     double eclClusterCellId(
const Particle* particle)
 
  245       const ECLCluster* cluster = particle->getECLCluster();
 
  247         return cluster->getMaxECellId();
 
  249       return std::numeric_limits<float>::quiet_NaN();
 
  253     const std::array<int, 69> lastCellIDperThetaID{48,   96,  160,  224,  288,  384,  480,  576,  672,  768,  864,
 
  254             1008, 1152, 1296, 1440, 1584, 1728, 1872, 2016, 2160, 2304, 2448,
 
  255             2592, 2736, 2880, 3024, 3168, 3312, 3456, 3600, 3744, 3888, 4032,
 
  256             4176, 4320, 4464, 4608, 4752, 4896, 5040, 5184, 5328, 5472, 5616,
 
  257             5760, 5904, 6048, 6192, 6336, 6480, 6624, 6768, 6912, 7056, 7200,
 
  258             7344, 7488, 7632, 7776, 7920, 8064, 8160, 8256, 8352, 8448, 8544,
 
  261     double eclClusterThetaId(
const Particle* particle)
 
  264       const ECLCluster* cluster = particle->getECLCluster();
 
  266         int cellID = cluster->getMaxECellId();
 
  267         return std::distance(lastCellIDperThetaID.begin(), std::lower_bound(lastCellIDperThetaID.begin(), lastCellIDperThetaID.end(),
 
  270       return std::numeric_limits<float>::quiet_NaN();
 
  273     double eclClusterPhiId(
const Particle* particle)
 
  276       const ECLCluster* cluster = particle->getECLCluster();
 
  278         int cellID = cluster->getMaxECellId();
 
  282           int closestinlist = lastCellIDperThetaID[std::distance(lastCellIDperThetaID.begin(), std::lower_bound(lastCellIDperThetaID.begin(),
 
  283                                                                                                              lastCellIDperThetaID.end(), cellID)) - 1];
 
  284           return cellID - closestinlist - 1;
 
  287       return std::numeric_limits<float>::quiet_NaN();
 
  290     double eclClusterTiming(
const Particle* particle)
 
  293       const ECLCluster* cluster = particle->getECLCluster();
 
  295         return cluster->getTime();
 
  297       return std::numeric_limits<float>::quiet_NaN();
 
  300     double eclClusterHasFailedTiming(
const Particle* particle)
 
  302       const ECLCluster* cluster = particle->getECLCluster();
 
  304         return cluster->hasFailedFitTime();
 
  306       return std::numeric_limits<float>::quiet_NaN();
 
  309     double eclClusterErrorTiming(
const Particle* particle)
 
  312       const ECLCluster* cluster = particle->getECLCluster();
 
  314         return cluster->getDeltaTime99();
 
  316       return std::numeric_limits<float>::quiet_NaN();
 
  319     double eclClusterHasFailedErrorTiming(
const Particle* particle)
 
  321       const ECLCluster* cluster = particle->getECLCluster();
 
  323         return cluster->hasFailedTimeResolution();
 
  325       return std::numeric_limits<float>::quiet_NaN();
 
  328     double eclClusterTheta(
const Particle* particle)
 
  332       const ECLCluster* cluster = particle->getECLCluster();
 
  335         ROOT::Math::PxPyPzEVector p4Cluster = clutls.Get4MomentumFromCluster(cluster, particle->getECLClusterEHypothesisBit());
 
  337         return frame.getMomentum(p4Cluster).Theta();
 
  339       return std::numeric_limits<float>::quiet_NaN();
 
  342     double eclClusterErrorTheta(
const Particle* particle)
 
  345       const ECLCluster* cluster  = particle->getECLCluster();
 
  347         return cluster->getUncertaintyTheta();
 
  349       return std::numeric_limits<float>::quiet_NaN();
 
  352     double eclClusterErrorPhi(
const Particle* particle)
 
  355       const ECLCluster* cluster = particle->getECLCluster();
 
  357         return cluster->getUncertaintyPhi();
 
  359       return std::numeric_limits<float>::quiet_NaN();
 
  362     double eclClusterPhi(
const Particle* particle)
 
  366       const ECLCluster* cluster = particle->getECLCluster();
 
  369         ROOT::Math::PxPyPzEVector p4Cluster = clutls.Get4MomentumFromCluster(cluster, particle->getECLClusterEHypothesisBit());
 
  371         return frame.getMomentum(p4Cluster).Phi();
 
  373       return std::numeric_limits<float>::quiet_NaN();
 
  376     double eclClusterR(
const Particle* particle)
 
  379       const ECLCluster* cluster = particle->getECLCluster();
 
  381         return cluster->getR();
 
  383       return std::numeric_limits<float>::quiet_NaN();
 
  386     double eclClusterE1E9(
const Particle* particle)
 
  389       const ECLCluster* cluster = particle->getECLCluster();
 
  391         return cluster->getE1oE9();
 
  393       return std::numeric_limits<float>::quiet_NaN();
 
  396     double eclClusterE9E21(
const Particle* particle)
 
  399       const ECLCluster* cluster = particle->getECLCluster();
 
  401         return cluster->getE9oE21();
 
  403       return std::numeric_limits<float>::quiet_NaN();
 
  406     double eclClusterAbsZernikeMoment40(
const Particle* particle)
 
  409       const ECLCluster* cluster = particle->getECLCluster();
 
  411         return cluster->getAbsZernike40();
 
  413       return std::numeric_limits<float>::quiet_NaN();
 
  416     double eclClusterAbsZernikeMoment51(
const Particle* particle)
 
  419       const ECLCluster* cluster = particle->getECLCluster();
 
  421         return cluster->getAbsZernike51();
 
  423       return std::numeric_limits<float>::quiet_NaN();
 
  426     double eclClusterZernikeMVA(
const Particle* particle)
 
  429       const ECLCluster* cluster = particle->getECLCluster();
 
  431         return cluster->getZernikeMVA();
 
  433       return std::numeric_limits<float>::quiet_NaN();
 
  436     double eclClusterSecondMoment(
const Particle* particle)
 
  439       const ECLCluster* cluster = particle->getECLCluster();
 
  441         return cluster->getSecondMoment();
 
  443       return std::numeric_limits<float>::quiet_NaN();
 
  446     double eclClusterLAT(
const Particle* particle)
 
  449       const ECLCluster* cluster = particle->getECLCluster();
 
  451         return cluster->getLAT();
 
  453       return std::numeric_limits<float>::quiet_NaN();
 
  456     double eclClusterNHits(
const Particle* particle)
 
  459       const ECLCluster* cluster = particle->getECLCluster();
 
  461         return cluster->getNumberOfCrystals();
 
  463       return std::numeric_limits<float>::quiet_NaN();
 
  466     double eclClusterTrackMatched(
const Particle* particle)
 
  469       const ECLCluster* cluster = particle->getECLCluster();
 
  471         const Track* track = cluster->getRelatedFrom<Track>();
 
  478       return std::numeric_limits<float>::quiet_NaN();
 
  481     double nECLClusterTrackMatches(
const Particle* particle)
 
  484       const ECLCluster* cluster = particle->getECLCluster();
 
  486         return std::numeric_limits<double>::quiet_NaN();
 
  489       size_t out = cluster->getRelationsFrom<Track>().size();
 
  493     double eclClusterConnectedRegionId(
const Particle* particle)
 
  495       const ECLCluster* cluster = particle->getECLCluster();
 
  497         return cluster->getConnectedRegionId();
 
  499       return std::numeric_limits<float>::quiet_NaN();
 
  502     double eclClusterId(
const Particle* particle)
 
  504       const ECLCluster* cluster = particle->getECLCluster();
 
  506         return cluster->getClusterId();
 
  508       return std::numeric_limits<float>::quiet_NaN();
 
  511     double eclClusterHasNPhotonsHypothesis(
const Particle* particle)
 
  513       const ECLCluster* cluster = particle->getECLCluster();
 
  517       return std::numeric_limits<float>::quiet_NaN();
 
  520     double eclClusterHasNeutralHadronHypothesis(
const Particle* particle)
 
  522       const ECLCluster* cluster = particle->getECLCluster();
 
  526       return std::numeric_limits<float>::quiet_NaN();
 
  529     double eclClusterHasPulseShapeDiscrimination(
const Particle* particle)
 
  531       const ECLCluster* cluster = particle->getECLCluster();
 
  533         return cluster->hasPulseShapeDiscrimination();
 
  535       return std::numeric_limits<float>::quiet_NaN();
 
  538     double eclExtTheta(
const Particle* particle)
 
  541       const Track* track = particle->getTrack();
 
  544         auto* eclinfo = track->getRelatedTo<ECLEnergyCloseToTrack>();
 
  547           return eclinfo->getExtTheta();
 
  549           B2WARNING(
"Relation to ECLEnergyCloseToTrack not found, did you forget to run ECLTrackCalDigitMatchModule?");
 
  550           return std::numeric_limits<float>::quiet_NaN();
 
  554       return std::numeric_limits<float>::quiet_NaN();
 
  557     double eclExtPhi(
const Particle* particle)
 
  560       const Track* track = particle->getTrack();
 
  563         auto* eclinfo = track->getRelatedTo<ECLEnergyCloseToTrack>();
 
  566           return eclinfo->getExtPhi();
 
  568           B2WARNING(
"Relation to ECLEnergyCloseToTrack not found, did you forget to run ECLTrackCalDigitMatchModule?");
 
  569           return std::numeric_limits<float>::quiet_NaN();
 
  573       return std::numeric_limits<float>::quiet_NaN();
 
  576     double eclExtPhiId(
const Particle* particle)
 
  578       const Track* track = particle->getTrack();
 
  581         auto* eclinfo = track->getRelatedTo<ECLEnergyCloseToTrack>();
 
  584           return eclinfo->getExtPhiId();
 
  586           B2WARNING(
"Relation to ECLEnergyCloseToTrack not found, did you forget to run ECLTrackCalDigitMatchModule?");
 
  587           return std::numeric_limits<float>::quiet_NaN();
 
  591       return std::numeric_limits<float>::quiet_NaN();
 
  594     double weightedAverageECLTime(
const Particle* particle)
 
  596       int nDaughters = particle->getNDaughters();
 
  597       if (nDaughters < 1) {
 
  598         B2WARNING(
"The provided particle has no daughters!");
 
  599         return std::numeric_limits<float>::quiet_NaN();
 
  602       double numer = 0, denom = 0;
 
  603       int numberOfClusterDaughters = 0;
 
  605       auto weightedECLTimeAverage = [&numer, &denom, &numberOfClusterDaughters](
const Particle * p) {
 
  606         const ECLCluster* cluster = p->getECLCluster();
 
  607         if (cluster and not cluster->hasFailedFitTime()) {
 
  608           numberOfClusterDaughters ++;
 
  610           double time = cluster->getTime();
 
  611           B2DEBUG(10, 
"time[" << numberOfClusterDaughters << 
"] = " << time);
 
  612           double deltatime = cluster->getDeltaTime99();
 
  613           B2DEBUG(10, 
"deltatime[" << numberOfClusterDaughters << 
"] = " << deltatime);
 
  614           numer += time / pow(deltatime, 2);
 
  615           B2DEBUG(11, 
"numer[" << numberOfClusterDaughters << 
"] = " << numer);
 
  616           denom += 1 / pow(deltatime, 2);
 
  617           B2DEBUG(11, 
"denom[" << numberOfClusterDaughters << 
"] = " << denom);
 
  622       particle->forEachDaughter(weightedECLTimeAverage, 
true, 
true);
 
  624       if (numberOfClusterDaughters < 1) {
 
  625         B2WARNING(
"There are no clusters or cluster matches amongst the daughters of the provided particle!");
 
  626         return std::numeric_limits<float>::quiet_NaN();
 
  630         B2WARNING(
"The denominator of the weighted mean is zero!");
 
  631         return std::numeric_limits<float>::quiet_NaN();
 
  633         B2DEBUG(10, 
"numer/denom = " << numer / denom);
 
  634         return numer / denom;
 
  638     double maxWeightedDistanceFromAverageECLTime(
const Particle* particle)
 
  640       int nDaughters = particle->getNDaughters();
 
  641       if (nDaughters < 1) {
 
  642         B2WARNING(
"The provided particle has no daughters!");
 
  643         return std::numeric_limits<float>::quiet_NaN();
 
  646       double maxTimeDiff = -DBL_MAX;
 
  647       int numberOfClusterDaughters = 0;
 
  649       double averageECLTime = weightedAverageECLTime(particle);
 
  651       auto maxTimeDifference = [&maxTimeDiff, &numberOfClusterDaughters, &averageECLTime](
const Particle * p) {
 
  653         const ECLCluster* cluster = p->getECLCluster();
 
  655           numberOfClusterDaughters ++;
 
  657           double time = cluster->getTime();
 
  658           B2DEBUG(10, 
"time[" << numberOfClusterDaughters << 
"] = " << time);
 
  659           double deltatime = cluster->getDeltaTime99();
 
  660           B2DEBUG(10, 
"deltatime[" << numberOfClusterDaughters << 
"] = " << deltatime);
 
  661           double maxTimeDiff_temp = fabs((time - averageECLTime) / deltatime);
 
  662           B2DEBUG(11, 
"maxTimeDiff_temp[" << numberOfClusterDaughters << 
"] = " << maxTimeDiff_temp);
 
  663           if (maxTimeDiff_temp > maxTimeDiff)
 
  664             maxTimeDiff = maxTimeDiff_temp;
 
  665           B2DEBUG(11, 
"maxTimeDiff[" << numberOfClusterDaughters << 
"] = " << maxTimeDiff);
 
  670       particle->forEachDaughter(maxTimeDifference, 
true, 
true);
 
  672       if (numberOfClusterDaughters < 1) {
 
  673         B2WARNING(
"There are no clusters or cluster matches amongst the daughters of the provided particle!");
 
  674         return std::numeric_limits<float>::quiet_NaN();
 
  677       if (maxTimeDiff < 0) {
 
  678         B2WARNING(
"The max time difference is negative!");
 
  679         return std::numeric_limits<float>::quiet_NaN();
 
  681         B2DEBUG(10, 
"maxTimeDiff = " << maxTimeDiff);
 
  686     double eclClusterMdstIndex(
const Particle* particle)
 
  688       const ECLCluster* cluster = particle->getECLCluster();
 
  690         return cluster->getArrayIndex();
 
  691       } 
else return std::numeric_limits<double>::quiet_NaN();
 
  693       return std::numeric_limits<double>::quiet_NaN();
 
  700     double nECLOutOfTimeCrystalsFWDEndcap(
const Particle*)
 
  702       StoreObjPtr<EventLevelClusteringInfo> elci;
 
  703       if (!elci) 
return std::numeric_limits<double>::quiet_NaN();
 
  704       return (
double)elci->getNECLCalDigitsOutOfTimeFWD();
 
  707     double nECLOutOfTimeCrystalsBarrel(
const Particle*)
 
  709       StoreObjPtr<EventLevelClusteringInfo> elci;
 
  710       if (!elci) 
return std::numeric_limits<double>::quiet_NaN();
 
  711       return (
double)elci->getNECLCalDigitsOutOfTimeBarrel();
 
  714     double nECLOutOfTimeCrystalsBWDEndcap(
const Particle*)
 
  716       StoreObjPtr<EventLevelClusteringInfo> elci;
 
  717       if (!elci) 
return std::numeric_limits<double>::quiet_NaN();
 
  718       return (
double)elci->getNECLCalDigitsOutOfTimeBWD();
 
  721     double nECLOutOfTimeCrystals(
const Particle*)
 
  723       StoreObjPtr<EventLevelClusteringInfo> elci;
 
  724       if (!elci) 
return std::numeric_limits<double>::quiet_NaN();
 
  725       return (
double)elci->getNECLCalDigitsOutOfTime();
 
  728     double nRejectedECLShowersFWDEndcap(
const Particle*)
 
  730       StoreObjPtr<EventLevelClusteringInfo> elci;
 
  731       if (!elci) 
return std::numeric_limits<double>::quiet_NaN();
 
  732       return (
double)elci->getNECLShowersRejectedFWD();
 
  735     double nRejectedECLShowersBarrel(
const Particle*)
 
  737       StoreObjPtr<EventLevelClusteringInfo> elci;
 
  738       if (!elci) 
return std::numeric_limits<double>::quiet_NaN();
 
  739       return (
double)elci->getNECLShowersRejectedBarrel();
 
  742     double nRejectedECLShowersBWDEndcap(
const Particle*)
 
  744       StoreObjPtr<EventLevelClusteringInfo> elci;
 
  745       if (!elci) 
return std::numeric_limits<double>::quiet_NaN();
 
  746       return (
double)elci->getNECLShowersRejectedBWD();
 
  749     double nRejectedECLShowers(
const Particle*)
 
  751       StoreObjPtr<EventLevelClusteringInfo> elci;
 
  752       if (!elci) 
return std::numeric_limits<double>::quiet_NaN();
 
  753       return (
double) elci->getNECLShowersRejected();
 
  756     double eclClusterEoP(
const Particle* part)
 
  758       double E = eclClusterE(part);
 
  759       if (part->hasExtraInfo(
"bremsCorrectedPhotonEnergy")) {
 
  760         E += part->getExtraInfo(
"bremsCorrectedPhotonEnergy");
 
  762       const double p =  part->getMomentumMagnitude();
 
  763       if (0 == p) { 
return std::numeric_limits<float>::quiet_NaN();}
 
  767     double eclClusterOnlyInvariantMass(
const Particle* part)
 
  769       int nDaughters = part->getNDaughters();
 
  770       ROOT::Math::PxPyPzEVector sum;
 
  772       if (nDaughters < 1) {
 
  773         return part->getMass();
 
  775         int nClusterDaughters = 0;
 
  776         std::stack<const Particle*> stacked;
 
  778         while (!stacked.empty()) {
 
  779           const Particle* current = stacked.top();
 
  782           const ECLCluster* cluster = current->getECLCluster();
 
  785             nClusterDaughters ++;
 
  787             ROOT::Math::PxPyPzEVector p4Cluster = clutls.Get4MomentumFromCluster(cluster, clusterBit);
 
  790             const std::vector<Particle*> daughters = current->getDaughters();
 
  791             nDaughters = current->getNDaughters();
 
  792             for (
int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
 
  793               stacked.push(daughters[iDaughter]);
 
  798         if (nClusterDaughters < 1) {
 
  799           B2WARNING(
"There are no clusters amongst the daughters of the provided particle!");
 
  800           return std::numeric_limits<double>::quiet_NaN();
 
  802         B2DEBUG(10, 
"Number of daughters with cluster associated = " << nClusterDaughters);
 
  809       std::string cutString = 
"";
 
  810       if (arguments.size() > 0) {
 
  811         cutString = arguments[0];
 
  815       std::string photonlistname = 
"gamma:all";
 
  816       if (arguments.size() > 1) {
 
  817         photonlistname = arguments[1];
 
  820       std::string tracklistname = 
"e-:all";
 
  821       if (arguments.size() > 2) {
 
  822         tracklistname = arguments[2];
 
  825       auto func = [cut, photonlistname, tracklistname](
const Particle * particle) -> 
double {
 
  829           B2WARNING(
"The variable photonHasOverlap is supposed to be calculated for photons. Returning NaN.");
 
  830           return std::numeric_limits<double>::quiet_NaN();
 
  833         StoreObjPtr<ParticleList> photonlist(photonlistname);
 
  834         if (!(photonlist.isValid()))
 
  836           B2WARNING(
"The provided particle list " << photonlistname << 
" does not exist." 
  837                     " Therefore, the variable photonHasOverlap can not be calculated. Returning NaN.");
 
  838           return std::numeric_limits<double>::quiet_NaN();
 
  842           B2WARNING(
"The list " << photonlistname << 
" does not contain photons." 
  843                     " Therefore, the variable photonHasOverlap can not be calculated reliably. Returning NaN.");
 
  844           return std::numeric_limits<double>::quiet_NaN();
 
  847         StoreObjPtr<ParticleList> tracklist(tracklistname);
 
  848         if (!(tracklist.isValid()))
 
  850           B2WARNING(
"The provided particle list " << tracklistname << 
" does not exist." 
  851                     " Therefore, the variable photonHasOverlap can not be calculated. Returning NaN.");
 
  852           return std::numeric_limits<double>::quiet_NaN();
 
  856           B2WARNING(
"The list " << tracklistname << 
" does not contain charged final state particles." 
  857                     " Therefore, the variable photonHasOverlap can not be calculated reliably. Returning NaN.");
 
  858           return std::numeric_limits<double>::quiet_NaN();
 
  861         double connectedRegionID = eclClusterConnectedRegionID(particle);
 
  862         unsigned mdstArrayIndex = particle->getMdstArrayIndex();
 
  864         for (
unsigned int i = 0; i < photonlist->getListSize(); i++)
 
  866           const Particle* part = photonlist->getParticle(i);
 
  869           if (part->getMdstArrayIndex() == mdstArrayIndex) {
 
  874           if (!cut->check(part)) {
 
  878           if (connectedRegionID == eclClusterConnectedRegionID(part)) {
 
  883         for (
unsigned int i = 0; i < tracklist->getListSize(); i++)
 
  885           const Particle* part = tracklist->getParticle(i);
 
  888           if (!cut->check(part)) {
 
  892           if (connectedRegionID == eclClusterConnectedRegionID(part)) {
 
  902     VARIABLE_GROUP(
"ECL Cluster related");
 
  903     REGISTER_VARIABLE(
"clusterEoP", eclClusterEoP, R
"DOC( 
  904 Returns ratio of uncorrelated energy E over momentum p, a convenience 
  905 alias for (clusterE / p). 
  907     REGISTER_VARIABLE("clusterReg", eclClusterDetectionRegion, R
"DOC( 
  908 Returns an integer code for the ECL region of a cluster. 
  910 - 1: forward, 2: barrel, 3: backward, 
  911 - 11: between FWD and barrel, 13: between BWD and barrel, 
  914     REGISTER_VARIABLE("clusterDeltaLTemp", eclClusterDeltaL, R
"DOC( 
  915 | Returns DeltaL for the shower shape. 
  916 | A cluster comprises the energy depositions of several crystals. All these crystals have slightly 
  917   different orientations in space. A shower direction can be constructed by calculating the weighted 
  918   average of these orientations using the corresponding energy depositions as weights. The intersection 
  919   (more precisely the point of closest approach) of the vector with this direction originating from the 
  920   cluster center and an extrapolated track can be used as reference for the calculation of the shower 
  921   depth. It is defined as the distance between this intersection and the cluster center. 
  924     This distance is calculated on the reconstructed level and is temporarily 
  925     included to the ECL cluster MDST data format for studying purposes. If it is found 
  926     that it is not crucial for physics analysis then this variable will be removed 
  928     Therefore, keep in mind that this variable might be removed in the future! 
  931     | Please read `this <importantNoteECL>` first. 
  932     | Lower limit: :math:`-250.0` 
  933     | Upper limit: :math:`250.0` 
  934     | Precision: :math:`10` bit 
  939     REGISTER_VARIABLE("minC2TDist", eclClusterIsolation, R"DOC( 
  940 Returns the distance between the ECL cluster and its nearest track.  
  942 For all tracks in the event, the distance between each of their extrapolated hits in the ECL and the ECL shower  
  943 position is calculated, and the overall smallest distance is returned. The track array index of the track that is  
  944 closest to the ECL cluster can be retrieved using `minC2TDistID`.  
  946 If the calculated distance is greater than :math:`250.0`, the returned distance will be capped at :math:`250.0`.  
  947 If there are no extrapolated hits found in the ECL for the event, NaN will be returned.  
  950     This distance is calculated on the reconstructed level. 
  953     | Please read `this <importantNoteECL>` first. 
  954     | Lower limit: :math:`0.0` 
  955     | Upper limit: :math:`250.0` 
  956     | Precision: :math:`10` bit 
  960     REGISTER_VARIABLE("minC2TDistID", eclClusterIsolationID, R"DOC( 
  961 Returns the track array index of the nearest track to the ECL cluster. The nearest track is calculcated  
  962 using the `minC2TDist` variable.  
  964     REGISTER_METAVARIABLE("minC2TDistVar(variable,particleList=pi-:all)", eclClusterIsolationVar, R
"DOC( 
  965 Returns the variable value of the nearest track to the given ECL cluster as calculated by `minC2TDist`. The  
  966 first argument is the variable name, e.g. `nCDCHits`, while the second (optional) argument is the particle list name which  
  967 will be used to pick up the nearest track in the calculation of `minC2TDist`. The default particle list used  
  969 )DOC", Manager::VariableDataType::c_double); 
  970     REGISTER_VARIABLE("clusterE", eclClusterE, R
"DOC( 
  971 Returns ECL cluster's energy corrected for leakage and background. 
  973 The raw photon energy is given by the weighted sum of all ECL crystal energies within the ECL cluster. 
  974 The weights per crystals are :math:`\leq 1` after cluster energy splitting in the case of overlapping 
  975 clusters. The number of crystals that are included in the sum depends on a initial energy estimation 
  976 and local beam background levels at the highest energy crystal position. It is optimized to minimize 
  977 the core width (resolution) of true photons. Photon energy distributions always show a low energy tail 
  978 due to unavoidable longitudinal and transverse leakage that can be further modified by the clustering 
  979 algorithm and beam backgrounds.The peak position of the photon energy distributions are corrected to 
  980 match the true photon energy in MC: 
  982 - Leakage correction: Using large MC samples of mono-energetic single photons, a correction factor 
  983   :math:`f` as function of reconstructed detector position, reconstructed photon energy and beam backgrounds 
  984   is determined via :math:`f = \frac{\text{peak_reconstructed}}{\text{energy_true}}`. 
  986 - Cluster energy calibration (data only): To reach the target precision of :math:`< 1.8\%` energy 
  987   resolution for high energetic photons, the remaining difference between MC and data must be calibrated 
  988   using kinematically fit muon pairs. This calibration is only applied to data and not to MC and will 
  989   take time to develop. 
  991 - Energy Bias Correction module, sub-percent correction, is NOT applied on clusterE, but on photon energy 
  992   and momentum. Only applied to data. 
  994 It is important to note that after perfect leakage correction and cluster energy calibration, 
  995 the :math:`\pi^{0}` mass peak will be shifted slightly to smaller values than the PDG average 
  996 due to the low energy tails of photons. The :math:`\pi^{0}` mass peak must not be corrected 
  997 to the PDG value by adjusting the reconstructed photon energies. Selection criteria based on 
  998 the mass for :math:`\pi^{0}` candidates must be based on the biased value. Most analysis 
  999 will used mass constrained :math:`\pi^{0}` s anyhow. 
 1002     We only store clusters with :math:`E > 20\,` MeV. 
 1005     | Please read `this <importantNoteECL>` first. 
 1006     | Lower limit: :math:`-5` (:math:`e^{-5} = 0.00674\,` GeV) 
 1007     | Upper limit: :math:`3.0` (:math:`e^3 = 20.08553\,` GeV) 
 1008     | Precision: :math:`18` bit 
 1009     | This value can be changed to a different reference frame with :b2:var:`useCMSFrame`. 
 1013     REGISTER_VARIABLE("clusterErrorE", eclClusterErrorE, R"DOC( 
 1014 Returns ECL cluster's uncertainty on energy 
 1015 (from background level and energy dependent tabulation). 
 1018     REGISTER_VARIABLE("clusterErrorPhi", eclClusterErrorPhi, R"DOC( 
 1019 Returns ECL cluster's uncertainty on :math:`\phi` 
 1020 (from background level and energy dependent tabulation). 
 1023     REGISTER_VARIABLE("clusterErrorTheta", eclClusterErrorTheta, R"DOC( 
 1024 Returns ECL cluster's uncertainty on :math:`\theta` 
 1025 (from background level and energy dependent tabulation). 
 1029     REGISTER_VARIABLE("clusterR", eclClusterR, R"DOC( 
 1030 Returns ECL cluster's centroid distance from :math:`(0,0,0)`. 
 1033     REGISTER_VARIABLE("clusterPhi", eclClusterPhi, R"DOC( 
 1034 Returns ECL cluster's azimuthal angle :math:`\phi` 
 1035 (this is not generally equal to a photon azimuthal angle). 
 1037 | The direction of a cluster is given by the connecting line of :math:`\,(0,0,0)\,` and 
 1038   cluster centroid position in the ECL. 
 1039 | The cluster centroid position is calculated using up to 21 crystals (5x5 excluding corners) 
 1040   after cluster energy splitting in the case of overlapping clusters. 
 1041 | The centroid position is the logarithmically weighted average of all crystals evaluated at 
 1042   the crystal centers. Cluster centroids are generally biased towards the centers of the 
 1043   highest energetic crystal. This effect is larger for low energetic photons. 
 1044 | Beam backgrounds slightly decrease the position resolution, mainly for low energetic photons. 
 1047     Radius of a cluster is almost constant in the barrel and should not be used directly in any selection. 
 1049 Unlike for charged tracks, the uncertainty (covariance) of the photon directions is not determined 
 1050 based on individual cluster properties but taken from on MC-based parametrizations of the resolution 
 1051 as function of true photon energy, true photon direction and beam background level. 
 1054     Users must use the actual particle direction (done automatically in the modularAnalysis using the average 
 1055     IP position (can be changed if needed)) and not the ECL Cluster direction (position in the ECL measured 
 1056     from :math:`(0,0,0)`) for particle kinematics. 
 1059     | Please read `this <importantNoteECL>` first. 
 1060     | Lower limit: :math:`-\pi` 
 1061     | Upper limit: :math:`\pi` 
 1062     | Precision: :math:`16` bit 
 1066     REGISTER_VARIABLE("clusterConnectedRegionID", eclClusterConnectedRegionID, R"DOC( 
 1067 Returns ECL cluster's connected region ID. 
 1069     REGISTER_VARIABLE("clusterTheta", eclClusterTheta, R
"DOC( 
 1070 Returns ECL cluster's polar angle :math:`\theta` 
 1071 (this is not generally equal to a photon polar angle). 
 1073 | The direction of a cluster is given by the connecting line of :math:`\,(0,0,0)\,` and 
 1074   cluster centroid position in the ECL. 
 1075 | The cluster centroid position is calculated using up to 21 crystals (5x5 excluding corners) 
 1076   after cluster energy splitting in the case of overlapping clusters. 
 1077 | The centroid position is the logarithmically weighted average of all crystals evaluated at 
 1078   the crystal centers. Cluster centroids are generally biased towards the centers of the 
 1079   highest energetic crystal. This effect is larger for low energetic photons. 
 1080 | Beam backgrounds slightly decrease the position resolution, mainly for low energetic photons. 
 1083     Radius of a cluster is almost constant in the barrel and should not be used directly in any selection. 
 1085 Unlike for charged tracks, the uncertainty (covariance) of the photon directions is not determined 
 1086 based on individual cluster properties but taken from on MC-based parametrizations of the resolution 
 1087 as function of true photon energy, true photon direction and beam background level. 
 1090     Users must use the actual particle direction (done automatically in the modularAnalysis using the average 
 1091     IP position (can be changed if needed)) and not the ECL Cluster direction (position in the ECL measured 
 1092     from :math:`(0,0,0)`) for particle kinematics. 
 1095     | Please read `this <importantNoteECL>` first. 
 1096     | Lower limit: :math:`0.0` 
 1097     | Upper limit: :math:`\pi` 
 1098     | Precision: :math:`16` bit 
 1102     REGISTER_VARIABLE("clusterTiming", eclClusterTiming, R"DOC( 
 1104 Returns the time of the ECL cluster. It is calculated as the Photon timing minus the Event t0. 
 1105 Photon timing is given by the fitted time of the recorded waveform of the highest energy crystal in the 
 1106 cluster. After all calibrations and corrections (including Time-Of-Flight), photons from the interaction 
 1107 point (IP) should have a Photon timing that corresponds to the Event t0, :math:`t_{0}`.  The Event t0 is the 
 1108 time of the event and may be measured by a different sub-detector (see Event t0 documentation). For an ECL 
 1109 cluster produced at the interaction point in time with the event, the cluster time should be consistent with zero 
 1110 within the uncertainties. Special values are returned if the fit for the Photon timing fails (see 
 1111 documentation for `clusterHasFailedTiming`). (For MC, the calibrations and corrections are not fully simulated). 
 1114     | Please read `this <importantNoteECL>` first. 
 1115     | Lower limit: :math:`-1000.0` 
 1116     | Upper limit: :math:`1000.0` 
 1117     | Precision: :math:`12` bit 
 1121 Returns the trigger cell (TC) time of the ECL cluster (photon). 
 1122 This information is available only in Belle data since experiment 31, and not available in Belle MC. 
 1123 Clusters produced at the interaction point in time with the event, have TC time in the range of 9000-11000 
 1124 Calculated based on the Appendix of Belle note 831. 
 1127     | In case this variable is obtained from Belle data that is stored in Belle II mdst/udst format, it will be truncated to: 
 1128     | Lower limit: :math:`-1000.0` 
 1129     | Upper limit: :math:`1000.0` 
 1130     | Precision: :math:`12` bit 
 1134     REGISTER_VARIABLE("clusterHasFailedTiming", eclClusterHasFailedTiming, R"DOC( 
 1135 Status bit for if the ECL cluster's timing fit failed. Photon timing is given by the fitted time 
 1136 of the recorded waveform of the highest energetic crystal in a cluster; however, that fit can fail and so 
 1137 this variable tells the user if that has happened. 
 1139     REGISTER_VARIABLE("clusterErrorTiming", eclClusterErrorTiming, R
"DOC( 
 1140 Returns ECL cluster's timing uncertainty that contains :math:`99\%` of true photons (dt99). 
 1142 The photon timing uncertainty is currently determined using MC. The resulting parametrization depends on 
 1143 the true energy deposition in the highest energetic crystal and the local beam background level in that crystal. 
 1144 The resulting timing distribution is non-Gaussian and for each photon the value dt99 is stored, 
 1145 where :math:`|\text{timing}| / \text{dt99} < 1` is designed to give a :math:`99\%` 
 1146 timing efficiency for true photons from the IP. 
 1147 The resulting efficiency is approximately flat in energy and independent of beam background levels. 
 1149 Very large values of dt99 are an indication of failed waveform fits in the ECL. 
 1150 We remove such clusters in most physics photon lists. 
 1153     | Please read `this <importantNoteECL>` first. 
 1154     | Lower limit: :math:`0.0` 
 1155     | Upper limit: :math:`1000.0` 
 1156     | Precision: :math:`12` bit 
 1159     In real data there will be a sizeable number of high energetic Bhabha events 
 1160     (from previous or later bunch collisions) that can easily be rejected by timing cuts. 
 1161     However, these events create large ECL clusters that can overlap with other ECL clusters 
 1162     and it is not clear that a simple rejection is the correction strategy. 
 1166     REGISTER_VARIABLE("clusterHasFailedErrorTiming", eclClusterHasFailedErrorTiming, R"DOC( 
 1167 Status bit for if the ECL cluster's timing uncertainty calculation failed. Photon timing is given by the fitted time 
 1168 of the recorded waveform of the highest energetic crystal in a cluster; however, that fit can fail and so 
 1169 this variable tells the user if that has happened. 
 1171     REGISTER_VARIABLE("clusterHighestE", eclClusterHighestE, R
"DOC( 
 1172 Returns energy of the highest energetic crystal in the ECL cluster after reweighting. 
 1175     This variable must be used carefully since it can bias shower selection 
 1176     towards photons that hit crystals in the center and hence have a large energy 
 1177     deposition in the highest energy crystal. 
 1180     | Please read `this <importantNoteECL>` first. 
 1181     | Lower limit: :math:`-5` (:math:`e^{-5} = 0.00674\,` GeV) 
 1182     | Upper limit: :math:`3.0` (:math:`e^3 = 20.08553\,` GeV) 
 1183     | Precision: :math:`18` bit 
 1187     REGISTER_VARIABLE("clusterCellID", eclClusterCellId, 
 1188                       "Returns cellId of the crystal with highest energy in the ECLCluster."); 
 1189     REGISTER_VARIABLE("clusterThetaID", eclClusterThetaId, 
 1190                       "Returns thetaId of the crystal with highest energy in the ECLCluster."); 
 1191     REGISTER_VARIABLE("clusterPhiID", eclClusterPhiId, 
 1192                       "Returns phiId of the crystal with highest energy in the ECLCluster."); 
 1193     REGISTER_VARIABLE("clusterE1E9", eclClusterE1E9, R"DOC( 
 1194 Returns ratio of energies of the central crystal, E1, and 3x3 crystals, E9, around the central crystal. 
 1195 Since :math:`E1 \leq E9`, this ratio is :math:`\leq 1` and tends towards larger values for photons 
 1196 and smaller values for hadrons. 
 1199     | Please read `this <importantNoteECL>` first. 
 1200     | Lower limit: :math:`0.0` 
 1201     | Upper limit: :math:`1.0` 
 1202     | Precision: :math:`10` bit 
 1204     REGISTER_VARIABLE("clusterE9E25", eclClusterE9E25, R
"DOC( 
 1205 Deprecated - kept for backwards compatibility - returns clusterE9E21. 
 1207     REGISTER_VARIABLE("clusterE9E21", eclClusterE9E21, R
"DOC( 
 1208 Returns ratio of energies in inner 3x3 crystals, E9, and 5x5 crystals around the central crystal without corners. 
 1209 Since :math:`E9 \leq E21`, this ratio is :math:`\leq 1` and tends towards larger values for photons 
 1210 and smaller values for hadrons. 
 1213     | Please read `this <importantNoteECL>` first. 
 1214     | Lower limit: :math:`0.0` 
 1215     | Upper limit: :math:`1.0` 
 1216     | Precision: :math:`10` bit 
 1218     REGISTER_VARIABLE("clusterAbsZernikeMoment40", eclClusterAbsZernikeMoment40, R
"DOC( 
 1219 Returns absolute value of Zernike moment 40 (:math:`|Z_{40}|`). (shower shape variable). 
 1222     | Please read `this <importantNoteECL>` first. 
 1223     | Lower limit: :math:`0.0` 
 1224     | Upper limit: :math:`1.7` 
 1225     | Precision: :math:`10` bit 
 1227     REGISTER_VARIABLE("clusterAbsZernikeMoment51", eclClusterAbsZernikeMoment51, R
"DOC( 
 1228 Returns absolute value of Zernike moment 51 (:math:`|Z_{51}|`). (shower shape variable). 
 1231     | Please read `this <importantNoteECL>` first. 
 1232     | Lower limit: :math:`0.0` 
 1233     | Upper limit: :math:`1.2` 
 1234     | Precision: :math:`10` bit 
 1236     REGISTER_VARIABLE("clusterZernikeMVA", eclClusterZernikeMVA, R
"DOC( 
 1237 Returns output of a MVA using eleven Zernike moments of the cluster. Zernike moments are calculated per 
 1238 shower in a plane perpendicular to the shower direction via 
 1241     |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| 
 1243 where n, m are the integers, :math:`i` runs over the crystals in the shower, 
 1244 :math:`E_{i}` is the energy of the i-th crystal in the shower, 
 1245 :math:`R_{nm}` is a polynomial of degree :math:`n`, 
 1246 :math:`\rho_{i}` is the radial distance of the :math:`i`-th crystal in the perpendicular plane, 
 1247 and :math:`\alpha_{i}` is the polar angle of the :math:`i`-th crystal in the perpendicular plane. 
 1248 As a crystal can be related to more than one shower, :math:`w_{i}` is the fraction of the 
 1249 energy of the :math:`i`-th crystal associated with the shower. 
 1251 More details about the implementation can be found in `BELLE2-NOTE-TE-2017-001 <https://docs.belle2.org/record/454?ln=en>`_ . 
 1253 More details about Zernike polynomials can be found in `Wikipedia <https://en.wikipedia.org/wiki/Zernike_polynomials>`_ . 
 1255 | For cluster with hypothesisId==N1: raw MVA output. 
 1256 | For cluster with hypothesisId==N2: 1 - prod{clusterZernikeMVA}, where the product is on all N1 showers 
 1257   belonging to the same connected region (shower shape variable). 
 1260     | Please read `this <importantNoteECL>` first. 
 1261     | Lower limit: :math:`0.0` 
 1262     | Upper limit: :math:`1.0` 
 1263     | Precision: :math:`10` bit 
 1265     REGISTER_VARIABLE("clusterSecondMoment", eclClusterSecondMoment, R
"DOC( 
 1266 Returns second moment :math:`S`. It is defined as: 
 1269     S = \frac{\sum_{i=0}^{n} w_{i} E_{i} r^2_{i}}{\sum_{i=0}^{n} w_{i} E_{i}} 
 1271 where :math:`E_{i} = (E_0, E_1, ...)` are the single crystal energies sorted by energy, :math:`w_{i}` is 
 1272 the crystal weight, and :math:`r_{i}` is the distance of the :math:`i`-th digit to the shower center projected 
 1273 to a plane perpendicular to the shower axis. 
 1276     | Please read `this <importantNoteECL>` first. 
 1277     | Lower limit: :math:`0.0` 
 1278     | Upper limit: :math:`40.0` 
 1279     | Precision: :math:`10` bit 
 1282 )DOC",":math:`\\text{cm}^2`"); 
 1283     REGISTER_VARIABLE("clusterLAT", eclClusterLAT, R"DOC( 
 1284 Returns lateral energy distribution (shower variable). It is defined as following: 
 1287     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}} 
 1289 where :math:`E_{i} = (E_{0}, E_{1}, ...)` are the single crystal energies sorted by energy 
 1290 (:math:`E_{0}` is the highest energy and :math:`E_{1}` the second highest), :math:`w_{i}` 
 1291 is the crystal weight, :math:`r_{i}` is the distance of the :math:`i`-th digit to the 
 1292 shower center projected to a plane perpendicular to the shower axis, 
 1293 and :math:`r_{0} \approx 5\,cm` is the distance between two crystals. 
 1295 clusterLAT peaks around 0.3 for radially symmetrical electromagnetic showers and is larger 
 1296 for hadronic events, and electrons with a close-by radiative or Bremsstrahlung photon. 
 1299     | Please read `this <importantNoteECL>` first. 
 1300     | Lower limit: :math:`0.0` 
 1301     | Upper limit: :math:`1.0` 
 1302     | Precision: :math:`10` bit 
 1304     REGISTER_VARIABLE("clusterNHits", eclClusterNHits, R
"DOC( 
 1305 Returns sum of weights :math:`w_{i}` (:math:`w_{i} \leq 1`) of all crystals in an ECL cluster. 
 1306 For non-overlapping clusters this is equal to the number of crystals in the cluster. 
 1307 In case of energy splitting among nearby clusters, this can be a non-integer value. 
 1310     | Please read `this <importantNoteECL>` first. 
 1311     | Lower limit: :math:`0.0` 
 1312     | Upper limit: :math:`200.0` 
 1313     | Precision: :math:`10` bit 
 1314     | If fractional weights are not of interest, this value should be cast to the nearest integer. 
 1316     REGISTER_VARIABLE("clusterTrackMatch", eclClusterTrackMatched, R
"DOC( 
 1317 Returns 1.0 if at least one reconstructed charged track is matched to the ECL cluster. 
 1319 Every reconstructed charged track is extrapolated into the ECL. 
 1320 Every ECL crystal that is crossed by the track extrapolation is marked. 
 1321 Each ECL cluster that contains any marked crystal is matched to the track. 
 1322 Multiple tracks can be matched to one cluster and multiple clusters can be matched to one track. 
 1323 It is conceptually correct to have two tracks matched to the same cluster. 
 1325     REGISTER_VARIABLE("nECLClusterTrackMatches", nECLClusterTrackMatches, R
"DOC( 
 1326 Returns number of charged tracks matched to this cluster. 
 1329     Sometimes (perfectly correctly) two tracks are extrapolated into the same cluster. 
 1331     - For charged particles, this should return at least 1 (but sometimes 2 or more). 
 1332     - For neutrals, this should always return 0. 
 1333     - Returns NaN if there is no cluster. 
 1335     REGISTER_VARIABLE("clusterHasPulseShapeDiscrimination", eclClusterHasPulseShapeDiscrimination, R
"DOC( 
 1336 Status bit to indicate if cluster has digits with waveforms that passed energy and :math:`\chi^2` 
 1337 thresholds for computing PSD variables. 
 1339     REGISTER_VARIABLE("beamBackgroundSuppression", beamBackgroundSuppression, R
"DOC( 
 1340 Returns the output of an MVA classifier that uses shower-related variables to distinguish true photon clusters from beam background clusters. 
 1341 Class 1 is for true photon clusters while class 0 is for beam background clusters. 
 1343 The MVA has been trained using MC and the features used are: 
 1346 - `clusterPulseShapeDiscriminationMVA` 
 1349 - `clusterZernikeMVA` 
 1351 Both run-dependent and run-independent weights are available. For more information on this, and for usage recommendations, please see 
 1352 the `Neutrals Performance Confluence Page <https://confluence.desy.de/display/BI/Neutrals+Performance>`_. 
 1354     REGISTER_VARIABLE("fakePhotonSuppression", fakePhotonSuppression, R
"DOC( 
 1355 Returns the output of an MVA classifier that uses shower-related variables to distinguish true photon clusters from fake photon clusters (e.g. split-offs, 
 1356 track-cluster matching failures etc.). Class 1 is for true photon clusters while class 0 is for fake photon clusters.  
 1358 The MVA has been trained using MC and the features are: 
 1360 - `clusterPulseShapeDiscriminationMVA` 
 1362 - `clusterZernikeMVA` 
 1367 This 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.  
 1369 Both run-dependent and run-independent weights are available. For more information on this, and for usage recommendations, please see 
 1370 the `Neutrals Performance Confluence Page <https://confluence.desy.de/display/BI/Neutrals+Performance>`_. 
 1372     REGISTER_VARIABLE("hadronicSplitOffSuppression", hadronicSplitOffSuppression, R
"DOC( 
 1373 Returns the output of an MVA classifier that uses shower-related variables to distinguish true photon clusters from hadronic splitoff clusters. 
 1376 - 1 for true photon clusters 
 1377 - 0 for hadronic splitoff clusters 
 1379 The 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):  
 1381 - `clusterPulseShapeDiscriminationMVA` 
 1383 - `clusterZernikeMVA` 
 1387 - `clusterSecondMoment` 
 1389     MAKE_DEPRECATED("hadronicSplitOffSuppression", 
false, 
"light-2302-genetta", R
"DOC( 
 1390                      Use the variable `fakePhotonSuppression` instead which is maintained and uses the latest weight files. 
 1392     REGISTER_VARIABLE("clusterKlId", eclClusterKlId, R
"DOC( 
 1393 Returns MVA classifier that uses ECL clusters variables to discriminate Klong clusters from em background. 
 1398     REGISTER_VARIABLE("clusterPulseShapeDiscriminationMVA", eclPulseShapeDiscriminationMVA, R
"DOC( 
 1399 Returns MVA classifier that uses pulse shape discrimination to identify electromagnetic vs hadronic showers. 
 1401 - 1 for electromagnetic showers 
 1402 - 0 for hadronic showers 
 1404     REGISTER_VARIABLE("clusterNumberOfHadronDigits", eclClusterNumberOfHadronDigits, R
"DOC( 
 1405 Returns ECL cluster's number of hadron digits in cluster (pulse shape discrimination variable). 
 1406 Weighted sum of digits in cluster with significant scintillation emission (:math:`> 3\,` MeV) 
 1407 in the hadronic scintillation component. 
 1408 Computed only using cluster digits with energy :math:`> 50\,` MeV and good offline waveform fit :math:`\chi^2`. 
 1411     | Please read `this <importantNoteECL>` first. 
 1412     | Lower limit: :math:`0.0` 
 1413     | Upper limit: :math:`255.0` 
 1414     | Precision: :math:`18` bit 
 1416     REGISTER_VARIABLE("clusterClusterID", eclClusterId, R
"DOC( 
 1417 Returns ECL cluster ID of this ECL cluster within the connected region (CR) to which it belongs to. 
 1419     REGISTER_VARIABLE("clusterHasNPhotons", eclClusterHasNPhotonsHypothesis, R
"DOC( 
 1420 Returns 1.0 if cluster has the 'N photons' hypothesis (historically called 'N1'), 
 1421 0.0 if not, and NaN if no cluster is associated to the particle. 
 1423     REGISTER_VARIABLE("clusterHasNeutralHadron", eclClusterHasNeutralHadronHypothesis, R
"DOC( 
 1424 Returns 1.0 if the cluster has the 'neutral hadrons' hypothesis (historically called 'N2'), 
 1425 0.0 if not, and NaN if no cluster is associated to the particle. 
 1427     REGISTER_VARIABLE("eclExtTheta", eclExtTheta, R
"DOC( 
 1428 Returns extrapolated :math:`\theta` of particle track associated to the cluster (if any). Requires module ECLTrackCalDigitMatch to be executed. 
 1431     REGISTER_VARIABLE("eclExtPhi", eclExtPhi, R"DOC( 
 1432 Returns extrapolated :math:`\phi` of particle track associated to the cluster (if any). Requires module ECLTrackCalDigitMatch to be executed.. 
 1435     REGISTER_VARIABLE("eclExtPhiId", eclExtPhiId, R"DOC( 
 1436 Returns extrapolated :math:`\phi` ID of particle track associated to the cluster (if any). Requires module ECLTrackCalDigitMatch to be executed.. 
 1438     REGISTER_VARIABLE("weightedAverageECLTime", weightedAverageECLTime, R
"DOC( 
 1439 Returns ECL weighted average time of all clusters (neutrals) and matched clusters (charged) of daughters 
 1440 (of any generation) of the provided particle. 
 1443     REGISTER_VARIABLE(
"maxWeightedDistanceFromAverageECLTime", maxWeightedDistanceFromAverageECLTime, R
"DOC( 
 1444 Returns maximum weighted distance between time of the cluster of a photon and the ECL average time, amongst 
 1445 the clusters (neutrals) and matched clusters (charged) of daughters (of all generations) of the provided particle. 
 1448     REGISTER_VARIABLE(
"clusterMdstIndex", eclClusterMdstIndex, R
"DOC( 
 1449 StoreArray index(0 - based) of the MDST ECLCluster (useful for track-based particles matched to a cluster). 
 1452     REGISTER_VARIABLE("nECLOutOfTimeCrystals", nECLOutOfTimeCrystals, R
"DOC( 
 1453 [Eventbased] Returns the number of crystals (ECLCalDigits) that are out of time. 
 1456     REGISTER_VARIABLE("nECLOutOfTimeCrystalsFWDEndcap", nECLOutOfTimeCrystalsFWDEndcap, R
"DOC( 
 1457 [Eventbased] Returns the number of crystals (ECLCalDigits) that are out of time in the forward endcap. 
 1460     REGISTER_VARIABLE("nECLOutOfTimeCrystalsBarrel", nECLOutOfTimeCrystalsBarrel, R
"DOC( 
 1461 [Eventbased] Returns the number of crystals (ECLCalDigits) that are out of time in the barrel. 
 1464     REGISTER_VARIABLE("nECLOutOfTimeCrystalsBWDEndcap", nECLOutOfTimeCrystalsBWDEndcap, R
"DOC( 
 1465 [Eventbased] Returns the number of crystals (ECLCalDigits) that are out of time in the backward endcap. 
 1468     REGISTER_VARIABLE("nRejectedECLShowers", nRejectedECLShowers, R
"DOC( 
 1469 [Eventbased] Returns the number of showers in the ECL that do not become clusters. 
 1472     REGISTER_VARIABLE("nRejectedECLShowersFWDEndcap", nRejectedECLShowersFWDEndcap, R
"DOC( 
 1473 [Eventbased] Returns the number of showers in the ECL that do not become clusters, from the forward endcap. 
 1474 If the number exceeds 255 (uint8_t maximum value) the variable is set to 255. 
 1477     REGISTER_VARIABLE("nRejectedECLShowersBarrel", nRejectedECLShowersBarrel, R
"DOC( 
 1478 [Eventbased] Returns the number of showers in the ECL that do not become clusters, from the barrel. 
 1479 If the number exceeds 255 (uint8_t maximum value) the variable is set to 255. 
 1482     REGISTER_VARIABLE("nRejectedECLShowersBWDEndcap", nRejectedECLShowersBWDEndcap, R
"DOC( 
 1483 [Eventbased] Returns the number of showers in the ECL that do not become clusters, from the backward endcap. 
 1484 If the number exceeds 255 (uint8_t maximum value) the variable is set to 255. 
 1487     REGISTER_VARIABLE("eclClusterOnlyInvariantMass", eclClusterOnlyInvariantMass, R
"DOC( 
 1488 [Expert] The invariant mass calculated from all ECLCluster daughters (i.e. photons) and 
 1489 cluster-matched tracks using the cluster 4-momenta. 
 1491 Used for ECL-based dark sector physics and debugging track-cluster matching. 
 1493 )DOC","GeV/:math:`\\text{c}^2`"); 
 1495     REGISTER_METAVARIABLE("photonHasOverlap(cutString, photonlistname, tracklistname)", photonHasOverlap, R"DOC( 
 1496       Returns true if the connected ECL region of the particle's cluster is shared by another particle's cluster. 
 1497       Neutral and charged cluster are considered. 
 1498       A cut string can be provided to ignore cluster that do not satisfy the given criteria. 
 1499       By default, the ParticleList ``gamma:all`` is used for the check of neutral ECL cluster, 
 1500       and the ParticleList ``e-:all`` for the check of charged ECL cluster. 
 1501       However, one can customize the name of the ParticleLists via additional arguments. 
 1502       If no argument or only a cut string is provided and ``gamma:all`` or ``e-:all`` does not exist 
 1503       or if the variable is requested for a particle that is not a photon, NaN is returned. 
 1504       )DOC", Manager::VariableDataType::c_double); 
 1506     REGISTER_VARIABLE("clusterUncorrE", eclClusterUncorrectedE, R
"DOC( 
 1507 [Expert] [Calibration] Returns ECL cluster's uncorrected energy. That is, before leakage corrections. 
 1508 This variable should only be used for study of the ECL. Please see :b2:var:`clusterE`. 
int getPDGCode() const
PDG code.
static const ParticleSet chargedStableSet
set of charged stable particles
static const ParticleType photon
photon particle
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
@ 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...
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.
#define MAKE_DEPRECATED(name, make_fatal, version, description)
Registers a variable as deprecated.
Abstract base class for different kinds of events.