Belle II Software development
MCInitialParticles.h
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#pragma once
10
11#include <framework/gearbox/Const.h>
12#include <framework/utilities/LabToCms.h>
13
14#include <Math/Boost.h>
15#include <Math/LorentzRotation.h>
16#include <Math/Vector3D.h>
17#include <Math/Vector4D.h>
18
19namespace Belle2 {
30 class MCInitialParticles: public TObject {
31
32 public:
33
37 c_generateCMS = 1 << 0,
45 c_smearVertex = 1 << 3,
48 };
49
51 MCInitialParticles(): TObject() {}
52
56
59 {
60 m_her = b.m_her; m_ler = b.m_ler; m_vertex = b.m_vertex;
61 m_validFlag = b.m_validFlag; m_generationFlags = b.m_generationFlags;
62 resetBoost();
63 return *this;
64 }
65
67 bool operator==(const MCInitialParticles& b) const
68 {
69 // FIXME: sin(x) returns slightly different values on
70 // different platforms in some cases so we cannot just do an equality
71 // comparison. We need to do this more elegantly, this is just for
72 // testing if it solves all problems
73#if defined(MCP_DBL_CMP) || defined(MCP_VEC3_CMP) || defined(MCP_VEC4_CMP)
74#error Macro already defined, cannot continue
75#endif
76#define MCP_DBL_CMP(a,b,x) ((a.X()==b.X())||(std::abs(a.X()-b.X())<1e-10))
77#define MCP_VEC3_CMP(a,b) (MCP_DBL_CMP(a,b,X) && MCP_DBL_CMP(a,b,Y) && MCP_DBL_CMP(a,b,Z))
78#define MCP_VEC4_CMP(a,b) (MCP_VEC3_CMP(a,b) && MCP_DBL_CMP(a,b,E))
79 return MCP_VEC4_CMP(m_her, b.m_her) && MCP_VEC4_CMP(m_ler, b.m_ler) && MCP_VEC3_CMP(m_vertex, b.m_vertex)
80 && m_validFlag == b.m_validFlag &&
81 m_generationFlags == b.m_generationFlags;
82#undef MCP_DBL_CMP
83#undef MCP_VEC3_CMP
84#undef MCP_VEC4_CMP
85 }
86
89
96 void set(const ROOT::Math::PxPyPzEVector& her, const ROOT::Math::PxPyPzEVector& ler, const ROOT::Math::XYZVector& vertex)
97 {
98 m_her = her;
99 m_ler = ler;
100 m_vertex = vertex;
101 m_validFlag = true;
102 resetBoost();
103 }
104
115 void setByLorentzTransformation(double Ecms, double bX, double bY, double bZ, double angleXZ, double angleYZ,
116 const ROOT::Math::XYZVector& vertex)
117 {
118 if (m_labToCMS) delete m_labToCMS;
119 if (m_CMSToLab) delete m_CMSToLab;
120
121 m_invariantMass = Ecms;
122 m_CMSToLab = new ROOT::Math::LorentzRotation();
123 m_labToCMS = new ROOT::Math::LorentzRotation();
124 *m_CMSToLab = cmsToLab(bX, bY, bZ, angleXZ, angleYZ);
125 *m_labToCMS = m_CMSToLab->Inverse();
126
127 const double me = Const::electron.getMass();
128 double p = sqrt(Ecms * Ecms / 4 - me * me);
129 m_her = (*m_CMSToLab) * ROOT::Math::PxPyPzEVector(0.0, 0.0, p, Ecms / 2);
130 m_ler = (*m_CMSToLab) * ROOT::Math::PxPyPzEVector(0.0, 0.0, -p, Ecms / 2);
131
132 m_vertex = vertex;
133 m_validFlag = true;
134 }
135
137 void setHER(const ROOT::Math::PxPyPzEVector& her)
138 {
139 m_her = her;
140 resetBoost();
141 }
142
144 void setLER(const ROOT::Math::PxPyPzEVector& ler)
145 {
146 m_ler = ler;
147 resetBoost();
148 }
149
151 void setVertex(const ROOT::Math::XYZVector& vertex)
152 {
153 m_vertex = vertex;
154 }
155
157 void setTime(double time) {m_time = time;}
158
160 virtual void setGenerationFlags(int flags) { m_generationFlags = flags; }
161
163 const ROOT::Math::PxPyPzEVector& getHER() const { return m_her; }
164
166 const ROOT::Math::PxPyPzEVector& getLER() const { return m_ler; }
167
169 const ROOT::Math::XYZVector& getVertex() const { return m_vertex; }
170
172 double getTime() const {return m_time;}
173
175 double getEnergy() const { return (m_her + m_ler).E(); }
176
178 double getMass() const { calculateBoost(); return m_invariantMass; }
179
181 ROOT::Math::XYZVector getBoostVector() const { calculateBoost(); return -(m_her + m_ler).BoostToCM(); }
182
184 ROOT::Math::PxPyPzEVector getBoostedHER() const { calculateBoost(); return m_boostedHER; }
185
187 const ROOT::Math::LorentzRotation& getLabToCMS() const
188 {
189 calculateBoost(); return *m_labToCMS;
190 }
191
193 const ROOT::Math::LorentzRotation& getCMSToLab() const
194 {
195 calculateBoost(); return *m_CMSToLab;
196 }
197
199 bool getValidFlag() const { return m_validFlag; }
200
202 int getGenerationFlags() const { return m_generationFlags; }
203
205 bool hasGenerationFlags(int flags) const { return (m_generationFlags & flags) == flags; }
206
209 std::string getGenerationFlagString(const std::string& separator = " ") const;
210
211
219 static ROOT::Math::LorentzRotation cmsToLab(double bX, double bY, double bZ, double angleXZ, double angleYZ);
220
221 private:
222
224 void calculateBoost() const;
226 void resetBoost();
228 ROOT::Math::PxPyPzEVector m_her;
230 ROOT::Math::PxPyPzEVector m_ler;
232 ROOT::Math::XYZVector m_vertex;
234 double m_time = 0;
236 mutable ROOT::Math::LorentzRotation* m_labToCMS{nullptr};
238 mutable ROOT::Math::LorentzRotation* m_CMSToLab{nullptr};
240 mutable double m_invariantMass{0.0};
242 mutable ROOT::Math::PxPyPzEVector m_boostedHER;
244 bool m_validFlag = false;
245
246 protected:
247
250
251 private:
252
255 };
256
258 {
259 if (m_labToCMS)
260 return;
261
262 ROOT::Math::PxPyPzEVector beam = m_her + m_ler;
263 // Save the invariant mass because it's used very often in analysis
264 m_invariantMass = beam.M();
265
266 // If we generate events in CMS we already are in CMS and there is no
267 // transformation so let's use the identity
269 m_labToCMS = new ROOT::Math::LorentzRotation();
270 m_CMSToLab = new ROOT::Math::LorentzRotation();
271 return;
272 }
273
274 m_boostedHER = ROOT::Math::LorentzRotation(ROOT::Math::Boost(beam.BoostToCM())) * m_her;
275 double cmsAngleXZ = atan(m_boostedHER.X() / m_boostedHER.Z());
276 double cmsAngleYZ = atan(m_boostedHER.Y() / m_boostedHER.Z());
277 m_labToCMS = new ROOT::Math::LorentzRotation(
278 LabToCms::rotateLabToCms(-beam.BoostToCM(), cmsAngleXZ, cmsAngleYZ)
279 );
280
281 //cache derived quantities
282 m_CMSToLab = new ROOT::Math::LorentzRotation(m_labToCMS->Inverse());
283 }
284
285
286
287
288
290 {
291 delete m_labToCMS;
292 delete m_CMSToLab;
293 m_labToCMS = nullptr;
294 m_CMSToLab = nullptr;
295 m_invariantMass = 0.0;
296 m_boostedHER = ROOT::Math::PxPyPzEVector();
297 }
298
300} //Belle2 namespace
double getMass() const
Particle mass.
Definition: UnitConst.cc:356
static const ChargedStable electron
electron particle
Definition: Const.h:659
static ROOT::Math::LorentzRotation rotateLabToCms(const ROOT::Math::XYZVector &boostVector, double cmsAngleXZ, double cmsAngleYZ)
Function takes 3D boostVector and angles of the HER momentum in the CM system obtained by pure boost.
Definition: LabToCms.h:43
This class contains the initial state for the given event.
int m_generationFlags
Flags to be used when generating events.
const ROOT::Math::PxPyPzEVector & getHER() const
Get 4vector of the high energy beam.
ROOT::Math::LorentzRotation * m_CMSToLab
transient
bool operator==(const MCInitialParticles &b) const
Equality operator.
int getGenerationFlags() const
Get the generation flags to be used for event generation (ORed combination of EGenerationFlags)
bool getValidFlag() const
Get the flag to check if a valid MCInitialParticles object was already generated and filled in an eve...
ROOT::Math::PxPyPzEVector m_boostedHER
transient
ROOT::Math::PxPyPzEVector getBoostedHER() const
Get the 4-vector of electron beam in CM system obtained by pure boost.
void setVertex(const ROOT::Math::XYZVector &vertex)
Set the vertex position.
double m_time
collision time
const ROOT::Math::PxPyPzEVector & getLER() const
Get 4vector of the low energy beam.
ROOT::Math::PxPyPzEVector m_ler
LER 4vector.
void setTime(double time)
Set collision time.
double getEnergy() const
Get the actual collision energy (in lab system)
ROOT::Math::PxPyPzEVector m_her
HER 4vector.
EGenerationFlags
Possible Flags for initial event generation.
@ c_smearBeam
smear the full beam momentum (energy and direction)
@ c_generateCMS
generate initial event in CMS instead of lab
@ c_smearBeamEnergy
smear energy of HER and LER (but not direction)
@ c_smearBeamDirection
smear direction of HER and LER (but not energy)
ROOT::Math::LorentzRotation * m_labToCMS
Boost from Lab into CMS.
const ROOT::Math::LorentzRotation & getLabToCMS() const
Return the LorentzRotation to convert from lab to CMS frame.
ROOT::Math::XYZVector m_vertex
collision position
const ROOT::Math::LorentzRotation & getCMSToLab() const
Return the LorentzRotation to convert from CMS to lab frame.
std::string getGenerationFlagString(const std::string &separator=" ") const
Return string representation of all active flags for printing.
MCInitialParticles & operator=(const MCInitialParticles &b)
Assignment operator.
virtual void setGenerationFlags(int flags)
Set the generation flags to be used for event generation (ORed combination of EGenerationFlags)
static ROOT::Math::LorentzRotation cmsToLab(double bX, double bY, double bZ, double angleXZ, double angleYZ)
Return the LorentzRotation from CMS to LAB based on the following parameters.
MCInitialParticles()
Default constructor.
void setByLorentzTransformation(double Ecms, double bX, double bY, double bZ, double angleXZ, double angleYZ, const ROOT::Math::XYZVector &vertex)
Initialize the event values from CMS energy and parameters of the Lorentz transformation between LAB ...
bool hasGenerationFlags(int flags) const
Check if a certain set of EGenerationFlags is set.
const ROOT::Math::XYZVector & getVertex() const
Get the position of the collision.
void setHER(const ROOT::Math::PxPyPzEVector &her)
Set the High Energy Beam 4-momentum.
MCInitialParticles(const MCInitialParticles &b)
Copy constructor.
void setLER(const ROOT::Math::PxPyPzEVector &ler)
Set the Low Energy Beam 4-momentum.
virtual ~MCInitialParticles()
Free memory of the LorentzRotation if it was created.
ROOT::Math::XYZVector getBoostVector() const
Get the boost vector (velocity of system produced in the collision)
double getMass() const
Get the invariant mass of the collision (= energy in CMS)
double getTime() const
Get collision time.
void set(const ROOT::Math::PxPyPzEVector &her, const ROOT::Math::PxPyPzEVector &ler, const ROOT::Math::XYZVector &vertex)
Set the initial event values, i.e.
ClassDef(MCInitialParticles, 5)
ROOT Dictionary.
void resetBoost()
Reset cached transformations after changing parameters.
void calculateBoost() const
Calculate the Lorentz transformations LAB->CMS & CMS->LAB if necessary.
static const double me
electron mass
Definition: beamHelpers.h:177
double atan(double a)
atan for double
Definition: beamHelpers.h:34
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.