Belle II Software light-2406-ragdoll
HelicityVariables.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 header.
10#include <analysis/variables/HelicityVariables.h>
11
12#include <analysis/VariableManager/Utility.h>
13
14#include <analysis/variables/EventVariables.h>
15
16#include <analysis/dataobjects/Particle.h>
17
18#include <analysis/utility/ReferenceFrame.h>
19#include <analysis/VariableManager/Manager.h>
20
21#include <framework/utilities/Conversion.h>
22#include <framework/gearbox/Const.h>
23
24#include <boost/algorithm/string.hpp>
25
26#include <Math/Boost.h>
27#include <Math/Vector4D.h>
28using namespace ROOT::Math;
29#include <cmath>
30
31namespace Belle2 {
36 namespace Variable {
37
38 double cosHelicityAngleMomentum(const Particle* part)
39 {
40
41 const auto& frame = ReferenceFrame::GetCurrent();
42 B2Vector3D motherBoost = frame.getMomentum(part).BoostToCM();
43 B2Vector3D motherMomentum = frame.getMomentum(part).Vect();
44 const auto& daughters = part -> getDaughters() ;
45
46 if (daughters.size() == 2) {
47
48 // Only for pi0 -> gamma gamma, gamma -> e+ e-
49 bool isOneConversion = false;
50 if (part->getPDGCode() == Const::pi0.getPDGCode()) {
51 for (auto& idaughter : daughters) {
52 // both daughter must be gamma
53 if (idaughter -> getPDGCode() != Const::photon.getPDGCode()) {
54 isOneConversion = false;
55 break;
56 }
57 // check if one of gammas has two daughters
58 if (idaughter -> getNDaughters() == 2) {
59 if (std::abs(idaughter -> getDaughters()[0]-> getPDGCode()) == Const::electron.getPDGCode()
60 && std::abs(idaughter -> getDaughters()[1]-> getPDGCode()) == Const::electron.getPDGCode()) { // e+ e-
61 isOneConversion = true;
62 }
63 }
64 }
65 }
66
67 if (isOneConversion) {
68 B2WARNING("cosHelicityAngleMomentum: Special treatment for pi0->gamma gamma, gamma -> e+ e-, is called. "
69 "This treatment is going to be deprecated and we recommend using ``cosHelicityAngleMomentumPi0Dalitz`` "
70 "If you find this message in another case, it must be a bug. Please report it to the software mailing list.");
71
72 //only for pi0 decay where one gamma converts
73 PxPyPzEVector pGamma;
74 for (auto& idaughter : daughters) {
75 if (idaughter -> getNDaughters() == 2) continue;
76 else pGamma = frame.getMomentum(idaughter);
77 }
78
79 pGamma = Boost(motherBoost) * pGamma;
80
81 return std::cos(motherMomentum.Angle(pGamma.Vect()));
82
83 } else {
84 PxPyPzEVector pDaughter1 = frame.getMomentum(daughters[0]);
85 PxPyPzEVector pDaughter2 = frame.getMomentum(daughters[1]);
86
87 pDaughter1 = Boost(motherBoost) * pDaughter1;
88 pDaughter2 = Boost(motherBoost) * pDaughter2;
89
90 B2Vector3D p12 = (pDaughter2 - pDaughter1).Vect();
91
92 return std::cos(motherMomentum.Angle(p12));
93 }
94
95 } else if (daughters.size() == 3) {
96
97 PxPyPzEVector pDaughter1 = frame.getMomentum(daughters[0]);
98 PxPyPzEVector pDaughter2 = frame.getMomentum(daughters[1]);
99 PxPyPzEVector pDaughter3 = frame.getMomentum(daughters[2]);
100
101 pDaughter1 = Boost(motherBoost) * pDaughter1;
102 pDaughter2 = Boost(motherBoost) * pDaughter2;
103 pDaughter3 = Boost(motherBoost) * pDaughter3;
104
105 B2Vector3D p12 = (pDaughter2 - pDaughter1).Vect();
106 B2Vector3D p13 = (pDaughter3 - pDaughter1).Vect();
107
108 B2Vector3D n = p12.Cross(p13);
109
110 return std::cos(motherMomentum.Angle(n));
111
112 } else return Const::doubleNaN;
113
114 }
115
116 double cosHelicityAngleMomentumPi0Dalitz(const Particle* part)
117 {
118
119 const auto& frame = ReferenceFrame::GetCurrent();
120 B2Vector3D motherBoost = frame.getMomentum(part).BoostToCM();
121 B2Vector3D motherMomentum = frame.getMomentum(part).Vect();
122 const auto& daughters = part -> getDaughters() ;
123
124
125 if (daughters.size() == 3) {
126
127 PxPyPzEVector pGamma;
128
129 for (auto& idaughter : daughters) {
130 if (std::abs(idaughter -> getPDGCode()) == Const::photon.getPDGCode()) {
131 pGamma = frame.getMomentum(idaughter);
132 break;
133 }
134 }
135 pGamma = Boost(motherBoost) * pGamma;
136
137 return std::cos(motherMomentum.Angle(pGamma.Vect()));
138
139 } else if (daughters.size() == 2) { // only for pi0 -> gamma gamma, gamma -> e+ e-
140
141 PxPyPzEVector pGamma;
142
143 // both daughters must be gamma
144 if (daughters[0] -> getPDGCode() != Const::photon.getPDGCode() or
145 daughters[1] -> getPDGCode() != Const::photon.getPDGCode())
146 return Const::doubleNaN;
147
148 if (daughters[0] -> getNDaughters() == 2 and daughters[1] -> getNDaughters() == 0) {
149 if (std::abs(daughters[0] -> getDaughters()[0]-> getPDGCode()) == Const::electron.getPDGCode()
150 && std::abs(daughters[0] -> getDaughters()[1]-> getPDGCode()) == Const::electron.getPDGCode()) { // e+ e-
151 pGamma = frame.getMomentum(daughters[1]);
152 } else {
153 return Const::doubleNaN;
154 }
155 } else if (daughters[0] -> getNDaughters() == 0 and daughters[1] -> getNDaughters() == 2) {
156 if (std::abs(daughters[1] -> getDaughters()[0]-> getPDGCode()) == Const::electron.getPDGCode()
157 && std::abs(daughters[1] -> getDaughters()[1]-> getPDGCode()) == Const::electron.getPDGCode()) { // e+ e-
158 pGamma = frame.getMomentum(daughters[0]);
159 } else {
160 return Const::doubleNaN;
161 }
162 } else {
163 return Const::doubleNaN;
164 }
165
166 pGamma = Boost(motherBoost) * pGamma;
167
168 return std::cos(motherMomentum.Angle(pGamma.Vect()));
169
170 } else return Const::doubleNaN;
171
172 }
173
174
175 double cosHelicityAngleBeamMomentum(const Particle* mother, const std::vector<double>& index)
176 {
177 if (index.size() != 1) {
178 B2FATAL("Wrong number of arguments for cosHelicityAngleIfCMSIsTheMother");
179 }
180
181 int idau = std::lround(index[0]);
182
183 const Particle* part = mother->getDaughter(idau);
184 if (!part) {
185 B2FATAL("Couldn't find the " << idau << "th daughter");
186 }
187
188 PxPyPzEVector beam4Vector(getBeamPx(nullptr), getBeamPy(nullptr), getBeamPz(nullptr), getBeamE(nullptr));
189 PxPyPzEVector part4Vector = part->get4Vector();
190 PxPyPzEVector mother4Vector = mother->get4Vector();
191
192 B2Vector3D motherBoost = mother4Vector.BoostToCM();
193
194 beam4Vector = Boost(motherBoost) * beam4Vector;
195 part4Vector = Boost(motherBoost) * part4Vector;
196
197 return - part4Vector.Vect().Dot(beam4Vector.Vect()) / part4Vector.P() / beam4Vector.P();
198 }
199
200
201 double cosHelicityAngle(const Particle* mother, const std::vector<double>& indices)
202 {
203 if (indices.size() != 2) {
204 B2FATAL("Wrong number of arguments for cosHelicityAngleIfRefFrameIsTheDaughter: two are needed.");
205 }
206
207 int iDau = std::lround(indices[0]);
208 int iGrandDau = std::lround(indices[1]);
209
210 const Particle* daughter = mother->getDaughter(iDau);
211 if (!daughter)
212 B2FATAL("Couldn't find the " << iDau << "th daughter.");
213
214 const Particle* grandDaughter = daughter->getDaughter(iGrandDau);
215 if (!grandDaughter)
216 B2FATAL("Couldn't find the " << iGrandDau << "th daughter of the " << iDau << "th daughter.");
217
218 PxPyPzEVector mother4Vector = mother->get4Vector();
219 PxPyPzEVector daughter4Vector = daughter->get4Vector();
220 PxPyPzEVector grandDaughter4Vector = grandDaughter->get4Vector();
221
222 B2Vector3D daughterBoost = daughter4Vector.BoostToCM();
223
224 // We boost the momentum of the mother and of the granddaughter to the reference frame of the daughter.
225 grandDaughter4Vector = Boost(daughterBoost) * grandDaughter4Vector;
226 mother4Vector = Boost(daughterBoost) * mother4Vector;
227
228 return - grandDaughter4Vector.Vect().Dot(mother4Vector.Vect()) / grandDaughter4Vector.P() / mother4Vector.P();
229 }
230
231 double cosAcoplanarityAngle(const Particle* mother, const std::vector<double>& granddaughters)
232 {
233 if (granddaughters.size() != 2) {
234 B2FATAL("Wrong number of arguments for cosAcoplanarityAngleIfRefFrameIsTheMother: two are needed.");
235 }
236
237 if (mother->getNDaughters() != 2)
238 B2FATAL("cosAcoplanarityAngleIfRefFrameIsTheMother: this variable works only for two-body decays.");
239
240 int iGrandDau1 = std::lround(granddaughters[0]);
241 int iGrandDau2 = std::lround(granddaughters[1]);
242
243 const Particle* daughter1 = mother->getDaughter(0);
244 const Particle* daughter2 = mother->getDaughter(1);
245
246 const Particle* grandDaughter1 = daughter1->getDaughter(iGrandDau1);
247 if (!grandDaughter1)
248 B2FATAL("Couldn't find the " << iGrandDau1 << "th daughter of the first daughter.");
249
250 const Particle* grandDaughter2 = daughter2->getDaughter(iGrandDau2);
251 if (!grandDaughter2)
252 B2FATAL("Couldn't find the " << iGrandDau2 << "th daughter of the second daughter.");
253
254 PxPyPzEVector mother4Vector = mother->get4Vector();
255 PxPyPzEVector daughter4Vector1 = daughter1->get4Vector();
256 PxPyPzEVector daughter4Vector2 = daughter2->get4Vector();
257 PxPyPzEVector grandDaughter4Vector1 = grandDaughter1->get4Vector();
258 PxPyPzEVector grandDaughter4Vector2 = grandDaughter2->get4Vector();
259
260 B2Vector3D motherBoost = mother4Vector.BoostToCM();
261 B2Vector3D daughter1Boost = daughter4Vector1.BoostToCM();
262 B2Vector3D daughter2Boost = daughter4Vector2.BoostToCM();
263
264 // Boosting daughters to reference frame of the mother
265 daughter4Vector1 = Boost(motherBoost) * daughter4Vector1;
266 daughter4Vector2 = Boost(motherBoost) * daughter4Vector2;
267
268 // Boosting each granddaughter to reference frame of its mother
269 grandDaughter4Vector1 = Boost(daughter1Boost) * grandDaughter4Vector1;
270 grandDaughter4Vector2 = Boost(daughter2Boost) * grandDaughter4Vector2;
271
272 // We calculate the normal vectors of the decay two planes
273 B2Vector3D normalVector1 = daughter4Vector1.Vect().Cross(grandDaughter4Vector1.Vect());
274 B2Vector3D normalVector2 = daughter4Vector2.Vect().Cross(grandDaughter4Vector2.Vect());
275
276 return std::cos(normalVector1.Angle(normalVector2));
277 }
278
279 double cosHelicityAnglePrimary(const Particle* part)
280 {
281 return part->getCosHelicity();
282 }
283
284 double cosHelicityAngleDaughter(const Particle* part, const std::vector<double>& indices)
285 {
286 if ((indices.size() == 0) || (indices.size() > 2)) {
287 B2FATAL("Wrong number of arguments for cosHelicityAngleDaughter: one or two are needed.");
288 }
289
290 int iDaughter = std::lround(indices[0]);
291 int iGrandDaughter = 0;
292 if (indices.size() == 2) {
293 iGrandDaughter = std::lround(indices[1]);
294 }
295
296 return part->getCosHelicityDaughter(iDaughter, iGrandDaughter);
297 }
298
299 double acoplanarityAngle(const Particle* part)
300 {
301 return part->getAcoplanarity();
302 }
303
304
305 double cosHelicityAngleForQuasiTwoBodyDecay(const Particle* mother, const std::vector<double>& indices)
306 {
307 if (indices.size() != 2) {
308 B2FATAL("Wrong number of arguments for cosHelicityAngleForQuasiTwoBodyDecay: two are needed.");
309 }
310
311 if (mother->getNDaughters() != 3)
312 return Const::doubleNaN;
313
314 int iDau = std::lround(indices[0]);
315 int jDau = std::lround(indices[1]);
316
317 const Particle* iDaughter = mother->getDaughter(iDau);
318 if (!iDaughter)
319 return Const::doubleNaN;
320
321 const Particle* jDaughter = mother->getDaughter(jDau);
322 if (!jDaughter)
323 return Const::doubleNaN;
324
325 PxPyPzEVector mother4Vector = mother->get4Vector();
326 PxPyPzEVector iDaughter4Vector = iDaughter->get4Vector();
327 PxPyPzEVector jDaughter4Vector = jDaughter->get4Vector();
328
329 PxPyPzEVector resonance4Vector = iDaughter4Vector + jDaughter4Vector;
330 B2Vector3D resonanceBoost = resonance4Vector.BoostToCM();
331
332 iDaughter4Vector = Boost(resonanceBoost) * iDaughter4Vector;
333 mother4Vector = Boost(resonanceBoost) * mother4Vector;
334
335 return - iDaughter4Vector.Vect().Dot(mother4Vector.Vect()) / iDaughter4Vector.P() / mother4Vector.P();
336 }
337
338 Manager::FunctionPtr momentaTripleProduct(const std::vector<std::string>& arguments)
339 {
340 if (arguments.size() != 3) {
341 B2FATAL("Wrong number of arguments for momentaTripleProduct: three (particles) are needed.");
342 }
343
344 // wrap with func and return it
345 auto func = [arguments](const Particle * mother) -> double {
346 auto iDau = arguments[0];
347 auto jDau = arguments[1];
348 auto kDau = arguments[2];
349
350 const Particle* iDaughter = mother->getParticleFromGeneralizedIndexString(iDau);
351 if (!iDaughter) B2FATAL("Couldn't find the " << iDau << "th daughter.");
352 const Particle* jDaughter = mother->getParticleFromGeneralizedIndexString(jDau);
353 if (!jDaughter) B2FATAL("Couldn't find the " << jDau << "th daughter.");
354 const Particle* kDaughter = mother->getParticleFromGeneralizedIndexString(kDau);
355 if (!kDaughter) B2FATAL("Couldn't find the " << kDau << "th daughter.");
356
357 PxPyPzEVector mother4Vector = mother->get4Vector();
358 PxPyPzEVector iDaughter4Vector = iDaughter->get4Vector();
359 PxPyPzEVector jDaughter4Vector = jDaughter->get4Vector();
360 PxPyPzEVector kDaughter4Vector = kDaughter->get4Vector();
361
362 B2Vector3D motherBoost = mother4Vector.BoostToCM();
363
364 // We boost the momenta of offspring to the reference frame of the mother.
365 iDaughter4Vector = Boost(motherBoost) * iDaughter4Vector;
366 jDaughter4Vector = Boost(motherBoost) * jDaughter4Vector;
367 kDaughter4Vector = Boost(motherBoost) * kDaughter4Vector;
368
369 // cross product: p_j x p_k
370 B2Vector3D jkDaughterCrossProduct = jDaughter4Vector.Vect().Cross(kDaughter4Vector.Vect());
371 // triple product: p_i * (p_j x p_k)
372 return iDaughter4Vector.Vect().Dot(jkDaughterCrossProduct) ;
373 };
374 return func;
375 }
376
377
378 VARIABLE_GROUP("Helicity variables");
379
380 REGISTER_VARIABLE("cosHelicityAngleMomentum", cosHelicityAngleMomentum, R"DOC(
381 If the given particle has two daughters: cosine of the angle between the line defined by the momentum difference of the two daughters
382 in the frame of the given particle (mother)
383 and the momentum of the given particle in the lab frame.
384
385 If the given particle has three daughters: cosine of the angle between the normal vector of the plane defined by
386 the momenta of the three daughters in the frame of the given particle (mother)
387 and the momentum of the given particle in the lab frame.
388
389 Otherwise, it returns 0.)DOC");
390
391 REGISTER_VARIABLE("cosHelicityAngleMomentumPi0Dalitz", cosHelicityAngleMomentumPi0Dalitz, R"DOC(
392 To be used for the decay :math:`\pi^0 \to e^+ e^- \gamma`:
393 cosine of the angle between the momentum of the gamma in the frame of the given particle (mother)
394 and the momentum of the given particle in the lab frame.
395
396 One can call the variable for the decay :math:`\pi^0 \to \gamma \gamma, \gamma \to e^+ e^-` as well.
397
398 Otherwise, it returns 0.)DOC");
399
400 REGISTER_VARIABLE("cosHelicityAngleBeamMomentum(i)", cosHelicityAngleBeamMomentum, R"DOC(
401 Cosine of the helicity angle of the :math:`i`-th daughter of the particle provided,
402 assuming that the mother of the provided particle corresponds to the centre-of-mass system, whose parameters are
403 automatically loaded by the function, given the accelerator's conditions.)DOC");
404
405 REGISTER_VARIABLE("cosHelicityAngle(i, j)", cosHelicityAngle, R"DOC(
406 Cosine of the helicity angle between the momentum of the selected granddaughter and the direction opposite to the momentum of the provided particle
407 in the reference frame of the selected daughter (:math:`\theta_1` and :math:`\theta_2` in the
408 `PDG <https://journals.aps.org/prd/abstract/10.1103/PhysRevD.98.030001>`_ 2018, p. 722).
409
410 This variable needs two integer arguments: the first one, ``i``, is the index of the daughter and the second one, ``j`` is the index of the granddaughter.
411
412 For example, in the decay :math:`B^0 \to \left(J/\psi \to \mu^+ \mu^-\right) \left(K^{*0} \to K^+ \pi^-\right)`,
413 if the provided particle is :math:`B^0` and the selected indices are (0, 0),
414 the variable will return the angle between the momentum of the :math:`\mu^+` and the direction opposite to the momentum of
415 the :math:`B^0`, both momenta in the rest frame of the :math:`J/\psi`.
416
417 This variable is needed for angular analyses of :math:`B`-meson decays into two vector particles.)DOC");
418
419 REGISTER_VARIABLE("cosAcoplanarityAngle(i, j)", cosAcoplanarityAngle, R"DOC(
420 Cosine of the acoplanarity angle (:math:`\Phi` in the `PDG Polarization Review <http://pdg.lbl.gov/2019/reviews/rpp2018-rev-b-decays-polarization.pdf>`_).
421 Given a two-body decay, the acoplanarity angle is defined as
422 the angle between the two decay planes in the reference frame of the mother.
423
424 We calculate the acoplanarity angle as the angle between the two
425 normal vectors of the decay planes. Each normal vector is the cross product of the momentum of one daughter (in the frame of the mother) and the
426 momentum of one of the granddaughters (in the reference frame of the daughter).
427
428 This variable needs two integer arguments: the first one, ``i`` is the index of the first granddaughter, and the second one, ``j`` the index of the
429 second granddaughter.
430
431 For example, in the decay :math:`B^0 \to \left(J/\psi \to \mu^+ \mu^-\right) \left(K^{*0} \to K^+ \pi^-\right)`, if the provided particle is :math:`B^0` and the selected indices are (0, 0),
432 the variable will return the acoplanarity using the :math:`\mu^+` and the :math:`K^+` granddaughters.)DOC");
433
434 REGISTER_VARIABLE("cosHelicityAnglePrimary", cosHelicityAnglePrimary, R"DOC(
435 Cosine of the helicity angle (see``Particle::getCosHelicity``) assuming the center of mass system as mother rest frame.
436 See `PDG Polarization Review <http://pdg.lbl.gov/2019/reviews/rpp2018-rev-b-decays-polarization.pdf>`_ for the definition of the helicity angle.)DOC");
437
438 REGISTER_VARIABLE("cosHelicityAngleDaughter(i [, j] )", cosHelicityAngleDaughter, R"DOC(
439 Cosine of the helicity angle of the i-th daughter (see ``Particle::getCosHelicityDaughter``).
440 The optional second argument is the index of the granddaughter that defines the angle, default is 0.
441
442 For example, in the decay: :math:`B^0 \to \left(J/\psi \to \mu^+ \mu^-\right) \left(K^{*0} \to K^+ \pi^-\right)`, if the provided particle is :math:`B^0` and the selected index is 0,
443 the variable will return the helicity angle of the :math:`\mu^+`.
444 If the selected index is 1 the variable will return the helicity angle of the :math:`K^+` (defined via the rest frame of the :math:`K^{*0}`).
445 In rare cases if one wanted the helicity angle of the second granddaughter, indices 1,1 would return the helicity angle of the :math:`\pi^-`).
446
447 See `PDG Polarization Review <http://pdg.lbl.gov/2019/reviews/rpp2018-rev-b-decays-polarization.pdf>`_ for the definition of the helicity angle.)DOC");
448
449 REGISTER_VARIABLE("acoplanarityAngle", acoplanarityAngle, R"DOC(
450Acoplanarity angle (see ``Particle::getAcoplanarity``) assuming a two body decay of the particle and its daughters.
451See `PDG Polarization Review <http://pdg.lbl.gov/2019/reviews/rpp2018-rev-b-decays-polarization.pdf>`_ for the definition of the acoplanarity angle.
452
453)DOC", "rad");
454
455 REGISTER_VARIABLE("cosHelicityAngleForQuasiTwoBodyDecay(i, j)", cosHelicityAngleForQuasiTwoBodyDecay, R"DOC(
456 Cosine of the helicity angle between the momentum of the provided particle and the momentum of the first selected
457 daughter (i-th) in the reference frame of the sum of two selected daughters (i-th + j-th).
458
459 The variable is supposed to be used for the analysis of a quasi-two-body decay. The number of daughters of the given
460 particle must be three. Otherwise, the variable returns NaN.
461
462 For example, in the decay :math:`\bar{B}^0 \to D^+ K^- K^{*0}`, if the provided particle is :math:`\bar{B}^0` and
463 the selected indices are (1, 2), the variable will return the angle between the momentum of the :math:`\bar{B}^0`
464 and the momentum of the :math:`K^-`, both momenta in the rest frame of the :math:`K^- K^{*0}`.)DOC");
465
466 REGISTER_METAVARIABLE("momentaTripleProduct(i,j,k)", momentaTripleProduct, R"DOC(
467a triple-product of three momenta of offspring in the mother rest frame: :math:`C_T=\vec{p}_i\cdot(\vec{p}_j\times\vec{p}_k)`. For examples,
468In a four-body decay M->D1D2D3D4, momentaTripleProduct(0,1,2) returns CT using the momenta of D1D2D3 particles.
469In other decays involving secondary decay, e.g. for M->(R->D1D2)D3D4, momentaTripleProduct(0:0,1,2) returns C_T using momenta of D1D3D4 particles.
470)DOC", Manager::VariableDataType::c_double);
471
472 }
474}
B2Vector3< DataType > Cross(const B2Vector3< DataType > &p) const
Cross product.
Definition: B2Vector3.h:296
DataType Dot(const B2Vector3< DataType > &p) const
Scalar product.
Definition: B2Vector3.h:290
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType pi0
neutral pion particle
Definition: Const.h:674
static const double doubleNaN
quiet_NaN
Definition: Const.h:703
static const ParticleType photon
photon particle
Definition: Const.h:673
static const ChargedStable electron
electron particle
Definition: Const.h:659
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
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
Abstract base class for different kinds of events.
Definition: ClusterUtils.h:24