Belle II Software  light-2403-persian
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 <TVector3.h>
15 #include <Math/Boost.h>
16 #include <Math/LorentzRotation.h>
17 #include <Math/Vector3D.h>
18 #include <Math/Vector4D.h>
19 
20 namespace Belle2 {
31  class MCInitialParticles: public TObject {
32 
33  public:
34 
38  c_generateCMS = 1 << 0,
46  c_smearVertex = 1 << 3,
49  };
50 
52  MCInitialParticles(): TObject() {}
53 
57 
60  {
61  m_her = b.m_her; m_ler = b.m_ler; m_vertex = b.m_vertex;
62  m_validFlag = b.m_validFlag; m_generationFlags = b.m_generationFlags;
63  resetBoost();
64  return *this;
65  }
66 
68  bool operator==(const MCInitialParticles& b) const
69  {
70  // FIXME: sin(x) returns slightly different values on
71  // different platforms in some cases so we cannot just do an equality
72  // comparison. We need to do this more elegantly, this is just for
73  // testing if it solves all problems
74 #if defined(MCP_DBL_CMP) || defined(MCP_VEC3_CMP) || defined(MCP_VEC4_CMP)
75 #error Macro already defined, cannot continue
76 #endif
77 #define MCP_DBL_CMP(a,b,x) ((a.X()==b.X())||(std::abs(a.X()-b.X())<1e-10))
78 #define MCP_VEC3_CMP(a,b) (MCP_DBL_CMP(a,b,X) && MCP_DBL_CMP(a,b,Y) && MCP_DBL_CMP(a,b,Z))
79 #define MCP_VEC4_CMP(a,b) (MCP_VEC3_CMP(a,b) && MCP_DBL_CMP(a,b,E))
80  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)
81  && m_validFlag == b.m_validFlag &&
82  m_generationFlags == b.m_generationFlags;
83 #undef MCP_DBL_CMP
84 #undef MCP_VEC3_CMP
85 #undef MCP_VEC4_CMP
86  }
87 
90 
97  void set(const ROOT::Math::PxPyPzEVector& her, const ROOT::Math::PxPyPzEVector& ler, const ROOT::Math::XYZVector& vertex)
98  {
99  m_her = her;
100  m_ler = ler;
101  m_vertex = vertex;
102  m_validFlag = true;
103  resetBoost();
104  }
105 
116  void setByLorentzTransformation(double Ecms, double bX, double bY, double bZ, double angleXZ, double angleYZ,
117  const ROOT::Math::XYZVector& vertex)
118  {
119  if (m_labToCMS) delete m_labToCMS;
120  if (m_CMSToLab) delete m_CMSToLab;
121 
122  m_invariantMass = Ecms;
123  m_CMSToLab = new ROOT::Math::LorentzRotation();
124  m_labToCMS = new ROOT::Math::LorentzRotation();
125  *m_CMSToLab = cmsToLab(bX, bY, bZ, angleXZ, angleYZ);
126  *m_labToCMS = m_CMSToLab->Inverse();
127 
128  const double me = Const::electron.getMass();
129  double p = sqrt(Ecms * Ecms / 4 - me * me);
130  m_her = (*m_CMSToLab) * ROOT::Math::PxPyPzEVector(0.0, 0.0, p, Ecms / 2);
131  m_ler = (*m_CMSToLab) * ROOT::Math::PxPyPzEVector(0.0, 0.0, -p, Ecms / 2);
132 
133  m_vertex = vertex;
134  m_validFlag = true;
135  }
136 
138  void setHER(const ROOT::Math::PxPyPzEVector& her)
139  {
140  m_her = her;
141  resetBoost();
142  }
143 
145  void setLER(const ROOT::Math::PxPyPzEVector& ler)
146  {
147  m_ler = ler;
148  resetBoost();
149  }
150 
152  void setVertex(const ROOT::Math::XYZVector& vertex)
153  {
154  m_vertex = vertex;
155  }
156 
158  void setTime(double time) {m_time = time;}
159 
161  virtual void setGenerationFlags(int flags) { m_generationFlags = flags; }
162 
164  const ROOT::Math::PxPyPzEVector& getHER() const { return m_her; }
165 
167  const ROOT::Math::PxPyPzEVector& getLER() const { return m_ler; }
168 
170  const ROOT::Math::XYZVector& getVertex() const { return m_vertex; }
171 
173  double getTime() const {return m_time;}
174 
176  double getEnergy() const { return (m_her + m_ler).E(); }
177 
179  double getMass() const { calculateBoost(); return m_invariantMass; }
180 
182  ROOT::Math::XYZVector getBoostVector() const { calculateBoost(); return -(m_her + m_ler).BoostToCM(); }
183 
185  ROOT::Math::PxPyPzEVector getBoostedHER() const { calculateBoost(); return m_boostedHER; }
186 
188  const ROOT::Math::LorentzRotation& getLabToCMS() const
189  {
190  calculateBoost(); return *m_labToCMS;
191  }
192 
194  const ROOT::Math::LorentzRotation& getCMSToLab() const
195  {
196  calculateBoost(); return *m_CMSToLab;
197  }
198 
200  bool getValidFlag() const { return m_validFlag; }
201 
203  int getGenerationFlags() const { return m_generationFlags; }
204 
206  bool hasGenerationFlags(int flags) const { return (m_generationFlags & flags) == flags; }
207 
210  std::string getGenerationFlagString(const std::string& separator = " ") const;
211 
212 
220  static ROOT::Math::LorentzRotation cmsToLab(double bX, double bY, double bZ, double angleXZ, double angleYZ);
221 
222  private:
223 
225  void calculateBoost() const;
227  void resetBoost();
229  ROOT::Math::PxPyPzEVector m_her;
231  ROOT::Math::PxPyPzEVector m_ler;
233  ROOT::Math::XYZVector m_vertex;
235  double m_time = 0;
237  mutable ROOT::Math::LorentzRotation* m_labToCMS{nullptr};
239  mutable ROOT::Math::LorentzRotation* m_CMSToLab{nullptr};
241  mutable double m_invariantMass{0.0};
243  mutable ROOT::Math::PxPyPzEVector m_boostedHER;
245  bool m_validFlag = false;
246 
247  protected:
248 
251 
252  private:
253 
256  };
257 
259  {
260  if (m_labToCMS)
261  return;
262 
263  ROOT::Math::PxPyPzEVector beam = m_her + m_ler;
264  // Save the invariant mass because it's used very often in analysis
265  m_invariantMass = beam.M();
266 
267  // If we generate events in CMS we already are in CMS and there is no
268  // transformation so let's use the identity
270  m_labToCMS = new ROOT::Math::LorentzRotation();
271  m_CMSToLab = new ROOT::Math::LorentzRotation();
272  return;
273  }
274 
275  m_boostedHER = ROOT::Math::LorentzRotation(ROOT::Math::Boost(beam.BoostToCM())) * m_her;
276  double cmsAngleXZ = atan(m_boostedHER.X() / m_boostedHER.Z());
277  double cmsAngleYZ = atan(m_boostedHER.Y() / m_boostedHER.Z());
278  m_labToCMS = new ROOT::Math::LorentzRotation(
279  LabToCms::rotateLabToCms(-beam.BoostToCM(), cmsAngleXZ, cmsAngleYZ)
280  );
281 
282  //cache derived quantities
283  m_CMSToLab = new ROOT::Math::LorentzRotation(m_labToCMS->Inverse());
284  }
285 
286 
287 
288 
289 
291  {
292  delete m_labToCMS;
293  delete m_CMSToLab;
294  m_labToCMS = nullptr;
295  m_CMSToLab = nullptr;
296  m_invariantMass = 0.0;
297  m_boostedHER = ROOT::Math::PxPyPzEVector();
298  }
299 
301 } //Belle2 namespace
double getMass() const
Particle mass.
Definition: UnitConst.cc:356
static const ChargedStable electron
electron particle
Definition: Const.h:650
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.
const ROOT::Math::PxPyPzEVector & getLER() const
Get 4vector of the low energy beam.
int m_generationFlags
Flags to be used when generating events.
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.
const ROOT::Math::PxPyPzEVector & getHER() const
Get 4vector of the high energy beam.
double m_time
collision time
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::XYZVector & getVertex() const
Get the position of the collision.
MCInitialParticles & operator=(const MCInitialParticles &b)
Assignment operator.
ROOT::Math::XYZVector m_vertex
collision position
std::string getGenerationFlagString(const std::string &separator=" ") const
Return string representation of all active flags for printing.
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.
const ROOT::Math::LorentzRotation & getLabToCMS() const
Return the LorentzRotation to convert from lab to CMS frame.
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.
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.
const ROOT::Math::LorentzRotation & getCMSToLab() const
Return the LorentzRotation to convert from CMS to lab frame.
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.
Abstract base class for different kinds of events.
Definition: ClusterUtils.h:24