CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ParticleState.h
Go to the documentation of this file.
1 #ifndef Candidate_ParticleState_h
2 #define Candidate_ParticleState_h
3 
18 #include "Rtypes.h"
19 
20 namespace reco {
21 
22  class ParticleState {
23  public:
25  typedef int Charge;
35  ParticleState() : vertex_(0, 0, 0), qx3_(0), pdgId_(0), status_(0) {}
36 
39  const PtEtaPhiMass& p4,
40  const Point& vertex = Point(0, 0, 0),
41  int pdgId = 0,
42  int status = 0,
43  bool integerCharge = true)
44  : vertex_(vertex),
45  p4Polar_(p4.pt(), p4.eta(), p4.phi(), p4.mass()),
47  qx3_(integerCharge ? q * 3 : q),
48  pdgId_(pdgId),
49  status_(status) {}
50 
53  const LorentzVector& p4,
54  const Point& vertex = Point(0, 0, 0),
55  int pdgId = 0,
56  int status = 0,
57  bool integerCharge = true)
58  : vertex_(vertex),
59  p4Polar_(p4),
60  p4Cartesian_(p4),
61  qx3_(integerCharge ? q * 3 : q),
62  pdgId_(pdgId),
63  status_(status) {}
64 
67  const PolarLorentzVector& p4,
68  const Point& vertex = Point(0, 0, 0),
69  int pdgId = 0,
70  int status = 0,
71  bool integerCharge = true)
72  : vertex_(vertex),
73  p4Polar_(p4),
74  p4Cartesian_(p4),
75  qx3_(integerCharge ? q * 3 : q),
76  pdgId_(pdgId),
77  status_(status) {}
78 
80  const GlobalVector& p3,
81  float iEnergy,
82  float imass,
83  const Point& vertex = Point(0, 0, 0),
84  int pdgId = 0,
85  int status = 0,
86  bool integerCharge = true)
87  : vertex_(vertex),
88  p4Polar_(p3.perp(), p3.eta(), p3.phi(), imass),
89  p4Cartesian_(p3.x(), p3.y(), p3.z(), iEnergy),
90  qx3_(integerCharge ? q * 3 : q),
91  pdgId_(pdgId),
92  status_(status) {}
93 
95  inline void setCartesian() { p4Cartesian_ = p4Polar_; }
96 
98  int charge() const { return qx3_ / 3; }
100  void setCharge(Charge q) { qx3_ = q * 3; }
102  int threeCharge() const { return qx3_; }
104  void setThreeCharge(Charge qx3) { qx3_ = qx3; }
106  const LorentzVector& p4() const { return p4Cartesian_; }
108  const PolarLorentzVector& polarP4() const { return p4Polar_; }
110  Vector momentum() const { return p4Cartesian_.Vect(); }
113  Vector boostToCM() const { return p4Cartesian_.BoostToCM(); }
115  double p() const { return p4Cartesian_.P(); }
117  double energy() const { return p4Cartesian_.E(); }
119  double et() const { return (pt() <= 0) ? 0 : p4Cartesian_.Et(); }
121  double et2() const { return (pt() <= 0) ? 0 : p4Cartesian_.Et2(); }
123  double mass() const { return p4Polar_.mass(); }
125  double massSqr() const { return mass() * mass(); }
127  double mt() const { return p4Polar_.Mt(); }
129  double mtSqr() const { return p4Polar_.Mt2(); }
131  double px() const { return p4Cartesian_.Px(); }
133  double py() const { return p4Cartesian_.Py(); }
135  double pz() const { return p4Cartesian_.Pz(); }
137  double pt() const { return p4Polar_.pt(); }
139  double phi() const { return p4Polar_.phi(); }
141  double theta() const { return p4Cartesian_.Theta(); }
143  double eta() const { return p4Polar_.eta(); }
145  double rapidity() const { return p4Polar_.Rapidity(); }
147  double y() const { return rapidity(); }
149  void setP4(const LorentzVector& p4) {
150  p4Cartesian_ = p4;
151  p4Polar_ = p4;
152  }
154  void setP4(const PolarLorentzVector& p4) {
155  p4Polar_ = p4;
156  p4Cartesian_ = p4;
157  }
159  void setMass(double m) {
160  p4Polar_.SetM(m);
161  setCartesian();
162  }
163  void setPz(double pz) {
164  p4Cartesian_.SetPz(pz);
166  }
168  const Point& vertex() const { return vertex_; }
170  double vx() const { return vertex_.X(); }
172  double vy() const { return vertex_.Y(); }
174  double vz() const { return vertex_.Z(); }
176  void setVertex(const Point& vertex) { vertex_ = vertex; }
178  int pdgId() const { return pdgId_; }
179  // set PDG identifier
180  void setPdgId(int pdgId) { pdgId_ = pdgId; }
182  int status() const { return status_; }
184  void setStatus(int status) { status_ = status; }
188  bool longLived() const { return status_ & longLivedTag; }
192  bool massConstraint() const { return status_ & massConstraintTag; }
193 
194  private:
195  static const unsigned int longLivedTag = 65536;
196  static const unsigned int massConstraintTag = 131072;
197 
198  private:
201 
206 
209 
211  int pdgId_;
213  int status_;
214  };
215 
216 } // namespace reco
217 
218 #endif
double rapidity() const
repidity
ParticleState(Charge q, const PolarLorentzVector &p4, const Point &vertex=Point(0, 0, 0), int pdgId=0, int status=0, bool integerCharge=true)
constructor from values
Definition: ParticleState.h:66
double px() const
x coordinate of momentum vector
int pdgId_
PDG identifier.
Point vertex_
vertex position
int pdgId() const
PDG identifier.
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: ParticleState.h:29
void setCartesian()
set internal cache
Definition: ParticleState.h:95
double theta() const
momentum polar angle
LorentzVector p4Cartesian_
internal cache for p4
double pz() const
z coordinate of momentum vector
int threeCharge() const
electric charge
void setP4(const PolarLorentzVector &p4)
set 4-momentum
static const unsigned int longLivedTag
int Charge
electric charge type
Definition: ParticleState.h:25
math::XYZPoint Point
point in the space
Definition: ParticleState.h:31
void setMassConstraint()
set mass constraint flag
double energy() const
energy
void setPdgId(int pdgId)
void setPz(double pz)
PolarLorentzVector p4Polar_
four-momentum Lorentz vector
double y() const
repidity
double phi() const
momentum azimuthal angle
bool massConstraint() const
do mass constraint?
const LorentzVector & p4() const
four-momentum Lorentz vector
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double pt() const
transverse momentum
void setStatus(int status)
set status word
ParticleState(Charge q, const LorentzVector &p4, const Point &vertex=Point(0, 0, 0), int pdgId=0, int status=0, bool integerCharge=true)
constructor from values
Definition: ParticleState.h:52
double mt() const
transverse mass
double mass() const
mass
math::XYZVector Vector
point in the space
Definition: ParticleState.h:33
bool longLived() const
is long lived?
void setVertex(const Point &vertex)
set vertex
double mtSqr() const
transverse mass squared
double massSqr() const
mass squared
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: ParticleState.h:27
Vector boostToCM() const
void setThreeCharge(Charge qx3)
set electric charge
double et2() const
transverse energy squared (use this for cuts)!
ParticleState(Charge q, const GlobalVector &p3, float iEnergy, float imass, const Point &vertex=Point(0, 0, 0), int pdgId=0, int status=0, bool integerCharge=true)
Definition: ParticleState.h:79
void setLongLived()
set long lived flag
const PolarLorentzVector & polarP4() const
four-momentum Lorentz vector
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:28
double et() const
transverse energy
double vy() const
y coordinate of vertex position
double vz() const
z coordinate of vertex position
int charge() const
electric charge
Definition: ParticleState.h:98
Charge qx3_
electric charge
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
int status() const
status word
const Point & vertex() const
vertex position
ParticleState(Charge q, const PtEtaPhiMass &p4, const Point &vertex=Point(0, 0, 0), int pdgId=0, int status=0, bool integerCharge=true)
constructor from values
Definition: ParticleState.h:38
T perp() const
Magnitude of transverse component.
static const unsigned int massConstraintTag
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
Vector momentum() const
spatial momentum vector
double eta() const
momentum pseudorapidity
void setP4(const LorentzVector &p4)
set 4-momentum
ParticleState()
default constructor
Definition: ParticleState.h:35
double p() const
magnitude of momentum vector
double py() const
y coordinate of momentum vector
void setCharge(Charge q)
set electric charge
double vx() const
x coordinate of vertex position
int status_
status word
void setMass(double m)
set particle mass