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 #include<memory>
21 
22 
31 
33 public:
34 
35  typedef ROOT::Math::AxisAngle Rotation;
36  typedef ROOT::Math::Rotation3D Rotation3D;
37  typedef ROOT::Math::RotationX RotationX;
38  typedef ROOT::Math::RotationY RotationY;
39  typedef ROOT::Math::RotationZ RotationZ;
40  typedef ROOT::Math::Boost Boost;
41 
42  RawParticle();
43 
44  virtual ~RawParticle();
45 
50 
54  RawParticle(const int id,
55  const XYZTLorentzVector& p);
56 
61  const XYZTLorentzVector& p);
62 
67  const XYZTLorentzVector& xStart);
68 
72  RawParticle(double px, double py, double pz, double e);
73 
75  RawParticle(const RawParticle &p);
76 
78  RawParticle& operator = (const RawParticle & rhs );
79 
80 public:
81 
86  void setID(const int id);
87 
92  void setID(const std::string name);
93 
98  void setStatus(int istat);
99 
101  void setMass(float m);
102 
104  void setCharge(float q);
105 
107  void setT(const double t);
108 
110  void setVertex(const XYZTLorentzVector& vtx);
111  void setVertex(double xv, double yv, double zv, double tv);
112 
113  /*** methods to be overloaded to include vertex ***/
114 
119  void boost(double bx, double by, double bz);
120  void boost(const Boost& b);
121 
122  // inline void boost(const Hep3Vector<double> &bv );
123 
129  void rotate(double rphi, const XYZVector& raxis);
130  void rotate(const Rotation& r);
131  void rotate(const Rotation3D& r);
132 
134  // void rotateUz(Hep3Vector &nuz);
135 
139  void rotateX(double rphi);
140  void rotate(const RotationX& r);
141 
146  void rotateY(double rphi);
147  void rotate(const RotationY& r);
152  void rotateZ(double rphi);
153  void rotate(const RotationZ& r);
154 
156  void translate(const XYZVector& t);
157 
158  // inline RawParticle & transform(const HepRotation &rot);
159  // inline RawParticle & transform(const HepLorentzRotation &rot);
160 
165  void chargeConjugate();
166 
167  int pid() const;
168 
169  int status() const;
170 
171  double charge() const;
172 
173  double PDGcharge() const;
174 
175  double mass() const;
176 
177  double PDGmass() const;
178 
179  double PDGcTau() const;
180 
182  std::string PDGname() const;
183 
187  double eta() const;
188 
190  double cos2Theta() const;
191  double cos2ThetaV() const;
192 
193  double et() const;
194 
195  double x() const;
196  double X() const;
197 
198  double y() const;
199  double Y() const;
200 
201  double z() const;
202  double Z() const;
203 
204  double t() const;
205  double T() const;
206 
207  double r() const;
208  double R() const;
209 
210  double r2() const;
211  double R2() const;
212 
213  const XYZTLorentzVector& vertex() const;
214 
215  const XYZTLorentzVector& momentum() const;
216 
222  void printName() const;
223 
228  void print() const;
229 
234  int isUsed() const {return myUsed;}
235 
238  void use() { myUsed = 1;}
239 
242  void reUse() { myUsed = 0;}
243 
244  private:
245 
246  void init();
247 
248  protected:
249 
251  int myId;
252  int myStatus;
253  int myUsed;
254  double myCharge;
255  double myMass;
257 
258  private:
260 };
261 
262 
263 std::ostream& operator <<(std::ostream& o , const RawParticle& p);
264 
265 inline int RawParticle::pid() const { return myId; }
266 inline int RawParticle::status() const { return myStatus; }
267 inline double RawParticle::eta() const { return -std::log(std::tan(this->theta()/2.)); }
268 inline double RawParticle::cos2Theta() const { return Pz()*Pz()/Vect().Mag2(); }
269 inline double RawParticle::cos2ThetaV() const { return Z()*Z()/myVertex.Vect().Mag2(); }
270 inline double RawParticle::x() const { return myVertex.Px(); }
271 inline double RawParticle::y() const { return myVertex.Py(); }
272 inline double RawParticle::z() const { return myVertex.Pz(); }
273 inline double RawParticle::t() const { return myVertex.E(); }
274 inline double RawParticle::X() const { return myVertex.Px(); }
275 inline double RawParticle::Y() const { return myVertex.Py(); }
276 inline double RawParticle::Z() const { return myVertex.Pz(); }
277 inline double RawParticle::T() const { return myVertex.E(); }
278 inline double RawParticle::R() const { return std::sqrt(R2()); }
279 inline double RawParticle::R2() const { return myVertex.Perp2(); }
280 inline double RawParticle::r() const { return std::sqrt(r2()); }
281 inline double RawParticle::r2() const { return myVertex.Perp2(); }
282 inline double RawParticle::charge() const { return myCharge; }
283 inline double RawParticle::mass() const { return myMass; }
284 
285 inline const XYZTLorentzVector& RawParticle::vertex() const { return myVertex ; }
286 inline const XYZTLorentzVector& RawParticle::momentum() const { return *this; }
287 
288 inline void RawParticle::setVertex(const XYZTLorentzVector& vtx) { myVertex = vtx; }
289 inline void RawParticle::setVertex(double a, double b, double c, double d) { myVertex.SetXYZT(a,b,c,d); }
290 
292  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
293 }
294 
296  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
297 }
298 
300  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
301 }
302 
304  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
305 }
306 
308  XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
309 }
310 
312  XYZTLorentzVector p ( b(momentum()) ); SetXYZT(p.Px(),p.Py(),p.Pz(),p.E());
313 }
314 
315 inline void RawParticle::translate(const XYZVector& tr) {
316  myVertex.SetXYZT(X()+tr.X(),Y()+tr.Y(),Z()+tr.Z(),T());
317 }
318 
319 
320 
321 
322 
323 #endif
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:139
void translate(const XYZVector &t)
Definition: RawParticle.h:315
const ParticleTable * tab
Definition: RawParticle.h:259
void boost(double bx, double by, double bz)
Definition: RawParticle.cc:183
void rotateZ(double rphi)
Definition: RawParticle.cc:176
void use()
Definition: RawParticle.h:238
XYZTLorentzVector myVertex
the four vector of the vertex
Definition: RawParticle.h:250
int status() const
get the particle status
Definition: RawParticle.h:266
void rotateY(double rphi)
Definition: RawParticle.cc:169
double y() const
y of vertex
Definition: RawParticle.h:271
int myId
the particle id number HEP-PID
Definition: RawParticle.h:251
ROOT::Math::Rotation3D Rotation3D
Definition: RawParticle.h:36
void printName() const
Definition: RawParticle.cc:200
double z() const
z of vertex
Definition: RawParticle.h:272
double PDGmass() const
get the THEORETICAL mass
Definition: RawParticle.cc:255
double PDGcTau() const
get the THEORETICAL lifetime
Definition: RawParticle.cc:264
double myMass
the RECONSTRUCTED mass
Definition: RawParticle.h:255
Geom::Theta< T > theta() const
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:265
void chargeConjugate()
Definition: RawParticle.cc:144
void setT(const double t)
set the time of creation
Definition: RawParticle.cc:150
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:187
double mass() const
get the MEASURED mass
Definition: RawParticle.h:283
ROOT::Math::Boost Boost
Definition: RawParticle.h:40
double x() const
x of vertex
Definition: RawParticle.h:270
ROOT::Math::RotationZ RotationZ
Definition: RawParticle.h:39
void setMass(float m)
set the RECONSTRUCTED mass
Definition: RawParticle.cc:134
int isUsed() const
Definition: RawParticle.h:234
int myStatus
the status code according to PYTHIA
Definition: RawParticle.h:252
double R() const
vertex radius
Definition: RawParticle.h:278
double PDGcharge() const
get the THEORETICAL charge
Definition: RawParticle.cc:246
double R2() const
vertex radius**2
Definition: RawParticle.h:279
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
ROOT::Math::RotationX RotationX
Definition: RawParticle.h:37
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:286
virtual ~RawParticle()
Definition: RawParticle.cc:69
double t() const
vertex time
Definition: RawParticle.h:273
double r2() const
vertex radius**2
Definition: RawParticle.h:281
math::XYZVector XYZVector
double myCharge
the MEASURED charge
Definition: RawParticle.h:254
T sqrt(T t)
Definition: SSEVec.h:48
void print() const
Definition: RawParticle.cc:212
double cos2ThetaV() const
Definition: RawParticle.h:269
double Y() const
y of vertex
Definition: RawParticle.h:275
void setID(const int id)
Definition: RawParticle.cc:102
double Z() const
z of vertex
Definition: RawParticle.h:276
ROOT::Math::RotationY RotationY
Definition: RawParticle.h:38
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
void rotate(double rphi, const XYZVector &raxis)
Definition: RawParticle.cc:155
double charge() const
get the MEASURED charge
Definition: RawParticle.h:282
HepPDT::ParticleData ParticleData
int myUsed
status of the locking
Definition: RawParticle.h:253
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
void setStatus(int istat)
Definition: RawParticle.cc:129
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:162
double r() const
vertex radius
Definition: RawParticle.h:280
double X() const
x of vertex
Definition: RawParticle.h:274
std::string PDGname() const
get the PDG name
Definition: RawParticle.cc:189
double a
Definition: hdecay.h:121
ROOT::Math::AxisAngle Rotation
Definition: RawParticle.h:35
double eta() const
Definition: RawParticle.h:267
double T() const
vertex time
Definition: RawParticle.h:277
const ParticleData * myInfo
The pointer to the PDG info.
Definition: RawParticle.h:256
RawParticle & operator=(const RawParticle &rhs)
Definition: RawParticle.cc:74
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:288
double cos2Theta() const
Cos**2(theta) is faster to determine than eta.
Definition: RawParticle.h:268
void init()
Definition: RawParticle.cc:91
double et() const
get the transverse energy
Definition: RawParticle.cc:307
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
void reUse()
Definition: RawParticle.h:242