10 #include <analysis/variables/MetaVariables.h>
12 #include <analysis/VariableManager/Utility.h>
13 #include <analysis/dataobjects/Particle.h>
14 #include <analysis/dataobjects/ParticleList.h>
15 #include <analysis/dataobjects/RestOfEvent.h>
16 #include <analysis/utility/PCmsLabTransform.h>
17 #include <analysis/utility/ReferenceFrame.h>
18 #include <analysis/utility/EvtPDLUtil.h>
19 #include <analysis/utility/ParticleCopy.h>
20 #include <analysis/utility/ValueIndexPairSorting.h>
21 #include <analysis/ClusterUtility/ClusterUtils.h>
22 #include <analysis/variables/VariableFormulaConstructor.h>
24 #include <framework/logging/Logger.h>
25 #include <framework/datastore/StoreArray.h>
26 #include <framework/datastore/StoreObjPtr.h>
27 #include <framework/dataobjects/EventExtraInfo.h>
28 #include <framework/utilities/Conversion.h>
29 #include <framework/utilities/MakeROOTCompatible.h>
30 #include <framework/gearbox/Const.h>
32 #include <mdst/dataobjects/Track.h>
33 #include <mdst/dataobjects/MCParticle.h>
34 #include <mdst/dataobjects/ECLCluster.h>
35 #include <mdst/dataobjects/TrackFitResult.h>
37 #include <boost/algorithm/string.hpp>
45 #include <TDatabasePDG.h>
46 #include <Math/Vector4D.h>
57 if (arguments.size() == 1) {
59 auto func = [var](
const Particle * particle) ->
double {
60 UseReferenceFrame<RestFrame> frame(particle);
61 double result = std::get<double>(var->function(particle));
66 B2FATAL(
"Wrong number of arguments for meta function useRestFrame");
72 if (arguments.size() == 1) {
74 auto func = [var](
const Particle * particle) ->
double {
75 UseReferenceFrame<CMSFrame> frame;
76 double result = std::get<double>(var->function(particle));
81 B2FATAL(
"Wrong number of arguments for meta function useCMSFrame");
87 if (arguments.size() == 1) {
89 auto func = [var](
const Particle * particle) ->
double {
90 UseReferenceFrame<LabFrame> frame;
91 double result = std::get<double>(var->function(particle));
96 B2FATAL(
"Wrong number of arguments for meta function useLabFrame");
102 if (arguments.size() == 2) {
105 int daughterIndexTagB = 0;
107 daughterIndexTagB = Belle2::convertString<int>(arguments[1]);
108 }
catch (std::invalid_argument&) {
109 B2FATAL(
"Second argument of useTagSideRecoilRestFrame meta function must be integer!");
112 auto func = [var, daughterIndexTagB](
const Particle * particle) ->
double {
113 if (particle->getPDGCode() != 300553)
115 B2ERROR(
"Variable should only be used on a Upsilon(4S) Particle List!");
120 ROOT::Math::PxPyPzEVector pSigB = T.getBeamFourMomentum() - particle->getDaughter(daughterIndexTagB)->get4Vector();
121 Particle tmp(pSigB, -particle->getDaughter(daughterIndexTagB)->getPDGCode());
123 UseReferenceFrame<RestFrame> frame(&tmp);
124 double result = std::get<double>(var->function(particle));
130 B2FATAL(
"Wrong number of arguments for meta function useTagSideRecoilRestFrame");
136 if (arguments.size() == 2) {
138 std::string listName = arguments[1];
139 auto func = [var, listName](
const Particle * particle) ->
double {
140 StoreObjPtr<ParticleList> list(listName);
141 unsigned listSize = list->getListSize();
145 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."
146 <<
LogVar(
"ParticleList", listName)
147 <<
LogVar(
"Number of candidates in the list", listSize));
148 const Particle* p = list->getParticle(0);
149 UseReferenceFrame<RestFrame> frame(p);
150 double result = std::get<double>(var->function(particle));
155 B2FATAL(
"Wrong number of arguments for meta function useParticleRestFrame.");
161 if (arguments.size() == 2) {
163 std::string listName = arguments[1];
164 auto func = [var, listName](
const Particle * particle) ->
double {
165 StoreObjPtr<ParticleList> list(listName);
166 unsigned listSize = list->getListSize();
170 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."
171 <<
LogVar(
"ParticleList", listName)
172 <<
LogVar(
"Number of candidates in the list", listSize));
173 const Particle* p = list->getParticle(0);
175 ROOT::Math::PxPyPzEVector recoil = T.getBeamFourMomentum() - p->get4Vector();
177 Particle pRecoil(recoil, 0);
178 pRecoil.setVertex(particle->getVertex());
179 UseReferenceFrame<RestFrame> frame(&pRecoil);
180 double result = std::get<double>(var->function(particle));
185 B2FATAL(
"Wrong number of arguments for meta function useParticleRestFrame.");
191 if (arguments.size() >= 2) {
193 auto func = [var, arguments](
const Particle * particle) ->
double {
196 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
198 for (
unsigned int i = 1; i < arguments.size(); i++)
200 auto generalizedIndex = arguments[i];
201 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
203 pSum += dauPart->get4Vector();
207 Particle tmp(pSum, 0);
208 UseReferenceFrame<RestFrame> frame(&tmp);
209 double result = std::get<double>(var->function(particle));
214 B2FATAL(
"Wrong number of arguments for meta function useDaughterRestFrame.");
220 if (arguments.size() == 1) {
221 auto extraInfoName = arguments[0];
222 auto func = [extraInfoName](
const Particle * particle) ->
double {
223 if (particle ==
nullptr)
225 B2WARNING(
"Returns NaN because the particle is nullptr! If you want EventExtraInfo variables, please use eventExtraInfo() instead");
228 if (particle->hasExtraInfo(extraInfoName))
230 return particle->getExtraInfo(extraInfoName);
238 B2FATAL(
"Wrong number of arguments for meta function extraInfo");
244 if (arguments.size() == 1) {
245 auto extraInfoName = arguments[0];
246 auto func = [extraInfoName](
const Particle*) ->
double {
247 StoreObjPtr<EventExtraInfo> eventExtraInfo;
248 if (not eventExtraInfo.isValid())
250 if (eventExtraInfo->hasExtraInfo(extraInfoName))
252 return eventExtraInfo->getExtraInfo(extraInfoName);
260 B2FATAL(
"Wrong number of arguments for meta function extraInfo");
266 if (arguments.size() == 1) {
269 auto func = [var, key](
const Particle*) ->
double {
271 StoreObjPtr<EventExtraInfo> eventExtraInfo;
272 if (not eventExtraInfo.isValid())
273 eventExtraInfo.create();
274 if (eventExtraInfo->hasExtraInfo(key))
276 return eventExtraInfo->getExtraInfo(key);
280 auto var_result = var->function(
nullptr);
281 if (std::holds_alternative<double>(var_result)) {
282 value = std::get<double>(var_result);
283 }
else if (std::holds_alternative<int>(var_result)) {
284 return std::get<int>(var_result);
285 }
else if (std::holds_alternative<bool>(var_result)) {
286 return std::get<bool>(var_result);
288 eventExtraInfo->addExtraInfo(key, value);
294 B2FATAL(
"Wrong number of arguments for meta function eventCached");
300 if (arguments.size() == 1) {
303 auto func = [var, key](
const Particle * particle) ->
double {
305 if (particle->hasExtraInfo(key))
307 return particle->getExtraInfo(key);
310 double value = std::get<double>(var->function(particle));
318 const_cast<Particle*
>(particle)->addExtraInfo(key, value);
324 B2FATAL(
"Wrong number of arguments for meta function particleCached");
333 if (arguments.size() != 1) B2FATAL(
"Wrong number of arguments for meta function formula");
334 FormulaParser<VariableFormulaConstructor> parser;
336 return parser.parse(arguments[0]);
337 }
catch (std::runtime_error& e) {
344 if (arguments.size() <= 1) {
346 std::string cutString;
347 if (arguments.size() == 1)
348 cutString = arguments[0];
350 auto func = [cut](
const Particle*) ->
int {
352 int number_of_tracks = 0;
353 StoreArray<Track> tracks;
354 for (
const auto& track : tracks)
356 const TrackFitResult* trackFit = track.getTrackFitResultWithClosestMass(
Const::pion);
357 if (trackFit->getChargeSign() == 0) {
361 if (cut->check(&particle))
366 return number_of_tracks;
371 B2FATAL(
"Wrong number of arguments for meta function nCleanedTracks");
377 if (arguments.size() <= 1) {
379 std::string cutString;
380 if (arguments.size() == 1)
381 cutString = arguments[0];
383 auto func = [cut](
const Particle*) ->
int {
385 int number_of_clusters = 0;
386 StoreArray<ECLCluster> clusters;
387 for (
const auto& cluster : clusters)
393 Particle particle(&cluster);
394 if (cut->check(&particle))
395 number_of_clusters++;
398 return number_of_clusters;
403 B2FATAL(
"Wrong number of arguments for meta function nCleanedECLClusters");
409 if (arguments.size() == 1) {
410 std::string cutString = arguments[0];
412 auto func = [cut](
const Particle * particle) ->
bool {
413 if (cut->check(particle))
420 B2FATAL(
"Wrong number of arguments for meta function passesCut");
426 if (arguments.size() == 1) {
427 std::string cutString = arguments[0];
429 auto func = [cut](
const Particle*) ->
bool {
430 if (cut->check(
nullptr))
437 B2FATAL(
"Wrong number of arguments for meta function passesEventCut");
443 if (arguments.size() == 2) {
446 pdgCode = Belle2::convertString<int>(arguments[0]);
447 }
catch (std::invalid_argument&) {
448 B2FATAL(
"The first argument of varFor meta function must be a positive integer!");
451 auto func = [pdgCode, var](
const Particle * particle) ->
double {
452 if (std::abs(particle->getPDGCode()) == std::abs(pdgCode))
454 auto var_result = var->function(particle);
455 if (std::holds_alternative<double>(var_result)) {
456 return std::get<double>(var_result);
457 }
else if (std::holds_alternative<int>(var_result)) {
458 return std::get<int>(var_result);
459 }
else if (std::holds_alternative<bool>(var_result)) {
460 return std::get<bool>(var_result);
466 B2FATAL(
"Wrong number of arguments for meta function varFor");
472 if (arguments.size() == 1) {
474 auto func = [var](
const Particle * particle) ->
double {
475 if (particle->getMCParticle())
480 auto var_result = var->function(particle);
481 if (std::holds_alternative<double>(var_result)) {
482 return std::get<double>(var_result);
483 }
else if (std::holds_alternative<int>(var_result)) {
484 return std::get<int>(var_result);
485 }
else if (std::holds_alternative<bool>(var_result)) {
486 return std::get<bool>(var_result);
493 B2FATAL(
"Wrong number of arguments for meta function varForMCGen");
499 if (arguments.size() == 1) {
500 std::string listName = arguments[0];
501 auto func = [listName](
const Particle * particle) ->
int {
504 StoreObjPtr<ParticleList> listOfParticles(listName);
506 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to nParticlesInList");
508 return listOfParticles->getListSize();
513 B2FATAL(
"Wrong number of arguments for meta function nParticlesInList");
520 if (arguments.size() != 1) {
521 B2FATAL(
"Wrong number of arguments for isInList");
523 auto listName = arguments[0];
525 auto func = [listName](
const Particle * particle) ->
bool {
528 StoreObjPtr<ParticleList> list(listName);
529 if (!(list.isValid()))
531 B2FATAL(
"Invalid Listname " << listName <<
" given to isInList");
535 return list->contains(particle);
544 if (arguments.size() != 1) {
545 B2FATAL(
"Wrong number of arguments for sourceObjectIsInList");
547 auto listName = arguments[0];
549 auto func = [listName](
const Particle * particle) ->
int {
552 StoreObjPtr<ParticleList> list(listName);
553 if (!(list.isValid()))
555 B2FATAL(
"Invalid Listname " << listName <<
" given to sourceObjectIsInList");
561 if (particlesource == Particle::EParticleSourceObject::c_Composite
562 or particlesource == Particle::EParticleSourceObject::c_Undefined)
568 for (
unsigned i = 0; i < list->getListSize(); ++i)
570 Particle* iparticle = list->getParticle(i);
571 if (particlesource == iparticle->getParticleSource())
572 if (particle->getMdstArrayIndex() == iparticle->getMdstArrayIndex())
584 if (arguments.size() != 1) {
585 B2FATAL(
"Wrong number of arguments for mcParticleIsInMCList");
587 auto listName = arguments[0];
589 auto func = [listName](
const Particle * particle) ->
bool {
592 StoreObjPtr<ParticleList> list(listName);
593 if (!(list.isValid()))
594 B2FATAL(
"Invalid Listname " << listName <<
" given to mcParticleIsInMCList");
597 const MCParticle* mcp = particle->getMCParticle();
598 if (mcp ==
nullptr)
return false;
601 for (
unsigned i = 0; i < list->getListSize(); ++i)
603 const MCParticle* imcp = list->getParticle(i)->
getMCParticle();
604 if ((imcp !=
nullptr) and (mcp->getArrayIndex() == imcp->getArrayIndex()))
614 B2WARNING(
"isDaughterOfList is outdated and replaced by isDescendantOfList.");
615 std::vector<std::string> new_arguments = arguments;
616 new_arguments.push_back(std::string(
"1"));
617 return isDescendantOfList(new_arguments);
622 B2WARNING(
"isGrandDaughterOfList is outdated and replaced by isDescendantOfList.");
623 std::vector<std::string> new_arguments = arguments;
624 new_arguments.push_back(std::string(
"2"));
625 return isDescendantOfList(new_arguments);
630 if (arguments.size() > 0) {
631 auto listNames = arguments;
632 auto func = [listNames](
const Particle * particle) ->
bool {
634 int generation_flag = -1;
637 generation_flag = Belle2::convertString<int>(listNames.back());
638 }
catch (std::exception& e) {}
640 for (
auto& iListName : listNames)
643 Belle2::convertString<int>(iListName);
645 }
catch (std::exception& e) {}
648 auto list_comparison = [](
auto&&
self,
const Particle * m,
const Particle * p,
int flag)->
bool {
650 for (
unsigned i = 0; i < m->getNDaughters(); ++i)
652 const Particle* daughter = m->getDaughter(i);
653 if ((flag == 1.) or (flag < 0)) {
654 if (p->isCopyOf(daughter)) {
660 if (daughter->getNDaughters() > 0) {
661 result =
self(
self, daughter, p, flag - 1);
671 StoreObjPtr<ParticleList> listOfParticles(iListName);
673 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << iListName <<
" given to isDescendantOfList");
675 for (
unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
676 Particle* iParticle = listOfParticles->getParticle(i);
677 output = list_comparison(list_comparison, iParticle, particle, generation_flag);
687 B2FATAL(
"Wrong number of arguments for meta function isDescendantOfList");
693 if (arguments.size() > 0) {
694 auto listNames = arguments;
695 auto func = [listNames](
const Particle * particle) ->
bool {
697 int generation_flag = -1;
700 generation_flag = Belle2::convertString<int>(listNames.back());
701 }
catch (std::exception& e) {}
703 if (particle->getMCParticle() ==
nullptr)
708 for (
auto& iListName : listNames)
711 std::stod(iListName);
713 }
catch (std::exception& e) {}
715 auto list_comparison = [](
auto&&
self,
const Particle * m,
const Particle * p,
int flag)->
bool {
717 for (
unsigned i = 0; i < m->getNDaughters(); ++i)
719 const Particle* daughter = m->getDaughter(i);
720 if ((flag == 1.) or (flag < 0)) {
721 if (daughter->getMCParticle() !=
nullptr) {
722 if (p->getMCParticle()->getArrayIndex() == daughter->getMCParticle()->getArrayIndex()) {
728 if (daughter->getNDaughters() > 0) {
729 result =
self(
self, daughter, p, flag - 1);
739 StoreObjPtr<ParticleList> listOfParticles(iListName);
741 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << iListName <<
" given to isMCDescendantOfList");
743 for (
unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
744 Particle* iParticle = listOfParticles->getParticle(i);
745 output = list_comparison(list_comparison, iParticle, particle, generation_flag);
755 B2FATAL(
"Wrong number of arguments for meta function isMCDescendantOfList");
761 if (arguments.size() == 1) {
763 auto func = [var](
const Particle * particle) ->
double {
764 double product = 1.0;
765 if (particle->getNDaughters() == 0)
769 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
771 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
772 product *= std::get<double>(var->function(particle->getDaughter(j)));
774 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
776 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
777 product *= std::get<int>(var->function(particle->getDaughter(j)));
784 B2FATAL(
"Wrong number of arguments for meta function daughterProductOf");
790 if (arguments.size() == 1) {
792 auto func = [var](
const Particle * particle) ->
double {
794 if (particle->getNDaughters() == 0)
798 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
800 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
801 sum += std::get<double>(var->function(particle->getDaughter(j)));
803 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
805 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
806 sum += std::get<int>(var->function(particle->getDaughter(j)));
813 B2FATAL(
"Wrong number of arguments for meta function daughterSumOf");
819 if (arguments.size() == 1) {
821 auto func = [var](
const Particle * particle) ->
double {
823 if (particle->getNDaughters() == 0)
827 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
829 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
830 double iValue = std::get<double>(var->function(particle->getDaughter(j)));
831 if (std::isnan(iValue))
continue;
832 if (std::isnan(min)) min = iValue;
833 if (iValue < min) min = iValue;
835 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
837 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
838 int iValue = std::get<int>(var->function(particle->getDaughter(j)));
839 if (std::isnan(min)) min = iValue;
840 if (iValue < min) min = iValue;
847 B2FATAL(
"Wrong number of arguments for meta function daughterLowest");
853 if (arguments.size() == 1) {
855 auto func = [var](
const Particle * particle) ->
double {
857 if (particle->getNDaughters() == 0)
861 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
863 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
864 double iValue = std::get<double>(var->function(particle->getDaughter(j)));
865 if (std::isnan(iValue))
continue;
866 if (std::isnan(max)) max = iValue;
867 if (iValue > max) max = iValue;
869 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
871 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
872 int iValue = std::get<int>(var->function(particle->getDaughter(j)));
873 if (std::isnan(max)) max = iValue;
874 if (iValue > max) max = iValue;
881 B2FATAL(
"Wrong number of arguments for meta function daughterHighest");
887 if (arguments.size() == 3) {
888 int iDaughterNumber = 0;
889 int jDaughterNumber = 0;
891 iDaughterNumber = Belle2::convertString<int>(arguments[0]);
892 jDaughterNumber = Belle2::convertString<int>(arguments[1]);
893 }
catch (std::invalid_argument&) {
894 B2FATAL(
"First two arguments of daughterDiffOf meta function must be integers!");
896 auto variablename = arguments[2];
897 auto func = [variablename, iDaughterNumber, jDaughterNumber](
const Particle * particle) ->
double {
898 if (particle ==
nullptr)
900 if (iDaughterNumber >=
int(particle->getNDaughters()) || jDaughterNumber >= int(particle->getNDaughters()))
905 auto result_j = var->function(particle->getDaughter(jDaughterNumber));
906 auto result_i = var->function(particle->getDaughter(iDaughterNumber));
908 if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i)) {
909 diff = std::get<double>(result_j) - std::get<double>(result_i);
910 }
else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i)) {
911 diff = std::get<int>(result_j) - std::get<int>(result_i);
913 throw std::runtime_error(
"Bad variant access");
915 if (variablename ==
"phi" or variablename ==
"clusterPhi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)"))
916 or std::regex_match(variablename, std::regex(
"use.*Frame\\(clusterPhi\\)"))) {
917 if (fabs(diff) > M_PI) {
919 diff = diff - 2 * M_PI;
921 diff = 2 * M_PI + diff;
930 B2FATAL(
"Wrong number of arguments for meta function daughterDiffOf");
936 if (arguments.size() == 3) {
937 int iDaughterNumber = 0;
938 int jDaughterNumber = 0;
940 iDaughterNumber = Belle2::convertString<int>(arguments[0]);
941 jDaughterNumber = Belle2::convertString<int>(arguments[1]);
942 }
catch (std::invalid_argument&) {
943 B2FATAL(
"First two arguments of mcDaughterDiffOf meta function must be integers!");
945 auto variablename = arguments[2];
946 auto func = [variablename, iDaughterNumber, jDaughterNumber](
const Particle * particle) ->
double {
947 if (particle ==
nullptr)
949 if (iDaughterNumber >=
int(particle->getNDaughters()) || jDaughterNumber >= int(particle->getNDaughters()))
951 if (particle->getDaughter(jDaughterNumber)->getMCParticle() ==
nullptr || particle->getDaughter(iDaughterNumber)->getMCParticle() ==
nullptr)
955 const MCParticle* iMcDaughter = particle->getDaughter(iDaughterNumber)->getMCParticle();
956 const MCParticle* jMcDaughter = particle->getDaughter(jDaughterNumber)->getMCParticle();
957 Particle iTmpPart(iMcDaughter);
958 Particle jTmpPart(jMcDaughter);
960 auto result_j = var->function(&jTmpPart);
961 auto result_i = var->function(&iTmpPart);
963 if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i)) {
964 diff = std::get<double>(result_j) - std::get<double>(result_i);
965 }
else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i)) {
966 diff = std::get<int>(result_j) - std::get<int>(result_i);
968 throw std::runtime_error(
"Bad variant access");
970 if (variablename ==
"phi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)"))) {
971 if (fabs(diff) > M_PI) {
973 diff = diff - 2 * M_PI;
975 diff = 2 * M_PI + diff;
984 B2FATAL(
"Wrong number of arguments for meta function mcDaughterDiffOf");
990 if (arguments.size() == 5) {
991 int iDaughterNumber = 0, jDaughterNumber = 0, agrandDaughterNumber = 0, bgrandDaughterNumber = 0;
993 iDaughterNumber = Belle2::convertString<int>(arguments[0]);
994 jDaughterNumber = Belle2::convertString<int>(arguments[1]);
995 agrandDaughterNumber = Belle2::convertString<int>(arguments[2]);
996 bgrandDaughterNumber = Belle2::convertString<int>(arguments[3]);
997 }
catch (std::invalid_argument&) {
998 B2FATAL(
"First four arguments of grandDaughterDiffOf meta function must be integers!");
1000 auto variablename = arguments[4];
1001 auto func = [variablename, iDaughterNumber, jDaughterNumber, agrandDaughterNumber,
1002 bgrandDaughterNumber](
const Particle * particle) ->
double {
1003 if (particle ==
nullptr)
1005 if (iDaughterNumber >=
int(particle->getNDaughters()) || jDaughterNumber >= int(particle->getNDaughters()))
1007 if (agrandDaughterNumber >=
int((particle->getDaughter(iDaughterNumber))->getNDaughters()) || bgrandDaughterNumber >= int((particle->getDaughter(jDaughterNumber))->getNDaughters()))
1012 double diff = std::get<double>(var->function((particle->getDaughter(jDaughterNumber))->getDaughter(
1013 bgrandDaughterNumber))) - std::get<double>(var->function((particle->getDaughter(iDaughterNumber))->getDaughter(
1014 agrandDaughterNumber)));
1015 if (variablename ==
"phi" or variablename ==
"clusterPhi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)"))
1016 or std::regex_match(variablename, std::regex(
"use.*Frame\\(clusterPhi\\)"))) {
1017 if (fabs(diff) > M_PI) {
1019 diff = diff - 2 * M_PI;
1021 diff = 2 * M_PI + diff;
1030 B2FATAL(
"Wrong number of arguments for meta function grandDaughterDiffOf");
1036 std::vector<std::string> new_arguments = arguments;
1037 new_arguments.push_back(std::string(
"phi"));
1038 return daughterDiffOf(new_arguments);
1043 std::vector<std::string> new_arguments = arguments;
1044 new_arguments.push_back(std::string(
"phi"));
1045 return mcDaughterDiffOf(new_arguments);
1050 std::vector<std::string> new_arguments = arguments;
1051 new_arguments.push_back(std::string(
"phi"));
1052 return grandDaughterDiffOf(new_arguments);
1057 std::vector<std::string> new_arguments = arguments;
1058 new_arguments.push_back(std::string(
"clusterPhi"));
1059 return daughterDiffOf(new_arguments);
1064 std::vector<std::string> new_arguments = arguments;
1065 new_arguments.push_back(std::string(
"clusterPhi"));
1066 return grandDaughterDiffOf(new_arguments);
1071 std::vector<std::string> new_arguments = arguments;
1072 new_arguments.push_back(std::string(
"useCMSFrame(phi)"));
1073 return daughterDiffOf(new_arguments);
1078 std::vector<std::string> new_arguments = arguments;
1079 new_arguments.push_back(std::string(
"useCMSFrame(phi)"));
1080 return mcDaughterDiffOf(new_arguments);
1085 std::vector<std::string> new_arguments = arguments;
1086 new_arguments.push_back(std::string(
"useCMSFrame(clusterPhi)"));
1087 return daughterDiffOf(new_arguments);
1092 if (arguments.size() == 3) {
1093 int iDaughterNumber = 0;
1094 int jDaughterNumber = 0;
1096 iDaughterNumber = Belle2::convertString<int>(arguments[0]);
1097 jDaughterNumber = Belle2::convertString<int>(arguments[1]);
1098 }
catch (std::invalid_argument&) {
1099 B2FATAL(
"First two arguments of daughterDiffOf meta function must be integers!");
1102 auto func = [var, iDaughterNumber, jDaughterNumber](
const Particle * particle) ->
double {
1103 if (particle ==
nullptr)
1105 if (iDaughterNumber >=
int(particle->getNDaughters()) || jDaughterNumber >= int(particle->getNDaughters()))
1109 double iValue, jValue;
1110 if (std::holds_alternative<double>(var->function(particle->getDaughter(jDaughterNumber)))) {
1111 iValue = std::get<double>(var->function(particle->getDaughter(iDaughterNumber)));
1112 jValue = std::get<double>(var->function(particle->getDaughter(jDaughterNumber)));
1113 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(jDaughterNumber)))) {
1114 iValue = std::get<int>(var->function(particle->getDaughter(iDaughterNumber)));
1115 jValue = std::get<int>(var->function(particle->getDaughter(jDaughterNumber)));
1117 return (jValue - iValue) / (jValue + iValue);
1122 B2FATAL(
"Wrong number of arguments for meta function daughterNormDiffOf");
1128 if (arguments.size() == 2) {
1129 int daughterNumber = 0;
1131 daughterNumber = Belle2::convertString<int>(arguments[0]);
1132 }
catch (std::invalid_argument&) {
1133 B2FATAL(
"First argument of daughterMotherDiffOf meta function must be integer!");
1135 auto variablename = arguments[1];
1136 auto func = [variablename, daughterNumber](
const Particle * particle) ->
double {
1137 if (particle ==
nullptr)
1139 if (daughterNumber >=
int(particle->getNDaughters()))
1144 auto result_mother = var->function(particle);
1145 auto result_daughter = var->function(particle->getDaughter(daughterNumber));
1147 if (std::holds_alternative<double>(result_mother) && std::holds_alternative<double>(result_daughter)) {
1148 diff = std::get<double>(result_mother) - std::get<double>(result_daughter);
1149 }
else if (std::holds_alternative<int>(result_mother) && std::holds_alternative<int>(result_daughter)) {
1150 diff = std::get<int>(result_mother) - std::get<int>(result_daughter);
1152 throw std::runtime_error(
"Bad variant access");
1155 if (variablename ==
"phi" or variablename ==
"useCMSFrame(phi)") {
1156 if (fabs(diff) > M_PI) {
1158 diff = diff - 2 * M_PI;
1160 diff = 2 * M_PI + diff;
1169 B2FATAL(
"Wrong number of arguments for meta function daughterMotherDiffOf");
1175 if (arguments.size() == 2) {
1176 int daughterNumber = 0;
1178 daughterNumber = Belle2::convertString<int>(arguments[0]);
1179 }
catch (std::invalid_argument&) {
1180 B2FATAL(
"First argument of daughterMotherDiffOf meta function must be integer!");
1183 auto func = [var, daughterNumber](
const Particle * particle) ->
double {
1184 if (particle ==
nullptr)
1186 if (daughterNumber >=
int(particle->getNDaughters()))
1190 double daughterValue = 0.0, motherValue = 0.0;
1191 if (std::holds_alternative<double>(var->function(particle))) {
1192 daughterValue = std::get<double>(var->function(particle->getDaughter(daughterNumber)));
1193 motherValue = std::get<double>(var->function(particle));
1194 }
else if (std::holds_alternative<int>(var->function(particle))) {
1195 daughterValue = std::get<int>(var->function(particle->getDaughter(daughterNumber)));
1196 motherValue = std::get<int>(var->function(particle));
1198 return (motherValue - daughterValue) / (motherValue + daughterValue);
1203 B2FATAL(
"Wrong number of arguments for meta function daughterMotherNormDiffOf");
1209 if (arguments.size() == 2 || arguments.size() == 3) {
1211 auto func = [arguments](
const Particle * particle) ->
double {
1212 if (particle ==
nullptr)
1215 std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1219 for (
auto& generalizedIndex : arguments)
1221 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1223 pDaus.push_back(frame.getMomentum(dauPart));
1225 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1231 if (pDaus.size() == 2)
1238 B2FATAL(
"Wrong number of arguments for meta function daughterAngle");
1242 double grandDaughterDecayAngle(
const Particle* particle,
const std::vector<double>& arguments)
1244 if (arguments.size() == 2) {
1249 int daughterIndex = std::lround(arguments[0]);
1250 if (daughterIndex >=
int(particle->getNDaughters()))
1252 const Particle* dau = particle->getDaughter(daughterIndex);
1254 int grandDaughterIndex = std::lround(arguments[1]);
1255 if (grandDaughterIndex >=
int(dau->getNDaughters()))
1258 B2Vector3D boost = dau->get4Vector().BoostToCM();
1260 ROOT::Math::PxPyPzEVector motherMomentum = - particle->get4Vector();
1261 motherMomentum = ROOT::Math::Boost(boost) * motherMomentum;
1263 ROOT::Math::PxPyPzEVector grandDaughterMomentum = dau->getDaughter(grandDaughterIndex)->get4Vector();
1264 grandDaughterMomentum = ROOT::Math::Boost(boost) * grandDaughterMomentum;
1269 B2FATAL(
"The variable grandDaughterDecayAngle needs exactly two integers as arguments!");
1275 if (arguments.size() == 2 || arguments.size() == 3) {
1277 auto func = [arguments](
const Particle * particle) ->
double {
1278 if (particle ==
nullptr)
1281 std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1285 if (particle->getParticleSource() == Particle::EParticleSourceObject::c_MCParticle)
1287 for (
auto& generalizedIndex : arguments) {
1288 const MCParticle* mcPart = particle->getMCParticle();
1289 if (mcPart ==
nullptr)
1292 if (dauMcPart ==
nullptr)
1295 pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1299 for (
auto& generalizedIndex : arguments) {
1300 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1301 if (dauPart ==
nullptr)
1305 if (dauMcPart ==
nullptr)
1308 pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1313 if (pDaus.size() == 2)
1320 B2FATAL(
"Wrong number of arguments for meta function mcDaughterAngle");
1324 double daughterClusterAngleInBetween(
const Particle* particle,
const std::vector<double>& daughterIndices)
1326 if (daughterIndices.size() == 2) {
1327 int daughterIndexi = std::lround(daughterIndices[0]);
1328 int daughterIndexj = std::lround(daughterIndices[1]);
1329 if (std::max(daughterIndexi, daughterIndexj) >=
int(particle->getNDaughters())) {
1332 const ECLCluster* clusteri = particle->getDaughter(daughterIndexi)->getECLCluster();
1333 const ECLCluster* clusterj = particle->getDaughter(daughterIndexj)->getECLCluster();
1334 if (clusteri and clusterj) {
1338 ClusterUtils clusutils;
1339 B2Vector3D pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit)).Vect();
1340 B2Vector3D pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit)).Vect();
1341 return pi.Angle(pj);
1345 }
else if (daughterIndices.size() == 3) {
1346 int daughterIndexi = std::lround(daughterIndices[0]);
1347 int daughterIndexj = std::lround(daughterIndices[1]);
1348 int daughterIndexk = std::lround(daughterIndices[2]);
1349 if (std::max(std::max(daughterIndexi, daughterIndexj), daughterIndexk) >=
int(particle->getNDaughters())) {
1352 const ECLCluster* clusteri = (particle->getDaughter(daughterIndices[0]))->getECLCluster();
1353 const ECLCluster* clusterj = (particle->getDaughter(daughterIndices[1]))->getECLCluster();
1354 const ECLCluster* clusterk = (particle->getDaughter(daughterIndices[2]))->getECLCluster();
1355 if (clusteri and clusterj and clusterk) {
1360 ClusterUtils clusutils;
1361 B2Vector3D pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit)).Vect();
1362 B2Vector3D pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit)).Vect();
1363 B2Vector3D pk = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterk, clusterkBit)).Vect();
1364 return pk.
Angle(pi + pj);
1369 B2FATAL(
"Wrong number of arguments for daughterClusterAngleInBetween!");
1375 if (arguments.size() > 1) {
1376 auto func = [arguments](
const Particle * particle) ->
double {
1378 ROOT::Math::PxPyPzEVector pSum;
1380 for (
auto& generalizedIndex : arguments)
1382 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1384 pSum += frame.getMomentum(dauPart);
1393 B2FATAL(
"Wrong number of arguments for meta function daughterInvM. At least two integers are needed.");
1399 if (arguments.size() == 2) {
1403 divideBy = Belle2::convertString<int>(arguments[1]);
1404 }
catch (std::invalid_argument&) {
1405 B2FATAL(
"Second argument of modulo meta function must be integer!");
1407 auto func = [var, divideBy](
const Particle * particle) ->
int {
1408 auto var_result = var->function(particle);
1409 if (std::holds_alternative<double>(var_result))
1411 return int(std::get<double>(var_result)) % divideBy;
1412 }
else if (std::holds_alternative<int>(var_result))
1414 return std::get<int>(var_result) % divideBy;
1415 }
else if (std::holds_alternative<bool>(var_result))
1417 return int(std::get<bool>(var_result)) % divideBy;
1422 B2FATAL(
"Wrong number of arguments for meta function modulo");
1428 if (arguments.size() == 1) {
1431 auto func = [var](
const Particle * particle) ->
bool {
return std::isnan(std::get<double>(var->function(particle))); };
1434 B2FATAL(
"Wrong number of arguments for meta function isNAN");
1440 if (arguments.size() == 2) {
1442 double defaultOutput;
1444 defaultOutput = Belle2::convertString<double>(arguments[1]);
1445 }
catch (std::invalid_argument&) {
1446 B2FATAL(
"The second argument of ifNANgiveX meta function must be a number!");
1448 auto func = [var, defaultOutput](
const Particle * particle) ->
double {
1449 double output = std::get<double>(var->function(particle));
1450 if (std::isnan(output))
return defaultOutput;
1455 B2FATAL(
"Wrong number of arguments for meta function ifNANgiveX");
1461 if (arguments.size() == 1) {
1464 auto func = [var](
const Particle * particle) ->
bool {
return std::isinf(std::get<double>(var->function(particle))); };
1467 B2FATAL(
"Wrong number of arguments for meta function isInfinity");
1473 if (arguments.size() >= 2) {
1479 for (
size_t i = 1; i < arguments.size(); ++i) {
1481 finalMask |= Belle2::convertString<int>(arguments[i]);
1482 }
catch (std::invalid_argument&) {
1483 B2FATAL(
"The input flags to meta function unmask() should be integer!");
1489 auto func = [var, finalMask](
const Particle * particle) ->
double {
1491 auto var_result = var->function(particle);
1492 if (std::holds_alternative<double>(var_result))
1495 if (std::isnan(std::get<double>(var_result))) {
1498 value = int(std::get<double>(var_result));
1499 }
else if (std::holds_alternative<int>(var_result))
1501 value = std::get<int>(var_result);
1505 value &= (~finalMask);
1512 B2FATAL(
"Meta function unmask needs at least two arguments!");
1518 if (arguments.size() == 3) {
1520 std::string cutString = arguments[0];
1526 auto func = [cut, variableIfTrue, variableIfFalse](
const Particle * particle) ->
double {
1527 if (particle ==
nullptr)
1529 if (cut->check(particle))
1531 auto var_result = variableIfTrue->function(particle);
1532 if (std::holds_alternative<double>(var_result)) {
1533 return std::get<double>(var_result);
1534 }
else if (std::holds_alternative<int>(var_result)) {
1535 return std::get<int>(var_result);
1536 }
else if (std::holds_alternative<bool>(var_result)) {
1537 return std::get<bool>(var_result);
1541 auto var_result = variableIfFalse->function(particle);
1542 if (std::holds_alternative<double>(var_result)) {
1543 return std::get<double>(var_result);
1544 }
else if (std::holds_alternative<int>(var_result)) {
1545 return std::get<int>(var_result);
1546 }
else if (std::holds_alternative<bool>(var_result)) {
1547 return std::get<bool>(var_result);
1554 B2FATAL(
"Wrong number of arguments for meta function conditionalVariableSelector");
1560 if (arguments.size() > 0) {
1561 std::vector<const Variable::Manager::Var*> variables;
1562 for (
auto& argument : arguments)
1565 auto func = [variables, arguments](
const Particle * particle) ->
double {
1566 double pValueProduct = 1.;
1567 for (
auto variable : variables)
1569 double pValue = std::get<double>(variable->function(particle));
1573 pValueProduct *= pValue;
1575 double pValueSum = 1.;
1576 double factorial = 1.;
1577 for (
unsigned int i = 1; i < arguments.size(); ++i)
1580 pValueSum += pow(-std::log(pValueProduct), i) / factorial;
1582 return pValueProduct * pValueSum;
1586 B2FATAL(
"Wrong number of arguments for meta function pValueCombination");
1592 if (arguments.size() == 1) {
1594 auto func = [var](
const Particle * particle) ->
double {
1595 auto var_result = var->function(particle);
1596 if (std::holds_alternative<double>(var_result))
1598 return std::abs(std::get<double>(var_result));
1599 }
else if (std::holds_alternative<int>(var_result))
1601 return std::abs(std::get<int>(var_result));
1606 B2FATAL(
"Wrong number of arguments for meta function abs");
1612 if (arguments.size() == 2) {
1617 B2FATAL(
"One or both of the used variables doesn't exist!");
1619 auto func = [var1, var2](
const Particle * particle) ->
double {
1621 auto var_result1 = var1->function(particle);
1622 auto var_result2 = var2->function(particle);
1623 if (std::holds_alternative<double>(var_result1))
1625 val1 = std::get<double>(var_result1);
1626 }
else if (std::holds_alternative<int>(var_result1))
1628 val1 = std::get<int>(var_result1);
1630 if (std::holds_alternative<double>(var_result2))
1632 val2 = std::get<double>(var_result2);
1633 }
else if (std::holds_alternative<int>(var_result2))
1635 val2 = std::get<int>(var_result2);
1637 return std::max(val1, val2);
1641 B2FATAL(
"Wrong number of arguments for meta function max");
1647 if (arguments.size() == 2) {
1652 B2FATAL(
"One or both of the used variables doesn't exist!");
1654 auto func = [var1, var2](
const Particle * particle) ->
double {
1656 auto var_result1 = var1->function(particle);
1657 auto var_result2 = var2->function(particle);
1658 if (std::holds_alternative<double>(var_result1))
1660 val1 = std::get<double>(var_result1);
1661 }
else if (std::holds_alternative<int>(var_result1))
1663 val1 = std::get<int>(var_result1);
1665 if (std::holds_alternative<double>(var_result2))
1667 val2 = std::get<double>(var_result2);
1668 }
else if (std::holds_alternative<int>(var_result2))
1670 val2 = std::get<int>(var_result2);
1672 return std::min(val1, val2);
1676 B2FATAL(
"Wrong number of arguments for meta function min");
1682 if (arguments.size() == 1) {
1684 auto func = [var](
const Particle * particle) ->
double {
1685 auto var_result = var->function(particle);
1686 if (std::holds_alternative<double>(var_result))
1687 return std::sin(std::get<double>(var_result));
1688 else if (std::holds_alternative<int>(var_result))
1689 return std::sin(std::get<int>(var_result));
1694 B2FATAL(
"Wrong number of arguments for meta function sin");
1700 if (arguments.size() == 1) {
1702 auto func = [var](
const Particle * particle) ->
double {
1703 auto var_result = var->function(particle);
1704 if (std::holds_alternative<double>(var_result))
1705 return std::asin(std::get<double>(var_result));
1706 else if (std::holds_alternative<int>(var_result))
1707 return std::asin(std::get<int>(var_result));
1712 B2FATAL(
"Wrong number of arguments for meta function asin");
1718 if (arguments.size() == 1) {
1720 auto func = [var](
const Particle * particle) ->
double {
1721 auto var_result = var->function(particle);
1722 if (std::holds_alternative<double>(var_result))
1723 return std::cos(std::get<double>(var_result));
1724 else if (std::holds_alternative<int>(var_result))
1725 return std::cos(std::get<int>(var_result));
1730 B2FATAL(
"Wrong number of arguments for meta function cos");
1736 if (arguments.size() == 1) {
1738 auto func = [var](
const Particle * particle) ->
double {
1739 auto var_result = var->function(particle);
1740 if (std::holds_alternative<double>(var_result))
1741 return std::acos(std::get<double>(var_result));
1742 else if (std::holds_alternative<int>(var_result))
1743 return std::acos(std::get<int>(var_result));
1748 B2FATAL(
"Wrong number of arguments for meta function acos");
1754 if (arguments.size() == 1) {
1756 auto func = [var](
const Particle * particle) ->
double {
return std::tan(std::get<double>(var->function(particle))); };
1759 B2FATAL(
"Wrong number of arguments for meta function tan");
1765 if (arguments.size() == 1) {
1767 auto func = [var](
const Particle * particle) ->
double {
return std::atan(std::get<double>(var->function(particle))); };
1770 B2FATAL(
"Wrong number of arguments for meta function atan");
1776 if (arguments.size() == 1) {
1778 auto func = [var](
const Particle * particle) ->
double {
1779 auto var_result = var->function(particle);
1780 if (std::holds_alternative<double>(var_result))
1781 return std::exp(std::get<double>(var_result));
1782 else if (std::holds_alternative<int>(var_result))
1783 return std::exp(std::get<int>(var_result));
1788 B2FATAL(
"Wrong number of arguments for meta function exp");
1794 if (arguments.size() == 1) {
1796 auto func = [var](
const Particle * particle) ->
double {
1797 auto var_result = var->function(particle);
1798 if (std::holds_alternative<double>(var_result))
1799 return std::log(std::get<double>(var_result));
1800 else if (std::holds_alternative<int>(var_result))
1801 return std::log(std::get<int>(var_result));
1806 B2FATAL(
"Wrong number of arguments for meta function log");
1812 if (arguments.size() == 1) {
1814 auto func = [var](
const Particle * particle) ->
double {
1815 auto var_result = var->function(particle);
1816 if (std::holds_alternative<double>(var_result))
1817 return std::log10(std::get<double>(var_result));
1818 else if (std::holds_alternative<int>(var_result))
1819 return std::log10(std::get<int>(var_result));
1824 B2FATAL(
"Wrong number of arguments for meta function log10");
1830 if (arguments.size() == 1) {
1832 auto func = [var](
const Particle * particle) ->
double {
1833 if (particle ==
nullptr)
1836 StoreArray<Particle> particles;
1837 if (!particle->hasExtraInfo(
"original_index"))
1840 auto originalParticle = particles[particle->getExtraInfo(
"original_index")];
1841 if (!originalParticle)
1843 auto var_result = var->function(originalParticle);
1844 if (std::holds_alternative<double>(var_result))
1846 return std::get<double>(var_result);
1847 }
else if (std::holds_alternative<int>(var_result))
1849 return std::get<int>(var_result);
1850 }
else if (std::holds_alternative<bool>(var_result))
1852 return std::get<bool>(var_result);
1857 B2FATAL(
"Wrong number of arguments for meta function originalParticle");
1863 if (arguments.size() == 2) {
1864 int daughterNumber = 0;
1866 daughterNumber = Belle2::convertString<int>(arguments[0]);
1867 }
catch (std::invalid_argument&) {
1868 B2FATAL(
"First argument of daughter meta function must be integer!");
1871 auto func = [var, daughterNumber](
const Particle * particle) ->
double {
1872 if (particle ==
nullptr)
1874 if (daughterNumber >=
int(particle->getNDaughters()))
1878 auto var_result = var->function(particle->getDaughter(daughterNumber));
1879 if (std::holds_alternative<double>(var_result)) {
1880 return std::get<double>(var_result);
1881 }
else if (std::holds_alternative<int>(var_result)) {
1882 return std::get<int>(var_result);
1883 }
else if (std::holds_alternative<bool>(var_result)) {
1884 return std::get<bool>(var_result);
1890 B2FATAL(
"Wrong number of arguments for meta function daughter");
1896 if (arguments.size() == 2) {
1897 int daughterNumber = 0;
1899 daughterNumber = Belle2::convertString<int>(arguments[0]);
1900 }
catch (std::invalid_argument&) {
1901 B2FATAL(
"First argument of daughter meta function must be integer!");
1904 auto func = [var, daughterNumber](
const Particle * particle) ->
double {
1905 if (particle ==
nullptr)
1907 if (daughterNumber >=
int(particle->getNDaughters()))
1911 StoreArray<Particle> particles;
1912 if (!particle->getDaughter(daughterNumber)->hasExtraInfo(
"original_index"))
1914 auto originalDaughter = particles[particle->getDaughter(daughterNumber)->getExtraInfo(
"original_index")];
1915 if (!originalDaughter)
1918 auto var_result = var->function(originalDaughter);
1919 if (std::holds_alternative<double>(var_result)) {
1920 return std::get<double>(var_result);
1921 }
else if (std::holds_alternative<int>(var_result)) {
1922 return std::get<int>(var_result);
1923 }
else if (std::holds_alternative<bool>(var_result)) {
1924 return std::get<bool>(var_result);
1930 B2FATAL(
"Wrong number of arguments for meta function daughter");
1936 if (arguments.size() == 2) {
1937 int daughterNumber = 0;
1939 daughterNumber = Belle2::convertString<int>(arguments[0]);
1940 }
catch (std::invalid_argument&) {
1941 B2FATAL(
"First argument of mcDaughter meta function must be integer!");
1944 auto func = [var, daughterNumber](
const Particle * particle) ->
double {
1945 if (particle ==
nullptr)
1947 if (particle->getMCParticle())
1949 if (daughterNumber >=
int(particle->getMCParticle()->getNDaughters())) {
1952 Particle tempParticle = Particle(particle->getMCParticle()->getDaughters().at(daughterNumber));
1953 auto var_result = var->function(&tempParticle);
1954 if (std::holds_alternative<double>(var_result)) {
1955 return std::get<double>(var_result);
1956 }
else if (std::holds_alternative<int>(var_result)) {
1957 return std::get<int>(var_result);
1958 }
else if (std::holds_alternative<bool>(var_result)) {
1959 return std::get<bool>(var_result);
1968 B2FATAL(
"Wrong number of arguments for meta function mcDaughter");
1974 if (arguments.size() == 1) {
1976 auto func = [var](
const Particle * particle) ->
double {
1977 if (particle ==
nullptr)
1979 if (particle->getMCParticle())
1981 if (particle->getMCParticle()->getMother() ==
nullptr) {
1984 Particle tempParticle = Particle(particle->getMCParticle()->getMother());
1985 auto var_result = var->function(&tempParticle);
1986 if (std::holds_alternative<double>(var_result)) {
1987 return std::get<double>(var_result);
1988 }
else if (std::holds_alternative<int>(var_result)) {
1989 return std::get<int>(var_result);
1990 }
else if (std::holds_alternative<bool>(var_result)) {
1991 return std::get<bool>(var_result);
2000 B2FATAL(
"Wrong number of arguments for meta function mcMother");
2006 if (arguments.size() == 2) {
2007 int particleNumber = 0;
2009 particleNumber = Belle2::convertString<int>(arguments[0]);
2010 }
catch (std::invalid_argument&) {
2011 B2FATAL(
"First argument of genParticle meta function must be integer!");
2015 auto func = [var, particleNumber](
const Particle*) ->
double {
2016 StoreArray<MCParticle> mcParticles(
"MCParticles");
2017 if (particleNumber >= mcParticles.getEntries())
2022 MCParticle* mcParticle = mcParticles[particleNumber];
2023 Particle part = Particle(mcParticle);
2024 auto var_result = var->function(&part);
2025 if (std::holds_alternative<double>(var_result))
2027 return std::get<double>(var_result);
2028 }
else if (std::holds_alternative<int>(var_result))
2030 return std::get<int>(var_result);
2031 }
else if (std::holds_alternative<bool>(var_result))
2033 return std::get<bool>(var_result);
2038 B2FATAL(
"Wrong number of arguments for meta function genParticle");
2044 if (arguments.size() == 1) {
2047 auto func = [var](
const Particle*) ->
double {
2048 StoreArray<MCParticle> mcParticles(
"MCParticles");
2049 if (mcParticles.getEntries() == 0)
2054 MCParticle* mcUpsilon4S = mcParticles[0];
2055 if (mcUpsilon4S->isInitial()) mcUpsilon4S = mcParticles[2];
2056 if (mcUpsilon4S->getPDG() != 300553)
2061 Particle upsilon4S = Particle(mcUpsilon4S);
2062 auto var_result = var->function(&upsilon4S);
2063 if (std::holds_alternative<double>(var_result))
2065 return std::get<double>(var_result);
2066 }
else if (std::holds_alternative<int>(var_result))
2068 return std::get<int>(var_result);
2069 }
else if (std::holds_alternative<bool>(var_result))
2071 return std::get<bool>(var_result);
2076 B2FATAL(
"Wrong number of arguments for meta function genUpsilon4S");
2082 if (arguments.size() == 4) {
2083 std::string listName = arguments[0];
2084 std::string rankedVariableName = arguments[1];
2085 std::string returnVariableName = arguments[2];
2086 std::string extraInfoName = rankedVariableName +
"_rank";
2089 rank = Belle2::convertString<int>(arguments[3]);
2090 }
catch (std::invalid_argument&) {
2091 B2ERROR(
"3rd argument of getVariableByRank meta function (Rank) must be an integer!");
2096 auto func = [var, rank, extraInfoName, listName](
const Particle*)->
double {
2097 StoreObjPtr<ParticleList> list(listName);
2099 const unsigned int numParticles = list->getListSize();
2100 for (
unsigned int i = 0; i < numParticles; i++)
2102 const Particle* p = list->getParticle(i);
2103 if (p->getExtraInfo(extraInfoName) == rank) {
2104 auto var_result = var->function(p);
2105 if (std::holds_alternative<double>(var_result)) {
2106 return std::get<double>(var_result);
2107 }
else if (std::holds_alternative<int>(var_result)) {
2108 return std::get<int>(var_result);
2109 }
else if (std::holds_alternative<bool>(var_result)) {
2110 return std::get<bool>(var_result);
2115 return std::numeric_limits<double>::signaling_NaN();
2119 B2FATAL(
"Wrong number of arguments for meta function getVariableByRank");
2125 if (arguments.size() == 1 or arguments.size() == 2) {
2127 std::string listName = arguments[0];
2128 std::string cutString =
"";
2130 if (arguments.size() == 2) {
2131 cutString = arguments[1];
2136 auto func = [listName, cut](
const Particle*) ->
int {
2138 StoreObjPtr<ParticleList> list(listName);
2140 for (
unsigned int i = 0; i < list->getListSize(); i++)
2142 const Particle* particle = list->getParticle(i);
2143 if (cut->check(particle)) {
2151 B2FATAL(
"Wrong number of arguments for meta function countInList");
2157 if (arguments.size() == 2 or arguments.size() == 3) {
2159 std::string roeListName = arguments[0];
2160 std::string cutString = arguments[1];
2162 if (arguments.size() == 2) {
2163 B2INFO(
"Use pdgCode of electron as default in meta variable veto, other arguments: " << roeListName <<
", " << cutString);
2166 pdgCode = Belle2::convertString<int>(arguments[2]);;
2167 }
catch (std::invalid_argument&) {
2168 B2FATAL(
"Third argument of veto meta function must be integer!");
2175 auto func = [roeListName, cut, pdgCode, flavourType](
const Particle * particle) ->
bool {
2176 StoreObjPtr<ParticleList> roeList(roeListName);
2177 ROOT::Math::PxPyPzEVector vec = particle->get4Vector();
2178 for (
unsigned int i = 0; i < roeList->getListSize(); i++)
2180 const Particle* roeParticle = roeList->getParticle(i);
2181 if (not particle->overlapsWith(roeParticle)) {
2182 ROOT::Math::PxPyPzEVector tempCombination = roeParticle->get4Vector() + vec;
2183 std::vector<int> indices = { particle->getArrayIndex(), roeParticle->getArrayIndex() };
2184 Particle tempParticle = Particle(tempCombination, pdgCode, flavourType, indices, particle->getArrayPointer());
2185 if (cut->check(&tempParticle)) {
2194 B2FATAL(
"Wrong number of arguments for meta function veto");
2200 if (arguments.size() == 1) {
2201 std::string cutString = arguments[0];
2203 auto func = [cut](
const Particle * particle) ->
int {
2205 for (
auto& daughter : particle->getDaughters())
2207 if (cut->check(daughter))
2214 B2FATAL(
"Wrong number of arguments for meta function countDaughters");
2220 if (arguments.size() == 1) {
2221 std::string cutString = arguments[0];
2223 auto func = [cut](
const Particle * particle) ->
int {
2225 std::vector<const Particle*> fspDaughters;
2226 particle->fillFSPDaughters(fspDaughters);
2229 for (
auto& daughter : fspDaughters)
2231 if (cut->check(daughter))
2238 B2FATAL(
"Wrong number of arguments for meta function countFSPDaughters");
2244 if (arguments.size() == 1) {
2245 std::string cutString = arguments[0];
2247 auto func = [cut](
const Particle * particle) ->
int {
2249 std::vector<const Particle*> allDaughters;
2250 particle->fillAllDaughters(allDaughters);
2253 for (
auto& daughter : allDaughters)
2255 if (cut->check(daughter))
2262 B2FATAL(
"Wrong number of arguments for meta function countDescendants");
2266 Manager::FunctionPtr numberOfNonOverlappingParticles(
const std::vector<std::string>& arguments)
2269 auto func = [arguments](
const Particle * particle) ->
int {
2271 int _numberOfNonOverlappingParticles = 0;
2272 for (
const auto& listName : arguments)
2274 StoreObjPtr<ParticleList> list(listName);
2275 if (not list.isValid()) {
2276 B2FATAL(
"Invalid list named " << listName <<
" encountered in numberOfNonOverlappingParticles.");
2278 for (
unsigned int i = 0; i < list->getListSize(); i++) {
2279 const Particle* p = list->getParticle(i);
2280 if (not particle->overlapsWith(p)) {
2281 _numberOfNonOverlappingParticles++;
2285 return _numberOfNonOverlappingParticles;
2294 if (arguments.size() == 1) {
2296 auto func = [var](
const Particle * particle) ->
double {
2297 const MCParticle* mcp = particle->getMCParticle();
2302 Particle tmpPart(mcp);
2303 auto var_result = var->function(&tmpPart);
2304 if (std::holds_alternative<double>(var_result))
2306 return std::get<double>(var_result);
2307 }
else if (std::holds_alternative<int>(var_result))
2309 return std::get<int>(var_result);
2310 }
else if (std::holds_alternative<bool>(var_result))
2312 return std::get<bool>(var_result);
2317 B2FATAL(
"Wrong number of arguments for meta function matchedMC");
2323 if (arguments.size() == 1) {
2326 auto func = [var](
const Particle * particle) ->
double {
2328 const ECLCluster* cluster = particle->getECLCluster();
2331 auto mcps = cluster->getRelationsTo<MCParticle>();
2334 std::vector<std::pair<double, int>> weightsAndIndices;
2335 for (
unsigned int i = 0; i < mcps.size(); ++i)
2336 weightsAndIndices.emplace_back(mcps.weight(i), i);
2339 std::sort(weightsAndIndices.begin(), weightsAndIndices.end(),
2340 ValueIndexPairSorting::higherPair<decltype(weightsAndIndices)::value_type>);
2343 const MCParticle* mcp = mcps.object(weightsAndIndices[0].second);
2345 Particle tmpPart(mcp);
2346 auto var_result = var->function(&tmpPart);
2347 if (std::holds_alternative<double>(var_result))
2349 return std::get<double>(var_result);
2350 }
else if (std::holds_alternative<int>(var_result))
2352 return std::get<int>(var_result);
2353 }
else if (std::holds_alternative<bool>(var_result))
2355 return std::get<bool>(var_result);
2364 B2FATAL(
"Wrong number of arguments for meta function clusterBestMatchedMCParticle");
2368 double matchedMCHasPDG(
const Particle* particle,
const std::vector<double>& pdgCode)
2370 if (pdgCode.size() != 1) {
2371 B2FATAL(
"Too many arguments provided to matchedMCHasPDG!");
2373 int inputPDG = std::lround(pdgCode[0]);
2375 const MCParticle* mcp = particle->getMCParticle();
2379 return std::abs(mcp->getPDG()) == inputPDG;
2384 if (arguments.size() == 1) {
2385 std::string listName = arguments[0];
2386 auto func = [listName](
const Particle * particle) ->
double {
2389 StoreObjPtr<ParticleList> listOfParticles(listName);
2391 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalEnergyOfParticlesInList");
2392 double totalEnergy = 0;
2393 int nParticles = listOfParticles->getListSize();
2394 for (
int i = 0; i < nParticles; i++)
2396 const Particle* part = listOfParticles->getParticle(i);
2398 totalEnergy += frame.getMomentum(part).E();
2405 B2FATAL(
"Wrong number of arguments for meta function totalEnergyOfParticlesInList");
2411 if (arguments.size() == 1) {
2412 std::string listName = arguments[0];
2413 auto func = [listName](
const Particle*) ->
double {
2414 StoreObjPtr<ParticleList> listOfParticles(listName);
2416 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPxOfParticlesInList");
2418 int nParticles = listOfParticles->getListSize();
2420 for (
int i = 0; i < nParticles; i++)
2422 const Particle* part = listOfParticles->getParticle(i);
2423 totalPx += frame.getMomentum(part).Px();
2429 B2FATAL(
"Wrong number of arguments for meta function totalPxOfParticlesInList");
2435 if (arguments.size() == 1) {
2436 std::string listName = arguments[0];
2437 auto func = [listName](
const Particle*) ->
double {
2438 StoreObjPtr<ParticleList> listOfParticles(listName);
2440 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPyOfParticlesInList");
2442 int nParticles = listOfParticles->getListSize();
2444 for (
int i = 0; i < nParticles; i++)
2446 const Particle* part = listOfParticles->getParticle(i);
2447 totalPy += frame.getMomentum(part).Py();
2453 B2FATAL(
"Wrong number of arguments for meta function totalPyOfParticlesInList");
2459 if (arguments.size() == 1) {
2460 std::string listName = arguments[0];
2461 auto func = [listName](
const Particle*) ->
double {
2462 StoreObjPtr<ParticleList> listOfParticles(listName);
2464 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPzOfParticlesInList");
2466 int nParticles = listOfParticles->getListSize();
2468 for (
int i = 0; i < nParticles; i++)
2470 const Particle* part = listOfParticles->getParticle(i);
2471 totalPz += frame.getMomentum(part).Pz();
2477 B2FATAL(
"Wrong number of arguments for meta function totalPzOfParticlesInList");
2483 if (arguments.size() > 0) {
2485 auto func = [arguments](
const Particle * particle) ->
double {
2487 ROOT::Math::PxPyPzEVector total4Vector;
2489 std::vector<Particle*> particlePool;
2492 for (
const auto& argument : arguments)
2494 StoreObjPtr <ParticleList> listOfParticles(argument);
2496 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << argument <<
" given to invMassInLists");
2497 int nParticles = listOfParticles->getListSize();
2498 for (
int i = 0; i < nParticles; i++) {
2499 bool overlaps =
false;
2500 Particle* part = listOfParticles->getParticle(i);
2501 for (
auto poolPart : particlePool) {
2502 if (part->overlapsWith(poolPart)) {
2508 total4Vector += part->get4Vector();
2509 particlePool.push_back(part);
2513 double invariantMass = total4Vector.M();
2514 return invariantMass;
2519 B2FATAL(
"Wrong number of arguments for meta function invMassInLists");
2523 Manager::FunctionPtr totalECLEnergyOfParticlesInList(
const std::vector<std::string>& arguments)
2525 if (arguments.size() == 1) {
2526 std::string listName = arguments[0];
2527 auto func = [listName](
const Particle * particle) ->
double {
2530 StoreObjPtr<ParticleList> listOfParticles(listName);
2532 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalEnergyOfParticlesInList");
2533 double totalEnergy = 0;
2534 int nParticles = listOfParticles->getListSize();
2535 for (
int i = 0; i < nParticles; i++)
2537 const Particle* part = listOfParticles->getParticle(i);
2538 const ECLCluster* cluster = part->getECLCluster();
2540 if (cluster !=
nullptr) {
2541 totalEnergy += cluster->getEnergy(clusterHypothesis);
2549 B2FATAL(
"Wrong number of arguments for meta function totalECLEnergyOfParticlesInList");
2555 if (arguments.size() == 1) {
2556 std::string listName = arguments[0];
2557 auto func = [listName](
const Particle*) ->
double {
2558 StoreObjPtr<ParticleList> listOfParticles(listName);
2560 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to maxPtInList");
2561 int nParticles = listOfParticles->getListSize();
2564 for (
int i = 0; i < nParticles; i++)
2566 const Particle* part = listOfParticles->getParticle(i);
2567 const double Pt = frame.getMomentum(part).Pt();
2568 if (Pt > maxPt) maxPt = Pt;
2574 B2FATAL(
"Wrong number of arguments for meta function maxPtInList");
2578 Manager::FunctionPtr eclClusterTrackMatchedWithCondition(
const std::vector<std::string>& arguments)
2580 if (arguments.size() <= 1) {
2582 std::string cutString;
2583 if (arguments.size() == 1)
2584 cutString = arguments[0];
2586 auto func = [cut](
const Particle * particle) ->
double {
2588 if (particle ==
nullptr)
2591 const ECLCluster* cluster = particle->getECLCluster();
2595 auto tracks = cluster->getRelationsFrom<Track>();
2597 for (
const auto& track : tracks) {
2600 if (cut->check(&trackParticle))
2609 B2FATAL(
"Wrong number of arguments for meta function eclClusterSpecialTrackMatched");
2615 if (arguments.size() == 2) {
2616 std::string listName = arguments[0];
2619 auto func = [listName, var](
const Particle*) ->
double {
2620 StoreObjPtr<ParticleList> listOfParticles(listName);
2622 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to averageValueInList");
2623 int nParticles = listOfParticles->getListSize();
2624 if (nParticles == 0)
2629 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2631 for (
int i = 0; i < nParticles; i++) {
2632 average += std::get<double>(var->function(listOfParticles->getParticle(i))) / nParticles;
2634 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2636 for (
int i = 0; i < nParticles; i++) {
2637 average += std::get<int>(var->function(listOfParticles->getParticle(i))) / nParticles;
2644 B2FATAL(
"Wrong number of arguments for meta function averageValueInList");
2650 if (arguments.size() == 2) {
2651 std::string listName = arguments[0];
2654 auto func = [listName, var](
const Particle*) ->
double {
2655 StoreObjPtr<ParticleList> listOfParticles(listName);
2657 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to medianValueInList");
2658 int nParticles = listOfParticles->getListSize();
2659 if (nParticles == 0)
2663 std::vector<double> valuesInList;
2664 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2666 for (
int i = 0; i < nParticles; i++) {
2667 valuesInList.push_back(std::get<double>(var->function(listOfParticles->getParticle(i))));
2669 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2671 for (
int i = 0; i < nParticles; i++) {
2672 valuesInList.push_back(std::get<int>(var->function(listOfParticles->getParticle(i))));
2675 std::sort(valuesInList.begin(), valuesInList.end());
2676 if (nParticles % 2 != 0)
2678 return valuesInList[nParticles / 2];
2681 return 0.5 * (valuesInList[nParticles / 2] + valuesInList[nParticles / 2 - 1]);
2686 B2FATAL(
"Wrong number of arguments for meta function medianValueInList");
2693 if (arguments.size() != 1)
2694 B2FATAL(
"Wrong number of arguments for meta function angleToClosestInList");
2696 std::string listname = arguments[0];
2698 auto func = [listname](
const Particle * particle) ->
double {
2700 StoreObjPtr<ParticleList> list(listname);
2701 if (not list.isValid())
2702 B2FATAL(
"Invalid particle list name " << listname <<
" given to angleToClosestInList");
2705 if (list->getListSize() == 0)
2710 const auto p_this =
B2Vector3D(frame.getMomentum(particle).Vect());
2713 double minAngle = 2 * M_PI;
2714 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2716 const Particle* compareme = list->getParticle(i);
2717 const auto p_compare =
B2Vector3D(frame.getMomentum(compareme).Vect());
2718 double angle = p_compare.Angle(p_this);
2719 if (minAngle > angle) minAngle = angle;
2729 if (arguments.size() != 2)
2730 B2FATAL(
"Wrong number of arguments for meta function closestInList");
2732 std::string listname = arguments[0];
2737 auto func = [listname, var](
const Particle * particle) ->
double {
2739 StoreObjPtr<ParticleList> list(listname);
2740 if (not list.isValid())
2741 B2FATAL(
"Invalid particle list name " << listname <<
" given to closestInList");
2745 const auto p_this =
B2Vector3D(frame.getMomentum(particle).Vect());
2748 double minAngle = 2 * M_PI;
2750 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2752 const Particle* compareme = list->getParticle(i);
2753 const auto p_compare =
B2Vector3D(frame.getMomentum(compareme).Vect());
2754 double angle = p_compare.Angle(p_this);
2755 if (minAngle > angle) {
2763 auto var_result = var->function(list->getParticle(iClosest));
2764 if (std::holds_alternative<double>(var_result))
2766 return std::get<double>(var_result);
2767 }
else if (std::holds_alternative<int>(var_result))
2769 return std::get<int>(var_result);
2770 }
else if (std::holds_alternative<bool>(var_result))
2772 return std::get<bool>(var_result);
2781 if (arguments.size() != 1)
2782 B2FATAL(
"Wrong number of arguments for meta function angleToMostB2BInList");
2784 std::string listname = arguments[0];
2786 auto func = [listname](
const Particle * particle) ->
double {
2788 StoreObjPtr<ParticleList> list(listname);
2789 if (not list.isValid())
2790 B2FATAL(
"Invalid particle list name " << listname <<
" given to angleToMostB2BInList");
2793 if (list->getListSize() == 0)
2798 const auto p_this =
B2Vector3D(frame.getMomentum(particle).Vect());
2802 double maxAngle = 0;
2803 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2805 const Particle* compareme = list->getParticle(i);
2806 const auto p_compare =
B2Vector3D(frame.getMomentum(compareme).Vect());
2807 double angle = p_compare.Angle(p_this);
2808 if (maxAngle < angle) maxAngle = angle;
2818 if (arguments.size() != 2)
2819 B2FATAL(
"Wrong number of arguments for meta function mostB2BInList");
2821 std::string listname = arguments[0];
2826 auto func = [listname, var](
const Particle * particle) ->
double {
2828 StoreObjPtr<ParticleList> list(listname);
2829 if (not list.isValid())
2830 B2FATAL(
"Invalid particle list name " << listname <<
" given to mostB2BInList");
2834 const auto p_this =
B2Vector3D(frame.getMomentum(particle).Vect());
2838 double maxAngle = -1.0;
2840 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2842 const Particle* compareme = list->getParticle(i);
2843 const auto p_compare =
B2Vector3D(frame.getMomentum(compareme).Vect());
2844 double angle = p_compare.Angle(p_this);
2845 if (maxAngle < angle) {
2853 auto var_result = var->function(list->getParticle(iMostB2B));
2854 if (std::holds_alternative<double>(var_result))
2856 return std::get<double>(var_result);
2857 }
else if (std::holds_alternative<int>(var_result))
2859 return std::get<int>(var_result);
2860 }
else if (std::holds_alternative<bool>(var_result))
2862 return std::get<bool>(var_result);
2870 if (arguments.size() == 1) {
2871 std::string listName = arguments[0];
2872 auto func = [listName](
const Particle*) ->
double {
2873 StoreObjPtr<ParticleList> listOfParticles(listName);
2875 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to maxOpeningAngleInList");
2876 int nParticles = listOfParticles->getListSize();
2881 double maxOpeningAngle = -1;
2882 for (
int i = 0; i < nParticles; i++)
2884 B2Vector3D v1 = frame.getMomentum(listOfParticles->getParticle(i)).Vect();
2885 for (
int j = i + 1; j < nParticles; j++) {
2886 B2Vector3D v2 = frame.getMomentum(listOfParticles->getParticle(j)).Vect();
2887 const double angle =
v1.Angle(v2);
2888 if (angle > maxOpeningAngle) maxOpeningAngle = angle;
2891 return maxOpeningAngle;
2895 B2FATAL(
"Wrong number of arguments for meta function maxOpeningAngleInList");
2902 if (arguments.size() >= 2) {
2907 auto func = [var, arguments](
const Particle * particle) ->
double {
2908 if (particle ==
nullptr)
2910 B2WARNING(
"Trying to access a daughter that does not exist. Skipping");
2916 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
2920 for (
unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++)
2922 auto generalizedIndex = arguments[iCoord];
2923 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
2925 pSum += frame.getMomentum(dauPart);
2927 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
2933 Particle sumOfDaughters(pSum, 100);
2935 auto var_result = var->function(&sumOfDaughters);
2937 if (std::holds_alternative<double>(var_result))
2939 return std::get<double>(var_result);
2940 }
else if (std::holds_alternative<int>(var_result))
2942 return std::get<int>(var_result);
2943 }
else if (std::holds_alternative<bool>(var_result))
2945 return std::get<bool>(var_result);
2950 B2FATAL(
"Wrong number of arguments for meta function daughterCombination");
2953 Manager::FunctionPtr useAlternativeDaughterHypothesis(
const std::vector<std::string>& arguments)
2966 if (arguments.size() >= 2) {
2977 std::unordered_map<unsigned int, int> mapOfReplacedDaughters;
2980 for (
unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++) {
2981 auto replacedDauString = arguments[iCoord];
2983 std::vector<std::string> indexAndMass;
2984 boost::split(indexAndMass, replacedDauString, boost::is_any_of(
":"));
2987 if (indexAndMass.size() > 2) {
2988 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 "
2989 << replacedDauString <<
", while a correct syntax looks like 0:K+.");
2993 if (indexAndMass.size() < 2) {
2994 B2WARNING(
"The string indicating which daughter's mass should be replaced contains only one colon-separated element instead of two. The offending string is "
2995 << replacedDauString <<
", while a correct syntax looks like 0:K+.");
3002 dauIndex = Belle2::convertString<int>(indexAndMass[0]);
3003 }
catch (std::invalid_argument&) {
3004 B2FATAL(
"Found the string " << indexAndMass[0] <<
"instead of a daughter index.");
3008 TParticlePDG* particlePDG = TDatabasePDG::Instance()->GetParticle(indexAndMass[1].c_str());
3010 B2WARNING(
"Particle not in evt.pdl file! " << indexAndMass[1]);
3015 int pdgCode = particlePDG->PdgCode();
3016 mapOfReplacedDaughters[dauIndex] = pdgCode;
3020 if (mapOfReplacedDaughters.size() != arguments.size() - 1)
3021 B2FATAL(
"Overlapped daughter's index is detected in the meta-variable useAlternativeDaughterHypothesis");
3029 auto func = [var, mapOfReplacedDaughters](
const Particle * particle) ->
double {
3030 if (particle ==
nullptr)
3032 B2WARNING(
"Trying to access a particle that does not exist. Skipping");
3042 ROOT::Math::PxPyPzMVector pSum(0, 0, 0, 0);
3044 for (
unsigned int iDau = 0; iDau < particle->getNDaughters(); iDau++)
3046 const Particle* dauPart = particle->getDaughter(iDau);
3048 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << iDau);
3052 ROOT::Math::PxPyPzMVector dauMom = ROOT::Math::PxPyPzMVector(frame.getMomentum(dauPart));
3056 pdgCode = mapOfReplacedDaughters.at(iDau);
3057 }
catch (std::out_of_range&) {
3064 double p_x = dauMom.Px();
3065 double p_y = dauMom.Py();
3066 double p_z = dauMom.Pz();
3067 dauMom.SetCoordinates(p_x, p_y, p_z, TDatabasePDG::Instance()->GetParticle(pdgCode)->Mass());
3068 const_cast<Particle*
>(dummy->getDaughter(iDau))->set4VectorDividingByMomentumScaling(ROOT::Math::PxPyPzEVector(dauMom));
3071 const int charge = dummy->getDaughter(iDau)->getCharge();
3072 if (TDatabasePDG::Instance()->GetParticle(pdgCode)->Charge() / 3.0 ==
charge)
3073 const_cast<Particle*
>(dummy->getDaughter(iDau))->setPDGCode(pdgCode);
3075 const_cast<Particle*
>(dummy->getDaughter(iDau))->setPDGCode(-1 * pdgCode);
3081 dummy->set4Vector(ROOT::Math::PxPyPzEVector(pSum));
3083 auto var_result = var->function(dummy);
3086 if (std::holds_alternative<double>(var_result))
3088 return std::get<double>(var_result);
3089 }
else if (std::holds_alternative<int>(var_result))
3091 return std::get<int>(var_result);
3092 }
else if (std::holds_alternative<bool>(var_result))
3094 return std::get<bool>(var_result);
3100 B2FATAL(
"Wrong number of arguments for meta function useAlternativeDaughterHypothesis");
3105 if (arguments.size() == 2) {
3107 std::string arg = arguments[0];
3109 TParticlePDG* part = TDatabasePDG::Instance()->GetParticle(arg.c_str());
3111 if (part !=
nullptr) {
3112 pdg_code = std::abs(part->PdgCode());
3115 pdg_code = Belle2::convertString<int>(arg);
3116 }
catch (std::exception& e) {}
3119 if (pdg_code == -1) {
3120 B2FATAL(
"Ancestor " + arg +
" is not recognised. Please provide valid PDG code or particle name.");
3123 auto func = [pdg_code, var](
const Particle * particle) ->
double {
3124 const Particle* p = particle;
3126 int ancestor_level = std::get<double>(
Manager::Instance().getVariable(
"hasAncestor(" + std::to_string(pdg_code) +
", 0)")->
function(p));
3127 if ((ancestor_level <= 0) or (std::isnan(ancestor_level)))
3132 const MCParticle* i_p = p->getMCParticle();
3134 for (
int a = 0; a < ancestor_level ; a = a + 1)
3136 i_p = i_p->getMother();
3139 auto var_result = var->function(&m_p);
3140 if (std::holds_alternative<double>(var_result))
3142 return std::get<double>(var_result);
3143 }
else if (std::holds_alternative<int>(var_result))
3145 return std::get<int>(var_result);
3146 }
else if (std::holds_alternative<bool>(var_result))
3148 return std::get<bool>(var_result);
3153 B2FATAL(
"Wrong number of arguments for meta function varForFirstMCAncestorOfType (expected 2: type and variable of interest)");
3157 VARIABLE_GROUP(
"MetaFunctions");
3158 REGISTER_METAVARIABLE(
"nCleanedECLClusters(cut)", nCleanedECLClusters,
3159 "[Eventbased] Returns the number of clean Clusters in the event\n"
3160 "Clean clusters are defined by the clusters which pass the given cut assuming a photon hypothesis.",
3161 Manager::VariableDataType::c_int);
3162 REGISTER_METAVARIABLE(
"nCleanedTracks(cut)", nCleanedTracks,
3163 "[Eventbased] Returns the number of clean Tracks in the event\n"
3164 "Clean tracks are defined by the tracks which pass the given cut assuming a pion hypothesis.", Manager::VariableDataType::c_int);
3165 REGISTER_METAVARIABLE(
"formula(v1 + v2 * [v3 - v4] / v5^v6)", formula, R
"DOCSTRING(
3166 Returns the result of the given formula, where v1 to vN are variables or floating
3167 point numbers. Currently the only supported operations are addition (``+``),
3168 subtraction (``-``), multiplication (``*``), division (``/``) and power (``^``
3169 or ``**``). Parenthesis can be in the form of square brackets ``[v1 * v2]``
3170 or normal brackets ``(v1 * v2)``. It will work also with variables taking
3171 arguments. Operator precedence is taken into account. For example ::
3173 (daughter(0, E) + daughter(1, E))**2 - p**2 + 0.138
3175 .. versionchanged:: release-03-00-00
3176 now both, ``[]`` and ``()`` can be used for grouping operations, ``**`` can
3177 be used for exponent and float literals are possible directly in the
3179 )DOCSTRING", Manager::VariableDataType::c_double);
3180 REGISTER_METAVARIABLE("useRestFrame(variable)", useRestFrame,
3181 "Returns the value of the variable using the rest frame of the given particle as current reference frame.\n"
3182 "E.g. ``useRestFrame(daughter(0, p))`` returns the total momentum of the first daughter in its mother's rest-frame", Manager::VariableDataType::c_double);
3183 REGISTER_METAVARIABLE(
"useCMSFrame(variable)", useCMSFrame,
3184 "Returns the value of the variable using the CMS frame as current reference frame.\n"
3185 "E.g. ``useCMSFrame(E)`` returns the energy of a particle in the CMS frame.", Manager::VariableDataType::c_double);
3186 REGISTER_METAVARIABLE(
"useLabFrame(variable)", useLabFrame, R
"DOC(
3187 Returns the value of ``variable`` in the *lab* frame.
3190 The lab frame is the default reference frame, usually you don't need to use this meta-variable.
3191 E.g. ``useLabFrame(E)`` returns the energy of a particle in the Lab frame, same as just ``E``.
3193 Specifying the lab frame is useful in some corner-cases. For example:
3194 ``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.
3195 )DOC", Manager::VariableDataType::c_double);
3196 REGISTER_METAVARIABLE("useTagSideRecoilRestFrame(variable, daughterIndexTagB)", useTagSideRecoilRestFrame,
3197 "Returns the value of the variable in the rest frame of the recoiling particle to the tag side B meson.\n"
3198 "The variable should only be applied to an Upsilon(4S) list.\n"
3199 "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);
3200 REGISTER_METAVARIABLE(
"useParticleRestFrame(variable, particleList)", useParticleRestFrame,
3201 "Returns the value of the variable in the rest frame of the first Particle contained in the given ParticleList.\n"
3202 "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3203 "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3204 "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);
3205 REGISTER_METAVARIABLE(
"useRecoilParticleRestFrame(variable, particleList)", useRecoilParticleRestFrame,
3206 "Returns the value of the variable in the rest frame of recoil system against the first Particle contained in the given ParticleList.\n"
3207 "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3208 "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3209 "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);
3210 REGISTER_METAVARIABLE(
"useDaughterRestFrame(variable, daughterIndex_1, [daughterIndex_2, ... daughterIndex_3])", useDaughterRestFrame,
3211 "Returns the value of the variable in the rest frame of the selected daughter particle.\n"
3212 "The daughter is identified via generalized daughter index, e.g. ``0:1`` identifies the second daughter (1) "
3213 "of the first daughter (0). If the daughter index is invalid, it returns NaN.\n"
3214 "If two or more indices are given, the rest frame of the sum of the daughters is used.",
3215 Manager::VariableDataType::c_double);
3216 REGISTER_METAVARIABLE(
"passesCut(cut)", passesCut,
3217 "Returns 1 if particle passes the cut otherwise 0.\n"
3218 "Useful if you want to write out if a particle would have passed a cut or not.", Manager::VariableDataType::c_bool);
3219 REGISTER_METAVARIABLE(
"passesEventCut(cut)", passesEventCut,
3220 "[Eventbased] Returns 1 if event passes the cut otherwise 0.\n"
3221 "Useful if you want to select events passing a cut without looping into particles, such as for skimming.\n", Manager::VariableDataType::c_bool);
3222 REGISTER_METAVARIABLE(
"countDaughters(cut)", countDaughters,
3223 "Returns number of direct daughters which satisfy the cut.\n"
3224 "Used by the skimming package (for what exactly?)", Manager::VariableDataType::c_int);
3225 REGISTER_METAVARIABLE(
"countFSPDaughters(cut)", countDescendants,
3226 "Returns number of final-state daughters which satisfy the cut.",
3227 Manager::VariableDataType::c_int);
3228 REGISTER_METAVARIABLE(
"countDescendants(cut)", countDescendants,
3229 "Returns number of descendants for all generations which satisfy the cut.",
3230 Manager::VariableDataType::c_int);
3231 REGISTER_METAVARIABLE(
"varFor(pdgCode, variable)", varFor,
3232 "Returns the value of the variable for the given particle if its abs(pdgCode) agrees with the given one.\n"
3233 "E.g. ``varFor(11, p)`` returns the momentum if the particle is an electron or a positron.", Manager::VariableDataType::c_double);
3234 REGISTER_METAVARIABLE(
"varForMCGen(variable)", varForMCGen,
3235 "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"
3236 "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"
3237 "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);
3238 REGISTER_METAVARIABLE(
"nParticlesInList(particleListName)", nParticlesInList,
3239 "[Eventbased] Returns number of particles in the given particle List.", Manager::VariableDataType::c_int);
3240 REGISTER_METAVARIABLE(
"isInList(particleListName)", isInList,
3241 "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);
3242 REGISTER_METAVARIABLE(
"isDaughterOfList(particleListNames)", isDaughterOfList,
3243 "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);
3244 REGISTER_METAVARIABLE(
"isDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isDescendantOfList, R
"DOC(
3245 Returns 1 if the given particle appears in the decay chain of the particles in the given ParticleLists.
3247 Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3249 * ``isDescendantOfList(<particle_list>,1)`` returns 1 if particle is a daughter of the list,
3250 * ``isDescendantOfList(<particle_list>,2)`` returns 1 if particle is a granddaughter of the list,
3251 * ``isDescendantOfList(<particle_list>,3)`` returns 1 if particle is a great-granddaughter of the list, etc.
3252 * Default value is ``-1`` that is inclusive for all generations.
3253 )DOC", Manager::VariableDataType::c_bool);
3254 REGISTER_METAVARIABLE("isMCDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isMCDescendantOfList, R
"DOC(
3255 Returns 1 if the given particle is linked to the same MC particle as any reconstructed daughter of the decay lists.
3257 Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3259 * ``isMCDescendantOfList(<particle_list>,1)`` returns 1 if particle is matched to the same particle as any daughter of the list,
3260 * ``isMCDescendantOfList(<particle_list>,2)`` returns 1 if particle is matched to the same particle as any granddaughter of the list,
3261 * ``isMCDescendantOfList(<particle_list>,3)`` returns 1 if particle is matched to the same particle as any great-granddaughter of the list, etc.
3262 * Default value is ``-1`` that is inclusive for all generations.
3264 It makes only sense for lists created with `fillParticleListFromMC` function with ``addDaughters=True`` argument.
3265 )DOC", Manager::VariableDataType::c_bool);
3267 REGISTER_METAVARIABLE("sourceObjectIsInList(particleListName)", sourceObjectIsInList, R
"DOC(
3268 Returns 1 if the underlying mdst object (e.g. track, or cluster) was used to create a particle in ``particleListName``, 0 if not.
3271 This only makes sense for particles that are not composite. Returns -1 for composite particles.
3272 )DOC", Manager::VariableDataType::c_int);
3274 REGISTER_METAVARIABLE("mcParticleIsInMCList(particleListName)", mcParticleIsInMCList, R
"DOC(
3275 Returns 1 if the particle's matched MC particle is also matched to a particle in ``particleListName``
3276 (or if either of the lists were filled from generator level `modularAnalysis.fillParticleListFromMC`.)
3278 .. seealso:: :b2:var:`isMCDescendantOfList` to check daughters.
3279 )DOC", Manager::VariableDataType::c_bool);
3281 REGISTER_METAVARIABLE("isGrandDaughterOfList(particleListNames)", isGrandDaughterOfList,
3282 "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);
3283 REGISTER_METAVARIABLE(
"originalParticle(variable)", originalParticle, R
"DOC(
3284 Returns value of variable for the original particle from which the given particle is copied.
3286 The copy of particle is created, for example, when the vertex fit updates the daughters and `modularAnalysis.copyParticles` is called.
3287 Returns NaN if the given particle is not copied and so there is no original particle.
3288 )DOC", Manager::VariableDataType::c_double);
3289 REGISTER_METAVARIABLE("daughter(i, variable)", daughter, R
"DOC(
3290 Returns value of variable for the i-th daughter. E.g.
3292 * ``daughter(0, p)`` returns the total momentum of the first daughter.
3293 * ``daughter(0, daughter(1, p)`` returns the total momentum of the second daughter of the first daughter.
3295 Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3296 )DOC", Manager::VariableDataType::c_double);
3297 REGISTER_METAVARIABLE("originalDaughter(i, variable)", originalDaughter, R
"DOC(
3298 Returns value of variable for the original particle from which the i-th daughter is copied.
3300 The copy of particle is created, for example, when the vertex fit updates the daughters and `modularAnalysis.copyParticles` is called.
3301 Returns NaN if the daughter is not copied and so there is no original daughter.
3303 Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3304 )DOC", Manager::VariableDataType::c_double);
3305 REGISTER_METAVARIABLE("mcDaughter(i, variable)", mcDaughter, R
"DOC(
3306 Returns the value of the requested variable for the i-th Monte Carlo daughter of the particle.
3308 Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3309 or if the i-th MC daughter does not exist.
3311 E.g. ``mcDaughter(0, PDG)`` will return the PDG code of the first MC daughter of the matched MC
3312 particle of the reconstructed particle the function is applied to.
3314 The meta variable can also be nested: ``mcDaughter(0, mcDaughter(1, PDG))``.
3315 )DOC", Manager::VariableDataType::c_double);
3316 REGISTER_METAVARIABLE("mcMother(variable)", mcMother, R
"DOC(
3317 Returns the value of the requested variable for the Monte Carlo mother of the particle.
3319 Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3320 or if the MC mother does not exist.
3322 E.g. ``mcMother(PDG)`` will return the PDG code of the MC mother of the matched MC
3323 particle of the reconstructed particle the function is applied to.
3325 The meta variable can also be nested: ``mcMother(mcMother(PDG))``.
3326 )DOC", Manager::VariableDataType::c_double);
3327 REGISTER_METAVARIABLE("genParticle(index, variable)", genParticle, R
"DOC(
3328 [Eventbased] Returns the ``variable`` for the ith generator particle.
3329 The arguments of the function must be the ``index`` of the particle in the MCParticle Array,
3330 and ``variable``, the name of the function or variable for that generator particle.
3331 If ``index`` goes beyond the length of the MCParticles array, NaN will be returned.
3333 E.g. ``genParticle(0, p)`` returns the total momentum of the first MCParticle, which in a generic decay up to MC15 is
3334 the Upsilon(4S) and for MC16 and beyond the initial electron.
3335 )DOC", Manager::VariableDataType::c_double);
3336 REGISTER_METAVARIABLE("genUpsilon4S(variable)", genUpsilon4S, R
"DOC(
3337 [Eventbased] Returns the ``variable`` evaluated for the generator-level :math:`\Upsilon(4S)`.
3338 If no generator level :math:`\Upsilon(4S)` exists for the event, NaN will be returned.
3340 E.g. ``genUpsilon4S(p)`` returns the total momentum of the :math:`\Upsilon(4S)` in a generic decay.
3341 ``genUpsilon4S(mcDaughter(1, p))`` returns the total momentum of the second daughter of the
3342 generator-level :math:`\Upsilon(4S)` (i.e. the momentum of the second B meson in a generic decay).
3343 )DOC", Manager::VariableDataType::c_double);
3344 REGISTER_METAVARIABLE("daughterProductOf(variable)", daughterProductOf,
3345 "Returns product of a variable over all daughters.\n"
3346 "E.g. ``daughterProductOf(extraInfo(SignalProbability))`` returns the product of the SignalProbabilitys of all daughters.", Manager::VariableDataType::c_double);
3347 REGISTER_METAVARIABLE(
"daughterSumOf(variable)", daughterSumOf,
3348 "Returns sum of a variable over all daughters.\n"
3349 "E.g. ``daughterSumOf(nDaughters)`` returns the number of grand-daughters.", Manager::VariableDataType::c_double);
3350 REGISTER_METAVARIABLE(
"daughterLowest(variable)", daughterLowest,
3351 "Returns the lowest value of the given variable among all daughters.\n"
3352 "E.g. ``useCMSFrame(daughterLowest(p))`` returns the lowest momentum in CMS frame.", Manager::VariableDataType::c_double);
3353 REGISTER_METAVARIABLE(
"daughterHighest(variable)", daughterHighest,
3354 "Returns the highest value of the given variable among all daughters.\n"
3355 "E.g. ``useCMSFrame(daughterHighest(p))`` returns the highest momentum in CMS frame.", Manager::VariableDataType::c_double);
3356 REGISTER_METAVARIABLE(
"daughterDiffOf(i, j, variable)", daughterDiffOf,
3357 "Returns the difference of a variable between the two given daughters.\n"
3358 "E.g. ``useRestFrame(daughterDiffOf(0, 1, p))`` returns the momentum difference between first and second daughter in the rest frame of the given particle.\n"
3359 "(That means that it returns :math:`p_j - p_i`)\n"
3360 "Nota Bene: for the particular case 'variable=phi' you should use the :b2:var:`daughterDiffOfPhi` function.", Manager::VariableDataType::c_double);
3361 REGISTER_METAVARIABLE(
"mcDaughterDiffOf(i, j, variable)", mcDaughterDiffOf,
3362 "MC matched version of the `daughterDiffOf` function.", Manager::VariableDataType::c_double);
3363 REGISTER_METAVARIABLE(
"grandDaughterDiffOf(i, j, variable)", grandDaughterDiffOf,
3364 "Returns the difference of a variable between the first daughters of the two given daughters.\n"
3365 "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"
3366 "(That means that it returns :math:`p_j - p_i`)\n"
3367 "Nota Bene: for the particular case 'variable=phi' you should use the :b2:var:`grandDaughterDiffOfPhi` function.", Manager::VariableDataType::c_double);
3368 REGISTER_METAVARIABLE(
"daughterDiffOfPhi(i, j)", daughterDiffOfPhi,
3369 "Returns the difference in :math:`\\phi` between the two given daughters. The unit of the angle is ``rad``.\n"
3370 "The difference is signed and takes account of the ordering of the given daughters.\n"
3371 "The function returns :math:`\\phi_j - \\phi_i`.", Manager::VariableDataType::c_double);
3372 MAKE_DEPRECATED(
"daughterDiffOfPhi",
false,
"release-06-00-00", R
"DOC(
3373 The difference of the azimuthal angle :math:`\\phi` of two daughters can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3374 REGISTER_METAVARIABLE("mcDaughterDiffOfPhi(i, j)", mcDaughterDiffOfPhi,
3375 "MC matched version of the `daughterDiffOfPhi` function. The unit of the angle is ``rad``", Manager::VariableDataType::c_double);
3376 MAKE_DEPRECATED(
"mcDaughterDiffOfPhi",
false,
"release-06-00-00", R
"DOC(
3377 The difference of the azimuthal angle :math:`\\phi` of the MC partners of two daughters can be calculated with the generic variable :b2:var:`mcDaughterDiffOf`.)DOC");
3378 REGISTER_METAVARIABLE("grandDaughterDiffOfPhi(i, j)", grandDaughterDiffOfPhi,
3379 "Returns the difference in :math:`\\phi` between the first daughters of the two given daughters. The unit of the angle is ``rad``.\n"
3380 "The difference is signed and takes account of the ordering of the given daughters.\n"
3381 "The function returns :math:`\\phi_j - \\phi_i`.\n", Manager::VariableDataType::c_double);
3382 MAKE_DEPRECATED(
"grandDaughterDiffOfPhi",
false,
"release-06-00-00", R
"DOC(
3383 The difference of the azimuthal angle :math:`\\phi` of two granddaughters can be calculated with the generic variable :b2:var:`grandDaughterDiffOf`.)DOC");
3384 REGISTER_METAVARIABLE("daughterDiffOfClusterPhi(i, j)", daughterDiffOfClusterPhi,
3385 "Returns the difference in :math:`\\phi` between the ECLClusters of two given daughters. The unit of the angle is ``rad``.\n"
3386 "The difference is signed and takes account of the ordering of the given daughters.\n"
3387 "The function returns :math:`\\phi_j - \\phi_i`.\n"
3388 "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.", Manager::VariableDataType::c_double);
3389 MAKE_DEPRECATED(
"daughterDiffOfClusterPhi",
false,
"release-06-00-00", R
"DOC(
3390 The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two daughters can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3391 REGISTER_METAVARIABLE("grandDaughterDiffOfClusterPhi(i, j)", grandDaughterDiffOfClusterPhi,
3392 "Returns the difference in :math:`\\phi` between the ECLClusters of the daughters of the two given daughters. The unit of the angle is ``rad``.\n"
3393 "The difference is signed and takes account of the ordering of the given daughters.\n"
3394 "The function returns :math:`\\phi_j - \\phi_i`.\n"
3395 "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.\n", Manager::VariableDataType::c_double);
3396 MAKE_DEPRECATED(
"grandDaughterDiffOfClusterPhi",
false,
"release-06-00-00", R
"DOC(
3397 The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two granddaughters can be calculated with the generic variable :b2:var:`grandDaughterDiffOf`.)DOC");
3398 REGISTER_METAVARIABLE("daughterDiffOfPhiCMS(i, j)", daughterDiffOfPhiCMS,
3399 "Returns the difference in :math:`\\phi` between the two given daughters in the CMS frame. The unit of the angle is ``rad``.\n"
3400 "The difference is signed and takes account of the ordering of the given daughters.\n"
3401 "The function returns :math:`\\phi_j - \\phi_i`.", Manager::VariableDataType::c_double);
3402 MAKE_DEPRECATED(
"daughterDiffOfPhiCMS",
false,
"release-06-00-00", R
"DOC(
3403 The difference of the azimuthal angle :math:`\\phi` of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3404 REGISTER_METAVARIABLE("mcDaughterDiffOfPhiCMS(i, j)", daughterDiffOfPhiCMS,
3405 "MC matched version of the `daughterDiffOfPhiCMS` function. The unit of the angle is ``rad``", Manager::VariableDataType::c_double);
3406 MAKE_DEPRECATED(
"mcDaughterDiffOfPhiCMS",
false,
"release-06-00-00", R
"DOC(
3407 The difference of the azimuthal angle :math:`\\phi` of the MC partners of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`mcDaughterDiffOf`.)DOC");
3408 REGISTER_METAVARIABLE("daughterDiffOfClusterPhiCMS(i, j)", daughterDiffOfClusterPhiCMS,
3409 "Returns the difference in :math:`\\phi` between the ECLClusters of two given daughters in the CMS frame. The unit of the angle is ``rad``.\n"
3410 "The difference is signed and takes account of the ordering of the given daughters.\n"
3411 "The function returns :math:`\\phi_j - \\phi_i``.\n"
3412 "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.", Manager::VariableDataType::c_double);
3413 MAKE_DEPRECATED(
"daughterDiffOfClusterPhiCMS",
false,
"release-06-00-00", R
"DOC(
3414 The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3415 REGISTER_METAVARIABLE("daughterNormDiffOf(i, j, variable)", daughterNormDiffOf,
3416 "Returns the normalized difference of a variable between the two given daughters.\n"
3417 "E.g. ``daughterNormDiffOf(0, 1, p)`` returns the normalized momentum difference between first and second daughter in the lab frame.", Manager::VariableDataType::c_double);
3418 REGISTER_METAVARIABLE(
"daughterMotherDiffOf(i, variable)", daughterMotherDiffOf,
3419 "Returns the difference of a variable between the given daughter and the mother particle itself.\n"
3420 "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);
3421 REGISTER_METAVARIABLE(
"daughterMotherNormDiffOf(i, variable)", daughterMotherNormDiffOf,
3422 "Returns the normalized difference of a variable between the given daughter and the mother particle itself.\n"
3423 "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);
3424 REGISTER_METAVARIABLE(
"daughterAngle(daughterIndex_1, daughterIndex_2[, daughterIndex_3])", daughterAngle, R
"DOC(
3425 Returns the angle in between any pair of particles belonging to the same decay tree.
3426 The unit of the angle is ``rad``.
3428 The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3429 daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3430 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3431 identifies the second daughter of the root particle.
3433 Both two and three generalized indexes can be given to ``daughterAngle``. If two indices are given, the
3434 variable returns the angle between the momenta of the two given particles. If three indices are given, the
3435 variable returns the angle between the momentum of the third particle and a vector which is the sum of the
3436 first two daughter momenta.
3439 ``daughterAngle(0, 3)`` will return the angle between the first and fourth daughter.
3440 ``daughterAngle(0, 1, 3)`` will return the angle between the fourth daughter and the sum of the first and
3442 ``daughterAngle(0:0, 3:0)`` will return the angle between the first daughter of the first daughter, and
3443 the first daughter of the fourth daughter.
3445 )DOC", Manager::VariableDataType::c_double);
3446 REGISTER_METAVARIABLE("mcDaughterAngle(daughterIndex_1, daughterIndex_2, [daughterIndex_3])", mcDaughterAngle,
3447 "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);
3448 REGISTER_VARIABLE(
"grandDaughterDecayAngle(i, j)", grandDaughterDecayAngle,
3449 "Returns the decay angle of the granddaughter in the daughter particle's rest frame.\n"
3450 "It is calculated with respect to the reverted momentum vector of the particle.\n"
3451 "Two arguments representing the daughter and granddaughter indices have to be provided as arguments.\n\n",
"rad");
3452 REGISTER_VARIABLE(
"daughterClusterAngleInBetween(i, j)", daughterClusterAngleInBetween,
3453 "Returns the angle between clusters associated to the two daughters."
3454 "If two indices given: returns the angle between the momenta of the clusters associated to the two given daughters."
3455 "If three indices given: returns the angle between the momentum of the third particle's cluster and a vector "
3456 "which is the sum of the first two daughter's cluster momenta."
3457 "Returns nan if any of the daughters specified don't have an associated cluster."
3458 "The arguments in the argument vector must be integers corresponding to the ith and jth (and kth) daughters.\n\n",
"rad");
3459 REGISTER_METAVARIABLE(
"daughterInvM(i[, j, ...])", daughterInvM, R
"DOC(
3460 Returns the invariant mass adding the Lorentz vectors of the given daughters. The unit of the invariant mass is GeV/:math:`\text{c}^2`
3461 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.
3463 Daughters from different generations of the decay tree can be combined using generalized daughter indexes,
3464 which are simply colon-separated daughter indexes for each generation, starting from the root particle. For
3465 example, ``0:1:3`` identifies the fourth daughter (3) of the second daughter (1) of the first daughter(0) of
3466 the mother particle.
3468 Returns NaN if the given daughter-index is out of bound (>= number of daughters))DOC", Manager::VariableDataType::c_double);
3469 REGISTER_METAVARIABLE("extraInfo(name)", extraInfo,
3470 "Returns extra info stored under the given name.\n"
3471 "The extraInfo has to be set by a module first.\n"
3472 "E.g. ``extraInfo(SignalProbability)`` returns the SignalProbability calculated by the ``MVAExpert`` module.\n"
3473 "If nothing is set under the given name or if the particle is a nullptr, NaN is returned.\n"
3474 "In the latter case please use `eventExtraInfo` if you want to access an EventExtraInfo variable.", Manager::VariableDataType::c_double);
3475 REGISTER_METAVARIABLE(
"eventExtraInfo(name)", eventExtraInfo,
3476 "[Eventbased] Returns extra info stored under the given name in the event extra info.\n"
3477 "The extraInfo has to be set first by another module like MVAExpert in event mode.\n"
3478 "If nothing is set under this name, NaN is returned.", Manager::VariableDataType::c_double);
3479 REGISTER_METAVARIABLE(
"eventCached(variable)", eventCached,
3480 "[Eventbased] Returns value of event-based variable and caches this value in the EventExtraInfo.\n"
3481 "The result of second call to this variable in the same event will be provided from the cache.\n"
3482 "It is recommended to use this variable in order to declare custom aliases as event-based. This is "
3483 "necessary if using the eventwise mode of variablesToNtuple).", Manager::VariableDataType::c_double);
3484 REGISTER_METAVARIABLE(
"particleCached(variable)", particleCached,
3485 "Returns value of given variable and caches this value in the ParticleExtraInfo of the provided particle.\n"
3486 "The result of second call to this variable on the same particle will be provided from the cache.", Manager::VariableDataType::c_double);
3487 REGISTER_METAVARIABLE(
"modulo(variable, n)", modulo,
3488 "Returns rest of division of variable by n.", Manager::VariableDataType::c_int);
3489 REGISTER_METAVARIABLE(
"abs(variable)", abs,
3490 "Returns absolute value of the given variable.\n"
3491 "E.g. abs(mcPDG) returns the absolute value of the mcPDG, which is often useful for cuts.", Manager::VariableDataType::c_double);
3492 REGISTER_METAVARIABLE(
"max(var1,var2)", max,
"Returns max value of two variables.\n", Manager::VariableDataType::c_double);
3493 REGISTER_METAVARIABLE(
"min(var1,var2)", min,
"Returns min value of two variables.\n", Manager::VariableDataType::c_double);
3494 REGISTER_METAVARIABLE(
"sin(variable)", sin,
"Returns sine value of the given variable.", Manager::VariableDataType::c_double);
3495 REGISTER_METAVARIABLE(
"asin(variable)", asin,
"Returns arcsine of the given variable. The unit of the asin() is ``rad``", Manager::VariableDataType::c_double);
3496 REGISTER_METAVARIABLE(
"cos(variable)", cos,
"Returns cosine value of the given variable.", Manager::VariableDataType::c_double);
3497 REGISTER_METAVARIABLE(
"acos(variable)", acos,
"Returns arccosine value of the given variable. The unit of the acos() is ``rad``", Manager::VariableDataType::c_double);
3498 REGISTER_METAVARIABLE(
"tan(variable)",
tan,
"Returns tangent value of the given variable.", Manager::VariableDataType::c_double);
3499 REGISTER_METAVARIABLE(
"atan(variable)",
atan,
"Returns arctangent value of the given variable. The unit of the atan() is ``rad``", Manager::VariableDataType::c_double);
3500 REGISTER_METAVARIABLE(
"exp(variable)", exp,
"Returns exponential evaluated for the given variable.", Manager::VariableDataType::c_double);
3501 REGISTER_METAVARIABLE(
"log(variable)", log,
"Returns natural logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3502 REGISTER_METAVARIABLE(
"log10(variable)", log10,
"Returns base-10 logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3503 REGISTER_METAVARIABLE(
"isNAN(variable)", isNAN,
3504 "Returns true if variable value evaluates to nan (determined via std::isnan(double)).\n"
3505 "Useful for debugging.", Manager::VariableDataType::c_bool);
3506 REGISTER_METAVARIABLE(
"ifNANgiveX(variable, x)", ifNANgiveX,
3507 "Returns x (has to be a number) if variable value is nan (determined via std::isnan(double)).\n"
3508 "Useful for technical purposes while training MVAs.", Manager::VariableDataType::c_double);
3509 REGISTER_METAVARIABLE(
"isInfinity(variable)", isInfinity,
3510 "Returns true if variable value evaluates to infinity (determined via std::isinf(double)).\n"
3511 "Useful for debugging.", Manager::VariableDataType::c_bool);
3512 REGISTER_METAVARIABLE(
"unmask(variable, flag1, flag2, ...)", unmask,
3513 "unmask(variable, flag1, flag2, ...) or unmask(variable, mask) sets certain bits in the variable to zero.\n"
3514 "For example, if you want to set the second, fourth and fifth bits to zero, you could call \n"
3515 "``unmask(variable, 2, 8, 16)`` or ``unmask(variable, 26)``.\n"
3516 "", Manager::VariableDataType::c_double);
3517 REGISTER_METAVARIABLE(
"conditionalVariableSelector(cut, variableIfTrue, variableIfFalse)", conditionalVariableSelector,
3518 "Returns one of the two supplied variables, depending on whether the particle passes the supplied cut.\n"
3519 "The first variable is returned if the particle passes the cut, and the second variable is returned otherwise.", Manager::VariableDataType::c_double);
3520 REGISTER_METAVARIABLE(
"pValueCombination(p1, p2, ...)", pValueCombination,
3521 "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"
3522 "If any of the p-values is invalid, i.e. smaller than zero, -1 is returned.", Manager::VariableDataType::c_double);
3523 REGISTER_METAVARIABLE(
"veto(particleList, cut, pdgCode = 11)", veto,
3524 "Combines current particle with particles from the given particle list and returns 1 if the combination passes the provided cut. \n"
3525 "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"
3526 "around the neutral Pion mass (e.g. ``0.130 < M < 0.140``). \n"
3527 "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);
3528 REGISTER_METAVARIABLE(
"matchedMC(variable)", matchedMC,
3529 "Returns variable output for the matched MCParticle by constructing a temporary Particle from it.\n"
3530 "This may not work too well if your variable requires accessing daughters of the particle.\n"
3531 "E.g. ``matchedMC(p)`` returns the total momentum of the related MCParticle.\n"
3532 "Returns NaN if no matched MCParticle exists.", Manager::VariableDataType::c_double);
3533 REGISTER_METAVARIABLE(
"clusterBestMatchedMCParticle(variable)", clusterBestMatchedMCParticle,
3534 "Returns variable output for the MCParticle that is best-matched with the ECLCluster of the given Particle.\n"
3535 "E.g. To get the energy of the MCParticle that matches best with an ECLCluster, one could use ``clusterBestMatchedMCParticle(E)``\n"
3536 "When the variable is called for ``gamma`` and if the ``gamma`` is matched with MCParticle, it works same as `matchedMC`.\n"
3537 "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"
3538 "Returns NaN if the particle is not matched to an ECLCluster, or if the ECLCluster has no matching MCParticles", Manager::VariableDataType::c_double);
3539 REGISTER_METAVARIABLE(
"countInList(particleList, cut='')", countInList,
"[Eventbased] "
3540 "Returns number of particle which pass given in cut in the specified particle list.\n"
3541 "Useful for creating statistics about the number of particles in a list.\n"
3542 "E.g. ``countInList(e+, isSignal == 1)`` returns the number of correctly reconstructed electrons in the event.\n"
3543 "The variable is event-based and does not need a valid particle pointer as input.", Manager::VariableDataType::c_int);
3544 REGISTER_METAVARIABLE(
"getVariableByRank(particleList, rankedVariableName, variableName, rank)", getVariableByRank, R
"DOC(
3545 [Eventbased] Returns the value of ``variableName`` for the candidate in the ``particleList`` with the requested ``rank``.
3548 The `BestCandidateSelection` module available via `rankByHighest` / `rankByLowest` has to be used before.
3551 The first candidate matching the given rank is used.
3552 Thus, it is not recommended to use this variable in conjunction with ``allowMultiRank`` in the `BestCandidateSelection` module.
3554 The suffix ``_rank`` is automatically added to the argument ``rankedVariableName``,
3555 which either has to be the name of the variable used to order the candidates or the selected outputVariable name without the ending ``_rank``.
3556 This means that your selected name for the rank variable has to end with ``_rank``.
3558 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>`_
3559 )DOC", Manager::VariableDataType::c_double);
3560 REGISTER_VARIABLE("matchedMCHasPDG(PDGCode)", matchedMCHasPDG,
3561 "Returns if the absolute value of the PDGCode of the MCParticle related to the Particle matches a given PDGCode."
3562 "Returns 0/NAN/1 if PDGCode does not match/is not available/ matches");
3563 REGISTER_METAVARIABLE(
"numberOfNonOverlappingParticles(pList1, pList2, ...)", numberOfNonOverlappingParticles,
3564 "Returns the number of non-overlapping particles in the given particle lists"
3565 "Useful to check if there is additional physics going on in the detector if one reconstructed the Y4S", Manager::VariableDataType::c_int);
3566 REGISTER_METAVARIABLE(
"totalEnergyOfParticlesInList(particleListName)", totalEnergyOfParticlesInList,
3567 "[Eventbased] Returns the total energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3568 REGISTER_METAVARIABLE(
"totalPxOfParticlesInList(particleListName)", totalPxOfParticlesInList,
3569 "[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);
3570 REGISTER_METAVARIABLE(
"totalPyOfParticlesInList(particleListName)", totalPyOfParticlesInList,
3571 "[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);
3572 REGISTER_METAVARIABLE(
"totalPzOfParticlesInList(particleListName)", totalPzOfParticlesInList,
3573 "[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);
3574 REGISTER_METAVARIABLE(
"invMassInLists(pList1, pList2, ...)", invMassInLists,
3575 "[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);
3576 REGISTER_METAVARIABLE(
"totalECLEnergyOfParticlesInList(particleListName)", totalECLEnergyOfParticlesInList,
3577 "[Eventbased] Returns the total ECL energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3578 REGISTER_METAVARIABLE(
"maxPtInList(particleListName)", maxPtInList,
3579 "[Eventbased] Returns maximum transverse momentum Pt in the given particle List. The unit of the transverse momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3580 REGISTER_METAVARIABLE(
"eclClusterSpecialTrackMatched(cut)", eclClusterTrackMatchedWithCondition,
3581 "Returns if at least one Track that satisfies the given condition is related to the ECLCluster of the Particle.", Manager::VariableDataType::c_double);
3582 REGISTER_METAVARIABLE(
"averageValueInList(particleListName, variable)", averageValueInList,
3583 "[Eventbased] Returns the arithmetic mean of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3584 REGISTER_METAVARIABLE(
"medianValueInList(particleListName, variable)", medianValueInList,
3585 "[Eventbased] Returns the median value of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3586 REGISTER_METAVARIABLE(
"angleToClosestInList(particleListName)", angleToClosestInList,
3587 "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);
3588 REGISTER_METAVARIABLE(
"closestInList(particleListName, variable)", closestInList,
3589 "Returns `variable` for the closest particle (smallest opening angle) in the list provided.", Manager::VariableDataType::c_double);
3590 REGISTER_METAVARIABLE(
"angleToMostB2BInList(particleListName)", angleToMostB2BInList,
3591 "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);
3592 REGISTER_METAVARIABLE(
"mostB2BInList(particleListName, variable)", mostB2BInList,
3593 "Returns `variable` for the most back-to-back particle (closest opening angle to 180) in the list provided.", Manager::VariableDataType::c_double);
3594 REGISTER_METAVARIABLE(
"maxOpeningAngleInList(particleListName)", maxOpeningAngleInList,
3595 "[Eventbased] Returns maximum opening angle in the given particle List. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
3596 REGISTER_METAVARIABLE(
"daughterCombination(variable, daughterIndex_1, daughterIndex_2 ... daughterIndex_n)", daughterCombination,R
"DOC(
3597 Returns a ``variable`` function only of the 4-momentum calculated on an arbitrary set of (grand)daughters.
3600 ``variable`` can only be a function of the daughters' 4-momenta.
3602 Daughters from different generations of the decay tree can be combined using generalized daughter indexes, which are simply colon-separated
3603 the list of daughter indexes, starting from the root particle: for example, ``0:1:3`` identifies the fourth
3604 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle.
3607 ``daughterCombination(M, 0, 3, 4)`` will return the invariant mass of the system made of the first, fourth and fifth daughter of particle.
3608 ``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.
3610 )DOC", Manager::VariableDataType::c_double);
3611 REGISTER_METAVARIABLE("useAlternativeDaughterHypothesis(variable, daughterIndex_1:newMassHyp_1, ..., daughterIndex_n:newMassHyp_n)", useAlternativeDaughterHypothesis,R
"DOC(
3612 Returns a ``variable`` calculated using new mass hypotheses for (some of) the particle's daughters.
3615 ``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.
3616 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.
3617 Also, the track fit is not performed again: the variable only re-calculates the 4-vectors using different mass assumptions.
3618 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).
3621 Generalized daughter indexes are not supported (yet!): this variable can be used only on first-generation daughters.
3624 ``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.
3625 ``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.
3627 )DOC", Manager::VariableDataType::c_double);
3628 REGISTER_METAVARIABLE("varForFirstMCAncestorOfType(type, variable)",varForFirstMCAncestorOfType,R
"DOC(Returns requested variable of the first ancestor of the given type.
3629 Ancestor type can be set up by PDG code or by particle name (check evt.pdl for valid particle names))DOC", Manager::VariableDataType::c_double)
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
int getPDGCode() const
PDG code.
static const ChargedStable pion
charged pion particle
static const double doubleNaN
quiet_NaN
static const ChargedStable electron
electron particle
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
@ c_nPhotons
CR is split into n photons (N1)
static std::unique_ptr< GeneralCut > compile(const std::string &cut)
Creates an instance of a cut and returns a unique_ptr to it, if you need a copy-able object instead y...
@ c_Initial
bit 5: Particle is initial such as e+ or e- and not going to Geant4
@ c_PrimaryParticle
bit 0: Particle is primary particle.
@ c_IsVirtual
bit 4: Particle is virtual and not going to Geant4.
static std::string makeROOTCompatible(std::string str)
Remove special characters that ROOT dislikes in branch names, e.g.
const MCParticle * getMCParticle() const
Returns the pointer to the MCParticle object that was used to create this Particle (ParticleType == c...
EParticleSourceObject
particle source enumerators
const Particle * getParticleFromGeneralizedIndexString(const std::string &generalizedIndex) const
Explores the decay tree of the particle and returns the (grand^n)daughter identified by a generalized...
@ c_Unflavored
Is its own antiparticle or we don't know whether it is a particle/antiparticle.
@ c_Flavored
Is either particle or antiparticle.
static const ReferenceFrame & GetCurrent()
Get current rest frame.
std::function< VarVariant(const Particle *)> FunctionPtr
functions stored take a const Particle* and return VarVariant.
const Var * getVariable(std::string name)
Get the variable belonging to the given key.
static Manager & Instance()
get singleton instance.
Class to store variables with their name which were sent to the logging service.
#define MAKE_DEPRECATED(name, make_fatal, version, description)
Registers a variable as deprecated.
B2Vector3< double > B2Vector3D
typedef for common usage with double
double atan(double a)
atan for double
double tan(double a)
tan for double
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
bool hasAntiParticle(int pdgCode)
Checks if the particle with given pdg code has an anti-particle or not.
Particle * copyParticle(const Particle *original)
Function takes argument Particle and creates a copy of it and copies of all its (grand-)^n-daughters.
Abstract base class for different kinds of events.
const std::vector< double > v2
MATLAB generated random vector.
const std::vector< double > v1
MATLAB generated random vector.