CMS 3D CMS Logo

CaloParticle.h
Go to the documentation of this file.
1 #ifndef SimDataFormats_CaloParticle_h
2 #define SimDataFormats_CaloParticle_h
3 
11 #include <vector>
12 
13 class SimTrack;
14 class EncodedEventId;
15 
16 class CaloParticle {
17  friend std::ostream &operator<<(std::ostream &s, CaloParticle const &tp);
18 
19 public:
20  typedef int Charge;
25 
28  typedef std::vector<SimTrack>::const_iterator g4t_iterator;
30 
31  CaloParticle();
32 
33  CaloParticle(const SimTrack &simtrk);
34  CaloParticle(EncodedEventId eventID, uint32_t particleID); // for PU
35 
36  // destructor
37  ~CaloParticle();
38 
43  int pdgId() const {
44  if (genParticles_.empty())
45  return g4Tracks_[0].type();
46  else
47  return (*genParticles_.begin())->pdgId();
48  }
49 
54  EncodedEventId eventId() const { return event_; }
55 
56  uint64_t particleId() const { return particleId_; }
57 
58  // Setters for G4 and reco::GenParticle
61  void addG4Track(const SimTrack &t) { g4Tracks_.push_back(t); }
65  g4t_iterator g4Track_begin() const { return g4Tracks_.begin(); }
66  g4t_iterator g4Track_end() const { return g4Tracks_.end(); }
69 
70  // Getters for Embd and Sim Tracks
72  const SimClusterRefVector &simClusters() const { return simClusters_; }
73  // Only for clusters from the signal vertex
74  const std::vector<SimTrack> &g4Tracks() const { return g4Tracks_; }
75 
77 
79  float charge() const { return g4Tracks_[0].charge(); }
81  int threeCharge() const { return lrintf(3.f * charge()); }
82 
85  const math::XYZTLorentzVectorF &p4() const { return theMomentum_; }
86 
88  math::XYZVectorF momentum() const { return p4().Vect(); }
89 
91  math::XYZVectorF boostToCM() const { return p4().BoostToCM(); }
92 
95  float p() const { return p4().P(); }
96 
98  float energy() const { return p4().E(); }
99 
101  float et() const { return p4().Et(); }
102 
104  float mass() const { return p4().M(); }
105 
107  float massSqr() const { return pow(mass(), 2); }
108 
110  float mt() const { return p4().Mt(); }
111 
114  float mtSqr() const { return p4().Mt2(); }
115 
118  float px() const { return p4().Px(); }
119 
122  float py() const { return p4().Py(); }
123 
126  float pz() const { return p4().Pz(); }
127 
130  float pt() const { return p4().Pt(); }
131 
134  float phi() const { return p4().Phi(); }
135 
138  float theta() const { return p4().Theta(); }
139 
142  float eta() const { return p4().Eta(); }
143 
146  float rapidity() const { return p4().Rapidity(); }
147 
149  float y() const { return rapidity(); }
150 
155  int status() const { return genParticles_.empty() ? -99 : (*genParticles_[0]).status(); }
156 
157  static const unsigned int longLivedTag;
158 
160  bool longLived() const { return status() & longLivedTag; }
161 
163  int numberOfSimHits() const { return nsimhits_; }
164 
166  int numberOfRecHits() const { return hits_.size(); }
167 
169  void addRecHitAndFraction(uint32_t hit, float fraction) {
170  hits_.emplace_back(hit);
171  fractions_.emplace_back(fraction);
172  }
173 
175  std::vector<std::pair<uint32_t, float>> hits_and_fractions() const {
176  std::vector<std::pair<uint32_t, float>> result;
177  for (size_t i = 0; i < hits_.size(); ++i) {
178  result.emplace_back(hits_[i], fractions_[i]);
179  }
180  return result;
181  }
182 
184  float simEnergy() const { return simhit_energy_; }
185 
187  void addSimHit(const PCaloHit &hit) {
188  simhit_energy_ += hit.energy();
189  ++nsimhits_;
190  }
191 
192 protected:
195 
196  uint32_t particleId_{0};
197  float simhit_energy_{0.f};
198  std::vector<uint32_t> hits_;
199  std::vector<float> fractions_;
200 
202 
204  std::vector<SimTrack> g4Tracks_;
206 
208 };
209 
210 #endif // SimDataFormats_CaloParticle_H
const reco::GenParticleRefVector & genParticles() const
Definition: CaloParticle.h:71
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
int Charge
electric charge type
Definition: CaloParticle.h:20
std::vector< uint32_t > hits_
Definition: CaloParticle.h:198
const std::vector< SimTrack > & g4Tracks() const
Definition: CaloParticle.h:74
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
uint64_t nsimhits_
Definition: CaloParticle.h:193
const SimClusterRefVector & simClusters() const
Definition: CaloParticle.h:72
bool longLived() const
is long lived?
Definition: CaloParticle.h:160
int status() const
Status word.
Definition: CaloParticle.h:155
SimClusterRefVector::iterator sc_iterator
Definition: CaloParticle.h:29
float px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: CaloParticle.h:118
void addG4Track(const SimTrack &t)
Definition: CaloParticle.h:61
void addRecHitAndFraction(uint32_t hit, float fraction)
add rechit with fraction
Definition: CaloParticle.h:169
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
Definition: CaloParticle.h:21
reco::GenParticleRefVector genParticles_
Definition: CaloParticle.h:205
std::vector< float > fractions_
Definition: CaloParticle.h:199
math::XYZVectorD Vector
point in the space
Definition: CaloParticle.h:24
float py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: CaloParticle.h:122
EncodedEventId event_
Definition: CaloParticle.h:194
float mt() const
Transverse mass. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:110
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: CaloParticle.h:22
SimClusterRefVector simClusters_
Definition: CaloParticle.h:207
std::vector< SimTrack >::const_iterator g4t_iterator
Definition: CaloParticle.h:28
math::XYZTLorentzVectorF theMomentum_
Definition: CaloParticle.h:201
math::XYZPointD Point
point in the space
Definition: CaloParticle.h:23
reco::GenParticleRefVector::iterator genp_iterator
reference to reco::GenParticle
Definition: CaloParticle.h:27
float rapidity() const
Rapidity. Note this is taken from the simtrack before the calorimeter.
Definition: CaloParticle.h:146
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:98
static const unsigned int longLivedTag
long lived flag
Definition: CaloParticle.h:157
uint32_t particleId_
Definition: CaloParticle.h:196
float pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: CaloParticle.h:126
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
int numberOfRecHits() const
Gives the total number of SimHits, in the cluster.
Definition: CaloParticle.h:166
void addGenParticle(const reco::GenParticleRef &ref)
Definition: CaloParticle.h:59
void addSimHit(const PCaloHit &hit)
add simhit&#39;s energy to cluster
Definition: CaloParticle.h:187
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:85
float pt() const
Transverse momentum. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:130
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: CaloParticle.h:184
math::XYZVectorF boostToCM() const
Vector to boost to the particle centre of mass frame.
Definition: CaloParticle.h:91
sc_iterator simCluster_end() const
Definition: CaloParticle.h:68
double f[11][100]
float theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:138
float simhit_energy_
Definition: CaloParticle.h:197
void addSimCluster(const SimClusterRef &ref)
Definition: CaloParticle.h:60
int threeCharge() const
Gives charge in unit of quark charge (should be 3 times "charge()")
Definition: CaloParticle.h:81
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
int pdgId() const
PDG ID.
Definition: CaloParticle.h:43
float phi() const
Momentum azimuthal angle. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:134
float et() const
Transverse energy. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:101
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
friend std::ostream & operator<<(std::ostream &s, CaloParticle const &tp)
Definition: CaloParticle.cc:30
unsigned long long uint64_t
Definition: Time.h:13
void clear()
Clear the vector.
Definition: RefVector.h:142
float p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:95
int numberOfSimHits() const
Gives the total number of SimHits, in the cluster.
Definition: CaloParticle.h:163
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
g4t_iterator g4Track_begin() const
Definition: CaloParticle.h:65
sc_iterator simCluster_begin() const
Definition: CaloParticle.h:67
genp_iterator genParticle_end() const
Definition: CaloParticle.h:64
float y() const
Same as rapidity().
Definition: CaloParticle.h:149
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: CaloParticle.h:204
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:79
float massSqr() const
Mass squared. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:107
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this CaloParticle.
Definition: CaloParticle.h:175
void clearSimClusters()
Definition: CaloParticle.h:76
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
g4t_iterator g4Track_end() const
Definition: CaloParticle.h:66
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
float mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:114
uint64_t particleId() const
Definition: CaloParticle.h:56
float eta() const
Momentum pseudorapidity. Note this is taken from the simtrack before the calorimeter.
Definition: CaloParticle.h:142
genp_iterator genParticle_begin() const
iterators
Definition: CaloParticle.h:63
math::XYZVectorF momentum() const
spatial momentum vector
Definition: CaloParticle.h:88
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
EncodedEventId eventId() const
Signal source, crossing number.
Definition: CaloParticle.h:54
float mass() const
Mass. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:104