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/EventKinematics.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>
45#include <TDatabasePDG.h>
46#include <Math/Vector4D.h>
47#include <Math/VectorUtil.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) {
105 auto daughterFunction = convertToDaughterIndex({arguments[1]});
106 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
107 int daughterIndexTagB = std::get<int>(daughterFunction(particle));
108 if (daughterIndexTagB < 0)
111 if (particle->getPDGCode() != 300553)
113 B2ERROR(
"Variable should only be used on a Upsilon(4S) Particle List!");
118 ROOT::Math::PxPyPzEVector pSigB = T.getBeamFourMomentum() - particle->getDaughter(daughterIndexTagB)->get4Vector();
119 Particle tmp(pSigB, -particle->getDaughter(daughterIndexTagB)->getPDGCode());
121 UseReferenceFrame<RestFrame> frame(&tmp);
122 double result = std::get<double>(var->function(particle));
128 B2FATAL(
"Wrong number of arguments for meta function useTagSideRecoilRestFrame");
134 if (arguments.size() == 2) {
136 std::string listName = arguments[1];
137 auto func = [var, listName](
const Particle * particle) ->
double {
138 StoreObjPtr<ParticleList> list(listName);
139 unsigned listSize = list->getListSize();
143 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."
144 << LogVar(
"ParticleList", listName)
145 << LogVar(
"Number of candidates in the list", listSize));
146 const Particle* p = list->getParticle(0);
147 UseReferenceFrame<RestFrame> frame(p);
148 double result = std::get<double>(var->function(particle));
153 B2FATAL(
"Wrong number of arguments for meta function useParticleRestFrame.");
159 if (arguments.size() == 2) {
161 std::string listName = arguments[1];
162 auto func = [var, listName](
const Particle * particle) ->
double {
163 StoreObjPtr<ParticleList> list(listName);
164 unsigned listSize = list->getListSize();
168 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."
169 << LogVar(
"ParticleList", listName)
170 << LogVar(
"Number of candidates in the list", listSize));
171 const Particle* p = list->getParticle(0);
173 ROOT::Math::PxPyPzEVector recoil = T.getBeamFourMomentum() - p->get4Vector();
175 Particle pRecoil(recoil, 0);
176 pRecoil.setVertex(particle->getVertex());
177 UseReferenceFrame<RestFrame> frame(&pRecoil);
178 double result = std::get<double>(var->function(particle));
183 B2FATAL(
"Wrong number of arguments for meta function useParticleRestFrame.");
189 if (arguments.size() >= 2) {
191 auto func = [var, arguments](
const Particle * particle) ->
double {
194 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
196 for (
unsigned int i = 1; i < arguments.size(); i++)
198 auto generalizedIndex = arguments[i];
199 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
201 pSum += dauPart->get4Vector();
205 Particle tmp(pSum, 0);
206 UseReferenceFrame<RestFrame> frame(&tmp);
207 double result = std::get<double>(var->function(particle));
212 B2FATAL(
"Wrong number of arguments for meta function useDaughterRestFrame.");
218 if (arguments.size() >= 2) {
220 auto func = [var, arguments](
const Particle * particle) ->
double {
223 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
225 for (
unsigned int i = 1; i < arguments.size(); i++)
227 auto generalizedIndex = arguments[i];
228 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
230 pSum += dauPart->get4Vector();
235 ROOT::Math::PxPyPzEVector recoil = T.getBeamFourMomentum() - pSum;
237 Particle pRecoil(recoil, 0);
238 UseReferenceFrame<RestFrame> frame(&pRecoil);
239 double result = std::get<double>(var->function(particle));
244 B2FATAL(
"Wrong number of arguments for meta function useDaughterRecoilRestFrame.");
250 if (arguments.size() == 1) {
252 auto func = [var](
const Particle * particle) ->
double {
253 int index = ancestorBIndex(particle);
255 StoreArray<MCParticle> mcparticles;
256 Particle temp(mcparticles[index]);
257 UseReferenceFrame<RestFrame> frame(&temp);
258 double result = std::get<double>(var->function(particle));
263 B2FATAL(
"Wrong number of arguments for meta function useMCancestorBRestFrame.");
269 if (arguments.size() == 1) {
270 auto extraInfoName = arguments[0];
271 auto func = [extraInfoName](
const Particle * particle) ->
double {
272 if (particle ==
nullptr)
274 B2WARNING(
"Returns NaN because the particle is nullptr! If you want EventExtraInfo variables, please use eventExtraInfo() instead");
277 if (particle->hasExtraInfo(extraInfoName))
279 return particle->getExtraInfo(extraInfoName);
287 B2FATAL(
"Wrong number of arguments for meta function extraInfo");
293 if (arguments.size() == 1) {
294 auto extraInfoName = arguments[0];
295 auto func = [extraInfoName](
const Particle*) ->
double {
296 StoreObjPtr<EventExtraInfo> eventExtraInfo;
297 if (not eventExtraInfo.isValid())
299 if (eventExtraInfo->hasExtraInfo(extraInfoName))
301 return eventExtraInfo->getExtraInfo(extraInfoName);
309 B2FATAL(
"Wrong number of arguments for meta function extraInfo");
315 if (arguments.size() == 1) {
318 auto func = [var, key](
const Particle*) ->
double {
320 StoreObjPtr<EventExtraInfo> eventExtraInfo;
321 if (not eventExtraInfo.isValid())
322 eventExtraInfo.create();
323 if (eventExtraInfo->hasExtraInfo(key))
325 return eventExtraInfo->getExtraInfo(key);
329 auto var_result = var->function(
nullptr);
330 if (std::holds_alternative<double>(var_result)) {
331 value = std::get<double>(var_result);
332 }
else if (std::holds_alternative<int>(var_result)) {
333 return std::get<int>(var_result);
334 }
else if (std::holds_alternative<bool>(var_result)) {
335 return std::get<bool>(var_result);
337 eventExtraInfo->addExtraInfo(key, value);
343 B2FATAL(
"Wrong number of arguments for meta function eventCached");
349 if (arguments.size() == 1) {
352 auto func = [var, key](
const Particle * particle) ->
double {
354 if (particle->hasExtraInfo(key))
356 return particle->getExtraInfo(key);
359 double value = std::get<double>(var->function(particle));
367 const_cast<Particle*
>(particle)->addExtraInfo(key, value);
373 B2FATAL(
"Wrong number of arguments for meta function particleCached");
382 if (arguments.size() != 1) B2FATAL(
"Wrong number of arguments for meta function formula");
383 FormulaParser<VariableFormulaConstructor> parser;
385 return parser.parse(arguments[0]);
386 }
catch (std::runtime_error& e) {
393 if (arguments.size() <= 1) {
395 std::string cutString;
396 if (arguments.size() == 1)
397 cutString = arguments[0];
399 auto func = [cut](
const Particle*) ->
int {
401 int number_of_tracks = 0;
402 StoreArray<Track> tracks;
403 for (
const auto& track : tracks)
405 const TrackFitResult* trackFit = track.getTrackFitResultWithClosestMass(
Const::pion);
406 if (!trackFit)
continue;
407 if (trackFit->getChargeSign() == 0) {
411 if (cut->check(&particle))
416 return number_of_tracks;
421 B2FATAL(
"Wrong number of arguments for meta function nCleanedTracks");
427 if (arguments.size() <= 1) {
429 std::string cutString;
430 if (arguments.size() == 1)
431 cutString = arguments[0];
433 auto func = [cut](
const Particle*) ->
int {
435 int number_of_clusters = 0;
436 StoreArray<ECLCluster> clusters;
437 for (
const auto& cluster : clusters)
443 Particle particle(&cluster);
444 if (cut->check(&particle))
445 number_of_clusters++;
448 return number_of_clusters;
453 B2FATAL(
"Wrong number of arguments for meta function nCleanedECLClusters");
459 if (arguments.size() == 1) {
460 std::string cutString = arguments[0];
462 auto func = [cut](
const Particle * particle) ->
bool {
463 if (cut->check(particle))
470 B2FATAL(
"Wrong number of arguments for meta function passesCut");
476 if (arguments.size() == 1) {
477 std::string cutString = arguments[0];
479 auto func = [cut](
const Particle*) ->
bool {
480 if (cut->check(
nullptr))
487 B2FATAL(
"Wrong number of arguments for meta function passesEventCut");
493 if (arguments.size() == 2) {
497 }
catch (std::invalid_argument&) {
498 B2FATAL(
"The first argument of varFor meta function must be a positive integer!");
501 auto func = [pdgCode, var](
const Particle * particle) ->
double {
502 if (std::abs(particle->getPDGCode()) == std::abs(pdgCode))
504 auto var_result = var->function(particle);
505 if (std::holds_alternative<double>(var_result)) {
506 return std::get<double>(var_result);
507 }
else if (std::holds_alternative<int>(var_result)) {
508 return std::get<int>(var_result);
509 }
else if (std::holds_alternative<bool>(var_result)) {
510 return std::get<bool>(var_result);
516 B2FATAL(
"Wrong number of arguments for meta function varFor");
522 if (arguments.size() == 1) {
524 auto func = [var](
const Particle * particle) ->
double {
525 if (particle->getMCParticle())
530 auto var_result = var->function(particle);
531 if (std::holds_alternative<double>(var_result)) {
532 return std::get<double>(var_result);
533 }
else if (std::holds_alternative<int>(var_result)) {
534 return std::get<int>(var_result);
535 }
else if (std::holds_alternative<bool>(var_result)) {
536 return std::get<bool>(var_result);
543 B2FATAL(
"Wrong number of arguments for meta function varForMCGen");
549 if (arguments.size() == 1) {
550 std::string listName = arguments[0];
551 auto func = [listName](
const Particle * particle) ->
int {
554 StoreObjPtr<ParticleList> listOfParticles(listName);
556 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to nParticlesInList");
558 return listOfParticles->getListSize();
563 B2FATAL(
"Wrong number of arguments for meta function nParticlesInList");
570 if (arguments.size() != 1) {
571 B2FATAL(
"Wrong number of arguments for isInList");
573 auto listName = arguments[0];
575 auto func = [listName](
const Particle * particle) ->
bool {
578 StoreObjPtr<ParticleList> list(listName);
579 if (!(list.isValid()))
581 B2FATAL(
"Invalid Listname " << listName <<
" given to isInList");
585 return list->contains(particle);
594 if (arguments.size() != 1) {
595 B2FATAL(
"Wrong number of arguments for sourceObjectIsInList");
597 auto listName = arguments[0];
599 auto func = [listName](
const Particle * particle) ->
int {
602 StoreObjPtr<ParticleList> list(listName);
603 if (!(list.isValid()))
605 B2FATAL(
"Invalid Listname " << listName <<
" given to sourceObjectIsInList");
611 if (particlesource == Particle::EParticleSourceObject::c_Composite
612 or particlesource == Particle::EParticleSourceObject::c_Undefined)
618 for (
unsigned i = 0; i < list->getListSize(); ++i)
620 Particle* iparticle = list->getParticle(i);
621 if (particle->getMdstSource() == iparticle->getMdstSource())
633 if (arguments.size() != 1) {
634 B2FATAL(
"Wrong number of arguments for mcParticleIsInMCList");
636 auto listName = arguments[0];
638 auto func = [listName](
const Particle * particle) ->
bool {
641 StoreObjPtr<ParticleList> list(listName);
642 if (!(list.isValid()))
643 B2FATAL(
"Invalid Listname " << listName <<
" given to mcParticleIsInMCList");
646 const MCParticle* mcp = particle->getMCParticle();
647 if (mcp ==
nullptr)
return false;
650 for (
unsigned i = 0; i < list->getListSize(); ++i)
652 const MCParticle* imcp = list->getParticle(i)->getMCParticle();
653 if ((imcp !=
nullptr) and (mcp->getArrayIndex() == imcp->getArrayIndex()))
663 B2WARNING(
"isDaughterOfList is outdated and replaced by isDescendantOfList.");
664 std::vector<std::string> new_arguments = arguments;
665 new_arguments.push_back(std::string(
"1"));
666 return isDescendantOfList(new_arguments);
671 B2WARNING(
"isGrandDaughterOfList is outdated and replaced by isDescendantOfList.");
672 std::vector<std::string> new_arguments = arguments;
673 new_arguments.push_back(std::string(
"2"));
674 return isDescendantOfList(new_arguments);
679 if (arguments.size() > 0) {
680 auto listNames = arguments;
681 auto func = [listNames](
const Particle * particle) ->
bool {
683 int generation_flag = -1;
687 }
catch (std::exception& e) {}
689 for (
auto& iListName : listNames)
694 }
catch (std::exception& e) {}
697 auto list_comparison = [](
auto&& self,
const Particle * m,
const Particle * p,
int flag)->
bool {
699 for (
unsigned i = 0; i < m->getNDaughters(); ++i)
701 const Particle* daughter = m->getDaughter(i);
702 if ((flag == 1.) or (flag < 0)) {
703 if (p->isCopyOf(daughter)) {
709 if (daughter->getNDaughters() > 0) {
710 result = self(self, daughter, p, flag - 1);
720 StoreObjPtr<ParticleList> listOfParticles(iListName);
722 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << iListName <<
" given to isDescendantOfList");
724 for (
unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
725 Particle* iParticle = listOfParticles->getParticle(i);
726 output = list_comparison(list_comparison, iParticle, particle, generation_flag);
736 B2FATAL(
"Wrong number of arguments for meta function isDescendantOfList");
742 if (arguments.size() > 0) {
743 auto listNames = arguments;
744 auto func = [listNames](
const Particle * particle) ->
bool {
746 int generation_flag = -1;
750 }
catch (std::exception& e) {}
752 if (particle->getMCParticle() ==
nullptr)
757 for (
auto& iListName : listNames)
760 std::stod(iListName);
762 }
catch (std::exception& e) {}
764 auto list_comparison = [](
auto&& self,
const Particle * m,
const Particle * p,
int flag)->
bool {
766 for (
unsigned i = 0; i < m->getNDaughters(); ++i)
768 const Particle* daughter = m->getDaughter(i);
769 if ((flag == 1.) or (flag < 0)) {
770 if (daughter->getMCParticle() !=
nullptr) {
771 if (p->getMCParticle()->getArrayIndex() == daughter->getMCParticle()->getArrayIndex()) {
777 if (daughter->getNDaughters() > 0) {
778 result = self(self, daughter, p, flag - 1);
788 StoreObjPtr<ParticleList> listOfParticles(iListName);
790 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << iListName <<
" given to isMCDescendantOfList");
792 for (
unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
793 Particle* iParticle = listOfParticles->getParticle(i);
794 output = list_comparison(list_comparison, iParticle, particle, generation_flag);
804 B2FATAL(
"Wrong number of arguments for meta function isMCDescendantOfList");
810 if (arguments.size() == 1) {
812 auto func = [var](
const Particle * particle) ->
double {
813 double product = 1.0;
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 product *= 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 product *= std::get<int>(var->function(particle->getDaughter(j)));
833 B2FATAL(
"Wrong number of arguments for meta function daughterProductOf");
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 sum += std::get<double>(var->function(particle->getDaughter(j)));
852 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
854 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
855 sum += std::get<int>(var->function(particle->getDaughter(j)));
862 B2FATAL(
"Wrong number of arguments for meta function daughterSumOf");
868 if (arguments.size() == 1) {
870 auto func = [var](
const Particle * particle) ->
double {
872 if (particle->getNDaughters() == 0)
876 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
878 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
879 double iValue = std::get<double>(var->function(particle->getDaughter(j)));
880 if (std::isnan(iValue))
continue;
881 if (std::isnan(min)) min = iValue;
882 if (iValue < min) min = iValue;
884 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
886 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
887 int iValue = std::get<int>(var->function(particle->getDaughter(j)));
888 if (std::isnan(min)) min = iValue;
889 if (iValue < min) min = iValue;
896 B2FATAL(
"Wrong number of arguments for meta function daughterLowest");
902 if (arguments.size() == 1) {
904 auto func = [var](
const Particle * particle) ->
double {
906 if (particle->getNDaughters() == 0)
910 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
912 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
913 double iValue = std::get<double>(var->function(particle->getDaughter(j)));
914 if (std::isnan(iValue))
continue;
915 if (std::isnan(max)) max = iValue;
916 if (iValue > max) max = iValue;
918 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
920 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
921 int iValue = std::get<int>(var->function(particle->getDaughter(j)));
922 if (std::isnan(max)) max = iValue;
923 if (iValue > max) max = iValue;
930 B2FATAL(
"Wrong number of arguments for meta function daughterHighest");
936 if (arguments.size() == 3) {
937 auto func = [arguments](
const Particle * particle) ->
double {
938 if (particle ==
nullptr)
940 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
941 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
942 auto variablename = arguments[2];
943 if (dau_i ==
nullptr || dau_j ==
nullptr)
945 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
949 auto result_j = var->function(dau_j);
950 auto result_i = var->function(dau_i);
952 if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i))
954 diff = std::get<double>(result_j) - std::get<double>(result_i);
955 }
else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i))
957 diff = std::get<int>(result_j) - std::get<int>(result_i);
960 throw std::runtime_error(
"Bad variant access");
962 if (variablename ==
"phi" or variablename ==
"clusterPhi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)"))
963 or std::regex_match(variablename, std::regex(
"use.*Frame\\(clusterPhi\\)")))
965 if (fabs(diff) > M_PI) {
967 diff = diff - 2 * M_PI;
969 diff = 2 * M_PI + diff;
977 B2FATAL(
"Wrong number of arguments for meta function daughterDiffOf");
983 if (arguments.size() == 3) {
984 auto func = [arguments](
const Particle * particle) ->
double {
985 if (particle ==
nullptr)
987 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
988 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
989 auto variablename = arguments[2];
990 if (dau_i ==
nullptr || dau_j ==
nullptr)
992 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
995 const MCParticle* iMcDaughter = dau_i->getMCParticle();
996 const MCParticle* jMcDaughter = dau_j->getMCParticle();
997 if (iMcDaughter ==
nullptr || jMcDaughter ==
nullptr)
999 Particle iTmpPart(iMcDaughter);
1000 Particle jTmpPart(jMcDaughter);
1002 auto result_j = var->function(&jTmpPart);
1003 auto result_i = var->function(&iTmpPart);
1005 if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i))
1007 diff = std::get<double>(result_j) - std::get<double>(result_i);
1008 }
else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i))
1010 diff = std::get<int>(result_j) - std::get<int>(result_i);
1013 throw std::runtime_error(
"Bad variant access");
1015 if (variablename ==
"phi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)")))
1017 if (fabs(diff) > M_PI) {
1019 diff = diff - 2 * M_PI;
1021 diff = 2 * M_PI + diff;
1029 B2FATAL(
"Wrong number of arguments for meta function mcDaughterDiffOf");
1035 if (arguments.size() == 5) {
1041 }
catch (std::invalid_argument&) {
1042 B2FATAL(
"First four arguments of grandDaughterDiffOf meta function must be integers!");
1044 std::vector<std::string> new_arguments;
1045 new_arguments.push_back(std::string(arguments[0] +
":" + arguments[2]));
1046 new_arguments.push_back(std::string(arguments[1] +
":" + arguments[3]));
1047 new_arguments.push_back(arguments[4]);
1048 return daughterDiffOf(new_arguments);
1050 B2FATAL(
"Wrong number of arguments for meta function grandDaughterDiffOf");
1056 if (arguments.size() == 3) {
1057 auto func = [arguments](
const Particle * particle) ->
double {
1058 if (particle ==
nullptr)
1060 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
1061 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
1062 if (!(dau_i && dau_j))
1064 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
1068 double iValue, jValue;
1069 if (std::holds_alternative<double>(var->function(dau_j)))
1071 iValue = std::get<double>(var->function(dau_i));
1072 jValue = std::get<double>(var->function(dau_j));
1073 }
else if (std::holds_alternative<int>(var->function(dau_j)))
1075 iValue = std::get<int>(var->function(dau_i));
1076 jValue = std::get<int>(var->function(dau_j));
1078 return (jValue - iValue) / (jValue + iValue);
1082 B2FATAL(
"Wrong number of arguments for meta function daughterNormDiffOf");
1088 if (arguments.size() == 2) {
1089 auto daughterFunction = convertToDaughterIndex({arguments[0]});
1090 std::string variableName = arguments[1];
1091 auto func = [daughterFunction, variableName](
const Particle * particle) ->
double {
1092 if (particle ==
nullptr)
1094 int daughterNumber = std::get<int>(daughterFunction(particle));
1095 if (daughterNumber >=
int(particle->getNDaughters()) or daughterNumber < 0)
1098 auto result_mother = var->function(particle);
1099 auto result_daughter = var->function(particle->getDaughter(daughterNumber));
1101 if (std::holds_alternative<double>(result_mother) && std::holds_alternative<double>(result_daughter))
1103 diff = std::get<double>(result_mother) - std::get<double>(result_daughter);
1104 }
else if (std::holds_alternative<int>(result_mother) && std::holds_alternative<int>(result_daughter))
1106 diff = std::get<int>(result_mother) - std::get<int>(result_daughter);
1109 throw std::runtime_error(
"Bad variant access");
1112 if (variableName ==
"phi" or variableName ==
"useCMSFrame(phi)")
1114 if (fabs(diff) > M_PI) {
1116 diff = diff - 2 * M_PI;
1118 diff = 2 * M_PI + diff;
1126 B2FATAL(
"Wrong number of arguments for meta function daughterMotherDiffOf");
1132 if (arguments.size() == 2) {
1133 auto daughterFunction = convertToDaughterIndex({arguments[0]});
1135 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
1136 if (particle ==
nullptr)
1138 int daughterNumber = std::get<int>(daughterFunction(particle));
1139 if (daughterNumber >=
int(particle->getNDaughters()) or daughterNumber < 0)
1141 double daughterValue = 0.0, motherValue = 0.0;
1142 if (std::holds_alternative<double>(var->function(particle)))
1144 daughterValue = std::get<double>(var->function(particle->getDaughter(daughterNumber)));
1145 motherValue = std::get<double>(var->function(particle));
1146 }
else if (std::holds_alternative<int>(var->function(particle)))
1148 daughterValue = std::get<int>(var->function(particle->getDaughter(daughterNumber)));
1149 motherValue = std::get<int>(var->function(particle));
1151 return (motherValue - daughterValue) / (motherValue + daughterValue);
1155 B2FATAL(
"Wrong number of arguments for meta function daughterMotherNormDiffOf");
1161 if (arguments.size() >= 1) {
1163 auto func = [arguments](
const Particle * particle) ->
double {
1164 if (particle ==
nullptr)
1169 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
1170 for (
auto& generalizedIndex : arguments)
1172 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1173 if (dauPart) pSum += frame.getMomentum(dauPart);
1175 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1181 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
1182 ROOT::Math::PxPyPzEVector pRecoil = frame.getMomentum(pIN - particle->get4Vector());
1184 return ROOT::Math::VectorUtil::Angle(pRecoil, pSum);
1188 B2FATAL(
"Wrong number of arguments for meta function angleBetweenDaughterAndRecoil");
1192 Manager::FunctionPtr angleBetweenDaughterAndMissingMomentum(
const std::vector<std::string>& arguments)
1194 if (arguments.size() >= 1) {
1195 auto func = [arguments](
const Particle * particle) ->
double {
1196 if (particle ==
nullptr)
1199 StoreObjPtr<EventKinematics> evtShape;
1202 B2WARNING(
"Cannot find missing momentum information, did you forget to run EventKinematicsModule?");
1205 ROOT::Math::XYZVector missingMomentumCMS = evtShape->getMissingMomentumCMS();
1206 ROOT::Math::PxPyPzEVector missingTotalMomentumCMS(missingMomentumCMS.X(),
1207 missingMomentumCMS.Y(),
1208 missingMomentumCMS.Z(),
1209 evtShape->getMissingEnergyCMS());
1211 ROOT::Math::PxPyPzEVector missingTotalMomentumLab = T.rotateCmsToLab() * missingTotalMomentumCMS;
1214 ROOT::Math::PxPyPzEVector pMiss = frame.getMomentum(missingTotalMomentumLab);
1216 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
1217 for (
auto& generalizedIndex : arguments)
1219 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1220 if (dauPart) pSum += frame.getMomentum(dauPart);
1222 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1227 return ROOT::Math::VectorUtil::Angle(pMiss, pSum);
1231 B2FATAL(
"Wrong number of arguments for meta function angleBetweenDaughterAndMissingMomentum");
1237 if (arguments.size() == 2 || arguments.size() == 3) {
1239 auto func = [arguments](
const Particle * particle) ->
double {
1240 if (particle ==
nullptr)
1243 std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1247 for (
auto& generalizedIndex : arguments)
1249 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1251 pDaus.push_back(frame.getMomentum(dauPart));
1253 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1259 if (pDaus.size() == 2)
1260 return ROOT::Math::VectorUtil::Angle(pDaus[0], pDaus[1]);
1262 return ROOT::Math::VectorUtil::Angle(pDaus[2], pDaus[0] + pDaus[1]);
1266 B2FATAL(
"Wrong number of arguments for meta function daughterAngle");
1270 double grandDaughterDecayAngle(
const Particle* particle,
const std::vector<double>& arguments)
1272 if (arguments.size() == 2) {
1277 int daughterIndex = std::lround(arguments[0]);
1278 if (daughterIndex >=
int(particle->getNDaughters()))
1280 const Particle* dau = particle->getDaughter(daughterIndex);
1282 int grandDaughterIndex = std::lround(arguments[1]);
1283 if (grandDaughterIndex >=
int(dau->getNDaughters()))
1286 ROOT::Math::XYZVector boost = dau->get4Vector().BoostToCM();
1288 ROOT::Math::PxPyPzEVector motherMomentum = - particle->get4Vector();
1289 motherMomentum = ROOT::Math::Boost(boost) * motherMomentum;
1291 ROOT::Math::PxPyPzEVector grandDaughterMomentum = dau->getDaughter(grandDaughterIndex)->get4Vector();
1292 grandDaughterMomentum = ROOT::Math::Boost(boost) * grandDaughterMomentum;
1294 return ROOT::Math::VectorUtil::Angle(motherMomentum, grandDaughterMomentum);
1297 B2FATAL(
"The variable grandDaughterDecayAngle needs exactly two integers as arguments!");
1303 if (arguments.size() == 2 || arguments.size() == 3) {
1305 auto func = [arguments](
const Particle * particle) ->
double {
1306 if (particle ==
nullptr)
1309 std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1313 if (particle->getParticleSource() == Particle::EParticleSourceObject::c_MCParticle)
1315 for (
auto& generalizedIndex : arguments) {
1316 const MCParticle* mcPart = particle->getMCParticle();
1317 if (mcPart ==
nullptr)
1319 const MCParticle* dauMcPart = mcPart->getParticleFromGeneralizedIndexString(generalizedIndex);
1320 if (dauMcPart ==
nullptr)
1323 pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1327 for (
auto& generalizedIndex : arguments) {
1328 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1329 if (dauPart ==
nullptr)
1332 const MCParticle* dauMcPart = dauPart->getMCParticle();
1333 if (dauMcPart ==
nullptr)
1336 pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1341 if (pDaus.size() == 2)
1342 return ROOT::Math::VectorUtil::Angle(pDaus[0], pDaus[1]);
1344 return ROOT::Math::VectorUtil::Angle(pDaus[2], pDaus[0] + pDaus[1]);
1348 B2FATAL(
"Wrong number of arguments for meta function mcDaughterAngle");
1352 double daughterClusterAngleInBetween(
const Particle* particle,
const std::vector<double>& daughterIndices)
1354 if (daughterIndices.size() == 2) {
1355 int daughterIndexi = std::lround(daughterIndices[0]);
1356 int daughterIndexj = std::lround(daughterIndices[1]);
1357 if (std::max(daughterIndexi, daughterIndexj) >=
int(particle->getNDaughters())) {
1360 const ECLCluster* clusteri = particle->getDaughter(daughterIndexi)->getECLCluster();
1361 const ECLCluster* clusterj = particle->getDaughter(daughterIndexj)->getECLCluster();
1362 if (clusteri and clusterj) {
1366 ClusterUtils clusutils;
1367 ROOT::Math::PxPyPzEVector pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit));
1368 ROOT::Math::PxPyPzEVector pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit));
1369 return ROOT::Math::VectorUtil::Angle(pi, pj);
1373 }
else if (daughterIndices.size() == 3) {
1374 int daughterIndexi = std::lround(daughterIndices[0]);
1375 int daughterIndexj = std::lround(daughterIndices[1]);
1376 int daughterIndexk = std::lround(daughterIndices[2]);
1377 if (std::max(std::max(daughterIndexi, daughterIndexj), daughterIndexk) >=
int(particle->getNDaughters())) {
1380 const ECLCluster* clusteri = (particle->getDaughter(daughterIndices[0]))->getECLCluster();
1381 const ECLCluster* clusterj = (particle->getDaughter(daughterIndices[1]))->getECLCluster();
1382 const ECLCluster* clusterk = (particle->getDaughter(daughterIndices[2]))->getECLCluster();
1383 if (clusteri and clusterj and clusterk) {
1388 ClusterUtils clusutils;
1389 ROOT::Math::PxPyPzEVector pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit));
1390 ROOT::Math::PxPyPzEVector pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit));
1391 ROOT::Math::PxPyPzEVector pk = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterk, clusterkBit));
1392 return ROOT::Math::VectorUtil::Angle(pk, pi + pj);
1397 B2FATAL(
"Wrong number of arguments for daughterClusterAngleInBetween!");
1403 if (arguments.size() > 1) {
1404 auto func = [arguments](
const Particle * particle) ->
double {
1406 ROOT::Math::PxPyPzEVector pSum;
1408 for (
auto& generalizedIndex : arguments)
1410 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1412 pSum += frame.getMomentum(dauPart);
1421 B2FATAL(
"Wrong number of arguments for meta function daughterInvM. At least two integers are needed.");
1427 if (arguments.size() == 2) {
1432 }
catch (std::invalid_argument&) {
1433 B2FATAL(
"Second argument of modulo meta function must be integer!");
1435 auto func = [var, divideBy](
const Particle * particle) ->
int {
1436 auto var_result = var->function(particle);
1437 if (std::holds_alternative<double>(var_result))
1439 return int(std::get<double>(var_result)) % divideBy;
1440 }
else if (std::holds_alternative<int>(var_result))
1442 return std::get<int>(var_result) % divideBy;
1443 }
else if (std::holds_alternative<bool>(var_result))
1445 return int(std::get<bool>(var_result)) % divideBy;
1450 B2FATAL(
"Wrong number of arguments for meta function modulo");
1456 if (arguments.size() == 1) {
1459 auto func = [var](
const Particle * particle) ->
bool {
return std::isnan(std::get<double>(var->function(particle))); };
1462 B2FATAL(
"Wrong number of arguments for meta function isNAN");
1468 if (arguments.size() == 2) {
1470 double defaultOutput;
1473 }
catch (std::invalid_argument&) {
1474 B2FATAL(
"The second argument of ifNANgiveX meta function must be a number!");
1476 auto func = [var, defaultOutput](
const Particle * particle) ->
double {
1477 double output = std::get<double>(var->function(particle));
1478 if (std::isnan(output))
return defaultOutput;
1483 B2FATAL(
"Wrong number of arguments for meta function ifNANgiveX");
1489 if (arguments.size() == 1) {
1492 auto func = [var](
const Particle * particle) ->
bool {
return std::isinf(std::get<double>(var->function(particle))); };
1495 B2FATAL(
"Wrong number of arguments for meta function isInfinity");
1501 if (arguments.size() >= 2) {
1507 for (
size_t i = 1; i < arguments.size(); ++i) {
1510 }
catch (std::invalid_argument&) {
1511 B2FATAL(
"The input flags to meta function unmask() should be integer!");
1517 auto func = [var, finalMask](
const Particle * particle) ->
double {
1519 auto var_result = var->function(particle);
1520 if (std::holds_alternative<double>(var_result))
1523 if (std::isnan(std::get<double>(var_result))) {
1526 value = int(std::get<double>(var_result));
1527 }
else if (std::holds_alternative<int>(var_result))
1529 value = std::get<int>(var_result);
1533 value &= (~finalMask);
1540 B2FATAL(
"Meta function unmask needs at least two arguments!");
1546 if (arguments.size() == 3) {
1548 std::string cutString = arguments[0];
1554 auto func = [cut, variableIfTrue, variableIfFalse](
const Particle * particle) ->
double {
1555 if (particle ==
nullptr)
1557 if (cut->check(particle))
1559 auto var_result = variableIfTrue->function(particle);
1560 if (std::holds_alternative<double>(var_result)) {
1561 return std::get<double>(var_result);
1562 }
else if (std::holds_alternative<int>(var_result)) {
1563 return std::get<int>(var_result);
1564 }
else if (std::holds_alternative<bool>(var_result)) {
1565 return std::get<bool>(var_result);
1569 auto var_result = variableIfFalse->function(particle);
1570 if (std::holds_alternative<double>(var_result)) {
1571 return std::get<double>(var_result);
1572 }
else if (std::holds_alternative<int>(var_result)) {
1573 return std::get<int>(var_result);
1574 }
else if (std::holds_alternative<bool>(var_result)) {
1575 return std::get<bool>(var_result);
1582 B2FATAL(
"Wrong number of arguments for meta function conditionalVariableSelector");
1588 if (arguments.size() > 0) {
1589 std::vector<const Variable::Manager::Var*> variables;
1590 for (
auto& argument : arguments)
1593 auto func = [variables, arguments](
const Particle * particle) ->
double {
1594 double pValueProduct = 1.;
1595 for (
auto variable : variables)
1597 double pValue = std::get<double>(variable->function(particle));
1601 pValueProduct *= pValue;
1603 double pValueSum = 1.;
1604 double factorial = 1.;
1605 for (
unsigned int i = 1; i < arguments.size(); ++i)
1608 pValueSum += pow(-std::log(pValueProduct), i) / factorial;
1610 return pValueProduct * pValueSum;
1614 B2FATAL(
"Wrong number of arguments for meta function pValueCombination");
1620 if (arguments.size() == 1) {
1622 auto func = [var](
const Particle * particle) ->
double {
1623 double pValueProduct = 1.;
1624 if (particle->getNDaughters() == 0)
1629 for (
unsigned j = 0; j < particle->getNDaughters(); ++j)
1631 double pValue = std::get<double>(var->function(particle->getDaughter(j)));
1632 if (pValue < 0)
return -1;
1633 else pValueProduct *= pValue;
1636 double pValueSum = 1.;
1637 double factorial = 1.;
1638 for (
unsigned int i = 1; i < particle->getNDaughters(); ++i)
1641 pValueSum += pow(-std::log(pValueProduct), i) / factorial;
1643 return pValueProduct * pValueSum;
1647 B2FATAL(
"Wrong number of arguments for meta function pValueCombinationOfDaughters");
1653 if (arguments.size() == 1) {
1655 auto func = [var](
const Particle * particle) ->
double {
1656 auto var_result = var->function(particle);
1657 if (std::holds_alternative<double>(var_result))
1659 return std::abs(std::get<double>(var_result));
1660 }
else if (std::holds_alternative<int>(var_result))
1662 return std::abs(std::get<int>(var_result));
1667 B2FATAL(
"Wrong number of arguments for meta function abs");
1673 if (arguments.size() == 2) {
1678 B2FATAL(
"One or both of the used variables doesn't exist!");
1680 auto func = [var1, var2](
const Particle * particle) ->
double {
1682 auto var_result1 = var1->function(particle);
1683 auto var_result2 = var2->function(particle);
1684 if (std::holds_alternative<double>(var_result1))
1686 val1 = std::get<double>(var_result1);
1687 }
else if (std::holds_alternative<int>(var_result1))
1689 val1 = std::get<int>(var_result1);
1691 if (std::holds_alternative<double>(var_result2))
1693 val2 = std::get<double>(var_result2);
1694 }
else if (std::holds_alternative<int>(var_result2))
1696 val2 = std::get<int>(var_result2);
1698 return std::max(val1, val2);
1702 B2FATAL(
"Wrong number of arguments for meta function max");
1708 if (arguments.size() == 2) {
1713 B2FATAL(
"One or both of the used variables doesn't exist!");
1715 auto func = [var1, var2](
const Particle * particle) ->
double {
1717 auto var_result1 = var1->function(particle);
1718 auto var_result2 = var2->function(particle);
1719 if (std::holds_alternative<double>(var_result1))
1721 val1 = std::get<double>(var_result1);
1722 }
else if (std::holds_alternative<int>(var_result1))
1724 val1 = std::get<int>(var_result1);
1726 if (std::holds_alternative<double>(var_result2))
1728 val2 = std::get<double>(var_result2);
1729 }
else if (std::holds_alternative<int>(var_result2))
1731 val2 = std::get<int>(var_result2);
1733 return std::min(val1, val2);
1737 B2FATAL(
"Wrong number of arguments for meta function min");
1743 if (arguments.size() == 1) {
1745 auto func = [var](
const Particle * particle) ->
double {
1746 auto var_result = var->function(particle);
1747 if (std::holds_alternative<double>(var_result))
1748 return std::sin(std::get<double>(var_result));
1749 else if (std::holds_alternative<int>(var_result))
1750 return std::sin(std::get<int>(var_result));
1755 B2FATAL(
"Wrong number of arguments for meta function sin");
1761 if (arguments.size() == 1) {
1763 auto func = [var](
const Particle * particle) ->
double {
1764 auto var_result = var->function(particle);
1765 if (std::holds_alternative<double>(var_result))
1766 return std::asin(std::get<double>(var_result));
1767 else if (std::holds_alternative<int>(var_result))
1768 return std::asin(std::get<int>(var_result));
1773 B2FATAL(
"Wrong number of arguments for meta function asin");
1779 if (arguments.size() == 1) {
1781 auto func = [var](
const Particle * particle) ->
double {
1782 auto var_result = var->function(particle);
1783 if (std::holds_alternative<double>(var_result))
1784 return std::cos(std::get<double>(var_result));
1785 else if (std::holds_alternative<int>(var_result))
1786 return std::cos(std::get<int>(var_result));
1791 B2FATAL(
"Wrong number of arguments for meta function cos");
1797 if (arguments.size() == 1) {
1799 auto func = [var](
const Particle * particle) ->
double {
1800 auto var_result = var->function(particle);
1801 if (std::holds_alternative<double>(var_result))
1802 return std::acos(std::get<double>(var_result));
1803 else if (std::holds_alternative<int>(var_result))
1804 return std::acos(std::get<int>(var_result));
1809 B2FATAL(
"Wrong number of arguments for meta function acos");
1815 if (arguments.size() == 1) {
1817 auto func = [var](
const Particle * particle) ->
double {
return std::tan(std::get<double>(var->function(particle))); };
1820 B2FATAL(
"Wrong number of arguments for meta function tan");
1826 if (arguments.size() == 1) {
1828 auto func = [var](
const Particle * particle) ->
double {
return std::atan(std::get<double>(var->function(particle))); };
1831 B2FATAL(
"Wrong number of arguments for meta function atan");
1837 if (arguments.size() == 1) {
1839 auto func = [var](
const Particle * particle) ->
double {
1840 auto var_result = var->function(particle);
1841 if (std::holds_alternative<double>(var_result))
1842 return std::exp(std::get<double>(var_result));
1843 else if (std::holds_alternative<int>(var_result))
1844 return std::exp(std::get<int>(var_result));
1849 B2FATAL(
"Wrong number of arguments for meta function exp");
1855 if (arguments.size() == 1) {
1857 auto func = [var](
const Particle * particle) ->
double {
1858 auto var_result = var->function(particle);
1859 if (std::holds_alternative<double>(var_result))
1860 return std::log(std::get<double>(var_result));
1861 else if (std::holds_alternative<int>(var_result))
1862 return std::log(std::get<int>(var_result));
1867 B2FATAL(
"Wrong number of arguments for meta function log");
1873 if (arguments.size() == 1) {
1875 auto func = [var](
const Particle * particle) ->
double {
1876 auto var_result = var->function(particle);
1877 if (std::holds_alternative<double>(var_result))
1878 return std::log10(std::get<double>(var_result));
1879 else if (std::holds_alternative<int>(var_result))
1880 return std::log10(std::get<int>(var_result));
1885 B2FATAL(
"Wrong number of arguments for meta function log10");
1891 if (arguments.size() == 1) {
1893 auto func = [var](
const Particle * particle) ->
double {
1894 if (particle ==
nullptr)
1897 StoreArray<Particle> particles;
1898 if (!particle->hasExtraInfo(
"original_index"))
1901 auto originalParticle = particles[particle->getExtraInfo(
"original_index")];
1902 if (!originalParticle)
1904 auto var_result = var->function(originalParticle);
1905 if (std::holds_alternative<double>(var_result))
1907 return std::get<double>(var_result);
1908 }
else if (std::holds_alternative<int>(var_result))
1910 return std::get<int>(var_result);
1911 }
else if (std::holds_alternative<bool>(var_result))
1913 return std::get<bool>(var_result);
1918 B2FATAL(
"Wrong number of arguments for meta function originalParticle");
1924 if (arguments.size() == 2) {
1925 auto daughterFunction = convertToDaughterIndex({arguments[0]});
1927 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
1928 if (particle ==
nullptr)
1930 int daughterNumber = std::get<int>(daughterFunction(particle));
1931 if (daughterNumber >=
int(particle->getNDaughters()) or daughterNumber < 0)
1933 auto var_result = var->function(particle->getDaughter(daughterNumber));
1934 if (std::holds_alternative<double>(var_result))
1936 return std::get<double>(var_result);
1937 }
else if (std::holds_alternative<int>(var_result))
1939 return std::get<int>(var_result);
1940 }
else if (std::holds_alternative<bool>(var_result))
1942 return std::get<bool>(var_result);
1947 B2FATAL(
"Wrong number of arguments for meta function daughter");
1953 if (arguments.size() == 2) {
1954 auto daughterFunction = convertToDaughterIndex({arguments[0]});
1956 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
1957 if (particle ==
nullptr)
1959 int daughterNumber = std::get<int>(daughterFunction(particle));
1960 if (daughterNumber >=
int(particle->getNDaughters()) or daughterNumber < 0)
1964 StoreArray<Particle> particles;
1965 if (!particle->getDaughter(daughterNumber)->hasExtraInfo(
"original_index"))
1967 auto originalDaughter = particles[particle->getDaughter(daughterNumber)->getExtraInfo(
"original_index")];
1968 if (!originalDaughter)
1971 auto var_result = var->function(originalDaughter);
1972 if (std::holds_alternative<double>(var_result)) {
1973 return std::get<double>(var_result);
1974 }
else if (std::holds_alternative<int>(var_result)) {
1975 return std::get<int>(var_result);
1976 }
else if (std::holds_alternative<bool>(var_result)) {
1977 return std::get<bool>(var_result);
1983 B2FATAL(
"Wrong number of arguments for meta function daughter");
1989 if (arguments.size() == 1) {
1990 std::string daughterString = arguments[0];
1991 auto func = [daughterString](
const Particle * particle) ->
int {
1992 if (particle ==
nullptr)
1994 int daughterNumber = 0;
1998 }
catch (std::invalid_argument&)
2000 auto daughterFunction = convertToInt({daughterString,
"-1"});
2001 auto daughterVarResult = daughterFunction(particle);
2002 daughterNumber = std::get<int>(daughterVarResult);
2004 return daughterNumber;
2008 B2FATAL(
"Wrong number of arguments for meta function convertToDaughterIndex");
2014 if (arguments.size() == 2) {
2015 auto daughterFunction = convertToDaughterIndex({arguments[0]});
2017 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
2018 if (particle ==
nullptr)
2020 if (particle->getMCParticle())
2022 int daughterNumber = std::get<int>(daughterFunction(particle));
2023 if (daughterNumber >=
int(particle->getMCParticle()->getNDaughters()) or daughterNumber < 0)
2025 Particle tempParticle = Particle(particle->getMCParticle()->getDaughters().at(daughterNumber));
2026 auto var_result = var->function(&tempParticle);
2027 if (std::holds_alternative<double>(var_result)) {
2028 return std::get<double>(var_result);
2029 }
else if (std::holds_alternative<int>(var_result)) {
2030 return std::get<int>(var_result);
2031 }
else if (std::holds_alternative<bool>(var_result)) {
2032 return std::get<bool>(var_result);
2043 B2FATAL(
"Wrong number of arguments for meta function mcDaughter");
2049 if (arguments.size() == 1) {
2051 auto func = [var](
const Particle * particle) ->
double {
2052 if (particle ==
nullptr)
2054 if (particle->getMCParticle())
2056 if (particle->getMCParticle()->getMother() ==
nullptr) {
2059 Particle tempParticle = Particle(particle->getMCParticle()->getMother());
2060 auto var_result = var->function(&tempParticle);
2061 if (std::holds_alternative<double>(var_result)) {
2062 return std::get<double>(var_result);
2063 }
else if (std::holds_alternative<int>(var_result)) {
2064 return std::get<int>(var_result);
2065 }
else if (std::holds_alternative<bool>(var_result)) {
2066 return std::get<bool>(var_result);
2075 B2FATAL(
"Wrong number of arguments for meta function mcMother");
2081 if (arguments.size() == 2) {
2082 int particleNumber = 0;
2085 }
catch (std::invalid_argument&) {
2086 B2FATAL(
"First argument of genParticle meta function must be integer!");
2090 auto func = [var, particleNumber](
const Particle*) ->
double {
2091 StoreArray<MCParticle> mcParticles(
"MCParticles");
2092 if (particleNumber >= mcParticles.getEntries())
2097 MCParticle* mcParticle = mcParticles[particleNumber];
2098 Particle part = Particle(mcParticle);
2099 auto var_result = var->function(&part);
2100 if (std::holds_alternative<double>(var_result))
2102 return std::get<double>(var_result);
2103 }
else if (std::holds_alternative<int>(var_result))
2105 return std::get<int>(var_result);
2106 }
else if (std::holds_alternative<bool>(var_result))
2108 return std::get<bool>(var_result);
2113 B2FATAL(
"Wrong number of arguments for meta function genParticle");
2119 if (arguments.size() == 1) {
2122 auto func = [var](
const Particle*) ->
double {
2123 StoreArray<MCParticle> mcParticles(
"MCParticles");
2124 if (mcParticles.getEntries() == 0)
2129 MCParticle* mcUpsilon4S = mcParticles[0];
2130 if (mcUpsilon4S->isInitial()) mcUpsilon4S = mcParticles[2];
2131 if (mcUpsilon4S->getPDG() != 300553)
2136 Particle upsilon4S = Particle(mcUpsilon4S);
2137 auto var_result = var->function(&upsilon4S);
2138 if (std::holds_alternative<double>(var_result))
2140 return std::get<double>(var_result);
2141 }
else if (std::holds_alternative<int>(var_result))
2143 return std::get<int>(var_result);
2144 }
else if (std::holds_alternative<bool>(var_result))
2146 return std::get<bool>(var_result);
2151 B2FATAL(
"Wrong number of arguments for meta function genUpsilon4S");
2157 if (arguments.size() == 4) {
2158 std::string listName = arguments[0];
2159 std::string rankedVariableName = arguments[1];
2160 std::string returnVariableName = arguments[2];
2161 std::string extraInfoName = rankedVariableName +
"_rank";
2165 }
catch (std::invalid_argument&) {
2166 B2ERROR(
"3rd argument of getVariableByRank meta function (Rank) must be an integer!");
2171 auto func = [var, rank, extraInfoName, listName](
const Particle*)->
double {
2172 StoreObjPtr<ParticleList> list(listName);
2174 const unsigned int numParticles = list->getListSize();
2175 for (
unsigned int i = 0; i < numParticles; i++)
2177 const Particle* p = list->getParticle(i);
2178 if (p->getExtraInfo(extraInfoName) == rank) {
2179 auto var_result = var->function(p);
2180 if (std::holds_alternative<double>(var_result)) {
2181 return std::get<double>(var_result);
2182 }
else if (std::holds_alternative<int>(var_result)) {
2183 return std::get<int>(var_result);
2184 }
else if (std::holds_alternative<bool>(var_result)) {
2185 return std::get<bool>(var_result);
2190 return std::numeric_limits<double>::signaling_NaN();
2194 B2FATAL(
"Wrong number of arguments for meta function getVariableByRank");
2200 if (arguments.size() == 1 or arguments.size() == 2) {
2202 std::string listName = arguments[0];
2203 std::string cutString =
"";
2205 if (arguments.size() == 2) {
2206 cutString = arguments[1];
2211 auto func = [listName, cut](
const Particle*) ->
int {
2213 StoreObjPtr<ParticleList> list(listName);
2215 for (
unsigned int i = 0; i < list->getListSize(); i++)
2217 const Particle* particle = list->getParticle(i);
2218 if (cut->check(particle)) {
2226 B2FATAL(
"Wrong number of arguments for meta function countInList");
2232 if (arguments.size() == 2 or arguments.size() == 3) {
2234 std::string roeListName = arguments[0];
2235 std::string cutString = arguments[1];
2237 if (arguments.size() == 2) {
2238 B2INFO(
"Use pdgCode of electron as default in meta variable veto, other arguments: " << roeListName <<
", " << cutString);
2242 }
catch (std::invalid_argument&) {
2243 B2FATAL(
"Third argument of veto meta function must be integer!");
2250 auto func = [roeListName, cut, pdgCode, flavourType](
const Particle * particle) ->
bool {
2251 StoreObjPtr<ParticleList> roeList(roeListName);
2252 ROOT::Math::PxPyPzEVector vec = particle->get4Vector();
2253 for (
unsigned int i = 0; i < roeList->getListSize(); i++)
2255 const Particle* roeParticle = roeList->getParticle(i);
2256 if (not particle->overlapsWith(roeParticle)) {
2257 ROOT::Math::PxPyPzEVector tempCombination = roeParticle->get4Vector() + vec;
2258 std::vector<int> indices = { particle->getArrayIndex(), roeParticle->getArrayIndex() };
2259 Particle tempParticle = Particle(tempCombination, pdgCode, flavourType, indices, particle->getArrayPointer());
2260 if (cut->check(&tempParticle)) {
2269 B2FATAL(
"Wrong number of arguments for meta function veto");
2275 if (arguments.size() == 1) {
2276 std::string cutString = arguments[0];
2278 auto func = [cut](
const Particle * particle) ->
int {
2280 for (
auto& daughter : particle->getDaughters())
2282 if (cut->check(daughter))
2289 B2FATAL(
"Wrong number of arguments for meta function countDaughters");
2295 if (arguments.size() == 1) {
2296 std::string cutString = arguments[0];
2298 auto func = [cut](
const Particle * particle) ->
int {
2300 std::vector<const Particle*> fspDaughters;
2301 particle->fillFSPDaughters(fspDaughters);
2304 for (
auto& daughter : fspDaughters)
2306 if (cut->check(daughter))
2313 B2FATAL(
"Wrong number of arguments for meta function countFSPDaughters");
2319 if (arguments.size() == 1) {
2320 std::string cutString = arguments[0];
2322 auto func = [cut](
const Particle * particle) ->
int {
2324 std::vector<const Particle*> allDaughters;
2325 particle->fillAllDaughters(allDaughters);
2328 for (
auto& daughter : allDaughters)
2330 if (cut->check(daughter))
2337 B2FATAL(
"Wrong number of arguments for meta function countDescendants");
2341 Manager::FunctionPtr numberOfNonOverlappingParticles(
const std::vector<std::string>& arguments)
2344 auto func = [arguments](
const Particle * particle) ->
int {
2346 int _numberOfNonOverlappingParticles = 0;
2347 for (
const auto& listName : arguments)
2349 StoreObjPtr<ParticleList> list(listName);
2350 if (not list.isValid()) {
2351 B2FATAL(
"Invalid list named " << listName <<
" encountered in numberOfNonOverlappingParticles.");
2353 for (
unsigned int i = 0; i < list->getListSize(); i++) {
2354 const Particle* p = list->getParticle(i);
2355 if (not particle->overlapsWith(p)) {
2356 _numberOfNonOverlappingParticles++;
2360 return _numberOfNonOverlappingParticles;
2367 void appendDaughtersRecursive(Particle* mother, StoreArray<Particle>& container)
2370 auto* mcmother = mother->getRelated<MCParticle>();
2375 for (
auto* mcdaughter : mcmother->getDaughters()) {
2377 Particle tmp_daughter(mcdaughter);
2378 Particle* new_daughter = container.appendNew(tmp_daughter);
2379 new_daughter->addRelationTo(mcdaughter);
2380 mother->appendDaughter(new_daughter,
false);
2382 if (mcdaughter->getNDaughters() > 0)
2383 appendDaughtersRecursive(new_daughter, container);
2389 if (arguments.size() == 1) {
2391 auto func = [var](
const Particle * particle) ->
double {
2392 const MCParticle* mcp = particle->getMCParticle();
2397 StoreArray<Particle> tempParticles(
"tempParticles");
2398 tempParticles.clear();
2399 Particle tmpPart(mcp);
2400 Particle* newPart = tempParticles.appendNew(tmpPart);
2401 newPart->addRelationTo(mcp);
2403 appendDaughtersRecursive(newPart, tempParticles);
2405 auto var_result = var->function(newPart);
2406 if (std::holds_alternative<double>(var_result))
2408 return std::get<double>(var_result);
2409 }
else if (std::holds_alternative<int>(var_result))
2411 return std::get<int>(var_result);
2412 }
else if (std::holds_alternative<bool>(var_result))
2414 return std::get<bool>(var_result);
2419 B2FATAL(
"Wrong number of arguments for meta function matchedMC");
2425 if (arguments.size() == 1) {
2428 auto func = [var](
const Particle * particle) ->
double {
2430 const ECLCluster* cluster = particle->getECLCluster();
2433 auto mcps = cluster->getRelationsTo<MCParticle>();
2436 std::vector<std::pair<double, int>> weightsAndIndices;
2437 for (
unsigned int i = 0; i < mcps.size(); ++i)
2438 weightsAndIndices.emplace_back(mcps.weight(i), i);
2441 std::sort(weightsAndIndices.begin(), weightsAndIndices.end(),
2442 ValueIndexPairSorting::higherPair<
decltype(weightsAndIndices)::value_type>);
2445 const MCParticle* mcp = mcps.object(weightsAndIndices[0].second);
2447 StoreArray<Particle> tempParticles(
"tempParticles");
2448 tempParticles.clear();
2449 Particle tmpPart(mcp);
2450 Particle* newPart = tempParticles.appendNew(tmpPart);
2451 newPart->addRelationTo(mcp);
2453 appendDaughtersRecursive(newPart, tempParticles);
2455 auto var_result = var->function(newPart);
2456 if (std::holds_alternative<double>(var_result))
2458 return std::get<double>(var_result);
2459 }
else if (std::holds_alternative<int>(var_result))
2461 return std::get<int>(var_result);
2462 }
else if (std::holds_alternative<bool>(var_result))
2464 return std::get<bool>(var_result);
2473 B2FATAL(
"Wrong number of arguments for meta function clusterBestMatchedMCParticle");
2479 if (arguments.size() == 1) {
2482 auto func = [var](
const Particle * particle) ->
double {
2484 const ECLCluster* cluster = particle->getECLCluster();
2487 auto mcps = cluster->getRelationsTo<MCParticle>();
2490 std::map<int, double> mapMCParticleIndxAndWeight;
2491 getKlongWeightMap(particle, mapMCParticleIndxAndWeight);
2494 if (mapMCParticleIndxAndWeight.size() == 0)
2498 auto maxMap = std::max_element(mapMCParticleIndxAndWeight.begin(), mapMCParticleIndxAndWeight.end(),
2499 [](
const auto & x,
const auto & y) { return x.second < y.second; }
2502 StoreArray<MCParticle> mcparticles;
2503 const MCParticle* mcKlong = mcparticles[maxMap->first];
2505 Particle tmpPart(mcKlong);
2506 auto var_result = var->function(&tmpPart);
2507 if (std::holds_alternative<double>(var_result))
2509 return std::get<double>(var_result);
2510 }
else if (std::holds_alternative<int>(var_result))
2512 return std::get<int>(var_result);
2513 }
else if (std::holds_alternative<bool>(var_result))
2515 return std::get<bool>(var_result);
2524 B2FATAL(
"Wrong number of arguments for meta function clusterBestMatchedMCKlong");
2528 double matchedMCHasPDG(
const Particle* particle,
const std::vector<double>& pdgCode)
2530 if (pdgCode.size() != 1) {
2531 B2FATAL(
"Too many arguments provided to matchedMCHasPDG!");
2533 int inputPDG = std::lround(pdgCode[0]);
2535 const MCParticle* mcp = particle->getMCParticle();
2539 return std::abs(mcp->getPDG()) == inputPDG;
2544 if (arguments.size() == 1) {
2545 std::string listName = arguments[0];
2546 auto func = [listName](
const Particle * particle) ->
double {
2549 StoreObjPtr<ParticleList> listOfParticles(listName);
2551 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalEnergyOfParticlesInList");
2552 double totalEnergy = 0;
2553 int nParticles = listOfParticles->getListSize();
2554 for (
int i = 0; i < nParticles; i++)
2556 const Particle* part = listOfParticles->getParticle(i);
2558 totalEnergy += frame.getMomentum(part).E();
2565 B2FATAL(
"Wrong number of arguments for meta function totalEnergyOfParticlesInList");
2571 if (arguments.size() == 1) {
2572 std::string listName = arguments[0];
2573 auto func = [listName](
const Particle*) ->
double {
2574 StoreObjPtr<ParticleList> listOfParticles(listName);
2576 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPxOfParticlesInList");
2578 int nParticles = listOfParticles->getListSize();
2580 for (
int i = 0; i < nParticles; i++)
2582 const Particle* part = listOfParticles->getParticle(i);
2583 totalPx += frame.getMomentum(part).Px();
2589 B2FATAL(
"Wrong number of arguments for meta function totalPxOfParticlesInList");
2595 if (arguments.size() == 1) {
2596 std::string listName = arguments[0];
2597 auto func = [listName](
const Particle*) ->
double {
2598 StoreObjPtr<ParticleList> listOfParticles(listName);
2600 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPyOfParticlesInList");
2602 int nParticles = listOfParticles->getListSize();
2604 for (
int i = 0; i < nParticles; i++)
2606 const Particle* part = listOfParticles->getParticle(i);
2607 totalPy += frame.getMomentum(part).Py();
2613 B2FATAL(
"Wrong number of arguments for meta function totalPyOfParticlesInList");
2619 if (arguments.size() == 1) {
2620 std::string listName = arguments[0];
2621 auto func = [listName](
const Particle*) ->
double {
2622 StoreObjPtr<ParticleList> listOfParticles(listName);
2624 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPzOfParticlesInList");
2626 int nParticles = listOfParticles->getListSize();
2628 for (
int i = 0; i < nParticles; i++)
2630 const Particle* part = listOfParticles->getParticle(i);
2631 totalPz += frame.getMomentum(part).Pz();
2637 B2FATAL(
"Wrong number of arguments for meta function totalPzOfParticlesInList");
2643 if (arguments.size() > 0) {
2645 auto func = [arguments](
const Particle * particle) ->
double {
2647 ROOT::Math::PxPyPzEVector total4Vector;
2649 std::vector<Particle*> particlePool;
2652 for (
const auto& argument : arguments)
2654 StoreObjPtr <ParticleList> listOfParticles(argument);
2656 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << argument <<
" given to invMassInLists");
2657 int nParticles = listOfParticles->getListSize();
2658 for (
int i = 0; i < nParticles; i++) {
2659 bool overlaps =
false;
2660 Particle* part = listOfParticles->getParticle(i);
2661 for (
auto poolPart : particlePool) {
2662 if (part->overlapsWith(poolPart)) {
2668 total4Vector += part->get4Vector();
2669 particlePool.push_back(part);
2673 double invariantMass = total4Vector.M();
2674 return invariantMass;
2679 B2FATAL(
"Wrong number of arguments for meta function invMassInLists");
2683 Manager::FunctionPtr totalECLEnergyOfParticlesInList(
const std::vector<std::string>& arguments)
2685 if (arguments.size() == 1) {
2686 std::string listName = arguments[0];
2687 auto func = [listName](
const Particle * particle) ->
double {
2690 StoreObjPtr<ParticleList> listOfParticles(listName);
2692 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalEnergyOfParticlesInList");
2693 double totalEnergy = 0;
2694 int nParticles = listOfParticles->getListSize();
2695 for (
int i = 0; i < nParticles; i++)
2697 const Particle* part = listOfParticles->getParticle(i);
2698 const ECLCluster* cluster = part->getECLCluster();
2700 if (cluster !=
nullptr) {
2701 totalEnergy += cluster->getEnergy(clusterHypothesis);
2709 B2FATAL(
"Wrong number of arguments for meta function totalECLEnergyOfParticlesInList");
2715 if (arguments.size() == 1) {
2716 std::string listName = arguments[0];
2717 auto func = [listName](
const Particle*) ->
double {
2718 StoreObjPtr<ParticleList> listOfParticles(listName);
2720 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to maxPtInList");
2721 int nParticles = listOfParticles->getListSize();
2724 for (
int i = 0; i < nParticles; i++)
2726 const Particle* part = listOfParticles->getParticle(i);
2727 const double Pt = frame.getMomentum(part).Pt();
2728 if (Pt > maxPt) maxPt = Pt;
2734 B2FATAL(
"Wrong number of arguments for meta function maxPtInList");
2738 Manager::FunctionPtr eclClusterTrackMatchedWithCondition(
const std::vector<std::string>& arguments)
2740 if (arguments.size() <= 1) {
2742 std::string cutString;
2743 if (arguments.size() == 1)
2744 cutString = arguments[0];
2746 auto func = [cut](
const Particle * particle) ->
double {
2748 if (particle ==
nullptr)
2751 const ECLCluster* cluster = particle->getECLCluster();
2755 auto tracks = cluster->getRelationsFrom<Track>();
2757 for (
const auto& track : tracks) {
2760 if (cut->check(&trackParticle))
2769 B2FATAL(
"Wrong number of arguments for meta function eclClusterSpecialTrackMatched");
2775 if (arguments.size() == 2) {
2776 std::string listName = arguments[0];
2779 auto func = [listName, var](
const Particle*) ->
double {
2780 StoreObjPtr<ParticleList> listOfParticles(listName);
2782 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to averageValueInList");
2783 int nParticles = listOfParticles->getListSize();
2784 if (nParticles == 0)
2789 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2791 for (
int i = 0; i < nParticles; i++) {
2792 average += std::get<double>(var->function(listOfParticles->getParticle(i))) / nParticles;
2794 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2796 for (
int i = 0; i < nParticles; i++) {
2797 average += std::get<int>(var->function(listOfParticles->getParticle(i))) / nParticles;
2804 B2FATAL(
"Wrong number of arguments for meta function averageValueInList");
2810 if (arguments.size() == 2) {
2811 std::string listName = arguments[0];
2814 auto func = [listName, var](
const Particle*) ->
double {
2815 StoreObjPtr<ParticleList> listOfParticles(listName);
2817 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to medianValueInList");
2818 int nParticles = listOfParticles->getListSize();
2819 if (nParticles == 0)
2823 std::vector<double> valuesInList;
2824 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2826 for (
int i = 0; i < nParticles; i++) {
2827 valuesInList.push_back(std::get<double>(var->function(listOfParticles->getParticle(i))));
2829 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2831 for (
int i = 0; i < nParticles; i++) {
2832 valuesInList.push_back(std::get<int>(var->function(listOfParticles->getParticle(i))));
2835 std::sort(valuesInList.begin(), valuesInList.end());
2836 if (nParticles % 2 != 0)
2838 return valuesInList[nParticles / 2];
2841 return 0.5 * (valuesInList[nParticles / 2] + valuesInList[nParticles / 2 - 1]);
2846 B2FATAL(
"Wrong number of arguments for meta function medianValueInList");
2852 if (arguments.size() == 2) {
2853 std::string listName = arguments[0];
2856 auto func = [listName, var](
const Particle*) ->
double {
2857 StoreObjPtr<ParticleList> listOfParticles(listName);
2859 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to sumValueInList");
2860 int nParticles = listOfParticles->getListSize();
2861 if (nParticles == 0)
2866 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2868 for (
int i = 0; i < nParticles; i++) {
2869 sum += std::get<double>(var->function(listOfParticles->getParticle(i)));
2871 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2873 for (
int i = 0; i < nParticles; i++) {
2874 sum += std::get<int>(var->function(listOfParticles->getParticle(i)));
2881 B2FATAL(
"Wrong number of arguments for meta function sumValueInList");
2887 if (arguments.size() == 2) {
2888 std::string listName = arguments[0];
2891 auto func = [listName, var](
const Particle*) ->
double {
2892 StoreObjPtr<ParticleList> listOfParticles(listName);
2894 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to productValueInList");
2895 int nParticles = listOfParticles->getListSize();
2896 if (nParticles == 0)
2901 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2903 for (
int i = 0; i < nParticles; i++) {
2904 product *= std::get<double>(var->function(listOfParticles->getParticle(i)));
2906 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2908 for (
int i = 0; i < nParticles; i++) {
2909 product *= std::get<int>(var->function(listOfParticles->getParticle(i)));
2916 B2FATAL(
"Wrong number of arguments for meta function productValueInList");
2923 if (arguments.size() != 1)
2924 B2FATAL(
"Wrong number of arguments for meta function angleToClosestInList");
2926 std::string listname = arguments[0];
2928 auto func = [listname](
const Particle * particle) ->
double {
2930 StoreObjPtr<ParticleList> list(listname);
2931 if (not list.isValid())
2932 B2FATAL(
"Invalid particle list name " << listname <<
" given to angleToClosestInList");
2935 if (list->getListSize() == 0)
2940 const auto p_this = frame.getMomentum(particle);
2943 double minAngle = 2 * M_PI;
2944 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2946 const Particle* compareme = list->getParticle(i);
2947 const auto p_compare = frame.getMomentum(compareme);
2948 double angle = ROOT::Math::VectorUtil::Angle(p_compare, p_this);
2949 if (minAngle > angle) minAngle = angle;
2959 if (arguments.size() != 2)
2960 B2FATAL(
"Wrong number of arguments for meta function closestInList");
2962 std::string listname = arguments[0];
2967 auto func = [listname, var](
const Particle * particle) ->
double {
2969 StoreObjPtr<ParticleList> list(listname);
2970 if (not list.isValid())
2971 B2FATAL(
"Invalid particle list name " << listname <<
" given to closestInList");
2975 const auto p_this = frame.getMomentum(particle);
2978 double minAngle = 2 * M_PI;
2980 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2982 const Particle* compareme = list->getParticle(i);
2983 const auto p_compare = frame.getMomentum(compareme);
2984 double angle = ROOT::Math::VectorUtil::Angle(p_compare, p_this);
2985 if (minAngle > angle) {
2993 auto var_result = var->function(list->getParticle(iClosest));
2994 if (std::holds_alternative<double>(var_result))
2996 return std::get<double>(var_result);
2997 }
else if (std::holds_alternative<int>(var_result))
2999 return std::get<int>(var_result);
3000 }
else if (std::holds_alternative<bool>(var_result))
3002 return std::get<bool>(var_result);
3011 if (arguments.size() != 1)
3012 B2FATAL(
"Wrong number of arguments for meta function angleToMostB2BInList");
3014 std::string listname = arguments[0];
3016 auto func = [listname](
const Particle * particle) ->
double {
3018 StoreObjPtr<ParticleList> list(listname);
3019 if (not list.isValid())
3020 B2FATAL(
"Invalid particle list name " << listname <<
" given to angleToMostB2BInList");
3023 if (list->getListSize() == 0)
3028 const auto p_this = frame.getMomentum(particle);
3032 double maxAngle = 0;
3033 for (
unsigned int i = 0; i < list->getListSize(); ++i)
3035 const Particle* compareme = list->getParticle(i);
3036 const auto p_compare = frame.getMomentum(compareme);
3037 double angle = ROOT::Math::VectorUtil::Angle(p_compare, p_this);
3038 if (maxAngle < angle) maxAngle = angle;
3048 if (arguments.size() != 1)
3049 B2FATAL(
"Wrong number of arguments for meta function deltaPhiToMostB2BPhiInList");
3051 std::string listname = arguments[0];
3053 auto func = [listname](
const Particle * particle) ->
double {
3055 StoreObjPtr<ParticleList> list(listname);
3056 if (not list.isValid())
3057 B2FATAL(
"Invalid particle list name " << listname <<
" given to deltaPhiToMostB2BPhiInList");
3060 if (list->getListSize() == 0)
3065 const auto phi_this = frame.getMomentum(particle).Phi();
3068 double maxAngle = 0;
3069 for (
unsigned int i = 0; i < list->getListSize(); ++i)
3071 const Particle* compareme = list->getParticle(i);
3072 const auto phi_compare = frame.getMomentum(compareme).Phi();
3073 double angle = std::abs(phi_compare - phi_this);
3074 if (angle > M_PI) {angle = 2 * M_PI - angle;}
3075 if (maxAngle < angle) maxAngle = angle;
3085 if (arguments.size() != 2)
3086 B2FATAL(
"Wrong number of arguments for meta function mostB2BInList");
3088 std::string listname = arguments[0];
3093 auto func = [listname, var](
const Particle * particle) ->
double {
3095 StoreObjPtr<ParticleList> list(listname);
3096 if (not list.isValid())
3097 B2FATAL(
"Invalid particle list name " << listname <<
" given to mostB2BInList");
3101 const auto p_this = frame.getMomentum(particle);
3105 double maxAngle = -1.0;
3107 for (
unsigned int i = 0; i < list->getListSize(); ++i)
3109 const Particle* compareme = list->getParticle(i);
3110 const auto p_compare = frame.getMomentum(compareme);
3111 double angle = ROOT::Math::VectorUtil::Angle(p_compare, p_this);
3112 if (maxAngle < angle) {
3120 auto var_result = var->function(list->getParticle(iMostB2B));
3121 if (std::holds_alternative<double>(var_result))
3123 return std::get<double>(var_result);
3124 }
else if (std::holds_alternative<int>(var_result))
3126 return std::get<int>(var_result);
3127 }
else if (std::holds_alternative<bool>(var_result))
3129 return std::get<bool>(var_result);
3137 if (arguments.size() == 1) {
3138 std::string listName = arguments[0];
3139 auto func = [listName](
const Particle*) ->
double {
3140 StoreObjPtr<ParticleList> listOfParticles(listName);
3142 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to maxOpeningAngleInList");
3143 int nParticles = listOfParticles->getListSize();
3148 double maxOpeningAngle = -1;
3149 for (
int i = 0; i < nParticles; i++)
3151 ROOT::Math::PxPyPzEVector v1 = frame.getMomentum(listOfParticles->getParticle(i));
3152 for (
int j = i + 1; j < nParticles; j++) {
3153 ROOT::Math::PxPyPzEVector v2 = frame.getMomentum(listOfParticles->getParticle(j));
3154 const double angle = ROOT::Math::VectorUtil::Angle(v1, v2);
3155 if (angle > maxOpeningAngle) maxOpeningAngle = angle;
3158 return maxOpeningAngle;
3162 B2FATAL(
"Wrong number of arguments for meta function maxOpeningAngleInList");
3169 if (arguments.size() >= 2) {
3174 auto func = [var, arguments](
const Particle * particle) ->
double {
3175 if (particle ==
nullptr)
3177 B2WARNING(
"Trying to access a daughter that does not exist. Skipping");
3183 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
3187 for (
unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++)
3189 auto generalizedIndex = arguments[iCoord];
3190 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
3192 pSum += frame.getMomentum(dauPart);
3194 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
3200 Particle sumOfDaughters(pSum, 100);
3202 auto var_result = var->function(&sumOfDaughters);
3204 if (std::holds_alternative<double>(var_result))
3206 return std::get<double>(var_result);
3207 }
else if (std::holds_alternative<int>(var_result))
3209 return std::get<int>(var_result);
3210 }
else if (std::holds_alternative<bool>(var_result))
3212 return std::get<bool>(var_result);
3217 B2FATAL(
"Wrong number of arguments for meta function daughterCombination");
3220 Manager::FunctionPtr useAlternativeDaughterHypothesis(
const std::vector<std::string>& arguments)
3233 if (arguments.size() >= 2) {
3244 std::unordered_map<unsigned int, int> mapOfReplacedDaughters;
3247 for (
unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++) {
3248 auto replacedDauString = arguments[iCoord];
3250 std::vector<std::string> indexAndMass;
3251 boost::split(indexAndMass, replacedDauString, boost::is_any_of(
":"));
3254 if (indexAndMass.size() > 2) {
3255 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 "
3256 << replacedDauString <<
", while a correct syntax looks like 0:K+.");
3260 if (indexAndMass.size() < 2) {
3261 B2WARNING(
"The string indicating which daughter's mass should be replaced contains only one colon-separated element instead of two. The offending string is "
3262 << replacedDauString <<
", while a correct syntax looks like 0:K+.");
3270 }
catch (std::invalid_argument&) {
3271 B2FATAL(
"Found the string " << indexAndMass[0] <<
"instead of a daughter index.");
3275 TParticlePDG* particlePDG = TDatabasePDG::Instance()->GetParticle(indexAndMass[1].c_str());
3277 B2WARNING(
"Particle not in evt.pdl file! " << indexAndMass[1]);
3282 int pdgCode = particlePDG->PdgCode();
3283 mapOfReplacedDaughters[dauIndex] = pdgCode;
3287 if (mapOfReplacedDaughters.size() != arguments.size() - 1)
3288 B2FATAL(
"Overlapped daughter's index is detected in the meta-variable useAlternativeDaughterHypothesis");
3296 auto func = [var, mapOfReplacedDaughters](
const Particle * particle) ->
double {
3297 if (particle ==
nullptr)
3299 B2WARNING(
"Trying to access a particle that does not exist. Skipping");
3309 ROOT::Math::PxPyPzMVector pSum(0, 0, 0, 0);
3311 for (
unsigned int iDau = 0; iDau < particle->getNDaughters(); iDau++)
3313 const Particle* dauPart = particle->getDaughter(iDau);
3315 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << iDau);
3319 ROOT::Math::PxPyPzMVector dauMom = ROOT::Math::PxPyPzMVector(frame.getMomentum(dauPart));
3323 pdgCode = mapOfReplacedDaughters.at(iDau);
3324 }
catch (std::out_of_range&) {
3331 double p_x = dauMom.Px();
3332 double p_y = dauMom.Py();
3333 double p_z = dauMom.Pz();
3334 dauMom.SetCoordinates(p_x, p_y, p_z, TDatabasePDG::Instance()->GetParticle(pdgCode)->Mass());
3335 const_cast<Particle*
>(dummy->getDaughter(iDau))->set4VectorDividingByMomentumScaling(ROOT::Math::PxPyPzEVector(dauMom));
3338 const int charge = dummy->getDaughter(iDau)->getCharge();
3339 if (TDatabasePDG::Instance()->GetParticle(pdgCode)->Charge() / 3.0 == charge)
3340 const_cast<Particle*
>(dummy->getDaughter(iDau))->setPDGCode(pdgCode);
3342 const_cast<Particle*
>(dummy->getDaughter(iDau))->setPDGCode(-1 * pdgCode);
3348 dummy->set4Vector(ROOT::Math::PxPyPzEVector(pSum));
3350 auto var_result = var->function(dummy);
3353 if (std::holds_alternative<double>(var_result))
3355 return std::get<double>(var_result);
3356 }
else if (std::holds_alternative<int>(var_result))
3358 return std::get<int>(var_result);
3359 }
else if (std::holds_alternative<bool>(var_result))
3361 return std::get<bool>(var_result);
3367 B2FATAL(
"Wrong number of arguments for meta function useAlternativeDaughterHypothesis");
3372 if (arguments.size() == 2) {
3374 std::string arg = arguments[0];
3376 TParticlePDG* part = TDatabasePDG::Instance()->GetParticle(arg.c_str());
3378 if (part !=
nullptr) {
3379 pdg_code = std::abs(part->PdgCode());
3383 }
catch (std::exception& e) {}
3386 if (pdg_code == -1) {
3387 B2FATAL(
"Ancestor " + arg +
" is not recognised. Please provide valid PDG code or particle name.");
3390 auto func = [pdg_code, var](
const Particle * particle) ->
double {
3391 const Particle* p = particle;
3393 int ancestor_level = std::get<double>(
Manager::Instance().getVariable(
"hasAncestor(" + std::to_string(pdg_code) +
", 0)")->function(p));
3394 if ((ancestor_level <= 0) or (std::isnan(ancestor_level)))
3399 const MCParticle* i_p = p->getMCParticle();
3401 for (
int a = 0; a < ancestor_level ; a = a + 1)
3403 i_p = i_p->getMother();
3406 StoreArray<Particle> tempParticles(
"tempParticles");
3407 tempParticles.clear();
3409 Particle* newPart = tempParticles.appendNew(m_p);
3410 newPart->addRelationTo(i_p);
3412 appendDaughtersRecursive(newPart, tempParticles);
3414 auto var_result = var->function(newPart);
3415 if (std::holds_alternative<double>(var_result))
3417 return std::get<double>(var_result);
3418 }
else if (std::holds_alternative<int>(var_result))
3420 return std::get<int>(var_result);
3421 }
else if (std::holds_alternative<bool>(var_result))
3423 return std::get<bool>(var_result);
3428 B2FATAL(
"Wrong number of arguments for meta function varForFirstMCAncestorOfType (expected 2: type and variable of interest)");
3434 if (arguments.size() != 1) {
3435 B2FATAL(
"Number of arguments for nTrackFitResults must be 1, particleType or PDGcode");
3438 std::string arg = arguments[0];
3439 TParticlePDG* part = TDatabasePDG::Instance()->GetParticle(arg.c_str());
3441 if (part !=
nullptr) {
3442 absPdg = std::abs(part->PdgCode());
3446 }
catch (std::exception& e) {}
3449 auto func = [absPdg](
const Particle*) ->
int {
3451 Const::ChargedStable type(absPdg);
3452 StoreArray<Track> tracks;
3454 int nTrackFitResults = 0;
3456 for (
const auto& track : tracks)
3458 const TrackFitResult* trackFit = track.getTrackFitResultWithClosestMass(type);
3460 if (!trackFit)
continue;
3461 if (trackFit->getChargeSign() == 0)
continue;
3466 return nTrackFitResults;
3475 if (arguments.size() == 2) {
3478 auto func = [var, default_val](
const Particle * particle) ->
int {
3479 auto var_result = var->function(particle);
3480 if (std::holds_alternative<double>(var_result))
3482 double value = std::get<double>(var_result);
3483 if (value > std::numeric_limits<int>::max())
3484 value = std::numeric_limits<int>::max();
3485 if (value < std::numeric_limits<int>::min())
3486 value = std::numeric_limits<int>::min();
3487 if (std::isnan(value))
3488 value = default_val;
3489 return static_cast<int>(value);
3490 }
else if (std::holds_alternative<int>(var_result))
3491 return std::get<int>(var_result);
3492 else if (std::holds_alternative<bool>(var_result))
3493 return static_cast<int>(std::get<bool>(var_result));
3494 else return default_val;
3498 B2FATAL(
"Wrong number of arguments for meta function int, please provide variable name and replacement value for NaN!");
3502 VARIABLE_GROUP(
"MetaFunctions");
3503 REGISTER_METAVARIABLE(
"nCleanedECLClusters(cut)", nCleanedECLClusters,
3504 "[Eventbased] Returns the number of clean Clusters in the event\n"
3505 "Clean clusters are defined by the clusters which pass the given cut assuming a photon hypothesis.",
3506 Manager::VariableDataType::c_int);
3507 REGISTER_METAVARIABLE(
"nCleanedTracks(cut)", nCleanedTracks,
3508 "[Eventbased] Returns the number of clean Tracks in the event\n"
3509 "Clean tracks are defined by the tracks which pass the given cut assuming a pion hypothesis.", Manager::VariableDataType::c_int);
3510 REGISTER_METAVARIABLE(
"formula(v1 + v2 * [v3 - v4] / v5^v6)", formula, R
"DOCSTRING(
3511Returns the result of the given formula, where v1 to vN are variables or floating
3512point numbers. Currently the only supported operations are addition (``+``),
3513subtraction (``-``), multiplication (``*``), division (``/``) and power (``^``
3514or ``**``). Parenthesis can be in the form of square brackets ``[v1 * v2]``
3515or normal brackets ``(v1 * v2)``. It will work also with variables taking
3516arguments. Operator precedence is taken into account. For example ::
3518 (daughter(0, E) + daughter(1, E))**2 - p**2 + 0.138
3520.. versionchanged:: release-03-00-00
3521 now both, ``[]`` and ``()`` can be used for grouping operations, ``**`` can
3522 be used for exponent and float literals are possible directly in the
3524)DOCSTRING", Manager::VariableDataType::c_double);
3525 REGISTER_METAVARIABLE("useRestFrame(variable)", useRestFrame,
3526 "Returns the value of the variable using the rest frame of the given particle as current reference frame.\n"
3527 "E.g. ``useRestFrame(daughter(0, p))`` returns the total momentum of the first daughter in its mother's rest-frame", Manager::VariableDataType::c_double);
3528 REGISTER_METAVARIABLE(
"useCMSFrame(variable)", useCMSFrame,
3529 "Returns the value of the variable using the CMS frame as current reference frame.\n"
3530 "E.g. ``useCMSFrame(E)`` returns the energy of a particle in the CMS frame.", Manager::VariableDataType::c_double);
3531 REGISTER_METAVARIABLE(
"useLabFrame(variable)", useLabFrame, R
"DOC(
3532Returns the value of ``variable`` in the *lab* frame.
3535 The lab frame is the default reference frame, usually you don't need to use this meta-variable.
3536 E.g. ``useLabFrame(E)`` returns the energy of a particle in the Lab frame, same as just ``E``.
3538Specifying the lab frame is useful in some corner-cases. For example:
3539``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.
3540)DOC", Manager::VariableDataType::c_double);
3541 REGISTER_METAVARIABLE("useTagSideRecoilRestFrame(variable, daughterIndexTagB)", useTagSideRecoilRestFrame,
3542 "Returns the value of the variable in the rest frame of the recoiling particle to the tag side B meson.\n"
3543 "The variable should only be applied to an Upsilon(4S) list.\n"
3544 "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);
3545 REGISTER_METAVARIABLE(
"useParticleRestFrame(variable, particleList)", useParticleRestFrame,
3546 "Returns the value of the variable in the rest frame of the first Particle contained in the given ParticleList.\n"
3547 "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3548 "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3549 "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);
3550 REGISTER_METAVARIABLE(
"useRecoilParticleRestFrame(variable, particleList)", useRecoilParticleRestFrame,
3551 "Returns the value of the variable in the rest frame of recoil system against the first Particle contained in the given ParticleList.\n"
3552 "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3553 "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3554 "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);
3555 REGISTER_METAVARIABLE(
"useDaughterRestFrame(variable, daughterIndex_1, [daughterIndex_2, ... daughterIndex_3])", useDaughterRestFrame,
3556 "Returns the value of the variable in the rest frame of the selected daughter particle.\n"
3557 "The daughter is identified via generalized daughter index, e.g. ``0:1`` identifies the second daughter (1) "
3558 "of the first daughter (0). If the daughter index is invalid, it returns NaN.\n"
3559 "If two or more indices are given, the rest frame of the sum of the daughters is used.",
3560 Manager::VariableDataType::c_double);
3561 REGISTER_METAVARIABLE(
"useDaughterRecoilRestFrame(variable, daughterIndex_1, [daughterIndex_2, ... daughterIndex_3])", useDaughterRecoilRestFrame,
3562 "Returns the value of the variable in the rest frame of the recoil of the selected daughter particle.\n"
3563 "The daughter is identified via generalized daughter index, e.g. ``0:1`` identifies the second daughter (1) "
3564 "of the first daughter (0). If the daughter index is invalid, it returns NaN.\n"
3565 "If two or more indices are given, the rest frame of the sum of the daughters is used.",
3566 Manager::VariableDataType::c_double);
3567 REGISTER_METAVARIABLE(
"useMCancestorBRestFrame(variable)", useMCancestorBRestFrame,
3568 "Returns the value of the variable in the rest frame of the ancestor B MC particle.\n"
3569 "If no B or no MC-matching is found, it returns NaN.", Manager::VariableDataType::c_double);
3570 REGISTER_METAVARIABLE(
"passesCut(cut)", passesCut,
3571 "Returns 1 if particle passes the cut otherwise 0.\n"
3572 "Useful if you want to write out if a particle would have passed a cut or not.", Manager::VariableDataType::c_bool);
3573 REGISTER_METAVARIABLE(
"passesEventCut(cut)", passesEventCut,
3574 "[Eventbased] Returns 1 if event passes the cut otherwise 0.\n"
3575 "Useful if you want to select events passing a cut without looping into particles, such as for skimming.\n", Manager::VariableDataType::c_bool);
3576 REGISTER_METAVARIABLE(
"countDaughters(cut)", countDaughters,
3577 "Returns number of direct daughters which satisfy the cut.\n"
3578 "Used by the skimming package (for what exactly?)", Manager::VariableDataType::c_int);
3579 REGISTER_METAVARIABLE(
"countFSPDaughters(cut)", countDescendants,
3580 "Returns number of final-state daughters which satisfy the cut.",
3581 Manager::VariableDataType::c_int);
3582 REGISTER_METAVARIABLE(
"countDescendants(cut)", countDescendants,
3583 "Returns number of descendants for all generations which satisfy the cut.",
3584 Manager::VariableDataType::c_int);
3585 REGISTER_METAVARIABLE(
"varFor(pdgCode, variable)", varFor,
3586 "Returns the value of the variable for the given particle if its abs(pdgCode) agrees with the given one.\n"
3587 "E.g. ``varFor(11, p)`` returns the momentum if the particle is an electron or a positron.", Manager::VariableDataType::c_double);
3588 REGISTER_METAVARIABLE(
"varForMCGen(variable)", varForMCGen,
3589 "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"
3590 "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"
3591 "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);
3592 REGISTER_METAVARIABLE(
"nParticlesInList(particleListName)", nParticlesInList,
3593 "[Eventbased] Returns number of particles in the given particle List.", Manager::VariableDataType::c_int);
3594 REGISTER_METAVARIABLE(
"isInList(particleListName)", isInList,
3595 "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);
3596 REGISTER_METAVARIABLE(
"isDaughterOfList(particleListNames)", isDaughterOfList,
3597 "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);
3598 REGISTER_METAVARIABLE(
"isDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isDescendantOfList, R
"DOC(
3599 Returns 1 if the given particle appears in the decay chain of the particles in the given ParticleLists.
3601 Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3603 * ``isDescendantOfList(<particle_list>,1)`` returns 1 if particle is a daughter of the list,
3604 * ``isDescendantOfList(<particle_list>,2)`` returns 1 if particle is a granddaughter of the list,
3605 * ``isDescendantOfList(<particle_list>,3)`` returns 1 if particle is a great-granddaughter of the list, etc.
3606 * Default value is ``-1`` that is inclusive for all generations.
3607 )DOC", Manager::VariableDataType::c_bool);
3608 REGISTER_METAVARIABLE("isMCDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isMCDescendantOfList, R
"DOC(
3609 Returns 1 if the given particle is linked to the same MC particle as any reconstructed daughter of the decay lists.
3611 Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3613 * ``isMCDescendantOfList(<particle_list>,1)`` returns 1 if particle is matched to the same particle as any daughter of the list,
3614 * ``isMCDescendantOfList(<particle_list>,2)`` returns 1 if particle is matched to the same particle as any granddaughter of the list,
3615 * ``isMCDescendantOfList(<particle_list>,3)`` returns 1 if particle is matched to the same particle as any great-granddaughter of the list, etc.
3616 * Default value is ``-1`` that is inclusive for all generations.
3618 It makes only sense for lists created with `fillParticleListFromMC` function with ``addDaughters=True`` argument.
3619 )DOC", Manager::VariableDataType::c_bool);
3621 REGISTER_METAVARIABLE("sourceObjectIsInList(particleListName)", sourceObjectIsInList, R
"DOC(
3622Returns 1 if the underlying mdst object (e.g. track, or cluster) was used to create a particle in ``particleListName``, 0 if not.
3625 This only makes sense for particles that are not composite. Returns -1 for composite particles.
3626)DOC", Manager::VariableDataType::c_int);
3628 REGISTER_METAVARIABLE("mcParticleIsInMCList(particleListName)", mcParticleIsInMCList, R
"DOC(
3629Returns 1 if the particle's matched MC particle is also matched to a particle in ``particleListName``
3630(or if either of the lists were filled from generator level `modularAnalysis.fillParticleListFromMC`.)
3632.. seealso:: :b2:var:`isMCDescendantOfList` to check daughters.
3633)DOC", Manager::VariableDataType::c_bool);
3635 REGISTER_METAVARIABLE("isGrandDaughterOfList(particleListNames)", isGrandDaughterOfList,
3636 "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);
3637 REGISTER_METAVARIABLE(
"originalParticle(variable)", originalParticle, R
"DOC(
3638 Returns value of variable for the original particle from which the given particle is copied.
3640 The copy of particle is created, for example, when the vertex fit updates the daughters and `modularAnalysis.copyParticles` is called.
3641 Returns NaN if the given particle is not copied and so there is no original particle.
3642 )DOC", Manager::VariableDataType::c_double);
3643 REGISTER_METAVARIABLE("daughter(i, variable)", daughter, R
"DOC(
3644 Returns value of variable for the i-th daughter. E.g.
3646 * ``daughter(0, p)`` returns the total momentum of the first daughter.
3647 * ``daughter(0, daughter(1, p)`` returns the total momentum of the second daughter of the first daughter.
3649 Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3650 )DOC", Manager::VariableDataType::c_double);
3651 REGISTER_METAVARIABLE("originalDaughter(i, variable)", originalDaughter, R
"DOC(
3652 Returns value of variable for the original particle from which the i-th daughter is copied.
3654 The copy of particle is created, for example, when the vertex fit updates the daughters and `modularAnalysis.copyParticles` is called.
3655 Returns NaN if the daughter is not copied and so there is no original daughter.
3657 Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3658 )DOC", Manager::VariableDataType::c_double);
3659 REGISTER_METAVARIABLE("mcDaughter(i, variable)", mcDaughter, R
"DOC(
3660 Returns the value of the requested variable for the i-th Monte Carlo daughter of the particle.
3662 Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3663 or if the i-th MC daughter does not exist.
3665 E.g. ``mcDaughter(0, PDG)`` will return the PDG code of the first MC daughter of the matched MC
3666 particle of the reconstructed particle the function is applied to.
3668 The meta variable can also be nested: ``mcDaughter(0, mcDaughter(1, PDG))``.
3669 )DOC", Manager::VariableDataType::c_double);
3670 REGISTER_METAVARIABLE("mcMother(variable)", mcMother, R
"DOC(
3671 Returns the value of the requested variable for the Monte Carlo mother of the particle.
3673 Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3674 or if the MC mother does not exist.
3676 E.g. ``mcMother(PDG)`` will return the PDG code of the MC mother of the matched MC
3677 particle of the reconstructed particle the function is applied to.
3679 The meta variable can also be nested: ``mcMother(mcMother(PDG))``.
3680 )DOC", Manager::VariableDataType::c_double);
3681 REGISTER_METAVARIABLE("genParticle(index, variable)", genParticle, R
"DOC(
3682[Eventbased] Returns the ``variable`` for the ith generator particle.
3683The arguments of the function must be the ``index`` of the particle in the MCParticle Array,
3684and ``variable``, the name of the function or variable for that generator particle.
3685If ``index`` goes beyond the length of the MCParticles array, NaN will be returned.
3687E.g. ``genParticle(0, p)`` returns the total momentum of the first MCParticle, which in a generic decay up to MC15 is
3688the Upsilon(4S) and for MC16 and beyond the initial electron.
3689)DOC", Manager::VariableDataType::c_double);
3690 REGISTER_METAVARIABLE("genUpsilon4S(variable)", genUpsilon4S, R
"DOC(
3691[Eventbased] Returns the ``variable`` evaluated for the generator-level :math:`\Upsilon(4S)`.
3692If no generator level :math:`\Upsilon(4S)` exists for the event, NaN will be returned.
3694E.g. ``genUpsilon4S(p)`` returns the total momentum of the :math:`\Upsilon(4S)` in a generic decay.
3695``genUpsilon4S(mcDaughter(1, p))`` returns the total momentum of the second daughter of the
3696generator-level :math:`\Upsilon(4S)` (i.e. the momentum of the second B meson in a generic decay).
3697)DOC", Manager::VariableDataType::c_double);
3698 REGISTER_METAVARIABLE("daughterProductOf(variable)", daughterProductOf,
3699 "Returns product of a variable over all daughters.\n"
3700 "E.g. ``daughterProductOf(extraInfo(SignalProbability))`` returns the product of the SignalProbabilitys of all daughters.", Manager::VariableDataType::c_double);
3701 REGISTER_METAVARIABLE(
"daughterSumOf(variable)", daughterSumOf,
3702 "Returns sum of a variable over all daughters.\n"
3703 "E.g. ``daughterSumOf(nDaughters)`` returns the number of grand-daughters.", Manager::VariableDataType::c_double);
3704 REGISTER_METAVARIABLE(
"daughterLowest(variable)", daughterLowest,
3705 "Returns the lowest value of the given variable among all daughters.\n"
3706 "E.g. ``useCMSFrame(daughterLowest(p))`` returns the lowest momentum in CMS frame.", Manager::VariableDataType::c_double);
3707 REGISTER_METAVARIABLE(
"daughterHighest(variable)", daughterHighest,
3708 "Returns the highest value of the given variable among all daughters.\n"
3709 "E.g. ``useCMSFrame(daughterHighest(p))`` returns the highest momentum in CMS frame.", Manager::VariableDataType::c_double);
3710 REGISTER_METAVARIABLE(
"daughterDiffOf(daughterIndex_i, daughterIndex_j, variable)", daughterDiffOf, R
"DOC(
3711 Returns the difference of a variable between the two given daughters.
3712 E.g. ``useRestFrame(daughterDiffOf(0, 1, p))`` returns the momentum difference between first and second daughter in the rest frame of the given particle.
3713 (That means that it returns :math:`p_j - p_i`)
3715 The daughters can be provided as generalized daughter indexes, which are simply colon-separated
3716 lists of daughter indexes, ordered starting from the root particle. For example, ``0:1``
3717 identifies the second daughter (1) of the first daughter (0) of the mother particle.
3719 )DOC", Manager::VariableDataType::c_double);
3720 REGISTER_METAVARIABLE("mcDaughterDiffOf(i, j, variable)", mcDaughterDiffOf,
3721 "MC matched version of the `daughterDiffOf` function.", Manager::VariableDataType::c_double);
3722 REGISTER_METAVARIABLE(
"grandDaughterDiffOf(i, j, variable)", grandDaughterDiffOf,
3723 "Returns the difference of a variable between the first daughters of the two given daughters.\n"
3724 "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"
3725 "(That means that it returns :math:`p_j - p_i`)", Manager::VariableDataType::c_double);
3726 MAKE_DEPRECATED(
"grandDaughterDiffOf",
false,
"light-2402-ocicat", R
"DOC(
3727 The difference between any combination of (grand-)daughters can be calculated with the more general variable :b2:var:`daughterDiffOf`
3728 by using generalized daughter indexes.)DOC");
3729 REGISTER_METAVARIABLE("daughterNormDiffOf(i, j, variable)", daughterNormDiffOf,
3730 "Returns the normalized difference of a variable between the two given daughters.\n"
3731 "E.g. ``daughterNormDiffOf(0, 1, p)`` returns the normalized momentum difference between first and second daughter in the lab frame.", Manager::VariableDataType::c_double);
3732 REGISTER_METAVARIABLE(
"daughterMotherDiffOf(i, variable)", daughterMotherDiffOf,
3733 "Returns the difference of a variable between the given daughter and the mother particle itself.\n"
3734 "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);
3735 REGISTER_METAVARIABLE(
"daughterMotherNormDiffOf(i, variable)", daughterMotherNormDiffOf,
3736 "Returns the normalized difference of a variable between the given daughter and the mother particle itself.\n"
3737 "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);
3738 REGISTER_METAVARIABLE(
"angleBetweenDaughterAndRecoil(daughterIndex_1, daughterIndex_2, ... )", angleBetweenDaughterAndRecoil, R
"DOC(
3739 Returns the angle between the momentum recoiling against the particle and the sum of the momenta of the given daughters.
3740 The unit of the angle is ``rad``.
3742 The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3743 daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3744 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3745 identifies the second daughter of the root particle.
3747 At least one generalized index has to be given to ``angleBetweenDaughterAndRecoil``.
3750 ``angleBetweenDaughterAndRecoil(0)`` will return the angle between pRecoil and the momentum of the first daughter.
3752 ``angleBetweenDaughterAndRecoil(0, 1)`` will return the angle between pRecoil and the sum of the momenta of the first and second daughter.
3754 ``angleBetweenDaughterAndRecoil(0:0, 3:0)`` will return the angle between pRecoil and the sum of the momenta of the: first daughter of the first daughter, and
3755 the first daughter of the fourth daughter.)DOC", Manager::VariableDataType::c_double);
3756 REGISTER_METAVARIABLE("angleBetweenDaughterAndMissingMomentum(daughterIndex_1, daughterIndex_2, ... )", angleBetweenDaughterAndMissingMomentum, R
"DOC(
3757 Returns the angle between the missing momentum in the event and the sum of the momenta of the given daughters.
3758 The unit of the angle is ``rad``. EventKinematics module has to be called to use this.
3760 The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3761 daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3762 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3763 identifies the second daughter of the root particle.
3765 At least one generalized index has to be given to ``angleBetweenDaughterAndMissingMomentum``.
3768 ``angleBetweenDaughterAndMissingMomentum(0)`` will return the angle between missMom and the momentum of the first daughter.
3770 ``angleBetweenDaughterAndMissingMomentum(0, 1)`` will return the angle between missMom and the sum of the momenta of the first and second daughter.
3772 ``angleBetweenDaughterAndMissingMomentum(0:0, 3:0)`` will return the angle between missMom and the sum of the momenta of the: first daughter of the first daughter, and
3773 the first daughter of the fourth daughter.)DOC", Manager::VariableDataType::c_double);
3774 REGISTER_METAVARIABLE("daughterAngle(daughterIndex_1, daughterIndex_2[, daughterIndex_3])", daughterAngle, R
"DOC(
3775 Returns the angle in between any pair of particles belonging to the same decay tree.
3776 The unit of the angle is ``rad``.
3778 The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3779 daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3780 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3781 identifies the second daughter of the root particle.
3783 Both two and three generalized indexes can be given to ``daughterAngle``. If two indices are given, the
3784 variable returns the angle between the momenta of the two given particles. If three indices are given, the
3785 variable returns the angle between the momentum of the third particle and a vector which is the sum of the
3786 first two daughter momenta.
3789 ``daughterAngle(0, 3)`` will return the angle between the first and fourth daughter.
3790 ``daughterAngle(0, 1, 3)`` will return the angle between the fourth daughter and the sum of the first and
3792 ``daughterAngle(0:0, 3:0)`` will return the angle between the first daughter of the first daughter, and
3793 the first daughter of the fourth daughter.
3795 )DOC", Manager::VariableDataType::c_double);
3796 REGISTER_METAVARIABLE("mcDaughterAngle(daughterIndex_1, daughterIndex_2, [daughterIndex_3])", mcDaughterAngle,
3797 "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);
3798 REGISTER_VARIABLE(
"grandDaughterDecayAngle(i, j)", grandDaughterDecayAngle,
3799 "Returns the decay angle of the granddaughter in the daughter particle's rest frame.\n"
3800 "It is calculated with respect to the reverted momentum vector of the particle.\n"
3801 "Two arguments representing the daughter and granddaughter indices have to be provided as arguments.\n\n",
"rad");
3802 REGISTER_VARIABLE(
"daughterClusterAngleInBetween(i, j)", daughterClusterAngleInBetween,
3803 "Returns the angle between clusters associated to the two daughters."
3804 "If two indices given: returns the angle between the momenta of the clusters associated to the two given daughters."
3805 "If three indices given: returns the angle between the momentum of the third particle's cluster and a vector "
3806 "which is the sum of the first two daughter's cluster momenta."
3807 "Returns nan if any of the daughters specified don't have an associated cluster."
3808 "The arguments in the argument vector must be integers corresponding to the ith and jth (and kth) daughters.\n\n",
"rad");
3809 REGISTER_METAVARIABLE(
"daughterInvM(i[, j, ...])", daughterInvM, R
"DOC(
3810 Returns the invariant mass adding the Lorentz vectors of the given daughters. The unit of the invariant mass is GeV/:math:`\text{c}^2`
3811 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.
3813 Daughters from different generations of the decay tree can be combined using generalized daughter indexes,
3814 which are simply colon-separated daughter indexes for each generation, starting from the root particle. For
3815 example, ``0:1:3`` identifies the fourth daughter (3) of the second daughter (1) of the first daughter(0) of
3816 the mother particle.
3818 Returns NaN if the given daughter-index is out of bound (>= number of daughters))DOC", Manager::VariableDataType::c_double);
3819 REGISTER_METAVARIABLE("extraInfo(name)", extraInfo,
3820 "Returns extra info stored under the given name.\n"
3821 "The extraInfo has to be set by a module first.\n"
3822 "E.g. ``extraInfo(SignalProbability)`` returns the SignalProbability calculated by the ``MVAExpert`` module.\n"
3823 "If nothing is set under the given name or if the particle is a nullptr, NaN is returned.\n"
3824 "In the latter case please use `eventExtraInfo` if you want to access an EventExtraInfo variable.", Manager::VariableDataType::c_double);
3825 REGISTER_METAVARIABLE(
"eventExtraInfo(name)", eventExtraInfo,
3826 "[Eventbased] Returns extra info stored under the given name in the event extra info.\n"
3827 "The extraInfo has to be set first by another module like MVAExpert in event mode.\n"
3828 "If nothing is set under this name, NaN is returned.", Manager::VariableDataType::c_double);
3829 REGISTER_METAVARIABLE(
"eventCached(variable)", eventCached,
3830 "[Eventbased] Returns value of event-based variable and caches this value in the EventExtraInfo.\n"
3831 "The result of second call to this variable in the same event will be provided from the cache.\n"
3832 "It is recommended to use this variable in order to declare custom aliases as event-based. This is "
3833 "necessary if using the eventwise mode of variablesToNtuple).", Manager::VariableDataType::c_double);
3834 REGISTER_METAVARIABLE(
"particleCached(variable)", particleCached,
3835 "Returns value of given variable and caches this value in the ParticleExtraInfo of the provided particle.\n"
3836 "The result of second call to this variable on the same particle will be provided from the cache.", Manager::VariableDataType::c_double);
3837 REGISTER_METAVARIABLE(
"modulo(variable, n)", modulo,
3838 "Returns rest of division of variable by n.", Manager::VariableDataType::c_int);
3839 REGISTER_METAVARIABLE(
"abs(variable)", abs,
3840 "Returns absolute value of the given variable.\n"
3841 "E.g. abs(mcPDG) returns the absolute value of the mcPDG, which is often useful for cuts.", Manager::VariableDataType::c_double);
3842 REGISTER_METAVARIABLE(
"max(var1,var2)", max,
"Returns max value of two variables.\n", Manager::VariableDataType::c_double);
3843 REGISTER_METAVARIABLE(
"min(var1,var2)", min,
"Returns min value of two variables.\n", Manager::VariableDataType::c_double);
3844 REGISTER_METAVARIABLE(
"sin(variable)", sin,
"Returns sine value of the given variable.", Manager::VariableDataType::c_double);
3845 REGISTER_METAVARIABLE(
"asin(variable)", asin,
"Returns arcsine of the given variable. The unit of the asin() is ``rad``", Manager::VariableDataType::c_double);
3846 REGISTER_METAVARIABLE(
"cos(variable)", cos,
"Returns cosine value of the given variable.", Manager::VariableDataType::c_double);
3847 REGISTER_METAVARIABLE(
"acos(variable)", acos,
"Returns arccosine value of the given variable. The unit of the acos() is ``rad``", Manager::VariableDataType::c_double);
3848 REGISTER_METAVARIABLE(
"tan(variable)", tan,
"Returns tangent value of the given variable.", Manager::VariableDataType::c_double);
3849 REGISTER_METAVARIABLE(
"atan(variable)", atan,
"Returns arctangent value of the given variable. The unit of the atan() is ``rad``", Manager::VariableDataType::c_double);
3850 REGISTER_METAVARIABLE(
"exp(variable)", exp,
"Returns exponential evaluated for the given variable.", Manager::VariableDataType::c_double);
3851 REGISTER_METAVARIABLE(
"log(variable)", log,
"Returns natural logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3852 REGISTER_METAVARIABLE(
"log10(variable)", log10,
"Returns base-10 logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3853 REGISTER_METAVARIABLE(
"int(variable, nan_replacement)", convertToInt, R
"DOC(
3854 Casts the output of the variable to an integer value.
3857 Overflow and underflow are clipped at maximum and minimum values, respectively. NaN values are replaced with the value of the 2nd argument.
3859 )DOC", Manager::VariableDataType::c_int);
3860 REGISTER_METAVARIABLE("isNAN(variable)", isNAN,
3861 "Returns true if variable value evaluates to nan (determined via std::isnan(double)).\n"
3862 "Useful for debugging.", Manager::VariableDataType::c_bool);
3863 REGISTER_METAVARIABLE(
"ifNANgiveX(variable, x)", ifNANgiveX,
3864 "Returns x (has to be a number) if variable value is nan (determined via std::isnan(double)).\n"
3865 "Useful for technical purposes while training MVAs.", Manager::VariableDataType::c_double);
3866 REGISTER_METAVARIABLE(
"isInfinity(variable)", isInfinity,
3867 "Returns true if variable value evaluates to infinity (determined via std::isinf(double)).\n"
3868 "Useful for debugging.", Manager::VariableDataType::c_bool);
3869 REGISTER_METAVARIABLE(
"unmask(variable, flag1, flag2, ...)", unmask,
3870 "unmask(variable, flag1, flag2, ...) or unmask(variable, mask) sets certain bits in the variable to zero.\n"
3871 "For example, if you want to set the second, fourth and fifth bits to zero, you could call \n"
3872 "``unmask(variable, 2, 8, 16)`` or ``unmask(variable, 26)``.\n"
3873 "", Manager::VariableDataType::c_double);
3874 REGISTER_METAVARIABLE(
"conditionalVariableSelector(cut, variableIfTrue, variableIfFalse)", conditionalVariableSelector,
3875 "Returns one of the two supplied variables, depending on whether the particle passes the supplied cut.\n"
3876 "The first variable is returned if the particle passes the cut, and the second variable is returned otherwise.", Manager::VariableDataType::c_double);
3877 REGISTER_METAVARIABLE(
"pValueCombination(p1, p2, ...)", pValueCombination,
3878 "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"
3879 "If any of the p-values is invalid, i.e. smaller than zero, -1 is returned.", Manager::VariableDataType::c_double);
3880 REGISTER_METAVARIABLE(
"pValueCombinationOfDaughters(variable)", pValueCombinationOfDaughters,
3881 "Returns the combined p-value of the daughter 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"
3882 "If any of the p-values is invalid, i.e. smaller than zero, -1 is returned.", Manager::VariableDataType::c_double);
3883 REGISTER_METAVARIABLE(
"veto(particleList, cut, pdgCode = 11)", veto,
3884 "Combines current particle with particles from the given particle list and returns 1 if the combination passes the provided cut. \n"
3885 "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"
3886 "around the neutral Pion mass (e.g. ``0.130 < M < 0.140``). \n"
3887 "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);
3888 REGISTER_METAVARIABLE(
"matchedMC(variable)", matchedMC,
3889 "Returns variable output for the matched MCParticle by constructing a temporary Particle from it.\n"
3890 "This may not work too well if your variable requires accessing daughters of the particle.\n"
3891 "E.g. ``matchedMC(p)`` returns the total momentum of the related MCParticle.\n"
3892 "Returns NaN if no matched MCParticle exists.", Manager::VariableDataType::c_double);
3893 REGISTER_METAVARIABLE(
"clusterBestMatchedMCParticle(variable)", clusterBestMatchedMCParticle,
3894 "Returns variable output for the MCParticle that is best-matched with the ECLCluster of the given Particle.\n"
3895 "E.g. To get the energy of the MCParticle that matches best with an ECLCluster, one could use ``clusterBestMatchedMCParticle(E)``\n"
3896 "When the variable is called for ``gamma`` and if the ``gamma`` is matched with MCParticle, it works same as `matchedMC`.\n"
3897 "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"
3898 "Returns NaN if the particle is not matched to an ECLCluster, or if the ECLCluster has no matching MCParticles", Manager::VariableDataType::c_double);
3899 REGISTER_METAVARIABLE(
"varForBestMatchedMCKlong(variable)", clusterBestMatchedMCKlong,
3900 "Returns variable output for the Klong MCParticle which has the best match with the ECLCluster of the given Particle.\n"
3901 "Returns NaN if the particle is not matched to an ECLCluster, or if the ECLCluster has no matching Klong MCParticle", Manager::VariableDataType::c_double);
3903 REGISTER_METAVARIABLE(
"countInList(particleList, cut='')", countInList,
"[Eventbased] "
3904 "Returns number of particle which pass given in cut in the specified particle list.\n"
3905 "Useful for creating statistics about the number of particles in a list.\n"
3906 "E.g. ``countInList(e+, isSignal == 1)`` returns the number of correctly reconstructed electrons in the event.\n"
3907 "The variable is event-based and does not need a valid particle pointer as input.", Manager::VariableDataType::c_int);
3908 REGISTER_METAVARIABLE(
"getVariableByRank(particleList, rankedVariableName, variableName, rank)", getVariableByRank, R
"DOC(
3909 [Eventbased] Returns the value of ``variableName`` for the candidate in the ``particleList`` with the requested ``rank``.
3912 The `BestCandidateSelection` module available via `rankByHighest` / `rankByLowest` has to be used before.
3915 The first candidate matching the given rank is used.
3916 Thus, it is not recommended to use this variable in conjunction with ``allowMultiRank`` in the `BestCandidateSelection` module.
3918 The suffix ``_rank`` is automatically added to the argument ``rankedVariableName``,
3919 which either has to be the name of the variable used to order the candidates or the selected outputVariable name without the ending ``_rank``.
3920 This means that your selected name for the rank variable has to end with ``_rank``.
3922 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>`_
3923 )DOC", Manager::VariableDataType::c_double);
3924 REGISTER_VARIABLE("matchedMCHasPDG(PDGCode)", matchedMCHasPDG,
3925 "Returns if the absolute value of the PDGCode of the MCParticle related to the Particle matches a given PDGCode."
3926 "Returns 0/NAN/1 if PDGCode does not match/is not available/ matches");
3927 REGISTER_METAVARIABLE(
"numberOfNonOverlappingParticles(pList1, pList2, ...)", numberOfNonOverlappingParticles,
3928 "Returns the number of non-overlapping particles in the given particle lists"
3929 "Useful to check if there is additional physics going on in the detector if one reconstructed the Y4S", Manager::VariableDataType::c_int);
3930 REGISTER_METAVARIABLE(
"totalEnergyOfParticlesInList(particleListName)", totalEnergyOfParticlesInList,
3931 "[Eventbased] Returns the total energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3932 REGISTER_METAVARIABLE(
"totalPxOfParticlesInList(particleListName)", totalPxOfParticlesInList,
3933 "[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);
3934 REGISTER_METAVARIABLE(
"totalPyOfParticlesInList(particleListName)", totalPyOfParticlesInList,
3935 "[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);
3936 REGISTER_METAVARIABLE(
"totalPzOfParticlesInList(particleListName)", totalPzOfParticlesInList,
3937 "[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);
3938 REGISTER_METAVARIABLE(
"invMassInLists(pList1, pList2, ...)", invMassInLists,
3939 "[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);
3940 REGISTER_METAVARIABLE(
"totalECLEnergyOfParticlesInList(particleListName)", totalECLEnergyOfParticlesInList,
3941 "[Eventbased] Returns the total ECL energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3942 REGISTER_METAVARIABLE(
"maxPtInList(particleListName)", maxPtInList,
3943 "[Eventbased] Returns maximum transverse momentum Pt in the given particle List. The unit of the transverse momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3944 REGISTER_METAVARIABLE(
"eclClusterSpecialTrackMatched(cut)", eclClusterTrackMatchedWithCondition,
3945 "Returns if at least one Track that satisfies the given condition is related to the ECLCluster of the Particle.", Manager::VariableDataType::c_double);
3946 REGISTER_METAVARIABLE(
"averageValueInList(particleListName, variable)", averageValueInList,
3947 "[Eventbased] Returns the arithmetic mean of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3948 REGISTER_METAVARIABLE(
"medianValueInList(particleListName, variable)", medianValueInList,
3949 "[Eventbased] Returns the median value of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3950 REGISTER_METAVARIABLE(
"sumValueInList(particleListName, variable)", sumValueInList,
3951 "[Eventbased] Returns the sum of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3952 REGISTER_METAVARIABLE(
"productValueInList(particleListName, variable)", productValueInList,
3953 "[Eventbased] Returns the product of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3954 REGISTER_METAVARIABLE(
"angleToClosestInList(particleListName)", angleToClosestInList,
3955 "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);
3956 REGISTER_METAVARIABLE(
"closestInList(particleListName, variable)", closestInList,
3957 "Returns `variable` for the closest particle (smallest opening angle) in the list provided.", Manager::VariableDataType::c_double);
3958 REGISTER_METAVARIABLE(
"angleToMostB2BInList(particleListName)", angleToMostB2BInList,
3959 "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);
3960 REGISTER_METAVARIABLE(
"deltaPhiToMostB2BPhiInList(particleListName)", deltaPhiToMostB2BPhiInList,
3961 "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);
3962 REGISTER_METAVARIABLE(
"mostB2BInList(particleListName, variable)", mostB2BInList,
3963 "Returns `variable` for the most back-to-back particle (closest opening angle to 180) in the list provided.", Manager::VariableDataType::c_double);
3964 REGISTER_METAVARIABLE(
"maxOpeningAngleInList(particleListName)", maxOpeningAngleInList,
3965 "[Eventbased] Returns maximum opening angle in the given particle List. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
3966 REGISTER_METAVARIABLE(
"daughterCombination(variable, daughterIndex_1, daughterIndex_2 ... daughterIndex_n)", daughterCombination,R
"DOC(
3967Returns a ``variable`` function only of the 4-momentum calculated on an arbitrary set of (grand)daughters.
3970 ``variable`` can only be a function of the daughters' 4-momenta.
3972Daughters from different generations of the decay tree can be combined using generalized daughter indexes, which are simply colon-separated
3973the list of daughter indexes, starting from the root particle: for example, ``0:1:3`` identifies the fourth
3974daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle.
3977 ``daughterCombination(M, 0, 3, 4)`` will return the invariant mass of the system made of the first, fourth and fifth daughter of particle.
3978 ``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.
3980)DOC", Manager::VariableDataType::c_double);
3981 REGISTER_METAVARIABLE("useAlternativeDaughterHypothesis(variable, daughterIndex_1:newMassHyp_1, ..., daughterIndex_n:newMassHyp_n)", useAlternativeDaughterHypothesis,R
"DOC(
3982Returns a ``variable`` calculated using new mass hypotheses for (some of) the particle's daughters.
3985 ``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.
3986 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.
3987 Also, the track fit is not performed again: the variable only re-calculates the 4-vectors using different mass assumptions.
3988 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).
3991 Generalized daughter indexes are not supported (yet!): this variable can be used only on first-generation daughters.
3994 ``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.
3995 ``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.
3997)DOC", Manager::VariableDataType::c_double);
3998 REGISTER_METAVARIABLE("varForFirstMCAncestorOfType(type, variable)",varForFirstMCAncestorOfType,R
"DOC(Returns requested variable of the first ancestor of the given type.
3999Ancestor type can be set up by PDG code or by particle name (check evt.pdl for valid particle names))DOC", Manager::VariableDataType::c_double);
4001 REGISTER_METAVARIABLE("nTrackFitResults(particleType)", nTrackFitResults,
4002 "[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).",
4003 Manager::VariableDataType::c_int);
4005 REGISTER_METAVARIABLE(
"convertToDaughterIndex(variable)", convertToDaughterIndex, R
"DOC(Converts the variable of the given particle into integer and returns it if it is a valid daughter index, else returns -1.)DOC", Manager::VariableDataType::c_int);
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)
@ 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.
EParticleSourceObject
particle source enumerators
@ 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.
#define MAKE_DEPRECATED(name, make_fatal, version, description)
Registers a variable as deprecated.
T convertString(const std::string &str)
Converts a string to type T (one of float, double, long double, int, long int, unsigned long int).
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.