CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RawParticle.h
Go to the documentation of this file.
1 #ifndef RAWPARTTICLE_H
2 #define RAWPARTTICLE_H
3 
7 
8 #include "Math/GenVector/RotationX.h"
9 #include "Math/GenVector/RotationY.h"
10 #include "Math/GenVector/RotationZ.h"
11 #include "Math/GenVector/Rotation3D.h"
12 #include "Math/GenVector/AxisAngle.h"
13 #include "Math/GenVector/Boost.h"
14 
16 
17 #include <string>
18 #include <iosfwd>
19 
20 
21 
30 
32 public:
33 
34  typedef ROOT::Math::AxisAngle Rotation;
35  typedef ROOT::Math::Rotation3D Rotation3D;
36  typedef ROOT::Math::RotationX RotationX;
37  typedef ROOT::Math::RotationY RotationY;
38  typedef ROOT::Math::RotationZ RotationZ;
39  typedef ROOT::Math::Boost Boost;
40 
41  RawParticle();
42 
43  virtual ~RawParticle();
44 
49 
53  RawParticle(const int id,
54  const XYZTLorentzVector& p);
55 
60  const XYZTLorentzVector& p);
61 
66  const XYZTLorentzVector& xStart);
67 
71  RawParticle(double px, double py, double pz, double e);
72 
74  RawParticle(const RawParticle &p);
75 
77  RawParticle& operator = (const RawParticle & rhs );
78 
79 public:
80 
85  void setID(const int id);
86 
91  void setID(const std::string name);
92 
97  void setStatus(int istat);
98 
100  void setMass(float m);
101 
103  void setCharge(float q);
104 
106  void setT(const double t);
107 
109  void setVertex(const XYZTLorentzVector& vtx);
110  void setVertex(double xv, double yv, double zv, double tv);
111 
112  /*** methods to be overloaded to include vertex ***/
113 
118  void boost(double bx, double by, double bz);
119  void boost(const Boost& b);
120 
121  // inline void boost(const Hep3Vector<double> &bv );
122 
128  void rotate(double rphi, const XYZVector& raxis);
129  void rotate(const Rotation& r);
130  void rotate(const Rotation3D& r);
131 
133  // void rotateUz(Hep3Vector &nuz);
134 
138  void rotateX(double rphi);
139  void rotate(const RotationX& r);
140 
145  void rotateY(double rphi);
146  void rotate(const RotationY& r);
151  void rotateZ(double rphi);
152  void rotate(const RotationZ& r);
153 
155  void translate(const XYZVector& t);
156 
157  // inline RawParticle & transform(const HepRotation &rot);
158  // inline RawParticle & transform(const HepLorentzRotation &rot);
159 
164  void chargeConjugate();
165 
166  int pid() const;
167 
168  int status() const;
169 
170  double charge() const;
171 
172  double PDGcharge() const;
173 
174  double mass() const;
175 
176  double PDGmass() const;
177 
178  double PDGcTau() const;
179 
181  std::string PDGname() const;
182 
186  double eta() const;
187 
189  double cos2Theta() const;
190  double cos2ThetaV() const;
191 
192  double et() const;
193 
194  double x() const;
195  double X() const;
196 
197  double y() const;
198  double Y() const;
199 
200  double z() const;
201  double Z() const;
202 
203  double t() const;
204  double T() const;
205 
206  double r() const;
207  double R() const;
208 
209  double r2() const;
210  double R2() const;
211 
212  const XYZTLorentzVector& vertex() const;
213 
214  const XYZTLorentzVector& momentum() const;
215 
221  void printName() const;
222 
227  void print() const;
228 
233  int isUsed() const {return myUsed;}
234 
237  void use() { myUsed = 1;}
238 
241  void reUse() { myUsed = 0;}
242 
243  private:
244 
245  void init();
246 
247  protected:
248 
250  int myId;
251  int myStatus;
252  int myUsed;
253  double myCharge;
254  double myMass;
256 
257  private:
259 };
260 
261 
262 std::ostream& operator <<(std::ostream& o , const RawParticle& p);
263 
264 inline int RawParticle::pid() const { return myId; }
265 inline int RawParticle::status() const { return myStatus; }
266 inline double RawParticle::eta() const { return -std::log(std::tan(this->theta()/2.)); }
267 inline double RawParticle::cos2Theta() const { return Pz()*Pz()/Vect().Mag2(); }
268 inline double RawParticle::cos2ThetaV() const { return Z()*Z()/myVertex.Vect().Mag2(); }
269 inline double RawParticle::x() const { return myVertex.Px(); }
270 inline double RawParticle::y() const { return myVertex.Py(); }
271 inline double RawParticle::z() const { return myVertex.Pz(); }
272 inline double RawParticle::t() const { return myVertex.E(); }
273 inline double RawParticle::X() const { return myVertex.Px(); }
274 inline double RawParticle::Y() const { return myVertex.Py(); }
275 inline double RawParticle::Z() const { return myVertex.Pz(); }
276 inline double RawParticle::T() const { return myVertex.E(); }
277 inline double RawParticle::R() const { return std::sqrt(R2()); }
278 inline double RawParticle::R2() const { return myVertex.Perp2(); }
279 inline double RawParticle::r() const { return std::sqrt(r2()); }
280 inline double RawParticle::r2() const { return myVertex.Perp2(); }
281 inline double RawParticle::charge() const { return myCharge; }
282 inline double RawParticle::mass() const { return myMass; }
283 
284 inline const XYZTLorentzVector& RawParticle::vertex() const { return myVertex ; }
285 inline const XYZTLorentzVector& RawParticle::momentum() const { return *this; }
286 
287 inline void RawParticle::setVertex(const XYZTLorentzVector& vtx) { myVertex = vtx; }
288 inline void RawParticle::setVertex(double a, double b, double c, double d) { myVertex.SetXYZT(a,b,c,d); }
289 
291  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
292 }
293 
295  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
296 }
297 
299  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
300 }
301 
303  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
304 }
305 
307  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
308 }
309 
311  XYZTLorentzVector p ( b(momentum()) ); SetXYZT(p.Px(),p.Py(),p.Pz(),p.E());
312 }
313 
314 inline void RawParticle::translate(const XYZVector& tr) {
315  myVertex.SetXYZT(X()+tr.X(),Y()+tr.Y(),Z()+tr.Z(),T());
316 }
317 
318 
319 
320 
321 
322 #endif
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:138
void translate(const XYZVector &t)
Definition: RawParticle.h:314
void boost(double bx, double by, double bz)
Definition: RawParticle.cc:182
void rotateZ(double rphi)
Definition: RawParticle.cc:175
void use()
Definition: RawParticle.h:237
XYZTLorentzVector myVertex
the four vector of the vertex
Definition: RawParticle.h:249
int status() const
get the particle status
Definition: RawParticle.h:265
void rotateY(double rphi)
Definition: RawParticle.cc:168
double y() const
y of vertex
Definition: RawParticle.h:270
int myId
the particle id number HEP-PID
Definition: RawParticle.h:250
ROOT::Math::Rotation3D Rotation3D
Definition: RawParticle.h:35
void printName() const
Definition: RawParticle.cc:199
double z() const
z of vertex
Definition: RawParticle.h:271
double PDGmass() const
get the THEORETICAL mass
Definition: RawParticle.cc:254
double PDGcTau() const
get the THEORETICAL lifetime
Definition: RawParticle.cc:263
double myMass
the RECONSTRUCTED mass
Definition: RawParticle.h:254
Geom::Theta< T > theta() const
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:264
void chargeConjugate()
Definition: RawParticle.cc:143
void setT(const double t)
set the time of creation
Definition: RawParticle.cc:149
ParticleTable * tab
Definition: RawParticle.h:258
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:187
double mass() const
get the MEASURED mass
Definition: RawParticle.h:282
ROOT::Math::Boost Boost
Definition: RawParticle.h:39
double x() const
x of vertex
Definition: RawParticle.h:269
ROOT::Math::RotationZ RotationZ
Definition: RawParticle.h:38
void setMass(float m)
set the RECONSTRUCTED mass
Definition: RawParticle.cc:133
int isUsed() const
Definition: RawParticle.h:233
int myStatus
the status code according to PYTHIA
Definition: RawParticle.h:251
double R() const
vertex radius
Definition: RawParticle.h:277
double PDGcharge() const
get the THEORETICAL charge
Definition: RawParticle.cc:245
double R2() const
vertex radius**2
Definition: RawParticle.h:278
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
ROOT::Math::RotationX RotationX
Definition: RawParticle.h:36
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:285
virtual ~RawParticle()
Definition: RawParticle.cc:68
double t() const
vertex time
Definition: RawParticle.h:272
double r2() const
vertex radius**2
Definition: RawParticle.h:280
math::XYZVector XYZVector
double myCharge
the MEASURED charge
Definition: RawParticle.h:253
T sqrt(T t)
Definition: SSEVec.h:48
void print() const
Definition: RawParticle.cc:211
double cos2ThetaV() const
Definition: RawParticle.h:268
double Y() const
y of vertex
Definition: RawParticle.h:274
void setID(const int id)
Definition: RawParticle.cc:101
double Z() const
z of vertex
Definition: RawParticle.h:275
ROOT::Math::RotationY RotationY
Definition: RawParticle.h:37
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
void rotate(double rphi, const XYZVector &raxis)
Definition: RawParticle.cc:154
double charge() const
get the MEASURED charge
Definition: RawParticle.h:281
HepPDT::ParticleData ParticleData
int myUsed
status of the locking
Definition: RawParticle.h:252
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:284
void setStatus(int istat)
Definition: RawParticle.cc:128
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
double b
Definition: hdecay.h:120
void rotateX(double rphi)
Definition: RawParticle.cc:161
double r() const
vertex radius
Definition: RawParticle.h:279
double X() const
x of vertex
Definition: RawParticle.h:273
std::string PDGname() const
get the PDG name
Definition: RawParticle.cc:188
double a
Definition: hdecay.h:121
ROOT::Math::AxisAngle Rotation
Definition: RawParticle.h:34
double eta() const
Definition: RawParticle.h:266
double T() const
vertex time
Definition: RawParticle.h:276
const ParticleData * myInfo
The pointer to the PDG info.
Definition: RawParticle.h:255
RawParticle & operator=(const RawParticle &rhs)
Definition: RawParticle.cc:73
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:287
double cos2Theta() const
Cos**2(theta) is faster to determine than eta.
Definition: RawParticle.h:267
void init()
Definition: RawParticle.cc:90
double et() const
get the transverse energy
Definition: RawParticle.cc:306
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
void reUse()
Definition: RawParticle.h:241