10 #include <analysis/variables/EventShapeVariables.h>
12 #include <analysis/dataobjects/Particle.h>
13 #include <analysis/dataobjects/EventShapeContainer.h>
15 #include <analysis/utility/ReferenceFrame.h>
17 #include <framework/logging/Logger.h>
18 #include <framework/datastore/StoreObjPtr.h>
19 #include <framework/utilities/Conversion.h>
23 #include <boost/algorithm/string.hpp>
32 double foxWolframR(
const Particle*,
const std::vector<double>& index)
34 if (index.size() != 1) {
35 B2ERROR(
"foxWolframR cannot be called without providing the moment order");
36 return std::numeric_limits<float>::quiet_NaN();
39 int order = std::lround(index[0]);
41 if (order < 0 || order > 8) {
42 B2ERROR(
"The Fox-Wolfram moment order must be within 0 and 8.");
43 return std::numeric_limits<float>::quiet_NaN();
46 StoreObjPtr<EventShapeContainer> evtShapeCont;
48 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
49 return std::numeric_limits<float>::quiet_NaN();
51 if (evtShapeCont->getFWMoment(0) == 0) {
52 B2ERROR(
"The 0th-order FoxWolfram moment is zero");
53 return std::numeric_limits<float>::quiet_NaN();
55 return evtShapeCont->getFWMoment(order) / evtShapeCont->getFWMoment(0);
58 double foxWolframH(
const Particle*,
const std::vector<double>& index)
60 if (index.size() != 1) {
61 B2ERROR(
"foxWolframH cannot be called without providing the moment order");
62 return std::numeric_limits<float>::quiet_NaN();
65 int order = std::lround(index[0]);
67 if (order < 0 || order > 8) {
68 B2ERROR(
"The Fox-Wolfram moment order must be within 0 and 8.");
69 return std::numeric_limits<float>::quiet_NaN();
72 StoreObjPtr<EventShapeContainer> evtShapeCont;
74 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
75 return std::numeric_limits<float>::quiet_NaN();
77 return evtShapeCont->getFWMoment(order);
82 if (arguments.size() != 2) {
83 B2ERROR(
"harmonicMoment requires two arguments: the harmonic order (0-8) and the reference axis name (thrust or collision)");
89 order = Belle2::convertString<int>(arguments[0]);
90 }
catch (std::invalid_argument&) {
91 B2ERROR(
"First argument of harmonicMoment must be an integer");
94 std::string axisName = arguments[1];
95 boost::to_lower(axisName);
97 if (order < 0 || order > 8) {
98 B2ERROR(
"The Fox-Wolfram moment order must be within 0 and 8.");
101 if (axisName !=
"thrust" && axisName !=
"collision") {
102 B2ERROR(
"Invalid axis name " << arguments[1] <<
". The valid options are thrust and collision");
106 auto func = [order, axisName](
const Particle*) ->
double{
107 StoreObjPtr<EventShapeContainer> evtShapeCont;
110 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
111 return std::numeric_limits<float>::quiet_NaN();
113 if (axisName ==
"thrust")
114 return evtShapeCont->getHarmonicMomentThrust(order);
116 return evtShapeCont->getHarmonicMomentCollision(order);
123 if (arguments.size() != 2) {
124 B2ERROR(
"cleoCone requires two arguments: the cone order (0-9) and the reference axis name (thrust or collision)");
130 order = Belle2::convertString<int>(arguments[0]);
131 }
catch (std::invalid_argument&) {
132 B2ERROR(
"Argument of cleoCone must be an integer");
135 std::string axisName = arguments[1];
136 boost::to_lower(axisName);
138 if (order < 0 || order > 8) {
139 B2ERROR(
"The CLEO cone order must be within 0 and 8.");
142 if (axisName !=
"thrust" && axisName !=
"collision") {
143 B2ERROR(
"Invalid axis name " << arguments[1] <<
". The valid options are thrust and collision");
147 auto func = [order, axisName](
const Particle*) ->
double{
148 StoreObjPtr<EventShapeContainer> evtShapeCont;
151 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
152 return std::numeric_limits<float>::quiet_NaN();
154 if (axisName ==
"thrust")
155 return evtShapeCont->getCleoConeThrust(order);
157 return evtShapeCont->getCleoConeCollision(order);
162 double foxWolframR1(
const Particle*)
164 StoreObjPtr<EventShapeContainer> evtShapeCont;
166 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
167 return std::numeric_limits<float>::quiet_NaN();
169 if (evtShapeCont->getFWMoment(0) == 0) {
170 B2ERROR(
"The 0th-order FoxWolfram moment is zero");
171 return std::numeric_limits<float>::quiet_NaN();
173 return evtShapeCont->getFWMoment(1) / evtShapeCont->getFWMoment(0);
176 double foxWolframR2(
const Particle*)
178 StoreObjPtr<EventShapeContainer> evtShapeCont;
180 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
181 return std::numeric_limits<float>::quiet_NaN();
183 if (evtShapeCont->getFWMoment(0) == 0) {
184 B2ERROR(
"The 0th-order FoxWolfram moment is zero");
185 return std::numeric_limits<float>::quiet_NaN();
187 return evtShapeCont->getFWMoment(2) / evtShapeCont->getFWMoment(0);
190 double foxWolframR3(
const Particle*)
192 StoreObjPtr<EventShapeContainer> evtShapeCont;
194 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
195 return std::numeric_limits<float>::quiet_NaN();
197 if (evtShapeCont->getFWMoment(0) == 0) {
198 B2ERROR(
"The 0th-order FoxWolfram moment is zero");
199 return std::numeric_limits<float>::quiet_NaN();
201 return evtShapeCont->getFWMoment(3) / evtShapeCont->getFWMoment(0);
204 double foxWolframR4(
const Particle*)
206 StoreObjPtr<EventShapeContainer> evtShapeCont;
208 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
209 return std::numeric_limits<float>::quiet_NaN();
211 if (evtShapeCont->getFWMoment(0) == 0) {
212 B2ERROR(
"The 0th-order FoxWolfram moment is zero");
213 return std::numeric_limits<float>::quiet_NaN();
215 return evtShapeCont->getFWMoment(4) / evtShapeCont->getFWMoment(0);
218 double harmonicMomentThrust0(
const Particle*)
220 StoreObjPtr<EventShapeContainer> evtShapeCont;
222 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
223 return std::numeric_limits<float>::quiet_NaN();
225 return evtShapeCont->getHarmonicMomentThrust(0);
228 double harmonicMomentThrust1(
const Particle*)
230 StoreObjPtr<EventShapeContainer> evtShapeCont;
232 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
233 return std::numeric_limits<float>::quiet_NaN();
235 return evtShapeCont->getHarmonicMomentThrust(1);
238 double harmonicMomentThrust2(
const Particle*)
240 StoreObjPtr<EventShapeContainer> evtShapeCont;
242 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
243 return std::numeric_limits<float>::quiet_NaN();
245 return evtShapeCont->getHarmonicMomentThrust(2);
248 double harmonicMomentThrust3(
const Particle*)
250 StoreObjPtr<EventShapeContainer> evtShapeCont;
252 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
253 return std::numeric_limits<float>::quiet_NaN();
255 return evtShapeCont->getHarmonicMomentThrust(3);
258 double harmonicMomentThrust4(
const Particle*)
260 StoreObjPtr<EventShapeContainer> evtShapeCont;
262 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
263 return std::numeric_limits<float>::quiet_NaN();
265 return evtShapeCont->getHarmonicMomentThrust(4);
268 double cleoConeThrust0(
const Particle*)
270 StoreObjPtr<EventShapeContainer> evtShapeCont;
272 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
273 return std::numeric_limits<float>::quiet_NaN();
275 return evtShapeCont->getCleoConeThrust(0);
278 double cleoConeThrust1(
const Particle*)
280 StoreObjPtr<EventShapeContainer> evtShapeCont;
282 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
283 return std::numeric_limits<float>::quiet_NaN();
285 return evtShapeCont->getCleoConeThrust(1);
288 double cleoConeThrust2(
const Particle*)
290 StoreObjPtr<EventShapeContainer> evtShapeCont;
292 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
293 return std::numeric_limits<float>::quiet_NaN();
295 return evtShapeCont->getCleoConeThrust(2);
298 double cleoConeThrust3(
const Particle*)
300 StoreObjPtr<EventShapeContainer> evtShapeCont;
302 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
303 return std::numeric_limits<float>::quiet_NaN();
305 return evtShapeCont->getCleoConeThrust(3);
308 double cleoConeThrust4(
const Particle*)
310 StoreObjPtr<EventShapeContainer> evtShapeCont;
312 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
313 return std::numeric_limits<float>::quiet_NaN();
315 return evtShapeCont->getCleoConeThrust(4);
318 double cleoConeThrust5(
const Particle*)
320 StoreObjPtr<EventShapeContainer> evtShapeCont;
322 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
323 return std::numeric_limits<float>::quiet_NaN();
325 return evtShapeCont->getCleoConeThrust(5);
328 double cleoConeThrust6(
const Particle*)
330 StoreObjPtr<EventShapeContainer> evtShapeCont;
332 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
333 return std::numeric_limits<float>::quiet_NaN();
335 return evtShapeCont->getCleoConeThrust(6);
338 double cleoConeThrust7(
const Particle*)
340 StoreObjPtr<EventShapeContainer> evtShapeCont;
342 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
343 return std::numeric_limits<float>::quiet_NaN();
345 return evtShapeCont->getCleoConeThrust(7);
348 double cleoConeThrust8(
const Particle*)
350 StoreObjPtr<EventShapeContainer> evtShapeCont;
352 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
353 return std::numeric_limits<float>::quiet_NaN();
355 return evtShapeCont->getCleoConeThrust(8);
358 double sphericity(
const Particle*)
360 StoreObjPtr<EventShapeContainer> evtShapeCont;
362 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
363 return std::numeric_limits<float>::quiet_NaN();
366 return 1.5 * (evtShapeCont->getSphericityEigenvalue(1) + evtShapeCont->getSphericityEigenvalue(2)) ;
369 double aplanarity(
const Particle*)
371 StoreObjPtr<EventShapeContainer> evtShapeCont;
373 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
374 return std::numeric_limits<float>::quiet_NaN();
377 return 1.5 * evtShapeCont->getSphericityEigenvalue(2);
380 double forwardHemisphereMass(
const Particle*)
382 StoreObjPtr<EventShapeContainer> evtShapeCont;
384 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
385 return std::numeric_limits<float>::quiet_NaN();
387 return evtShapeCont->getForwardHemisphere4Momentum().M();
390 double forwardHemisphereX(
const Particle*)
392 StoreObjPtr<EventShapeContainer> evtShapeCont;
394 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
395 return std::numeric_limits<float>::quiet_NaN();
397 return evtShapeCont->getForwardHemisphere4Momentum().px();
400 double forwardHemisphereY(
const Particle*)
402 StoreObjPtr<EventShapeContainer> evtShapeCont;
404 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
405 return std::numeric_limits<float>::quiet_NaN();
407 return evtShapeCont->getForwardHemisphere4Momentum().py();
410 double forwardHemisphereZ(
const Particle*)
412 StoreObjPtr<EventShapeContainer> evtShapeCont;
414 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
415 return std::numeric_limits<float>::quiet_NaN();
417 return evtShapeCont->getForwardHemisphere4Momentum().pz();
420 double forwardHemisphereMomentum(
const Particle*)
422 StoreObjPtr<EventShapeContainer> evtShapeCont;
424 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
425 return std::numeric_limits<float>::quiet_NaN();
427 return evtShapeCont->getForwardHemisphere4Momentum().P();
430 double forwardHemisphereEnergy(
const Particle*)
432 StoreObjPtr<EventShapeContainer> evtShapeCont;
434 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
435 return std::numeric_limits<float>::quiet_NaN();
437 return evtShapeCont->getForwardHemisphere4Momentum().E();
440 double backwardHemisphereMass(
const Particle*)
442 StoreObjPtr<EventShapeContainer> evtShapeCont;
444 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
445 return std::numeric_limits<float>::quiet_NaN();
447 return evtShapeCont->getBackwardHemisphere4Momentum().M();
450 double backwardHemisphereX(
const Particle*)
452 StoreObjPtr<EventShapeContainer> evtShapeCont;
454 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
455 return std::numeric_limits<float>::quiet_NaN();
457 return evtShapeCont->getBackwardHemisphere4Momentum().px();
460 double backwardHemisphereY(
const Particle*)
462 StoreObjPtr<EventShapeContainer> evtShapeCont;
464 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
465 return std::numeric_limits<float>::quiet_NaN();
467 return evtShapeCont->getBackwardHemisphere4Momentum().py();
470 double backwardHemisphereZ(
const Particle*)
472 StoreObjPtr<EventShapeContainer> evtShapeCont;
474 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
475 return std::numeric_limits<float>::quiet_NaN();
477 return evtShapeCont->getBackwardHemisphere4Momentum().pz();
480 double backwardHemisphereMomentum(
const Particle*)
482 StoreObjPtr<EventShapeContainer> evtShapeCont;
484 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
485 return std::numeric_limits<float>::quiet_NaN();
487 return evtShapeCont->getBackwardHemisphere4Momentum().P();
490 double backwardHemisphereEnergy(
const Particle*)
492 StoreObjPtr<EventShapeContainer> evtShapeCont;
494 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
495 return std::numeric_limits<float>::quiet_NaN();
497 return evtShapeCont->getBackwardHemisphere4Momentum().E();
500 double thrust(
const Particle*)
502 StoreObjPtr<EventShapeContainer> evtShapeCont;
504 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
505 return std::numeric_limits<float>::quiet_NaN();
507 return evtShapeCont->getThrust();
510 double thrustAxisX(
const Particle*)
512 StoreObjPtr<EventShapeContainer> evtShapeCont;
514 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
515 return std::numeric_limits<float>::quiet_NaN();
517 return evtShapeCont->getThrustAxis().X();
520 double thrustAxisY(
const Particle*)
522 StoreObjPtr<EventShapeContainer> evtShapeCont;
524 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
525 return std::numeric_limits<float>::quiet_NaN();
527 return evtShapeCont->getThrustAxis().Y();
530 double thrustAxisZ(
const Particle*)
532 StoreObjPtr<EventShapeContainer> evtShapeCont;
534 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
535 return std::numeric_limits<float>::quiet_NaN();
537 return evtShapeCont->getThrustAxis().Z();
540 double thrustAxisCosTheta(
const Particle*)
542 StoreObjPtr<EventShapeContainer> evtShapeCont;
544 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
545 return std::numeric_limits<float>::quiet_NaN();
547 return cos(evtShapeCont->getThrustAxis().Theta());
552 if (arguments.size() == 1) {
553 auto variableName = arguments[0];
557 auto func = [var](
const Particle * particle) ->
double {
558 StoreObjPtr<EventShapeContainer> evtShapeCont;
561 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
562 return std::numeric_limits<float>::quiet_NaN();
565 ROOT::Math::XYZVector newZ = evtShapeCont->getThrustAxis();
566 ROOT::Math::XYZVector newY(0, 0, 0);
567 if (newZ.Z() == 0 and newZ.Y() == 0)
572 newY.SetZ(-newZ.Y());
574 ROOT::Math::XYZVector newX = newY.Cross(newZ);
576 UseReferenceFrame<CMSRotationFrame> signalframe(newX, newY, newZ);
578 return std::get<double>(var->function(particle));
582 B2FATAL(
"Wrong number of arguments for meta function useThrustFrame. It only takes one argument, the variable name.");
586 VARIABLE_GROUP(
"EventShape");
588 REGISTER_VARIABLE(
"foxWolframR(i)", foxWolframR, R
"DOC(
589 [Eventbased] Ratio of the i-th to the 0-th order Fox Wolfram moments. The order ``i`` can go from 0 up to 8th.
591 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
592 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
594 REGISTER_VARIABLE("foxWolframH(i)", foxWolframH, R
"DOC(
595 [Eventbased] Returns i-th order Fox Wolfram moment. The order ``i`` can go from 0 up to 8th."
597 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
598 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
599 )DOC",":math:`\\text{GeV}^2`/:math:`\\text{c}^2`");
600 REGISTER_METAVARIABLE("harmonicMoment(i, axisName)", harmonicMoment, R"DOC(
601 [Eventbased] Returns i-th order harmonic moment, calculated with respect to the axis ``axisName``.
602 The order ``i`` can go from 0 up to 8th, the ``axisName`` can be either 'thrust' or 'collision'.
603 The unit of the harmonic moment is ``GeV/c``.
605 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
606 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
607 )DOC", Manager::VariableDataType::c_double);
608 REGISTER_METAVARIABLE("cleoCone(i, axisName)", cleoCone, R
"DOC(
609 [Eventbased] Returns i-th order Cleo cone, calculated with respect to the axis ``axisName``. The order ``i`` can go from 0 up to 8th, the ``axisName`` can be either 'thrust' or 'collision'.
611 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
612 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
613 )DOC", Manager::VariableDataType::c_double);
614 REGISTER_METAVARIABLE("useThrustFrame(variable)", useThrustFrame, R
"DOC(
615 Evaluates a variable value in the thrust reference frame.
617 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
618 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
619 )DOC", Manager::VariableDataType::c_double);
622 REGISTER_VARIABLE("foxWolframR1", foxWolframR1, R
"DOC(
623 [Eventbased] Ratio of the 1-st to the 0-th order Fox Wolfram moments. This is just an alias of foxWolframR(1) defined for the user's convenience.
625 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
626 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
628 REGISTER_VARIABLE("foxWolframR2", foxWolframR2, R
"DOC(
629 [Eventbased] Ratio of the 2-nd to the 0-th order Fox Wolfram moments. This is just an alias of foxWolframR(2) defined for the user's convenience.
631 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
632 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
634 REGISTER_VARIABLE("foxWolframR3", foxWolframR3, R
"DOC(
635 [Eventbased] Ratio of the 3-rd to the 0-th order Fox Wolfram moments. This is just an alias of foxWolframR(3) defined for the user's convenience.
637 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
638 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
640 REGISTER_VARIABLE("foxWolframR4", foxWolframR4, R
"DOC(
641 [Eventbased] Ratio of the 4-th to the 0-th order Fox Wolfram moments. This is just an alias of foxWolframR(4) defined for the user's convenience.
643 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
644 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
647 REGISTER_VARIABLE("harmonicMomentThrust0", harmonicMomentThrust0, R
"DOC(
648 [Eventbased] Harmonic moment of the 0th order calculated with respect to the thrust axis.
650 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
651 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
653 REGISTER_VARIABLE("harmonicMomentThrust1", harmonicMomentThrust1, R"DOC(
654 [Eventbased] Harmonic moment of the 1st order calculated with respect to the thrust axis.
656 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
657 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
659 REGISTER_VARIABLE("harmonicMomentThrust2", harmonicMomentThrust2, R"DOC(
660 [Eventbased] Harmonic moment of the 2nd order calculated with respect to the thrust axis.
662 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
663 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
665 REGISTER_VARIABLE("harmonicMomentThrust3", harmonicMomentThrust3, R"DOC(
666 [Eventbased] Harmonic moment of the 3rd order calculated with respect to the thrust axis.
668 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
669 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
671 REGISTER_VARIABLE("harmonicMomentThrust4", harmonicMomentThrust4, R"DOC(
672 [Eventbased] Harmonic moment of the 4th order calculated with respect to the thrust axis.
674 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
675 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
678 REGISTER_VARIABLE("cleoConeThrust0", cleoConeThrust0, R"DOC(
679 [Eventbased] 0th Cleo cone calculated with respect to the thrust axis.
681 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
682 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
684 REGISTER_VARIABLE("cleoConeThrust1", cleoConeThrust1, R
"DOC(
685 [Eventbased] 1st Cleo cone calculated with respect to the thrust axis.
687 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
688 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
690 REGISTER_VARIABLE("cleoConeThrust2", cleoConeThrust2, R
"DOC(
691 [Eventbased] 2nd Cleo cone calculated with respect to the thrust axis.
693 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
694 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
696 REGISTER_VARIABLE("cleoConeThrust3", cleoConeThrust3, R
"DOC(
697 [Eventbased] 3rd Cleo cone calculated with respect to the thrust axis.
699 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
700 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
702 REGISTER_VARIABLE("cleoConeThrust4", cleoConeThrust4, R
"DOC(
703 [Eventbased] 4th Cleo cone calculated with respect to the thrust axis.
705 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
706 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
708 REGISTER_VARIABLE("cleoConeThrust5", cleoConeThrust5, R
"DOC(
709 [Eventbased] 5th Cleo cone calculated with respect to the thrust axis.
711 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
712 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
714 REGISTER_VARIABLE("cleoConeThrust6", cleoConeThrust6, R
"DOC(
715 [Eventbased] 6th Cleo cone calculated with respect to the thrust axis.
717 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
718 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
720 REGISTER_VARIABLE("cleoConeThrust7", cleoConeThrust7, R
"DOC(
721 [Eventbased] 7th Cleo cone calculated with respect to the thrust axis.
723 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
724 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
726 REGISTER_VARIABLE("cleoConeThrust8", cleoConeThrust8, R
"DOC(
727 [Eventbased] 8th Cleo cone calculated with respect to the thrust axis.
729 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
730 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
734 REGISTER_VARIABLE("sphericity", sphericity, R
"DOC(
735 [Eventbased] Event sphericity, defined as the linear combination of the sphericity eigenvalues :math:`\\lambda_i`: :math:`S = (3/2)(\\lambda_2+\\lambda_3)`.
737 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
738 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
740 REGISTER_VARIABLE("aplanarity", aplanarity, R
"DOC(
741 [Eventbased] Event aplanarity, defined as the 3/2 of the third sphericity eigenvalue.
743 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
744 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
747 REGISTER_VARIABLE("thrust", thrust, R
"DOC(
748 [Eventbased] Event thrust.
750 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
751 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
753 REGISTER_VARIABLE("thrustAxisX", thrustAxisX, R
"DOC(
754 [Eventbased] X component of the thrust axis.
756 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
757 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
759 REGISTER_VARIABLE("thrustAxisY", thrustAxisY, R
"DOC(
760 [Eventbased] Y component of the thrust axis.
762 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
763 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
765 REGISTER_VARIABLE("thrustAxisZ", thrustAxisZ, R
"DOC(
766 [Eventbased] Z component of the thrust axis.
768 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
769 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
771 REGISTER_VARIABLE("thrustAxisCosTheta", thrustAxisCosTheta, R
"DOC(
772 [Eventbased] Cosine of the polar angle component of the thrust axis.
774 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
775 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
778 REGISTER_VARIABLE("forwardHemisphereMass", forwardHemisphereMass, R
"DOC(
779 [Eventbased] Invariant mass of the particles flying in the same direction as the thrust axis.
781 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
782 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
783 )DOC","GeV/:math:`\\text{c}^2`");
784 REGISTER_VARIABLE("forwardHemisphereX", forwardHemisphereX, R"DOC(
785 [Eventbased] X component of the total momentum of the particles flying in the same direction as the thrust axis.
787 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
788 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
790 REGISTER_VARIABLE("forwardHemisphereY", forwardHemisphereY, R"DOC(
791 [Eventbased] Y component of the total momentum of the particles flying in the same direction as the thrust axis.
793 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
794 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
796 REGISTER_VARIABLE("forwardHemisphereZ", forwardHemisphereZ, R"DOC(
797 [Eventbased] Z component of the total momentum of the particles flying in the same direction of the thrust axis.
799 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
800 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
802 REGISTER_VARIABLE("forwardHemisphereMomentum", forwardHemisphereMomentum, R"DOC(
803 [Eventbased] Total momentum of the particles flying in the same direction as the thrust axis.
805 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
806 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
808 REGISTER_VARIABLE("forwardHemisphereEnergy", forwardHemisphereEnergy, R"DOC(
809 [Eventbased] Total energy of the particles flying in the same direction as the thrust axis.
811 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
812 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
814 REGISTER_VARIABLE("backwardHemisphereMass", backwardHemisphereMass, R"DOC(
815 [Eventbased] Invariant mass of the particles flying in the direction opposite to the thrust axis.
817 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
818 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
819 )DOC","GeV/:math:`\\text{c}^2`");
820 REGISTER_VARIABLE("backwardHemisphereX", backwardHemisphereX, R"DOC(
821 [Eventbased] X component of the total momentum of the particles flying in the direction opposite to the thrust axis.
823 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
824 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
826 REGISTER_VARIABLE("backwardHemisphereY", backwardHemisphereY, R"DOC(
827 [Eventbased] Y component of the total momentum of the particles flying in the direction opposite to the thrust axis.
829 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
830 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
832 REGISTER_VARIABLE("backwardHemisphereZ", backwardHemisphereZ, R"DOC(
833 [Eventbased] Z component of the total momentum of the particles flying in the direction opposite to the thrust axis.
835 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
836 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
838 REGISTER_VARIABLE("backwardHemisphereMomentum", backwardHemisphereMomentum, R"DOC(
839 [Eventbased] Total momentum of the particles flying in the direction opposite to the thrust axis.
841 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
842 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
844 REGISTER_VARIABLE("backwardHemisphereEnergy", backwardHemisphereEnergy, R"DOC(
845 [Eventbased] Total energy of the particles flying in the direction opposite to the thrust axis.
847 .. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
848 .. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
std::function< VarVariant(const Particle *)> FunctionPtr
functions stored take a const Particle* and return VarVariant.
const Var * getVariable(std::string name)
Get the variable belonging to the given key.
static Manager & Instance()
get singleton instance.
Abstract base class for different kinds of events.