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 B2FATAL(
"foxWolframR cannot be called without providing the moment order");
37 int order = std::lround(index[0]);
39 if (order < 0 || order > 8)
40 B2FATAL(
"The Fox-Wolfram moment order must be within 0 and 8.");
42 StoreObjPtr<EventShapeContainer> evtShapeCont;
44 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
47 if (evtShapeCont->getFWMoment(0) == 0) {
48 B2INFO(
"The 0th-order FoxWolfram moment is zero");
51 return evtShapeCont->getFWMoment(order) / evtShapeCont->getFWMoment(0);
54 double foxWolframH(
const Particle*,
const std::vector<double>& index)
56 if (index.size() != 1)
57 B2FATAL(
"foxWolframH cannot be called without providing the moment order");
59 int order = std::lround(index[0]);
61 if (order < 0 || order > 8)
62 B2FATAL(
"The Fox-Wolfram moment order must be within 0 and 8.");
64 StoreObjPtr<EventShapeContainer> evtShapeCont;
66 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
69 return evtShapeCont->getFWMoment(order);
74 if (arguments.size() != 2)
75 B2FATAL(
"harmonicMoment requires two arguments: the harmonic order (0-8) and the reference axis name (thrust or collision)");
79 order = Belle2::convertString<int>(arguments[0]);
80 }
catch (std::invalid_argument&) {
81 B2FATAL(
"First argument of harmonicMoment must be an integer");
83 std::string axisName = arguments[1];
84 boost::to_lower(axisName);
86 if (order < 0 || order > 8)
87 B2FATAL(
"The Fox-Wolfram moment order must be within 0 and 8.");
89 if (axisName !=
"thrust" && axisName !=
"collision")
90 B2FATAL(
"Invalid axis name " << arguments[1] <<
". The valid options are thrust and collision");
92 auto func = [order, axisName](
const Particle*) ->
double{
93 StoreObjPtr<EventShapeContainer> evtShapeCont;
96 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
99 if (axisName ==
"thrust")
100 return evtShapeCont->getHarmonicMomentThrust(order);
102 return evtShapeCont->getHarmonicMomentCollision(order);
109 if (arguments.size() != 2)
110 B2FATAL(
"cleoCone requires two arguments: the cone order (0-9) and the reference axis name (thrust or collision)");
114 order = Belle2::convertString<int>(arguments[0]);
115 }
catch (std::invalid_argument&) {
116 B2FATAL(
"Argument of cleoCone must be an integer");
118 std::string axisName = arguments[1];
119 boost::to_lower(axisName);
121 if (order < 0 || order > 8)
122 B2FATAL(
"The CLEO cone order must be within 0 and 8.");
124 if (axisName !=
"thrust" && axisName !=
"collision")
125 B2FATAL(
"Invalid axis name " << arguments[1] <<
". The valid options are thrust and collision");
127 auto func = [order, axisName](
const Particle*) ->
double{
128 StoreObjPtr<EventShapeContainer> evtShapeCont;
131 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
134 if (axisName ==
"thrust")
135 return evtShapeCont->getCleoConeThrust(order);
137 return evtShapeCont->getCleoConeCollision(order);
142 double foxWolframR1(
const Particle*)
144 StoreObjPtr<EventShapeContainer> evtShapeCont;
146 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
149 if (evtShapeCont->getFWMoment(0) == 0) {
150 B2INFO(
"The 0th-order FoxWolfram moment is zero");
153 return evtShapeCont->getFWMoment(1) / evtShapeCont->getFWMoment(0);
156 double foxWolframR2(
const Particle*)
158 StoreObjPtr<EventShapeContainer> evtShapeCont;
160 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
163 if (evtShapeCont->getFWMoment(0) == 0) {
164 B2INFO(
"The 0th-order FoxWolfram moment is zero");
167 return evtShapeCont->getFWMoment(2) / evtShapeCont->getFWMoment(0);
170 double foxWolframR3(
const Particle*)
172 StoreObjPtr<EventShapeContainer> evtShapeCont;
174 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
177 if (evtShapeCont->getFWMoment(0) == 0) {
178 B2INFO(
"The 0th-order FoxWolfram moment is zero");
181 return evtShapeCont->getFWMoment(3) / evtShapeCont->getFWMoment(0);
184 double foxWolframR4(
const Particle*)
186 StoreObjPtr<EventShapeContainer> evtShapeCont;
188 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
191 if (evtShapeCont->getFWMoment(0) == 0) {
192 B2INFO(
"The 0th-order FoxWolfram moment is zero");
195 return evtShapeCont->getFWMoment(4) / evtShapeCont->getFWMoment(0);
198 double harmonicMomentThrust0(
const Particle*)
200 StoreObjPtr<EventShapeContainer> evtShapeCont;
202 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
205 return evtShapeCont->getHarmonicMomentThrust(0);
208 double harmonicMomentThrust1(
const Particle*)
210 StoreObjPtr<EventShapeContainer> evtShapeCont;
212 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
215 return evtShapeCont->getHarmonicMomentThrust(1);
218 double harmonicMomentThrust2(
const Particle*)
220 StoreObjPtr<EventShapeContainer> evtShapeCont;
222 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
225 return evtShapeCont->getHarmonicMomentThrust(2);
228 double harmonicMomentThrust3(
const Particle*)
230 StoreObjPtr<EventShapeContainer> evtShapeCont;
232 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
235 return evtShapeCont->getHarmonicMomentThrust(3);
238 double harmonicMomentThrust4(
const Particle*)
240 StoreObjPtr<EventShapeContainer> evtShapeCont;
242 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
245 return evtShapeCont->getHarmonicMomentThrust(4);
248 double cleoConeThrust0(
const Particle*)
250 StoreObjPtr<EventShapeContainer> evtShapeCont;
252 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
255 return evtShapeCont->getCleoConeThrust(0);
258 double cleoConeThrust1(
const Particle*)
260 StoreObjPtr<EventShapeContainer> evtShapeCont;
262 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
265 return evtShapeCont->getCleoConeThrust(1);
268 double cleoConeThrust2(
const Particle*)
270 StoreObjPtr<EventShapeContainer> evtShapeCont;
272 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
275 return evtShapeCont->getCleoConeThrust(2);
278 double cleoConeThrust3(
const Particle*)
280 StoreObjPtr<EventShapeContainer> evtShapeCont;
282 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
285 return evtShapeCont->getCleoConeThrust(3);
288 double cleoConeThrust4(
const Particle*)
290 StoreObjPtr<EventShapeContainer> evtShapeCont;
292 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
295 return evtShapeCont->getCleoConeThrust(4);
298 double cleoConeThrust5(
const Particle*)
300 StoreObjPtr<EventShapeContainer> evtShapeCont;
302 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
305 return evtShapeCont->getCleoConeThrust(5);
308 double cleoConeThrust6(
const Particle*)
310 StoreObjPtr<EventShapeContainer> evtShapeCont;
312 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
315 return evtShapeCont->getCleoConeThrust(6);
318 double cleoConeThrust7(
const Particle*)
320 StoreObjPtr<EventShapeContainer> evtShapeCont;
322 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
325 return evtShapeCont->getCleoConeThrust(7);
328 double cleoConeThrust8(
const Particle*)
330 StoreObjPtr<EventShapeContainer> evtShapeCont;
332 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
335 return evtShapeCont->getCleoConeThrust(8);
338 double sphericity(
const Particle*)
340 StoreObjPtr<EventShapeContainer> evtShapeCont;
342 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
346 return 1.5 * (evtShapeCont->getSphericityEigenvalue(1) + evtShapeCont->getSphericityEigenvalue(2)) ;
349 double aplanarity(
const Particle*)
351 StoreObjPtr<EventShapeContainer> evtShapeCont;
353 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
357 return 1.5 * evtShapeCont->getSphericityEigenvalue(2);
360 double forwardHemisphereMass(
const Particle*)
362 StoreObjPtr<EventShapeContainer> evtShapeCont;
364 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
367 return evtShapeCont->getForwardHemisphere4Momentum().M();
370 double forwardHemisphereX(
const Particle*)
372 StoreObjPtr<EventShapeContainer> evtShapeCont;
374 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
377 return evtShapeCont->getForwardHemisphere4Momentum().px();
380 double forwardHemisphereY(
const Particle*)
382 StoreObjPtr<EventShapeContainer> evtShapeCont;
384 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
387 return evtShapeCont->getForwardHemisphere4Momentum().py();
390 double forwardHemisphereZ(
const Particle*)
392 StoreObjPtr<EventShapeContainer> evtShapeCont;
394 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
397 return evtShapeCont->getForwardHemisphere4Momentum().pz();
400 double forwardHemisphereMomentum(
const Particle*)
402 StoreObjPtr<EventShapeContainer> evtShapeCont;
404 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
407 return evtShapeCont->getForwardHemisphere4Momentum().P();
410 double forwardHemisphereEnergy(
const Particle*)
412 StoreObjPtr<EventShapeContainer> evtShapeCont;
414 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
417 return evtShapeCont->getForwardHemisphere4Momentum().E();
420 double backwardHemisphereMass(
const Particle*)
422 StoreObjPtr<EventShapeContainer> evtShapeCont;
424 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
427 return evtShapeCont->getBackwardHemisphere4Momentum().M();
430 double backwardHemisphereX(
const Particle*)
432 StoreObjPtr<EventShapeContainer> evtShapeCont;
434 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
437 return evtShapeCont->getBackwardHemisphere4Momentum().px();
440 double backwardHemisphereY(
const Particle*)
442 StoreObjPtr<EventShapeContainer> evtShapeCont;
444 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
447 return evtShapeCont->getBackwardHemisphere4Momentum().py();
450 double backwardHemisphereZ(
const Particle*)
452 StoreObjPtr<EventShapeContainer> evtShapeCont;
454 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
457 return evtShapeCont->getBackwardHemisphere4Momentum().pz();
460 double backwardHemisphereMomentum(
const Particle*)
462 StoreObjPtr<EventShapeContainer> evtShapeCont;
464 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
467 return evtShapeCont->getBackwardHemisphere4Momentum().P();
470 double backwardHemisphereEnergy(
const Particle*)
472 StoreObjPtr<EventShapeContainer> evtShapeCont;
474 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
477 return evtShapeCont->getBackwardHemisphere4Momentum().E();
480 double thrust(
const Particle*)
482 StoreObjPtr<EventShapeContainer> evtShapeCont;
484 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
487 return evtShapeCont->getThrust();
490 double thrustAxisX(
const Particle*)
492 StoreObjPtr<EventShapeContainer> evtShapeCont;
494 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
497 return evtShapeCont->getThrustAxis().X();
500 double thrustAxisY(
const Particle*)
502 StoreObjPtr<EventShapeContainer> evtShapeCont;
504 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
507 return evtShapeCont->getThrustAxis().Y();
510 double thrustAxisZ(
const Particle*)
512 StoreObjPtr<EventShapeContainer> evtShapeCont;
514 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
517 return evtShapeCont->getThrustAxis().Z();
520 double thrustAxisCosTheta(
const Particle*)
522 StoreObjPtr<EventShapeContainer> evtShapeCont;
524 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
527 return cos(evtShapeCont->getThrustAxis().Theta());
532 if (arguments.size() == 1) {
533 auto variableName = arguments[0];
537 auto func = [var](
const Particle * particle) ->
double {
538 StoreObjPtr<EventShapeContainer> evtShapeCont;
541 B2ERROR(
"No EventShapeContainer object has been found in the datastore");
545 ROOT::Math::XYZVector newZ = evtShapeCont->getThrustAxis();
546 ROOT::Math::XYZVector newY(0, 0, 0);
547 if (newZ.Z() == 0 and newZ.Y() == 0)
552 newY.SetZ(-newZ.Y());
554 ROOT::Math::XYZVector newX = newY.Cross(newZ);
556 UseReferenceFrame<CMSRotationFrame> signalframe(newX, newY, newZ);
558 return std::get<double>(var->function(particle));
562 B2FATAL(
"Wrong number of arguments for meta function useThrustFrame. It only takes one argument, the variable name.");
566 VARIABLE_GROUP(
"EventShape");
568 REGISTER_VARIABLE(
"foxWolframR(i)", foxWolframR, R
"DOC(
569[Eventbased] Ratio of the i-th to the 0-th order Fox Wolfram moments. The order ``i`` can go from 0 up to 8th.
571.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
572.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
574 REGISTER_VARIABLE("foxWolframH(i)", foxWolframH, R
"DOC(
575[Eventbased] Returns i-th order Fox Wolfram moment. The order ``i`` can go from 0 up to 8th."
577.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
578.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
580)DOC",":math:`\\text{GeV}^2`/:math:`\\text{c}^2`");
581 REGISTER_METAVARIABLE("harmonicMoment(i, axisName)", harmonicMoment, R"DOC(
582[Eventbased] Returns i-th order harmonic moment, calculated with respect to the axis ``axisName``.
583The order ``i`` can go from 0 up to 8th, the ``axisName`` can be either 'thrust' or 'collision'.
584The unit of the harmonic moment is ``GeV/c``.
586.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
587.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
588)DOC", Manager::VariableDataType::c_double);
589 REGISTER_METAVARIABLE("cleoCone(i, axisName)", cleoCone, R
"DOC(
590[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'.
592.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
593.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
594)DOC", Manager::VariableDataType::c_double);
595 REGISTER_METAVARIABLE("useThrustFrame(variable)", useThrustFrame, R
"DOC(
596Evaluates a variable value in the thrust reference frame.
598.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
599.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
600)DOC", Manager::VariableDataType::c_double);
603 REGISTER_VARIABLE("foxWolframR1", foxWolframR1, R
"DOC(
604[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.
606.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
607.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
609 REGISTER_VARIABLE("foxWolframR2", foxWolframR2, R
"DOC(
610[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.
612.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
613.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
615 REGISTER_VARIABLE("foxWolframR3", foxWolframR3, R
"DOC(
616[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.
618.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
619.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
621 REGISTER_VARIABLE("foxWolframR4", foxWolframR4, R
"DOC(
622[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.
624.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
625.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
628 REGISTER_VARIABLE("harmonicMomentThrust0", harmonicMomentThrust0, R
"DOC(
629[Eventbased] Harmonic moment of the 0th order calculated with respect to the thrust axis.
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`.
635 REGISTER_VARIABLE("harmonicMomentThrust1", harmonicMomentThrust1, R"DOC(
636[Eventbased] Harmonic moment of the 1st order calculated with respect to the thrust axis.
638.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
639.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
642 REGISTER_VARIABLE("harmonicMomentThrust2", harmonicMomentThrust2, R"DOC(
643[Eventbased] Harmonic moment of the 2nd order calculated with respect to the thrust axis.
645.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
646.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
649 REGISTER_VARIABLE("harmonicMomentThrust3", harmonicMomentThrust3, R"DOC(
650[Eventbased] Harmonic moment of the 3rd order calculated with respect to the thrust axis.
652.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
653.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
656 REGISTER_VARIABLE("harmonicMomentThrust4", harmonicMomentThrust4, R"DOC(
657[Eventbased] Harmonic moment of the 4th order calculated with respect to the thrust axis.
659.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
660.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
664 REGISTER_VARIABLE("cleoConeThrust0", cleoConeThrust0, R"DOC(
665[Eventbased] 0th Cleo cone calculated with respect to the thrust axis.
667.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
668.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
670 REGISTER_VARIABLE("cleoConeThrust1", cleoConeThrust1, R
"DOC(
671[Eventbased] 1st Cleo cone calculated with respect to the thrust axis.
673.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
674.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
676 REGISTER_VARIABLE("cleoConeThrust2", cleoConeThrust2, R
"DOC(
677[Eventbased] 2nd Cleo cone calculated with respect to the thrust axis.
679.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
680.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
682 REGISTER_VARIABLE("cleoConeThrust3", cleoConeThrust3, R
"DOC(
683[Eventbased] 3rd Cleo cone calculated with respect to the thrust axis.
685.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
686.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
688 REGISTER_VARIABLE("cleoConeThrust4", cleoConeThrust4, R
"DOC(
689[Eventbased] 4th Cleo cone calculated with respect to the thrust axis.
691.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
692.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
694 REGISTER_VARIABLE("cleoConeThrust5", cleoConeThrust5, R
"DOC(
695[Eventbased] 5th Cleo cone calculated with respect to the thrust axis.
697.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
698.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
700 REGISTER_VARIABLE("cleoConeThrust6", cleoConeThrust6, R
"DOC(
701[Eventbased] 6th Cleo cone calculated with respect to the thrust axis.
703.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
704.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
706 REGISTER_VARIABLE("cleoConeThrust7", cleoConeThrust7, R
"DOC(
707[Eventbased] 7th Cleo cone calculated with respect to the thrust axis.
709.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
710.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
712 REGISTER_VARIABLE("cleoConeThrust8", cleoConeThrust8, R
"DOC(
713[Eventbased] 8th Cleo cone calculated with respect to the thrust axis.
715.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
716.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
720 REGISTER_VARIABLE("sphericity", sphericity, R
"DOC(
721[Eventbased] Event sphericity, defined as the linear combination of the sphericity eigenvalues :math:`\lambda_i`: :math:`S = (3/2)(\lambda_2+\lambda_3)`.
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("aplanarity", aplanarity, R
"DOC(
727[Eventbased] Event aplanarity, defined as the 3/2 of the third sphericity eigenvalue.
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`.
733 REGISTER_VARIABLE("thrust", thrust, R
"DOC(
734[Eventbased] Event thrust.
736.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
737.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
739 REGISTER_VARIABLE("thrustAxisX", thrustAxisX, R
"DOC(
740[Eventbased] X component of the thrust axis.
742.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
743.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
745 REGISTER_VARIABLE("thrustAxisY", thrustAxisY, R
"DOC(
746[Eventbased] Y component of the thrust axis.
748.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
749.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
751 REGISTER_VARIABLE("thrustAxisZ", thrustAxisZ, R
"DOC(
752[Eventbased] Z component of the thrust axis.
754.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
755.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
757 REGISTER_VARIABLE("thrustAxisCosTheta", thrustAxisCosTheta, R
"DOC(
758[Eventbased] Cosine of the polar angle component of the thrust axis.
760.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
761.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
764 REGISTER_VARIABLE("forwardHemisphereMass", forwardHemisphereMass, R
"DOC(
765[Eventbased] Invariant mass of the particles flying in the same direction as the thrust axis.
767.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
768.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
770)DOC","GeV/:math:`\\text{c}^2`");
771 REGISTER_VARIABLE("forwardHemisphereX", forwardHemisphereX, R"DOC(
772[Eventbased] X component of the total momentum of the particles flying in the same direction as 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("forwardHemisphereY", forwardHemisphereY, R"DOC(
779[Eventbased] Y component of the total momentum 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`.
785 REGISTER_VARIABLE("forwardHemisphereZ", forwardHemisphereZ, R"DOC(
786[Eventbased] Z component of the total momentum of the particles flying in the same direction of the thrust axis.
788.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
789.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
792 REGISTER_VARIABLE("forwardHemisphereMomentum", forwardHemisphereMomentum, R"DOC(
793[Eventbased] Total momentum of the particles flying in the same direction as the thrust axis.
795.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
796.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
799 REGISTER_VARIABLE("forwardHemisphereEnergy", forwardHemisphereEnergy, R"DOC(
800[Eventbased] Total energy of the particles flying in the same direction as the thrust axis.
802.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
803.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
806 REGISTER_VARIABLE("backwardHemisphereMass", backwardHemisphereMass, R"DOC(
807[Eventbased] Invariant mass of the particles flying in the direction opposite to the thrust axis.
809.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
810.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
812)DOC","GeV/:math:`\\text{c}^2`");
813 REGISTER_VARIABLE("backwardHemisphereX", backwardHemisphereX, R"DOC(
814[Eventbased] X component of the total momentum of the particles flying in the direction opposite to the thrust axis.
816.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
817.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
820 REGISTER_VARIABLE("backwardHemisphereY", backwardHemisphereY, R"DOC(
821[Eventbased] Y 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`.
827 REGISTER_VARIABLE("backwardHemisphereZ", backwardHemisphereZ, R"DOC(
828[Eventbased] Z component of the total momentum of the particles flying in the direction opposite to the thrust axis.
830.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
831.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
834 REGISTER_VARIABLE("backwardHemisphereMomentum", backwardHemisphereMomentum, R"DOC(
835[Eventbased] Total momentum of the particles flying in the direction opposite to the thrust axis.
837.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
838.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
841 REGISTER_VARIABLE("backwardHemisphereEnergy", backwardHemisphereEnergy, R"DOC(
842[Eventbased] Total energy of the particles flying in the direction opposite to the thrust axis.
844.. warning:: You have to run the Event Shape builder module for this variable to be meaningful.
845.. seealso:: :ref:`analysis_eventshape` and `modularAnalysis.buildEventShape`.
static const double doubleNaN
quiet_NaN
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.