Belle II Software  light-2212-foldex
MetaVariables.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 // Own include
10 #include <analysis/variables/MetaVariables.h>
11 
12 #include <analysis/VariableManager/Utility.h>
13 #include <framework/dataobjects/EventExtraInfo.h>
14 #include <analysis/dataobjects/Particle.h>
15 #include <analysis/dataobjects/ParticleList.h>
16 #include <analysis/dataobjects/RestOfEvent.h>
17 #include <analysis/utility/PCmsLabTransform.h>
18 #include <analysis/utility/ReferenceFrame.h>
19 #include <analysis/utility/EvtPDLUtil.h>
20 #include <analysis/utility/ParticleCopy.h>
21 #include <analysis/ClusterUtility/ClusterUtils.h>
22 #include <analysis/variables/VariableFormulaConstructor.h>
23 
24 #include <framework/logging/Logger.h>
25 #include <framework/datastore/StoreArray.h>
26 #include <framework/datastore/StoreObjPtr.h>
27 #include <framework/utilities/Conversion.h>
28 #include <framework/utilities/MakeROOTCompatible.h>
29 #include <framework/gearbox/Const.h>
30 
31 #include <mdst/dataobjects/Track.h>
32 #include <mdst/dataobjects/MCParticle.h>
33 #include <mdst/dataobjects/ECLCluster.h>
34 #include <mdst/dataobjects/TrackFitResult.h>
35 
36 #include <boost/algorithm/string.hpp>
37 #include <limits>
38 
39 #include <cmath>
40 #include <stdexcept>
41 #include <memory>
42 #include <regex>
43 
44 #include <TDatabasePDG.h>
45 #include <Math/Vector4D.h>
46 
47 namespace Belle2 {
52  namespace Variable {
53 
54  Manager::FunctionPtr useRestFrame(const std::vector<std::string>& arguments)
55  {
56  if (arguments.size() == 1) {
57  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
58  auto func = [var](const Particle * particle) -> double {
59  UseReferenceFrame<RestFrame> frame(particle);
60  double result = std::get<double>(var->function(particle));
61  return result;
62  };
63  return func;
64  } else {
65  B2FATAL("Wrong number of arguments for meta function useRestFrame");
66  }
67  }
68 
69  Manager::FunctionPtr useCMSFrame(const std::vector<std::string>& arguments)
70  {
71  if (arguments.size() == 1) {
72  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
73  auto func = [var](const Particle * particle) -> double {
74  UseReferenceFrame<CMSFrame> frame;
75  double result = std::get<double>(var->function(particle));
76  return result;
77  };
78  return func;
79  } else {
80  B2FATAL("Wrong number of arguments for meta function useCMSFrame");
81  }
82  }
83 
84  Manager::FunctionPtr useLabFrame(const std::vector<std::string>& arguments)
85  {
86  if (arguments.size() == 1) {
87  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
88  auto func = [var](const Particle * particle) -> double {
89  UseReferenceFrame<LabFrame> frame;
90  double result = std::get<double>(var->function(particle));
91  return result;
92  };
93  return func;
94  } else {
95  B2FATAL("Wrong number of arguments for meta function useLabFrame");
96  }
97  }
98 
99  Manager::FunctionPtr useTagSideRecoilRestFrame(const std::vector<std::string>& arguments)
100  {
101  if (arguments.size() == 2) {
102  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
103 
104  int daughterIndexTagB = 0;
105  try {
106  daughterIndexTagB = Belle2::convertString<int>(arguments[1]);
107  } catch (std::invalid_argument&) {
108  B2FATAL("Second argument of useTagSideRecoilRestFrame meta function must be integer!");
109  }
110 
111  auto func = [var, daughterIndexTagB](const Particle * particle) -> double {
112  if (particle->getPDGCode() != 300553)
113  {
114  B2ERROR("Variable should only be used on a Upsilon(4S) Particle List!");
115  return std::numeric_limits<float>::quiet_NaN();
116  }
117 
118  PCmsLabTransform T;
119  ROOT::Math::PxPyPzEVector pSigB = T.getBeamFourMomentum() - particle->getDaughter(daughterIndexTagB)->get4Vector();
120  Particle tmp(pSigB, -particle->getDaughter(daughterIndexTagB)->getPDGCode());
121 
122  UseReferenceFrame<RestFrame> frame(&tmp);
123  double result = std::get<double>(var->function(particle));
124  return result;
125  };
126 
127  return func;
128  } else {
129  B2FATAL("Wrong number of arguments for meta function useTagSideRecoilRestFrame");
130  }
131  }
132 
133  Manager::FunctionPtr useParticleRestFrame(const std::vector<std::string>& arguments)
134  {
135  if (arguments.size() == 2) {
136  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
137  std::string listName = arguments[1];
138  auto func = [var, listName](const Particle * particle) -> double {
139  StoreObjPtr<ParticleList> list(listName);
140  unsigned listSize = list->getListSize();
141  if (listSize == 0)
142  return std::numeric_limits<float>::quiet_NaN();
143  if (listSize > 1)
144  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."
145  << LogVar("ParticleList", listName)
146  << LogVar("Number of candidates in the list", listSize));
147  const Particle* p = list->getParticle(0);
148  UseReferenceFrame<RestFrame> frame(p);
149  double result = std::get<double>(var->function(particle));
150  return result;
151  };
152  return func;
153  } else {
154  B2FATAL("Wrong number of arguments for meta function useParticleRestFrame.");
155  }
156  }
157 
158  Manager::FunctionPtr useRecoilParticleRestFrame(const std::vector<std::string>& arguments)
159  {
160  if (arguments.size() == 2) {
161  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
162  std::string listName = arguments[1];
163  auto func = [var, listName](const Particle * particle) -> double {
164  StoreObjPtr<ParticleList> list(listName);
165  unsigned listSize = list->getListSize();
166  if (listSize == 0)
167  return std::numeric_limits<float>::quiet_NaN();
168  if (listSize > 1)
169  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."
170  << LogVar("ParticleList", listName)
171  << LogVar("Number of candidates in the list", listSize));
172  const Particle* p = list->getParticle(0);
173  PCmsLabTransform T;
174  ROOT::Math::PxPyPzEVector recoil = T.getBeamFourMomentum() - p->get4Vector();
175  /* Let's use 0 as PDG code to avoid wrong assumptions. */
176  Particle pRecoil(recoil, 0);
177  pRecoil.setVertex(particle->getVertex());
178  UseReferenceFrame<RestFrame> frame(&pRecoil);
179  double result = std::get<double>(var->function(particle));
180  return result;
181  };
182  return func;
183  } else {
184  B2FATAL("Wrong number of arguments for meta function useParticleRestFrame.");
185  }
186  }
187 
188  Manager::FunctionPtr useDaughterRestFrame(const std::vector<std::string>& arguments)
189  {
190  if (arguments.size() >= 2) {
191  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
192  auto func = [var, arguments](const Particle * particle) -> double {
193 
194  // Sum of the 4-momenta of all the selected daughters
195  ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
196 
197  for (unsigned int i = 1; i < arguments.size(); i++)
198  {
199  auto generalizedIndex = arguments[i];
200  const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
201  if (dauPart)
202  pSum += dauPart->get4Vector();
203  else
204  return std::numeric_limits<float>::quiet_NaN();
205  }
206  Particle tmp(pSum, 0);
207  UseReferenceFrame<RestFrame> frame(&tmp);
208  double result = std::get<double>(var->function(particle));
209  return result;
210  };
211  return func;
212  } else {
213  B2FATAL("Wrong number of arguments for meta function useDaughterRestFrame.");
214  }
215  }
216 
217  Manager::FunctionPtr extraInfo(const std::vector<std::string>& arguments)
218  {
219  if (arguments.size() == 1) {
220  auto extraInfoName = arguments[0];
221  auto func = [extraInfoName](const Particle * particle) -> double {
222  if (particle == nullptr)
223  {
224  B2WARNING("Returns NaN because the particle is nullptr! If you want EventExtraInfo variables, please use eventExtraInfo() instead");
225  return std::numeric_limits<float>::quiet_NaN();
226  }
227  if (particle->hasExtraInfo(extraInfoName))
228  {
229  return particle->getExtraInfo(extraInfoName);
230  } else
231  {
232  return std::numeric_limits<float>::quiet_NaN();
233  }
234  };
235  return func;
236  } else {
237  B2FATAL("Wrong number of arguments for meta function extraInfo");
238  }
239  }
240 
241  Manager::FunctionPtr eventExtraInfo(const std::vector<std::string>& arguments)
242  {
243  if (arguments.size() == 1) {
244  auto extraInfoName = arguments[0];
245  auto func = [extraInfoName](const Particle*) -> double {
246  StoreObjPtr<EventExtraInfo> eventExtraInfo;
247  if (not eventExtraInfo.isValid())
248  return std::numeric_limits<float>::quiet_NaN();
249  if (eventExtraInfo->hasExtraInfo(extraInfoName))
250  {
251  return eventExtraInfo->getExtraInfo(extraInfoName);
252  } else
253  {
254  return std::numeric_limits<float>::quiet_NaN();
255  }
256  };
257  return func;
258  } else {
259  B2FATAL("Wrong number of arguments for meta function extraInfo");
260  }
261  }
262 
263  Manager::FunctionPtr eventCached(const std::vector<std::string>& arguments)
264  {
265  if (arguments.size() == 1) {
266  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
267  std::string key = std::string("__") + MakeROOTCompatible::makeROOTCompatible(var->name);
268  auto func = [var, key](const Particle*) -> double {
269 
270  StoreObjPtr<EventExtraInfo> eventExtraInfo;
271  if (not eventExtraInfo.isValid())
272  eventExtraInfo.create();
273  if (eventExtraInfo->hasExtraInfo(key))
274  {
275  return eventExtraInfo->getExtraInfo(key);
276  } else
277  {
278  double value = std::numeric_limits<double>::quiet_NaN();
279  auto var_result = var->function(nullptr);
280  if (std::holds_alternative<double>(var_result)) {
281  value = std::get<double>(var_result);
282  } else if (std::holds_alternative<int>(var_result)) {
283  return std::get<int>(var_result);
284  } else if (std::holds_alternative<bool>(var_result)) {
285  return std::get<bool>(var_result);
286  }
287  eventExtraInfo->addExtraInfo(key, value);
288  return value;
289  }
290  };
291  return func;
292  } else {
293  B2FATAL("Wrong number of arguments for meta function eventCached");
294  }
295  }
296 
297  Manager::FunctionPtr particleCached(const std::vector<std::string>& arguments)
298  {
299  if (arguments.size() == 1) {
300  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
301  std::string key = std::string("__") + MakeROOTCompatible::makeROOTCompatible(var->name);
302  auto func = [var, key](const Particle * particle) -> double {
303 
304  if (particle->hasExtraInfo(key))
305  {
306  return particle->getExtraInfo(key);
307  } else
308  {
309  double value = std::get<double>(var->function(particle));
310  // Remove constness from Particle pointer.
311  // The extra-info is used as a cache in our case,
312  // indicated by the double-underscore in front of the key.
313  // One could implement the cache as a separate property of the particle object
314  // and mark it as mutable, however, this would only lead to code duplication
315  // and an increased size of the particle object.
316  // Thus, we decided to use the extra-info field and cast away the const in this case.
317  const_cast<Particle*>(particle)->addExtraInfo(key, value);
318  return value;
319  }
320  };
321  return func;
322  } else {
323  B2FATAL("Wrong number of arguments for meta function particleCached");
324  }
325  }
326 
327  // Formula of other variables, going to require a space between all operators and operations.
328  // Later can add some check for : (colon) trailing + or - to distinguish between particle lists
329  // and operations, but for now cbf.
330  Manager::FunctionPtr formula(const std::vector<std::string>& arguments)
331  {
332  if (arguments.size() != 1) B2FATAL("Wrong number of arguments for meta function formula");
333  FormulaParser<VariableFormulaConstructor> parser;
334  try {
335  return parser.parse(arguments[0]);
336  } catch (std::runtime_error& e) {
337  B2FATAL(e.what());
338  }
339  }
340 
341  Manager::FunctionPtr nCleanedTracks(const std::vector<std::string>& arguments)
342  {
343  if (arguments.size() <= 1) {
344 
345  std::string cutString;
346  if (arguments.size() == 1)
347  cutString = arguments[0];
348  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
349  auto func = [cut](const Particle*) -> int {
350 
351  int number_of_tracks = 0;
352  StoreArray<Track> tracks;
353  for (const auto& track : tracks)
354  {
355  const TrackFitResult* trackFit = track.getTrackFitResultWithClosestMass(Const::pion);
356  if (trackFit->getChargeSign() == 0) {
357  // Ignore track
358  } else {
359  Particle particle(&track, Const::pion);
360  if (cut->check(&particle))
361  number_of_tracks++;
362  }
363  }
364 
365  return number_of_tracks;
366 
367  };
368  return func;
369  } else {
370  B2FATAL("Wrong number of arguments for meta function nCleanedTracks");
371  }
372  }
373 
374  Manager::FunctionPtr nCleanedECLClusters(const std::vector<std::string>& arguments)
375  {
376  if (arguments.size() <= 1) {
377 
378  std::string cutString;
379  if (arguments.size() == 1)
380  cutString = arguments[0];
381  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
382  auto func = [cut](const Particle*) -> int {
383 
384  int number_of_clusters = 0;
385  StoreArray<ECLCluster> clusters;
386  for (const auto& cluster : clusters)
387  {
388  // look only at momentum of N1 (n photons) ECLClusters
389  if (!cluster.hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons))
390  continue;
391 
392  Particle particle(&cluster);
393  if (cut->check(&particle))
394  number_of_clusters++;
395  }
396 
397  return number_of_clusters;
398 
399  };
400  return func;
401  } else {
402  B2FATAL("Wrong number of arguments for meta function nCleanedECLClusters");
403  }
404  }
405 
406  Manager::FunctionPtr passesCut(const std::vector<std::string>& arguments)
407  {
408  if (arguments.size() == 1) {
409  std::string cutString = arguments[0];
410  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
411  auto func = [cut](const Particle * particle) -> bool {
412  if (cut->check(particle))
413  return 1;
414  else
415  return 0;
416  };
417  return func;
418  } else {
419  B2FATAL("Wrong number of arguments for meta function passesCut");
420  }
421  }
422 
423  Manager::FunctionPtr passesEventCut(const std::vector<std::string>& arguments)
424  {
425  if (arguments.size() == 1) {
426  std::string cutString = arguments[0];
427  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
428  auto func = [cut](const Particle*) -> bool {
429  if (cut->check(nullptr))
430  return 1;
431  else
432  return 0;
433  };
434  return func;
435  } else {
436  B2FATAL("Wrong number of arguments for meta function passesEventCut");
437  }
438  }
439 
440  Manager::FunctionPtr varFor(const std::vector<std::string>& arguments)
441  {
442  if (arguments.size() == 2) {
443  int pdgCode = 0;
444  try {
445  pdgCode = Belle2::convertString<int>(arguments[0]);
446  } catch (std::invalid_argument&) {
447  B2FATAL("The first argument of varFor meta function must be a positive integer!");
448  }
449  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
450  auto func = [pdgCode, var](const Particle * particle) -> double {
451  if (std::abs(particle->getPDGCode()) == std::abs(pdgCode))
452  {
453  auto var_result = var->function(particle);
454  if (std::holds_alternative<double>(var_result)) {
455  return std::get<double>(var_result);
456  } else if (std::holds_alternative<int>(var_result)) {
457  return std::get<int>(var_result);
458  } else if (std::holds_alternative<bool>(var_result)) {
459  return std::get<bool>(var_result);
460  } else return std::numeric_limits<double>::quiet_NaN();
461  } else return std::numeric_limits<double>::quiet_NaN();
462  };
463  return func;
464  } else {
465  B2FATAL("Wrong number of arguments for meta function varFor");
466  }
467  }
468 
469  Manager::FunctionPtr varForMCGen(const std::vector<std::string>& arguments)
470  {
471  if (arguments.size() == 1) {
472  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
473  auto func = [var](const Particle * particle) -> double {
474  if (particle->getMCParticle())
475  {
476  if (particle->getMCParticle()->getStatus(MCParticle::c_PrimaryParticle)
477  && (! particle->getMCParticle()->getStatus(MCParticle::c_IsVirtual))
478  && (! particle->getMCParticle()->getStatus(MCParticle::c_Initial))) {
479  auto var_result = var->function(particle);
480  if (std::holds_alternative<double>(var_result)) {
481  return std::get<double>(var_result);
482  } else if (std::holds_alternative<int>(var_result)) {
483  return std::get<int>(var_result);
484  } else if (std::holds_alternative<bool>(var_result)) {
485  return std::get<bool>(var_result);
486  } else return std::numeric_limits<double>::quiet_NaN();
487  } else return std::numeric_limits<double>::quiet_NaN();
488  } else return std::numeric_limits<double>::quiet_NaN();
489  };
490  return func;
491  } else {
492  B2FATAL("Wrong number of arguments for meta function varForMCGen");
493  }
494  }
495 
496  Manager::FunctionPtr nParticlesInList(const std::vector<std::string>& arguments)
497  {
498  if (arguments.size() == 1) {
499  std::string listName = arguments[0];
500  auto func = [listName](const Particle * particle) -> int {
501 
502  (void) particle;
503  StoreObjPtr<ParticleList> listOfParticles(listName);
504 
505  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << listName << " given to nParticlesInList");
506 
507  return listOfParticles->getListSize();
508 
509  };
510  return func;
511  } else {
512  B2FATAL("Wrong number of arguments for meta function nParticlesInList");
513  }
514  }
515 
516  Manager::FunctionPtr isInList(const std::vector<std::string>& arguments)
517  {
518  // unpack arguments, there should be only one: the name of the list we're checking
519  if (arguments.size() != 1) {
520  B2FATAL("Wrong number of arguments for isInList");
521  }
522  auto listName = arguments[0];
523 
524  auto func = [listName](const Particle * particle) -> bool {
525 
526  // check the list exists
527  StoreObjPtr<ParticleList> list(listName);
528  if (!(list.isValid()))
529  {
530  B2FATAL("Invalid Listname " << listName << " given to isInList");
531  }
532 
533  // is the particle in the list?
534  return list->contains(particle);
535 
536  };
537  return func;
538  }
539 
540  Manager::FunctionPtr sourceObjectIsInList(const std::vector<std::string>& arguments)
541  {
542  // unpack arguments, there should be only one: the name of the list we're checking
543  if (arguments.size() != 1) {
544  B2FATAL("Wrong number of arguments for sourceObjectIsInList");
545  }
546  auto listName = arguments[0];
547 
548  auto func = [listName](const Particle * particle) -> int {
549 
550  // check the list exists
551  StoreObjPtr<ParticleList> list(listName);
552  if (!(list.isValid()))
553  {
554  B2FATAL("Invalid Listname " << listName << " given to sourceObjectIsInList");
555  }
556 
557  // this only makes sense for particles that are *not* composite and come
558  // from some mdst object (tracks, clusters..)
559  Particle::EParticleSourceObject particlesource = particle->getParticleSource();
560  if (particlesource == Particle::EParticleSourceObject::c_Composite
561  or particlesource == Particle::EParticleSourceObject::c_Undefined)
562  return -1;
563 
564  // it *is* possible to have a particle list from different sources (like
565  // hadrons from the ECL and KLM) so we have to check each particle in
566  // the list individually
567  for (unsigned i = 0; i < list->getListSize(); ++i)
568  {
569  Particle* iparticle = list->getParticle(i);
570  if (particlesource == iparticle->getParticleSource())
571  if (particle->getMdstArrayIndex() == iparticle->getMdstArrayIndex())
572  return 1;
573  }
574  return 0;
575 
576  };
577  return func;
578  }
579 
580  Manager::FunctionPtr mcParticleIsInMCList(const std::vector<std::string>& arguments)
581  {
582  // unpack arguments, there should be only one: the name of the list we're checking
583  if (arguments.size() != 1) {
584  B2FATAL("Wrong number of arguments for mcParticleIsInMCList");
585  }
586  auto listName = arguments[0];
587 
588  auto func = [listName](const Particle * particle) -> bool {
589 
590  // check the list exists
591  StoreObjPtr<ParticleList> list(listName);
592  if (!(list.isValid()))
593  B2FATAL("Invalid Listname " << listName << " given to mcParticleIsInMCList");
594 
595  // this can only be true for mc-matched particles or particles are created from MCParticles
596  const MCParticle* mcp = particle->getMCParticle();
597  if (mcp == nullptr) return false;
598 
599  // check every particle in the input list is not matched to (or created from) the same MCParticle
600  for (unsigned i = 0; i < list->getListSize(); ++i)
601  {
602  const MCParticle* imcp = list->getParticle(i)->getMCParticle();
603  if ((imcp != nullptr) and (mcp->getArrayIndex() == imcp->getArrayIndex()))
604  return true;
605  }
606  return false;
607  };
608  return func;
609  }
610 
611  Manager::FunctionPtr isDaughterOfList(const std::vector<std::string>& arguments)
612  {
613  B2WARNING("isDaughterOfList is outdated and replaced by isDescendantOfList.");
614  std::vector<std::string> new_arguments = arguments;
615  new_arguments.push_back(std::string("1"));
616  return isDescendantOfList(new_arguments);
617  }
618 
619  Manager::FunctionPtr isGrandDaughterOfList(const std::vector<std::string>& arguments)
620  {
621  B2WARNING("isGrandDaughterOfList is outdated and replaced by isDescendantOfList.");
622  std::vector<std::string> new_arguments = arguments;
623  new_arguments.push_back(std::string("2"));
624  return isDescendantOfList(new_arguments);
625  }
626 
627  Manager::FunctionPtr isDescendantOfList(const std::vector<std::string>& arguments)
628  {
629  if (arguments.size() > 0) {
630  auto listNames = arguments;
631  auto func = [listNames](const Particle * particle) -> bool {
632  bool output = false;
633  int generation_flag = -1;
634  try
635  {
636  generation_flag = Belle2::convertString<int>(listNames.back());
637  } catch (std::exception& e) {}
638 
639  for (auto& iListName : listNames)
640  {
641  try {
642  Belle2::convertString<int>(iListName);
643  continue;
644  } catch (std::exception& e) {}
645 
646  // Creating recursive lambda
647  auto list_comparison = [](auto&& self, const Particle * m, const Particle * p, int flag)-> bool {
648  bool result = false;
649  for (unsigned i = 0; i < m->getNDaughters(); ++i)
650  {
651  const Particle* daughter = m->getDaughter(i);
652  if ((flag == 1.) or (flag < 0)) {
653  if (p->isCopyOf(daughter)) {
654  return true;
655  }
656  }
657 
658  if (flag != 1.) {
659  if (daughter->getNDaughters() > 0) {
660  result = self(self, daughter, p, flag - 1);
661  if (result == 1) {
662  return true;
663  }
664  }
665  }
666  }
667  return result;
668  };
669 
670  StoreObjPtr<ParticleList> listOfParticles(iListName);
671 
672  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << iListName << " given to isDescendantOfList");
673 
674  for (unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
675  Particle* iParticle = listOfParticles->getParticle(i);
676  output = list_comparison(list_comparison, iParticle, particle, generation_flag);
677  if (output) {
678  return output;
679  }
680  }
681  }
682  return output;
683  };
684  return func;
685  } else {
686  B2FATAL("Wrong number of arguments for meta function isDescendantOfList");
687  }
688  }
689 
690  Manager::FunctionPtr isMCDescendantOfList(const std::vector<std::string>& arguments)
691  {
692  if (arguments.size() > 0) {
693  auto listNames = arguments;
694  auto func = [listNames](const Particle * particle) -> bool {
695  bool output = false;
696  int generation_flag = -1;
697  try
698  {
699  generation_flag = Belle2::convertString<int>(listNames.back());
700  } catch (std::exception& e) {}
701 
702  if (particle->getMCParticle() == nullptr)
703  {
704  return false;
705  }
706 
707  for (auto& iListName : listNames)
708  {
709  try {
710  std::stod(iListName);
711  continue;
712  } catch (std::exception& e) {}
713  // Creating recursive lambda
714  auto list_comparison = [](auto&& self, const Particle * m, const Particle * p, int flag)-> bool {
715  bool result = false;
716  for (unsigned i = 0; i < m->getNDaughters(); ++i)
717  {
718  const Particle* daughter = m->getDaughter(i);
719  if ((flag == 1.) or (flag < 0)) {
720  if (daughter->getMCParticle() != nullptr) {
721  if (p->getMCParticle()->getArrayIndex() == daughter->getMCParticle()->getArrayIndex()) {
722  return true;
723  }
724  }
725  }
726  if (flag != 1.) {
727  if (daughter->getNDaughters() > 0) {
728  result = self(self, daughter, p, flag - 1);
729  if (result) {
730  return true;
731  }
732  }
733  }
734  }
735  return result;
736  };
737 
738  StoreObjPtr<ParticleList> listOfParticles(iListName);
739 
740  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << iListName << " given to isMCDescendantOfList");
741 
742  for (unsigned i = 0; i < listOfParticles->getListSize(); ++i) {
743  Particle* iParticle = listOfParticles->getParticle(i);
744  output = list_comparison(list_comparison, iParticle, particle, generation_flag);
745  if (output) {
746  return output;
747  }
748  }
749  }
750  return output;
751  };
752  return func;
753  } else {
754  B2FATAL("Wrong number of arguments for meta function isMCDescendantOfList");
755  }
756  }
757 
758  Manager::FunctionPtr daughterProductOf(const std::vector<std::string>& arguments)
759  {
760  if (arguments.size() == 1) {
761  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
762  auto func = [var](const Particle * particle) -> double {
763  double product = 1.0;
764  if (particle->getNDaughters() == 0)
765  {
766  return std::numeric_limits<double>::quiet_NaN();
767  }
768  if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
769  {
770  for (unsigned j = 0; j < particle->getNDaughters(); ++j) {
771  product *= std::get<double>(var->function(particle->getDaughter(j)));
772  }
773  } else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
774  {
775  for (unsigned j = 0; j < particle->getNDaughters(); ++j) {
776  product *= std::get<int>(var->function(particle->getDaughter(j)));
777  }
778  } else return std::numeric_limits<double>::quiet_NaN();
779  return product;
780  };
781  return func;
782  } else {
783  B2FATAL("Wrong number of arguments for meta function daughterProductOf");
784  }
785  }
786 
787  Manager::FunctionPtr daughterSumOf(const std::vector<std::string>& arguments)
788  {
789  if (arguments.size() == 1) {
790  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
791  auto func = [var](const Particle * particle) -> double {
792  double sum = 0.0;
793  if (particle->getNDaughters() == 0)
794  {
795  return std::numeric_limits<double>::quiet_NaN();
796  }
797  if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
798  {
799  for (unsigned j = 0; j < particle->getNDaughters(); ++j) {
800  sum += std::get<double>(var->function(particle->getDaughter(j)));
801  }
802  } else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
803  {
804  for (unsigned j = 0; j < particle->getNDaughters(); ++j) {
805  sum += std::get<int>(var->function(particle->getDaughter(j)));
806  }
807  } else return std::numeric_limits<double>::quiet_NaN();
808  return sum;
809  };
810  return func;
811  } else {
812  B2FATAL("Wrong number of arguments for meta function daughterSumOf");
813  }
814  }
815 
816  Manager::FunctionPtr daughterLowest(const std::vector<std::string>& arguments)
817  {
818  if (arguments.size() == 1) {
819  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
820  auto func = [var](const Particle * particle) -> double {
821  double min = std::numeric_limits<double>::quiet_NaN();
822  if (particle->getNDaughters() == 0)
823  {
824  return std::numeric_limits<double>::quiet_NaN();
825  }
826  if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
827  {
828  for (unsigned j = 0; j < particle->getNDaughters(); ++j) {
829  double iValue = std::get<double>(var->function(particle->getDaughter(j)));
830  if (std::isnan(iValue)) continue;
831  if (std::isnan(min)) min = iValue;
832  if (iValue < min) min = iValue;
833  }
834  } else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
835  {
836  for (unsigned j = 0; j < particle->getNDaughters(); ++j) {
837  int iValue = std::get<int>(var->function(particle->getDaughter(j)));
838  if (std::isnan(min)) min = iValue;
839  if (iValue < min) min = iValue;
840  }
841  }
842  return min;
843  };
844  return func;
845  } else {
846  B2FATAL("Wrong number of arguments for meta function daughterLowest");
847  }
848  }
849 
850  Manager::FunctionPtr daughterHighest(const std::vector<std::string>& arguments)
851  {
852  if (arguments.size() == 1) {
853  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
854  auto func = [var](const Particle * particle) -> double {
855  double max = std::numeric_limits<double>::quiet_NaN();
856  if (particle->getNDaughters() == 0)
857  {
858  return std::numeric_limits<double>::quiet_NaN();
859  }
860  if (std::holds_alternative<double>(var->function(particle->getDaughter(0))))
861  {
862  for (unsigned j = 0; j < particle->getNDaughters(); ++j) {
863  double iValue = std::get<double>(var->function(particle->getDaughter(j)));
864  if (std::isnan(iValue)) continue;
865  if (std::isnan(max)) max = iValue;
866  if (iValue > max) max = iValue;
867  }
868  } else if (std::holds_alternative<int>(var->function(particle->getDaughter(0))))
869  {
870  for (unsigned j = 0; j < particle->getNDaughters(); ++j) {
871  int iValue = std::get<int>(var->function(particle->getDaughter(j)));
872  if (std::isnan(max)) max = iValue;
873  if (iValue > max) max = iValue;
874  }
875  }
876  return max;
877  };
878  return func;
879  } else {
880  B2FATAL("Wrong number of arguments for meta function daughterHighest");
881  }
882  }
883 
884  Manager::FunctionPtr daughterDiffOf(const std::vector<std::string>& arguments)
885  {
886  if (arguments.size() == 3) {
887  int iDaughterNumber = 0;
888  int jDaughterNumber = 0;
889  try {
890  iDaughterNumber = Belle2::convertString<int>(arguments[0]);
891  jDaughterNumber = Belle2::convertString<int>(arguments[1]);
892  } catch (std::invalid_argument&) {
893  B2FATAL("First two arguments of daughterDiffOf meta function must be integers!");
894  }
895  auto variablename = arguments[2];
896  auto func = [variablename, iDaughterNumber, jDaughterNumber](const Particle * particle) -> double {
897  if (particle == nullptr)
898  return std::numeric_limits<double>::quiet_NaN();
899  if (iDaughterNumber >= int(particle->getNDaughters()) || jDaughterNumber >= int(particle->getNDaughters()))
900  return std::numeric_limits<double>::quiet_NaN();
901  else
902  {
903  const Variable::Manager::Var* var = Manager::Instance().getVariable(variablename);
904  auto result_j = var->function(particle->getDaughter(jDaughterNumber));
905  auto result_i = var->function(particle->getDaughter(iDaughterNumber));
906  double diff = std::numeric_limits<double>::quiet_NaN();
907  if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i)) {
908  diff = std::get<double>(result_j) - std::get<double>(result_i);
909  } else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i)) {
910  diff = std::get<int>(result_j) - std::get<int>(result_i);
911  } else {
912  throw std::runtime_error("Bad variant access");
913  }
914  if (variablename == "phi" or variablename == "clusterPhi" or std::regex_match(variablename, std::regex("use.*Frame\\(phi\\)"))
915  or std::regex_match(variablename, std::regex("use.*Frame\\(clusterPhi\\)"))) {
916  if (fabs(diff) > M_PI) {
917  if (diff > M_PI) {
918  diff = diff - 2 * M_PI;
919  } else {
920  diff = 2 * M_PI + diff;
921  }
922  }
923  }
924  return diff;
925  }
926  };
927  return func;
928  } else {
929  B2FATAL("Wrong number of arguments for meta function daughterDiffOf");
930  }
931  }
932 
933  Manager::FunctionPtr mcDaughterDiffOf(const std::vector<std::string>& arguments)
934  {
935  if (arguments.size() == 3) {
936  int iDaughterNumber = 0;
937  int jDaughterNumber = 0;
938  try {
939  iDaughterNumber = Belle2::convertString<int>(arguments[0]);
940  jDaughterNumber = Belle2::convertString<int>(arguments[1]);
941  } catch (std::invalid_argument&) {
942  B2FATAL("First two arguments of mcDaughterDiffOf meta function must be integers!");
943  }
944  auto variablename = arguments[2];
945  auto func = [variablename, iDaughterNumber, jDaughterNumber](const Particle * particle) -> double {
946  if (particle == nullptr)
947  return std::numeric_limits<double>::quiet_NaN();
948  if (iDaughterNumber >= int(particle->getNDaughters()) || jDaughterNumber >= int(particle->getNDaughters()))
949  return std::numeric_limits<double>::quiet_NaN();
950  if (particle->getDaughter(jDaughterNumber)->getMCParticle() == nullptr || particle->getDaughter(iDaughterNumber)->getMCParticle() == nullptr)
951  return std::numeric_limits<double>::quiet_NaN();
952  else
953  {
954  const MCParticle* iMcDaughter = particle->getDaughter(iDaughterNumber)->getMCParticle();
955  const MCParticle* jMcDaughter = particle->getDaughter(jDaughterNumber)->getMCParticle();
956  Particle iTmpPart(iMcDaughter);
957  Particle jTmpPart(jMcDaughter);
958  const Variable::Manager::Var* var = Manager::Instance().getVariable(variablename);
959  auto result_j = var->function(&jTmpPart);
960  auto result_i = var->function(&iTmpPart);
961  double diff = std::numeric_limits<double>::quiet_NaN();
962  if (std::holds_alternative<double>(result_j) && std::holds_alternative<double>(result_i)) {
963  diff = std::get<double>(result_j) - std::get<double>(result_i);
964  } else if (std::holds_alternative<int>(result_j) && std::holds_alternative<int>(result_i)) {
965  diff = std::get<int>(result_j) - std::get<int>(result_i);
966  } else {
967  throw std::runtime_error("Bad variant access");
968  }
969  if (variablename == "phi" or std::regex_match(variablename, std::regex("use.*Frame\\(phi\\)"))) {
970  if (fabs(diff) > M_PI) {
971  if (diff > M_PI) {
972  diff = diff - 2 * M_PI;
973  } else {
974  diff = 2 * M_PI + diff;
975  }
976  }
977  }
978  return diff;
979  }
980  };
981  return func;
982  } else {
983  B2FATAL("Wrong number of arguments for meta function mcDaughterDiffOf");
984  }
985  }
986 
987  Manager::FunctionPtr grandDaughterDiffOf(const std::vector<std::string>& arguments)
988  {
989  if (arguments.size() == 5) {
990  int iDaughterNumber = 0, jDaughterNumber = 0, agrandDaughterNumber = 0, bgrandDaughterNumber = 0;
991  try {
992  iDaughterNumber = Belle2::convertString<int>(arguments[0]);
993  jDaughterNumber = Belle2::convertString<int>(arguments[1]);
994  agrandDaughterNumber = Belle2::convertString<int>(arguments[2]);
995  bgrandDaughterNumber = Belle2::convertString<int>(arguments[3]);
996  } catch (std::invalid_argument&) {
997  B2FATAL("First four arguments of grandDaughterDiffOf meta function must be integers!");
998  }
999  auto variablename = arguments[4];
1000  auto func = [variablename, iDaughterNumber, jDaughterNumber, agrandDaughterNumber,
1001  bgrandDaughterNumber](const Particle * particle) -> double {
1002  if (particle == nullptr)
1003  return std::numeric_limits<double>::quiet_NaN();
1004  if (iDaughterNumber >= int(particle->getNDaughters()) || jDaughterNumber >= int(particle->getNDaughters()))
1005  return std::numeric_limits<double>::quiet_NaN();
1006  if (agrandDaughterNumber >= int((particle->getDaughter(iDaughterNumber))->getNDaughters()) || bgrandDaughterNumber >= int((particle->getDaughter(jDaughterNumber))->getNDaughters()))
1007  return std::numeric_limits<double>::quiet_NaN();
1008  else
1009  {
1010  const Variable::Manager::Var* var = Manager::Instance().getVariable(variablename);
1011  double diff = std::get<double>(var->function((particle->getDaughter(jDaughterNumber))->getDaughter(
1012  bgrandDaughterNumber))) - std::get<double>(var->function((particle->getDaughter(iDaughterNumber))->getDaughter(
1013  agrandDaughterNumber)));
1014  if (variablename == "phi" or variablename == "clusterPhi" or std::regex_match(variablename, std::regex("use.*Frame\\(phi\\)"))
1015  or std::regex_match(variablename, std::regex("use.*Frame\\(clusterPhi\\)"))) {
1016  if (fabs(diff) > M_PI) {
1017  if (diff > M_PI) {
1018  diff = diff - 2 * M_PI;
1019  } else {
1020  diff = 2 * M_PI + diff;
1021  }
1022  }
1023  }
1024  return diff;
1025  }
1026  };
1027  return func;
1028  } else {
1029  B2FATAL("Wrong number of arguments for meta function grandDaughterDiffOf");
1030  }
1031  }
1032 
1033  Manager::FunctionPtr daughterDiffOfPhi(const std::vector<std::string>& arguments)
1034  {
1035  std::vector<std::string> new_arguments = arguments;
1036  new_arguments.push_back(std::string("phi"));
1037  return daughterDiffOf(new_arguments);
1038  }
1039 
1040  Manager::FunctionPtr mcDaughterDiffOfPhi(const std::vector<std::string>& arguments)
1041  {
1042  std::vector<std::string> new_arguments = arguments;
1043  new_arguments.push_back(std::string("phi"));
1044  return mcDaughterDiffOf(new_arguments);
1045  }
1046 
1047  Manager::FunctionPtr grandDaughterDiffOfPhi(const std::vector<std::string>& arguments)
1048  {
1049  std::vector<std::string> new_arguments = arguments;
1050  new_arguments.push_back(std::string("phi"));
1051  return grandDaughterDiffOf(new_arguments);
1052  }
1053 
1054  Manager::FunctionPtr daughterDiffOfClusterPhi(const std::vector<std::string>& arguments)
1055  {
1056  std::vector<std::string> new_arguments = arguments;
1057  new_arguments.push_back(std::string("clusterPhi"));
1058  return daughterDiffOf(new_arguments);
1059  }
1060 
1061  Manager::FunctionPtr grandDaughterDiffOfClusterPhi(const std::vector<std::string>& arguments)
1062  {
1063  std::vector<std::string> new_arguments = arguments;
1064  new_arguments.push_back(std::string("clusterPhi"));
1065  return grandDaughterDiffOf(new_arguments);
1066  }
1067 
1068  Manager::FunctionPtr daughterDiffOfPhiCMS(const std::vector<std::string>& arguments)
1069  {
1070  std::vector<std::string> new_arguments = arguments;
1071  new_arguments.push_back(std::string("useCMSFrame(phi)"));
1072  return daughterDiffOf(new_arguments);
1073  }
1074 
1075  Manager::FunctionPtr mcDaughterDiffOfPhiCMS(const std::vector<std::string>& arguments)
1076  {
1077  std::vector<std::string> new_arguments = arguments;
1078  new_arguments.push_back(std::string("useCMSFrame(phi)"));
1079  return mcDaughterDiffOf(new_arguments);
1080  }
1081 
1082  Manager::FunctionPtr daughterDiffOfClusterPhiCMS(const std::vector<std::string>& arguments)
1083  {
1084  std::vector<std::string> new_arguments = arguments;
1085  new_arguments.push_back(std::string("useCMSFrame(clusterPhi)"));
1086  return daughterDiffOf(new_arguments);
1087  }
1088 
1089  Manager::FunctionPtr daughterNormDiffOf(const std::vector<std::string>& arguments)
1090  {
1091  if (arguments.size() == 3) {
1092  int iDaughterNumber = 0;
1093  int jDaughterNumber = 0;
1094  try {
1095  iDaughterNumber = Belle2::convertString<int>(arguments[0]);
1096  jDaughterNumber = Belle2::convertString<int>(arguments[1]);
1097  } catch (std::invalid_argument&) {
1098  B2FATAL("First two arguments of daughterDiffOf meta function must be integers!");
1099  }
1100  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[2]);
1101  auto func = [var, iDaughterNumber, jDaughterNumber](const Particle * particle) -> double {
1102  if (particle == nullptr)
1103  return std::numeric_limits<double>::quiet_NaN();
1104  if (iDaughterNumber >= int(particle->getNDaughters()) || jDaughterNumber >= int(particle->getNDaughters()))
1105  return std::numeric_limits<double>::quiet_NaN();
1106  else
1107  {
1108  double iValue, jValue;
1109  if (std::holds_alternative<double>(var->function(particle->getDaughter(jDaughterNumber)))) {
1110  iValue = std::get<double>(var->function(particle->getDaughter(iDaughterNumber)));
1111  jValue = std::get<double>(var->function(particle->getDaughter(jDaughterNumber)));
1112  } else if (std::holds_alternative<int>(var->function(particle->getDaughter(jDaughterNumber)))) {
1113  iValue = std::get<int>(var->function(particle->getDaughter(iDaughterNumber)));
1114  jValue = std::get<int>(var->function(particle->getDaughter(jDaughterNumber)));
1115  } else return std::numeric_limits<double>::quiet_NaN();
1116  return (jValue - iValue) / (jValue + iValue);
1117  }
1118  };
1119  return func;
1120  } else {
1121  B2FATAL("Wrong number of arguments for meta function daughterNormDiffOf");
1122  }
1123  }
1124 
1125  Manager::FunctionPtr daughterMotherDiffOf(const std::vector<std::string>& arguments)
1126  {
1127  if (arguments.size() == 2) {
1128  int daughterNumber = 0;
1129  try {
1130  daughterNumber = Belle2::convertString<int>(arguments[0]);
1131  } catch (std::invalid_argument&) {
1132  B2FATAL("First argument of daughterMotherDiffOf meta function must be integer!");
1133  }
1134  auto variablename = arguments[1];
1135  auto func = [variablename, daughterNumber](const Particle * particle) -> double {
1136  if (particle == nullptr)
1137  return std::numeric_limits<double>::quiet_NaN();
1138  if (daughterNumber >= int(particle->getNDaughters()))
1139  return std::numeric_limits<double>::quiet_NaN();
1140  else
1141  {
1142  const Variable::Manager::Var* var = Manager::Instance().getVariable(variablename);
1143  auto result_mother = var->function(particle);
1144  auto result_daughter = var->function(particle->getDaughter(daughterNumber));
1145  double diff = std::numeric_limits<double>::quiet_NaN();
1146  if (std::holds_alternative<double>(result_mother) && std::holds_alternative<double>(result_daughter)) {
1147  diff = std::get<double>(result_mother) - std::get<double>(result_daughter);
1148  } else if (std::holds_alternative<int>(result_mother) && std::holds_alternative<int>(result_daughter)) {
1149  diff = std::get<int>(result_mother) - std::get<int>(result_daughter);
1150  } else {
1151  throw std::runtime_error("Bad variant access");
1152  }
1153 
1154  if (variablename == "phi" or variablename == "useCMSFrame(phi)") {
1155  if (fabs(diff) > M_PI) {
1156  if (diff > M_PI) {
1157  diff = diff - 2 * M_PI;
1158  } else {
1159  diff = 2 * M_PI + diff;
1160  }
1161  }
1162  }
1163  return diff;
1164  }
1165  };
1166  return func;
1167  } else {
1168  B2FATAL("Wrong number of arguments for meta function daughterMotherDiffOf");
1169  }
1170  }
1171 
1172  Manager::FunctionPtr daughterMotherNormDiffOf(const std::vector<std::string>& arguments)
1173  {
1174  if (arguments.size() == 2) {
1175  int daughterNumber = 0;
1176  try {
1177  daughterNumber = Belle2::convertString<int>(arguments[0]);
1178  } catch (std::invalid_argument&) {
1179  B2FATAL("First argument of daughterMotherDiffOf meta function must be integer!");
1180  }
1181  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
1182  auto func = [var, daughterNumber](const Particle * particle) -> double {
1183  if (particle == nullptr)
1184  return std::numeric_limits<double>::quiet_NaN();
1185  if (daughterNumber >= int(particle->getNDaughters()))
1186  return std::numeric_limits<double>::quiet_NaN();
1187  else
1188  {
1189  double daughterValue = 0.0, motherValue = 0.0;
1190  if (std::holds_alternative<double>(var->function(particle))) {
1191  daughterValue = std::get<double>(var->function(particle->getDaughter(daughterNumber)));
1192  motherValue = std::get<double>(var->function(particle));
1193  } else if (std::holds_alternative<int>(var->function(particle))) {
1194  daughterValue = std::get<int>(var->function(particle->getDaughter(daughterNumber)));
1195  motherValue = std::get<int>(var->function(particle));
1196  }
1197  return (motherValue - daughterValue) / (motherValue + daughterValue);
1198  }
1199  };
1200  return func;
1201  } else {
1202  B2FATAL("Wrong number of arguments for meta function daughterMotherNormDiffOf");
1203  }
1204  }
1205 
1206  Manager::FunctionPtr daughterAngle(const std::vector<std::string>& arguments)
1207  {
1208  if (arguments.size() == 2 || arguments.size() == 3) {
1209 
1210  auto func = [arguments](const Particle * particle) -> double {
1211  if (particle == nullptr)
1212  return std::numeric_limits<double>::quiet_NaN();
1213 
1214  std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1215  const auto& frame = ReferenceFrame::GetCurrent();
1216 
1217  // Parses the generalized indexes and fetches the 4-momenta of the particles of interest
1218  for (auto& generalizedIndex : arguments)
1219  {
1220  const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1221  if (dauPart)
1222  pDaus.push_back(frame.getMomentum(dauPart));
1223  else {
1224  B2WARNING("Trying to access a daughter that does not exist. Index = " << generalizedIndex);
1225  return std::numeric_limits<double>::quiet_NaN();
1226  }
1227  }
1228 
1229  // Calculates the angle between the selected particles
1230  if (pDaus.size() == 2)
1231  return B2Vector3D(pDaus[0].Vect()).Angle(B2Vector3D(pDaus[1].Vect()));
1232  else
1233  return B2Vector3D(pDaus[2].Vect()).Angle(B2Vector3D(pDaus[0].Vect() + pDaus[1].Vect()));
1234  };
1235  return func;
1236  } else {
1237  B2FATAL("Wrong number of arguments for meta function daughterAngle");
1238  }
1239  }
1240 
1241  double grandDaughterDecayAngle(const Particle* particle, const std::vector<double>& arguments)
1242  {
1243  if (arguments.size() == 2) {
1244 
1245  if (!particle)
1246  return std::numeric_limits<double>::quiet_NaN();
1247 
1248  int daughterIndex = std::lround(arguments[0]);
1249  if (daughterIndex >= int(particle->getNDaughters()))
1250  return std::numeric_limits<float>::quiet_NaN();
1251  const Particle* dau = particle->getDaughter(daughterIndex);
1252 
1253  int grandDaughterIndex = std::lround(arguments[1]);
1254  if (grandDaughterIndex >= int(dau->getNDaughters()))
1255  return std::numeric_limits<float>::quiet_NaN();
1256 
1257  B2Vector3D boost = dau->get4Vector().BoostToCM();
1258 
1259  ROOT::Math::PxPyPzEVector motherMomentum = - particle->get4Vector();
1260  motherMomentum = ROOT::Math::Boost(boost) * motherMomentum;
1261 
1262  ROOT::Math::PxPyPzEVector grandDaughterMomentum = dau->getDaughter(grandDaughterIndex)->get4Vector();
1263  grandDaughterMomentum = ROOT::Math::Boost(boost) * grandDaughterMomentum;
1264 
1265  return B2Vector3D(motherMomentum.Vect()).Angle(B2Vector3D(grandDaughterMomentum.Vect()));
1266 
1267  } else {
1268  B2FATAL("The variable grandDaughterDecayAngle needs exactly two integers as arguments!");
1269  }
1270  }
1271 
1272  Manager::FunctionPtr mcDaughterAngle(const std::vector<std::string>& arguments)
1273  {
1274  if (arguments.size() == 2 || arguments.size() == 3) {
1275 
1276  auto func = [arguments](const Particle * particle) -> double {
1277  if (particle == nullptr)
1278  return std::numeric_limits<double>::quiet_NaN();
1279 
1280  std::vector<ROOT::Math::PxPyPzEVector> pDaus;
1281  const auto& frame = ReferenceFrame::GetCurrent();
1282 
1283  // Parses the generalized indexes and fetches the 4-momenta of the particles of interest
1284  if (particle->getParticleSource() == Particle::EParticleSourceObject::c_MCParticle) // Check if MCParticle
1285  {
1286  for (auto& generalizedIndex : arguments) {
1287  const MCParticle* mcPart = particle->getMCParticle();
1288  if (mcPart == nullptr)
1289  return std::numeric_limits<double>::quiet_NaN();
1290  const MCParticle* dauMcPart = mcPart->getParticleFromGeneralizedIndexString(generalizedIndex);
1291  if (dauMcPart == nullptr)
1292  return std::numeric_limits<double>::quiet_NaN();
1293 
1294  pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1295  }
1296  } else
1297  {
1298  for (auto& generalizedIndex : arguments) {
1299  const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1300  if (dauPart == nullptr)
1301  return std::numeric_limits<double>::quiet_NaN();
1302 
1303  const MCParticle* dauMcPart = dauPart->getMCParticle();
1304  if (dauMcPart == nullptr)
1305  return std::numeric_limits<double>::quiet_NaN();
1306 
1307  pDaus.push_back(frame.getMomentum(dauMcPart->get4Vector()));
1308  }
1309  }
1310 
1311  // Calculates the angle between the selected particles
1312  if (pDaus.size() == 2)
1313  return B2Vector3D(pDaus[0].Vect()).Angle(B2Vector3D(pDaus[1].Vect()));
1314  else
1315  return B2Vector3D(pDaus[2].Vect()).Angle(B2Vector3D(pDaus[0].Vect() + pDaus[1].Vect()));
1316  };
1317  return func;
1318  } else {
1319  B2FATAL("Wrong number of arguments for meta function mcDaughterAngle");
1320  }
1321  }
1322 
1323  double daughterClusterAngleInBetween(const Particle* particle, const std::vector<double>& daughterIndices)
1324  {
1325  if (daughterIndices.size() == 2) {
1326  int daughterIndexi = std::lround(daughterIndices[0]);
1327  int daughterIndexj = std::lround(daughterIndices[1]);
1328  if (std::max(daughterIndexi, daughterIndexj) >= int(particle->getNDaughters())) {
1329  return std::numeric_limits<double>::quiet_NaN();
1330  } else {
1331  const ECLCluster* clusteri = particle->getDaughter(daughterIndexi)->getECLCluster();
1332  const ECLCluster* clusterj = particle->getDaughter(daughterIndexj)->getECLCluster();
1333  if (clusteri and clusterj) {
1334  const auto& frame = ReferenceFrame::GetCurrent();
1335  const ECLCluster::EHypothesisBit clusteriBit = (particle->getDaughter(daughterIndexi))->getECLClusterEHypothesisBit();
1336  const ECLCluster::EHypothesisBit clusterjBit = (particle->getDaughter(daughterIndexj))->getECLClusterEHypothesisBit();
1337  ClusterUtils clusutils;
1338  B2Vector3D pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit)).Vect();
1339  B2Vector3D pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit)).Vect();
1340  return pi.Angle(pj);
1341  }
1342  return std::numeric_limits<float>::quiet_NaN();
1343  }
1344  } else if (daughterIndices.size() == 3) {
1345  int daughterIndexi = std::lround(daughterIndices[0]);
1346  int daughterIndexj = std::lround(daughterIndices[1]);
1347  int daughterIndexk = std::lround(daughterIndices[2]);
1348  if (std::max(std::max(daughterIndexi, daughterIndexj), daughterIndexk) >= int(particle->getNDaughters())) {
1349  return std::numeric_limits<double>::quiet_NaN();
1350  } else {
1351  const ECLCluster* clusteri = (particle->getDaughter(daughterIndices[0]))->getECLCluster();
1352  const ECLCluster* clusterj = (particle->getDaughter(daughterIndices[1]))->getECLCluster();
1353  const ECLCluster* clusterk = (particle->getDaughter(daughterIndices[2]))->getECLCluster();
1354  if (clusteri and clusterj and clusterk) {
1355  const auto& frame = ReferenceFrame::GetCurrent();
1356  const ECLCluster::EHypothesisBit clusteriBit = (particle->getDaughter(daughterIndices[0]))->getECLClusterEHypothesisBit();
1357  const ECLCluster::EHypothesisBit clusterjBit = (particle->getDaughter(daughterIndices[1]))->getECLClusterEHypothesisBit();
1358  const ECLCluster::EHypothesisBit clusterkBit = (particle->getDaughter(daughterIndices[2]))->getECLClusterEHypothesisBit();
1359  ClusterUtils clusutils;
1360  B2Vector3D pi = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusteri, clusteriBit)).Vect();
1361  B2Vector3D pj = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterj, clusterjBit)).Vect();
1362  B2Vector3D pk = frame.getMomentum(clusutils.Get4MomentumFromCluster(clusterk, clusterkBit)).Vect();
1363  return pk.Angle(pi + pj);
1364  }
1365  return std::numeric_limits<float>::quiet_NaN();
1366  }
1367  } else {
1368  B2FATAL("Wrong number of arguments for daughterClusterAngleInBetween!");
1369  }
1370  }
1371 
1372  Manager::FunctionPtr daughterInvM(const std::vector<std::string>& arguments)
1373  {
1374  if (arguments.size() > 1) {
1375  auto func = [arguments](const Particle * particle) -> double {
1376  const auto& frame = ReferenceFrame::GetCurrent();
1377  ROOT::Math::PxPyPzEVector pSum;
1378 
1379  for (auto& generalizedIndex : arguments)
1380  {
1381  const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
1382  if (dauPart)
1383  pSum += frame.getMomentum(dauPart);
1384  else {
1385  return std::numeric_limits<float>::quiet_NaN();
1386  }
1387  }
1388  return pSum.M();
1389  };
1390  return func;
1391  } else {
1392  B2FATAL("Wrong number of arguments for meta function daughterInvM. At least two integers are needed.");
1393  }
1394  }
1395 
1396  Manager::FunctionPtr modulo(const std::vector<std::string>& arguments)
1397  {
1398  if (arguments.size() == 2) {
1399  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1400  int divideBy = 1;
1401  try {
1402  divideBy = Belle2::convertString<int>(arguments[1]);
1403  } catch (std::invalid_argument&) {
1404  B2FATAL("Second argument of modulo meta function must be integer!");
1405  }
1406  auto func = [var, divideBy](const Particle * particle) -> int {
1407  auto var_result = var->function(particle);
1408  if (std::holds_alternative<double>(var_result))
1409  {
1410  return int(std::get<double>(var_result)) % divideBy;
1411  } else if (std::holds_alternative<int>(var_result))
1412  {
1413  return std::get<int>(var_result) % divideBy;
1414  } else if (std::holds_alternative<bool>(var_result))
1415  {
1416  return int(std::get<bool>(var_result)) % divideBy;
1417  } else return 0;
1418  };
1419  return func;
1420  } else {
1421  B2FATAL("Wrong number of arguments for meta function modulo");
1422  }
1423  }
1424 
1425  Manager::FunctionPtr isNAN(const std::vector<std::string>& arguments)
1426  {
1427  if (arguments.size() == 1) {
1428  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1429 
1430  auto func = [var](const Particle * particle) -> bool { return std::isnan(std::get<double>(var->function(particle))); };
1431  return func;
1432  } else {
1433  B2FATAL("Wrong number of arguments for meta function isNAN");
1434  }
1435  }
1436 
1437  Manager::FunctionPtr ifNANgiveX(const std::vector<std::string>& arguments)
1438  {
1439  if (arguments.size() == 2) {
1440  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1441  double defaultOutput;
1442  try {
1443  defaultOutput = Belle2::convertString<double>(arguments[1]);
1444  } catch (std::invalid_argument&) {
1445  B2FATAL("The second argument of ifNANgiveX meta function must be a number!");
1446  }
1447  auto func = [var, defaultOutput](const Particle * particle) -> double {
1448  double output = std::get<double>(var->function(particle));
1449  if (std::isnan(output)) return defaultOutput;
1450  else return output;
1451  };
1452  return func;
1453  } else {
1454  B2FATAL("Wrong number of arguments for meta function ifNANgiveX");
1455  }
1456  }
1457 
1458  Manager::FunctionPtr isInfinity(const std::vector<std::string>& arguments)
1459  {
1460  if (arguments.size() == 1) {
1461  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1462 
1463  auto func = [var](const Particle * particle) -> bool { return std::isinf(std::get<double>(var->function(particle))); };
1464  return func;
1465  } else {
1466  B2FATAL("Wrong number of arguments for meta function isInfinity");
1467  }
1468  }
1469 
1470  Manager::FunctionPtr unmask(const std::vector<std::string>& arguments)
1471  {
1472  if (arguments.size() >= 2) {
1473  // get the function pointer of variable to be unmasked
1474  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1475 
1476  // get the final mask which summarize all the input masks
1477  int finalMask = 0;
1478  for (size_t i = 1; i < arguments.size(); ++i) {
1479  try {
1480  finalMask |= Belle2::convertString<int>(arguments[i]);
1481  } catch (std::invalid_argument&) {
1482  B2FATAL("The input flags to meta function unmask() should be integer!");
1483  return nullptr;
1484  }
1485  }
1486 
1487  // unmask the variable
1488  auto func = [var, finalMask](const Particle * particle) -> double {
1489  int value = 0;
1490  auto var_result = var->function(particle);
1491  if (std::holds_alternative<double>(var_result))
1492  {
1493  // judge if the value is nan before unmasking
1494  if (std::isnan(std::get<double>(var_result))) {
1495  return std::numeric_limits<double>::quiet_NaN();
1496  }
1497  value = int(std::get<double>(var_result));
1498  } else if (std::holds_alternative<int>(var_result))
1499  {
1500  value = std::get<int>(var_result);
1501  }
1502 
1503  // apply the final mask
1504  value &= (~finalMask);
1505 
1506  return value;
1507  };
1508  return func;
1509 
1510  } else {
1511  B2FATAL("Meta function unmask needs at least two arguments!");
1512  }
1513  }
1514 
1515  Manager::FunctionPtr conditionalVariableSelector(const std::vector<std::string>& arguments)
1516  {
1517  if (arguments.size() == 3) {
1518 
1519  std::string cutString = arguments[0];
1520  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
1521 
1522  const Variable::Manager::Var* variableIfTrue = Manager::Instance().getVariable(arguments[1]);
1523  const Variable::Manager::Var* variableIfFalse = Manager::Instance().getVariable(arguments[2]);
1524 
1525  auto func = [cut, variableIfTrue, variableIfFalse](const Particle * particle) -> double {
1526  if (particle == nullptr)
1527  return std::numeric_limits<float>::quiet_NaN();
1528  if (cut->check(particle))
1529  {
1530  auto var_result = variableIfTrue->function(particle);
1531  if (std::holds_alternative<double>(var_result)) {
1532  return std::get<double>(var_result);
1533  } else if (std::holds_alternative<int>(var_result)) {
1534  return std::get<int>(var_result);
1535  } else if (std::holds_alternative<bool>(var_result)) {
1536  return std::get<bool>(var_result);
1537  } else return std::numeric_limits<double>::quiet_NaN();
1538  } else
1539  {
1540  auto var_result = variableIfFalse->function(particle);
1541  if (std::holds_alternative<double>(var_result)) {
1542  return std::get<double>(var_result);
1543  } else if (std::holds_alternative<int>(var_result)) {
1544  return std::get<int>(var_result);
1545  } else if (std::holds_alternative<bool>(var_result)) {
1546  return std::get<bool>(var_result);
1547  } else return std::numeric_limits<double>::quiet_NaN();
1548  }
1549  };
1550  return func;
1551 
1552  } else {
1553  B2FATAL("Wrong number of arguments for meta function conditionalVariableSelector");
1554  }
1555  }
1556 
1557  Manager::FunctionPtr pValueCombination(const std::vector<std::string>& arguments)
1558  {
1559  if (arguments.size() > 0) {
1560  std::vector<const Variable::Manager::Var*> variables;
1561  for (auto& argument : arguments)
1562  variables.push_back(Manager::Instance().getVariable(argument));
1563 
1564  auto func = [variables, arguments](const Particle * particle) -> double {
1565  double pValueProduct = 1.;
1566  for (auto variable : variables)
1567  {
1568  double pValue = std::get<double>(variable->function(particle));
1569  if (pValue < 0)
1570  return -1;
1571  else
1572  pValueProduct *= pValue;
1573  }
1574  double pValueSum = 1.;
1575  double factorial = 1.;
1576  for (unsigned int i = 1; i < arguments.size(); ++i)
1577  {
1578  factorial *= i;
1579  pValueSum += pow(-std::log(pValueProduct), i) / factorial;
1580  }
1581  return pValueProduct * pValueSum;
1582  };
1583  return func;
1584  } else {
1585  B2FATAL("Wrong number of arguments for meta function pValueCombination");
1586  }
1587  }
1588 
1589  Manager::FunctionPtr abs(const std::vector<std::string>& arguments)
1590  {
1591  if (arguments.size() == 1) {
1592  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1593  auto func = [var](const Particle * particle) -> double {
1594  auto var_result = var->function(particle);
1595  if (std::holds_alternative<double>(var_result))
1596  {
1597  return std::abs(std::get<double>(var_result));
1598  } else if (std::holds_alternative<int>(var_result))
1599  {
1600  return std::abs(std::get<int>(var_result));
1601  } else return std::numeric_limits<double>::quiet_NaN();
1602  };
1603  return func;
1604  } else {
1605  B2FATAL("Wrong number of arguments for meta function abs");
1606  }
1607  }
1608 
1609  Manager::FunctionPtr max(const std::vector<std::string>& arguments)
1610  {
1611  if (arguments.size() == 2) {
1612  const Variable::Manager::Var* var1 = Manager::Instance().getVariable(arguments[0]);
1613  const Variable::Manager::Var* var2 = Manager::Instance().getVariable(arguments[1]);
1614 
1615  if (!var1 or !var2)
1616  B2FATAL("One or both of the used variables doesn't exist!");
1617 
1618  auto func = [var1, var2](const Particle * particle) -> double {
1619  double val1, val2;
1620  auto var_result1 = var1->function(particle);
1621  auto var_result2 = var2->function(particle);
1622  if (std::holds_alternative<double>(var_result1))
1623  {
1624  val1 = std::get<double>(var_result1);
1625  } else if (std::holds_alternative<int>(var_result1))
1626  {
1627  val1 = std::get<int>(var_result1);
1628  }
1629  if (std::holds_alternative<double>(var_result2))
1630  {
1631  val2 = std::get<double>(var_result2);
1632  } else if (std::holds_alternative<int>(var_result2))
1633  {
1634  val2 = std::get<int>(var_result2);
1635  }
1636  return std::max(val1, val2);
1637  };
1638  return func;
1639  } else {
1640  B2FATAL("Wrong number of arguments for meta function max");
1641  }
1642  }
1643 
1644  Manager::FunctionPtr min(const std::vector<std::string>& arguments)
1645  {
1646  if (arguments.size() == 2) {
1647  const Variable::Manager::Var* var1 = Manager::Instance().getVariable(arguments[0]);
1648  const Variable::Manager::Var* var2 = Manager::Instance().getVariable(arguments[1]);
1649 
1650  if (!var1 or !var2)
1651  B2FATAL("One or both of the used variables doesn't exist!");
1652 
1653  auto func = [var1, var2](const Particle * particle) -> double {
1654  double val1, val2;
1655  auto var_result1 = var1->function(particle);
1656  auto var_result2 = var2->function(particle);
1657  if (std::holds_alternative<double>(var_result1))
1658  {
1659  val1 = std::get<double>(var_result1);
1660  } else if (std::holds_alternative<int>(var_result1))
1661  {
1662  val1 = std::get<int>(var_result1);
1663  }
1664  if (std::holds_alternative<double>(var_result2))
1665  {
1666  val2 = std::get<double>(var_result2);
1667  } else if (std::holds_alternative<int>(var_result2))
1668  {
1669  val2 = std::get<int>(var_result2);
1670  }
1671  return std::min(val1, val2);
1672  };
1673  return func;
1674  } else {
1675  B2FATAL("Wrong number of arguments for meta function min");
1676  }
1677  }
1678 
1679  Manager::FunctionPtr sin(const std::vector<std::string>& arguments)
1680  {
1681  if (arguments.size() == 1) {
1682  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1683  auto func = [var](const Particle * particle) -> double {
1684  auto var_result = var->function(particle);
1685  if (std::holds_alternative<double>(var_result))
1686  return std::sin(std::get<double>(var_result));
1687  else if (std::holds_alternative<int>(var_result))
1688  return std::sin(std::get<int>(var_result));
1689  else return std::numeric_limits<double>::quiet_NaN();
1690  };
1691  return func;
1692  } else {
1693  B2FATAL("Wrong number of arguments for meta function sin");
1694  }
1695  }
1696 
1697  Manager::FunctionPtr asin(const std::vector<std::string>& arguments)
1698  {
1699  if (arguments.size() == 1) {
1700  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1701  auto func = [var](const Particle * particle) -> double {
1702  auto var_result = var->function(particle);
1703  if (std::holds_alternative<double>(var_result))
1704  return std::asin(std::get<double>(var_result));
1705  else if (std::holds_alternative<int>(var_result))
1706  return std::asin(std::get<int>(var_result));
1707  else return std::numeric_limits<double>::quiet_NaN();
1708  };
1709  return func;
1710  } else {
1711  B2FATAL("Wrong number of arguments for meta function asin");
1712  }
1713  }
1714 
1715  Manager::FunctionPtr cos(const std::vector<std::string>& arguments)
1716  {
1717  if (arguments.size() == 1) {
1718  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1719  auto func = [var](const Particle * particle) -> double {
1720  auto var_result = var->function(particle);
1721  if (std::holds_alternative<double>(var_result))
1722  return std::cos(std::get<double>(var_result));
1723  else if (std::holds_alternative<int>(var_result))
1724  return std::cos(std::get<int>(var_result));
1725  else return std::numeric_limits<double>::quiet_NaN();
1726  };
1727  return func;
1728  } else {
1729  B2FATAL("Wrong number of arguments for meta function cos");
1730  }
1731  }
1732 
1733  Manager::FunctionPtr acos(const std::vector<std::string>& arguments)
1734  {
1735  if (arguments.size() == 1) {
1736  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1737  auto func = [var](const Particle * particle) -> double {
1738  auto var_result = var->function(particle);
1739  if (std::holds_alternative<double>(var_result))
1740  return std::acos(std::get<double>(var_result));
1741  else if (std::holds_alternative<int>(var_result))
1742  return std::acos(std::get<int>(var_result));
1743  else return std::numeric_limits<double>::quiet_NaN();
1744  };
1745  return func;
1746  } else {
1747  B2FATAL("Wrong number of arguments for meta function acos");
1748  }
1749  }
1750 
1751  Manager::FunctionPtr tan(const std::vector<std::string>& arguments)
1752  {
1753  if (arguments.size() == 1) {
1754  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1755  auto func = [var](const Particle * particle) -> double { return std::tan(std::get<double>(var->function(particle))); };
1756  return func;
1757  } else {
1758  B2FATAL("Wrong number of arguments for meta function tan");
1759  }
1760  }
1761 
1762  Manager::FunctionPtr atan(const std::vector<std::string>& arguments)
1763  {
1764  if (arguments.size() == 1) {
1765  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1766  auto func = [var](const Particle * particle) -> double { return std::atan(std::get<double>(var->function(particle))); };
1767  return func;
1768  } else {
1769  B2FATAL("Wrong number of arguments for meta function atan");
1770  }
1771  }
1772 
1773  Manager::FunctionPtr exp(const std::vector<std::string>& arguments)
1774  {
1775  if (arguments.size() == 1) {
1776  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1777  auto func = [var](const Particle * particle) -> double {
1778  auto var_result = var->function(particle);
1779  if (std::holds_alternative<double>(var_result))
1780  return std::exp(std::get<double>(var_result));
1781  else if (std::holds_alternative<int>(var_result))
1782  return std::exp(std::get<int>(var_result));
1783  else return std::numeric_limits<double>::quiet_NaN();
1784  };
1785  return func;
1786  } else {
1787  B2FATAL("Wrong number of arguments for meta function exp");
1788  }
1789  }
1790 
1791  Manager::FunctionPtr log(const std::vector<std::string>& arguments)
1792  {
1793  if (arguments.size() == 1) {
1794  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1795  auto func = [var](const Particle * particle) -> double {
1796  auto var_result = var->function(particle);
1797  if (std::holds_alternative<double>(var_result))
1798  return std::log(std::get<double>(var_result));
1799  else if (std::holds_alternative<int>(var_result))
1800  return std::log(std::get<int>(var_result));
1801  else return std::numeric_limits<double>::quiet_NaN();
1802  };
1803  return func;
1804  } else {
1805  B2FATAL("Wrong number of arguments for meta function log");
1806  }
1807  }
1808 
1809  Manager::FunctionPtr log10(const std::vector<std::string>& arguments)
1810  {
1811  if (arguments.size() == 1) {
1812  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1813  auto func = [var](const Particle * particle) -> double {
1814  auto var_result = var->function(particle);
1815  if (std::holds_alternative<double>(var_result))
1816  return std::log10(std::get<double>(var_result));
1817  else if (std::holds_alternative<int>(var_result))
1818  return std::log10(std::get<int>(var_result));
1819  else return std::numeric_limits<double>::quiet_NaN();
1820  };
1821  return func;
1822  } else {
1823  B2FATAL("Wrong number of arguments for meta function log10");
1824  }
1825  }
1826 
1827  Manager::FunctionPtr daughter(const std::vector<std::string>& arguments)
1828  {
1829  if (arguments.size() == 2) {
1830  int daughterNumber = 0;
1831  try {
1832  daughterNumber = Belle2::convertString<int>(arguments[0]);
1833  } catch (std::invalid_argument&) {
1834  B2FATAL("First argument of daughter meta function must be integer!");
1835  }
1836  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
1837  auto func = [var, daughterNumber](const Particle * particle) -> double {
1838  if (particle == nullptr)
1839  return std::numeric_limits<float>::quiet_NaN();
1840  if (daughterNumber >= int(particle->getNDaughters()))
1841  return std::numeric_limits<float>::quiet_NaN();
1842  else
1843  {
1844  auto var_result = var->function(particle->getDaughter(daughterNumber));
1845  if (std::holds_alternative<double>(var_result)) {
1846  return std::get<double>(var_result);
1847  } else if (std::holds_alternative<int>(var_result)) {
1848  return std::get<int>(var_result);
1849  } else if (std::holds_alternative<bool>(var_result)) {
1850  return std::get<bool>(var_result);
1851  } else return std::numeric_limits<double>::quiet_NaN();
1852  }
1853  };
1854  return func;
1855  } else {
1856  B2FATAL("Wrong number of arguments for meta function daughter");
1857  }
1858  }
1859 
1860  Manager::FunctionPtr mcDaughter(const std::vector<std::string>& arguments)
1861  {
1862  if (arguments.size() == 2) {
1863  int daughterNumber = 0;
1864  try {
1865  daughterNumber = Belle2::convertString<int>(arguments[0]);
1866  } catch (std::invalid_argument&) {
1867  B2FATAL("First argument of mcDaughter meta function must be integer!");
1868  }
1869  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
1870  auto func = [var, daughterNumber](const Particle * particle) -> double {
1871  if (particle == nullptr)
1872  return std::numeric_limits<float>::quiet_NaN();
1873  if (particle->getMCParticle()) // has MC match or is MCParticle
1874  {
1875  if (daughterNumber >= int(particle->getMCParticle()->getNDaughters())) {
1876  return std::numeric_limits<float>::quiet_NaN();
1877  }
1878  Particle tempParticle = Particle(particle->getMCParticle()->getDaughters().at(daughterNumber));
1879  auto var_result = var->function(&tempParticle);
1880  if (std::holds_alternative<double>(var_result)) {
1881  return std::get<double>(var_result);
1882  } else if (std::holds_alternative<int>(var_result)) {
1883  return std::get<int>(var_result);
1884  } else if (std::holds_alternative<bool>(var_result)) {
1885  return std::get<bool>(var_result);
1886  } else return std::numeric_limits<double>::quiet_NaN();
1887  } else
1888  {
1889  return std::numeric_limits<float>::quiet_NaN();
1890  }
1891  };
1892  return func;
1893  } else {
1894  B2FATAL("Wrong number of arguments for meta function mcDaughter");
1895  }
1896  }
1897 
1898  Manager::FunctionPtr mcMother(const std::vector<std::string>& arguments)
1899  {
1900  if (arguments.size() == 1) {
1901  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1902  auto func = [var](const Particle * particle) -> double {
1903  if (particle == nullptr)
1904  return std::numeric_limits<float>::quiet_NaN();
1905  if (particle->getMCParticle()) // has MC match or is MCParticle
1906  {
1907  if (particle->getMCParticle()->getMother() == nullptr) {
1908  return std::numeric_limits<float>::quiet_NaN();
1909  }
1910  Particle tempParticle = Particle(particle->getMCParticle()->getMother());
1911  auto var_result = var->function(&tempParticle);
1912  if (std::holds_alternative<double>(var_result)) {
1913  return std::get<double>(var_result);
1914  } else if (std::holds_alternative<int>(var_result)) {
1915  return std::get<int>(var_result);
1916  } else if (std::holds_alternative<bool>(var_result)) {
1917  return std::get<bool>(var_result);
1918  } else return std::numeric_limits<double>::quiet_NaN();
1919  } else
1920  {
1921  return std::numeric_limits<float>::quiet_NaN();
1922  }
1923  };
1924  return func;
1925  } else {
1926  B2FATAL("Wrong number of arguments for meta function mcMother");
1927  }
1928  }
1929 
1930  Manager::FunctionPtr genParticle(const std::vector<std::string>& arguments)
1931  {
1932  if (arguments.size() == 2) {
1933  int particleNumber = 0;
1934  try {
1935  particleNumber = Belle2::convertString<int>(arguments[0]);
1936  } catch (std::invalid_argument&) {
1937  B2FATAL("First argument of genParticle meta function must be integer!");
1938  }
1939  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
1940 
1941  auto func = [var, particleNumber](const Particle*) -> double {
1942  StoreArray<MCParticle> mcParticles("MCParticles");
1943  if (particleNumber >= mcParticles.getEntries())
1944  {
1945  return std::numeric_limits<float>::quiet_NaN();
1946  }
1947 
1948  MCParticle* mcParticle = mcParticles[particleNumber];
1949  Particle part = Particle(mcParticle);
1950  auto var_result = var->function(&part);
1951  if (std::holds_alternative<double>(var_result))
1952  {
1953  return std::get<double>(var_result);
1954  } else if (std::holds_alternative<int>(var_result))
1955  {
1956  return std::get<int>(var_result);
1957  } else if (std::holds_alternative<bool>(var_result))
1958  {
1959  return std::get<bool>(var_result);
1960  } else return std::numeric_limits<double>::quiet_NaN();
1961  };
1962  return func;
1963  } else {
1964  B2FATAL("Wrong number of arguments for meta function genParticle");
1965  }
1966  }
1967 
1968  Manager::FunctionPtr genUpsilon4S(const std::vector<std::string>& arguments)
1969  {
1970  if (arguments.size() == 1) {
1971  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
1972 
1973  auto func = [var](const Particle*) -> double {
1974  StoreArray<MCParticle> mcParticles("MCParticles");
1975  if (mcParticles.getEntries() == 0)
1976  {
1977  return std::numeric_limits<float>::quiet_NaN();
1978  }
1979 
1980  MCParticle* mcUpsilon4S = mcParticles[0];
1981  if (mcUpsilon4S->getPDG() != 300553)
1982  {
1983  return std::numeric_limits<float>::quiet_NaN();
1984  }
1985 
1986  Particle upsilon4S = Particle(mcUpsilon4S);
1987  auto var_result = var->function(&upsilon4S);
1988  if (std::holds_alternative<double>(var_result))
1989  {
1990  return std::get<double>(var_result);
1991  } else if (std::holds_alternative<int>(var_result))
1992  {
1993  return std::get<int>(var_result);
1994  } else if (std::holds_alternative<bool>(var_result))
1995  {
1996  return std::get<bool>(var_result);
1997  } else return std::numeric_limits<double>::quiet_NaN();
1998  };
1999  return func;
2000  } else {
2001  B2FATAL("Wrong number of arguments for meta function genUpsilon4S");
2002  }
2003  }
2004 
2005  Manager::FunctionPtr getVariableByRank(const std::vector<std::string>& arguments)
2006  {
2007  if (arguments.size() == 4) {
2008  std::string listName = arguments[0];
2009  std::string rankedVariableName = arguments[1];
2010  std::string returnVariableName = arguments[2];
2011  std::string extraInfoName = rankedVariableName + "_rank";
2012  int rank = 1;
2013  try {
2014  rank = Belle2::convertString<int>(arguments[3]);
2015  } catch (std::invalid_argument&) {
2016  B2ERROR("3rd argument of getVariableByRank meta function (Rank) must be an integer!");
2017  return nullptr;
2018  }
2019 
2020  const Variable::Manager::Var* var = Manager::Instance().getVariable(returnVariableName);
2021  auto func = [var, rank, extraInfoName, listName](const Particle*)-> double {
2022  StoreObjPtr<ParticleList> list(listName);
2023 
2024  const unsigned int numParticles = list->getListSize();
2025  for (unsigned int i = 0; i < numParticles; i++)
2026  {
2027  const Particle* p = list->getParticle(i);
2028  if (p->getExtraInfo(extraInfoName) == rank) {
2029  auto var_result = var->function(p);
2030  if (std::holds_alternative<double>(var_result)) {
2031  return std::get<double>(var_result);
2032  } else if (std::holds_alternative<int>(var_result)) {
2033  return std::get<int>(var_result);
2034  } else if (std::holds_alternative<bool>(var_result)) {
2035  return std::get<bool>(var_result);
2036  } else return std::numeric_limits<double>::quiet_NaN();
2037  }
2038  }
2039  // return 0;
2040  return std::numeric_limits<double>::signaling_NaN();
2041  };
2042  return func;
2043  } else {
2044  B2FATAL("Wrong number of arguments for meta function getVariableByRank");
2045  }
2046  }
2047 
2048  Manager::FunctionPtr countInList(const std::vector<std::string>& arguments)
2049  {
2050  if (arguments.size() == 1 or arguments.size() == 2) {
2051 
2052  std::string listName = arguments[0];
2053  std::string cutString = "";
2054 
2055  if (arguments.size() == 2) {
2056  cutString = arguments[1];
2057  }
2058 
2059  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
2060 
2061  auto func = [listName, cut](const Particle*) -> int {
2062 
2063  StoreObjPtr<ParticleList> list(listName);
2064  int sum = 0;
2065  for (unsigned int i = 0; i < list->getListSize(); i++)
2066  {
2067  const Particle* particle = list->getParticle(i);
2068  if (cut->check(particle)) {
2069  sum++;
2070  }
2071  }
2072  return sum;
2073  };
2074  return func;
2075  } else {
2076  B2FATAL("Wrong number of arguments for meta function countInList");
2077  }
2078  }
2079 
2080  Manager::FunctionPtr veto(const std::vector<std::string>& arguments)
2081  {
2082  if (arguments.size() == 2 or arguments.size() == 3) {
2083 
2084  std::string roeListName = arguments[0];
2085  std::string cutString = arguments[1];
2086  int pdgCode = Const::electron.getPDGCode();
2087  if (arguments.size() == 2) {
2088  B2INFO("Use pdgCode of electron as default in meta variable veto, other arguments: " << roeListName << ", " << cutString);
2089  } else {
2090  try {
2091  pdgCode = Belle2::convertString<int>(arguments[2]);;
2092  } catch (std::invalid_argument&) {
2093  B2FATAL("Third argument of veto meta function must be integer!");
2094  }
2095  }
2096 
2098  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
2099 
2100  auto func = [roeListName, cut, pdgCode, flavourType](const Particle * particle) -> bool {
2101  StoreObjPtr<ParticleList> roeList(roeListName);
2102  ROOT::Math::PxPyPzEVector vec = particle->get4Vector();
2103  for (unsigned int i = 0; i < roeList->getListSize(); i++)
2104  {
2105  const Particle* roeParticle = roeList->getParticle(i);
2106  if (not particle->overlapsWith(roeParticle)) {
2107  ROOT::Math::PxPyPzEVector tempCombination = roeParticle->get4Vector() + vec;
2108  std::vector<int> indices = { particle->getArrayIndex(), roeParticle->getArrayIndex() };
2109  Particle tempParticle = Particle(tempCombination, pdgCode, flavourType, indices, particle->getArrayPointer());
2110  if (cut->check(&tempParticle)) {
2111  return 1;
2112  }
2113  }
2114  }
2115  return 0;
2116  };
2117  return func;
2118  } else {
2119  B2FATAL("Wrong number of arguments for meta function veto");
2120  }
2121  }
2122 
2123  Manager::FunctionPtr countDaughters(const std::vector<std::string>& arguments)
2124  {
2125  if (arguments.size() == 1) {
2126  std::string cutString = arguments[0];
2127  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
2128  auto func = [cut](const Particle * particle) -> int {
2129  int n = 0;
2130  for (auto& daughter : particle->getDaughters())
2131  {
2132  if (cut->check(daughter))
2133  ++n;
2134  }
2135  return n;
2136  };
2137  return func;
2138  } else {
2139  B2FATAL("Wrong number of arguments for meta function countDaughters");
2140  }
2141  }
2142 
2143  Manager::FunctionPtr numberOfNonOverlappingParticles(const std::vector<std::string>& arguments)
2144  {
2145 
2146  auto func = [arguments](const Particle * particle) -> int {
2147 
2148  int _numberOfNonOverlappingParticles = 0;
2149  for (const auto& listName : arguments)
2150  {
2151  StoreObjPtr<ParticleList> list(listName);
2152  if (not list.isValid()) {
2153  B2FATAL("Invalid list named " << listName << " encountered in numberOfNonOverlappingParticles.");
2154  }
2155  for (unsigned int i = 0; i < list->getListSize(); i++) {
2156  const Particle* p = list->getParticle(i);
2157  if (not particle->overlapsWith(p)) {
2158  _numberOfNonOverlappingParticles++;
2159  }
2160  }
2161  }
2162  return _numberOfNonOverlappingParticles;
2163  };
2164 
2165  return func;
2166 
2167  }
2168 
2169  Manager::FunctionPtr matchedMC(const std::vector<std::string>& arguments)
2170  {
2171  if (arguments.size() == 1) {
2172  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
2173  auto func = [var](const Particle * particle) -> double {
2174  const MCParticle* mcp = particle->getMCParticle();
2175  if (!mcp) // Has no MC match and is no MCParticle
2176  {
2177  return std::numeric_limits<float>::quiet_NaN();
2178  }
2179  Particle tmpPart(mcp);
2180  auto var_result = var->function(&tmpPart);
2181  if (std::holds_alternative<double>(var_result))
2182  {
2183  return std::get<double>(var_result);
2184  } else if (std::holds_alternative<int>(var_result))
2185  {
2186  return std::get<int>(var_result);
2187  } else if (std::holds_alternative<bool>(var_result))
2188  {
2189  return std::get<bool>(var_result);
2190  } else return std::numeric_limits<double>::quiet_NaN();
2191  };
2192  return func;
2193  } else {
2194  B2FATAL("Wrong number of arguments for meta function matchedMC");
2195  }
2196  }
2197 
2198  double matchedMCHasPDG(const Particle* particle, const std::vector<double>& pdgCode)
2199  {
2200  if (pdgCode.size() != 1) {
2201  B2FATAL("Too many arguments provided to matchedMCHasPDG!");
2202  }
2203  int inputPDG = std::lround(pdgCode[0]);
2204 
2205  const MCParticle* mcp = particle->getMCParticle();
2206  if (!mcp)
2207  return std::numeric_limits<double>::quiet_NaN();
2208 
2209  return std::abs(mcp->getPDG()) == inputPDG;
2210  }
2211 
2212  Manager::FunctionPtr totalEnergyOfParticlesInList(const std::vector<std::string>& arguments)
2213  {
2214  if (arguments.size() == 1) {
2215  std::string listName = arguments[0];
2216  auto func = [listName](const Particle * particle) -> double {
2217 
2218  (void) particle;
2219  StoreObjPtr<ParticleList> listOfParticles(listName);
2220 
2221  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << listName << " given to totalEnergyOfParticlesInList");
2222  double totalEnergy = 0;
2223  int nParticles = listOfParticles->getListSize();
2224  for (int i = 0; i < nParticles; i++)
2225  {
2226  const Particle* part = listOfParticles->getParticle(i);
2227  const auto& frame = ReferenceFrame::GetCurrent();
2228  totalEnergy += frame.getMomentum(part).E();
2229  }
2230  return totalEnergy;
2231 
2232  };
2233  return func;
2234  } else {
2235  B2FATAL("Wrong number of arguments for meta function totalEnergyOfParticlesInList");
2236  }
2237  }
2238 
2239  Manager::FunctionPtr totalPxOfParticlesInList(const std::vector<std::string>& arguments)
2240  {
2241  if (arguments.size() == 1) {
2242  std::string listName = arguments[0];
2243  auto func = [listName](const Particle*) -> double {
2244  StoreObjPtr<ParticleList> listOfParticles(listName);
2245 
2246  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << listName << " given to totalPxOfParticlesInList");
2247  double totalPx = 0;
2248  int nParticles = listOfParticles->getListSize();
2249  const auto& frame = ReferenceFrame::GetCurrent();
2250  for (int i = 0; i < nParticles; i++)
2251  {
2252  const Particle* part = listOfParticles->getParticle(i);
2253  totalPx += frame.getMomentum(part).Px();
2254  }
2255  return totalPx;
2256  };
2257  return func;
2258  } else {
2259  B2FATAL("Wrong number of arguments for meta function totalPxOfParticlesInList");
2260  }
2261  }
2262 
2263  Manager::FunctionPtr totalPyOfParticlesInList(const std::vector<std::string>& arguments)
2264  {
2265  if (arguments.size() == 1) {
2266  std::string listName = arguments[0];
2267  auto func = [listName](const Particle*) -> double {
2268  StoreObjPtr<ParticleList> listOfParticles(listName);
2269 
2270  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << listName << " given to totalPyOfParticlesInList");
2271  double totalPy = 0;
2272  int nParticles = listOfParticles->getListSize();
2273  const auto& frame = ReferenceFrame::GetCurrent();
2274  for (int i = 0; i < nParticles; i++)
2275  {
2276  const Particle* part = listOfParticles->getParticle(i);
2277  totalPy += frame.getMomentum(part).Py();
2278  }
2279  return totalPy;
2280  };
2281  return func;
2282  } else {
2283  B2FATAL("Wrong number of arguments for meta function totalPyOfParticlesInList");
2284  }
2285  }
2286 
2287  Manager::FunctionPtr totalPzOfParticlesInList(const std::vector<std::string>& arguments)
2288  {
2289  if (arguments.size() == 1) {
2290  std::string listName = arguments[0];
2291  auto func = [listName](const Particle*) -> double {
2292  StoreObjPtr<ParticleList> listOfParticles(listName);
2293 
2294  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << listName << " given to totalPzOfParticlesInList");
2295  double totalPz = 0;
2296  int nParticles = listOfParticles->getListSize();
2297  const auto& frame = ReferenceFrame::GetCurrent();
2298  for (int i = 0; i < nParticles; i++)
2299  {
2300  const Particle* part = listOfParticles->getParticle(i);
2301  totalPz += frame.getMomentum(part).Pz();
2302  }
2303  return totalPz;
2304  };
2305  return func;
2306  } else {
2307  B2FATAL("Wrong number of arguments for meta function totalPzOfParticlesInList");
2308  }
2309  }
2310 
2311  Manager::FunctionPtr invMassInLists(const std::vector<std::string>& arguments)
2312  {
2313  if (arguments.size() > 0) {
2314 
2315  auto func = [arguments](const Particle * particle) -> double {
2316 
2317  ROOT::Math::PxPyPzEVector total4Vector;
2318  // To make sure particles in particlesList don't overlap.
2319  std::vector<Particle*> particlePool;
2320 
2321  (void) particle;
2322  for (const auto& argument : arguments)
2323  {
2324  StoreObjPtr <ParticleList> listOfParticles(argument);
2325 
2326  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << argument << " given to invMassInLists");
2327  int nParticles = listOfParticles->getListSize();
2328  for (int i = 0; i < nParticles; i++) {
2329  bool overlaps = false;
2330  Particle* part = listOfParticles->getParticle(i);
2331  for (auto poolPart : particlePool) {
2332  if (part->overlapsWith(poolPart)) {
2333  overlaps = true;
2334  break;
2335  }
2336  }
2337  if (!overlaps) {
2338  total4Vector += part->get4Vector();
2339  particlePool.push_back(part);
2340  }
2341  }
2342  }
2343  double invariantMass = total4Vector.M();
2344  return invariantMass;
2345 
2346  };
2347  return func;
2348  } else {
2349  B2FATAL("Wrong number of arguments for meta function invMassInLists");
2350  }
2351  }
2352 
2353  Manager::FunctionPtr totalECLEnergyOfParticlesInList(const std::vector<std::string>& arguments)
2354  {
2355  if (arguments.size() == 1) {
2356  std::string listName = arguments[0];
2357  auto func = [listName](const Particle * particle) -> double {
2358 
2359  (void) particle;
2360  StoreObjPtr<ParticleList> listOfParticles(listName);
2361 
2362  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << listName << " given to totalEnergyOfParticlesInList");
2363  double totalEnergy = 0;
2364  int nParticles = listOfParticles->getListSize();
2365  for (int i = 0; i < nParticles; i++)
2366  {
2367  const Particle* part = listOfParticles->getParticle(i);
2368  const ECLCluster* cluster = part->getECLCluster();
2369  const ECLCluster::EHypothesisBit clusterHypothesis = part->getECLClusterEHypothesisBit();
2370  if (cluster != nullptr) {
2371  totalEnergy += cluster->getEnergy(clusterHypothesis);
2372  }
2373  }
2374  return totalEnergy;
2375 
2376  };
2377  return func;
2378  } else {
2379  B2FATAL("Wrong number of arguments for meta function totalECLEnergyOfParticlesInList");
2380  }
2381  }
2382 
2383  Manager::FunctionPtr maxPtInList(const std::vector<std::string>& arguments)
2384  {
2385  if (arguments.size() == 1) {
2386  std::string listName = arguments[0];
2387  auto func = [listName](const Particle*) -> double {
2388  StoreObjPtr<ParticleList> listOfParticles(listName);
2389 
2390  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << listName << " given to maxPtInList");
2391  int nParticles = listOfParticles->getListSize();
2392  const auto& frame = ReferenceFrame::GetCurrent();
2393  double maxPt = 0;
2394  for (int i = 0; i < nParticles; i++)
2395  {
2396  const Particle* part = listOfParticles->getParticle(i);
2397  const double Pt = frame.getMomentum(part).Pt();
2398  if (Pt > maxPt) maxPt = Pt;
2399  }
2400  return maxPt;
2401  };
2402  return func;
2403  } else {
2404  B2FATAL("Wrong number of arguments for meta function maxPtInList");
2405  }
2406  }
2407 
2408  Manager::FunctionPtr eclClusterTrackMatchedWithCondition(const std::vector<std::string>& arguments)
2409  {
2410  if (arguments.size() <= 1) {
2411 
2412  std::string cutString;
2413  if (arguments.size() == 1)
2414  cutString = arguments[0];
2415  std::shared_ptr<Variable::Cut> cut = std::shared_ptr<Variable::Cut>(Variable::Cut::compile(cutString));
2416  auto func = [cut](const Particle * particle) -> double {
2417 
2418  if (particle == nullptr)
2419  return std::numeric_limits<double>::quiet_NaN();
2420 
2421  const ECLCluster* cluster = particle->getECLCluster();
2422 
2423  if (cluster)
2424  {
2425  auto tracks = cluster->getRelationsFrom<Track>();
2426 
2427  for (const auto& track : tracks) {
2428  Particle trackParticle(&track, Belle2::Const::pion);
2429 
2430  if (cut->check(&trackParticle))
2431  return 1;
2432  }
2433  return 0;
2434  }
2435  return std::numeric_limits<double>::quiet_NaN();
2436  };
2437  return func;
2438  } else {
2439  B2FATAL("Wrong number of arguments for meta function eclClusterSpecialTrackMatched");
2440  }
2441  }
2442 
2443  Manager::FunctionPtr averageValueInList(const std::vector<std::string>& arguments)
2444  {
2445  if (arguments.size() == 2) {
2446  std::string listName = arguments[0];
2447  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
2448 
2449  auto func = [listName, var](const Particle*) -> double {
2450  StoreObjPtr<ParticleList> listOfParticles(listName);
2451 
2452  if (!(listOfParticles.isValid())) B2FATAL("Invalid list name " << listName << " given to averageValueInList");
2453  int nParticles = listOfParticles->getListSize();
2454  if (nParticles == 0)
2455  {
2456  return std::numeric_limits<double>::quiet_NaN();
2457  }
2458  double average = 0;
2459  if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2460  {
2461  for (int i = 0; i < nParticles; i++) {
2462  average += std::get<double>(var->function(listOfParticles->getParticle(i))) / nParticles;
2463  }
2464  } else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2465  {
2466  for (int i = 0; i < nParticles; i++) {
2467  average += std::get<int>(var->function(listOfParticles->getParticle(i))) / nParticles;
2468  }
2469  } else return std::numeric_limits<double>::quiet_NaN();
2470  return average;
2471  };
2472  return func;
2473  } else {
2474  B2FATAL("Wrong number of arguments for meta function averageValueInList");
2475  }
2476  }
2477 
2478  Manager::FunctionPtr medianValueInList(const std::vector<std::string>& arguments)
2479  {
2480  if (arguments.size() == 2) {
2481  std::string listName = arguments[0];
2482  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
2483 
2484  auto func = [listName, var](const Particle*) -> double {
2485  StoreObjPtr<ParticleList> listOfParticles(listName);
2486 
2487  if (!(listOfParticles.isValid())) B2FATAL("Invalid list name " << listName << " given to medianValueInList");
2488  int nParticles = listOfParticles->getListSize();
2489  if (nParticles == 0)
2490  {
2491  return std::numeric_limits<double>::quiet_NaN();
2492  }
2493  std::vector<double> valuesInList;
2494  if (std::holds_alternative<double>(var->function(listOfParticles->getParticle(0))))
2495  {
2496  for (int i = 0; i < nParticles; i++) {
2497  valuesInList.push_back(std::get<double>(var->function(listOfParticles->getParticle(i))));
2498  }
2499  } else if (std::holds_alternative<int>(var->function(listOfParticles->getParticle(0))))
2500  {
2501  for (int i = 0; i < nParticles; i++) {
2502  valuesInList.push_back(std::get<int>(var->function(listOfParticles->getParticle(i))));
2503  }
2504  } else return std::numeric_limits<double>::quiet_NaN();
2505  std::sort(valuesInList.begin(), valuesInList.end());
2506  if (nParticles % 2 != 0)
2507  {
2508  return valuesInList[nParticles / 2];
2509  } else
2510  {
2511  return 0.5 * (valuesInList[nParticles / 2] + valuesInList[nParticles / 2 - 1]);
2512  }
2513  };
2514  return func;
2515  } else {
2516  B2FATAL("Wrong number of arguments for meta function medianValueInList");
2517  }
2518  }
2519 
2520  Manager::FunctionPtr angleToClosestInList(const std::vector<std::string>& arguments)
2521  {
2522  // expecting the list name
2523  if (arguments.size() != 1)
2524  B2FATAL("Wrong number of arguments for meta function angleToClosestInList");
2525 
2526  std::string listname = arguments[0];
2527 
2528  auto func = [listname](const Particle * particle) -> double {
2529  // get the list and check it's valid
2530  StoreObjPtr<ParticleList> list(listname);
2531  if (not list.isValid())
2532  B2FATAL("Invalid particle list name " << listname << " given to angleToClosestInList");
2533 
2534  // check the list isn't empty
2535  if (list->getListSize() == 0)
2536  return std::numeric_limits<double>::quiet_NaN();
2537 
2538  // respect the current frame and get the momentum of our input
2539  const auto& frame = ReferenceFrame::GetCurrent();
2540  const auto p_this = B2Vector3D(frame.getMomentum(particle).Vect());
2541 
2542  // find the particle index with the smallest opening angle
2543  double minAngle = 2 * M_PI;
2544  for (unsigned int i = 0; i < list->getListSize(); ++i)
2545  {
2546  const Particle* compareme = list->getParticle(i);
2547  const auto p_compare = B2Vector3D(frame.getMomentum(compareme).Vect());
2548  double angle = p_compare.Angle(p_this);
2549  if (minAngle > angle) minAngle = angle;
2550  }
2551  return minAngle;
2552  };
2553  return func;
2554  }
2555 
2556  Manager::FunctionPtr closestInList(const std::vector<std::string>& arguments)
2557  {
2558  // expecting the list name and a variable name
2559  if (arguments.size() != 2)
2560  B2FATAL("Wrong number of arguments for meta function closestInList");
2561 
2562  std::string listname = arguments[0];
2563 
2564  // the requested variable and check it exists
2565  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
2566 
2567  auto func = [listname, var](const Particle * particle) -> double {
2568  // get the list and check it's valid
2569  StoreObjPtr<ParticleList> list(listname);
2570  if (not list.isValid())
2571  B2FATAL("Invalid particle list name " << listname << " given to closestInList");
2572 
2573  // respect the current frame and get the momentum of our input
2574  const auto& frame = ReferenceFrame::GetCurrent();
2575  const auto p_this = B2Vector3D(frame.getMomentum(particle).Vect());
2576 
2577  // find the particle index with the smallest opening angle
2578  double minAngle = 2 * M_PI;
2579  int iClosest = -1;
2580  for (unsigned int i = 0; i < list->getListSize(); ++i)
2581  {
2582  const Particle* compareme = list->getParticle(i);
2583  const auto p_compare = B2Vector3D(frame.getMomentum(compareme).Vect());
2584  double angle = p_compare.Angle(p_this);
2585  if (minAngle > angle) {
2586  minAngle = angle;
2587  iClosest = i;
2588  }
2589  }
2590 
2591  // final check that the list wasn't empty (or some other problem)
2592  if (iClosest == -1) return std::numeric_limits<double>::quiet_NaN();
2593  auto var_result = var->function(list->getParticle(iClosest));
2594  if (std::holds_alternative<double>(var_result))
2595  {
2596  return std::get<double>(var_result);
2597  } else if (std::holds_alternative<int>(var_result))
2598  {
2599  return std::get<int>(var_result);
2600  } else if (std::holds_alternative<bool>(var_result))
2601  {
2602  return std::get<bool>(var_result);
2603  } else return std::numeric_limits<double>::quiet_NaN();
2604  };
2605  return func;
2606  }
2607 
2608  Manager::FunctionPtr angleToMostB2BInList(const std::vector<std::string>& arguments)
2609  {
2610  // expecting the list name
2611  if (arguments.size() != 1)
2612  B2FATAL("Wrong number of arguments for meta function angleToMostB2BInList");
2613 
2614  std::string listname = arguments[0];
2615 
2616  auto func = [listname](const Particle * particle) -> double {
2617  // get the list and check it's valid
2618  StoreObjPtr<ParticleList> list(listname);
2619  if (not list.isValid())
2620  B2FATAL("Invalid particle list name " << listname << " given to angleToMostB2BInList");
2621 
2622  // check the list isn't empty
2623  if (list->getListSize() == 0)
2624  return std::numeric_limits<double>::quiet_NaN();
2625 
2626  // respect the current frame and get the momentum of our input
2627  const auto& frame = ReferenceFrame::GetCurrent();
2628  const auto p_this = B2Vector3D(frame.getMomentum(particle).Vect());
2629 
2630  // find the most back-to-back (the largest opening angle before they
2631  // start getting smaller again!)
2632  double maxAngle = 0;
2633  for (unsigned int i = 0; i < list->getListSize(); ++i)
2634  {
2635  const Particle* compareme = list->getParticle(i);
2636  const auto p_compare = B2Vector3D(frame.getMomentum(compareme).Vect());
2637  double angle = p_compare.Angle(p_this);
2638  if (maxAngle < angle) maxAngle = angle;
2639  }
2640  return maxAngle;
2641  };
2642  return func;
2643  }
2644 
2645  Manager::FunctionPtr mostB2BInList(const std::vector<std::string>& arguments)
2646  {
2647  // expecting the list name and a variable name
2648  if (arguments.size() != 2)
2649  B2FATAL("Wrong number of arguments for meta function mostB2BInList");
2650 
2651  std::string listname = arguments[0];
2652 
2653  // the requested variable and check it exists
2654  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
2655 
2656  auto func = [listname, var](const Particle * particle) -> double {
2657  // get the list and check it's valid
2658  StoreObjPtr<ParticleList> list(listname);
2659  if (not list.isValid())
2660  B2FATAL("Invalid particle list name " << listname << " given to mostB2BInList");
2661 
2662  // respect the current frame and get the momentum of our input
2663  const auto& frame = ReferenceFrame::GetCurrent();
2664  const auto p_this = B2Vector3D(frame.getMomentum(particle).Vect());
2665 
2666  // find the most back-to-back (the largest opening angle before they
2667  // start getting smaller again!)
2668  double maxAngle = -1.0;
2669  int iMostB2B = -1;
2670  for (unsigned int i = 0; i < list->getListSize(); ++i)
2671  {
2672  const Particle* compareme = list->getParticle(i);
2673  const auto p_compare = B2Vector3D(frame.getMomentum(compareme).Vect());
2674  double angle = p_compare.Angle(p_this);
2675  if (maxAngle < angle) {
2676  maxAngle = angle;
2677  iMostB2B = i;
2678  }
2679  }
2680 
2681  // final check that the list wasn't empty (or some other problem)
2682  if (iMostB2B == -1) return std::numeric_limits<double>::quiet_NaN();
2683  auto var_result = var->function(list->getParticle(iMostB2B));
2684  if (std::holds_alternative<double>(var_result))
2685  {
2686  return std::get<double>(var_result);
2687  } else if (std::holds_alternative<int>(var_result))
2688  {
2689  return std::get<int>(var_result);
2690  } else if (std::holds_alternative<bool>(var_result))
2691  {
2692  return std::get<bool>(var_result);
2693  } else return std::numeric_limits<double>::quiet_NaN();
2694  };
2695  return func;
2696  }
2697 
2698  Manager::FunctionPtr maxOpeningAngleInList(const std::vector<std::string>& arguments)
2699  {
2700  if (arguments.size() == 1) {
2701  std::string listName = arguments[0];
2702  auto func = [listName](const Particle*) -> double {
2703  StoreObjPtr<ParticleList> listOfParticles(listName);
2704 
2705  if (!(listOfParticles.isValid())) B2FATAL("Invalid Listname " << listName << " given to maxOpeningAngleInList");
2706  int nParticles = listOfParticles->getListSize();
2707  // return NaN if number of particles is less than 2
2708  if (nParticles < 2) return std::numeric_limits<double>::quiet_NaN();
2709 
2710  const auto& frame = ReferenceFrame::GetCurrent();
2711  double maxOpeningAngle = -1;
2712  for (int i = 0; i < nParticles; i++)
2713  {
2714  B2Vector3D v1 = frame.getMomentum(listOfParticles->getParticle(i)).Vect();
2715  for (int j = i + 1; j < nParticles; j++) {
2716  B2Vector3D v2 = frame.getMomentum(listOfParticles->getParticle(j)).Vect();
2717  const double angle = v1.Angle(v2);
2718  if (angle > maxOpeningAngle) maxOpeningAngle = angle;
2719  }
2720  }
2721  return maxOpeningAngle;
2722  };
2723  return func;
2724  } else {
2725  B2FATAL("Wrong number of arguments for meta function maxOpeningAngleInList");
2726  }
2727  }
2728 
2729  Manager::FunctionPtr daughterCombination(const std::vector<std::string>& arguments)
2730  {
2731  // Expect 2 or more arguments.
2732  if (arguments.size() >= 2) {
2733  // First argument is the variable name
2734  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
2735 
2736  // Core function: calculates a variable combining an arbitrary number of particles
2737  auto func = [var, arguments](const Particle * particle) -> double {
2738  if (particle == nullptr)
2739  {
2740  B2WARNING("Trying to access a daughter that does not exist. Skipping");
2741  return std::numeric_limits<float>::quiet_NaN();
2742  }
2743  const auto& frame = ReferenceFrame::GetCurrent();
2744 
2745  // Sum of the 4-momenta of all the selected daughters
2746  ROOT::Math::PxPyPzEVector pSum(0, 0, 0, 0);
2747 
2748  // Loop over the arguments. Each one of them is a generalizedIndex,
2749  // pointing to a particle in the decay tree.
2750  for (unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++)
2751  {
2752  auto generalizedIndex = arguments[iCoord];
2753  const Particle* dauPart = particle->getParticleFromGeneralizedIndexString(generalizedIndex);
2754  if (dauPart)
2755  pSum += frame.getMomentum(dauPart);
2756  else {
2757  B2WARNING("Trying to access a daughter that does not exist. Index = " << generalizedIndex);
2758  return std::numeric_limits<float>::quiet_NaN();
2759  }
2760  }
2761 
2762  // Make a dummy particle out of the sum of the 4-momenta of the selected daughters
2763  Particle sumOfDaughters(pSum, 100); // 100 is one of the special numbers
2764 
2765  auto var_result = var->function(&sumOfDaughters);
2766  // Calculate the variable on the dummy particle
2767  if (std::holds_alternative<double>(var_result))
2768  {
2769  return std::get<double>(var_result);
2770  } else if (std::holds_alternative<int>(var_result))
2771  {
2772  return std::get<int>(var_result);
2773  } else if (std::holds_alternative<bool>(var_result))
2774  {
2775  return std::get<bool>(var_result);
2776  } else return std::numeric_limits<double>::quiet_NaN();
2777  };
2778  return func;
2779  } else
2780  B2FATAL("Wrong number of arguments for meta function daughterCombination");
2781  }
2782 
2783  Manager::FunctionPtr useAlternativeDaughterHypothesis(const std::vector<std::string>& arguments)
2784  {
2785  /*
2786  `arguments` contains the variable to calculate and a list of colon-separated index-particle pairs.
2787  Overall, it looks like {"M", "0:K+", "1:p+", "3:e-"}.
2788  The code is thus divided in two parts:
2789  1) Parsing. A loop over the elements of `arguments` that first separates the variable from the rest, and then splits all the index:particle
2790  pairs, filling a std::vector with the indexes and another one with the new mass values.
2791  2) Replacing: A loop over the particle's daughters. We take the 4-momentum of each of them, recalculating it with a new mass if needed, and then we calculate
2792  the variable value using the sum of all the 4-momenta, both updated and non-updated ones.
2793  */
2794 
2795  // Expect 2 or more arguments.
2796  if (arguments.size() >= 2) {
2797 
2798  //----
2799  // 1) parsing
2800  //----
2801 
2802  // First argument is the variable name
2803  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[0]);
2804 
2805  // Parses the other arguments, which are in the form of index:particleName pairs,
2806  // and stores indexes and pdgs in std::unordered_map
2807  std::unordered_map<unsigned int, int> mapOfReplacedDaughters;
2808 
2809  // Loop over the arguments to parse them
2810  for (unsigned int iCoord = 1; iCoord < arguments.size(); iCoord++) {
2811  auto replacedDauString = arguments[iCoord];
2812  // Split the string in index and new mass
2813  std::vector<std::string> indexAndMass;
2814  boost::split(indexAndMass, replacedDauString, boost::is_any_of(":"));
2815 
2816  // Checks that the index:particleName pair is properly formatted.
2817  if (indexAndMass.size() > 2) {
2818  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 "
2819  << replacedDauString << ", while a correct syntax looks like 0:K+.");
2820  return nullptr;
2821  }
2822 
2823  if (indexAndMass.size() < 2) {
2824  B2WARNING("The string indicating which daughter's mass should be replaced contains only one colon-separated element instead of two. The offending string is "
2825  << replacedDauString << ", while a correct syntax looks like 0:K+.");
2826  return nullptr;
2827  }
2828 
2829  // indexAndMass[0] is the daughter index as string. Try to convert it
2830  int dauIndex = 0;
2831  try {
2832  dauIndex = Belle2::convertString<int>(indexAndMass[0]);
2833  } catch (std::invalid_argument&) {
2834  B2FATAL("Found the string " << indexAndMass[0] << "instead of a daughter index.");
2835  }
2836 
2837  // Determine PDG code corresponding to indexAndMass[1] using the particle names defined in evt.pdl
2838  TParticlePDG* particlePDG = TDatabasePDG::Instance()->GetParticle(indexAndMass[1].c_str());
2839  if (!particlePDG) {
2840  B2WARNING("Particle not in evt.pdl file! " << indexAndMass[1]);
2841  return nullptr;
2842  }
2843 
2844  // Stores the indexes and the pdgs in the map that will be passed to the lambda function
2845  int pdgCode = particlePDG->PdgCode();
2846  mapOfReplacedDaughters[dauIndex] = pdgCode;
2847  } // End of parsing
2848 
2849  // Check the size of mapOfReplacedDaughters
2850  if (mapOfReplacedDaughters.size() != arguments.size() - 1)
2851  B2FATAL("Overlapped daughter's index is detected in the meta-variable useAlternativeDaughterHypothesis");
2852 
2853  //----
2854  // 2) replacing
2855  //----
2856 
2857  // Core function: creates a new particle from the original one changing
2858  // some of the daughters' masses
2859  auto func = [var, mapOfReplacedDaughters](const Particle * particle) -> double {
2860  if (particle == nullptr)
2861  {
2862  B2WARNING("Trying to access a particle that does not exist. Skipping");
2863  return std::numeric_limits<float>::quiet_NaN();
2864  }
2865 
2866  const auto& frame = ReferenceFrame::GetCurrent();
2867 
2868  // Create a dummy particle from the given particle to overwrite its kinematics
2869  Particle* dummy = ParticleCopy::copyParticle(particle);
2870 
2871  // Sum of the 4-momenta of all the daughters with the new mass assumptions
2872  ROOT::Math::PxPyPzMVector pSum(0, 0, 0, 0);
2873 
2874  for (unsigned int iDau = 0; iDau < particle->getNDaughters(); iDau++)
2875  {
2876  const Particle* dauPart = particle->getDaughter(iDau);
2877  if (not dauPart) {
2878  B2WARNING("Trying to access a daughter that does not exist. Index = " << iDau);
2879  return std::numeric_limits<float>::quiet_NaN();
2880  }
2881 
2882  ROOT::Math::PxPyPzMVector dauMom = ROOT::Math::PxPyPzMVector(frame.getMomentum(dauPart));
2883 
2884  int pdgCode;
2885  try {
2886  pdgCode = mapOfReplacedDaughters.at(iDau);
2887  } catch (std::out_of_range&) {
2888  // iDau is not in mapOfReplacedDaughters
2889  pSum += dauMom;
2890  continue;
2891  }
2892 
2893  // overwrite the daughter's kinematics
2894  double p_x = dauMom.Px();
2895  double p_y = dauMom.Py();
2896  double p_z = dauMom.Pz();
2897  dauMom.SetCoordinates(p_x, p_y, p_z, TDatabasePDG::Instance()->GetParticle(pdgCode)->Mass());
2898  const_cast<Particle*>(dummy->getDaughter(iDau))->set4VectorDividingByMomentumScaling(ROOT::Math::PxPyPzEVector(dauMom));
2899 
2900  // overwrite the daughter's pdg
2901  const int charge = dummy->getDaughter(iDau)->getCharge();
2902  if (TDatabasePDG::Instance()->GetParticle(pdgCode)->Charge() / 3.0 == charge)
2903  const_cast<Particle*>(dummy->getDaughter(iDau))->setPDGCode(pdgCode);
2904  else
2905  const_cast<Particle*>(dummy->getDaughter(iDau))->setPDGCode(-1 * pdgCode);
2906 
2907  pSum += dauMom;
2908  } // End of loop over number of daughter
2909 
2910  // overwrite the particle's kinematics
2911  dummy->set4Vector(ROOT::Math::PxPyPzEVector(pSum));
2912 
2913  auto var_result = var->function(dummy);
2914 
2915  // Calculate the variable on the dummy particle
2916  if (std::holds_alternative<double>(var_result))
2917  {
2918  return std::get<double>(var_result);
2919  } else if (std::holds_alternative<int>(var_result))
2920  {
2921  return std::get<int>(var_result);
2922  } else if (std::holds_alternative<bool>(var_result))
2923  {
2924  return std::get<bool>(var_result);
2925  } else return std::numeric_limits<double>::quiet_NaN();
2926  }; // end of lambda function
2927  return func;
2928  }// end of check on number of arguments
2929  else
2930  B2FATAL("Wrong number of arguments for meta function useAlternativeDaughterHypothesis");
2931  }
2932 
2933  Manager::FunctionPtr varForFirstMCAncestorOfType(const std::vector<std::string>& arguments)
2934  {
2935  if (arguments.size() == 2) {
2936  int pdg_code = -1;
2937  std::string arg = arguments[0];
2938  const Variable::Manager::Var* var = Manager::Instance().getVariable(arguments[1]);
2939  TParticlePDG* part = TDatabasePDG::Instance()->GetParticle(arg.c_str());
2940 
2941  if (part != nullptr) {
2942  pdg_code = std::abs(part->PdgCode());
2943  } else {
2944  try {
2945  pdg_code = Belle2::convertString<int>(arg);
2946  } catch (std::exception& e) {}
2947  }
2948 
2949  if (pdg_code == -1) {
2950  B2FATAL("Ancestor " + arg + " is not recognised. Please provide valid PDG code or particle name.");
2951  }
2952 
2953  auto func = [pdg_code, var](const Particle * particle) -> double {
2954  const Particle* p = particle;
2955 
2956  int ancestor_level = std::get<double>(Manager::Instance().getVariable("hasAncestor(" + std::to_string(pdg_code) + ", 0)")->function(p));
2957  if ((ancestor_level <= 0) or (std::isnan(ancestor_level)))
2958  {
2959  return std::numeric_limits<float>::quiet_NaN();
2960  }
2961 
2962  const MCParticle* i_p = p->getMCParticle();
2963 
2964  for (int a = 0; a < ancestor_level ; a = a + 1)
2965  {
2966  i_p = i_p->getMother();
2967  }
2968  Particle m_p(i_p);
2969  auto var_result = var->function(&m_p);
2970  if (std::holds_alternative<double>(var_result))
2971  {
2972  return std::get<double>(var_result);
2973  } else if (std::holds_alternative<int>(var_result))
2974  {
2975  return std::get<int>(var_result);
2976  } else if (std::holds_alternative<bool>(var_result))
2977  {
2978  return std::get<bool>(var_result);
2979  } else return std::numeric_limits<double>::quiet_NaN();
2980  };
2981  return func;
2982  } else {
2983  B2FATAL("Wrong number of arguments for meta function varForFirstMCAncestorOfType (expected 2: type and variable of interest)");
2984  }
2985  }
2986 
2987  VARIABLE_GROUP("MetaFunctions");
2988  REGISTER_METAVARIABLE("nCleanedECLClusters(cut)", nCleanedECLClusters,
2989  "[Eventbased] Returns the number of clean Clusters in the event\n"
2990  "Clean clusters are defined by the clusters which pass the given cut assuming a photon hypothesis.",
2991  Manager::VariableDataType::c_int);
2992  REGISTER_METAVARIABLE("nCleanedTracks(cut)", nCleanedTracks,
2993  "[Eventbased] Returns the number of clean Tracks in the event\n"
2994  "Clean tracks are defined by the tracks which pass the given cut assuming a pion hypothesis.", Manager::VariableDataType::c_int);
2995  REGISTER_METAVARIABLE("formula(v1 + v2 * [v3 - v4] / v5^v6)", formula, R"DOCSTRING(
2996 Returns the result of the given formula, where v1 to vN are variables or floating
2997 point numbers. Currently the only supported operations are addition (``+``),
2998 subtraction (``-``), multiplication (``*``), division (``/``) and power (``^``
2999 or ``**``). Parenthesis can be in the form of square brackets ``[v1 * v2]``
3000 or normal brackets ``(v1 * v2)``. It will work also with variables taking
3001 arguments. Operator precedence is taken into account. For example ::
3002 
3003  (daughter(0, E) + daughter(1, E))**2 - p**2 + 0.138
3004 
3005 .. versionchanged:: release-03-00-00
3006  now both, ``[]`` and ``()`` can be used for grouping operations, ``**`` can
3007  be used for exponent and float literals are possible directly in the
3008  formula.
3009 )DOCSTRING", Manager::VariableDataType::c_double);
3010  REGISTER_METAVARIABLE("useRestFrame(variable)", useRestFrame,
3011  "Returns the value of the variable using the rest frame of the given particle as current reference frame.\n"
3012  "E.g. ``useRestFrame(daughter(0, p))`` returns the total momentum of the first daughter in its mother's rest-frame", Manager::VariableDataType::c_double);
3013  REGISTER_METAVARIABLE("useCMSFrame(variable)", useCMSFrame,
3014  "Returns the value of the variable using the CMS frame as current reference frame.\n"
3015  "E.g. ``useCMSFrame(E)`` returns the energy of a particle in the CMS frame.", Manager::VariableDataType::c_double);
3016  REGISTER_METAVARIABLE("useLabFrame(variable)", useLabFrame, R"DOC(
3017 Returns the value of ``variable`` in the *lab* frame.
3018 
3019 .. tip::
3020  The lab frame is the default reference frame, usually you don't need to use this meta-variable.
3021  E.g. ``useLabFrame(E)`` returns the energy of a particle in the Lab frame, same as just ``E``.
3022 
3023 Specifying the lab frame is useful in some corner-cases. For example:
3024 ``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.
3025 )DOC", Manager::VariableDataType::c_double);
3026  REGISTER_METAVARIABLE("useTagSideRecoilRestFrame(variable, daughterIndexTagB)", useTagSideRecoilRestFrame,
3027  "Returns the value of the variable in the rest frame of the recoiling particle to the tag side B meson.\n"
3028  "The variable should only be applied to an Upsilon(4S) list.\n"
3029  "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);
3030  REGISTER_METAVARIABLE("useParticleRestFrame(variable, particleList)", useParticleRestFrame,
3031  "Returns the value of the variable in the rest frame of the first Particle contained in the given ParticleList.\n"
3032  "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3033  "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3034  "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);
3035  REGISTER_METAVARIABLE("useRecoilParticleRestFrame(variable, particleList)", useRecoilParticleRestFrame,
3036  "Returns the value of the variable in the rest frame of recoil system against the first Particle contained in the given ParticleList.\n"
3037  "It is strongly recommended to pass a ParticleList that contains at most only one Particle in each event. "
3038  "When more than one Particle is present in the ParticleList, only the first Particle in the list is used for "
3039  "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);
3040  REGISTER_METAVARIABLE("useDaughterRestFrame(variable, daughterIndex_1, [daughterIndex_2, ... daughterIndex_3])", useDaughterRestFrame,
3041  "Returns the value of the variable in the rest frame of the selected daughter particle.\n"
3042  "The daughter is identified via generalized daughter index, e.g. ``0:1`` identifies the second daughter (1) "
3043  "of the first daughter (0). If the daughter index is invalid, it returns NaN.\n"
3044  "If two or more indices are given, the rest frame of the sum of the daughters is used.",
3045  Manager::VariableDataType::c_double);
3046  REGISTER_METAVARIABLE("passesCut(cut)", passesCut,
3047  "Returns 1 if particle passes the cut otherwise 0.\n"
3048  "Useful if you want to write out if a particle would have passed a cut or not.", Manager::VariableDataType::c_bool);
3049  REGISTER_METAVARIABLE("passesEventCut(cut)", passesEventCut,
3050  "[Eventbased] Returns 1 if event passes the cut otherwise 0.\n"
3051  "Useful if you want to select events passing a cut without looping into particles, such as for skimming.\n", Manager::VariableDataType::c_bool);
3052  REGISTER_METAVARIABLE("countDaughters(cut)", countDaughters,
3053  "Returns number of direct daughters which satisfy the cut.\n"
3054  "Used by the skimming package (for what exactly?)", Manager::VariableDataType::c_int);
3055  REGISTER_METAVARIABLE("varFor(pdgCode, variable)", varFor,
3056  "Returns the value of the variable for the given particle if its abs(pdgCode) agrees with the given one.\n"
3057  "E.g. ``varFor(11, p)`` returns the momentum if the particle is an electron or a positron.", Manager::VariableDataType::c_double);
3058  REGISTER_METAVARIABLE("varForMCGen(variable)", varForMCGen,
3059  "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"
3060  "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"
3061  "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);
3062  REGISTER_METAVARIABLE("nParticlesInList(particleListName)", nParticlesInList,
3063  "[Eventbased] Returns number of particles in the given particle List.", Manager::VariableDataType::c_int);
3064  REGISTER_METAVARIABLE("isInList(particleListName)", isInList,
3065  "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);
3066  REGISTER_METAVARIABLE("isDaughterOfList(particleListNames)", isDaughterOfList,
3067  "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);
3068  REGISTER_METAVARIABLE("isDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isDescendantOfList, R"DOC(
3069  Returns 1 if the given particle appears in the decay chain of the particles in the given ParticleLists.
3070 
3071  Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3072 
3073  * ``isDescendantOfList(<particle_list>,1)`` returns 1 if particle is a daughter of the list,
3074  * ``isDescendantOfList(<particle_list>,2)`` returns 1 if particle is a granddaughter of the list,
3075  * ``isDescendantOfList(<particle_list>,3)`` returns 1 if particle is a great-granddaughter of the list, etc.
3076  * Default value is ``-1`` that is inclusive for all generations.
3077  )DOC", Manager::VariableDataType::c_bool);
3078  REGISTER_METAVARIABLE("isMCDescendantOfList(particleListName[, anotherParticleListName][, generationFlag = -1])", isMCDescendantOfList, R"DOC(
3079  Returns 1 if the given particle is linked to the same MC particle as any reconstructed daughter of the decay lists.
3080 
3081  Passing an integer as the last argument, allows to check if the particle belongs to the specific generation:
3082 
3083  * ``isMCDescendantOfList(<particle_list>,1)`` returns 1 if particle is matched to the same particle as any daughter of the list,
3084  * ``isMCDescendantOfList(<particle_list>,2)`` returns 1 if particle is matched to the same particle as any granddaughter of the list,
3085  * ``isMCDescendantOfList(<particle_list>,3)`` returns 1 if particle is matched to the same particle as any great-granddaughter of the list, etc.
3086  * Default value is ``-1`` that is inclusive for all generations.
3087 
3088  It makes only sense for lists created with `fillParticleListFromMC` function with ``addDaughters=True`` argument.
3089  )DOC", Manager::VariableDataType::c_bool);
3090 
3091  REGISTER_METAVARIABLE("sourceObjectIsInList(particleListName)", sourceObjectIsInList, R"DOC(
3092 Returns 1 if the underlying mdst object (e.g. track, or cluster) was used to create a particle in ``particleListName``, 0 if not.
3093 
3094 .. note::
3095  This only makes sense for particles that are not composite. Returns -1 for composite particles.
3096 )DOC", Manager::VariableDataType::c_int);
3097 
3098  REGISTER_METAVARIABLE("mcParticleIsInMCList(particleListName)", mcParticleIsInMCList, R"DOC(
3099 Returns 1 if the particle's matched MC particle is also matched to a particle in ``particleListName``
3100 (or if either of the lists were filled from generator level `modularAnalysis.fillParticleListFromMC`.)
3101 
3102 .. seealso:: :b2:var:`isMCDescendantOfList` to check daughters.
3103 )DOC", Manager::VariableDataType::c_bool);
3104 
3105  REGISTER_METAVARIABLE("isGrandDaughterOfList(particleListNames)", isGrandDaughterOfList,
3106  "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);
3107  REGISTER_METAVARIABLE("daughter(i, variable)", daughter, R"DOC(
3108  Returns value of variable for the i-th daughter. E.g.
3109 
3110  * ``daughter(0, p)`` returns the total momentum of the first daughter.
3111  * ``daughter(0, daughter(1, p)`` returns the total momentum of the second daughter of the first daughter.
3112 
3113  Returns NaN if particle is nullptr or if the given daughter-index is out of bound (>= amount of daughters).
3114  )DOC", Manager::VariableDataType::c_double);
3115  REGISTER_METAVARIABLE("mcDaughter(i, variable)", mcDaughter, R"DOC(
3116  Returns the value of the requested variable for the i-th Monte Carlo daughter of the particle.
3117 
3118  Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3119  or if the i-th MC daughter does not exist.
3120 
3121  E.g. ``mcDaughter(0, PDG)`` will return the PDG code of the first MC daughter of the matched MC
3122  particle of the reconstructed particle the function is applied to.
3123 
3124  The meta variable can also be nested: ``mcDaughter(0, mcDaughter(1, PDG))``.
3125  )DOC", Manager::VariableDataType::c_double);
3126  REGISTER_METAVARIABLE("mcMother(variable)", mcMother, R"DOC(
3127  Returns the value of the requested variable for the Monte Carlo mother of the particle.
3128 
3129  Returns NaN if the particle is nullptr, if the particle is not matched to an MC particle,
3130  or if the MC mother does not exist.
3131 
3132  E.g. ``mcMother(PDG)`` will return the PDG code of the MC mother of the matched MC
3133  particle of the reconstructed particle the function is applied to.
3134 
3135  The meta variable can also be nested: ``mcMother(mcMother(PDG))``.
3136  )DOC", Manager::VariableDataType::c_double);
3137  REGISTER_METAVARIABLE("genParticle(index, variable)", genParticle, R"DOC(
3138 [Eventbased] Returns the ``variable`` for the ith generator particle.
3139 The arguments of the function must be the ``index`` of the particle in the MCParticle Array,
3140 and ``variable``, the name of the function or variable for that generator particle.
3141 If ``index`` goes beyond the length of the MCParticles array, NaN will be returned.
3142 
3143 E.g. ``genParticle(0, p)`` returns the total momentum of the first MCParticle, which is
3144 the Upsilon(4S) in a generic decay.
3145 ``genParticle(0, mcDaughter(1, p)`` returns the total momentum of the second daughter of
3146 the first MC Particle, which is the momentum of the second B meson in a generic decay.
3147 )DOC", Manager::VariableDataType::c_double);
3148  REGISTER_METAVARIABLE("genUpsilon4S(variable)", genUpsilon4S, R"DOC(
3149 [Eventbased] Returns the ``variable`` evaluated for the generator-level :math:`\Upsilon(4S)`.
3150 If no generator level :math:`\Upsilon(4S)` exists for the event, NaN will be returned.
3151 
3152 E.g. ``genUpsilon4S(p)`` returns the total momentum of the :math:`\Upsilon(4S)` in a generic decay.
3153 ``genUpsilon4S(mcDaughter(1, p)`` returns the total momentum of the second daughter of the
3154 generator-level :math:`\Upsilon(4S)` (i.e. the momentum of the second B meson in a generic decay.
3155 )DOC", Manager::VariableDataType::c_double);
3156  REGISTER_METAVARIABLE("daughterProductOf(variable)", daughterProductOf,
3157  "Returns product of a variable over all daughters.\n"
3158  "E.g. ``daughterProductOf(extraInfo(SignalProbability))`` returns the product of the SignalProbabilitys of all daughters.", Manager::VariableDataType::c_double);
3159  REGISTER_METAVARIABLE("daughterSumOf(variable)", daughterSumOf,
3160  "Returns sum of a variable over all daughters.\n"
3161  "E.g. ``daughterSumOf(nDaughters)`` returns the number of grand-daughters.", Manager::VariableDataType::c_double);
3162  REGISTER_METAVARIABLE("daughterLowest(variable)", daughterLowest,
3163  "Returns the lowest value of the given variable among all daughters.\n"
3164  "E.g. ``useCMSFrame(daughterLowest(p))`` returns the lowest momentum in CMS frame.", Manager::VariableDataType::c_double);
3165  REGISTER_METAVARIABLE("daughterHighest(variable)", daughterHighest,
3166  "Returns the highest value of the given variable among all daughters.\n"
3167  "E.g. ``useCMSFrame(daughterHighest(p))`` returns the highest momentum in CMS frame.", Manager::VariableDataType::c_double);
3168  REGISTER_METAVARIABLE("daughterDiffOf(i, j, variable)", daughterDiffOf,
3169  "Returns the difference of a variable between the two given daughters.\n"
3170  "E.g. ``useRestFrame(daughterDiffOf(0, 1, p))`` returns the momentum difference between first and second daughter in the rest frame of the given particle.\n"
3171  "(That means that it returns :math:`p_j - p_i`)\n"
3172  "Nota Bene: for the particular case 'variable=phi' you should use the :b2:var:`daughterDiffOfPhi` function.", Manager::VariableDataType::c_double);
3173  REGISTER_METAVARIABLE("mcDaughterDiffOf(i, j, variable)", mcDaughterDiffOf,
3174  "MC matched version of the `daughterDiffOf` function.", Manager::VariableDataType::c_double);
3175  REGISTER_METAVARIABLE("grandDaughterDiffOf(i, j, variable)", grandDaughterDiffOf,
3176  "Returns the difference of a variable between the first daughters of the two given daughters.\n"
3177  "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"
3178  "(That means that it returns :math:`p_j - p_i`)\n"
3179  "Nota Bene: for the particular case 'variable=phi' you should use the :b2:var:`grandDaughterDiffOfPhi` function.", Manager::VariableDataType::c_double);
3180  REGISTER_METAVARIABLE("daughterDiffOfPhi(i, j)", daughterDiffOfPhi,
3181  "Returns the difference in :math:`\\phi` between the two given daughters. The unit of the angle is ``rad``.\n"
3182  "The difference is signed and takes account of the ordering of the given daughters.\n"
3183  "The function returns :math:`\\phi_j - \\phi_i`.", Manager::VariableDataType::c_double);
3184  MAKE_DEPRECATED("daughterDiffOfPhi", false, "release-06-00-00", R"DOC(
3185  The difference of the azimuthal angle :math:`\\phi` of two daughters can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3186  REGISTER_METAVARIABLE("mcDaughterDiffOfPhi(i, j)", mcDaughterDiffOfPhi,
3187  "MC matched version of the `daughterDiffOfPhi` function. The unit of the angle is ``rad``", Manager::VariableDataType::c_double);
3188  MAKE_DEPRECATED("mcDaughterDiffOfPhi", false, "release-06-00-00", R"DOC(
3189  The difference of the azimuthal angle :math:`\\phi` of the MC partners of two daughters can be calculated with the generic variable :b2:var:`mcDaughterDiffOf`.)DOC");
3190  REGISTER_METAVARIABLE("grandDaughterDiffOfPhi(i, j)", grandDaughterDiffOfPhi,
3191  "Returns the difference in :math:`\\phi` between the first daughters of the two given daughters. The unit of the angle is ``rad``.\n"
3192  "The difference is signed and takes account of the ordering of the given daughters.\n"
3193  "The function returns :math:`\\phi_j - \\phi_i`.\n", Manager::VariableDataType::c_double);
3194  MAKE_DEPRECATED("grandDaughterDiffOfPhi", false, "release-06-00-00", R"DOC(
3195  The difference of the azimuthal angle :math:`\\phi` of two granddaughters can be calculated with the generic variable :b2:var:`grandDaughterDiffOf`.)DOC");
3196  REGISTER_METAVARIABLE("daughterDiffOfClusterPhi(i, j)", daughterDiffOfClusterPhi,
3197  "Returns the difference in :math:`\\phi` between the ECLClusters of two given daughters. The unit of the angle is ``rad``.\n"
3198  "The difference is signed and takes account of the ordering of the given daughters.\n"
3199  "The function returns :math:`\\phi_j - \\phi_i`.\n"
3200  "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.", Manager::VariableDataType::c_double);
3201  MAKE_DEPRECATED("daughterDiffOfClusterPhi", false, "release-06-00-00", R"DOC(
3202  The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two daughters can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3203  REGISTER_METAVARIABLE("grandDaughterDiffOfClusterPhi(i, j)", grandDaughterDiffOfClusterPhi,
3204  "Returns the difference in :math:`\\phi` between the ECLClusters of the daughters of the two given daughters. The unit of the angle is ``rad``.\n"
3205  "The difference is signed and takes account of the ordering of the given daughters.\n"
3206  "The function returns :math:`\\phi_j - \\phi_i`.\n"
3207  "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.\n", Manager::VariableDataType::c_double);
3208  MAKE_DEPRECATED("grandDaughterDiffOfClusterPhi", false, "release-06-00-00", R"DOC(
3209  The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two granddaughters can be calculated with the generic variable :b2:var:`grandDaughterDiffOf`.)DOC");
3210  REGISTER_METAVARIABLE("daughterDiffOfPhiCMS(i, j)", daughterDiffOfPhiCMS,
3211  "Returns the difference in :math:`\\phi` between the two given daughters in the CMS frame. The unit of the angle is ``rad``.\n"
3212  "The difference is signed and takes account of the ordering of the given daughters.\n"
3213  "The function returns :math:`\\phi_j - \\phi_i`.", Manager::VariableDataType::c_double);
3214  MAKE_DEPRECATED("daughterDiffOfPhiCMS", false, "release-06-00-00", R"DOC(
3215  The difference of the azimuthal angle :math:`\\phi` of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3216  REGISTER_METAVARIABLE("mcDaughterDiffOfPhiCMS(i, j)", daughterDiffOfPhiCMS,
3217  "MC matched version of the `daughterDiffOfPhiCMS` function. The unit of the angle is ``rad``", Manager::VariableDataType::c_double);
3218  MAKE_DEPRECATED("mcDaughterDiffOfPhiCMS", false, "release-06-00-00", R"DOC(
3219  The difference of the azimuthal angle :math:`\\phi` of the MC partners of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`mcDaughterDiffOf`.)DOC");
3220  REGISTER_METAVARIABLE("daughterDiffOfClusterPhiCMS(i, j)", daughterDiffOfClusterPhiCMS,
3221  "Returns the difference in :math:`\\phi` between the ECLClusters of two given daughters in the CMS frame. The unit of the angle is ``rad``.\n"
3222  "The difference is signed and takes account of the ordering of the given daughters.\n"
3223  "The function returns :math:`\\phi_j - \\phi_i``.\n"
3224  "The function returns NaN if at least one of the daughters is not matched to or not based on an ECLCluster.", Manager::VariableDataType::c_double);
3225  MAKE_DEPRECATED("daughterDiffOfClusterPhiCMS", false, "release-06-00-00", R"DOC(
3226  The difference of the azimuthal angle :math:`\\phi` of the related ECL clusters of two daughters in the CMS frame can be calculated with the generic variable :b2:var:`daughterDiffOf`.)DOC");
3227  REGISTER_METAVARIABLE("daughterNormDiffOf(i, j, variable)", daughterNormDiffOf,
3228  "Returns the normalized difference of a variable between the two given daughters.\n"
3229  "E.g. ``daughterNormDiffOf(0, 1, p)`` returns the normalized momentum difference between first and second daughter in the lab frame.", Manager::VariableDataType::c_double);
3230  REGISTER_METAVARIABLE("daughterMotherDiffOf(i, variable)", daughterMotherDiffOf,
3231  "Returns the difference of a variable between the given daughter and the mother particle itself.\n"
3232  "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);
3233  REGISTER_METAVARIABLE("daughterMotherNormDiffOf(i, variable)", daughterMotherNormDiffOf,
3234  "Returns the normalized difference of a variable between the given daughter and the mother particle itself.\n"
3235  "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);
3236  REGISTER_METAVARIABLE("daughterAngle(daughterIndex_1, daughterIndex_2[, daughterIndex_3])", daughterAngle, R"DOC(
3237  Returns the angle in between any pair of particles belonging to the same decay tree.
3238  The unit of the angle is ``rad``.
3239 
3240  The particles are identified via generalized daughter indexes, which are simply colon-separated lists of
3241  daughter indexes, ordered starting from the root particle. For example, ``0:1:3`` identifies the fourth
3242  daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle. ``1`` simply
3243  identifies the second daughter of the root particle.
3244 
3245  Both two and three generalized indexes can be given to ``daughterAngle``. If two indices are given, the
3246  variable returns the angle between the momenta of the two given particles. If three indices are given, the
3247  variable returns the angle between the momentum of the third particle and a vector which is the sum of the
3248  first two daughter momenta.
3249 
3250  .. tip::
3251  ``daughterAngle(0, 3)`` will return the angle between the first and fourth daughter.
3252  ``daughterAngle(0, 1, 3)`` will return the angle between the fourth daughter and the sum of the first and
3253  second daughter.
3254  ``daughterAngle(0:0, 3:0)`` will return the angle between the first daughter of the first daughter, and
3255  the first daughter of the fourth daughter.
3256 
3257  )DOC", Manager::VariableDataType::c_double);
3258  REGISTER_METAVARIABLE("mcDaughterAngle(daughterIndex_1, daughterIndex_2, [daughterIndex_3])", mcDaughterAngle,
3259  "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);
3260  REGISTER_VARIABLE("grandDaughterDecayAngle(i, j)", grandDaughterDecayAngle,
3261  "Returns the decay angle of the granddaughter in the daughter particle's rest frame.\n"
3262  "It is calculated with respect to the reverted momentum vector of the particle.\n"
3263  "Two arguments representing the daughter and granddaughter indices have to be provided as arguments.", "rad");
3264  REGISTER_VARIABLE("daughterClusterAngleInBetween(i, j)", daughterClusterAngleInBetween,
3265  "Returns the angle between clusters associated to the two daughters."
3266  "If two indices given: returns the angle between the momenta of the clusters associated to the two given daughters."
3267  "If three indices given: returns the angle between the momentum of the third particle's cluster and a vector "
3268  "which is the sum of the first two daughter's cluster momenta."
3269  "Returns nan if any of the daughters specified don't have an associated cluster."
3270  "The arguments in the argument vector must be integers corresponding to the ith and jth (and kth) daughters.", "rad");
3271  REGISTER_METAVARIABLE("daughterInvM(i[, j, ...])", daughterInvM, R"DOC(
3272  Returns the invariant mass adding the Lorentz vectors of the given daughters. The unit of the invariant mass is GeV/:math:`\text{c}^2`
3273  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.
3274 
3275  Daughters from different generations of the decay tree can be combined using generalized daughter indexes,
3276  which are simply colon-separated daughter indexes for each generation, starting from the root particle. For
3277  example, ``0:1:3`` identifies the fourth daughter (3) of the second daughter (1) of the first daughter(0) of
3278  the mother particle.
3279 
3280  Returns NaN if the given daughter-index is out of bound (>= number of daughters))DOC", Manager::VariableDataType::c_double);
3281  REGISTER_METAVARIABLE("extraInfo(name)", extraInfo,
3282  "Returns extra info stored under the given name.\n"
3283  "The extraInfo has to be set by a module first.\n"
3284  "E.g. ``extraInfo(SignalProbability)`` returns the SignalProbability calculated by the ``MVAExpert`` module.\n"
3285  "If nothing is set under the given name or if the particle is a nullptr, NaN is returned.\n"
3286  "In the latter case please use `eventExtraInfo` if you want to access an EventExtraInfo variable.", Manager::VariableDataType::c_double);
3287  REGISTER_METAVARIABLE("eventExtraInfo(name)", eventExtraInfo,
3288  "[Eventbased] Returns extra info stored under the given name in the event extra info.\n"
3289  "The extraInfo has to be set first by another module like MVAExpert in event mode.\n"
3290  "If nothing is set under this name, NaN is returned.", Manager::VariableDataType::c_double);
3291  REGISTER_METAVARIABLE("eventCached(variable)", eventCached,
3292  "[Eventbased] Returns value of event-based variable and caches this value in the EventExtraInfo.\n"
3293  "The result of second call to this variable in the same event will be provided from the cache.\n"
3294  "It is recommended to use this variable in order to declare custom aliases as event-based. This is "
3295  "necessary if using the eventwise mode of variablesToNtuple).", Manager::VariableDataType::c_double);
3296  REGISTER_METAVARIABLE("particleCached(variable)", particleCached,
3297  "Returns value of given variable and caches this value in the ParticleExtraInfo of the provided particle.\n"
3298  "The result of second call to this variable on the same particle will be provided from the cache.", Manager::VariableDataType::c_double);
3299  REGISTER_METAVARIABLE("modulo(variable, n)", modulo,
3300  "Returns rest of division of variable by n.", Manager::VariableDataType::c_int);
3301  REGISTER_METAVARIABLE("abs(variable)", abs,
3302  "Returns absolute value of the given variable.\n"
3303  "E.g. abs(mcPDG) returns the absolute value of the mcPDG, which is often useful for cuts.", Manager::VariableDataType::c_double);
3304  REGISTER_METAVARIABLE("max(var1,var2)", max, "Returns max value of two variables.\n", Manager::VariableDataType::c_double);
3305  REGISTER_METAVARIABLE("min(var1,var2)", min, "Returns min value of two variables.\n", Manager::VariableDataType::c_double);
3306  REGISTER_METAVARIABLE("sin(variable)", sin, "Returns sine value of the given variable.", Manager::VariableDataType::c_double);
3307  REGISTER_METAVARIABLE("asin(variable)", asin, "Returns arcsine of the given variable. The unit of the asin() is ``rad``", Manager::VariableDataType::c_double);
3308  REGISTER_METAVARIABLE("cos(variable)", cos, "Returns cosine value of the given variable.", Manager::VariableDataType::c_double);
3309  REGISTER_METAVARIABLE("acos(variable)", acos, "Returns arccosine value of the given variable. The unit of the acos() is ``rad``", Manager::VariableDataType::c_double);
3310  REGISTER_METAVARIABLE("tan(variable)", tan, "Returns tangent value of the given variable.", Manager::VariableDataType::c_double);
3311  REGISTER_METAVARIABLE("atan(variable)", atan, "Returns arctangent value of the given variable. The unit of the atan() is ``rad``", Manager::VariableDataType::c_double);
3312  REGISTER_METAVARIABLE("exp(variable)", exp, "Returns exponential evaluated for the given variable.", Manager::VariableDataType::c_double);
3313  REGISTER_METAVARIABLE("log(variable)", log, "Returns natural logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3314  REGISTER_METAVARIABLE("log10(variable)", log10, "Returns base-10 logarithm evaluated for the given variable.", Manager::VariableDataType::c_double);
3315  REGISTER_METAVARIABLE("isNAN(variable)", isNAN,
3316  "Returns true if variable value evaluates to nan (determined via std::isnan(double)).\n"
3317  "Useful for debugging.", Manager::VariableDataType::c_bool);
3318  REGISTER_METAVARIABLE("ifNANgiveX(variable, x)", ifNANgiveX,
3319  "Returns x (has to be a number) if variable value is nan (determined via std::isnan(double)).\n"
3320  "Useful for technical purposes while training MVAs.", Manager::VariableDataType::c_double);
3321  REGISTER_METAVARIABLE("isInfinity(variable)", isInfinity,
3322  "Returns true if variable value evaluates to infinity (determined via std::isinf(double)).\n"
3323  "Useful for debugging.", Manager::VariableDataType::c_bool);
3324  REGISTER_METAVARIABLE("unmask(variable, flag1, flag2, ...)", unmask,
3325  "unmask(variable, flag1, flag2, ...) or unmask(variable, mask) sets certain bits in the variable to zero.\n"
3326  "For example, if you want to set the second, fourth and fifth bits to zero, you could call \n"
3327  "``unmask(variable, 2, 8, 16)`` or ``unmask(variable, 26)``.\n"
3328  "", Manager::VariableDataType::c_double);
3329  REGISTER_METAVARIABLE("conditionalVariableSelector(cut, variableIfTrue, variableIfFalse)", conditionalVariableSelector,
3330  "Returns one of the two supplied variables, depending on whether the particle passes the supplied cut.\n"
3331  "The first variable is returned if the particle passes the cut, and the second variable is returned otherwise.", Manager::VariableDataType::c_double);
3332  REGISTER_METAVARIABLE("pValueCombination(p1, p2, ...)", pValueCombination,
3333  "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"
3334  "If any of the p-values is invalid, i.e. smaller than zero, -1 is returned.", Manager::VariableDataType::c_double);
3335  REGISTER_METAVARIABLE("veto(particleList, cut, pdgCode = 11)", veto,
3336  "Combines current particle with particles from the given particle list and returns 1 if the combination passes the provided cut. \n"
3337  "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"
3338  "around the neutral Pion mass (e.g. ``0.130 < M < 0.140``). \n"
3339  "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);
3340  REGISTER_METAVARIABLE("matchedMC(variable)", matchedMC,
3341  "Returns variable output for the matched MCParticle by constructing a temporary Particle from it.\n"
3342  "This may not work too well if your variable requires accessing daughters of the particle.\n"
3343  "E.g. ``matchedMC(p)`` returns the total momentum of the related MCParticle.\n"
3344  "Returns NaN if no matched MCParticle exists.", Manager::VariableDataType::c_double);
3345  REGISTER_METAVARIABLE("countInList(particleList, cut='')", countInList, "[Eventbased] "
3346  "Returns number of particle which pass given in cut in the specified particle list.\n"
3347  "Useful for creating statistics about the number of particles in a list.\n"
3348  "E.g. ``countInList(e+, isSignal == 1)`` returns the number of correctly reconstructed electrons in the event.\n"
3349  "The variable is event-based and does not need a valid particle pointer as input.", Manager::VariableDataType::c_int);
3350  REGISTER_METAVARIABLE("getVariableByRank(particleList, rankedVariableName, variableName, rank)", getVariableByRank, R"DOC(
3351  [Eventbased] Returns the value of ``variableName`` for the candidate in the ``particleList`` with the requested ``rank``.
3352 
3353  .. note::
3354  The `BestCandidateSelection` module available via `rankByHighest` / `rankByLowest` has to be used before.
3355 
3356  .. warning::
3357  The first candidate matching the given rank is used.
3358  Thus, it is not recommended to use this variable in conjunction with ``allowMultiRank`` in the `BestCandidateSelection` module.
3359 
3360  The suffix ``_rank`` is automatically added to the argument ``rankedVariableName``,
3361  which either has to be the name of the variable used to order the candidates or the selected outputVariable name without the ending ``_rank``.
3362  This means that your selected name for the rank variable has to end with ``_rank``.
3363 
3364  An example of this variable's usage is given in the tutorial `B2A602-BestCandidateSelection <https://stash.desy.de/projects/B2/repos/basf2/browse/analysis/examples/tutorials/B2A602-BestCandidateSelection.py>`_
3365  )DOC", Manager::VariableDataType::c_double);
3366  REGISTER_VARIABLE("matchedMCHasPDG(PDGCode)", matchedMCHasPDG,
3367  "Returns if the absolute value of the PDGCode of the MCParticle related to the Particle matches a given PDGCode."
3368  "Returns 0/NAN/1 if PDGCode does not match/is not available/ matches");
3369  REGISTER_METAVARIABLE("numberOfNonOverlappingParticles(pList1, pList2, ...)", numberOfNonOverlappingParticles,
3370  "Returns the number of non-overlapping particles in the given particle lists"
3371  "Useful to check if there is additional physics going on in the detector if one reconstructed the Y4S", Manager::VariableDataType::c_int);
3372  REGISTER_METAVARIABLE("totalEnergyOfParticlesInList(particleListName)", totalEnergyOfParticlesInList,
3373  "Returns the total energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3374  REGISTER_METAVARIABLE("totalPxOfParticlesInList(particleListName)", totalPxOfParticlesInList,
3375  "Returns the total momentum Px of particles in the given particle List. The unit of the momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3376  REGISTER_METAVARIABLE("totalPyOfParticlesInList(particleListName)", totalPyOfParticlesInList,
3377  "Returns the total momentum Py of particles in the given particle List. The unit of the momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3378  REGISTER_METAVARIABLE("totalPzOfParticlesInList(particleListName)", totalPzOfParticlesInList,
3379  "Returns the total momentum Pz of particles in the given particle List. The unit of the momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3380  REGISTER_METAVARIABLE("invMassInLists(pList1, pList2, ...)", invMassInLists,
3381  "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);
3382  REGISTER_METAVARIABLE("totalECLEnergyOfParticlesInList(particleListName)", totalECLEnergyOfParticlesInList,
3383  "Returns the total ECL energy of particles in the given particle List. The unit of the energy is ``GeV``", Manager::VariableDataType::c_double);
3384  REGISTER_METAVARIABLE("maxPtInList(particleListName)", maxPtInList,
3385  "Returns maximum transverse momentum Pt in the given particle List. The unit of the transverse momentum is ``GeV/c``", Manager::VariableDataType::c_double);
3386  REGISTER_METAVARIABLE("eclClusterSpecialTrackMatched(cut)", eclClusterTrackMatchedWithCondition,
3387  "Returns if at least one Track that satisfies the given condition is related to the ECLCluster of the Particle.", Manager::VariableDataType::c_double);
3388  REGISTER_METAVARIABLE("averageValueInList(particleListName, variable)", averageValueInList,
3389  "Returns the arithmetic mean of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3390  REGISTER_METAVARIABLE("medianValueInList(particleListName, variable)", medianValueInList,
3391  "Returns the median value of the given variable of the particles in the given particle list.", Manager::VariableDataType::c_double);
3392  REGISTER_METAVARIABLE("angleToClosestInList(particleListName)", angleToClosestInList,
3393  "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);
3394  REGISTER_METAVARIABLE("closestInList(particleListName, variable)", closestInList,
3395  "Returns `variable` for the closest particle (smallest opening angle) in the list provided.", Manager::VariableDataType::c_double);
3396  REGISTER_METAVARIABLE("angleToMostB2BInList(particleListName)", angleToMostB2BInList,
3397  "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);
3398  REGISTER_METAVARIABLE("mostB2BInList(particleListName, variable)", mostB2BInList,
3399  "Returns `variable` for the most back-to-back particle (closest opening angle to 180) in the list provided.", Manager::VariableDataType::c_double);
3400  REGISTER_METAVARIABLE("maxOpeningAngleInList(particleListName)", maxOpeningAngleInList,
3401  "Returns maximum opening angle in the given particle List. The unit of the angle is ``rad`` ", Manager::VariableDataType::c_double);
3402  REGISTER_METAVARIABLE("daughterCombination(variable, daughterIndex_1, daughterIndex_2 ... daughterIndex_n)", daughterCombination,R"DOC(
3403 Returns a ``variable`` function only of the 4-momentum calculated on an arbitrary set of (grand)daughters.
3404 
3405 .. warning::
3406  ``variable`` can only be a function of the daughters' 4-momenta.
3407 
3408 Daughters from different generations of the decay tree can be combined using generalized daughter indexes, which are simply colon-separated
3409 the list of daughter indexes, starting from the root particle: for example, ``0:1:3`` identifies the fourth
3410 daughter (3) of the second daughter (1) of the first daughter (0) of the mother particle.
3411 
3412 .. tip::
3413  ``daughterCombination(M, 0, 3, 4)`` will return the invariant mass of the system made of the first, fourth and fifth daughter of particle.
3414  ``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.
3415 
3416 )DOC", Manager::VariableDataType::c_double);
3417  REGISTER_METAVARIABLE("useAlternativeDaughterHypothesis(variable, daughterIndex_1:newMassHyp_1, ..., daughterIndex_n:newMassHyp_n)", useAlternativeDaughterHypothesis,R"DOC(
3418 Returns a ``variable`` calculated using new mass hypotheses for (some of) the particle's daughters.
3419 
3420 .. warning::
3421  ``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.
3422  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.
3423  Also, the track fit is not performed again: the variable only re-calculates the 4-vectors using different mass assumptions.
3424  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).
3425 
3426 .. warning::
3427  Generalized daughter indexes are not supported (yet!): this variable can be used only on first-generation daughters.
3428 
3429 .. tip::
3430  ``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.
3431  ``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.
3432 
3433 )DOC", Manager::VariableDataType::c_double);
3434  REGISTER_METAVARIABLE("varForFirstMCAncestorOfType(type, variable)",varForFirstMCAncestorOfType,R"DOC(Returns requested variable of the first ancestor of the given type.
3435 Ancestor type can be set up by PDG code or by particle name (check evt.pdl for valid particle names))DOC", Manager::VariableDataType::c_double)
3436 
3437  }
3439 }
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
Definition: B2Vector3.h:302
int getPDGCode() const
PDG code.
Definition: Const.h:464
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
static const ChargedStable electron
electron particle
Definition: Const.h:650
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
Definition: ECLCluster.h:31
@ c_nPhotons
CR is split into n photons (N1)
static std::unique_ptr< GeneralCut > compile(const std::string &cut)
Creates an instance of a cut and returns a unique_ptr to it, if you need a copy-able object instead y...
Definition: GeneralCut.h:84
@ c_Initial
bit 5: Particle is initial such as e+ or e- and not going to Geant4
Definition: MCParticle.h:57
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
@ c_IsVirtual
bit 4: Particle is virtual and not going to Geant4.
Definition: MCParticle.h:55
static std::string makeROOTCompatible(std::string str)
Remove special characters that ROOT dislikes in branch names, e.g.
const MCParticle * getMCParticle() const
Returns the pointer to the MCParticle object that was used to create this Particle (ParticleType == c...
Definition: Particle.cc:972
EParticleSourceObject
particle source enumerators
Definition: Particle.h:81
const Particle * getParticleFromGeneralizedIndexString(const std::string &generalizedIndex) const
Explores the decay tree of the particle and returns the (grand^n)daughter identified by a generalized...
Definition: Particle.cc:986
@ c_Unflavored
Is its own antiparticle or we don't know whether it is a particle/antiparticle.
Definition: Particle.h:94
@ c_Flavored
Is either particle or antiparticle.
Definition: Particle.h:95
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:663
static const ReferenceFrame & GetCurrent()
Get current rest frame.
std::function< VarVariant(const Particle *)> FunctionPtr
functions stored take a const Particle* and return VarVariant.
Definition: Manager.h:113
const Var * getVariable(std::string name)
Get the variable belonging to the given key.
Definition: Manager.cc:57
static Manager & Instance()
get singleton instance.
Definition: Manager.cc:25
Class to store variables with their name which were sent to the logging service.
#define MAKE_DEPRECATED(name, make_fatal, version, description)
Registers a variable as deprecated.
Definition: Manager.h:443
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
Definition: EvtPDLUtil.cc:44
bool hasAntiParticle(int pdgCode)
Checks if the particle with given pdg code has an anti-particle or not.
Definition: EvtPDLUtil.cc:12
Particle * copyParticle(const Particle *original)
Function takes argument Particle and creates a copy of it and copies of all its (grand-)^n-daughters.
Definition: ParticleCopy.cc:18
Abstract base class for different kinds of events.
Definition: ClusterUtils.h:23