10#include <analysis/variables/MetaVariables.h>
11#include <analysis/variables/MCTruthVariables.h>
13#include <analysis/VariableManager/Utility.h>
14#include <analysis/dataobjects/Particle.h>
15#include <analysis/dataobjects/ParticleList.h>
16#include <analysis/dataobjects/EventKinematics.h>
17#include <analysis/utility/PCmsLabTransform.h>
18#include <analysis/utility/ReferenceFrame.h>
19#include <analysis/utility/EvtPDLUtil.h>
20#include <analysis/utility/ParticleCopy.h>
21#include <analysis/utility/ValueIndexPairSorting.h>
22#include <analysis/ClusterUtility/ClusterUtils.h>
23#include <analysis/variables/VariableFormulaConstructor.h>
25#include <framework/logging/Logger.h>
26#include <framework/datastore/StoreArray.h>
27#include <framework/datastore/StoreObjPtr.h>
28#include <framework/dataobjects/EventExtraInfo.h>
29#include <framework/utilities/Conversion.h>
30#include <framework/utilities/MakeROOTCompatible.h>
31#include <framework/gearbox/Const.h>
33#include <mdst/dataobjects/Track.h>
34#include <mdst/dataobjects/MCParticle.h>
35#include <mdst/dataobjects/ECLCluster.h>
36#include <mdst/dataobjects/TrackFitResult.h>
38#include <boost/algorithm/string.hpp>
45#include <TDatabasePDG.h>
46#include <Math/Vector4D.h>
47#include <Math/VectorUtil.h>
55 double requireDoubleForFrameVariable(
const Variable::Manager::Var* var,
57 const std::string& frameFunction)
59 if (std::holds_alternative<double>(value)) {
60 return std::get<double>(value);
63 const char* returnedType = std::holds_alternative<int>(value) ?
"int" :
"bool";
64 B2ERROR(
"Meta function " << frameFunction <<
" expects a double variable, but '" << var->name
65 <<
"' returned " << returnedType <<
". Returning NaN.");
71 if (arguments.size() == 1) {
73 auto func = [var](
const Particle * particle) ->
double {
74 UseReferenceFrame<RestFrame> frame(particle);
75 return requireDoubleForFrameVariable(var, var->function(particle),
"useRestFrame");
79 B2FATAL(
"Wrong number of arguments for meta function useRestFrame");
85 if (arguments.size() == 1) {
87 auto func = [var](
const Particle * particle) ->
double {
88 UseReferenceFrame<CMSFrame> frame;
89 return requireDoubleForFrameVariable(var, var->function(particle),
"useCMSFrame");
93 B2FATAL(
"Wrong number of arguments for meta function useCMSFrame");
99 if (arguments.size() == 1) {
101 auto func = [var](
const Particle * particle) ->
double {
102 UseReferenceFrame<LabFrame> frame;
103 return requireDoubleForFrameVariable(var, var->function(particle),
"useLabFrame");
107 B2FATAL(
"Wrong number of arguments for meta function useLabFrame");
113 if (arguments.size() == 2) {
115 auto daughterFunction = convertToDaughterIndex({arguments[1]});
116 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
117 int daughterIndexTagB = std::get<int>(daughterFunction(particle));
118 if (daughterIndexTagB < 0)
121 if (particle->getPDGCode() != 300553)
123 B2ERROR(
"Variable should only be used on a Upsilon(4S) Particle List!");
128 ROOT::Math::PxPyPzEVector pSigB = T.getBeamFourMomentum() - particle->getDaughter(daughterIndexTagB)->get4Vector();
129 Particle tmp(pSigB, -particle->getDaughter(daughterIndexTagB)->getPDGCode());
131 UseReferenceFrame<RestFrame> frame(&tmp);
132 return requireDoubleForFrameVariable(var, var->function(particle),
"useTagSideRecoilRestFrame");
137 B2FATAL(
"Wrong number of arguments for meta function useTagSideRecoilRestFrame");
143 if (arguments.size() == 2) {
145 std::string listName = arguments[1];
146 auto func = [var, listName](
const Particle * particle) ->
double {
147 StoreObjPtr<ParticleList> list(listName);
148 unsigned listSize = list->getListSize();
152 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."
153 << LogVar(
"ParticleList", listName)
154 << LogVar(
"Number of candidates in the list", listSize));
155 const Particle* p = list->getParticle(0);
156 UseReferenceFrame<RestFrame> frame(p);
157 return requireDoubleForFrameVariable(var, var->function(particle),
"useParticleRestFrame");
161 B2FATAL(
"Wrong number of arguments for meta function useParticleRestFrame.");
167 if (arguments.size() == 2) {
169 std::string listName = arguments[1];
170 auto func = [var, listName](
const Particle * particle) ->
double {
171 StoreObjPtr<ParticleList> list(listName);
172 unsigned listSize = list->getListSize();
176 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."
177 << LogVar(
"ParticleList", listName)
178 << LogVar(
"Number of candidates in the list", listSize));
179 const Particle* p = list->getParticle(0);
181 ROOT::Math::PxPyPzEVector recoil = T.getBeamFourMomentum() - p->get4Vector();
183 Particle pRecoil(recoil, 0);
184 pRecoil.setVertex(particle->getVertex());
185 UseReferenceFrame<RestFrame> frame(&pRecoil);
186 return requireDoubleForFrameVariable(var, var->function(particle),
"useRecoilParticleRestFrame");
190 B2FATAL(
"Wrong number of arguments for meta function useParticleRestFrame.");
196 if (arguments.size() >= 2) {
198 auto func = [var, arguments](
const Particle * particle) ->
double {
201 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
203 for (
unsigned int i = 1; i < arguments.size(); i++)
205 auto generalizedIndex = arguments[i];
206 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
208 pSum += dauPart->get4Vector();
212 Particle tmp(pSum, 0);
213 UseReferenceFrame<RestFrame> frame(&tmp);
214 return requireDoubleForFrameVariable(var, var->function(particle),
"useDaughterRestFrame");
218 B2FATAL(
"Wrong number of arguments for meta function useDaughterRestFrame.");
224 if (arguments.size() >= 2) {
226 auto func = [var, arguments](
const Particle * particle) ->
double {
229 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
231 for (
unsigned int i = 1; i < arguments.size(); i++)
233 auto generalizedIndex = arguments[i];
234 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
236 pSum += dauPart->get4Vector();
241 ROOT::Math::PxPyPzEVector recoil = T.getBeamFourMomentum() - pSum;
243 Particle pRecoil(recoil, 0);
244 UseReferenceFrame<RestFrame> frame(&pRecoil);
245 return requireDoubleForFrameVariable(var, var->function(particle),
"useDaughterRecoilRestFrame");
249 B2FATAL(
"Wrong number of arguments for meta function useDaughterRecoilRestFrame.");
255 if (arguments.size() == 1) {
257 auto func = [var](
const Particle * particle) ->
double {
258 int index = ancestorBIndex(particle);
260 StoreArray<MCParticle> mcparticles;
261 Particle temp(mcparticles[index]);
262 UseReferenceFrame<RestFrame> frame(&temp);
263 return requireDoubleForFrameVariable(var, var->function(particle),
"useMCancestorBRestFrame");
267 B2FATAL(
"Wrong number of arguments for meta function useMCancestorBRestFrame.");
273 if (arguments.size() == 1) {
274 auto extraInfoName = arguments[0];
275 auto func = [extraInfoName](
const Particle * particle) ->
double {
276 if (particle ==
nullptr)
278 B2WARNING(
"Returns NaN because the particle is nullptr! If you want EventExtraInfo variables, please use eventExtraInfo() instead");
281 if (particle->hasExtraInfo(extraInfoName))
283 return particle->getExtraInfo(extraInfoName);
291 B2FATAL(
"Wrong number of arguments for meta function extraInfo");
297 if (arguments.size() == 1) {
298 auto extraInfoName = arguments[0];
299 auto func = [extraInfoName](
const Particle*) ->
double {
300 StoreObjPtr<EventExtraInfo> eventExtraInfo;
301 if (not eventExtraInfo.isValid())
303 if (eventExtraInfo->hasExtraInfo(extraInfoName))
305 return eventExtraInfo->getExtraInfo(extraInfoName);
313 B2FATAL(
"Wrong number of arguments for meta function extraInfo");
319 if (arguments.size() == 1) {
322 auto func = [var, key](
const Particle*) ->
double {
324 StoreObjPtr<EventExtraInfo> eventExtraInfo;
325 if (not eventExtraInfo.isValid())
326 eventExtraInfo.create();
327 if (eventExtraInfo->hasExtraInfo(key))
329 return eventExtraInfo->getExtraInfo(key);
333 auto var_result = var->function(
nullptr);
334 if (std::holds_alternative<double>(var_result)) {
335 value = std::get<double>(var_result);
336 }
else if (std::holds_alternative<int>(var_result)) {
337 return std::get<int>(var_result);
338 }
else if (std::holds_alternative<bool>(var_result)) {
339 return std::get<bool>(var_result);
341 eventExtraInfo->addExtraInfo(key, value);
347 B2FATAL(
"Wrong number of arguments for meta function eventCached");
353 if (arguments.size() == 1) {
356 auto func = [var, key](
const Particle * particle) ->
double {
358 if (particle->hasExtraInfo(key))
360 return particle->getExtraInfo(key);
363 double value = std::get<double>(var->function(particle));
371 const_cast<Particle*
>(particle)->addExtraInfo(key, value);
377 B2FATAL(
"Wrong number of arguments for meta function particleCached");
386 if (arguments.size() != 1) B2FATAL(
"Wrong number of arguments for meta function formula");
387 FormulaParser<VariableFormulaConstructor> parser;
389 return parser.parse(arguments[0]);
390 }
catch (std::runtime_error& e) {
397 if (arguments.size() <= 1) {
399 std::string cutString;
400 if (arguments.size() == 1)
401 cutString = arguments[0];
403 auto func = [cut](
const Particle*) ->
int {
405 int number_of_tracks = 0;
406 StoreArray<Track> tracks;
407 for (
const auto& track : tracks)
409 const TrackFitResult* trackFit = track.getTrackFitResultWithClosestMass(
Const::pion);
410 if (!trackFit)
continue;
411 if (trackFit->getChargeSign() == 0) {
415 if (cut->check(&particle))
420 return number_of_tracks;
425 B2FATAL(
"Wrong number of arguments for meta function nCleanedTracks");
431 if (arguments.size() <= 1) {
433 std::string cutString;
434 if (arguments.size() == 1)
435 cutString = arguments[0];
437 auto func = [cut](
const Particle*) ->
int {
439 int number_of_clusters = 0;
440 StoreArray<ECLCluster> clusters;
441 for (
const auto& cluster : clusters)
447 Particle particle(&cluster);
448 if (cut->check(&particle))
449 number_of_clusters++;
452 return number_of_clusters;
457 B2FATAL(
"Wrong number of arguments for meta function nCleanedECLClusters");
463 if (arguments.size() == 1) {
464 std::string cutString = arguments[0];
466 auto func = [cut](
const Particle * particle) ->
bool {
467 if (cut->check(particle))
474 B2FATAL(
"Wrong number of arguments for meta function passesCut");
480 if (arguments.size() == 1) {
481 std::string cutString = arguments[0];
483 auto func = [cut](
const Particle*) ->
bool {
484 if (cut->check(
nullptr))
491 B2FATAL(
"Wrong number of arguments for meta function passesEventCut");
497 if (arguments.size() == 2) {
501 }
catch (std::invalid_argument&) {
502 B2FATAL(
"The first argument of varFor meta function must be a positive integer!");
505 auto func = [pdgCode, var](
const Particle * particle) ->
double {
506 if (std::abs(particle->getPDGCode()) == std::abs(pdgCode))
508 auto var_result = var->function(particle);
509 if (std::holds_alternative<double>(var_result)) {
510 return std::get<double>(var_result);
511 }
else if (std::holds_alternative<int>(var_result)) {
512 return std::get<int>(var_result);
513 }
else if (std::holds_alternative<bool>(var_result)) {
514 return std::get<bool>(var_result);
520 B2FATAL(
"Wrong number of arguments for meta function varFor");
526 if (arguments.size() == 1) {
528 auto func = [var](
const Particle * particle) ->
double {
529 if (particle->getMCParticle())
534 auto var_result = var->function(particle);
535 if (std::holds_alternative<double>(var_result)) {
536 return std::get<double>(var_result);
537 }
else if (std::holds_alternative<int>(var_result)) {
538 return std::get<int>(var_result);
539 }
else if (std::holds_alternative<bool>(var_result)) {
540 return std::get<bool>(var_result);
547 B2FATAL(
"Wrong number of arguments for meta function varForMCGen");
553 if (arguments.size() == 1) {
554 std::string listName = arguments[0];
555 auto func = [listName](
const Particle * particle) ->
int {
558 StoreObjPtr<ParticleList> listOfParticles(listName);
560 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to nParticlesInList");
562 return listOfParticles->getListSize();
567 B2FATAL(
"Wrong number of arguments for meta function nParticlesInList");
574 if (arguments.size() != 1) {
575 B2FATAL(
"Wrong number of arguments for isInList");
577 auto listName = arguments[0];
579 auto func = [listName](
const Particle * particle) ->
bool {
582 StoreObjPtr<ParticleList> list(listName);
583 if (!(list.isValid()))
585 B2FATAL(
"Invalid Listname " << listName <<
" given to isInList");
589 return list->contains(particle);
598 if (arguments.size() != 1) {
599 B2FATAL(
"Wrong number of arguments for sourceObjectIsInList");
601 auto listName = arguments[0];
603 auto func = [listName](
const Particle * particle) ->
int {
606 StoreObjPtr<ParticleList> list(listName);
607 if (!(list.isValid()))
609 B2FATAL(
"Invalid Listname " << listName <<
" given to sourceObjectIsInList");
615 if (particlesource == Particle::EParticleSourceObject::c_Composite
616 or particlesource == Particle::EParticleSourceObject::c_Undefined)
622 for (
unsigned i = 0; i < list->getListSize(); ++i)
624 Particle* iparticle = list->getParticle(i);
625 if (particle->getMdstSource() == iparticle->getMdstSource())
637 if (arguments.size() != 1) {
638 B2FATAL(
"Wrong number of arguments for mcParticleIsInMCList");
640 auto listName = arguments[0];
642 auto func = [listName](
const Particle * particle) ->
bool {
645 StoreObjPtr<ParticleList> list(listName);
646 if (!(list.isValid()))
647 B2FATAL(
"Invalid Listname " << listName <<
" given to mcParticleIsInMCList");
650 const MCParticle* mcp = particle->getMCParticle();
651 if (mcp ==
nullptr)
return false;
654 for (
unsigned i = 0; i < list->getListSize(); ++i)
656 const MCParticle* imcp = list->getParticle(i)->getMCParticle();
657 if ((imcp !=
nullptr) and (mcp->getArrayIndex() == imcp->getArrayIndex()))
667 B2WARNING(
"isDaughterOfList is outdated and replaced by isDescendantOfList.");
668 std::vector<std::string> new_arguments = arguments;
669 new_arguments.push_back(std::string(
"1"));
670 return isDescendantOfList(new_arguments);
675 B2WARNING(
"isGrandDaughterOfList is outdated and replaced by isDescendantOfList.");
676 std::vector<std::string> new_arguments = arguments;
677 new_arguments.push_back(std::string(
"2"));
678 return isDescendantOfList(new_arguments);
683 if (arguments.size() > 0) {
684 auto listNames = arguments;
685 auto func = [listNames](
const Particle * particle) ->
bool {
687 int generation_flag = -1;
691 }
catch (std::exception& e) {}
693 for (
auto& iListName : listNames)
698 }
catch (std::exception& e) {}
701 auto list_comparison = [](
auto&& self,
const Particle * m,
const Particle * p,
int flag)->
bool {
703 for (
unsigned i = 0; i < m->getNDaughters(); ++i)
705 const Particle* daughter = m->getDaughter(i);
706 if ((flag == 1.) or (flag < 0)) {
707 if (p->isCopyOf(daughter)) {
713 if (daughter->getNDaughters() > 0) {
714 result = self(self, daughter, p, flag - 1);
724 StoreObjPtr<ParticleList> listOfParticles(iListName);
726 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << iListName <<
" given to isDescendantOfList");
728 for (
unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
729 Particle* iParticle = listOfParticles->getParticle(i);
730 output = list_comparison(list_comparison, iParticle, particle, generation_flag);
740 B2FATAL(
"Wrong number of arguments for meta function isDescendantOfList");
746 if (arguments.size() > 0) {
747 auto listNames = arguments;
748 auto func = [listNames](
const Particle * particle) ->
bool {
750 int generation_flag = -1;
754 }
catch (std::exception& e) {}
756 if (particle->getMCParticle() ==
nullptr)
761 for (
auto& iListName : listNames)
764 std::stod(iListName);
766 }
catch (std::exception& e) {}
768 auto list_comparison = [](
auto&& self,
const Particle * m,
const Particle * p,
int flag)->
bool {
770 for (
unsigned i = 0; i < m->getNDaughters(); ++i)
772 const Particle* daughter = m->getDaughter(i);
773 if ((flag == 1.) or (flag < 0)) {
774 if (daughter->getMCParticle() !=
nullptr) {
775 if (p->getMCParticle()->getArrayIndex() == daughter->getMCParticle()->getArrayIndex()) {
781 if (daughter->getNDaughters() > 0) {
782 result = self(self, daughter, p, flag - 1);
792 StoreObjPtr<ParticleList> listOfParticles(iListName);
794 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << iListName <<
" given to isMCDescendantOfList");
796 for (
unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
797 Particle* iParticle = listOfParticles->getParticle(i);
798 output = list_comparison(list_comparison, iParticle, particle, generation_flag);
808 B2FATAL(
"Wrong number of arguments for meta function isMCDescendantOfList");
814 if (arguments.size() == 1) {
816 auto func = [var](
const Particle * particle) ->
double {
817 double product = 1.0;
818 if (particle->getNDaughters() == 0)
822 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
824 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
825 product *= std::get<double>(var->function(particle->getDaughter(j)));
827 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
829 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
830 product *= std::get<int>(var->function(particle->getDaughter(j)));
837 B2FATAL(
"Wrong number of arguments for meta function daughterProductOf");
843 if (arguments.size() == 1) {
845 auto func = [var](
const Particle * particle) ->
double {
847 if (particle->getNDaughters() == 0)
851 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
853 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
854 sum += std::get<double>(var->function(particle->getDaughter(j)));
856 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
858 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
859 sum += std::get<int>(var->function(particle->getDaughter(j)));
866 B2FATAL(
"Wrong number of arguments for meta function daughterSumOf");
872 if (arguments.size() == 1) {
874 auto func = [var](
const Particle * particle) ->
double {
876 if (particle->getNDaughters() == 0)
880 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
882 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
883 double iValue = std::get<double>(var->function(particle->getDaughter(j)));
884 if (std::isnan(iValue))
continue;
885 if (std::isnan(min)) min = iValue;
886 if (iValue < min) min = iValue;
888 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
890 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
891 int iValue = std::get<int>(var->function(particle->getDaughter(j)));
892 if (std::isnan(min)) min = iValue;
893 if (iValue < min) min = iValue;
900 B2FATAL(
"Wrong number of arguments for meta function daughterLowest");
906 if (arguments.size() == 1) {
908 auto func = [var](
const Particle * particle) ->
double {
910 if (particle->getNDaughters() == 0)
914 if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
916 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
917 double iValue = std::get<double>(var->function(particle->getDaughter(j)));
918 if (std::isnan(iValue))
continue;
919 if (std::isnan(max)) max = iValue;
920 if (iValue > max) max = iValue;
922 }
else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
924 for (
unsigned j = 0; j < particle->getNDaughters(); ++j) {
925 int iValue = std::get<int>(var->function(particle->getDaughter(j)));
926 if (std::isnan(max)) max = iValue;
927 if (iValue > max) max = iValue;
934 B2FATAL(
"Wrong number of arguments for meta function daughterHighest");
940 if (arguments.size() == 3) {
941 auto func = [arguments](
const Particle * particle) ->
double {
942 if (particle ==
nullptr)
944 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
945 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
946 auto variablename = arguments[2];
947 if (dau_i ==
nullptr || dau_j ==
nullptr)
949 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
953 auto result_j = var->function(dau_j);
954 auto result_i = var->function(dau_i);
956 if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i))
958 diff = std::get<double>(result_j) - std::get<double>(result_i);
959 }
else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i))
961 diff = std::get<int>(result_j) - std::get<int>(result_i);
964 throw std::runtime_error(
"Bad variant access");
966 if (variablename ==
"phi" or variablename ==
"clusterPhi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)"))
967 or std::regex_match(variablename, std::regex(
"use.*Frame\\(clusterPhi\\)")))
969 if (fabs(diff) > M_PI) {
971 diff = diff - 2 * M_PI;
973 diff = 2 * M_PI + diff;
981 B2FATAL(
"Wrong number of arguments for meta function daughterDiffOf");
987 if (arguments.size() == 3) {
988 auto func = [arguments](
const Particle * particle) ->
double {
989 if (particle ==
nullptr)
991 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
992 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
993 auto variablename = arguments[2];
994 if (dau_i ==
nullptr || dau_j ==
nullptr)
996 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
999 const MCParticle* iMcDaughter = dau_i->getMCParticle();
1000 const MCParticle* jMcDaughter = dau_j->getMCParticle();
1001 if (iMcDaughter ==
nullptr || jMcDaughter ==
nullptr)
1003 Particle iTmpPart(iMcDaughter);
1004 Particle jTmpPart(jMcDaughter);
1006 auto result_j = var->function(&jTmpPart);
1007 auto result_i = var->function(&iTmpPart);
1009 if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i))
1011 diff = std::get<double>(result_j) - std::get<double>(result_i);
1012 }
else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i))
1014 diff = std::get<int>(result_j) - std::get<int>(result_i);
1017 throw std::runtime_error(
"Bad variant access");
1019 if (variablename ==
"phi" or std::regex_match(variablename, std::regex(
"use.*Frame\\(phi\\)")))
1021 if (fabs(diff) > M_PI) {
1023 diff = diff - 2 * M_PI;
1025 diff = 2 * M_PI + diff;
1033 B2FATAL(
"Wrong number of arguments for meta function mcDaughterDiffOf");
1039 if (arguments.size() == 5) {
1045 }
catch (std::invalid_argument&) {
1046 B2FATAL(
"First four arguments of grandDaughterDiffOf meta function must be integers!");
1048 std::vector<std::string> new_arguments;
1049 new_arguments.push_back(std::string(arguments[0] +
":" + arguments[2]));
1050 new_arguments.push_back(std::string(arguments[1] +
":" + arguments[3]));
1051 new_arguments.push_back(arguments[4]);
1052 return daughterDiffOf(new_arguments);
1054 B2FATAL(
"Wrong number of arguments for meta function grandDaughterDiffOf");
1060 if (arguments.size() == 3) {
1061 auto func = [arguments](
const Particle * particle) ->
double {
1062 if (particle ==
nullptr)
1064 const Particle* dau_i = particle->getParticleFromGeneralizedIndexString(arguments[0]);
1065 const Particle* dau_j = particle->getParticleFromGeneralizedIndexString(arguments[1]);
1066 if (!(dau_i && dau_j))
1068 B2ERROR(
"One of the first two arguments doesn't specify a valid (grand-)daughter!");
1072 double iValue, jValue;
1073 if (std::holds_alternative<double>(var->function(dau_j)))
1075 iValue = std::get<double>(var->function(dau_i));
1076 jValue = std::get<double>(var->function(dau_j));
1077 }
else if (std::holds_alternative<int>(var->function(dau_j)))
1079 iValue = std::get<int>(var->function(dau_i));
1080 jValue = std::get<int>(var->function(dau_j));
1082 return (jValue - iValue) / (jValue + iValue);
1086 B2FATAL(
"Wrong number of arguments for meta function daughterNormDiffOf");
1092 if (arguments.size() == 2) {
1093 auto daughterFunction = convertToDaughterIndex({arguments[0]});
1094 std::string variableName = arguments[1];
1095 auto func = [daughterFunction, variableName](
const Particle * particle) ->
double {
1096 if (particle ==
nullptr)
1098 int daughterNumber = std::get<int>(daughterFunction(particle));
1099 if (daughterNumber >=
int(particle->getNDaughters()) or daughterNumber < 0)
1102 auto result_mother = var->function(particle);
1103 auto result_daughter = var->function(particle->getDaughter(daughterNumber));
1105 if (std::holds_alternative<double>(result_mother) && std::holds_alternative<double>(result_daughter))
1107 diff = std::get<double>(result_mother) - std::get<double>(result_daughter);
1108 }
else if (std::holds_alternative<int>(result_mother) && std::holds_alternative<int>(result_daughter))
1110 diff = std::get<int>(result_mother) - std::get<int>(result_daughter);
1113 throw std::runtime_error(
"Bad variant access");
1116 if (variableName ==
"phi" or variableName ==
"useCMSFrame(phi)")
1118 if (fabs(diff) > M_PI) {
1120 diff = diff - 2 * M_PI;
1122 diff = 2 * M_PI + diff;
1130 B2FATAL(
"Wrong number of arguments for meta function daughterMotherDiffOf");
1136 if (arguments.size() == 2) {
1137 auto daughterFunction = convertToDaughterIndex({arguments[0]});
1139 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
1140 if (particle ==
nullptr)
1142 int daughterNumber = std::get<int>(daughterFunction(particle));
1143 if (daughterNumber >=
int(particle->getNDaughters()) or daughterNumber < 0)
1145 double daughterValue = 0.0, motherValue = 0.0;
1146 if (std::holds_alternative<double>(var->function(particle)))
1148 daughterValue = std::get<double>(var->function(particle->getDaughter(daughterNumber)));
1149 motherValue = std::get<double>(var->function(particle));
1150 }
else if (std::holds_alternative<int>(var->function(particle)))
1152 daughterValue = std::get<int>(var->function(particle->getDaughter(daughterNumber)));
1153 motherValue = std::get<int>(var->function(particle));
1155 return (motherValue - daughterValue) / (motherValue + daughterValue);
1159 B2FATAL(
"Wrong number of arguments for meta function daughterMotherNormDiffOf");
1165 if (arguments.size() >= 1) {
1167 auto func = [arguments](
const Particle * particle) ->
double {
1168 if (particle ==
nullptr)
1173 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
1174 for (
auto& generalizedIndex : arguments)
1176 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1177 if (dauPart) pSum += frame.getMomentum(dauPart);
1179 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1185 ROOT::Math::PxPyPzEVector pIN = T.getBeamFourMomentum();
1186 ROOT::Math::PxPyPzEVector pRecoil = frame.getMomentum(pIN - particle->get4Vector());
1188 return ROOT::Math::VectorUtil::Angle(pRecoil, pSum);
1192 B2FATAL(
"Wrong number of arguments for meta function angleBetweenDaughterAndRecoil");
1196 Manager::FunctionPtr angleBetweenDaughterAndMissingMomentum(
const std::vector<std::string>& arguments)
1198 if (arguments.size() >= 1) {
1199 auto func = [arguments](
const Particle * particle) ->
double {
1200 if (particle ==
nullptr)
1203 StoreObjPtr<EventKinematics> evtShape;
1206 B2WARNING(
"Cannot find missing momentum information, did you forget to run EventKinematicsModule?");
1209 ROOT::Math::XYZVector missingMomentumCMS = evtShape->getMissingMomentumCMS();
1210 ROOT::Math::PxPyPzEVector missingTotalMomentumCMS(missingMomentumCMS.X(),
1211 missingMomentumCMS.Y(),
1212 missingMomentumCMS.Z(),
1213 evtShape->getMissingEnergyCMS());
1215 ROOT::Math::PxPyPzEVector missingTotalMomentumLab = T.rotateCmsToLab() * missingTotalMomentumCMS;
1218 ROOT::Math::PxPyPzEVector pMiss = frame.getMomentum(missingTotalMomentumLab);
1220 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
1221 for (
auto& generalizedIndex : arguments)
1223 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1224 if (dauPart) pSum += frame.getMomentum(dauPart);
1226 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1231 return ROOT::Math::VectorUtil::Angle(pMiss, pSum);
1235 B2FATAL(
"Wrong number of arguments for meta function angleBetweenDaughterAndMissingMomentum");
1241 if (arguments.size() == 2 || arguments.size() == 3) {
1243 auto func = [arguments](
const Particle * particle) ->
double {
1244 if (particle ==
nullptr)
1247 std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1251 for (
auto& generalizedIndex : arguments)
1253 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1255 pDaus.push_back(frame.getMomentum(dauPart));
1257 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1263 if (pDaus.size() == 2)
1264 return ROOT::Math::VectorUtil::Angle(pDaus[0], pDaus[1]);
1266 return ROOT::Math::VectorUtil::Angle(pDaus[2], pDaus[0] + pDaus[1]);
1270 B2FATAL(
"Wrong number of arguments for meta function daughterAngle");
1274 double grandDaughterDecayAngle(
const Particle* particle,
const std::vector<double>& arguments)
1276 if (arguments.size() == 2) {
1281 int daughterIndex = std::lround(arguments[0]);
1282 if (daughterIndex >=
int(particle->getNDaughters()))
1284 const Particle* dau = particle->getDaughter(daughterIndex);
1286 int grandDaughterIndex = std::lround(arguments[1]);
1287 if (grandDaughterIndex >=
int(dau->getNDaughters()))
1290 ROOT::Math::XYZVector boost = dau->get4Vector().BoostToCM();
1292 ROOT::Math::PxPyPzEVector motherMomentum = - particle->get4Vector();
1293 motherMomentum = ROOT::Math::Boost(boost) * motherMomentum;
1295 ROOT::Math::PxPyPzEVector grandDaughterMomentum = dau->getDaughter(grandDaughterIndex)->get4Vector();
1296 grandDaughterMomentum = ROOT::Math::Boost(boost) * grandDaughterMomentum;
1298 return ROOT::Math::VectorUtil::Angle(motherMomentum, grandDaughterMomentum);
1301 B2FATAL(
"The variable grandDaughterDecayAngle needs exactly two integers as arguments!");
1307 if (arguments.size() == 2 || arguments.size() == 3) {
1309 auto func = [arguments](
const Particle * particle) ->
double {
1310 if (particle ==
nullptr)
1313 std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1317 if (particle->getParticleSource() == Particle::EParticleSourceObject::c_MCParticle)
1319 for (
auto& generalizedIndex : arguments) {
1320 const MCParticle* mcPart = particle->getMCParticle();
1321 if (mcPart ==
nullptr)
1323 const MCParticle* dauMcPart = mcPart->getParticleFromGeneralizedIndexString(generalizedIndex);
1324 if (dauMcPart ==
nullptr)
1327 pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1331 for (
auto& generalizedIndex : arguments) {
1332 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1333 if (dauPart ==
nullptr)
1336 const MCParticle* dauMcPart = dauPart->getMCParticle();
1337 if (dauMcPart ==
nullptr)
1340 pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1345 if (pDaus.size() == 2)
1346 return ROOT::Math::VectorUtil::Angle(pDaus[0], pDaus[1]);
1348 return ROOT::Math::VectorUtil::Angle(pDaus[2], pDaus[0] + pDaus[1]);
1352 B2FATAL(
"Wrong number of arguments for meta function mcDaughterAngle");
1356 double daughterClusterAngleInBetween(
const Particle* particle,
const std::vector<double>& daughterIndices)
1358 if (daughterIndices.size() == 2) {
1359 int daughterIndexi = std::lround(daughterIndices[0]);
1360 int daughterIndexj = std::lround(daughterIndices[1]);
1361 if (std::max(daughterIndexi, daughterIndexj) >=
int(particle->getNDaughters())) {
1364 const ECLCluster* clusteri = particle->getDaughter(daughterIndexi)->getECLCluster();
1365 const ECLCluster* clusterj = particle->getDaughter(daughterIndexj)->getECLCluster();
1366 if (clusteri and clusterj) {
1370 ClusterUtils clusutils;
1371 ROOT::Math::PxPyPzEVector pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit));
1372 ROOT::Math::PxPyPzEVector pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit));
1373 return ROOT::Math::VectorUtil::Angle(pi, pj);
1377 }
else if (daughterIndices.size() == 3) {
1378 int daughterIndexi = std::lround(daughterIndices[0]);
1379 int daughterIndexj = std::lround(daughterIndices[1]);
1380 int daughterIndexk = std::lround(daughterIndices[2]);
1381 if (std::max(std::max(daughterIndexi, daughterIndexj), daughterIndexk) >=
int(particle->getNDaughters())) {
1384 const ECLCluster* clusteri = (particle->getDaughter(daughterIndices[0]))->getECLCluster();
1385 const ECLCluster* clusterj = (particle->getDaughter(daughterIndices[1]))->getECLCluster();
1386 const ECLCluster* clusterk = (particle->getDaughter(daughterIndices[2]))->getECLCluster();
1387 if (clusteri and clusterj and clusterk) {
1392 ClusterUtils clusutils;
1393 ROOT::Math::PxPyPzEVector pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit));
1394 ROOT::Math::PxPyPzEVector pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit));
1395 ROOT::Math::PxPyPzEVector pk = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterk, clusterkBit));
1396 return ROOT::Math::VectorUtil::Angle(pk, pi + pj);
1401 B2FATAL(
"Wrong number of arguments for daughterClusterAngleInBetween!");
1407 if (arguments.size() > 1) {
1408 auto func = [arguments](
const Particle * particle) ->
double {
1410 ROOT::Math::PxPyPzEVector pSum;
1412 for (
auto& generalizedIndex : arguments)
1414 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1416 pSum += frame.getMomentum(dauPart);
1425 B2FATAL(
"Wrong number of arguments for meta function daughterInvM. At least two integers are needed.");
1431 if (arguments.size() == 2) {
1436 }
catch (std::invalid_argument&) {
1437 B2FATAL(
"Second argument of modulo meta function must be integer!");
1439 auto func = [var, divideBy](
const Particle * particle) ->
int {
1440 auto var_result = var->function(particle);
1441 if (std::holds_alternative<double>(var_result))
1443 return int(std::get<double>(var_result)) % divideBy;
1444 }
else if (std::holds_alternative<int>(var_result))
1446 return std::get<int>(var_result) % divideBy;
1447 }
else if (std::holds_alternative<bool>(var_result))
1449 return int(std::get<bool>(var_result)) % divideBy;
1454 B2FATAL(
"Wrong number of arguments for meta function modulo");
1460 if (arguments.size() == 1) {
1463 auto func = [var](
const Particle * particle) ->
bool {
return std::isnan(std::get<double>(var->function(particle))); };
1466 B2FATAL(
"Wrong number of arguments for meta function isNAN");
1472 if (arguments.size() == 2) {
1474 double defaultOutput;
1477 }
catch (std::invalid_argument&) {
1478 B2FATAL(
"The second argument of ifNANgiveX meta function must be a number!");
1480 auto func = [var, defaultOutput](
const Particle * particle) ->
double {
1481 double output = std::get<double>(var->function(particle));
1482 if (std::isnan(output))
return defaultOutput;
1487 B2FATAL(
"Wrong number of arguments for meta function ifNANgiveX");
1493 if (arguments.size() == 1) {
1496 auto func = [var](
const Particle * particle) ->
bool {
return std::isinf(std::get<double>(var->function(particle))); };
1499 B2FATAL(
"Wrong number of arguments for meta function isInfinity");
1505 if (arguments.size() >= 2) {
1511 for (
size_t i = 1; i < arguments.size(); ++i) {
1514 }
catch (std::invalid_argument&) {
1515 B2FATAL(
"The input flags to meta function unmask() should be integer!");
1521 auto func = [var, finalMask](
const Particle * particle) ->
double {
1523 auto var_result = var->function(particle);
1524 if (std::holds_alternative<double>(var_result))
1527 if (std::isnan(std::get<double>(var_result))) {
1530 value = int(std::get<double>(var_result));
1531 }
else if (std::holds_alternative<int>(var_result))
1533 value = std::get<int>(var_result);
1537 value &= (~finalMask);
1544 B2FATAL(
"Meta function unmask needs at least two arguments!");
1550 if (arguments.size() == 3) {
1552 std::string cutString = arguments[0];
1558 auto func = [cut, variableIfTrue, variableIfFalse](
const Particle * particle) ->
double {
1559 if (particle ==
nullptr)
1561 if (cut->check(particle))
1563 auto var_result = variableIfTrue->function(particle);
1564 if (std::holds_alternative<double>(var_result)) {
1565 return std::get<double>(var_result);
1566 }
else if (std::holds_alternative<int>(var_result)) {
1567 return std::get<int>(var_result);
1568 }
else if (std::holds_alternative<bool>(var_result)) {
1569 return std::get<bool>(var_result);
1573 auto var_result = variableIfFalse->function(particle);
1574 if (std::holds_alternative<double>(var_result)) {
1575 return std::get<double>(var_result);
1576 }
else if (std::holds_alternative<int>(var_result)) {
1577 return std::get<int>(var_result);
1578 }
else if (std::holds_alternative<bool>(var_result)) {
1579 return std::get<bool>(var_result);
1586 B2FATAL(
"Wrong number of arguments for meta function conditionalVariableSelector");
1592 if (arguments.size() > 0) {
1593 std::vector<const Variable::Manager::Var*> variables;
1594 for (
auto& argument : arguments)
1597 auto func = [variables, arguments](
const Particle * particle) ->
double {
1598 double pValueProduct = 1.;
1599 for (
auto variable : variables)
1601 double pValue = std::get<double>(variable->function(particle));
1605 pValueProduct *= pValue;
1607 double pValueSum = 1.;
1608 double factorial = 1.;
1609 for (
unsigned int i = 1; i < arguments.size(); ++i)
1612 pValueSum += pow(-std::log(pValueProduct), i) / factorial;
1614 return pValueProduct * pValueSum;
1618 B2FATAL(
"Wrong number of arguments for meta function pValueCombination");
1624 if (arguments.size() == 1) {
1626 auto func = [var](
const Particle * particle) ->
double {
1627 double pValueProduct = 1.;
1628 if (particle->getNDaughters() == 0)
1633 for (
unsigned j = 0; j < particle->getNDaughters(); ++j)
1635 double pValue = std::get<double>(var->function(particle->getDaughter(j)));
1636 if (pValue < 0)
return -1;
1637 else pValueProduct *= pValue;
1640 double pValueSum = 1.;
1641 double factorial = 1.;
1642 for (
unsigned int i = 1; i < particle->getNDaughters(); ++i)
1645 pValueSum += pow(-std::log(pValueProduct), i) / factorial;
1647 return pValueProduct * pValueSum;
1651 B2FATAL(
"Wrong number of arguments for meta function pValueCombinationOfDaughters");
1657 if (arguments.size() == 1) {
1659 auto func = [var](
const Particle * particle) ->
double {
1660 auto var_result = var->function(particle);
1661 if (std::holds_alternative<double>(var_result))
1663 return std::abs(std::get<double>(var_result));
1664 }
else if (std::holds_alternative<int>(var_result))
1666 return std::abs(std::get<int>(var_result));
1671 B2FATAL(
"Wrong number of arguments for meta function abs");
1677 if (arguments.size() == 2) {
1682 B2FATAL(
"One or both of the used variables doesn't exist!");
1684 auto func = [var1, var2](
const Particle * particle) ->
double {
1686 auto var_result1 = var1->function(particle);
1687 auto var_result2 = var2->function(particle);
1688 if (std::holds_alternative<double>(var_result1))
1690 val1 = std::get<double>(var_result1);
1691 }
else if (std::holds_alternative<int>(var_result1))
1693 val1 = std::get<int>(var_result1);
1694 }
else if (std::holds_alternative<bool>(var_result1))
1696 val1 = std::get<bool>(var_result1);
1699 B2FATAL(
"A variable in meta function max holds no double, int or bool values");
1701 if (std::holds_alternative<double>(var_result2))
1703 val2 = std::get<double>(var_result2);
1704 }
else if (std::holds_alternative<int>(var_result2))
1706 val2 = std::get<int>(var_result2);
1707 }
else if (std::holds_alternative<bool>(var_result2))
1709 val2 = std::get<bool>(var_result2);
1712 B2FATAL(
"A variable in meta function max holds no double, int or bool values");
1714 return std::max(val1, val2);
1718 B2FATAL(
"Wrong number of arguments for meta function max");
1724 if (arguments.size() == 2) {
1729 B2FATAL(
"One or both of the used variables doesn't exist!");
1731 auto func = [var1, var2](
const Particle * particle) ->
double {
1733 auto var_result1 = var1->function(particle);
1734 auto var_result2 = var2->function(particle);
1735 if (std::holds_alternative<double>(var_result1))
1737 val1 = std::get<double>(var_result1);
1738 }
else if (std::holds_alternative<int>(var_result1))
1740 val1 = std::get<int>(var_result1);
1741 }
else if (std::holds_alternative<bool>(var_result1))
1743 val1 = std::get<bool>(var_result1);
1746 B2FATAL(
"A variable in meta function min holds no double, int or bool values");
1748 if (std::holds_alternative<double>(var_result2))
1750 val2 = std::get<double>(var_result2);
1751 }
else if (std::holds_alternative<int>(var_result2))
1753 val2 = std::get<int>(var_result2);
1754 }
else if (std::holds_alternative<bool>(var_result2))
1756 val2 = std::get<bool>(var_result2);
1759 B2FATAL(
"A variable in meta function min holds no double, int or bool values");
1761 return std::min(val1, val2);
1765 B2FATAL(
"Wrong number of arguments for meta function min");
1771 if (arguments.size() == 1) {
1773 auto func = [var](
const Particle * particle) ->
double {
1774 auto var_result = var->function(particle);
1775 if (std::holds_alternative<double>(var_result))
1776 return std::sin(std::get<double>(var_result));
1777 else if (std::holds_alternative<int>(var_result))
1778 return std::sin(std::get<int>(var_result));
1783 B2FATAL(
"Wrong number of arguments for meta function sin");
1789 if (arguments.size() == 1) {
1791 auto func = [var](
const Particle * particle) ->
double {
1792 auto var_result = var->function(particle);
1793 if (std::holds_alternative<double>(var_result))
1794 return std::asin(std::get<double>(var_result));
1795 else if (std::holds_alternative<int>(var_result))
1796 return std::asin(std::get<int>(var_result));
1801 B2FATAL(
"Wrong number of arguments for meta function asin");
1807 if (arguments.size() == 1) {
1809 auto func = [var](
const Particle * particle) ->
double {
1810 auto var_result = var->function(particle);
1811 if (std::holds_alternative<double>(var_result))
1812 return std::cos(std::get<double>(var_result));
1813 else if (std::holds_alternative<int>(var_result))
1814 return std::cos(std::get<int>(var_result));
1819 B2FATAL(
"Wrong number of arguments for meta function cos");
1825 if (arguments.size() == 1) {
1827 auto func = [var](
const Particle * particle) ->
double {
1828 auto var_result = var->function(particle);
1829 if (std::holds_alternative<double>(var_result))
1830 return std::acos(std::get<double>(var_result));
1831 else if (std::holds_alternative<int>(var_result))
1832 return std::acos(std::get<int>(var_result));
1837 B2FATAL(
"Wrong number of arguments for meta function acos");
1843 if (arguments.size() == 1) {
1845 auto func = [var](
const Particle * particle) ->
double {
return std::tan(std::get<double>(var->function(particle))); };
1848 B2FATAL(
"Wrong number of arguments for meta function tan");
1854 if (arguments.size() == 1) {
1856 auto func = [var](
const Particle * particle) ->
double {
return std::atan(std::get<double>(var->function(particle))); };
1859 B2FATAL(
"Wrong number of arguments for meta function atan");
1865 if (arguments.size() == 1) {
1867 auto func = [var](
const Particle * particle) ->
double {
1868 auto var_result = var->function(particle);
1869 if (std::holds_alternative<double>(var_result))
1870 return std::exp(std::get<double>(var_result));
1871 else if (std::holds_alternative<int>(var_result))
1872 return std::exp(std::get<int>(var_result));
1877 B2FATAL(
"Wrong number of arguments for meta function exp");
1883 if (arguments.size() == 1) {
1885 auto func = [var](
const Particle * particle) ->
double {
1886 auto var_result = var->function(particle);
1887 if (std::holds_alternative<double>(var_result))
1888 return std::log(std::get<double>(var_result));
1889 else if (std::holds_alternative<int>(var_result))
1890 return std::log(std::get<int>(var_result));
1895 B2FATAL(
"Wrong number of arguments for meta function log");
1901 if (arguments.size() == 1) {
1903 auto func = [var](
const Particle * particle) ->
double {
1904 auto var_result = var->function(particle);
1905 if (std::holds_alternative<double>(var_result))
1906 return std::log10(std::get<double>(var_result));
1907 else if (std::holds_alternative<int>(var_result))
1908 return std::log10(std::get<int>(var_result));
1913 B2FATAL(
"Wrong number of arguments for meta function log10");
1919 if (arguments.size() == 1) {
1921 auto func = [var](
const Particle * particle) ->
double {
1922 if (particle ==
nullptr)
1925 StoreArray<Particle> particles;
1926 if (!particle->hasExtraInfo(
"original_index"))
1929 auto originalParticle = particles[particle->getExtraInfo(
"original_index")];
1930 if (!originalParticle)
1932 auto var_result = var->function(originalParticle);
1933 if (std::holds_alternative<double>(var_result))
1935 return std::get<double>(var_result);
1936 }
else if (std::holds_alternative<int>(var_result))
1938 return std::get<int>(var_result);
1939 }
else if (std::holds_alternative<bool>(var_result))
1941 return std::get<bool>(var_result);
1946 B2FATAL(
"Wrong number of arguments for meta function originalParticle");
1952 if (arguments.size() == 2) {
1953 auto daughterFunction = convertToDaughterIndex({arguments[0]});
1955 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
1956 if (particle ==
nullptr)
1958 int daughterNumber = std::get<int>(daughterFunction(particle));
1959 if (daughterNumber >=
int(particle->getNDaughters()) or daughterNumber < 0)
1961 auto var_result = var->function(particle->getDaughter(daughterNumber));
1962 if (std::holds_alternative<double>(var_result))
1964 return std::get<double>(var_result);
1965 }
else if (std::holds_alternative<int>(var_result))
1967 return std::get<int>(var_result);
1968 }
else if (std::holds_alternative<bool>(var_result))
1970 return std::get<bool>(var_result);
1975 B2FATAL(
"Wrong number of arguments for meta function daughter");
1981 if (arguments.size() == 2) {
1982 auto daughterFunction = convertToDaughterIndex({arguments[0]});
1984 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
1985 if (particle ==
nullptr)
1987 int daughterNumber = std::get<int>(daughterFunction(particle));
1988 if (daughterNumber >=
int(particle->getNDaughters()) or daughterNumber < 0)
1992 StoreArray<Particle> particles;
1993 if (!particle->getDaughter(daughterNumber)->hasExtraInfo(
"original_index"))
1995 auto originalDaughter = particles[particle->getDaughter(daughterNumber)->getExtraInfo(
"original_index")];
1996 if (!originalDaughter)
1999 auto var_result = var->function(originalDaughter);
2000 if (std::holds_alternative<double>(var_result)) {
2001 return std::get<double>(var_result);
2002 }
else if (std::holds_alternative<int>(var_result)) {
2003 return std::get<int>(var_result);
2004 }
else if (std::holds_alternative<bool>(var_result)) {
2005 return std::get<bool>(var_result);
2011 B2FATAL(
"Wrong number of arguments for meta function daughter");
2017 if (arguments.size() == 1) {
2018 std::string daughterString = arguments[0];
2019 auto func = [daughterString](
const Particle * particle) ->
int {
2020 if (particle ==
nullptr)
2022 int daughterNumber = 0;
2026 }
catch (std::invalid_argument&)
2028 auto daughterFunction = convertToInt({daughterString,
"-1"});
2029 auto daughterVarResult = daughterFunction(particle);
2030 daughterNumber = std::get<int>(daughterVarResult);
2032 return daughterNumber;
2036 B2FATAL(
"Wrong number of arguments for meta function convertToDaughterIndex");
2042 if (arguments.size() == 2) {
2043 auto daughterFunction = convertToDaughterIndex({arguments[0]});
2045 auto func = [var, daughterFunction](
const Particle * particle) ->
double {
2046 if (particle ==
nullptr)
2048 if (particle->getMCParticle())
2050 int daughterNumber = std::get<int>(daughterFunction(particle));
2051 if (daughterNumber >=
int(particle->getMCParticle()->getNDaughters()) or daughterNumber < 0)
2053 Particle tempParticle = Particle(particle->getMCParticle()->getDaughters().at(daughterNumber));
2054 auto var_result = var->function(&tempParticle);
2055 if (std::holds_alternative<double>(var_result)) {
2056 return std::get<double>(var_result);
2057 }
else if (std::holds_alternative<int>(var_result)) {
2058 return std::get<int>(var_result);
2059 }
else if (std::holds_alternative<bool>(var_result)) {
2060 return std::get<bool>(var_result);
2071 B2FATAL(
"Wrong number of arguments for meta function mcDaughter");
2077 if (arguments.size() == 1) {
2079 auto func = [var](
const Particle * particle) ->
double {
2080 if (particle ==
nullptr)
2082 if (particle->getMCParticle())
2084 if (particle->getMCParticle()->getMother() ==
nullptr) {
2087 Particle tempParticle = Particle(particle->getMCParticle()->getMother());
2088 auto var_result = var->function(&tempParticle);
2089 if (std::holds_alternative<double>(var_result)) {
2090 return std::get<double>(var_result);
2091 }
else if (std::holds_alternative<int>(var_result)) {
2092 return std::get<int>(var_result);
2093 }
else if (std::holds_alternative<bool>(var_result)) {
2094 return std::get<bool>(var_result);
2103 B2FATAL(
"Wrong number of arguments for meta function mcMother");
2109 if (arguments.size() == 2) {
2110 std::string indexString = arguments[0];
2113 auto func = [var, indexString](
const Particle * particle) ->
double {
2115 int particleNumber = 0;
2119 }
catch (std::invalid_argument&)
2121 auto indexFunction = convertToInt({indexString,
"-1"});
2122 auto indexVarResult = indexFunction(particle);
2123 particleNumber = std::get<int>(indexVarResult);
2126 StoreArray<MCParticle> mcParticles(
"MCParticles");
2127 if (particleNumber >= mcParticles.getEntries())
2132 MCParticle* mcParticle = mcParticles[particleNumber];
2133 Particle part = Particle(mcParticle);
2134 auto var_result = var->function(&part);
2135 if (std::holds_alternative<double>(var_result))
2137 return std::get<double>(var_result);
2138 }
else if (std::holds_alternative<int>(var_result))
2140 return std::get<int>(var_result);
2141 }
else if (std::holds_alternative<bool>(var_result))
2143 return std::get<bool>(var_result);
2148 B2FATAL(
"Wrong number of arguments for meta function genParticle");
2154 if (arguments.size() == 1) {
2157 auto func = [var](
const Particle*) ->
double {
2158 StoreArray<MCParticle> mcParticles(
"MCParticles");
2159 if (mcParticles.getEntries() == 0)
2164 MCParticle* mcUpsilon4S = mcParticles[0];
2165 if (mcUpsilon4S->isInitial()) mcUpsilon4S = mcParticles[2];
2166 if (mcUpsilon4S->getPDG() != 300553)
2171 Particle upsilon4S = Particle(mcUpsilon4S);
2172 auto var_result = var->function(&upsilon4S);
2173 if (std::holds_alternative<double>(var_result))
2175 return std::get<double>(var_result);
2176 }
else if (std::holds_alternative<int>(var_result))
2178 return std::get<int>(var_result);
2179 }
else if (std::holds_alternative<bool>(var_result))
2181 return std::get<bool>(var_result);
2186 B2FATAL(
"Wrong number of arguments for meta function genUpsilon4S");
2192 if (arguments.size() == 4) {
2193 std::string listName = arguments[0];
2194 std::string rankedVariableName = arguments[1];
2195 std::string returnVariableName = arguments[2];
2196 std::string extraInfoName = rankedVariableName +
"_rank";
2200 }
catch (std::invalid_argument&) {
2201 B2ERROR(
"3rd argument of getVariableByRank meta function (Rank) must be an integer!");
2206 auto func = [var, rank, extraInfoName, listName](
const Particle*)->
double {
2207 StoreObjPtr<ParticleList> list(listName);
2209 const unsigned int numParticles = list->getListSize();
2210 for (
unsigned int i = 0; i < numParticles; i++)
2212 const Particle* p = list->getParticle(i);
2213 if (p->getExtraInfo(extraInfoName) == rank) {
2214 auto var_result = var->function(p);
2215 if (std::holds_alternative<double>(var_result)) {
2216 return std::get<double>(var_result);
2217 }
else if (std::holds_alternative<int>(var_result)) {
2218 return std::get<int>(var_result);
2219 }
else if (std::holds_alternative<bool>(var_result)) {
2220 return std::get<bool>(var_result);
2225 return std::numeric_limits<double>::signaling_NaN();
2229 B2FATAL(
"Wrong number of arguments for meta function getVariableByRank");
2235 if (arguments.size() == 1 or arguments.size() == 2) {
2237 std::string listName = arguments[0];
2238 std::string cutString =
"";
2240 if (arguments.size() == 2) {
2241 cutString = arguments[1];
2246 auto func = [listName, cut](
const Particle*) ->
int {
2248 StoreObjPtr<ParticleList> list(listName);
2250 for (
unsigned int i = 0; i < list->getListSize(); i++)
2252 const Particle* particle = list->getParticle(i);
2253 if (cut->check(particle)) {
2261 B2FATAL(
"Wrong number of arguments for meta function countInList");
2267 if (arguments.size() == 2 or arguments.size() == 3) {
2269 std::string roeListName = arguments[0];
2270 std::string cutString = arguments[1];
2272 if (arguments.size() == 2) {
2273 B2INFO(
"Use pdgCode of electron as default in meta variable veto, other arguments: " << roeListName <<
", " << cutString);
2277 }
catch (std::invalid_argument&) {
2278 B2FATAL(
"Third argument of veto meta function must be integer!");
2285 auto func = [roeListName, cut, pdgCode, flavourType](
const Particle * particle) ->
bool {
2286 StoreObjPtr<ParticleList> roeList(roeListName);
2287 ROOT::Math::PxPyPzEVector vec = particle->get4Vector();
2288 for (
unsigned int i = 0; i < roeList->getListSize(); i++)
2290 const Particle* roeParticle = roeList->getParticle(i);
2291 if (not particle->overlapsWith(roeParticle)) {
2292 ROOT::Math::PxPyPzEVector tempCombination = roeParticle->get4Vector() + vec;
2293 std::vector<int> indices = { particle->getArrayIndex(), roeParticle->getArrayIndex() };
2294 Particle tempParticle = Particle(tempCombination, pdgCode, flavourType, indices, particle->getArrayPointer());
2295 if (cut->check(&tempParticle)) {
2304 B2FATAL(
"Wrong number of arguments for meta function veto");
2310 if (arguments.size() == 1) {
2311 std::string cutString = arguments[0];
2313 auto func = [cut](
const Particle * particle) ->
int {
2315 for (
auto& daughter : particle->getDaughters())
2317 if (cut->check(daughter))
2324 B2FATAL(
"Wrong number of arguments for meta function countDaughters");
2330 if (arguments.size() == 1) {
2331 std::string cutString = arguments[0];
2333 auto func = [cut](
const Particle * particle) ->
int {
2335 std::vector<const Particle*> fspDaughters;
2336 particle->fillFSPDaughters(fspDaughters);
2339 for (
auto& daughter : fspDaughters)
2341 if (cut->check(daughter))
2348 B2FATAL(
"Wrong number of arguments for meta function countFSPDaughters");
2354 if (arguments.size() == 1) {
2355 std::string cutString = arguments[0];
2357 auto func = [cut](
const Particle * particle) ->
int {
2359 std::vector<const Particle*> allDaughters;
2360 particle->fillAllDaughters(allDaughters);
2363 for (
auto& daughter : allDaughters)
2365 if (cut->check(daughter))
2372 B2FATAL(
"Wrong number of arguments for meta function countDescendants");
2376 Manager::FunctionPtr numberOfNonOverlappingParticles(
const std::vector<std::string>& arguments)
2379 auto func = [arguments](
const Particle * particle) ->
int {
2381 int _numberOfNonOverlappingParticles = 0;
2382 for (
const auto& listName : arguments)
2384 StoreObjPtr<ParticleList> list(listName);
2385 if (not list.isValid()) {
2386 B2FATAL(
"Invalid list named " << listName <<
" encountered in numberOfNonOverlappingParticles.");
2388 for (
unsigned int i = 0; i < list->getListSize(); i++) {
2389 const Particle* p = list->getParticle(i);
2390 if (not particle->overlapsWith(p)) {
2391 _numberOfNonOverlappingParticles++;
2395 return _numberOfNonOverlappingParticles;
2402 void appendDaughtersRecursive(Particle* mother, StoreArray<Particle>& container)
2405 auto* mcmother = mother->getRelated<MCParticle>();
2410 for (
auto* mcdaughter : mcmother->getDaughters()) {
2412 Particle tmp_daughter(mcdaughter);
2413 Particle* new_daughter = container.appendNew(tmp_daughter);
2414 new_daughter->addRelationTo(mcdaughter);
2415 mother->appendDaughter(new_daughter,
false);
2417 if (mcdaughter->getNDaughters() > 0)
2418 appendDaughtersRecursive(new_daughter, container);
2424 if (arguments.size() == 1) {
2426 auto func = [var](
const Particle * particle) ->
double {
2427 const MCParticle* mcp = particle->getMCParticle();
2432 StoreArray<Particle> tempParticles(
"tempParticles");
2433 tempParticles.clear();
2434 Particle tmpPart(mcp);
2435 Particle* newPart = tempParticles.appendNew(tmpPart);
2436 newPart->addRelationTo(mcp);
2438 appendDaughtersRecursive(newPart, tempParticles);
2440 auto var_result = var->function(newPart);
2441 if (std::holds_alternative<double>(var_result))
2443 return std::get<double>(var_result);
2444 }
else if (std::holds_alternative<int>(var_result))
2446 return std::get<int>(var_result);
2447 }
else if (std::holds_alternative<bool>(var_result))
2449 return std::get<bool>(var_result);
2454 B2FATAL(
"Wrong number of arguments for meta function matchedMC");
2460 if (arguments.size() == 1) {
2463 auto func = [var](
const Particle * particle) ->
double {
2465 const ECLCluster* cluster = particle->getECLCluster();
2468 auto mcps = cluster->getRelationsTo<MCParticle>();
2471 std::vector<std::pair<double, int>> weightsAndIndices;
2472 for (
unsigned int i = 0; i < mcps.size(); ++i)
2473 weightsAndIndices.emplace_back(mcps.weight(i), i);
2476 std::sort(weightsAndIndices.begin(), weightsAndIndices.end(),
2477 ValueIndexPairSorting::higherPair<
decltype(weightsAndIndices)::value_type>);
2480 const MCParticle* mcp = mcps.object(weightsAndIndices[0].second);
2482 StoreArray<Particle> tempParticles(
"tempParticles");
2483 tempParticles.clear();
2484 Particle tmpPart(mcp);
2485 Particle* newPart = tempParticles.appendNew(tmpPart);
2486 newPart->addRelationTo(mcp);
2488 appendDaughtersRecursive(newPart, tempParticles);
2490 auto var_result = var->function(newPart);
2491 if (std::holds_alternative<double>(var_result))
2493 return std::get<double>(var_result);
2494 }
else if (std::holds_alternative<int>(var_result))
2496 return std::get<int>(var_result);
2497 }
else if (std::holds_alternative<bool>(var_result))
2499 return std::get<bool>(var_result);
2508 B2FATAL(
"Wrong number of arguments for meta function clusterBestMatchedMCParticle");
2514 if (arguments.size() == 1) {
2517 auto func = [var](
const Particle * particle) ->
double {
2519 const ECLCluster* cluster = particle->getECLCluster();
2522 auto mcps = cluster->getRelationsTo<MCParticle>();
2525 std::map<int, double> mapMCParticleIndxAndWeight;
2526 getKlongWeightMap(particle, mapMCParticleIndxAndWeight);
2529 if (mapMCParticleIndxAndWeight.size() == 0)
2533 auto maxMap = std::max_element(mapMCParticleIndxAndWeight.begin(), mapMCParticleIndxAndWeight.end(),
2534 [](
const auto & x,
const auto & y) { return x.second < y.second; }
2537 StoreArray<MCParticle> mcparticles;
2538 const MCParticle* mcKlong = mcparticles[maxMap->first];
2540 Particle tmpPart(mcKlong);
2541 auto var_result = var->function(&tmpPart);
2542 if (std::holds_alternative<double>(var_result))
2544 return std::get<double>(var_result);
2545 }
else if (std::holds_alternative<int>(var_result))
2547 return std::get<int>(var_result);
2548 }
else if (std::holds_alternative<bool>(var_result))
2550 return std::get<bool>(var_result);
2559 B2FATAL(
"Wrong number of arguments for meta function clusterBestMatchedMCKlong");
2563 double matchedMCHasPDG(
const Particle* particle,
const std::vector<double>& pdgCode)
2565 if (pdgCode.size() != 1) {
2566 B2FATAL(
"Too many arguments provided to matchedMCHasPDG!");
2568 int inputPDG = std::lround(pdgCode[0]);
2570 const MCParticle* mcp = particle->getMCParticle();
2574 return std::abs(mcp->getPDG()) == inputPDG;
2579 if (arguments.size() == 1) {
2580 std::string listName = arguments[0];
2581 auto func = [listName](
const Particle * particle) ->
double {
2584 StoreObjPtr<ParticleList> listOfParticles(listName);
2586 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalEnergyOfParticlesInList");
2587 double totalEnergy = 0;
2588 int nParticles = listOfParticles->getListSize();
2589 for (
int i = 0; i < nParticles; i++)
2591 const Particle* part = listOfParticles->getParticle(i);
2593 totalEnergy += frame.getMomentum(part).E();
2600 B2FATAL(
"Wrong number of arguments for meta function totalEnergyOfParticlesInList");
2606 if (arguments.size() == 1) {
2607 std::string listName = arguments[0];
2608 auto func = [listName](
const Particle*) ->
double {
2609 StoreObjPtr<ParticleList> listOfParticles(listName);
2611 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPxOfParticlesInList");
2613 int nParticles = listOfParticles->getListSize();
2615 for (
int i = 0; i < nParticles; i++)
2617 const Particle* part = listOfParticles->getParticle(i);
2618 totalPx += frame.getMomentum(part).Px();
2624 B2FATAL(
"Wrong number of arguments for meta function totalPxOfParticlesInList");
2630 if (arguments.size() == 1) {
2631 std::string listName = arguments[0];
2632 auto func = [listName](
const Particle*) ->
double {
2633 StoreObjPtr<ParticleList> listOfParticles(listName);
2635 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPyOfParticlesInList");
2637 int nParticles = listOfParticles->getListSize();
2639 for (
int i = 0; i < nParticles; i++)
2641 const Particle* part = listOfParticles->getParticle(i);
2642 totalPy += frame.getMomentum(part).Py();
2648 B2FATAL(
"Wrong number of arguments for meta function totalPyOfParticlesInList");
2654 if (arguments.size() == 1) {
2655 std::string listName = arguments[0];
2656 auto func = [listName](
const Particle*) ->
double {
2657 StoreObjPtr<ParticleList> listOfParticles(listName);
2659 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalPzOfParticlesInList");
2661 int nParticles = listOfParticles->getListSize();
2663 for (
int i = 0; i < nParticles; i++)
2665 const Particle* part = listOfParticles->getParticle(i);
2666 totalPz += frame.getMomentum(part).Pz();
2672 B2FATAL(
"Wrong number of arguments for meta function totalPzOfParticlesInList");
2678 if (arguments.size() > 0) {
2680 auto func = [arguments](
const Particle * particle) ->
double {
2682 ROOT::Math::PxPyPzEVector total4Vector;
2684 std::vector<Particle*> particlePool;
2687 for (
const auto& argument : arguments)
2689 StoreObjPtr <ParticleList> listOfParticles(argument);
2691 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << argument <<
" given to invMassInLists");
2692 int nParticles = listOfParticles->getListSize();
2693 for (
int i = 0; i < nParticles; i++) {
2694 bool overlaps =
false;
2695 Particle* part = listOfParticles->getParticle(i);
2696 for (
auto poolPart : particlePool) {
2697 if (part->overlapsWith(poolPart)) {
2703 total4Vector += part->get4Vector();
2704 particlePool.push_back(part);
2708 double invariantMass = total4Vector.M();
2709 return invariantMass;
2714 B2FATAL(
"Wrong number of arguments for meta function invMassInLists");
2718 Manager::FunctionPtr totalECLEnergyOfParticlesInList(
const std::vector<std::string>& arguments)
2720 if (arguments.size() == 1) {
2721 std::string listName = arguments[0];
2722 auto func = [listName](
const Particle * particle) ->
double {
2725 StoreObjPtr<ParticleList> listOfParticles(listName);
2727 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to totalEnergyOfParticlesInList");
2728 double totalEnergy = 0;
2729 int nParticles = listOfParticles->getListSize();
2730 for (
int i = 0; i < nParticles; i++)
2732 const Particle* part = listOfParticles->getParticle(i);
2733 const ECLCluster* cluster = part->getECLCluster();
2735 if (cluster !=
nullptr) {
2736 totalEnergy += cluster->getEnergy(clusterHypothesis);
2744 B2FATAL(
"Wrong number of arguments for meta function totalECLEnergyOfParticlesInList");
2750 if (arguments.size() == 1) {
2751 std::string listName = arguments[0];
2752 auto func = [listName](
const Particle*) ->
double {
2753 StoreObjPtr<ParticleList> listOfParticles(listName);
2755 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to maxPtInList");
2756 int nParticles = listOfParticles->getListSize();
2759 for (
int i = 0; i < nParticles; i++)
2761 const Particle* part = listOfParticles->getParticle(i);
2762 const double Pt = frame.getMomentum(part).Pt();
2763 if (Pt > maxPt) maxPt = Pt;
2769 B2FATAL(
"Wrong number of arguments for meta function maxPtInList");
2773 Manager::FunctionPtr eclClusterTrackMatchedWithCondition(
const std::vector<std::string>& arguments)
2775 if (arguments.size() <= 1) {
2777 std::string cutString;
2778 if (arguments.size() == 1)
2779 cutString = arguments[0];
2781 auto func = [cut](
const Particle * particle) ->
double {
2783 if (particle ==
nullptr)
2786 const ECLCluster* cluster = particle->getECLCluster();
2790 auto tracks = cluster->getRelationsFrom<Track>();
2792 for (
const auto& track : tracks) {
2795 if (cut->check(&trackParticle))
2804 B2FATAL(
"Wrong number of arguments for meta function eclClusterSpecialTrackMatched");
2810 if (arguments.size() == 2) {
2811 std::string listName = arguments[0];
2814 auto func = [listName, var](
const Particle*) ->
double {
2815 StoreObjPtr<ParticleList> listOfParticles(listName);
2817 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to averageValueInList");
2818 int nParticles = listOfParticles->getListSize();
2819 if (nParticles == 0)
2824 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2826 for (
int i = 0; i < nParticles; i++) {
2827 average += std::get<double>(var->function(listOfParticles->getParticle(i))) / nParticles;
2829 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2831 for (
int i = 0; i < nParticles; i++) {
2832 average += std::get<int>(var->function(listOfParticles->getParticle(i))) / nParticles;
2839 B2FATAL(
"Wrong number of arguments for meta function averageValueInList");
2845 if (arguments.size() == 2) {
2846 std::string listName = arguments[0];
2849 auto func = [listName, var](
const Particle*) ->
double {
2850 StoreObjPtr<ParticleList> listOfParticles(listName);
2852 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to medianValueInList");
2853 int nParticles = listOfParticles->getListSize();
2854 if (nParticles == 0)
2858 std::vector<double> valuesInList;
2859 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2861 for (
int i = 0; i < nParticles; i++) {
2862 valuesInList.push_back(std::get<double>(var->function(listOfParticles->getParticle(i))));
2864 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2866 for (
int i = 0; i < nParticles; i++) {
2867 valuesInList.push_back(std::get<int>(var->function(listOfParticles->getParticle(i))));
2870 std::sort(valuesInList.begin(), valuesInList.end());
2871 if (nParticles % 2 != 0)
2873 return valuesInList[nParticles / 2];
2876 return 0.5 * (valuesInList[nParticles / 2] + valuesInList[nParticles / 2 - 1]);
2881 B2FATAL(
"Wrong number of arguments for meta function medianValueInList");
2887 if (arguments.size() == 2) {
2888 std::string listName = arguments[0];
2891 auto func = [listName, var](
const Particle*) ->
double {
2892 StoreObjPtr<ParticleList> listOfParticles(listName);
2894 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to sumValueInList");
2895 int nParticles = listOfParticles->getListSize();
2896 if (nParticles == 0)
2901 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2903 for (
int i = 0; i < nParticles; i++) {
2904 sum += std::get<double>(var->function(listOfParticles->getParticle(i)));
2906 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2908 for (
int i = 0; i < nParticles; i++) {
2909 sum += std::get<int>(var->function(listOfParticles->getParticle(i)));
2916 B2FATAL(
"Wrong number of arguments for meta function sumValueInList");
2922 if (arguments.size() == 2) {
2923 std::string listName = arguments[0];
2926 auto func = [listName, var](
const Particle*) ->
double {
2927 StoreObjPtr<ParticleList> listOfParticles(listName);
2929 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid list name " << listName <<
" given to productValueInList");
2930 int nParticles = listOfParticles->getListSize();
2931 if (nParticles == 0)
2936 if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2938 for (
int i = 0; i < nParticles; i++) {
2939 product *= std::get<double>(var->function(listOfParticles->getParticle(i)));
2941 }
else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2943 for (
int i = 0; i < nParticles; i++) {
2944 product *= std::get<int>(var->function(listOfParticles->getParticle(i)));
2951 B2FATAL(
"Wrong number of arguments for meta function productValueInList");
2958 if (arguments.size() != 1)
2959 B2FATAL(
"Wrong number of arguments for meta function angleToClosestInList");
2961 std::string listname = arguments[0];
2963 auto func = [listname](
const Particle * particle) ->
double {
2965 StoreObjPtr<ParticleList> list(listname);
2966 if (not list.isValid())
2967 B2FATAL(
"Invalid particle list name " << listname <<
" given to angleToClosestInList");
2970 if (list->getListSize() == 0)
2975 const auto p_this = frame.getMomentum(particle);
2978 double minAngle = 2 * M_PI;
2979 for (
unsigned int i = 0; i < list->getListSize(); ++i)
2981 const Particle* compareme = list->getParticle(i);
2982 const auto p_compare = frame.getMomentum(compareme);
2983 double angle = ROOT::Math::VectorUtil::Angle(p_compare, p_this);
2984 if (minAngle > angle) minAngle = angle;
2994 if (arguments.size() != 2)
2995 B2FATAL(
"Wrong number of arguments for meta function closestInList");
2997 std::string listname = arguments[0];
3002 auto func = [listname, var](
const Particle * particle) ->
double {
3004 StoreObjPtr<ParticleList> list(listname);
3005 if (not list.isValid())
3006 B2FATAL(
"Invalid particle list name " << listname <<
" given to closestInList");
3010 const auto p_this = frame.getMomentum(particle);
3013 double minAngle = 2 * M_PI;
3015 for (
unsigned int i = 0; i < list->getListSize(); ++i)
3017 const Particle* compareme = list->getParticle(i);
3018 const auto p_compare = frame.getMomentum(compareme);
3019 double angle = ROOT::Math::VectorUtil::Angle(p_compare, p_this);
3020 if (minAngle > angle) {
3028 auto var_result = var->function(list->getParticle(iClosest));
3029 if (std::holds_alternative<double>(var_result))
3031 return std::get<double>(var_result);
3032 }
else if (std::holds_alternative<int>(var_result))
3034 return std::get<int>(var_result);
3035 }
else if (std::holds_alternative<bool>(var_result))
3037 return std::get<bool>(var_result);
3046 if (arguments.size() != 1)
3047 B2FATAL(
"Wrong number of arguments for meta function angleToMostB2BInList");
3049 std::string listname = arguments[0];
3051 auto func = [listname](
const Particle * particle) ->
double {
3053 StoreObjPtr<ParticleList> list(listname);
3054 if (not list.isValid())
3055 B2FATAL(
"Invalid particle list name " << listname <<
" given to angleToMostB2BInList");
3058 if (list->getListSize() == 0)
3063 const auto p_this = frame.getMomentum(particle);
3067 double maxAngle = 0;
3068 for (
unsigned int i = 0; i < list->getListSize(); ++i)
3070 const Particle* compareme = list->getParticle(i);
3071 const auto p_compare = frame.getMomentum(compareme);
3072 double angle = ROOT::Math::VectorUtil::Angle(p_compare, p_this);
3073 if (maxAngle < angle) maxAngle = angle;
3083 if (arguments.size() != 1)
3084 B2FATAL(
"Wrong number of arguments for meta function deltaPhiToMostB2BPhiInList");
3086 std::string listname = arguments[0];
3088 auto func = [listname](
const Particle * particle) ->
double {
3090 StoreObjPtr<ParticleList> list(listname);
3091 if (not list.isValid())
3092 B2FATAL(
"Invalid particle list name " << listname <<
" given to deltaPhiToMostB2BPhiInList");
3095 if (list->getListSize() == 0)
3100 const auto phi_this = frame.getMomentum(particle).Phi();
3103 double maxAngle = 0;
3104 for (
unsigned int i = 0; i < list->getListSize(); ++i)
3106 const Particle* compareme = list->getParticle(i);
3107 const auto phi_compare = frame.getMomentum(compareme).Phi();
3108 double angle = std::abs(phi_compare - phi_this);
3109 if (angle > M_PI) {angle = 2 * M_PI - angle;}
3110 if (maxAngle < angle) maxAngle = angle;
3120 if (arguments.size() != 2)
3121 B2FATAL(
"Wrong number of arguments for meta function mostB2BInList");
3123 std::string listname = arguments[0];
3128 auto func = [listname, var](
const Particle * particle) ->
double {
3130 StoreObjPtr<ParticleList> list(listname);
3131 if (not list.isValid())
3132 B2FATAL(
"Invalid particle list name " << listname <<
" given to mostB2BInList");
3136 const auto p_this = frame.getMomentum(particle);
3140 double maxAngle = -1.0;
3142 for (
unsigned int i = 0; i < list->getListSize(); ++i)
3144 const Particle* compareme = list->getParticle(i);
3145 const auto p_compare = frame.getMomentum(compareme);
3146 double angle = ROOT::Math::VectorUtil::Angle(p_compare, p_this);
3147 if (maxAngle < angle) {
3155 auto var_result = var->function(list->getParticle(iMostB2B));
3156 if (std::holds_alternative<double>(var_result))
3158 return std::get<double>(var_result);
3159 }
else if (std::holds_alternative<int>(var_result))
3161 return std::get<int>(var_result);
3162 }
else if (std::holds_alternative<bool>(var_result))
3164 return std::get<bool>(var_result);
3172 if (arguments.size() == 1) {
3173 std::string listName = arguments[0];
3174 auto func = [listName](
const Particle*) ->
double {
3175 StoreObjPtr<ParticleList> listOfParticles(listName);
3177 if (!(listOfParticles.isValid())) B2FATAL(
"Invalid Listname " << listName <<
" given to maxOpeningAngleInList");
3178 int nParticles = listOfParticles->getListSize();
3183 double maxOpeningAngle = -1;
3184 for (
int i = 0; i < nParticles; i++)
3186 ROOT::Math::PxPyPzEVector v1 = frame.getMomentum(listOfParticles->getParticle(i));
3187 for (
int j = i + 1; j < nParticles; j++) {
3188 ROOT::Math::PxPyPzEVector v2 = frame.getMomentum(listOfParticles->getParticle(j));
3189 const double angle = ROOT::Math::VectorUtil::Angle(v1, v2);
3190 if (angle > maxOpeningAngle) maxOpeningAngle = angle;
3193 return maxOpeningAngle;
3197 B2FATAL(
"Wrong number of arguments for meta function maxOpeningAngleInList");
3204 if (arguments.size() >= 2) {
3209 auto func = [var, arguments](
const Particle * particle) ->
double {
3210 if (particle ==
nullptr)
3212 B2WARNING(
"Trying to access a daughter that does not exist. Skipping");
3218 ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
3222 for (
unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++)
3224 auto generalizedIndex = arguments[iCoord];
3225 const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
3227 pSum += frame.getMomentum(dauPart);
3229 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << generalizedIndex);
3235 Particle sumOfDaughters(pSum, 100);
3237 auto var_result = var->function(&sumOfDaughters);
3239 if (std::holds_alternative<double>(var_result))
3241 return std::get<double>(var_result);
3242 }
else if (std::holds_alternative<int>(var_result))
3244 return std::get<int>(var_result);
3245 }
else if (std::holds_alternative<bool>(var_result))
3247 return std::get<bool>(var_result);
3252 B2FATAL(
"Wrong number of arguments for meta function daughterCombination");
3255 Manager::FunctionPtr useAlternativeDaughterHypothesis(
const std::vector<std::string>& arguments)
3268 if (arguments.size() >= 2) {
3279 std::unordered_map<unsigned int, int> mapOfReplacedDaughters;
3282 for (
unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++) {
3283 auto replacedDauString = arguments[iCoord];
3285 std::vector<std::string> indexAndMass;
3286 boost::split(indexAndMass, replacedDauString, boost::is_any_of(
":"));
3289 if (indexAndMass.size() > 2) {
3290 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 "
3291 << replacedDauString <<
", while a correct syntax looks like 0:K+.");
3295 if (indexAndMass.size() < 2) {
3296 B2WARNING(
"The string indicating which daughter's mass should be replaced contains only one colon-separated element instead of two. The offending string is "
3297 << replacedDauString <<
", while a correct syntax looks like 0:K+.");
3305 }
catch (std::invalid_argument&) {
3306 B2FATAL(
"Found the string " << indexAndMass[0] <<
"instead of a daughter index.");
3310 TParticlePDG* particlePDG = TDatabasePDG::Instance()->GetParticle(indexAndMass[1].c_str());
3312 B2WARNING(
"Particle not in evt.pdl file! " << indexAndMass[1]);
3317 int pdgCode = particlePDG->PdgCode();
3318 mapOfReplacedDaughters[dauIndex] = pdgCode;
3322 if (mapOfReplacedDaughters.size() != arguments.size() - 1)
3323 B2FATAL(
"Overlapped daughter's index is detected in the meta-variable useAlternativeDaughterHypothesis");
3331 auto func = [var, mapOfReplacedDaughters](
const Particle * particle) ->
double {
3332 if (particle ==
nullptr)
3334 B2WARNING(
"Trying to access a particle that does not exist. Skipping");
3344 ROOT::Math::PxPyPzMVector pSum(0, 0, 0, 0);
3346 for (
unsigned int iDau = 0; iDau < particle->getNDaughters(); iDau++)
3348 const Particle* dauPart = particle->getDaughter(iDau);
3350 B2WARNING(
"Trying to access a daughter that does not exist. Index = " << iDau);
3354 ROOT::Math::PxPyPzMVector dauMom = ROOT::Math::PxPyPzMVector(frame.getMomentum(dauPart));
3358 pdgCode = mapOfReplacedDaughters.at(iDau);
3359 }
catch (std::out_of_range&) {
3366 double p_x = dauMom.Px();
3367 double p_y = dauMom.Py();
3368 double p_z = dauMom.Pz();
3369 dauMom.SetCoordinates(p_x, p_y, p_z, TDatabasePDG::Instance()->GetParticle(pdgCode)->Mass());
3370 const_cast<Particle*
>(dummy->getDaughter(iDau))->set4VectorDividingByMomentumScaling(ROOT::Math::PxPyPzEVector(dauMom));
3373 const int charge = dummy->getDaughter(iDau)->getCharge();
3374 if (TDatabasePDG::Instance()->GetParticle(pdgCode)->Charge() / 3.0 == charge)
3375 const_cast<Particle*
>(dummy->getDaughter(iDau))->setPDGCode(pdgCode);
3377 const_cast<Particle*
>(dummy->getDaughter(iDau))->setPDGCode(-1 * pdgCode);
3383 dummy->set4Vector(ROOT::Math::PxPyPzEVector(pSum));
3385 auto var_result = var->function(dummy);
3388 if (std::holds_alternative<double>(var_result))
3390 return std::get<double>(var_result);
3391 }
else if (std::holds_alternative<int>(var_result))
3393 return std::get<int>(var_result);
3394 }
else if (std::holds_alternative<bool>(var_result))
3396 return std::get<bool>(var_result);
3402 B2FATAL(
"Wrong number of arguments for meta function useAlternativeDaughterHypothesis");
3407 if (arguments.size() == 2) {
3409 std::string arg = arguments[0];
3411 TParticlePDG* part = TDatabasePDG::Instance()->GetParticle(arg.c_str());
3413 if (part !=
nullptr) {
3414 pdg_code = std::abs(part->PdgCode());
3418 }
catch (std::exception& e) {}
3421 if (pdg_code == -1) {
3422 B2FATAL(
"Ancestor " + arg +
" is not recognised. Please provide valid PDG code or particle name.");
3425 auto func = [pdg_code, var](
const Particle * particle) ->
double {
3426 const Particle* p = particle;
3428 int ancestor_level = std::get<double>(
Manager::Instance().getVariable(
"hasAncestor(" + std::to_string(pdg_code) +
", 0)")->function(p));
3429 if ((ancestor_level <= 0) or (std::isnan(ancestor_level)))
3434 const MCParticle* i_p = p->getMCParticle();
3436 for (
int a = 0; a < ancestor_level ; a = a + 1)
3438 i_p = i_p->getMother();
3441 StoreArray<Particle> tempParticles(
"tempParticles");
3442 tempParticles.clear();
3444 Particle* newPart = tempParticles.appendNew(m_p);
3445 newPart->addRelationTo(i_p);
3447 appendDaughtersRecursive(newPart, tempParticles);
3449 auto var_result = var->function(newPart);
3450 if (std::holds_alternative<double>(var_result))
3452 return std::get<double>(var_result);
3453 }
else if (std::holds_alternative<int>(var_result))
3455 return std::get<int>(var_result);
3456 }
else if (std::holds_alternative<bool>(var_result))
3458 return std::get<bool>(var_result);
3463 B2FATAL(
"Wrong number of arguments for meta function varForFirstMCAncestorOfType (expected 2: type and variable of interest)");
3469 if (arguments.size() != 1) {
3470 B2FATAL(
"Number of arguments for nTrackFitResults must be 1, particleType or PDGcode");
3473 std::string arg = arguments[0];
3474 TDatabasePDG* pdgDatabase = TDatabasePDG::Instance();
3475 TParticlePDG* part = pdgDatabase->GetParticle(arg.c_str());
3477 if (part !=
nullptr) {
3478 absPdg = std::abs(part->PdgCode());
3482 }
catch (
const std::exception&) {
3486 if (absPdg == 0 || pdgDatabase->GetParticle(absPdg) ==
nullptr) {
3487 B2FATAL(
"nTrackFitResults: argument '" << arg <<
"' is neither a valid particle name nor a PDG code");
3491 auto func = [absPdg](
const Particle*) ->
int {
3493 Const::ChargedStable type(absPdg);
3494 StoreArray<Track> tracks;
3496 int nTrackFitResults = 0;
3498 for (
const auto& track : tracks)
3500 const TrackFitResult* trackFit = track.getTrackFitResultWithClosestMass(type);
3502 if (!trackFit)
continue;
3503 if (trackFit->getChargeSign() == 0)
continue;
3508 return nTrackFitResults;
3517 if (arguments.size() == 2) {
3520 auto func = [var, default_val](
const Particle * particle) ->
int {
3521 auto var_result = var->function(particle);
3522 if (std::holds_alternative<double>(var_result))
3524 double value = std::get<double>(var_result);
3525 if (value > std::numeric_limits<int>::max())
3526 value = std::numeric_limits<int>::max();
3527 if (value < std::numeric_limits<int>::min())
3528 value = std::numeric_limits<int>::min();
3529 if (std::isnan(value))
3530 value = default_val;
3531 return static_cast<int>(value);
3532 }
else if (std::holds_alternative<int>(var_result))
3533 return std::get<int>(var_result);
3534 else if (std::holds_alternative<bool>(var_result))
3535 return static_cast<int>(std::get<bool>(var_result));
3536 else return default_val;
3540 B2FATAL(
"Wrong number of arguments for meta function int, please provide variable name and replacement value for NaN!");
3544 VARIABLE_GROUP(
"MetaFunctions");
3545 REGISTER_METAVARIABLE(
"nCleanedECLClusters(cut)", nCleanedECLClusters,
3546 "[Eventbased] Returns the number of clean Clusters in the event\n"
3547 "Clean clusters are defined by the clusters which pass the given cut assuming a photon hypothesis.",
3548 Manager::VariableDataType::c_int);
3549 REGISTER_METAVARIABLE(
"nCleanedTracks(cut)", nCleanedTracks,
3550 "[Eventbased] Returns the number of clean Tracks in the event\n"
3551 "Clean tracks are defined by the tracks which pass the given cut assuming a pion hypothesis.", Manager::VariableDataType::c_int);
3552 REGISTER_METAVARIABLE(
"formula(v1 + v2 * [v3 - v4] / v5^v6)", formula, R
"DOCSTRING(
3553Returns the result of the given formula, where v1 to vN are variables or floating
3554point numbers. Currently the only supported operations are addition (``+``),
3555subtraction (``-``), multiplication (``*``), division (``/``) and power (``^``
3556or ``**``). Parenthesis can be in the form of square brackets ``[v1 * v2]``
3557or normal brackets ``(v1 * v2)``. It will work also with variables taking
3558arguments. Operator precedence is taken into account. For example ::
3560 (daughter(0, E) + daughter(1, E))**2 - p**2 + 0.138
3562.. versionchanged:: release-03-00-00
3563 now both, ``[]`` and ``()`` can be used for grouping operations, ``**`` can
3564 be used for exponent and float literals are possible directly in the
3566)DOCSTRING", Manager::VariableDataType::c_double);
3567 REGISTER_METAVARIABLE("useRestFrame(variable)", useRestFrame,
3568 "Returns the value of the variable using the rest frame of the given particle as current reference frame.\n"
3569 "E.g. ``useRestFrame(daughter(0, p))`` returns the total momentum of the first daughter in its mother's rest-frame", Manager::VariableDataType::c_double);
3570 REGISTER_METAVARIABLE(
"useCMSFrame(variable)", useCMSFrame,
3571 "Returns the value of the variable using the CMS frame as current reference frame.\n"
3572 "E.g. ``useCMSFrame(E)`` returns the energy of a particle in the CMS frame.", Manager::VariableDataType::c_double);
3573 REGISTER_METAVARIABLE(
"useLabFrame(variable)", useLabFrame, R
"DOC(
3574Returns the value of ``variable`` in the *lab* frame.
3577 The lab frame is the default reference frame, usually you don't need to use this meta-variable.
3578 E.g. ``useLabFrame(E)`` returns the energy of a particle in the Lab frame, same as just ``E``.
3580Specifying the lab frame is useful in some corner-cases. For example:
3581``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.
3582)DOC", Manager::VariableDataType::c_double);
3583 REGISTER_METAVARIABLE("useTagSideRecoilRestFrame(variable, daughterIndexTagB)", useTagSideRecoilRestFrame,
3584 "Returns the value of the variable in the rest frame of the recoiling particle to the tag side B meson.\n"
3585 "The variable should only be applied to an Upsilon(4S) list.\n"
3586 "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);
3587 REGISTER_METAVARIABLE(
"useParticleRestFrame(variable, particleList)", useParticleRestFrame,
3588 "Returns the value of the variable in the rest frame of the first Particle contained in the given ParticleList.\n"
3589 "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3590 "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3591 "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);
3592 REGISTER_METAVARIABLE(
"useRecoilParticleRestFrame(variable, particleList)", useRecoilParticleRestFrame,
3593 "Returns the value of the variable in the rest frame of recoil system against the first Particle contained in the given ParticleList.\n"
3594 "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3595 "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3596 "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);
3597 REGISTER_METAVARIABLE(
"useDaughterRestFrame(variable, daughterIndex_1, [daughterIndex_2, ... daughterIndex_3])", useDaughterRestFrame,
3598 "Returns the value of the variable in the rest frame of the selected daughter particle.\n"
3599 "The daughter is identified via generalized daughter index, e.g. ``0:1`` identifies the second daughter (1) "
3600 "of the first daughter (0). If the daughter index is invalid, it returns NaN.\n"
3601 "If two or more indices are given, the rest frame of the sum of the daughters is used.",
3602 Manager::VariableDataType::c_double);
3603 REGISTER_METAVARIABLE(
"useDaughterRecoilRestFrame(variable, daughterIndex_1, [daughterIndex_2, ... daughterIndex_3])", useDaughterRecoilRestFrame,
3604 "Returns the value of the variable in the rest frame of the recoil of the selected daughter particle.\n"
3605 "The daughter is identified via generalized daughter index, e.g. ``0:1`` identifies the second daughter (1) "
3606 "of the first daughter (0). If the daughter index is invalid, it returns NaN.\n"
3607 "If two or more indices are given, the rest frame of the sum of the daughters is used.",
3608 Manager::VariableDataType::c_double);
3609 REGISTER_METAVARIABLE(
"useMCancestorBRestFrame(variable)", useMCancestorBRestFrame,
3610 "Returns the value of the variable in the rest frame of the ancestor B MC particle.\n"
3611 "If no B or no MC-matching is found, it returns NaN.", Manager::VariableDataType::c_double);
3612 REGISTER_METAVARIABLE(
"passesCut(cut)", passesCut,
3613 "Returns 1 if particle passes the cut otherwise 0.\n"
3614 "Useful if you want to write out if a particle would have passed a cut or not.", Manager::VariableDataType::c_bool);
3615 REGISTER_METAVARIABLE(
"passesEventCut(cut)", passesEventCut,
3616 "[Eventbased] Returns 1 if event passes the cut otherwise 0.\n"
3617 "Useful if you want to select events passing a cut without looping into particles, such as for skimming.\n", Manager::VariableDataType::c_bool);
3618 REGISTER_METAVARIABLE(
"countDaughters(cut)", countDaughters,
3619 "Returns number of direct daughters which satisfy the cut.\n"
3620 "Used by the skimming package (for what exactly?)", Manager::VariableDataType::c_int);
3621 REGISTER_METAVARIABLE(
"countFSPDaughters(cut)", countDescendants,
3622 "Returns number of final-state daughters which satisfy the cut.",
3623 Manager::VariableDataType::c_int);
3624 REGISTER_METAVARIABLE(
"countDescendants(cut)", countDescendants,
3625 "Returns number of descendants for all generations which satisfy the cut.",
3626 Manager::VariableDataType::c_int);
3627 REGISTER_METAVARIABLE(
"varFor(pdgCode, variable)", varFor,
3628 "Returns the value of the variable for the given particle if its abs(pdgCode) agrees with the given one.\n"
3629 "E.g. ``varFor(11, p)`` returns the momentum if the particle is an electron or a positron.", Manager::VariableDataType::c_double);
3630 REGISTER_METAVARIABLE(
"varForMCGen(variable)", varForMCGen,
3631 "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"
3632 "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"
3633 "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);
3634 REGISTER_METAVARIABLE(
"nParticlesInList(particleListName)", nParticlesInList,
3635 "[Eventbased] Returns number of particles in the given particle List.", Manager::VariableDataType::c_int);
3636 REGISTER_METAVARIABLE(
"isInList(particleListName)", isInList,
3637 "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);
3638 REGISTER_METAVARIABLE(
"isDaughterOfList(particleListNames)", isDaughterOfList,
3639 "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);
3640 REGISTER_METAVARIABLE(
"isDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isDescendantOfList, R
"DOC(
3641 Returns 1 if the given particle appears in the decay chain of the particles in the given ParticleLists.
3643 Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3645 * ``isDescendantOfList(<particle_list>,1)`` returns 1 if particle is a daughter of the list,
3646 * ``isDescendantOfList(<particle_list>,2)`` returns 1 if particle is a granddaughter of the list,
3647 * ``isDescendantOfList(<particle_list>,3)`` returns 1 if particle is a great-granddaughter of the list, etc.
3648 * Default value is ``-1`` that is inclusive for all generations.
3649 )DOC", Manager::VariableDataType::c_bool);
3650 REGISTER_METAVARIABLE("isMCDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isMCDescendantOfList, R
"DOC(
3651 Returns 1 if the given particle is linked to the same MC particle as any reconstructed daughter of the decay lists.
3653 Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3655 * ``isMCDescendantOfList(<particle_list>,1)`` returns 1 if particle is matched to the same particle as any daughter of the list,
3656 * ``isMCDescendantOfList(<particle_list>,2)`` returns 1 if particle is matched to the same particle as any granddaughter of the list,
3657 * ``isMCDescendantOfList(<particle_list>,3)`` returns 1 if particle is matched to the same particle as any great-granddaughter of the list, etc.
3658 * Default value is ``-1`` that is inclusive for all generations.
3660 It makes only sense for lists created with `fillParticleListFromMC` function with ``addDaughters=True`` argument.
3661 )DOC", Manager::VariableDataType::c_bool);
3663 REGISTER_METAVARIABLE("sourceObjectIsInList(particleListName)", sourceObjectIsInList, R
"DOC(
3664Returns 1 if the underlying mdst object (e.g. track, or cluster) was used to create a particle in ``particleListName``, 0 if not.
3667 This only makes sense for particles that are not composite. Returns -1 for composite particles.
3668)DOC", Manager::VariableDataType::c_int);
3670 REGISTER_METAVARIABLE("mcParticleIsInMCList(particleListName)", mcParticleIsInMCList, R
"DOC(
3671Returns 1 if the particle's matched MC particle is also matched to a particle in ``particleListName``
3672(or if either of the lists were filled from generator level `modularAnalysis.fillParticleListFromMC`.)
3674.. seealso:: :b2:var:`isMCDescendantOfList` to check daughters.
3675)DOC", Manager::VariableDataType::c_bool);
3677 REGISTER_METAVARIABLE("isGrandDaughterOfList(particleListNames)", isGrandDaughterOfList,
3678 "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);
3679 REGISTER_METAVARIABLE(
"originalParticle(variable)", originalParticle, R
"DOC(
3680 Returns value of variable for the original particle from which the given particle is copied.
3682 The copy of particle is created, for example, when the vertex fit updates the daughters and `modularAnalysis.copyParticles` is called.
3683 Returns NaN if the given particle is not copied and so there is no original particle.
3684 )DOC", Manager::VariableDataType::c_double);
3685 REGISTER_METAVARIABLE("daughter(i, variable)", daughter, R
"DOC(
3686 Returns value of variable for the i-th daughter. E.g.
3688 * ``daughter(0, p)`` returns the total momentum of the first daughter.
3689 * ``daughter(0, daughter(1, p)`` returns the total momentum of the second daughter of the first daughter.
3691 Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3692 )DOC", Manager::VariableDataType::c_double);
3693 REGISTER_METAVARIABLE("originalDaughter(i, variable)", originalDaughter, R
"DOC(
3694 Returns value of variable for the original particle from which the i-th daughter is copied.
3696 The copy of particle is created, for example, when the vertex fit updates the daughters and `modularAnalysis.copyParticles` is called.
3697 Returns NaN if the daughter is not copied and so there is no original daughter.
3699 Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3700 )DOC", Manager::VariableDataType::c_double);
3701 REGISTER_METAVARIABLE("mcDaughter(i, variable)", mcDaughter, R
"DOC(
3702 Returns the value of the requested variable for the i-th Monte Carlo daughter of the particle.
3704 Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3705 or if the i-th MC daughter does not exist.
3707 E.g. ``mcDaughter(0, PDG)`` will return the PDG code of the first MC daughter of the matched MC
3708 particle of the reconstructed particle the function is applied to.
3710 The meta variable can also be nested: ``mcDaughter(0, mcDaughter(1, PDG))``.
3711 )DOC", Manager::VariableDataType::c_double);
3712 REGISTER_METAVARIABLE("mcMother(variable)", mcMother, R
"DOC(
3713 Returns the value of the requested variable for the Monte Carlo mother of the particle.
3715 Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3716 or if the MC mother does not exist.
3718 E.g. ``mcMother(PDG)`` will return the PDG code of the MC mother of the matched MC
3719 particle of the reconstructed particle the function is applied to.
3721 The meta variable can also be nested: ``mcMother(mcMother(PDG))``.
3722 )DOC", Manager::VariableDataType::c_double);
3723 REGISTER_METAVARIABLE("genParticle(index, variable)", genParticle, R
"DOC(
3724[Eventbased] Returns the ``variable`` for the ith generator particle.
3725The arguments of the function must be the ``index`` of the particle in the MCParticle Array,
3726and ``variable``, the name of the function or variable for that generator particle.
3727If ``index`` goes beyond the length of the MCParticles array, NaN will be returned.
3729E.g. ``genParticle(0, p)`` returns the total momentum of the first MCParticle, which in a generic decay up to MC15 is
3730the Upsilon(4S) and for MC16 and beyond the initial electron.
3731)DOC", Manager::VariableDataType::c_double);
3732 REGISTER_METAVARIABLE("genUpsilon4S(variable)", genUpsilon4S, R
"DOC(
3733[Eventbased] Returns the ``variable`` evaluated for the generator-level :math:`\Upsilon(4S)`.
3734If no generator level :math:`\Upsilon(4S)` exists for the event, NaN will be returned.
3736E.g. ``genUpsilon4S(p)`` returns the total momentum of the :math:`\Upsilon(4S)` in a generic decay.
3737``genUpsilon4S(mcDaughter(1, p))`` returns the total momentum of the second daughter of the
3738generator-level :math:`\Upsilon(4S)` (i.e. the momentum of the second B meson in a generic decay).
3739)DOC", Manager::VariableDataType::c_double);
3740 REGISTER_METAVARIABLE("daughterProductOf(variable)", daughterProductOf,
3741 "Returns product of a variable over all daughters.\n"
3742 "E.g. ``daughterProductOf(extraInfo(SignalProbability))`` returns the product of the SignalProbabilitys of all daughters.", Manager::VariableDataType::c_double);
3743 REGISTER_METAVARIABLE(
"daughterSumOf(variable)", daughterSumOf,
3744 "Returns sum of a variable over all daughters.\n"
3745 "E.g. ``daughterSumOf(nDaughters)`` returns the number of grand-daughters.", Manager::VariableDataType::c_double);
3746 REGISTER_METAVARIABLE(
"daughterLowest(variable)", daughterLowest,
3747 "Returns the lowest value of the given variable among all daughters.\n"
3748 "E.g. ``useCMSFrame(daughterLowest(p))`` returns the lowest momentum in CMS frame.", Manager::VariableDataType::c_double);
3749 REGISTER_METAVARIABLE(
"daughterHighest(variable)", daughterHighest,
3750 "Returns the highest value of the given variable among all daughters.\n"
3751 "E.g. ``useCMSFrame(daughterHighest(p))`` returns the highest momentum in CMS frame.", Manager::VariableDataType::c_double);
3752 REGISTER_METAVARIABLE(
"daughterDiffOf(daughterIndex_i, daughterIndex_j, variable)", daughterDiffOf, R
"DOC(
3753 Returns the difference of a variable between the two given daughters.
3754 E.g. ``useRestFrame(daughterDiffOf(0, 1, p))`` returns the momentum difference between first and second daughter in the rest frame of the given particle.
3755 (That means that it returns :math:`p_j - p_i`)
3757 The daughters can be provided as generalized daughter indexes, which are simply colon-separated
3758 lists of daughter indexes, ordered starting from the root particle. For example, ``0:1``
3759 identifies the second daughter (1) of the first daughter (0) of the mother particle.
3761 )DOC", Manager::VariableDataType::c_double);
3762 REGISTER_METAVARIABLE("mcDaughterDiffOf(i, j, variable)", mcDaughterDiffOf,
3763 "MC matched version of the `daughterDiffOf` function.", Manager::VariableDataType::c_double);
3764 REGISTER_METAVARIABLE(
"grandDaughterDiffOf(i, j, variable)", grandDaughterDiffOf,
3765 "Returns the difference of a variable between the first daughters of the two given daughters.\n"
3766 "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"
3767 "(That means that it returns :math:`p_j - p_i`)", Manager::VariableDataType::c_double);
3768 MAKE_DEPRECATED(
"grandDaughterDiffOf",
false,
"light-2402-ocicat", R
"DOC(
3769 The difference between any combination of (grand-)daughters can be calculated with the more general variable :b2:var:`daughterDiffOf`
3770 by using generalized daughter indexes.)DOC");
3771 REGISTER_METAVARIABLE("daughterNormDiffOf(i, j, variable)", daughterNormDiffOf,
3772 "Returns the normalized difference of a variable between the two given daughters.\n"
3773 "E.g. ``daughterNormDiffOf(0, 1, p)`` returns the normalized momentum difference between first and second daughter in the lab frame.", Manager::VariableDataType::c_double);
3774 REGISTER_METAVARIABLE(
"daughterMotherDiffOf(i, variable)", daughterMotherDiffOf,
3775 "Returns the difference of a variable between the given daughter and the mother particle itself.\n"
3776 "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);
3777 REGISTER_METAVARIABLE(
"daughterMotherNormDiffOf(i, variable)", daughterMotherNormDiffOf,
3778 "Returns the normalized difference of a variable between the given daughter and the mother particle itself.\n"
3779 "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);
3780 REGISTER_METAVARIABLE(
"angleBetweenDaughterAndRecoil(daughterIndex_1, daughterIndex_2, ... )", angleBetweenDaughterAndRecoil, R
"DOC(
3781 Returns the angle between the momentum recoiling against the particle and the sum of the momenta of the given daughters.
3782 The unit of the angle is ``rad``.
3784 The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3785 daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3786 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3787 identifies the second daughter of the root particle.
3789 At least one generalized index has to be given to ``angleBetweenDaughterAndRecoil``.
3792 ``angleBetweenDaughterAndRecoil(0)`` will return the angle between pRecoil and the momentum of the first daughter.
3794 ``angleBetweenDaughterAndRecoil(0, 1)`` will return the angle between pRecoil and the sum of the momenta of the first and second daughter.
3796 ``angleBetweenDaughterAndRecoil(0:0, 3:0)`` will return the angle between pRecoil and the sum of the momenta of the: first daughter of the first daughter, and
3797 the first daughter of the fourth daughter.)DOC", Manager::VariableDataType::c_double);
3798 REGISTER_METAVARIABLE("angleBetweenDaughterAndMissingMomentum(daughterIndex_1, daughterIndex_2, ... )", angleBetweenDaughterAndMissingMomentum, R
"DOC(
3799 Returns the angle between the missing momentum in the event and the sum of the momenta of the given daughters.
3800 The unit of the angle is ``rad``. EventKinematics module has to be called to use this.
3802 The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3803 daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3804 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3805 identifies the second daughter of the root particle.
3807 At least one generalized index has to be given to ``angleBetweenDaughterAndMissingMomentum``.
3810 ``angleBetweenDaughterAndMissingMomentum(0)`` will return the angle between missMom and the momentum of the first daughter.
3812 ``angleBetweenDaughterAndMissingMomentum(0, 1)`` will return the angle between missMom and the sum of the momenta of the first and second daughter.
3814 ``angleBetweenDaughterAndMissingMomentum(0:0, 3:0)`` will return the angle between missMom and the sum of the momenta of the: first daughter of the first daughter, and
3815 the first daughter of the fourth daughter.)DOC", Manager::VariableDataType::c_double);
3816 REGISTER_METAVARIABLE("daughterAngle(daughterIndex_1, daughterIndex_2[, daughterIndex_3])", daughterAngle, R
"DOC(
3817 Returns the angle in between any pair of particles belonging to the same decay tree.
3818 The unit of the angle is ``rad``.
3820 The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3821 daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3822 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3823 identifies the second daughter of the root particle.
3825 Both two and three generalized indexes can be given to ``daughterAngle``. If two indices are given, the
3826 variable returns the angle between the momenta of the two given particles. If three indices are given, the
3827 variable returns the angle between the momentum of the third particle and a vector which is the sum of the
3828 first two daughter momenta.
3831 ``daughterAngle(0, 3)`` will return the angle between the first and fourth daughter.
3832 ``daughterAngle(0, 1, 3)`` will return the angle between the fourth daughter and the sum of the first and
3834 ``daughterAngle(0:0, 3:0)`` will return the angle between the first daughter of the first daughter, and
3835 the first daughter of the fourth daughter.
3837 )DOC", Manager::VariableDataType::c_double);
3838 REGISTER_METAVARIABLE("mcDaughterAngle(daughterIndex_1, daughterIndex_2, [daughterIndex_3])", mcDaughterAngle,
3839 "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);
3840 REGISTER_VARIABLE(
"grandDaughterDecayAngle(i, j)", grandDaughterDecayAngle,
3841 "Returns the decay angle of the granddaughter in the daughter particle's rest frame.\n"
3842 "It is calculated with respect to the reverted momentum vector of the particle.\n"
3843 "Two arguments representing the daughter and granddaughter indices have to be provided as arguments.\n\n",
"rad");
3844 REGISTER_VARIABLE(
"daughterClusterAngleInBetween(i, j)", daughterClusterAngleInBetween,
3845 "Returns the angle between clusters associated to the two daughters."
3846 "If two indices given: returns the angle between the momenta of the clusters associated to the two given daughters."
3847 "If three indices given: returns the angle between the momentum of the third particle's cluster and a vector "
3848 "which is the sum of the first two daughter's cluster momenta."
3849 "Returns nan if any of the daughters specified don't have an associated cluster."
3850 "The arguments in the argument vector must be integers corresponding to the ith and jth (and kth) daughters.\n\n",
"rad");
3851 REGISTER_METAVARIABLE(
"daughterInvM(i[, j, ...])", daughterInvM, R
"DOC(
3852 Returns the invariant mass adding the Lorentz vectors of the given daughters. The unit of the invariant mass is GeV/:math:`\text{c}^2`
3853 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.
3855 Daughters from different generations of the decay tree can be combined using generalized daughter indexes,
3856 which are simply colon-separated daughter indexes for each generation, starting from the root particle. For
3857 example, ``0:1:3`` identifies the fourth daughter (3) of the second daughter (1) of the first daughter(0) of
3858 the mother particle.
3860 Returns NaN if the given daughter-index is out of bound (>= number of daughters))DOC", Manager::VariableDataType::c_double);
3861 REGISTER_METAVARIABLE("extraInfo(name)", extraInfo,
3862 "Returns extra info stored under the given name.\n"
3863 "The extraInfo has to be set by a module first.\n"
3864 "E.g. ``extraInfo(SignalProbability)`` returns the SignalProbability calculated by the ``MVAExpert`` module.\n"
3865 "If nothing is set under the given name or if the particle is a nullptr, NaN is returned.\n"
3866 "In the latter case please use `eventExtraInfo` if you want to access an EventExtraInfo variable.", Manager::VariableDataType::c_double);
3867 REGISTER_METAVARIABLE(
"eventExtraInfo(name)", eventExtraInfo,
3868 "[Eventbased] Returns extra info stored under the given name in the event extra info.\n"
3869 "The extraInfo has to be set first by another module like MVAExpert in event mode.\n"
3870 "If nothing is set under this name, NaN is returned.", Manager::VariableDataType::c_double);
3871 REGISTER_METAVARIABLE(
"eventCached(variable)", eventCached,
3872 "[Eventbased] Returns value of event-based variable and caches this value in the EventExtraInfo.\n"
3873 "The result of second call to this variable in the same event will be provided from the cache.\n"
3874 "It is recommended to use this variable in order to declare custom aliases as event-based. This is "
3875 "necessary if using the eventwise mode of variablesToNtuple).", Manager::VariableDataType::c_double);
3876 REGISTER_METAVARIABLE(
"particleCached(variable)", particleCached,
3877 "Returns value of given variable and caches this value in the ParticleExtraInfo of the provided particle.\n"
3878 "The result of second call to this variable on the same particle will be provided from the cache.", Manager::VariableDataType::c_double);
3879 REGISTER_METAVARIABLE(
"modulo(variable, n)", modulo,
3880 "Returns rest of division of variable by n.", Manager::VariableDataType::c_int);
3881 REGISTER_METAVARIABLE(
"abs(variable)", abs,
3882 "Returns absolute value of the given variable.\n"
3883 "E.g. abs(mcPDG) returns the absolute value of the mcPDG, which is often useful for cuts.", Manager::VariableDataType::c_double);
3884 REGISTER_METAVARIABLE(
"max(var1,var2)", max,
"Returns max value of two variables.\n", Manager::VariableDataType::c_double);
3885 REGISTER_METAVARIABLE(
"min(var1,var2)", min,
"Returns min value of two variables.\n", Manager::VariableDataType::c_double);
3886 REGISTER_METAVARIABLE(
"sin(variable)", sin,
"Returns sine value of the given variable.", Manager::VariableDataType::c_double);
3887 REGISTER_METAVARIABLE(
"asin(variable)", asin,
"Returns arcsine of the given variable. The unit of the asin() is ``rad``", Manager::VariableDataType::c_double);
3888 REGISTER_METAVARIABLE(
"cos(variable)", cos,
"Returns cosine value of the given variable.", Manager::VariableDataType::c_double);
3889 REGISTER_METAVARIABLE(
"acos(variable)", acos,
"Returns arccosine value of the given variable. The unit of the acos() is ``rad``", Manager::VariableDataType::c_double);
3890 REGISTER_METAVARIABLE(
"tan(variable)", tan,
"Returns tangent value of the given variable.", Manager::VariableDataType::c_double);
3891 REGISTER_METAVARIABLE(
"atan(variable)", atan,
"Returns arctangent value of the given variable. The unit of the atan() is ``rad``", Manager::VariableDataType::c_double);
3892 REGISTER_METAVARIABLE(
"exp(variable)", exp,
"Returns exponential evaluated for the given variable.", Manager::VariableDataType::c_double);
3893 REGISTER_METAVARIABLE(
"log(variable)", log,
"Returns natural logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3894 REGISTER_METAVARIABLE(
"log10(variable)", log10,
"Returns base-10 logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3895 REGISTER_METAVARIABLE(
"int(variable, nan_replacement)", convertToInt, R
"DOC(
3896 Casts the output of the variable to an integer value.
3899 Overflow and underflow are clipped at maximum and minimum values, respectively. NaN values are replaced with the value of the 2nd argument.
3901 )DOC", Manager::VariableDataType::c_int);
3902 REGISTER_METAVARIABLE("isNAN(variable)", isNAN,
3903 "Returns true if variable value evaluates to nan (determined via std::isnan(double)).\n"
3904 "Useful for debugging.", Manager::VariableDataType::c_bool);
3905 REGISTER_METAVARIABLE(
"ifNANgiveX(variable, x)", ifNANgiveX,
3906 "Returns x (has to be a number) if variable value is nan (determined via std::isnan(double)).\n"
3907 "Useful for technical purposes while training MVAs.", Manager::VariableDataType::c_double);
3908 REGISTER_METAVARIABLE(
"isInfinity(variable)", isInfinity,
3909 "Returns true if variable value evaluates to infinity (determined via std::isinf(double)).\n"
3910 "Useful for debugging.", Manager::VariableDataType::c_bool);
3911 REGISTER_METAVARIABLE(
"unmask(variable, flag1, flag2, ...)", unmask,
3912 "unmask(variable, flag1, flag2, ...) or unmask(variable, mask) sets certain bits in the variable to zero.\n"
3913 "For example, if you want to set the second, fourth and fifth bits to zero, you could call \n"
3914 "``unmask(variable, 2, 8, 16)`` or ``unmask(variable, 26)``.\n"
3915 "", Manager::VariableDataType::c_double);
3916 REGISTER_METAVARIABLE(
"conditionalVariableSelector(cut, variableIfTrue, variableIfFalse)", conditionalVariableSelector,
3917 "Returns one of the two supplied variables, depending on whether the particle passes the supplied cut.\n"
3918 "The first variable is returned if the particle passes the cut, and the second variable is returned otherwise.", Manager::VariableDataType::c_double);
3919 REGISTER_METAVARIABLE(
"pValueCombination(p1, p2, ...)", pValueCombination,
3920 "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"
3921 "If any of the p-values is invalid, i.e. smaller than zero, -1 is returned.", Manager::VariableDataType::c_double);
3922 REGISTER_METAVARIABLE(
"pValueCombinationOfDaughters(variable)", pValueCombinationOfDaughters,
3923 "Returns the combined p-value of the daughter p-values according to the formula given in `Nucl. Instr. and Meth. A 411 (1998) 449 <https://doi.org/10.1016/S0168-9002(98)00293-9>`_ .\n"
3924 "If any of the p-values is invalid, i.e. smaller than zero, -1 is returned.", Manager::VariableDataType::c_double);
3925 REGISTER_METAVARIABLE(
"veto(particleList, cut, pdgCode = 11)", veto,
3926 "Combines current particle with particles from the given particle list and returns 1 if the combination passes the provided cut. \n"
3927 "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"
3928 "around the neutral Pion mass (e.g. ``0.130 < M < 0.140``). \n"
3929 "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);
3930 REGISTER_METAVARIABLE(
"matchedMC(variable)", matchedMC,
3931 "Returns variable output for the matched MCParticle by constructing a temporary Particle from it.\n"
3932 "This may not work too well if your variable requires accessing daughters of the particle.\n"
3933 "E.g. ``matchedMC(p)`` returns the total momentum of the related MCParticle.\n"
3934 "Returns NaN if no matched MCParticle exists.", Manager::VariableDataType::c_double);
3935 REGISTER_METAVARIABLE(
"clusterBestMatchedMCParticle(variable)", clusterBestMatchedMCParticle,
3936 "Returns variable output for the MCParticle that is best-matched with the ECLCluster of the given Particle.\n"
3937 "E.g. To get the energy of the MCParticle that matches best with an ECLCluster, one could use ``clusterBestMatchedMCParticle(E)``\n"
3938 "When the variable is called for ``gamma`` and if the ``gamma`` is matched with MCParticle, it works same as `matchedMC`.\n"
3939 "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"
3940 "Returns NaN if the particle is not matched to an ECLCluster, or if the ECLCluster has no matching MCParticles", Manager::VariableDataType::c_double);
3941 REGISTER_METAVARIABLE(
"varForBestMatchedMCKlong(variable)", clusterBestMatchedMCKlong,
3942 "Returns variable output for the Klong MCParticle which has the best match with the ECLCluster of the given Particle.\n"
3943 "Returns NaN if the particle is not matched to an ECLCluster, or if the ECLCluster has no matching Klong MCParticle", Manager::VariableDataType::c_double);
3945 REGISTER_METAVARIABLE(
"countInList(particleList, cut='')", countInList,
"[Eventbased] "
3946 "Returns number of particle which pass given in cut in the specified particle list.\n"
3947 "Useful for creating statistics about the number of particles in a list.\n"
3948 "E.g. ``countInList(e+, isSignal == 1)`` returns the number of correctly reconstructed electrons in the event.\n"
3949 "The variable is event-based and does not need a valid particle pointer as input.", Manager::VariableDataType::c_int);
3950 REGISTER_METAVARIABLE(
"getVariableByRank(particleList, rankedVariableName, variableName, rank)", getVariableByRank, R
"DOC(
3951 [Eventbased] Returns the value of ``variableName`` for the candidate in the ``particleList`` with the requested ``rank``.
3954 The `BestCandidateSelection` module available via `rankByHighest` / `rankByLowest` has to be used before.
3957 The first candidate matching the given rank is used.
3958 Thus, it is not recommended to use this variable in conjunction with ``allowMultiRank`` in the `BestCandidateSelection` module.
3960 The suffix ``_rank`` is automatically added to the argument ``rankedVariableName``,
3961 which either has to be the name of the variable used to order the candidates or the selected outputVariable name without the ending ``_rank``.
3962 This means that your selected name for the rank variable has to end with ``_rank``.
3964 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>`_
3965 )DOC", Manager::VariableDataType::c_double);
3966 REGISTER_VARIABLE("matchedMCHasPDG(PDGCode)", matchedMCHasPDG,
3967 "Returns if the absolute value of the PDGCode of the MCParticle related to the Particle matches a given PDGCode."
3968 "Returns 0/NAN/1 if PDGCode does not match/is not available/ matches");
3969 REGISTER_METAVARIABLE(
"numberOfNonOverlappingParticles(pList1, pList2, ...)", numberOfNonOverlappingParticles,
3970 "Returns the number of non-overlapping particles in the given particle lists"
3971 "Useful to check if there is additional physics going on in the detector if one reconstructed the Y4S", Manager::VariableDataType::c_int);
3972 REGISTER_METAVARIABLE(
"totalEnergyOfParticlesInList(particleListName)", totalEnergyOfParticlesInList,
3973 "[Eventbased] Returns the total energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3974 REGISTER_METAVARIABLE(
"totalPxOfParticlesInList(particleListName)", totalPxOfParticlesInList,
3975 "[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);
3976 REGISTER_METAVARIABLE(
"totalPyOfParticlesInList(particleListName)", totalPyOfParticlesInList,
3977 "[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);
3978 REGISTER_METAVARIABLE(
"totalPzOfParticlesInList(particleListName)", totalPzOfParticlesInList,
3979 "[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);
3980 REGISTER_METAVARIABLE(
"invMassInLists(pList1, pList2, ...)", invMassInLists,
3981 "[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);
3982 REGISTER_METAVARIABLE(
"totalECLEnergyOfParticlesInList(particleListName)", totalECLEnergyOfParticlesInList,
3983 "[Eventbased] Returns the total ECL energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3984 REGISTER_METAVARIABLE(
"maxPtInList(particleListName)", maxPtInList,
3985 "[Eventbased] Returns maximum transverse momentum Pt in the given particle List. The unit of the transverse momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3986 REGISTER_METAVARIABLE(
"eclClusterSpecialTrackMatched(cut)", eclClusterTrackMatchedWithCondition,
3987 "Returns if at least one Track that satisfies the given condition is related to the ECLCluster of the Particle.", Manager::VariableDataType::c_double);
3988 REGISTER_METAVARIABLE(
"averageValueInList(particleListName, variable)", averageValueInList,
3989 "[Eventbased] Returns the arithmetic mean of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3990 REGISTER_METAVARIABLE(
"medianValueInList(particleListName, variable)", medianValueInList,
3991 "[Eventbased] Returns the median value of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3992 REGISTER_METAVARIABLE(
"sumValueInList(particleListName, variable)", sumValueInList,
3993 "[Eventbased] Returns the sum of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3994 REGISTER_METAVARIABLE(
"productValueInList(particleListName, variable)", productValueInList,
3995 "[Eventbased] Returns the product of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3996 REGISTER_METAVARIABLE(
"angleToClosestInList(particleListName)", angleToClosestInList,
3997 "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);
3998 REGISTER_METAVARIABLE(
"closestInList(particleListName, variable)", closestInList,
3999 "Returns `variable` for the closest particle (smallest opening angle) in the list provided.", Manager::VariableDataType::c_double);
4000 REGISTER_METAVARIABLE(
"angleToMostB2BInList(particleListName)", angleToMostB2BInList,
4001 "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);
4002 REGISTER_METAVARIABLE(
"deltaPhiToMostB2BPhiInList(particleListName)", deltaPhiToMostB2BPhiInList,
4003 "Returns the abs(delta phi) between this particle and the most back-to-back particle in phi (closest opening angle to 180) in the list provided. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
4004 REGISTER_METAVARIABLE(
"mostB2BInList(particleListName, variable)", mostB2BInList,
4005 "Returns `variable` for the most back-to-back particle (closest opening angle to 180) in the list provided.", Manager::VariableDataType::c_double);
4006 REGISTER_METAVARIABLE(
"maxOpeningAngleInList(particleListName)", maxOpeningAngleInList,
4007 "[Eventbased] Returns maximum opening angle in the given particle List. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
4008 REGISTER_METAVARIABLE(
"daughterCombination(variable, daughterIndex_1, daughterIndex_2 ... daughterIndex_n)", daughterCombination,R
"DOC(
4009Returns a ``variable`` function only of the 4-momentum calculated on an arbitrary set of (grand)daughters.
4012 ``variable`` can only be a function of the daughters' 4-momenta.
4014Daughters from different generations of the decay tree can be combined using generalized daughter indexes, which are simply colon-separated
4015the list of daughter indexes, starting from the root particle: for example, ``0:1:3`` identifies the fourth
4016daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle.
4019 ``daughterCombination(M, 0, 3, 4)`` will return the invariant mass of the system made of the first, fourth and fifth daughter of particle.
4020 ``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.
4022)DOC", Manager::VariableDataType::c_double);
4023 REGISTER_METAVARIABLE("useAlternativeDaughterHypothesis(variable, daughterIndex_1:newMassHyp_1, ..., daughterIndex_n:newMassHyp_n)", useAlternativeDaughterHypothesis,R
"DOC(
4024Returns a ``variable`` calculated using new mass hypotheses for (some of) the particle's daughters.
4027 ``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.
4028 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.
4029 Also, the track fit is not performed again: the variable only re-calculates the 4-vectors using different mass assumptions.
4030 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).
4033 Generalized daughter indexes are not supported (yet!): this variable can be used only on first-generation daughters.
4036 ``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.
4037 ``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.
4039)DOC", Manager::VariableDataType::c_double);
4040 REGISTER_METAVARIABLE("varForFirstMCAncestorOfType(type, variable)",varForFirstMCAncestorOfType,R
"DOC(Returns requested variable of the first ancestor of the given type.
4041Ancestor type can be set up by PDG code or by particle name (check evt.pdl for valid particle names))DOC", Manager::VariableDataType::c_double);
4043 REGISTER_METAVARIABLE("nTrackFitResults(particleType)", nTrackFitResults,
4044 "[Eventbased] Returns the total number of TrackFitResults for a given particleType. The argument can be the name of particle (e.g. pi+) or PDG code (e.g. 211).",
4045 Manager::VariableDataType::c_int);
4047 REGISTER_METAVARIABLE(
"convertToDaughterIndex(variable)", convertToDaughterIndex, R
"DOC(Converts the variable of the given particle into integer and returns it if it is a valid daughter index, else returns -1.)DOC", Manager::VariableDataType::c_int);
int getPDGCode() const
PDG code.
static const ChargedStable pion
charged pion particle
static const double doubleNaN
quiet_NaN
static const ChargedStable electron
electron particle
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
@ c_nPhotons
CR is split into n photons (N1)
static std::unique_ptr< GeneralCut > compile(const std::string &cut)
@ c_Initial
bit 5: Particle is initial such as e+ or e- and not going to Geant4
@ c_PrimaryParticle
bit 0: Particle is primary particle.
@ c_IsVirtual
bit 4: Particle is virtual and not going to Geant4.
static std::string makeROOTCompatible(std::string str)
Remove special characters that ROOT dislikes in branch names, e.g.
EParticleSourceObject
particle source enumerators
@ c_Flavored
Is either particle or antiparticle.
static const ReferenceFrame & GetCurrent()
Get current rest frame.
std::function< VarVariant(const Particle *)> FunctionPtr
functions stored take a const Particle* and return VarVariant.
const Var * getVariable(std::string name)
Get the variable belonging to the given key.
std::variant< double, int, bool > VarVariant
NOTE: the python interface is documented manually in analysis/doc/Variables.rst (because we use ROOT ...
static Manager & Instance()
get singleton instance.
#define MAKE_DEPRECATED(name, make_fatal, version, description)
Registers a variable as deprecated.
T convertString(const std::string &str)
Converts a string to type T (one of float, double, long double, int, long int, unsigned long int).
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.