10#include <analysis/variables/MetaVariables.h>
11#include <analysis/variables/MCTruthVariables.h>
13#include <analysis/VariableManager/Utility.h>
14#include <analysis/dataobjects/Particle.h>
15#include <analysis/dataobjects/ParticleList.h>
16#include <analysis/dataobjects/RestOfEvent.h>
17#include <analysis/utility/PCmsLabTransform.h>
18#include <analysis/utility/ReferenceFrame.h>
19#include <analysis/utility/EvtPDLUtil.h>
20#include <analysis/utility/ParticleCopy.h>
21#include <analysis/utility/ValueIndexPairSorting.h>
22#include <analysis/ClusterUtility/ClusterUtils.h>
23#include <analysis/variables/VariableFormulaConstructor.h>
25#include <framework/logging/Logger.h>
26#include <framework/datastore/StoreArray.h>
27#include <framework/datastore/StoreObjPtr.h>
28#include <framework/dataobjects/EventExtraInfo.h>
29#include <framework/utilities/Conversion.h>
30#include <framework/utilities/MakeROOTCompatible.h>
31#include <framework/gearbox/Const.h>
33#include <mdst/dataobjects/Track.h>
34#include <mdst/dataobjects/MCParticle.h>
35#include <mdst/dataobjects/ECLCluster.h>
36#include <mdst/dataobjects/TrackFitResult.h>
38#include <boost/algorithm/string.hpp>
46#include <TDatabasePDG.h>
47#include <Math/Vector4D.h>
58 if (arguments.size() == 1) {
60 auto func = [var](
const Particle * particle) ->
double {
61 UseReferenceFrame<RestFrame> frame(particle);
62 double result = std::get<double>(var->function(particle));
67 B2FATAL(
"Wrong number of arguments for meta function useRestFrame");
73 if (arguments.size() == 1) {
75 auto func = [var](
const Particle * particle) ->
double {
76 UseReferenceFrame<CMSFrame> frame;
77 double result = std::get<double>(var->function(particle));
82 B2FATAL(
"Wrong number of arguments for meta function useCMSFrame");
88 if (arguments.size() == 1) {
90 auto func = [var](
const Particle * particle) ->
double {
91 UseReferenceFrame<LabFrame> frame;
92 double result = std::get<double>(var->function(particle));
97 B2FATAL(
"Wrong number of arguments for meta function useLabFrame");
103 if (arguments.size() == 2) {
106 int daughterIndexTagB = 0;
108 daughterIndexTagB = Belle2::convertString<int>(arguments[1]);
109 }
catch (std::invalid_argument&) {
110 B2FATAL(
"Second argument of useTagSideRecoilRestFrame meta function must be integer!");
113 auto func = [var, daughterIndexTagB](
const Particle * particle) ->
double {
114 if (particle->getPDGCode() != 300553)
116 B2ERROR(
"Variable should only be used on a Upsilon(4S) Particle List!");
121 ROOT::Math::PxPyPzEVector pSigB = T.getBeamFourMomentum() - particle->getDaughter(daughterIndexTagB)->get4Vector();
122 Particle tmp(pSigB, -particle->getDaughter(daughterIndexTagB)->getPDGCode());
124 UseReferenceFrame<RestFrame> frame(&tmp);
125 double result = std::get<double>(var->function(particle));
131 B2FATAL(
"Wrong number of arguments for meta function useTagSideRecoilRestFrame");
137 if (arguments.size() == 2) {
139 std::string listName = arguments[1];
140 auto func = [var, listName](
const Particle * particle) ->
double {
141 StoreObjPtr<ParticleList> list(listName);
142 unsigned listSize = list->getListSize();
146 B2WARNING(
"The selected ParticleList contains more than 1 Particles in this event. The variable useParticleRestFrame will use only the first candidate, and the result may not be the expected one."
147 <<
LogVar(
"ParticleList", listName)
148 <<
LogVar(
"Number of candidates in the list", listSize));
149 const Particle* p = list->getParticle(0);
150 UseReferenceFrame<RestFrame> frame(p);
151 double result = std::get<double>(var->function(particle));
156 B2FATAL(
"Wrong number of arguments for meta function useParticleRestFrame.");
162 if (arguments.size() == 2) {
164 std::string listName = arguments[1];
165 auto func = [var, listName](
const Particle * particle) ->
double {
166 StoreObjPtr<ParticleList> list(listName);
167 unsigned listSize = list->getListSize();
171 B2WARNING(
"The selected ParticleList contains more than 1 Particles in this event. The variable useParticleRestFrame will use only the first candidate, and the result may not be the expected one."
172 <<
LogVar(
"ParticleList", listName)
173 <<
LogVar(
"Number of candidates in the list", listSize));
174 const Particle* p = list->getParticle(0);
176 ROOT::Math::PxPyPzEVector recoil = T.getBeamFourMomentum() - p->get4Vector();
178 Particle pRecoil(recoil, 0);
179 pRecoil.setVertex(particle->getVertex());
180 UseReferenceFrame<RestFrame> frame(&pRecoil);
181 double result = std::get<double>(var->function(particle));
186 B2FATAL(
"Wrong number of arguments for meta function useParticleRestFrame.");
192 if (arguments.size() >= 2) {
194 auto func = [var, arguments](
const Particle * particle) ->
double {
197 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
199 for (
unsigned int i = 1; i < arguments.size(); i++)
201 auto generalizedIndex = arguments[i];
202 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
204 pSum += dauPart->get4Vector();
208 Particle tmp(pSum, 0);
209 UseReferenceFrame<RestFrame> frame(&tmp);
210 double result = std::get<double>(var->function(particle));
215 B2FATAL(
"Wrong number of arguments for meta function useDaughterRestFrame.");
221 if (arguments.size() == 1) {
223 auto func = [var](
const Particle * particle) ->
double {
224 int index = ancestorBIndex(particle);
226 StoreArray<MCParticle> mcparticles;
227 Particle temp(mcparticles[index]);
228 UseReferenceFrame<RestFrame> frame(&temp);
229 double result = std::get<double>(var->function(particle));
234 B2FATAL(
"Wrong number of arguments for meta function useMCancestorBRestFrame.");
240 if (arguments.size() == 1) {
241 auto extraInfoName = arguments[0];
242 auto func = [extraInfoName](
const Particle * particle) ->
double {
243 if (particle ==
nullptr)
245 B2WARNING(
"Returns NaN because the particle is nullptr! If you want EventExtraInfo variables, please use eventExtraInfo() instead");
248 if (particle->hasExtraInfo(extraInfoName))
250 return particle->getExtraInfo(extraInfoName);
258 B2FATAL(
"Wrong number of arguments for meta function extraInfo");
264 if (arguments.size() == 1) {
265 auto extraInfoName = arguments[0];
266 auto func = [extraInfoName](
const Particle*) ->
double {
267 StoreObjPtr<EventExtraInfo> eventExtraInfo;
268 if (not eventExtraInfo.isValid())
270 if (eventExtraInfo->hasExtraInfo(extraInfoName))
272 return eventExtraInfo->getExtraInfo(extraInfoName);
280 B2FATAL(
"Wrong number of arguments for meta function extraInfo");
286 if (arguments.size() == 1) {
289 auto func = [var, key](
const Particle*) ->
double {
291 StoreObjPtr<EventExtraInfo> eventExtraInfo;
292 if (not eventExtraInfo.isValid())
293 eventExtraInfo.create();
294 if (eventExtraInfo->hasExtraInfo(key))
296 return eventExtraInfo->getExtraInfo(key);
300 auto var_result = var->function(
nullptr);
301 if (std::holds_alternative<double>(var_result)) {
302 value = std::get<double>(var_result);
303 }
else if (std::holds_alternative<int>(var_result)) {
304 return std::get<int>(var_result);
305 }
else if (std::holds_alternative<bool>(var_result)) {
306 return std::get<bool>(var_result);
308 eventExtraInfo->addExtraInfo(key, value);
314 B2FATAL(
"Wrong number of arguments for meta function eventCached");
320 if (arguments.size() == 1) {
323 auto func = [var, key](
const Particle * particle) ->
double {
325 if (particle->hasExtraInfo(key))
327 return particle->getExtraInfo(key);
330 double value = std::get<double>(var->function(particle));
338 const_cast<Particle*
>(particle)->addExtraInfo(key, value);
344 B2FATAL(
"Wrong number of arguments for meta function particleCached");
353 if (arguments.size() != 1) B2FATAL(
"Wrong number of arguments for meta function formula");
354 FormulaParser<VariableFormulaConstructor> parser;
356 return parser.parse(arguments[0]);
357 }
catch (std::runtime_error& e) {
364 if (arguments.size() <= 1) {
366 std::string cutString;
367 if (arguments.size() == 1)
368 cutString = arguments[0];
370 auto func = [cut](
const Particle*) ->
int {
372 int number_of_tracks = 0;
373 StoreArray<Track> tracks;
374 for (
const auto& track : tracks)
376 const TrackFitResult* trackFit = track.getTrackFitResultWithClosestMass(
Const::pion);
377 if (trackFit->getChargeSign() == 0) {
381 if (cut->check(&particle))
386 return number_of_tracks;
391 B2FATAL(
"Wrong number of arguments for meta function nCleanedTracks");
397 if (arguments.size() <= 1) {
399 std::string cutString;
400 if (arguments.size() == 1)
401 cutString = arguments[0];
403 auto func = [cut](
const Particle*) ->
int {
405 int number_of_clusters = 0;
406 StoreArray<ECLCluster> clusters;
407 for (
const auto& cluster : clusters)
413 Particle particle(&cluster);
414 if (cut->check(&particle))
415 number_of_clusters++;
418 return number_of_clusters;
423 B2FATAL(
"Wrong number of arguments for meta function nCleanedECLClusters");
429 if (arguments.size() == 1) {
430 std::string cutString = arguments[0];
432 auto func = [cut](
const Particle * particle) ->
bool {
433 if (cut->check(particle))
440 B2FATAL(
"Wrong number of arguments for meta function passesCut");
446 if (arguments.size() == 1) {
447 std::string cutString = arguments[0];
449 auto func = [cut](
const Particle*) ->
bool {
450 if (cut->check(
nullptr))
457 B2FATAL(
"Wrong number of arguments for meta function passesEventCut");
463 if (arguments.size() == 2) {
466 pdgCode = Belle2::convertString<int>(arguments[0]);
467 }
catch (std::invalid_argument&) {
468 B2FATAL(
"The first argument of varFor meta function must be a positive integer!");
471 auto func = [pdgCode, var](
const Particle * particle) ->
double {
472 if (std::abs(particle->getPDGCode()) == std::abs(pdgCode))
474 auto var_result = var->function(particle);
475 if (std::holds_alternative<double>(var_result)) {
476 return std::get<double>(var_result);
477 }
else if (std::holds_alternative<int>(var_result)) {
478 return std::get<int>(var_result);
479 }
else if (std::holds_alternative<bool>(var_result)) {
480 return std::get<bool>(var_result);
486 B2FATAL(
"Wrong number of arguments for meta function varFor");
492 if (arguments.size() == 1) {
494 auto func = [var](
const Particle * particle) ->
double {
495 if (particle->getMCParticle())
500 auto var_result = var->function(particle);
501 if (std::holds_alternative<double>(var_result)) {
502 return std::get<double>(var_result);
503 }
else if (std::holds_alternative<int>(var_result)) {
504 return std::get<int>(var_result);
505 }
else if (std::holds_alternative<bool>(var_result)) {
506 return std::get<bool>(var_result);
513 B2FATAL(
"Wrong number of arguments for meta function varForMCGen");
519 if (arguments.size() == 1) {
520 std::string listName = arguments[0];
521 auto func = [listName](
const Particle * particle) ->
int {
524 StoreObjPtr<ParticleList> listOfParticles(listName);
526 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to nParticlesInList");
528 return listOfParticles->getListSize();
533 B2FATAL(
"Wrong number of arguments for meta function nParticlesInList");
540 if (arguments.size() != 1) {
541 B2FATAL(
"Wrong number of arguments for isInList");
543 auto listName = arguments[0];
545 auto func = [listName](
const Particle * particle) ->
bool {
548 StoreObjPtr<ParticleList> list(listName);
549 if (!(list.isValid()))
551 B2FATAL(
"Invalid Listname " << listName <<
" given to isInList");
555 return list->contains(particle);
564 if (arguments.size() != 1) {
565 B2FATAL(
"Wrong number of arguments for sourceObjectIsInList");
567 auto listName = arguments[0];
569 auto func = [listName](
const Particle * particle) ->
int {
572 StoreObjPtr<ParticleList> list(listName);
573 if (!(list.isValid()))
575 B2FATAL(
"Invalid Listname " << listName <<
" given to sourceObjectIsInList");
581 if (particlesource == Particle::EParticleSourceObject::c_Composite
582 or particlesource == Particle::EParticleSourceObject::c_Undefined)
588 for (
unsigned i = 0; i < list->getListSize(); ++i)
590 Particle* iparticle = list->getParticle(i);
591 if (particlesource == iparticle->getParticleSource())
592 if (particle->getMdstArrayIndex() == iparticle->getMdstArrayIndex())
604 if (arguments.size() != 1) {
605 B2FATAL(
"Wrong number of arguments for mcParticleIsInMCList");
607 auto listName = arguments[0];
609 auto func = [listName](
const Particle * particle) ->
bool {
612 StoreObjPtr<ParticleList> list(listName);
613 if (!(list.isValid()))
614 B2FATAL(
"Invalid Listname " << listName <<
" given to mcParticleIsInMCList");
617 const MCParticle* mcp = particle->getMCParticle();
618 if (mcp ==
nullptr)
return false;
621 for (
unsigned i = 0; i < list->getListSize(); ++i)
623 const MCParticle* imcp = list->getParticle(i)->
getMCParticle();
624 if ((imcp !=
nullptr) and (mcp->getArrayIndex() == imcp->getArrayIndex()))
634 B2WARNING(
"isDaughterOfList is outdated and replaced by isDescendantOfList.");
635 std::vector<std::string> new_arguments = arguments;
636 new_arguments.push_back(std::string(
"1"));
637 return isDescendantOfList(new_arguments);
642 B2WARNING(
"isGrandDaughterOfList is outdated and replaced by isDescendantOfList.");
643 std::vector<std::string> new_arguments = arguments;
644 new_arguments.push_back(std::string(
"2"));
645 return isDescendantOfList(new_arguments);
650 if (arguments.size() > 0) {
651 auto listNames = arguments;
652 auto func = [listNames](
const Particle * particle) ->
bool {
654 int generation_flag = -1;
657 generation_flag = Belle2::convertString<int>(listNames.back());
658 }
catch (std::exception& e) {}
660 for (
auto& iListName : listNames)
663 Belle2::convertString<int>(iListName);
665 }
catch (std::exception& e) {}
668 auto list_comparison = [](
auto&& self,
const Particle * m,
const Particle * p,
int flag)->
bool {
670 for (
unsigned i = 0; i < m->getNDaughters(); ++i)
672 const Particle* daughter = m->getDaughter(i);
673 if ((flag == 1.) or (flag < 0)) {
674 if (p->isCopyOf(daughter)) {
680 if (daughter->getNDaughters() > 0) {
681 result = self(self, daughter, p, flag - 1);
691 StoreObjPtr<ParticleList> listOfParticles(iListName);
693 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << iListName <<
" given to isDescendantOfList");
695 for (
unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
696 Particle* iParticle = listOfParticles->getParticle(i);
697 output = list_comparison(list_comparison, iParticle, particle, generation_flag);
707 B2FATAL(
"Wrong number of arguments for meta function isDescendantOfList");
713 if (arguments.size() > 0) {
714 auto listNames = arguments;
715 auto func = [listNames](
const Particle * particle) ->
bool {
717 int generation_flag = -1;
720 generation_flag = Belle2::convertString<int>(listNames.back());
721 }
catch (std::exception& e) {}
723 if (particle->getMCParticle() ==
nullptr)
728 for (
auto& iListName : listNames)
731 std::stod(iListName);
733 }
catch (std::exception& e) {}
735 auto list_comparison = [](
auto&& self,
const Particle * m,
const Particle * p,
int flag)->
bool {
737 for (
unsigned i = 0; i < m->getNDaughters(); ++i)
739 const Particle* daughter = m->getDaughter(i);
740 if ((flag == 1.) or (flag < 0)) {
741 if (daughter->getMCParticle() !=
nullptr) {
742 if (p->getMCParticle()->getArrayIndex() == daughter->getMCParticle()->getArrayIndex()) {
748 if (daughter->getNDaughters() > 0) {
749 result = self(self, daughter, p, flag - 1);
759 StoreObjPtr<ParticleList> listOfParticles(iListName);
761 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << iListName <<
" given to isMCDescendantOfList");
763 for (
unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
764 Particle* iParticle = listOfParticles->getParticle(i);
765 output = list_comparison(list_comparison, iParticle, particle, generation_flag);
775 B2FATAL(
"Wrong number of arguments for meta function isMCDescendantOfList");
781 if (arguments.size() == 1) {
783 auto func = [var](
const Particle * particle) ->
double {
784 double product = 1.0;
785 if (particle->getNDaughters() == 0)
789 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
791 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
792 product *= std::get<double>(var->function(particle->getDaughter(j)));
794 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
796 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
797 product *= std::get<int>(var->function(particle->getDaughter(j)));
804 B2FATAL(
"Wrong number of arguments for meta function daughterProductOf");
810 if (arguments.size() == 1) {
812 auto func = [var](
const Particle * particle) ->
double {
814 if (particle->getNDaughters() == 0)
818 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
820 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
821 sum += std::get<double>(var->function(particle->getDaughter(j)));
823 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
825 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
826 sum += std::get<int>(var->function(particle->getDaughter(j)));
833 B2FATAL(
"Wrong number of arguments for meta function daughterSumOf");
839 if (arguments.size() == 1) {
841 auto func = [var](
const Particle * particle) ->
double {
843 if (particle->getNDaughters() == 0)
847 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
849 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
850 double iValue = std::get<double>(var->function(particle->getDaughter(j)));
851 if (std::isnan(iValue))
continue;
852 if (std::isnan(min)) min = iValue;
853 if (iValue < min) min = iValue;
855 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
857 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
858 int iValue = std::get<int>(var->function(particle->getDaughter(j)));
859 if (std::isnan(min)) min = iValue;
860 if (iValue < min) min = iValue;
867 B2FATAL(
"Wrong number of arguments for meta function daughterLowest");
873 if (arguments.size() == 1) {
875 auto func = [var](
const Particle * particle) ->
double {
877 if (particle->getNDaughters() == 0)
881 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
883 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
884 double iValue = std::get<double>(var->function(particle->getDaughter(j)));
885 if (std::isnan(iValue))
continue;
886 if (std::isnan(max)) max = iValue;
887 if (iValue > max) max = iValue;
889 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
891 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
892 int iValue = std::get<int>(var->function(particle->getDaughter(j)));
893 if (std::isnan(max)) max = iValue;
894 if (iValue > max) max = iValue;
901 B2FATAL(
"Wrong number of arguments for meta function daughterHighest");
907 if (arguments.size() == 3) {
908 auto func = [arguments](
const Particle * particle) ->
double {
909 if (particle ==
nullptr)
911 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
912 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
913 auto variablename = arguments[2];
914 if (dau_i ==
nullptr || dau_j ==
nullptr)
916 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
920 auto result_j = var->function(dau_j);
921 auto result_i = var->function(dau_i);
923 if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i))
925 diff = std::get<double>(result_j) - std::get<double>(result_i);
926 }
else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i))
928 diff = std::get<int>(result_j) - std::get<int>(result_i);
931 throw std::runtime_error(
"Bad variant access");
933 if (variablename ==
"phi" or variablename ==
"clusterPhi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)"))
934 or std::regex_match(variablename, std::regex(
"use.*Frame\\(clusterPhi\\)")))
936 if (fabs(diff) > M_PI) {
938 diff = diff - 2 * M_PI;
940 diff = 2 * M_PI + diff;
948 B2FATAL(
"Wrong number of arguments for meta function daughterDiffOf");
954 if (arguments.size() == 3) {
955 auto func = [arguments](
const Particle * particle) ->
double {
956 if (particle ==
nullptr)
958 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
959 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
960 auto variablename = arguments[2];
961 if (dau_i ==
nullptr || dau_j ==
nullptr)
963 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
968 if (iMcDaughter ==
nullptr || jMcDaughter ==
nullptr)
970 Particle iTmpPart(iMcDaughter);
971 Particle jTmpPart(jMcDaughter);
973 auto result_j = var->function(&jTmpPart);
974 auto result_i = var->function(&iTmpPart);
976 if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i))
978 diff = std::get<double>(result_j) - std::get<double>(result_i);
979 }
else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i))
981 diff = std::get<int>(result_j) - std::get<int>(result_i);
984 throw std::runtime_error(
"Bad variant access");
986 if (variablename ==
"phi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)")))
988 if (fabs(diff) > M_PI) {
990 diff = diff - 2 * M_PI;
992 diff = 2 * M_PI + diff;
1000 B2FATAL(
"Wrong number of arguments for meta function mcDaughterDiffOf");
1006 if (arguments.size() == 5) {
1008 Belle2::convertString<int>(arguments[0]);
1009 Belle2::convertString<int>(arguments[1]);
1010 Belle2::convertString<int>(arguments[2]);
1011 Belle2::convertString<int>(arguments[3]);
1012 }
catch (std::invalid_argument&) {
1013 B2FATAL(
"First four arguments of grandDaughterDiffOf meta function must be integers!");
1015 std::vector<std::string> new_arguments;
1016 new_arguments.push_back(std::string(arguments[0] +
":" + arguments[2]));
1017 new_arguments.push_back(std::string(arguments[1] +
":" + arguments[3]));
1018 new_arguments.push_back(arguments[4]);
1019 return daughterDiffOf(new_arguments);
1021 B2FATAL(
"Wrong number of arguments for meta function grandDaughterDiffOf");
1027 std::vector<std::string> new_arguments = arguments;
1028 new_arguments.push_back(std::string(
"phi"));
1029 return daughterDiffOf(new_arguments);
1034 std::vector<std::string> new_arguments = arguments;
1035 new_arguments.push_back(std::string(
"phi"));
1036 return mcDaughterDiffOf(new_arguments);
1041 std::vector<std::string> new_arguments = arguments;
1042 new_arguments.push_back(std::string(
"phi"));
1043 return grandDaughterDiffOf(new_arguments);
1048 std::vector<std::string> new_arguments = arguments;
1049 new_arguments.push_back(std::string(
"clusterPhi"));
1050 return daughterDiffOf(new_arguments);
1055 std::vector<std::string> new_arguments = arguments;
1056 new_arguments.push_back(std::string(
"clusterPhi"));
1057 return grandDaughterDiffOf(new_arguments);
1062 std::vector<std::string> new_arguments = arguments;
1063 new_arguments.push_back(std::string(
"useCMSFrame(phi)"));
1064 return daughterDiffOf(new_arguments);
1069 std::vector<std::string> new_arguments = arguments;
1070 new_arguments.push_back(std::string(
"useCMSFrame(phi)"));
1071 return mcDaughterDiffOf(new_arguments);
1076 std::vector<std::string> new_arguments = arguments;
1077 new_arguments.push_back(std::string(
"useCMSFrame(clusterPhi)"));
1078 return daughterDiffOf(new_arguments);
1083 if (arguments.size() == 3) {
1084 auto func = [arguments](
const Particle * particle) ->
double {
1085 if (particle ==
nullptr)
1087 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
1088 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
1089 if (!(dau_i && dau_j))
1091 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
1095 double iValue, jValue;
1096 if (std::holds_alternative<double>(var->function(dau_j)))
1098 iValue = std::get<double>(var->function(dau_i));
1099 jValue = std::get<double>(var->function(dau_j));
1100 }
else if (std::holds_alternative<int>(var->function(dau_j)))
1102 iValue = std::get<int>(var->function(dau_i));
1103 jValue = std::get<int>(var->function(dau_j));
1105 return (jValue - iValue) / (jValue + iValue);
1109 B2FATAL(
"Wrong number of arguments for meta function daughterNormDiffOf");
1115 if (arguments.size() == 2) {
1116 int daughterNumber = 0;
1118 daughterNumber = Belle2::convertString<int>(arguments[0]);
1119 }
catch (std::invalid_argument&) {
1120 B2FATAL(
"First argument of daughterMotherDiffOf meta function must be integer!");
1122 auto variablename = arguments[1];
1123 auto func = [variablename, daughterNumber](
const Particle * particle) ->
double {
1124 if (particle ==
nullptr)
1126 if (daughterNumber >=
int(particle->getNDaughters()))
1131 auto result_mother = var->function(particle);
1132 auto result_daughter = var->function(particle->getDaughter(daughterNumber));
1134 if (std::holds_alternative<double>(result_mother) && std::holds_alternative<double>(result_daughter)) {
1135 diff = std::get<double>(result_mother) - std::get<double>(result_daughter);
1136 }
else if (std::holds_alternative<int>(result_mother) && std::holds_alternative<int>(result_daughter)) {
1137 diff = std::get<int>(result_mother) - std::get<int>(result_daughter);
1139 throw std::runtime_error(
"Bad variant access");
1142 if (variablename ==
"phi" or variablename ==
"useCMSFrame(phi)") {
1143 if (fabs(diff) > M_PI) {
1145 diff = diff - 2 * M_PI;
1147 diff = 2 * M_PI + diff;
1156 B2FATAL(
"Wrong number of arguments for meta function daughterMotherDiffOf");
1162 if (arguments.size() == 2) {
1163 int daughterNumber = 0;
1165 daughterNumber = Belle2::convertString<int>(arguments[0]);
1166 }
catch (std::invalid_argument&) {
1167 B2FATAL(
"First argument of daughterMotherDiffOf meta function must be integer!");
1170 auto func = [var, daughterNumber](
const Particle * particle) ->
double {
1171 if (particle ==
nullptr)
1173 if (daughterNumber >=
int(particle->getNDaughters()))
1177 double daughterValue = 0.0, motherValue = 0.0;
1178 if (std::holds_alternative<double>(var->function(particle))) {
1179 daughterValue = std::get<double>(var->function(particle->getDaughter(daughterNumber)));
1180 motherValue = std::get<double>(var->function(particle));
1181 }
else if (std::holds_alternative<int>(var->function(particle))) {
1182 daughterValue = std::get<int>(var->function(particle->getDaughter(daughterNumber)));
1183 motherValue = std::get<int>(var->function(particle));
1185 return (motherValue - daughterValue) / (motherValue + daughterValue);
1190 B2FATAL(
"Wrong number of arguments for meta function daughterMotherNormDiffOf");
1196 if (arguments.size() == 2 || arguments.size() == 3) {
1198 auto func = [arguments](
const Particle * particle) ->
double {
1199 if (particle ==
nullptr)
1202 std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1206 for (
auto& generalizedIndex : arguments)
1208 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1210 pDaus.push_back(frame.getMomentum(dauPart));
1212 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1218 if (pDaus.size() == 2)
1225 B2FATAL(
"Wrong number of arguments for meta function daughterAngle");
1229 double grandDaughterDecayAngle(
const Particle* particle,
const std::vector<double>& arguments)
1231 if (arguments.size() == 2) {
1236 int daughterIndex = std::lround(arguments[0]);
1237 if (daughterIndex >=
int(particle->getNDaughters()))
1239 const Particle* dau = particle->getDaughter(daughterIndex);
1241 int grandDaughterIndex = std::lround(arguments[1]);
1242 if (grandDaughterIndex >=
int(dau->getNDaughters()))
1245 B2Vector3D boost = dau->get4Vector().BoostToCM();
1247 ROOT::Math::PxPyPzEVector motherMomentum = - particle->get4Vector();
1248 motherMomentum = ROOT::Math::Boost(boost) * motherMomentum;
1250 ROOT::Math::PxPyPzEVector grandDaughterMomentum = dau->getDaughter(grandDaughterIndex)->get4Vector();
1251 grandDaughterMomentum = ROOT::Math::Boost(boost) * grandDaughterMomentum;
1256 B2FATAL(
"The variable grandDaughterDecayAngle needs exactly two integers as arguments!");
1262 if (arguments.size() == 2 || arguments.size() == 3) {
1264 auto func = [arguments](
const Particle * particle) ->
double {
1265 if (particle ==
nullptr)
1268 std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1272 if (particle->getParticleSource() == Particle::EParticleSourceObject::c_MCParticle)
1274 for (
auto& generalizedIndex : arguments) {
1275 const MCParticle* mcPart = particle->getMCParticle();
1276 if (mcPart ==
nullptr)
1279 if (dauMcPart ==
nullptr)
1282 pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1286 for (
auto& generalizedIndex : arguments) {
1287 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1288 if (dauPart ==
nullptr)
1292 if (dauMcPart ==
nullptr)
1295 pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1300 if (pDaus.size() == 2)
1307 B2FATAL(
"Wrong number of arguments for meta function mcDaughterAngle");
1311 double daughterClusterAngleInBetween(
const Particle* particle,
const std::vector<double>& daughterIndices)
1313 if (daughterIndices.size() == 2) {
1314 int daughterIndexi = std::lround(daughterIndices[0]);
1315 int daughterIndexj = std::lround(daughterIndices[1]);
1316 if (std::max(daughterIndexi, daughterIndexj) >=
int(particle->getNDaughters())) {
1319 const ECLCluster* clusteri = particle->getDaughter(daughterIndexi)->getECLCluster();
1320 const ECLCluster* clusterj = particle->getDaughter(daughterIndexj)->getECLCluster();
1321 if (clusteri and clusterj) {
1325 ClusterUtils clusutils;
1326 B2Vector3D pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit)).Vect();
1327 B2Vector3D pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit)).Vect();
1328 return pi.Angle(pj);
1332 }
else if (daughterIndices.size() == 3) {
1333 int daughterIndexi = std::lround(daughterIndices[0]);
1334 int daughterIndexj = std::lround(daughterIndices[1]);
1335 int daughterIndexk = std::lround(daughterIndices[2]);
1336 if (std::max(std::max(daughterIndexi, daughterIndexj), daughterIndexk) >=
int(particle->getNDaughters())) {
1339 const ECLCluster* clusteri = (particle->getDaughter(daughterIndices[0]))->getECLCluster();
1340 const ECLCluster* clusterj = (particle->getDaughter(daughterIndices[1]))->getECLCluster();
1341 const ECLCluster* clusterk = (particle->getDaughter(daughterIndices[2]))->getECLCluster();
1342 if (clusteri and clusterj and clusterk) {
1347 ClusterUtils clusutils;
1348 B2Vector3D pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit)).Vect();
1349 B2Vector3D pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit)).Vect();
1350 B2Vector3D pk = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterk, clusterkBit)).Vect();
1351 return pk.
Angle(pi + pj);
1356 B2FATAL(
"Wrong number of arguments for daughterClusterAngleInBetween!");
1362 if (arguments.size() > 1) {
1363 auto func = [arguments](
const Particle * particle) ->
double {
1365 ROOT::Math::PxPyPzEVector pSum;
1367 for (
auto& generalizedIndex : arguments)
1369 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1371 pSum += frame.getMomentum(dauPart);
1380 B2FATAL(
"Wrong number of arguments for meta function daughterInvM. At least two integers are needed.");
1386 if (arguments.size() == 2) {
1390 divideBy = Belle2::convertString<int>(arguments[1]);
1391 }
catch (std::invalid_argument&) {
1392 B2FATAL(
"Second argument of modulo meta function must be integer!");
1394 auto func = [var, divideBy](
const Particle * particle) ->
int {
1395 auto var_result = var->function(particle);
1396 if (std::holds_alternative<double>(var_result))
1398 return int(std::get<double>(var_result)) % divideBy;
1399 }
else if (std::holds_alternative<int>(var_result))
1401 return std::get<int>(var_result) % divideBy;
1402 }
else if (std::holds_alternative<bool>(var_result))
1404 return int(std::get<bool>(var_result)) % divideBy;
1409 B2FATAL(
"Wrong number of arguments for meta function modulo");
1415 if (arguments.size() == 1) {
1418 auto func = [var](
const Particle * particle) ->
bool {
return std::isnan(std::get<double>(var->function(particle))); };
1421 B2FATAL(
"Wrong number of arguments for meta function isNAN");
1427 if (arguments.size() == 2) {
1429 double defaultOutput;
1431 defaultOutput = Belle2::convertString<double>(arguments[1]);
1432 }
catch (std::invalid_argument&) {
1433 B2FATAL(
"The second argument of ifNANgiveX meta function must be a number!");
1435 auto func = [var, defaultOutput](
const Particle * particle) ->
double {
1436 double output = std::get<double>(var->function(particle));
1437 if (std::isnan(output))
return defaultOutput;
1442 B2FATAL(
"Wrong number of arguments for meta function ifNANgiveX");
1448 if (arguments.size() == 1) {
1451 auto func = [var](
const Particle * particle) ->
bool {
return std::isinf(std::get<double>(var->function(particle))); };
1454 B2FATAL(
"Wrong number of arguments for meta function isInfinity");
1460 if (arguments.size() >= 2) {
1466 for (
size_t i = 1; i < arguments.size(); ++i) {
1468 finalMask |= Belle2::convertString<int>(arguments[i]);
1469 }
catch (std::invalid_argument&) {
1470 B2FATAL(
"The input flags to meta function unmask() should be integer!");
1476 auto func = [var, finalMask](
const Particle * particle) ->
double {
1478 auto var_result = var->function(particle);
1479 if (std::holds_alternative<double>(var_result))
1482 if (std::isnan(std::get<double>(var_result))) {
1485 value = int(std::get<double>(var_result));
1486 }
else if (std::holds_alternative<int>(var_result))
1488 value = std::get<int>(var_result);
1492 value &= (~finalMask);
1499 B2FATAL(
"Meta function unmask needs at least two arguments!");
1505 if (arguments.size() == 3) {
1507 std::string cutString = arguments[0];
1513 auto func = [cut, variableIfTrue, variableIfFalse](
const Particle * particle) ->
double {
1514 if (particle ==
nullptr)
1516 if (cut->check(particle))
1518 auto var_result = variableIfTrue->function(particle);
1519 if (std::holds_alternative<double>(var_result)) {
1520 return std::get<double>(var_result);
1521 }
else if (std::holds_alternative<int>(var_result)) {
1522 return std::get<int>(var_result);
1523 }
else if (std::holds_alternative<bool>(var_result)) {
1524 return std::get<bool>(var_result);
1528 auto var_result = variableIfFalse->function(particle);
1529 if (std::holds_alternative<double>(var_result)) {
1530 return std::get<double>(var_result);
1531 }
else if (std::holds_alternative<int>(var_result)) {
1532 return std::get<int>(var_result);
1533 }
else if (std::holds_alternative<bool>(var_result)) {
1534 return std::get<bool>(var_result);
1541 B2FATAL(
"Wrong number of arguments for meta function conditionalVariableSelector");
1547 if (arguments.size() > 0) {
1548 std::vector<const Variable::Manager::Var*> variables;
1549 for (
auto& argument : arguments)
1552 auto func = [variables, arguments](
const Particle * particle) ->
double {
1553 double pValueProduct = 1.;
1554 for (
auto variable : variables)
1556 double pValue = std::get<double>(variable->function(particle));
1560 pValueProduct *= pValue;
1562 double pValueSum = 1.;
1563 double factorial = 1.;
1564 for (
unsigned int i = 1; i < arguments.size(); ++i)
1567 pValueSum += pow(-std::log(pValueProduct), i) / factorial;
1569 return pValueProduct * pValueSum;
1573 B2FATAL(
"Wrong number of arguments for meta function pValueCombination");
1579 if (arguments.size() == 1) {
1581 auto func = [var](
const Particle * particle) ->
double {
1582 auto var_result = var->function(particle);
1583 if (std::holds_alternative<double>(var_result))
1585 return std::abs(std::get<double>(var_result));
1586 }
else if (std::holds_alternative<int>(var_result))
1588 return std::abs(std::get<int>(var_result));
1593 B2FATAL(
"Wrong number of arguments for meta function abs");
1599 if (arguments.size() == 2) {
1604 B2FATAL(
"One or both of the used variables doesn't exist!");
1606 auto func = [var1, var2](
const Particle * particle) ->
double {
1608 auto var_result1 = var1->function(particle);
1609 auto var_result2 = var2->function(particle);
1610 if (std::holds_alternative<double>(var_result1))
1612 val1 = std::get<double>(var_result1);
1613 }
else if (std::holds_alternative<int>(var_result1))
1615 val1 = std::get<int>(var_result1);
1617 if (std::holds_alternative<double>(var_result2))
1619 val2 = std::get<double>(var_result2);
1620 }
else if (std::holds_alternative<int>(var_result2))
1622 val2 = std::get<int>(var_result2);
1624 return std::max(val1, val2);
1628 B2FATAL(
"Wrong number of arguments for meta function max");
1634 if (arguments.size() == 2) {
1639 B2FATAL(
"One or both of the used variables doesn't exist!");
1641 auto func = [var1, var2](
const Particle * particle) ->
double {
1643 auto var_result1 = var1->function(particle);
1644 auto var_result2 = var2->function(particle);
1645 if (std::holds_alternative<double>(var_result1))
1647 val1 = std::get<double>(var_result1);
1648 }
else if (std::holds_alternative<int>(var_result1))
1650 val1 = std::get<int>(var_result1);
1652 if (std::holds_alternative<double>(var_result2))
1654 val2 = std::get<double>(var_result2);
1655 }
else if (std::holds_alternative<int>(var_result2))
1657 val2 = std::get<int>(var_result2);
1659 return std::min(val1, val2);
1663 B2FATAL(
"Wrong number of arguments for meta function min");
1669 if (arguments.size() == 1) {
1671 auto func = [var](
const Particle * particle) ->
double {
1672 auto var_result = var->function(particle);
1673 if (std::holds_alternative<double>(var_result))
1674 return std::sin(std::get<double>(var_result));
1675 else if (std::holds_alternative<int>(var_result))
1676 return std::sin(std::get<int>(var_result));
1681 B2FATAL(
"Wrong number of arguments for meta function sin");
1687 if (arguments.size() == 1) {
1689 auto func = [var](
const Particle * particle) ->
double {
1690 auto var_result = var->function(particle);
1691 if (std::holds_alternative<double>(var_result))
1692 return std::asin(std::get<double>(var_result));
1693 else if (std::holds_alternative<int>(var_result))
1694 return std::asin(std::get<int>(var_result));
1699 B2FATAL(
"Wrong number of arguments for meta function asin");
1705 if (arguments.size() == 1) {
1707 auto func = [var](
const Particle * particle) ->
double {
1708 auto var_result = var->function(particle);
1709 if (std::holds_alternative<double>(var_result))
1710 return std::cos(std::get<double>(var_result));
1711 else if (std::holds_alternative<int>(var_result))
1712 return std::cos(std::get<int>(var_result));
1717 B2FATAL(
"Wrong number of arguments for meta function cos");
1723 if (arguments.size() == 1) {
1725 auto func = [var](
const Particle * particle) ->
double {
1726 auto var_result = var->function(particle);
1727 if (std::holds_alternative<double>(var_result))
1728 return std::acos(std::get<double>(var_result));
1729 else if (std::holds_alternative<int>(var_result))
1730 return std::acos(std::get<int>(var_result));
1735 B2FATAL(
"Wrong number of arguments for meta function acos");
1741 if (arguments.size() == 1) {
1743 auto func = [var](
const Particle * particle) ->
double {
return std::tan(std::get<double>(var->function(particle))); };
1746 B2FATAL(
"Wrong number of arguments for meta function tan");
1752 if (arguments.size() == 1) {
1754 auto func = [var](
const Particle * particle) ->
double {
return std::atan(std::get<double>(var->function(particle))); };
1757 B2FATAL(
"Wrong number of arguments for meta function atan");
1763 if (arguments.size() == 1) {
1765 auto func = [var](
const Particle * particle) ->
double {
1766 auto var_result = var->function(particle);
1767 if (std::holds_alternative<double>(var_result))
1768 return std::exp(std::get<double>(var_result));
1769 else if (std::holds_alternative<int>(var_result))
1770 return std::exp(std::get<int>(var_result));
1775 B2FATAL(
"Wrong number of arguments for meta function exp");
1781 if (arguments.size() == 1) {
1783 auto func = [var](
const Particle * particle) ->
double {
1784 auto var_result = var->function(particle);
1785 if (std::holds_alternative<double>(var_result))
1786 return std::log(std::get<double>(var_result));
1787 else if (std::holds_alternative<int>(var_result))
1788 return std::log(std::get<int>(var_result));
1793 B2FATAL(
"Wrong number of arguments for meta function log");
1799 if (arguments.size() == 1) {
1801 auto func = [var](
const Particle * particle) ->
double {
1802 auto var_result = var->function(particle);
1803 if (std::holds_alternative<double>(var_result))
1804 return std::log10(std::get<double>(var_result));
1805 else if (std::holds_alternative<int>(var_result))
1806 return std::log10(std::get<int>(var_result));
1811 B2FATAL(
"Wrong number of arguments for meta function log10");
1817 if (arguments.size() == 1) {
1819 auto func = [var](
const Particle * particle) ->
double {
1820 if (particle ==
nullptr)
1823 StoreArray<Particle> particles;
1824 if (!particle->hasExtraInfo(
"original_index"))
1827 auto originalParticle = particles[particle->getExtraInfo(
"original_index")];
1828 if (!originalParticle)
1830 auto var_result = var->function(originalParticle);
1831 if (std::holds_alternative<double>(var_result))
1833 return std::get<double>(var_result);
1834 }
else if (std::holds_alternative<int>(var_result))
1836 return std::get<int>(var_result);
1837 }
else if (std::holds_alternative<bool>(var_result))
1839 return std::get<bool>(var_result);
1844 B2FATAL(
"Wrong number of arguments for meta function originalParticle");
1850 if (arguments.size() == 2) {
1851 int daughterNumber = 0;
1853 daughterNumber = Belle2::convertString<int>(arguments[0]);
1854 }
catch (std::invalid_argument&) {
1855 B2FATAL(
"First argument of daughter meta function must be integer!");
1858 auto func = [var, daughterNumber](
const Particle * particle) ->
double {
1859 if (particle ==
nullptr)
1861 if (daughterNumber >=
int(particle->getNDaughters()))
1865 auto var_result = var->function(particle->getDaughter(daughterNumber));
1866 if (std::holds_alternative<double>(var_result)) {
1867 return std::get<double>(var_result);
1868 }
else if (std::holds_alternative<int>(var_result)) {
1869 return std::get<int>(var_result);
1870 }
else if (std::holds_alternative<bool>(var_result)) {
1871 return std::get<bool>(var_result);
1877 B2FATAL(
"Wrong number of arguments for meta function daughter");
1883 if (arguments.size() == 2) {
1884 int daughterNumber = 0;
1886 daughterNumber = Belle2::convertString<int>(arguments[0]);
1887 }
catch (std::invalid_argument&) {
1888 B2FATAL(
"First argument of daughter meta function must be integer!");
1891 auto func = [var, daughterNumber](
const Particle * particle) ->
double {
1892 if (particle ==
nullptr)
1894 if (daughterNumber >=
int(particle->getNDaughters()))
1898 StoreArray<Particle> particles;
1899 if (!particle->getDaughter(daughterNumber)->hasExtraInfo(
"original_index"))
1901 auto originalDaughter = particles[particle->getDaughter(daughterNumber)->getExtraInfo(
"original_index")];
1902 if (!originalDaughter)
1905 auto var_result = var->function(originalDaughter);
1906 if (std::holds_alternative<double>(var_result)) {
1907 return std::get<double>(var_result);
1908 }
else if (std::holds_alternative<int>(var_result)) {
1909 return std::get<int>(var_result);
1910 }
else if (std::holds_alternative<bool>(var_result)) {
1911 return std::get<bool>(var_result);
1917 B2FATAL(
"Wrong number of arguments for meta function daughter");
1923 if (arguments.size() == 2) {
1924 int daughterNumber = 0;
1926 daughterNumber = Belle2::convertString<int>(arguments[0]);
1927 }
catch (std::invalid_argument&) {
1928 B2FATAL(
"First argument of mcDaughter meta function must be integer!");
1931 auto func = [var, daughterNumber](
const Particle * particle) ->
double {
1932 if (particle ==
nullptr)
1934 if (particle->getMCParticle())
1936 if (daughterNumber >=
int(particle->getMCParticle()->getNDaughters())) {
1939 Particle tempParticle = Particle(particle->getMCParticle()->getDaughters().at(daughterNumber));
1940 auto var_result = var->function(&tempParticle);
1941 if (std::holds_alternative<double>(var_result)) {
1942 return std::get<double>(var_result);
1943 }
else if (std::holds_alternative<int>(var_result)) {
1944 return std::get<int>(var_result);
1945 }
else if (std::holds_alternative<bool>(var_result)) {
1946 return std::get<bool>(var_result);
1955 B2FATAL(
"Wrong number of arguments for meta function mcDaughter");
1961 if (arguments.size() == 1) {
1963 auto func = [var](
const Particle * particle) ->
double {
1964 if (particle ==
nullptr)
1966 if (particle->getMCParticle())
1968 if (particle->getMCParticle()->getMother() ==
nullptr) {
1971 Particle tempParticle = Particle(particle->getMCParticle()->getMother());
1972 auto var_result = var->function(&tempParticle);
1973 if (std::holds_alternative<double>(var_result)) {
1974 return std::get<double>(var_result);
1975 }
else if (std::holds_alternative<int>(var_result)) {
1976 return std::get<int>(var_result);
1977 }
else if (std::holds_alternative<bool>(var_result)) {
1978 return std::get<bool>(var_result);
1987 B2FATAL(
"Wrong number of arguments for meta function mcMother");
1993 if (arguments.size() == 2) {
1994 int particleNumber = 0;
1996 particleNumber = Belle2::convertString<int>(arguments[0]);
1997 }
catch (std::invalid_argument&) {
1998 B2FATAL(
"First argument of genParticle meta function must be integer!");
2002 auto func = [var, particleNumber](
const Particle*) ->
double {
2003 StoreArray<MCParticle> mcParticles(
"MCParticles");
2004 if (particleNumber >= mcParticles.getEntries())
2009 MCParticle* mcParticle = mcParticles[particleNumber];
2010 Particle part = Particle(mcParticle);
2011 auto var_result = var->function(&part);
2012 if (std::holds_alternative<double>(var_result))
2014 return std::get<double>(var_result);
2015 }
else if (std::holds_alternative<int>(var_result))
2017 return std::get<int>(var_result);
2018 }
else if (std::holds_alternative<bool>(var_result))
2020 return std::get<bool>(var_result);
2025 B2FATAL(
"Wrong number of arguments for meta function genParticle");
2031 if (arguments.size() == 1) {
2034 auto func = [var](
const Particle*) ->
double {
2035 StoreArray<MCParticle> mcParticles(
"MCParticles");
2036 if (mcParticles.getEntries() == 0)
2041 MCParticle* mcUpsilon4S = mcParticles[0];
2042 if (mcUpsilon4S->isInitial()) mcUpsilon4S = mcParticles[2];
2043 if (mcUpsilon4S->getPDG() != 300553)
2048 Particle upsilon4S = Particle(mcUpsilon4S);
2049 auto var_result = var->function(&upsilon4S);
2050 if (std::holds_alternative<double>(var_result))
2052 return std::get<double>(var_result);
2053 }
else if (std::holds_alternative<int>(var_result))
2055 return std::get<int>(var_result);
2056 }
else if (std::holds_alternative<bool>(var_result))
2058 return std::get<bool>(var_result);
2063 B2FATAL(
"Wrong number of arguments for meta function genUpsilon4S");
2069 if (arguments.size() == 4) {
2070 std::string listName = arguments[0];
2071 std::string rankedVariableName = arguments[1];
2072 std::string returnVariableName = arguments[2];
2073 std::string extraInfoName = rankedVariableName +
"_rank";
2076 rank = Belle2::convertString<int>(arguments[3]);
2077 }
catch (std::invalid_argument&) {
2078 B2ERROR(
"3rd argument of getVariableByRank meta function (Rank) must be an integer!");
2083 auto func = [var, rank, extraInfoName, listName](
const Particle*)->
double {
2084 StoreObjPtr<ParticleList> list(listName);
2086 const unsigned int numParticles = list->getListSize();
2087 for (
unsigned int i = 0; i < numParticles; i++)
2089 const Particle* p = list->getParticle(i);
2090 if (p->getExtraInfo(extraInfoName) == rank) {
2091 auto var_result = var->function(p);
2092 if (std::holds_alternative<double>(var_result)) {
2093 return std::get<double>(var_result);
2094 }
else if (std::holds_alternative<int>(var_result)) {
2095 return std::get<int>(var_result);
2096 }
else if (std::holds_alternative<bool>(var_result)) {
2097 return std::get<bool>(var_result);
2102 return std::numeric_limits<double>::signaling_NaN();
2106 B2FATAL(
"Wrong number of arguments for meta function getVariableByRank");
2112 if (arguments.size() == 1 or arguments.size() == 2) {
2114 std::string listName = arguments[0];
2115 std::string cutString =
"";
2117 if (arguments.size() == 2) {
2118 cutString = arguments[1];
2123 auto func = [listName, cut](
const Particle*) ->
int {
2125 StoreObjPtr<ParticleList> list(listName);
2127 for (
unsigned int i = 0; i < list->getListSize(); i++)
2129 const Particle* particle = list->getParticle(i);
2130 if (cut->check(particle)) {
2138 B2FATAL(
"Wrong number of arguments for meta function countInList");
2144 if (arguments.size() == 2 or arguments.size() == 3) {
2146 std::string roeListName = arguments[0];
2147 std::string cutString = arguments[1];
2149 if (arguments.size() == 2) {
2150 B2INFO(
"Use pdgCode of electron as default in meta variable veto, other arguments: " << roeListName <<
", " << cutString);
2153 pdgCode = Belle2::convertString<int>(arguments[2]);;
2154 }
catch (std::invalid_argument&) {
2155 B2FATAL(
"Third argument of veto meta function must be integer!");
2162 auto func = [roeListName, cut, pdgCode, flavourType](
const Particle * particle) ->
bool {
2163 StoreObjPtr<ParticleList> roeList(roeListName);
2164 ROOT::Math::PxPyPzEVector vec = particle->get4Vector();
2165 for (
unsigned int i = 0; i < roeList->getListSize(); i++)
2167 const Particle* roeParticle = roeList->getParticle(i);
2168 if (not particle->overlapsWith(roeParticle)) {
2169 ROOT::Math::PxPyPzEVector tempCombination = roeParticle->get4Vector() + vec;
2170 std::vector<int> indices = { particle->getArrayIndex(), roeParticle->getArrayIndex() };
2171 Particle tempParticle = Particle(tempCombination, pdgCode, flavourType, indices, particle->getArrayPointer());
2172 if (cut->check(&tempParticle)) {
2181 B2FATAL(
"Wrong number of arguments for meta function veto");
2187 if (arguments.size() == 1) {
2188 std::string cutString = arguments[0];
2190 auto func = [cut](
const Particle * particle) ->
int {
2192 for (
auto& daughter : particle->getDaughters())
2194 if (cut->check(daughter))
2201 B2FATAL(
"Wrong number of arguments for meta function countDaughters");
2207 if (arguments.size() == 1) {
2208 std::string cutString = arguments[0];
2210 auto func = [cut](
const Particle * particle) ->
int {
2212 std::vector<const Particle*> fspDaughters;
2213 particle->fillFSPDaughters(fspDaughters);
2216 for (
auto& daughter : fspDaughters)
2218 if (cut->check(daughter))
2225 B2FATAL(
"Wrong number of arguments for meta function countFSPDaughters");
2231 if (arguments.size() == 1) {
2232 std::string cutString = arguments[0];
2234 auto func = [cut](
const Particle * particle) ->
int {
2236 std::vector<const Particle*> allDaughters;
2237 particle->fillAllDaughters(allDaughters);
2240 for (
auto& daughter : allDaughters)
2242 if (cut->check(daughter))
2249 B2FATAL(
"Wrong number of arguments for meta function countDescendants");
2253 Manager::FunctionPtr numberOfNonOverlappingParticles(
const std::vector<std::string>& arguments)
2256 auto func = [arguments](
const Particle * particle) ->
int {
2258 int _numberOfNonOverlappingParticles = 0;
2259 for (
const auto& listName : arguments)
2261 StoreObjPtr<ParticleList> list(listName);
2262 if (not list.isValid()) {
2263 B2FATAL(
"Invalid list named " << listName <<
" encountered in numberOfNonOverlappingParticles.");
2265 for (
unsigned int i = 0; i < list->getListSize(); i++) {
2266 const Particle* p = list->getParticle(i);
2267 if (not particle->overlapsWith(p)) {
2268 _numberOfNonOverlappingParticles++;
2272 return _numberOfNonOverlappingParticles;
2279 void appendDaughtersRecursive(Particle* mother)
2282 auto* mcmother = mother->getRelated<MCParticle>();
2287 std::vector<MCParticle*> mcdaughters = mcmother->
getDaughters();
2288 StoreArray<Particle> particles;
2290 for (
auto& mcdaughter : mcdaughters) {
2292 Particle tmp_daughter(mcdaughter);
2293 Particle* new_daughter = particles.appendNew(tmp_daughter);
2294 new_daughter->addRelationTo(mcdaughter);
2295 mother->appendDaughter(new_daughter,
false);
2297 if (mcdaughter->getNDaughters() > 0)
2298 appendDaughtersRecursive(new_daughter);
2305 if (arguments.size() == 1) {
2307 auto func = [var](
const Particle * particle) ->
double {
2308 const MCParticle* mcp = particle->getMCParticle();
2313 Particle tmpPart(mcp);
2314 StoreArray<Particle> particles;
2315 Particle* newPart = particles.appendNew(tmpPart);
2316 newPart->addRelationTo(mcp);
2318 appendDaughtersRecursive(newPart);
2320 auto var_result = var->function(newPart);
2321 if (std::holds_alternative<double>(var_result))
2323 return std::get<double>(var_result);
2324 }
else if (std::holds_alternative<int>(var_result))
2326 return std::get<int>(var_result);
2327 }
else if (std::holds_alternative<bool>(var_result))
2329 return std::get<bool>(var_result);
2334 B2FATAL(
"Wrong number of arguments for meta function matchedMC");
2340 if (arguments.size() == 1) {
2343 auto func = [var](
const Particle * particle) ->
double {
2345 const ECLCluster* cluster = particle->getECLCluster();
2348 auto mcps = cluster->getRelationsTo<MCParticle>();
2351 std::vector<std::pair<double, int>> weightsAndIndices;
2352 for (
unsigned int i = 0; i < mcps.size(); ++i)
2353 weightsAndIndices.emplace_back(mcps.weight(i), i);
2356 std::sort(weightsAndIndices.begin(), weightsAndIndices.end(),
2357 ValueIndexPairSorting::higherPair<
decltype(weightsAndIndices)::value_type>);
2360 const MCParticle* mcp = mcps.object(weightsAndIndices[0].second);
2362 Particle tmpPart(mcp);
2363 StoreArray<Particle> particles;
2364 Particle* newPart = particles.appendNew(tmpPart);
2365 newPart->addRelationTo(mcp);
2367 appendDaughtersRecursive(newPart);
2369 auto var_result = var->function(newPart);
2370 if (std::holds_alternative<double>(var_result))
2372 return std::get<double>(var_result);
2373 }
else if (std::holds_alternative<int>(var_result))
2375 return std::get<int>(var_result);
2376 }
else if (std::holds_alternative<bool>(var_result))
2378 return std::get<bool>(var_result);
2387 B2FATAL(
"Wrong number of arguments for meta function clusterBestMatchedMCParticle");
2393 if (arguments.size() == 1) {
2396 auto func = [var](
const Particle * particle) ->
double {
2398 const ECLCluster* cluster = particle->getECLCluster();
2401 auto mcps = cluster->getRelationsTo<MCParticle>();
2404 std::map<int, double> mapMCParticleIndxAndWeight;
2405 getKlongWeightMap(particle, mapMCParticleIndxAndWeight);
2408 if (mapMCParticleIndxAndWeight.size() == 0)
2412 auto maxMap = std::max_element(mapMCParticleIndxAndWeight.begin(), mapMCParticleIndxAndWeight.end(),
2413 [](
const auto & x,
const auto & y) { return x.second < y.second; }
2416 StoreArray<MCParticle> mcparticles;
2417 const MCParticle* mcKlong = mcparticles[maxMap->first];
2419 Particle tmpPart(mcKlong);
2420 auto var_result = var->function(&tmpPart);
2421 if (std::holds_alternative<double>(var_result))
2423 return std::get<double>(var_result);
2424 }
else if (std::holds_alternative<int>(var_result))
2426 return std::get<int>(var_result);
2427 }
else if (std::holds_alternative<bool>(var_result))
2429 return std::get<bool>(var_result);
2438 B2FATAL(
"Wrong number of arguments for meta function clusterBestMatchedMCKlong");
2442 double matchedMCHasPDG(
const Particle* particle,
const std::vector<double>& pdgCode)
2444 if (pdgCode.size() != 1) {
2445 B2FATAL(
"Too many arguments provided to matchedMCHasPDG!");
2447 int inputPDG = std::lround(pdgCode[0]);
2449 const MCParticle* mcp = particle->getMCParticle();
2453 return std::abs(mcp->getPDG()) == inputPDG;
2458 if (arguments.size() == 1) {
2459 std::string listName = arguments[0];
2460 auto func = [listName](
const Particle * particle) ->
double {
2463 StoreObjPtr<ParticleList> listOfParticles(listName);
2465 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalEnergyOfParticlesInList");
2466 double totalEnergy = 0;
2467 int nParticles = listOfParticles->getListSize();
2468 for (
int i = 0; i < nParticles; i++)
2470 const Particle* part = listOfParticles->getParticle(i);
2472 totalEnergy += frame.getMomentum(part).E();
2479 B2FATAL(
"Wrong number of arguments for meta function totalEnergyOfParticlesInList");
2485 if (arguments.size() == 1) {
2486 std::string listName = arguments[0];
2487 auto func = [listName](
const Particle*) ->
double {
2488 StoreObjPtr<ParticleList> listOfParticles(listName);
2490 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPxOfParticlesInList");
2492 int nParticles = listOfParticles->getListSize();
2494 for (
int i = 0; i < nParticles; i++)
2496 const Particle* part = listOfParticles->getParticle(i);
2497 totalPx += frame.getMomentum(part).Px();
2503 B2FATAL(
"Wrong number of arguments for meta function totalPxOfParticlesInList");
2509 if (arguments.size() == 1) {
2510 std::string listName = arguments[0];
2511 auto func = [listName](
const Particle*) ->
double {
2512 StoreObjPtr<ParticleList> listOfParticles(listName);
2514 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPyOfParticlesInList");
2516 int nParticles = listOfParticles->getListSize();
2518 for (
int i = 0; i < nParticles; i++)
2520 const Particle* part = listOfParticles->getParticle(i);
2521 totalPy += frame.getMomentum(part).Py();
2527 B2FATAL(
"Wrong number of arguments for meta function totalPyOfParticlesInList");
2533 if (arguments.size() == 1) {
2534 std::string listName = arguments[0];
2535 auto func = [listName](
const Particle*) ->
double {
2536 StoreObjPtr<ParticleList> listOfParticles(listName);
2538 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPzOfParticlesInList");
2540 int nParticles = listOfParticles->getListSize();
2542 for (
int i = 0; i < nParticles; i++)
2544 const Particle* part = listOfParticles->getParticle(i);
2545 totalPz += frame.getMomentum(part).Pz();
2551 B2FATAL(
"Wrong number of arguments for meta function totalPzOfParticlesInList");
2557 if (arguments.size() > 0) {
2559 auto func = [arguments](
const Particle * particle) ->
double {
2561 ROOT::Math::PxPyPzEVector total4Vector;
2563 std::vector<Particle*> particlePool;
2566 for (
const auto& argument : arguments)
2568 StoreObjPtr <ParticleList> listOfParticles(argument);
2570 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << argument <<
" given to invMassInLists");
2571 int nParticles = listOfParticles->getListSize();
2572 for (
int i = 0; i < nParticles; i++) {
2573 bool overlaps =
false;
2574 Particle* part = listOfParticles->getParticle(i);
2575 for (
auto poolPart : particlePool) {
2576 if (part->overlapsWith(poolPart)) {
2582 total4Vector += part->get4Vector();
2583 particlePool.push_back(part);
2587 double invariantMass = total4Vector.M();
2588 return invariantMass;
2593 B2FATAL(
"Wrong number of arguments for meta function invMassInLists");
2597 Manager::FunctionPtr totalECLEnergyOfParticlesInList(
const std::vector<std::string>& arguments)
2599 if (arguments.size() == 1) {
2600 std::string listName = arguments[0];
2601 auto func = [listName](
const Particle * particle) ->
double {
2604 StoreObjPtr<ParticleList> listOfParticles(listName);
2606 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalEnergyOfParticlesInList");
2607 double totalEnergy = 0;
2608 int nParticles = listOfParticles->getListSize();
2609 for (
int i = 0; i < nParticles; i++)
2611 const Particle* part = listOfParticles->getParticle(i);
2612 const ECLCluster* cluster = part->getECLCluster();
2614 if (cluster !=
nullptr) {
2615 totalEnergy += cluster->getEnergy(clusterHypothesis);
2623 B2FATAL(
"Wrong number of arguments for meta function totalECLEnergyOfParticlesInList");
2629 if (arguments.size() == 1) {
2630 std::string listName = arguments[0];
2631 auto func = [listName](
const Particle*) ->
double {
2632 StoreObjPtr<ParticleList> listOfParticles(listName);
2634 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to maxPtInList");
2635 int nParticles = listOfParticles->getListSize();
2638 for (
int i = 0; i < nParticles; i++)
2640 const Particle* part = listOfParticles->getParticle(i);
2641 const double Pt = frame.getMomentum(part).Pt();
2642 if (Pt > maxPt) maxPt = Pt;
2648 B2FATAL(
"Wrong number of arguments for meta function maxPtInList");
2652 Manager::FunctionPtr eclClusterTrackMatchedWithCondition(
const std::vector<std::string>& arguments)
2654 if (arguments.size() <= 1) {
2656 std::string cutString;
2657 if (arguments.size() == 1)
2658 cutString = arguments[0];
2660 auto func = [cut](
const Particle * particle) ->
double {
2662 if (particle ==
nullptr)
2665 const ECLCluster* cluster = particle->getECLCluster();
2669 auto tracks = cluster->getRelationsFrom<Track>();
2671 for (
const auto& track : tracks) {
2674 if (cut->check(&trackParticle))
2683 B2FATAL(
"Wrong number of arguments for meta function eclClusterSpecialTrackMatched");
2689 if (arguments.size() == 2) {
2690 std::string listName = arguments[0];
2693 auto func = [listName, var](
const Particle*) ->
double {
2694 StoreObjPtr<ParticleList> listOfParticles(listName);
2696 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to averageValueInList");
2697 int nParticles = listOfParticles->getListSize();
2698 if (nParticles == 0)
2703 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2705 for (
int i = 0; i < nParticles; i++) {
2706 average += std::get<double>(var->function(listOfParticles->getParticle(i))) / nParticles;
2708 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2710 for (
int i = 0; i < nParticles; i++) {
2711 average += std::get<int>(var->function(listOfParticles->getParticle(i))) / nParticles;
2718 B2FATAL(
"Wrong number of arguments for meta function averageValueInList");
2724 if (arguments.size() == 2) {
2725 std::string listName = arguments[0];
2728 auto func = [listName, var](
const Particle*) ->
double {
2729 StoreObjPtr<ParticleList> listOfParticles(listName);
2731 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to medianValueInList");
2732 int nParticles = listOfParticles->getListSize();
2733 if (nParticles == 0)
2737 std::vector<double> valuesInList;
2738 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2740 for (
int i = 0; i < nParticles; i++) {
2741 valuesInList.push_back(std::get<double>(var->function(listOfParticles->getParticle(i))));
2743 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2745 for (
int i = 0; i < nParticles; i++) {
2746 valuesInList.push_back(std::get<int>(var->function(listOfParticles->getParticle(i))));
2749 std::sort(valuesInList.begin(), valuesInList.end());
2750 if (nParticles % 2 != 0)
2752 return valuesInList[nParticles / 2];
2755 return 0.5 * (valuesInList[nParticles / 2] + valuesInList[nParticles / 2 - 1]);
2760 B2FATAL(
"Wrong number of arguments for meta function medianValueInList");
2766 if (arguments.size() == 2) {
2767 std::string listName = arguments[0];
2770 auto func = [listName, var](
const Particle*) ->
double {
2771 StoreObjPtr<ParticleList> listOfParticles(listName);
2773 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to sumValueInList");
2774 int nParticles = listOfParticles->getListSize();
2775 if (nParticles == 0)
2780 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2782 for (
int i = 0; i < nParticles; i++) {
2783 sum += std::get<double>(var->function(listOfParticles->getParticle(i)));
2785 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2787 for (
int i = 0; i < nParticles; i++) {
2788 sum += std::get<int>(var->function(listOfParticles->getParticle(i)));
2795 B2FATAL(
"Wrong number of arguments for meta function sumValueInList");
2801 if (arguments.size() == 2) {
2802 std::string listName = arguments[0];
2805 auto func = [listName, var](
const Particle*) ->
double {
2806 StoreObjPtr<ParticleList> listOfParticles(listName);
2808 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to productValueInList");
2809 int nParticles = listOfParticles->getListSize();
2810 if (nParticles == 0)
2815 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2817 for (
int i = 0; i < nParticles; i++) {
2818 product *= std::get<double>(var->function(listOfParticles->getParticle(i)));
2820 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2822 for (
int i = 0; i < nParticles; i++) {
2823 product *= std::get<int>(var->function(listOfParticles->getParticle(i)));
2830 B2FATAL(
"Wrong number of arguments for meta function productValueInList");
2837 if (arguments.size() != 1)
2838 B2FATAL(
"Wrong number of arguments for meta function angleToClosestInList");
2840 std::string listname = arguments[0];
2842 auto func = [listname](
const Particle * particle) ->
double {
2844 StoreObjPtr<ParticleList> list(listname);
2845 if (not list.isValid())
2846 B2FATAL(
"Invalid particle list name " << listname <<
" given to angleToClosestInList");
2849 if (list->getListSize() == 0)
2854 const auto p_this =
B2Vector3D(frame.getMomentum(particle).Vect());
2857 double minAngle = 2 * M_PI;
2858 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2860 const Particle* compareme = list->getParticle(i);
2861 const auto p_compare =
B2Vector3D(frame.getMomentum(compareme).Vect());
2862 double angle = p_compare.Angle(p_this);
2863 if (minAngle > angle) minAngle = angle;
2873 if (arguments.size() != 2)
2874 B2FATAL(
"Wrong number of arguments for meta function closestInList");
2876 std::string listname = arguments[0];
2881 auto func = [listname, var](
const Particle * particle) ->
double {
2883 StoreObjPtr<ParticleList> list(listname);
2884 if (not list.isValid())
2885 B2FATAL(
"Invalid particle list name " << listname <<
" given to closestInList");
2889 const auto p_this =
B2Vector3D(frame.getMomentum(particle).Vect());
2892 double minAngle = 2 * M_PI;
2894 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2896 const Particle* compareme = list->getParticle(i);
2897 const auto p_compare =
B2Vector3D(frame.getMomentum(compareme).Vect());
2898 double angle = p_compare.Angle(p_this);
2899 if (minAngle > angle) {
2907 auto var_result = var->function(list->getParticle(iClosest));
2908 if (std::holds_alternative<double>(var_result))
2910 return std::get<double>(var_result);
2911 }
else if (std::holds_alternative<int>(var_result))
2913 return std::get<int>(var_result);
2914 }
else if (std::holds_alternative<bool>(var_result))
2916 return std::get<bool>(var_result);
2925 if (arguments.size() != 1)
2926 B2FATAL(
"Wrong number of arguments for meta function angleToMostB2BInList");
2928 std::string listname = arguments[0];
2930 auto func = [listname](
const Particle * particle) ->
double {
2932 StoreObjPtr<ParticleList> list(listname);
2933 if (not list.isValid())
2934 B2FATAL(
"Invalid particle list name " << listname <<
" given to angleToMostB2BInList");
2937 if (list->getListSize() == 0)
2942 const auto p_this =
B2Vector3D(frame.getMomentum(particle).Vect());
2946 double maxAngle = 0;
2947 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2949 const Particle* compareme = list->getParticle(i);
2950 const auto p_compare =
B2Vector3D(frame.getMomentum(compareme).Vect());
2951 double angle = p_compare.Angle(p_this);
2952 if (maxAngle < angle) maxAngle = angle;
2962 if (arguments.size() != 1)
2963 B2FATAL(
"Wrong number of arguments for meta function deltaPhiToMostB2BPhiInList");
2965 std::string listname = arguments[0];
2967 auto func = [listname](
const Particle * particle) ->
double {
2969 StoreObjPtr<ParticleList> list(listname);
2970 if (not list.isValid())
2971 B2FATAL(
"Invalid particle list name " << listname <<
" given to deltaPhiToMostB2BPhiInList");
2974 if (list->getListSize() == 0)
2979 const auto phi_this = frame.getMomentum(particle).Phi();
2982 double maxAngle = 0;
2983 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2985 const Particle* compareme = list->getParticle(i);
2986 const auto phi_compare = frame.getMomentum(compareme).Phi();
2987 double angle = std::abs(phi_compare - phi_this);
2988 if (angle > M_PI) {angle = 2 * M_PI - angle;}
2989 if (maxAngle < angle) maxAngle = angle;
2999 if (arguments.size() != 2)
3000 B2FATAL(
"Wrong number of arguments for meta function mostB2BInList");
3002 std::string listname = arguments[0];
3007 auto func = [listname, var](
const Particle * particle) ->
double {
3009 StoreObjPtr<ParticleList> list(listname);
3010 if (not list.isValid())
3011 B2FATAL(
"Invalid particle list name " << listname <<
" given to mostB2BInList");
3015 const auto p_this =
B2Vector3D(frame.getMomentum(particle).Vect());
3019 double maxAngle = -1.0;
3021 for (
unsigned int i = 0; i < list->getListSize(); ++i)
3023 const Particle* compareme = list->getParticle(i);
3024 const auto p_compare =
B2Vector3D(frame.getMomentum(compareme).Vect());
3025 double angle = p_compare.Angle(p_this);
3026 if (maxAngle < angle) {
3034 auto var_result = var->function(list->getParticle(iMostB2B));
3035 if (std::holds_alternative<double>(var_result))
3037 return std::get<double>(var_result);
3038 }
else if (std::holds_alternative<int>(var_result))
3040 return std::get<int>(var_result);
3041 }
else if (std::holds_alternative<bool>(var_result))
3043 return std::get<bool>(var_result);
3051 if (arguments.size() == 1) {
3052 std::string listName = arguments[0];
3053 auto func = [listName](
const Particle*) ->
double {
3054 StoreObjPtr<ParticleList> listOfParticles(listName);
3056 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to maxOpeningAngleInList");
3057 int nParticles = listOfParticles->getListSize();
3062 double maxOpeningAngle = -1;
3063 for (
int i = 0; i < nParticles; i++)
3065 B2Vector3D v1 = frame.getMomentum(listOfParticles->getParticle(i)).Vect();
3066 for (
int j = i + 1; j < nParticles; j++) {
3067 B2Vector3D v2 = frame.getMomentum(listOfParticles->getParticle(j)).Vect();
3068 const double angle =
v1.Angle(v2);
3069 if (angle > maxOpeningAngle) maxOpeningAngle = angle;
3072 return maxOpeningAngle;
3076 B2FATAL(
"Wrong number of arguments for meta function maxOpeningAngleInList");
3083 if (arguments.size() >= 2) {
3088 auto func = [var, arguments](
const Particle * particle) ->
double {
3089 if (particle ==
nullptr)
3091 B2WARNING(
"Trying to access a daughter that does not exist. Skipping");
3097 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
3101 for (
unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++)
3103 auto generalizedIndex = arguments[iCoord];
3104 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
3106 pSum += frame.getMomentum(dauPart);
3108 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
3114 Particle sumOfDaughters(pSum, 100);
3116 auto var_result = var->function(&sumOfDaughters);
3118 if (std::holds_alternative<double>(var_result))
3120 return std::get<double>(var_result);
3121 }
else if (std::holds_alternative<int>(var_result))
3123 return std::get<int>(var_result);
3124 }
else if (std::holds_alternative<bool>(var_result))
3126 return std::get<bool>(var_result);
3131 B2FATAL(
"Wrong number of arguments for meta function daughterCombination");
3134 Manager::FunctionPtr useAlternativeDaughterHypothesis(
const std::vector<std::string>& arguments)
3147 if (arguments.size() >= 2) {
3158 std::unordered_map<unsigned int, int> mapOfReplacedDaughters;
3161 for (
unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++) {
3162 auto replacedDauString = arguments[iCoord];
3164 std::vector<std::string> indexAndMass;
3165 boost::split(indexAndMass, replacedDauString, boost::is_any_of(
":"));
3168 if (indexAndMass.size() > 2) {
3169 B2WARNING(
"The string indicating which daughter's mass should be replaced contains more than two elements separated by a colon. Perhaps you tried to pass a generalized index, which is not supported yet for this variable. The offending string is "
3170 << replacedDauString <<
", while a correct syntax looks like 0:K+.");
3174 if (indexAndMass.size() < 2) {
3175 B2WARNING(
"The string indicating which daughter's mass should be replaced contains only one colon-separated element instead of two. The offending string is "
3176 << replacedDauString <<
", while a correct syntax looks like 0:K+.");
3183 dauIndex = Belle2::convertString<int>(indexAndMass[0]);
3184 }
catch (std::invalid_argument&) {
3185 B2FATAL(
"Found the string " << indexAndMass[0] <<
"instead of a daughter index.");
3189 TParticlePDG* particlePDG = TDatabasePDG::Instance()->GetParticle(indexAndMass[1].c_str());
3191 B2WARNING(
"Particle not in evt.pdl file! " << indexAndMass[1]);
3196 int pdgCode = particlePDG->PdgCode();
3197 mapOfReplacedDaughters[dauIndex] = pdgCode;
3201 if (mapOfReplacedDaughters.size() != arguments.size() - 1)
3202 B2FATAL(
"Overlapped daughter's index is detected in the meta-variable useAlternativeDaughterHypothesis");
3210 auto func = [var, mapOfReplacedDaughters](
const Particle * particle) ->
double {
3211 if (particle ==
nullptr)
3213 B2WARNING(
"Trying to access a particle that does not exist. Skipping");
3223 ROOT::Math::PxPyPzMVector pSum(0, 0, 0, 0);
3225 for (
unsigned int iDau = 0; iDau < particle->getNDaughters(); iDau++)
3227 const Particle* dauPart = particle->getDaughter(iDau);
3229 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << iDau);
3233 ROOT::Math::PxPyPzMVector dauMom = ROOT::Math::PxPyPzMVector(frame.getMomentum(dauPart));
3237 pdgCode = mapOfReplacedDaughters.at(iDau);
3238 }
catch (std::out_of_range&) {
3245 double p_x = dauMom.Px();
3246 double p_y = dauMom.Py();
3247 double p_z = dauMom.Pz();
3248 dauMom.SetCoordinates(p_x, p_y, p_z, TDatabasePDG::Instance()->GetParticle(pdgCode)->Mass());
3249 const_cast<Particle*
>(dummy->getDaughter(iDau))->set4VectorDividingByMomentumScaling(ROOT::Math::PxPyPzEVector(dauMom));
3252 const int charge = dummy->getDaughter(iDau)->getCharge();
3253 if (TDatabasePDG::Instance()->GetParticle(pdgCode)->Charge() / 3.0 ==
charge)
3254 const_cast<Particle*
>(dummy->getDaughter(iDau))->setPDGCode(pdgCode);
3256 const_cast<Particle*
>(dummy->getDaughter(iDau))->setPDGCode(-1 * pdgCode);
3262 dummy->set4Vector(ROOT::Math::PxPyPzEVector(pSum));
3264 auto var_result = var->function(dummy);
3267 if (std::holds_alternative<double>(var_result))
3269 return std::get<double>(var_result);
3270 }
else if (std::holds_alternative<int>(var_result))
3272 return std::get<int>(var_result);
3273 }
else if (std::holds_alternative<bool>(var_result))
3275 return std::get<bool>(var_result);
3281 B2FATAL(
"Wrong number of arguments for meta function useAlternativeDaughterHypothesis");
3286 if (arguments.size() == 2) {
3288 std::string arg = arguments[0];
3290 TParticlePDG* part = TDatabasePDG::Instance()->GetParticle(arg.c_str());
3292 if (part !=
nullptr) {
3293 pdg_code = std::abs(part->PdgCode());
3296 pdg_code = Belle2::convertString<int>(arg);
3297 }
catch (std::exception& e) {}
3300 if (pdg_code == -1) {
3301 B2FATAL(
"Ancestor " + arg +
" is not recognised. Please provide valid PDG code or particle name.");
3304 auto func = [pdg_code, var](
const Particle * particle) ->
double {
3305 const Particle* p = particle;
3307 int ancestor_level = std::get<double>(
Manager::Instance().getVariable(
"hasAncestor(" + std::to_string(pdg_code) +
", 0)")->function(p));
3308 if ((ancestor_level <= 0) or (std::isnan(ancestor_level)))
3313 const MCParticle* i_p = p->getMCParticle();
3315 for (
int a = 0; a < ancestor_level ; a = a + 1)
3317 i_p = i_p->getMother();
3321 StoreArray<Particle> particles;
3322 Particle* newPart = particles.appendNew(m_p);
3323 newPart->addRelationTo(i_p);
3325 appendDaughtersRecursive(newPart);
3327 auto var_result = var->function(newPart);
3328 if (std::holds_alternative<double>(var_result))
3330 return std::get<double>(var_result);
3331 }
else if (std::holds_alternative<int>(var_result))
3333 return std::get<int>(var_result);
3334 }
else if (std::holds_alternative<bool>(var_result))
3336 return std::get<bool>(var_result);
3341 B2FATAL(
"Wrong number of arguments for meta function varForFirstMCAncestorOfType (expected 2: type and variable of interest)");
3347 if (arguments.size() != 1) {
3348 B2FATAL(
"Number of arguments for nTrackFitResults must be 1, particleType or PDGcode");
3351 std::string arg = arguments[0];
3352 TParticlePDG* part = TDatabasePDG::Instance()->GetParticle(arg.c_str());
3354 if (part !=
nullptr) {
3355 absPdg = std::abs(part->PdgCode());
3358 absPdg = Belle2::convertString<int>(arg);
3359 }
catch (std::exception& e) {}
3362 auto func = [absPdg](
const Particle*) ->
int {
3364 Const::ChargedStable type(absPdg);
3365 StoreArray<Track> tracks;
3367 int nTrackFitResults = 0;
3369 for (
int i = 0; i < tracks.getEntries(); i++)
3371 const Track* track = tracks[i];
3372 const TrackFitResult* trackFit = track->getTrackFitResultWithClosestMass(type);
3374 if (!trackFit)
continue;
3375 if (trackFit->getChargeSign() == 0)
continue;
3380 return nTrackFitResults;
3386 VARIABLE_GROUP(
"MetaFunctions");
3387 REGISTER_METAVARIABLE(
"nCleanedECLClusters(cut)", nCleanedECLClusters,
3388 "[Eventbased] Returns the number of clean Clusters in the event\n"
3389 "Clean clusters are defined by the clusters which pass the given cut assuming a photon hypothesis.",
3390 Manager::VariableDataType::c_int);
3391 REGISTER_METAVARIABLE(
"nCleanedTracks(cut)", nCleanedTracks,
3392 "[Eventbased] Returns the number of clean Tracks in the event\n"
3393 "Clean tracks are defined by the tracks which pass the given cut assuming a pion hypothesis.", Manager::VariableDataType::c_int);
3394 REGISTER_METAVARIABLE(
"formula(v1 + v2 * [v3 - v4] / v5^v6)", formula, R
"DOCSTRING(
3395Returns the result of the given formula, where v1 to vN are variables or floating
3396point numbers. Currently the only supported operations are addition (``+``),
3397subtraction (``-``), multiplication (``*``), division (``/``) and power (``^``
3398or ``**``). Parenthesis can be in the form of square brackets ``[v1 * v2]``
3399or normal brackets ``(v1 * v2)``. It will work also with variables taking
3400arguments. Operator precedence is taken into account. For example ::
3402 (daughter(0, E) + daughter(1, E))**2 - p**2 + 0.138
3404.. versionchanged:: release-03-00-00
3405 now both, ``[]`` and ``()`` can be used for grouping operations, ``**`` can
3406 be used for exponent and float literals are possible directly in the
3408)DOCSTRING", Manager::VariableDataType::c_double);
3409 REGISTER_METAVARIABLE("useRestFrame(variable)", useRestFrame,
3410 "Returns the value of the variable using the rest frame of the given particle as current reference frame.\n"
3411 "E.g. ``useRestFrame(daughter(0, p))`` returns the total momentum of the first daughter in its mother's rest-frame", Manager::VariableDataType::c_double);
3412 REGISTER_METAVARIABLE(
"useCMSFrame(variable)", useCMSFrame,
3413 "Returns the value of the variable using the CMS frame as current reference frame.\n"
3414 "E.g. ``useCMSFrame(E)`` returns the energy of a particle in the CMS frame.", Manager::VariableDataType::c_double);
3415 REGISTER_METAVARIABLE(
"useLabFrame(variable)", useLabFrame, R
"DOC(
3416Returns the value of ``variable`` in the *lab* frame.
3419 The lab frame is the default reference frame, usually you don't need to use this meta-variable.
3420 E.g. ``useLabFrame(E)`` returns the energy of a particle in the Lab frame, same as just ``E``.
3422Specifying the lab frame is useful in some corner-cases. For example:
3423``useRestFrame(daughter(0, formula(E - useLabFrame(E))))`` which is the difference of the first daughter's energy in the rest frame of the mother (current particle) with the same daughter's lab-frame energy.
3424)DOC", Manager::VariableDataType::c_double);
3425 REGISTER_METAVARIABLE("useTagSideRecoilRestFrame(variable, daughterIndexTagB)", useTagSideRecoilRestFrame,
3426 "Returns the value of the variable in the rest frame of the recoiling particle to the tag side B meson.\n"
3427 "The variable should only be applied to an Upsilon(4S) list.\n"
3428 "E.g. ``useTagSideRecoilRestFrame(daughter(1, daughter(1, p)), 0)`` applied on a Upsilon(4S) list (``Upsilon(4S)->B+:tag B-:sig``) returns the momentum of the second daughter of the signal B meson in the signal B meson rest frame.", Manager::VariableDataType::c_double);
3429 REGISTER_METAVARIABLE(
"useParticleRestFrame(variable, particleList)", useParticleRestFrame,
3430 "Returns the value of the variable in the rest frame of the first Particle contained in the given ParticleList.\n"
3431 "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3432 "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3433 "computing the rest frame and a warning is thrown. If the given ParticleList is empty in an event, it returns NaN.", Manager::VariableDataType::c_double);
3434 REGISTER_METAVARIABLE(
"useRecoilParticleRestFrame(variable, particleList)", useRecoilParticleRestFrame,
3435 "Returns the value of the variable in the rest frame of recoil system against the first Particle contained in the given ParticleList.\n"
3436 "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3437 "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3438 "computing the rest frame and a warning is thrown. If the given ParticleList is empty in an event, it returns NaN.", Manager::VariableDataType::c_double);
3439 REGISTER_METAVARIABLE(
"useDaughterRestFrame(variable, daughterIndex_1, [daughterIndex_2, ... daughterIndex_3])", useDaughterRestFrame,
3440 "Returns the value of the variable in the rest frame of the selected daughter particle.\n"
3441 "The daughter is identified via generalized daughter index, e.g. ``0:1`` identifies the second daughter (1) "
3442 "of the first daughter (0). If the daughter index is invalid, it returns NaN.\n"
3443 "If two or more indices are given, the rest frame of the sum of the daughters is used.",
3444 Manager::VariableDataType::c_double);
3445 REGISTER_METAVARIABLE(
"useMCancestorBRestFrame(variable)", useMCancestorBRestFrame,
3446 "Returns the value of the variable in the rest frame of the ancestor B MC particle.\n"
3447 "If no B or no MC-matching is found, it returns NaN.", Manager::VariableDataType::c_double);
3448 REGISTER_METAVARIABLE(
"passesCut(cut)", passesCut,
3449 "Returns 1 if particle passes the cut otherwise 0.\n"
3450 "Useful if you want to write out if a particle would have passed a cut or not.", Manager::VariableDataType::c_bool);
3451 REGISTER_METAVARIABLE(
"passesEventCut(cut)", passesEventCut,
3452 "[Eventbased] Returns 1 if event passes the cut otherwise 0.\n"
3453 "Useful if you want to select events passing a cut without looping into particles, such as for skimming.\n", Manager::VariableDataType::c_bool);
3454 REGISTER_METAVARIABLE(
"countDaughters(cut)", countDaughters,
3455 "Returns number of direct daughters which satisfy the cut.\n"
3456 "Used by the skimming package (for what exactly?)", Manager::VariableDataType::c_int);
3457 REGISTER_METAVARIABLE(
"countFSPDaughters(cut)", countDescendants,
3458 "Returns number of final-state daughters which satisfy the cut.",
3459 Manager::VariableDataType::c_int);
3460 REGISTER_METAVARIABLE(
"countDescendants(cut)", countDescendants,
3461 "Returns number of descendants for all generations which satisfy the cut.",
3462 Manager::VariableDataType::c_int);
3463 REGISTER_METAVARIABLE(
"varFor(pdgCode, variable)", varFor,
3464 "Returns the value of the variable for the given particle if its abs(pdgCode) agrees with the given one.\n"
3465 "E.g. ``varFor(11, p)`` returns the momentum if the particle is an electron or a positron.", Manager::VariableDataType::c_double);
3466 REGISTER_METAVARIABLE(
"varForMCGen(variable)", varForMCGen,
3467 "Returns the value of the variable for the given particle if the MC particle related to it is primary, not virtual, and not initial.\n"
3468 "If no MC particle is related to the given particle, or the MC particle is not primary, virtual, or initial, NaN will be returned.\n"
3469 "E.g. ``varForMCGen(PDG)`` returns the PDG code of the MC particle related to the given particle if it is primary, not virtual, and not initial.", Manager::VariableDataType::c_double);
3470 REGISTER_METAVARIABLE(
"nParticlesInList(particleListName)", nParticlesInList,
3471 "[Eventbased] Returns number of particles in the given particle List.", Manager::VariableDataType::c_int);
3472 REGISTER_METAVARIABLE(
"isInList(particleListName)", isInList,
3473 "Returns 1 if the particle is in the list provided, 0 if not. Note that this only checks the particle given. For daughters of composite particles, please see :b2:var:`isDaughterOfList`.", Manager::VariableDataType::c_bool);
3474 REGISTER_METAVARIABLE(
"isDaughterOfList(particleListNames)", isDaughterOfList,
3475 "Returns 1 if the given particle is a daughter of at least one of the particles in the given particle Lists.", Manager::VariableDataType::c_bool);
3476 REGISTER_METAVARIABLE(
"isDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isDescendantOfList, R
"DOC(
3477 Returns 1 if the given particle appears in the decay chain of the particles in the given ParticleLists.
3479 Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3481 * ``isDescendantOfList(<particle_list>,1)`` returns 1 if particle is a daughter of the list,
3482 * ``isDescendantOfList(<particle_list>,2)`` returns 1 if particle is a granddaughter of the list,
3483 * ``isDescendantOfList(<particle_list>,3)`` returns 1 if particle is a great-granddaughter of the list, etc.
3484 * Default value is ``-1`` that is inclusive for all generations.
3485 )DOC", Manager::VariableDataType::c_bool);
3486 REGISTER_METAVARIABLE("isMCDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isMCDescendantOfList, R
"DOC(
3487 Returns 1 if the given particle is linked to the same MC particle as any reconstructed daughter of the decay lists.
3489 Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3491 * ``isMCDescendantOfList(<particle_list>,1)`` returns 1 if particle is matched to the same particle as any daughter of the list,
3492 * ``isMCDescendantOfList(<particle_list>,2)`` returns 1 if particle is matched to the same particle as any granddaughter of the list,
3493 * ``isMCDescendantOfList(<particle_list>,3)`` returns 1 if particle is matched to the same particle as any great-granddaughter of the list, etc.
3494 * Default value is ``-1`` that is inclusive for all generations.
3496 It makes only sense for lists created with `fillParticleListFromMC` function with ``addDaughters=True`` argument.
3497 )DOC", Manager::VariableDataType::c_bool);
3499 REGISTER_METAVARIABLE("sourceObjectIsInList(particleListName)", sourceObjectIsInList, R
"DOC(
3500Returns 1 if the underlying mdst object (e.g. track, or cluster) was used to create a particle in ``particleListName``, 0 if not.
3503 This only makes sense for particles that are not composite. Returns -1 for composite particles.
3504)DOC", Manager::VariableDataType::c_int);
3506 REGISTER_METAVARIABLE("mcParticleIsInMCList(particleListName)", mcParticleIsInMCList, R
"DOC(
3507Returns 1 if the particle's matched MC particle is also matched to a particle in ``particleListName``
3508(or if either of the lists were filled from generator level `modularAnalysis.fillParticleListFromMC`.)
3510.. seealso:: :b2:var:`isMCDescendantOfList` to check daughters.
3511)DOC", Manager::VariableDataType::c_bool);
3513 REGISTER_METAVARIABLE("isGrandDaughterOfList(particleListNames)", isGrandDaughterOfList,
3514 "Returns 1 if the given particle is a grand daughter of at least one of the particles in the given particle Lists.", Manager::VariableDataType::c_bool);
3515 REGISTER_METAVARIABLE(
"originalParticle(variable)", originalParticle, R
"DOC(
3516 Returns value of variable for the original particle from which the given particle is copied.
3518 The copy of particle is created, for example, when the vertex fit updates the daughters and `modularAnalysis.copyParticles` is called.
3519 Returns NaN if the given particle is not copied and so there is no original particle.
3520 )DOC", Manager::VariableDataType::c_double);
3521 REGISTER_METAVARIABLE("daughter(i, variable)", daughter, R
"DOC(
3522 Returns value of variable for the i-th daughter. E.g.
3524 * ``daughter(0, p)`` returns the total momentum of the first daughter.
3525 * ``daughter(0, daughter(1, p)`` returns the total momentum of the second daughter of the first daughter.
3527 Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3528 )DOC", Manager::VariableDataType::c_double);
3529 REGISTER_METAVARIABLE("originalDaughter(i, variable)", originalDaughter, R
"DOC(
3530 Returns value of variable for the original particle from which the i-th daughter is copied.
3532 The copy of particle is created, for example, when the vertex fit updates the daughters and `modularAnalysis.copyParticles` is called.
3533 Returns NaN if the daughter is not copied and so there is no original daughter.
3535 Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3536 )DOC", Manager::VariableDataType::c_double);
3537 REGISTER_METAVARIABLE("mcDaughter(i, variable)", mcDaughter, R
"DOC(
3538 Returns the value of the requested variable for the i-th Monte Carlo daughter of the particle.
3540 Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3541 or if the i-th MC daughter does not exist.
3543 E.g. ``mcDaughter(0, PDG)`` will return the PDG code of the first MC daughter of the matched MC
3544 particle of the reconstructed particle the function is applied to.
3546 The meta variable can also be nested: ``mcDaughter(0, mcDaughter(1, PDG))``.
3547 )DOC", Manager::VariableDataType::c_double);
3548 REGISTER_METAVARIABLE("mcMother(variable)", mcMother, R
"DOC(
3549 Returns the value of the requested variable for the Monte Carlo mother of the particle.
3551 Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3552 or if the MC mother does not exist.
3554 E.g. ``mcMother(PDG)`` will return the PDG code of the MC mother of the matched MC
3555 particle of the reconstructed particle the function is applied to.
3557 The meta variable can also be nested: ``mcMother(mcMother(PDG))``.
3558 )DOC", Manager::VariableDataType::c_double);
3559 REGISTER_METAVARIABLE("genParticle(index, variable)", genParticle, R
"DOC(
3560[Eventbased] Returns the ``variable`` for the ith generator particle.
3561The arguments of the function must be the ``index`` of the particle in the MCParticle Array,
3562and ``variable``, the name of the function or variable for that generator particle.
3563If ``index`` goes beyond the length of the MCParticles array, NaN will be returned.
3565E.g. ``genParticle(0, p)`` returns the total momentum of the first MCParticle, which in a generic decay up to MC15 is
3566the Upsilon(4S) and for MC16 and beyond the initial electron.
3567)DOC", Manager::VariableDataType::c_double);
3568 REGISTER_METAVARIABLE("genUpsilon4S(variable)", genUpsilon4S, R
"DOC(
3569[Eventbased] Returns the ``variable`` evaluated for the generator-level :math:`\Upsilon(4S)`.
3570If no generator level :math:`\Upsilon(4S)` exists for the event, NaN will be returned.
3572E.g. ``genUpsilon4S(p)`` returns the total momentum of the :math:`\Upsilon(4S)` in a generic decay.
3573``genUpsilon4S(mcDaughter(1, p))`` returns the total momentum of the second daughter of the
3574generator-level :math:`\Upsilon(4S)` (i.e. the momentum of the second B meson in a generic decay).
3575)DOC", Manager::VariableDataType::c_double);
3576 REGISTER_METAVARIABLE("daughterProductOf(variable)", daughterProductOf,
3577 "Returns product of a variable over all daughters.\n"
3578 "E.g. ``daughterProductOf(extraInfo(SignalProbability))`` returns the product of the SignalProbabilitys of all daughters.", Manager::VariableDataType::c_double);
3579 REGISTER_METAVARIABLE(
"daughterSumOf(variable)", daughterSumOf,
3580 "Returns sum of a variable over all daughters.\n"
3581 "E.g. ``daughterSumOf(nDaughters)`` returns the number of grand-daughters.", Manager::VariableDataType::c_double);
3582 REGISTER_METAVARIABLE(
"daughterLowest(variable)", daughterLowest,
3583 "Returns the lowest value of the given variable among all daughters.\n"
3584 "E.g. ``useCMSFrame(daughterLowest(p))`` returns the lowest momentum in CMS frame.", Manager::VariableDataType::c_double);
3585 REGISTER_METAVARIABLE(
"daughterHighest(variable)", daughterHighest,
3586 "Returns the highest value of the given variable among all daughters.\n"
3587 "E.g. ``useCMSFrame(daughterHighest(p))`` returns the highest momentum in CMS frame.", Manager::VariableDataType::c_double);
3588 REGISTER_METAVARIABLE(
"daughterDiffOf(daughterIndex_i, daughterIndex_j, variable)", daughterDiffOf, R
"DOC(
3589 Returns the difference of a variable between the two given daughters.
3590 E.g. ``useRestFrame(daughterDiffOf(0, 1, p))`` returns the momentum difference between first and second daughter in the rest frame of the given particle.
3591 (That means that it returns :math:`p_j - p_i`)
3593 The daughters can be provided as generalized daughter indexes, which are simply colon-separated
3594 lists of daughter indexes, ordered starting from the root particle. For example, ``0:1``
3595 identifies the second daughter (1) of the first daughter (0) of the mother particle.
3597 )DOC", Manager::VariableDataType::c_double);
3598 REGISTER_METAVARIABLE("mcDaughterDiffOf(i, j, variable)", mcDaughterDiffOf,
3599 "MC matched version of the `daughterDiffOf` function.", Manager::VariableDataType::c_double);
3600 REGISTER_METAVARIABLE(
"grandDaughterDiffOf(i, j, variable)", grandDaughterDiffOf,
3601 "Returns the difference of a variable between the first daughters of the two given daughters.\n"
3602 "E.g. ``useRestFrame(grandDaughterDiffOf(0, 1, p))`` returns the momentum difference between the first daughters of the first and second daughter in the rest frame of the given particle.\n"
3603 "(That means that it returns :math:`p_j - p_i`)", Manager::VariableDataType::c_double);
3604 MAKE_DEPRECATED(
"grandDaughterDiffOf",
false,
"light-2402-ocicat", R
"DOC(
3605 The difference between any combination of (grand-)daughters can be calculated with the more general variable :b2:var:`daughterDiffOf`
3606 by using generalized daughter indexes.)DOC");
3607 REGISTER_METAVARIABLE("daughterDiffOfPhi(i, j)", daughterDiffOfPhi,
3608 "Returns the difference in :math:`\\phi` between the two given daughters. The unit of the angle is ``rad``.\n"
3609 "The difference is signed and takes account of the ordering of the given daughters.\n"
3610 "The function returns :math:`\\phi_j - \\phi_i`.", Manager::VariableDataType::c_double);
3611 MAKE_DEPRECATED(
"daughterDiffOfPhi",
false,
"release-06-00-00", R
"DOC(
3612 The difference of the azimuthal angle :math:`\\phi` of two daughters can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3613 REGISTER_METAVARIABLE("mcDaughterDiffOfPhi(i, j)", mcDaughterDiffOfPhi,
3614 "MC matched version of the `daughterDiffOfPhi` function. The unit of the angle is ``rad``", Manager::VariableDataType::c_double);
3615 MAKE_DEPRECATED(
"mcDaughterDiffOfPhi",
false,
"release-06-00-00", R
"DOC(
3616 The difference of the azimuthal angle :math:`\\phi` of the MC partners of two daughters can be calculated with the generic variable :b2:var:`mcDaughterDiffOf`.)DOC");
3617 REGISTER_METAVARIABLE("grandDaughterDiffOfPhi(i, j)", grandDaughterDiffOfPhi,
3618 "Returns the difference in :math:`\\phi` between the first daughters of the two given daughters. The unit of the angle is ``rad``.\n"
3619 "The difference is signed and takes account of the ordering of the given daughters.\n"
3620 "The function returns :math:`\\phi_j - \\phi_i`.\n", Manager::VariableDataType::c_double);
3621 MAKE_DEPRECATED(
"grandDaughterDiffOfPhi",
false,
"release-06-00-00", R
"DOC(
3622 The difference of the azimuthal angle :math:`\\phi` of two granddaughters can be calculated with the generic variable :b2:var:`grandDaughterDiffOf`.)DOC");
3623 REGISTER_METAVARIABLE("daughterDiffOfClusterPhi(i, j)", daughterDiffOfClusterPhi,
3624 "Returns the difference in :math:`\\phi` between the ECLClusters of two given daughters. The unit of the angle is ``rad``.\n"
3625 "The difference is signed and takes account of the ordering of the given daughters.\n"
3626 "The function returns :math:`\\phi_j - \\phi_i`.\n"
3627 "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.", Manager::VariableDataType::c_double);
3628 MAKE_DEPRECATED(
"daughterDiffOfClusterPhi",
false,
"release-06-00-00", R
"DOC(
3629 The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two daughters can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3630 REGISTER_METAVARIABLE("grandDaughterDiffOfClusterPhi(i, j)", grandDaughterDiffOfClusterPhi,
3631 "Returns the difference in :math:`\\phi` between the ECLClusters of the daughters of the two given daughters. The unit of the angle is ``rad``.\n"
3632 "The difference is signed and takes account of the ordering of the given daughters.\n"
3633 "The function returns :math:`\\phi_j - \\phi_i`.\n"
3634 "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.\n", Manager::VariableDataType::c_double);
3635 MAKE_DEPRECATED(
"grandDaughterDiffOfClusterPhi",
false,
"release-06-00-00", R
"DOC(
3636 The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two granddaughters can be calculated with the generic variable :b2:var:`grandDaughterDiffOf`.)DOC");
3637 REGISTER_METAVARIABLE("daughterDiffOfPhiCMS(i, j)", daughterDiffOfPhiCMS,
3638 "Returns the difference in :math:`\\phi` between the two given daughters in the CMS frame. The unit of the angle is ``rad``.\n"
3639 "The difference is signed and takes account of the ordering of the given daughters.\n"
3640 "The function returns :math:`\\phi_j - \\phi_i`.", Manager::VariableDataType::c_double);
3641 MAKE_DEPRECATED(
"daughterDiffOfPhiCMS",
false,
"release-06-00-00", R
"DOC(
3642 The difference of the azimuthal angle :math:`\\phi` of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3643 REGISTER_METAVARIABLE("mcDaughterDiffOfPhiCMS(i, j)", daughterDiffOfPhiCMS,
3644 "MC matched version of the `daughterDiffOfPhiCMS` function. The unit of the angle is ``rad``", Manager::VariableDataType::c_double);
3645 MAKE_DEPRECATED(
"mcDaughterDiffOfPhiCMS",
false,
"release-06-00-00", R
"DOC(
3646 The difference of the azimuthal angle :math:`\\phi` of the MC partners of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`mcDaughterDiffOf`.)DOC");
3647 REGISTER_METAVARIABLE("daughterDiffOfClusterPhiCMS(i, j)", daughterDiffOfClusterPhiCMS,
3648 "Returns the difference in :math:`\\phi` between the ECLClusters of two given daughters in the CMS frame. The unit of the angle is ``rad``.\n"
3649 "The difference is signed and takes account of the ordering of the given daughters.\n"
3650 "The function returns :math:`\\phi_j - \\phi_i``.\n"
3651 "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.", Manager::VariableDataType::c_double);
3652 MAKE_DEPRECATED(
"daughterDiffOfClusterPhiCMS",
false,
"release-06-00-00", R
"DOC(
3653 The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3654 REGISTER_METAVARIABLE("daughterNormDiffOf(i, j, variable)", daughterNormDiffOf,
3655 "Returns the normalized difference of a variable between the two given daughters.\n"
3656 "E.g. ``daughterNormDiffOf(0, 1, p)`` returns the normalized momentum difference between first and second daughter in the lab frame.", Manager::VariableDataType::c_double);
3657 REGISTER_METAVARIABLE(
"daughterMotherDiffOf(i, variable)", daughterMotherDiffOf,
3658 "Returns the difference of a variable between the given daughter and the mother particle itself.\n"
3659 "E.g. ``useRestFrame(daughterMotherDiffOf(0, p))`` returns the momentum difference between the given particle and its first daughter in the rest frame of the mother.", Manager::VariableDataType::c_double);
3660 REGISTER_METAVARIABLE(
"daughterMotherNormDiffOf(i, variable)", daughterMotherNormDiffOf,
3661 "Returns the normalized difference of a variable between the given daughter and the mother particle itself.\n"
3662 "E.g. ``daughterMotherNormDiffOf(1, p)`` returns the normalized momentum difference between the given particle and its second daughter in the lab frame.", Manager::VariableDataType::c_double);
3663 REGISTER_METAVARIABLE(
"daughterAngle(daughterIndex_1, daughterIndex_2[, daughterIndex_3])", daughterAngle, R
"DOC(
3664 Returns the angle in between any pair of particles belonging to the same decay tree.
3665 The unit of the angle is ``rad``.
3667 The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3668 daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3669 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3670 identifies the second daughter of the root particle.
3672 Both two and three generalized indexes can be given to ``daughterAngle``. If two indices are given, the
3673 variable returns the angle between the momenta of the two given particles. If three indices are given, the
3674 variable returns the angle between the momentum of the third particle and a vector which is the sum of the
3675 first two daughter momenta.
3678 ``daughterAngle(0, 3)`` will return the angle between the first and fourth daughter.
3679 ``daughterAngle(0, 1, 3)`` will return the angle between the fourth daughter and the sum of the first and
3681 ``daughterAngle(0:0, 3:0)`` will return the angle between the first daughter of the first daughter, and
3682 the first daughter of the fourth daughter.
3684 )DOC", Manager::VariableDataType::c_double);
3685 REGISTER_METAVARIABLE("mcDaughterAngle(daughterIndex_1, daughterIndex_2, [daughterIndex_3])", mcDaughterAngle,
3686 "MC matched version of the `daughterAngle` function. Also works if applied directly to MC particles. The unit of the angle is ``rad``", Manager::VariableDataType::c_double);
3687 REGISTER_VARIABLE(
"grandDaughterDecayAngle(i, j)", grandDaughterDecayAngle,
3688 "Returns the decay angle of the granddaughter in the daughter particle's rest frame.\n"
3689 "It is calculated with respect to the reverted momentum vector of the particle.\n"
3690 "Two arguments representing the daughter and granddaughter indices have to be provided as arguments.\n\n",
"rad");
3691 REGISTER_VARIABLE(
"daughterClusterAngleInBetween(i, j)", daughterClusterAngleInBetween,
3692 "Returns the angle between clusters associated to the two daughters."
3693 "If two indices given: returns the angle between the momenta of the clusters associated to the two given daughters."
3694 "If three indices given: returns the angle between the momentum of the third particle's cluster and a vector "
3695 "which is the sum of the first two daughter's cluster momenta."
3696 "Returns nan if any of the daughters specified don't have an associated cluster."
3697 "The arguments in the argument vector must be integers corresponding to the ith and jth (and kth) daughters.\n\n",
"rad");
3698 REGISTER_METAVARIABLE(
"daughterInvM(i[, j, ...])", daughterInvM, R
"DOC(
3699 Returns the invariant mass adding the Lorentz vectors of the given daughters. The unit of the invariant mass is GeV/:math:`\text{c}^2`
3700 E.g. ``daughterInvM(0, 1, 2)`` returns the invariant Mass :math:`m = \sqrt{(p_0 + p_1 + p_2)^2}` of the first, second and third daughter.
3702 Daughters from different generations of the decay tree can be combined using generalized daughter indexes,
3703 which are simply colon-separated daughter indexes for each generation, starting from the root particle. For
3704 example, ``0:1:3`` identifies the fourth daughter (3) of the second daughter (1) of the first daughter(0) of
3705 the mother particle.
3707 Returns NaN if the given daughter-index is out of bound (>= number of daughters))DOC", Manager::VariableDataType::c_double);
3708 REGISTER_METAVARIABLE("extraInfo(name)", extraInfo,
3709 "Returns extra info stored under the given name.\n"
3710 "The extraInfo has to be set by a module first.\n"
3711 "E.g. ``extraInfo(SignalProbability)`` returns the SignalProbability calculated by the ``MVAExpert`` module.\n"
3712 "If nothing is set under the given name or if the particle is a nullptr, NaN is returned.\n"
3713 "In the latter case please use `eventExtraInfo` if you want to access an EventExtraInfo variable.", Manager::VariableDataType::c_double);
3714 REGISTER_METAVARIABLE(
"eventExtraInfo(name)", eventExtraInfo,
3715 "[Eventbased] Returns extra info stored under the given name in the event extra info.\n"
3716 "The extraInfo has to be set first by another module like MVAExpert in event mode.\n"
3717 "If nothing is set under this name, NaN is returned.", Manager::VariableDataType::c_double);
3718 REGISTER_METAVARIABLE(
"eventCached(variable)", eventCached,
3719 "[Eventbased] Returns value of event-based variable and caches this value in the EventExtraInfo.\n"
3720 "The result of second call to this variable in the same event will be provided from the cache.\n"
3721 "It is recommended to use this variable in order to declare custom aliases as event-based. This is "
3722 "necessary if using the eventwise mode of variablesToNtuple).", Manager::VariableDataType::c_double);
3723 REGISTER_METAVARIABLE(
"particleCached(variable)", particleCached,
3724 "Returns value of given variable and caches this value in the ParticleExtraInfo of the provided particle.\n"
3725 "The result of second call to this variable on the same particle will be provided from the cache.", Manager::VariableDataType::c_double);
3726 REGISTER_METAVARIABLE(
"modulo(variable, n)", modulo,
3727 "Returns rest of division of variable by n.", Manager::VariableDataType::c_int);
3728 REGISTER_METAVARIABLE(
"abs(variable)", abs,
3729 "Returns absolute value of the given variable.\n"
3730 "E.g. abs(mcPDG) returns the absolute value of the mcPDG, which is often useful for cuts.", Manager::VariableDataType::c_double);
3731 REGISTER_METAVARIABLE(
"max(var1,var2)", max,
"Returns max value of two variables.\n", Manager::VariableDataType::c_double);
3732 REGISTER_METAVARIABLE(
"min(var1,var2)", min,
"Returns min value of two variables.\n", Manager::VariableDataType::c_double);
3733 REGISTER_METAVARIABLE(
"sin(variable)", sin,
"Returns sine value of the given variable.", Manager::VariableDataType::c_double);
3734 REGISTER_METAVARIABLE(
"asin(variable)", asin,
"Returns arcsine of the given variable. The unit of the asin() is ``rad``", Manager::VariableDataType::c_double);
3735 REGISTER_METAVARIABLE(
"cos(variable)", cos,
"Returns cosine value of the given variable.", Manager::VariableDataType::c_double);
3736 REGISTER_METAVARIABLE(
"acos(variable)", acos,
"Returns arccosine value of the given variable. The unit of the acos() is ``rad``", Manager::VariableDataType::c_double);
3737 REGISTER_METAVARIABLE(
"tan(variable)", tan,
"Returns tangent value of the given variable.", Manager::VariableDataType::c_double);
3738 REGISTER_METAVARIABLE(
"atan(variable)", atan,
"Returns arctangent value of the given variable. The unit of the atan() is ``rad``", Manager::VariableDataType::c_double);
3739 REGISTER_METAVARIABLE(
"exp(variable)", exp,
"Returns exponential evaluated for the given variable.", Manager::VariableDataType::c_double);
3740 REGISTER_METAVARIABLE(
"log(variable)", log,
"Returns natural logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3741 REGISTER_METAVARIABLE(
"log10(variable)", log10,
"Returns base-10 logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3742 REGISTER_METAVARIABLE(
"isNAN(variable)", isNAN,
3743 "Returns true if variable value evaluates to nan (determined via std::isnan(double)).\n"
3744 "Useful for debugging.", Manager::VariableDataType::c_bool);
3745 REGISTER_METAVARIABLE(
"ifNANgiveX(variable, x)", ifNANgiveX,
3746 "Returns x (has to be a number) if variable value is nan (determined via std::isnan(double)).\n"
3747 "Useful for technical purposes while training MVAs.", Manager::VariableDataType::c_double);
3748 REGISTER_METAVARIABLE(
"isInfinity(variable)", isInfinity,
3749 "Returns true if variable value evaluates to infinity (determined via std::isinf(double)).\n"
3750 "Useful for debugging.", Manager::VariableDataType::c_bool);
3751 REGISTER_METAVARIABLE(
"unmask(variable, flag1, flag2, ...)", unmask,
3752 "unmask(variable, flag1, flag2, ...) or unmask(variable, mask) sets certain bits in the variable to zero.\n"
3753 "For example, if you want to set the second, fourth and fifth bits to zero, you could call \n"
3754 "``unmask(variable, 2, 8, 16)`` or ``unmask(variable, 26)``.\n"
3755 "", Manager::VariableDataType::c_double);
3756 REGISTER_METAVARIABLE(
"conditionalVariableSelector(cut, variableIfTrue, variableIfFalse)", conditionalVariableSelector,
3757 "Returns one of the two supplied variables, depending on whether the particle passes the supplied cut.\n"
3758 "The first variable is returned if the particle passes the cut, and the second variable is returned otherwise.", Manager::VariableDataType::c_double);
3759 REGISTER_METAVARIABLE(
"pValueCombination(p1, p2, ...)", pValueCombination,
3760 "Returns the combined p-value of the provided p-values according to the formula given in `Nucl. Instr. and Meth. A 411 (1998) 449 <https://doi.org/10.1016/S0168-9002(98)00293-9>`_ .\n"
3761 "If any of the p-values is invalid, i.e. smaller than zero, -1 is returned.", Manager::VariableDataType::c_double);
3762 REGISTER_METAVARIABLE(
"veto(particleList, cut, pdgCode = 11)", veto,
3763 "Combines current particle with particles from the given particle list and returns 1 if the combination passes the provided cut. \n"
3764 "For instance one can apply this function on a signal Photon and provide a list of all photons in the rest of event and a cut \n"
3765 "around the neutral Pion mass (e.g. ``0.130 < M < 0.140``). \n"
3766 "If a combination of the signal Photon with a ROE photon fits this criteria, hence looks like a neutral pion, the veto-Metavariable will return 1", Manager::VariableDataType::c_bool);
3767 REGISTER_METAVARIABLE(
"matchedMC(variable)", matchedMC,
3768 "Returns variable output for the matched MCParticle by constructing a temporary Particle from it.\n"
3769 "This may not work too well if your variable requires accessing daughters of the particle.\n"
3770 "E.g. ``matchedMC(p)`` returns the total momentum of the related MCParticle.\n"
3771 "Returns NaN if no matched MCParticle exists.", Manager::VariableDataType::c_double);
3772 REGISTER_METAVARIABLE(
"clusterBestMatchedMCParticle(variable)", clusterBestMatchedMCParticle,
3773 "Returns variable output for the MCParticle that is best-matched with the ECLCluster of the given Particle.\n"
3774 "E.g. To get the energy of the MCParticle that matches best with an ECLCluster, one could use ``clusterBestMatchedMCParticle(E)``\n"
3775 "When the variable is called for ``gamma`` and if the ``gamma`` is matched with MCParticle, it works same as `matchedMC`.\n"
3776 "If the variable is called for ``gamma`` that fails to match with an MCParticle, it provides the mdst-level MCMatching information abouth the ECLCluster.\n"
3777 "Returns NaN if the particle is not matched to an ECLCluster, or if the ECLCluster has no matching MCParticles", Manager::VariableDataType::c_double);
3778 REGISTER_METAVARIABLE(
"varForBestMatchedMCKlong(variable)", clusterBestMatchedMCKlong,
3779 "Returns variable output for the Klong MCParticle which has the best match with the ECLCluster of the given Particle.\n"
3780 "Returns NaN if the particle is not matched to an ECLCluster, or if the ECLCluster has no matching Klong MCParticle", Manager::VariableDataType::c_double);
3782 REGISTER_METAVARIABLE(
"countInList(particleList, cut='')", countInList,
"[Eventbased] "
3783 "Returns number of particle which pass given in cut in the specified particle list.\n"
3784 "Useful for creating statistics about the number of particles in a list.\n"
3785 "E.g. ``countInList(e+, isSignal == 1)`` returns the number of correctly reconstructed electrons in the event.\n"
3786 "The variable is event-based and does not need a valid particle pointer as input.", Manager::VariableDataType::c_int);
3787 REGISTER_METAVARIABLE(
"getVariableByRank(particleList, rankedVariableName, variableName, rank)", getVariableByRank, R
"DOC(
3788 [Eventbased] Returns the value of ``variableName`` for the candidate in the ``particleList`` with the requested ``rank``.
3791 The `BestCandidateSelection` module available via `rankByHighest` / `rankByLowest` has to be used before.
3794 The first candidate matching the given rank is used.
3795 Thus, it is not recommended to use this variable in conjunction with ``allowMultiRank`` in the `BestCandidateSelection` module.
3797 The suffix ``_rank`` is automatically added to the argument ``rankedVariableName``,
3798 which either has to be the name of the variable used to order the candidates or the selected outputVariable name without the ending ``_rank``.
3799 This means that your selected name for the rank variable has to end with ``_rank``.
3801 An example of this variable's usage is given in the tutorial `B2A602-BestCandidateSelection <https://gitlab.desy.de/belle2/software/basf2/-/tree/main/analysis/examples/tutorials/B2A602-BestCandidateSelection.py>`_
3802 )DOC", Manager::VariableDataType::c_double);
3803 REGISTER_VARIABLE("matchedMCHasPDG(PDGCode)", matchedMCHasPDG,
3804 "Returns if the absolute value of the PDGCode of the MCParticle related to the Particle matches a given PDGCode."
3805 "Returns 0/NAN/1 if PDGCode does not match/is not available/ matches");
3806 REGISTER_METAVARIABLE(
"numberOfNonOverlappingParticles(pList1, pList2, ...)", numberOfNonOverlappingParticles,
3807 "Returns the number of non-overlapping particles in the given particle lists"
3808 "Useful to check if there is additional physics going on in the detector if one reconstructed the Y4S", Manager::VariableDataType::c_int);
3809 REGISTER_METAVARIABLE(
"totalEnergyOfParticlesInList(particleListName)", totalEnergyOfParticlesInList,
3810 "[Eventbased] Returns the total energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3811 REGISTER_METAVARIABLE(
"totalPxOfParticlesInList(particleListName)", totalPxOfParticlesInList,
3812 "[Eventbased] Returns the total momentum Px of particles in the given particle List. The unit of the momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3813 REGISTER_METAVARIABLE(
"totalPyOfParticlesInList(particleListName)", totalPyOfParticlesInList,
3814 "[Eventbased] Returns the total momentum Py of particles in the given particle List. The unit of the momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3815 REGISTER_METAVARIABLE(
"totalPzOfParticlesInList(particleListName)", totalPzOfParticlesInList,
3816 "[Eventbased] Returns the total momentum Pz of particles in the given particle List. The unit of the momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3817 REGISTER_METAVARIABLE(
"invMassInLists(pList1, pList2, ...)", invMassInLists,
3818 "[Eventbased] Returns the invariant mass of the combination of particles in the given particle lists. The unit of the invariant mass is GeV/:math:`\\text{c}^2` ", Manager::VariableDataType::c_double);
3819 REGISTER_METAVARIABLE(
"totalECLEnergyOfParticlesInList(particleListName)", totalECLEnergyOfParticlesInList,
3820 "[Eventbased] Returns the total ECL energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3821 REGISTER_METAVARIABLE(
"maxPtInList(particleListName)", maxPtInList,
3822 "[Eventbased] Returns maximum transverse momentum Pt in the given particle List. The unit of the transverse momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3823 REGISTER_METAVARIABLE(
"eclClusterSpecialTrackMatched(cut)", eclClusterTrackMatchedWithCondition,
3824 "Returns if at least one Track that satisfies the given condition is related to the ECLCluster of the Particle.", Manager::VariableDataType::c_double);
3825 REGISTER_METAVARIABLE(
"averageValueInList(particleListName, variable)", averageValueInList,
3826 "[Eventbased] Returns the arithmetic mean of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3827 REGISTER_METAVARIABLE(
"medianValueInList(particleListName, variable)", medianValueInList,
3828 "[Eventbased] Returns the median value of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3829 REGISTER_METAVARIABLE(
"sumValueInList(particleListName, variable)", sumValueInList,
3830 "[Eventbased] Returns the sum of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3831 REGISTER_METAVARIABLE(
"productValueInList(particleListName, variable)", productValueInList,
3832 "[Eventbased] Returns the product of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3833 REGISTER_METAVARIABLE(
"angleToClosestInList(particleListName)", angleToClosestInList,
3834 "Returns the angle between this particle and the closest particle (smallest opening angle) in the list provided. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
3835 REGISTER_METAVARIABLE(
"closestInList(particleListName, variable)", closestInList,
3836 "Returns `variable` for the closest particle (smallest opening angle) in the list provided.", Manager::VariableDataType::c_double);
3837 REGISTER_METAVARIABLE(
"angleToMostB2BInList(particleListName)", angleToMostB2BInList,
3838 "Returns the angle between this particle and the most back-to-back particle (closest opening angle to 180) in the list provided. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
3839 REGISTER_METAVARIABLE(
"deltaPhiToMostB2BPhiInList(particleListName)", deltaPhiToMostB2BPhiInList,
3840 "Returns the abs(delta phi) between this particle and the most back-to-back particle in phi (closest opening angle to 180) in the list provided. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
3841 REGISTER_METAVARIABLE(
"mostB2BInList(particleListName, variable)", mostB2BInList,
3842 "Returns `variable` for the most back-to-back particle (closest opening angle to 180) in the list provided.", Manager::VariableDataType::c_double);
3843 REGISTER_METAVARIABLE(
"maxOpeningAngleInList(particleListName)", maxOpeningAngleInList,
3844 "[Eventbased] Returns maximum opening angle in the given particle List. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
3845 REGISTER_METAVARIABLE(
"daughterCombination(variable, daughterIndex_1, daughterIndex_2 ... daughterIndex_n)", daughterCombination,R
"DOC(
3846Returns a ``variable`` function only of the 4-momentum calculated on an arbitrary set of (grand)daughters.
3849 ``variable`` can only be a function of the daughters' 4-momenta.
3851Daughters from different generations of the decay tree can be combined using generalized daughter indexes, which are simply colon-separated
3852the list of daughter indexes, starting from the root particle: for example, ``0:1:3`` identifies the fourth
3853daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle.
3856 ``daughterCombination(M, 0, 3, 4)`` will return the invariant mass of the system made of the first, fourth and fifth daughter of particle.
3857 ``daughterCombination(M, 0:0, 3:0)`` will return the invariant mass of the system made of the first daughter of the first daughter and the first daughter of the fourth daughter.
3859)DOC", Manager::VariableDataType::c_double);
3860 REGISTER_METAVARIABLE("useAlternativeDaughterHypothesis(variable, daughterIndex_1:newMassHyp_1, ..., daughterIndex_n:newMassHyp_n)", useAlternativeDaughterHypothesis,R
"DOC(
3861Returns a ``variable`` calculated using new mass hypotheses for (some of) the particle's daughters.
3864 ``variable`` can only be a function of the particle 4-momentum, which is re-calculated as the sum of the daughters' 4-momenta, and the daughters' 4-momentum.
3865 This means that if you made a kinematic fit without updating the daughters' momenta, the result of this variable will not reflect the effect of the kinematic fit.
3866 Also, the track fit is not performed again: the variable only re-calculates the 4-vectors using different mass assumptions.
3867 In the variable, a copy of the given particle is created with daughters' alternative mass assumption (i.e. the original particle and daughters are not changed).
3870 Generalized daughter indexes are not supported (yet!): this variable can be used only on first-generation daughters.
3873 ``useAlternativeDaughterHypothesis(M, 0:K+, 2:pi-)`` will return the invariant mass of the particle assuming that the first daughter is a kaon and the third is a pion, instead of whatever was used in reconstructing the decay.
3874 ``useAlternativeDaughterHypothesis(mRecoil, 1:p+)`` will return the recoil mass of the particle assuming that the second daughter is a proton instead of whatever was used in reconstructing the decay.
3876)DOC", Manager::VariableDataType::c_double);
3877 REGISTER_METAVARIABLE("varForFirstMCAncestorOfType(type, variable)",varForFirstMCAncestorOfType,R
"DOC(Returns requested variable of the first ancestor of the given type.
3878Ancestor type can be set up by PDG code or by particle name (check evt.pdl for valid particle names))DOC", Manager::VariableDataType::c_double);
3880 REGISTER_METAVARIABLE("nTrackFitResults(particleType)", nTrackFitResults,
3881 "[Eventbased] Returns the total number of TrackFitResults for a given particleType. The argument can be the name of particle (e.g. pi+) or PDG code (e.g. 211).",
3882 Manager::VariableDataType::c_int);
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
int getPDGCode() const
PDG code.
static const ChargedStable pion
charged pion particle
static const double doubleNaN
quiet_NaN
static const ChargedStable electron
electron 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)
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...
@ c_Initial
bit 5: Particle is initial such as e+ or e- and not going to Geant4
@ c_PrimaryParticle
bit 0: Particle is primary particle.
@ c_IsVirtual
bit 4: Particle is virtual and not going to Geant4.
static std::string makeROOTCompatible(std::string str)
Remove special characters that ROOT dislikes in branch names, e.g.
const MCParticle * getMCParticle() const
Returns the pointer to the MCParticle object that was used to create this Particle (ParticleType == c...
EParticleSourceObject
particle source enumerators
const Particle * getParticleFromGeneralizedIndexString(const std::string &generalizedIndex) const
Explores the decay tree of the particle and returns the (grand^n)daughter identified by a generalized...
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
@ c_Unflavored
Is its own antiparticle or we don't know whether it is a particle/antiparticle.
@ c_Flavored
Is either particle or antiparticle.
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
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
bool hasAntiParticle(int pdgCode)
Checks if the particle with given pdg code has an anti-particle or not.
Particle * copyParticle(const Particle *original)
Function takes argument Particle and creates a copy of it and copies of all its (grand-)^n-daughters.
Abstract base class for different kinds of events.
const std::vector< double > v2
MATLAB generated random vector.
const std::vector< double > v1
MATLAB generated random vector.