CMS 3D CMS Logo

CaloParticle.h
Go to the documentation of this file.
1 #ifndef SimDataFormats_CaloParticle_h
2 #define SimDataFormats_CaloParticle_h
3 
4 #include <vector>
12 
13 class SimTrack;
14 class EncodedEventId;
15 
17 {
18  friend std::ostream& operator<< (std::ostream& s, CaloParticle const& tp);
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() ) return g4Tracks_[0].type();
45  else return (*genParticles_.begin())->pdgId();
46  }
47 
53  return event_;
54  }
55 
56  uint64_t particleId() const {
57  return particleId_;
58  }
59 
60  // Setters for G4 and reco::GenParticle
62  void addSimCluster( const SimClusterRef& ref) { simClusters_.push_back( ref ); }
63  void addG4Track( const SimTrack& t) { g4Tracks_.push_back(t); }
65  genp_iterator genParticle_begin() const { return genParticles_.begin(); }
66  genp_iterator genParticle_end() const { return genParticles_.end(); }
67  g4t_iterator g4Track_begin() const { return g4Tracks_.begin(); }
68  g4t_iterator g4Track_end() const { return g4Tracks_.end(); }
69  sc_iterator simCluster_begin() const { return simClusters_.begin(); }
70  sc_iterator simCluster_end() const { return simClusters_.end(); }
71 
72  // Getters for Embd and Sim Tracks
74  const SimClusterRefVector& simClusters() const { return simClusters_; }
75  // Only for clusters from the signal vertex
76  const std::vector<SimTrack>& g4Tracks() const { return g4Tracks_; }
77 
79 
81  float charge() const { return g4Tracks_[0].charge(); }
83  int threeCharge() const { return lrintf(3.f*charge()); }
84 
86  const math::XYZTLorentzVectorF& p4() const {
87  return theMomentum_;
88  }
89 
92  return p4().Vect();
93  }
94 
97  return p4().BoostToCM();
98  }
99 
101  float p() const {
102  return p4().P();
103  }
104 
106  float energy() const {
107  return p4().E();
108  }
109 
111  float et() const {
112  return p4().Et();
113  }
114 
116  float mass() const {
117  return p4().M();
118  }
119 
121  float massSqr() const {
122  return pow( mass(), 2 );
123  }
124 
126  float mt() const {
127  return p4().Mt();
128  }
129 
131  float mtSqr() const {
132  return p4().Mt2();
133  }
134 
136  float px() const {
137  return p4().Px();
138  }
139 
141  float py() const {
142  return p4().Py();
143  }
144 
146  float pz() const {
147  return p4().Pz();
148  }
149 
151  float pt() const {
152  return p4().Pt();
153  }
154 
156  float phi() const {
157  return p4().Phi();
158  }
159 
161  float theta() const {
162  return p4().Theta();
163  }
164 
166  float eta() const {
167  return p4().Eta();
168  }
169 
171  float rapidity() const {
172  return p4().Rapidity();
173  }
174 
176  float y() const {
177  return rapidity();
178  }
179 
183  int status() const {
184  return genParticles_.empty() ? -99 : (*genParticles_[0]).status();
185  }
186 
187  static const unsigned int longLivedTag;
188 
190  bool longLived() const { return status()&longLivedTag;}
191 
193  int numberOfSimHits() const {return nsimhits_;}
194 
196  int numberOfRecHits() const {return hits_.size();}
197 
199  void addRecHitAndFraction(uint32_t hit, float fraction) {
200  hits_.emplace_back(hit);
201  fractions_.emplace_back(fraction);
202  }
203 
205  std::vector<std::pair<uint32_t,float> > hits_and_fractions() const {
206  std::vector<std::pair<uint32_t,float> > result;
207  for(size_t i = 0; i < hits_.size(); ++i) {
208  result.emplace_back(hits_[i],fractions_[i]);
209  }
210  return result;
211  }
212 
214  float simEnergy() const { return simhit_energy_; }
215 
217  void addSimHit(const PCaloHit& hit) { simhit_energy_ += hit.energy(); }
218 
219 private:
222 
223  uint32_t particleId_;
225  std::vector<uint32_t> hits_;
226  std::vector<float> fractions_;
227 
229 
231  std::vector<SimTrack> g4Tracks_;
233 
235 };
236 
237 #endif // SimDataFormats_CaloParticle_H
float pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: CaloParticle.h:146
type
Definition: HCALResponse.h:21
int Charge
electric charge type
Definition: CaloParticle.h:20
std::vector< uint32_t > hits_
Definition: CaloParticle.h:225
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: CaloParticle.h:214
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:220
double energy() const
Definition: PCaloHit.h:29
SimClusterRefVector::iterator sc_iterator
Definition: CaloParticle.h:29
void addG4Track(const SimTrack &t)
Definition: CaloParticle.h:63
void addRecHitAndFraction(uint32_t hit, float fraction)
add rechit with fraction
Definition: CaloParticle.h:199
EncodedEventId eventId() const
Signal source, crossing number.
Definition: CaloParticle.h:52
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:81
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
Definition: CaloParticle.h:21
float p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:101
reco::GenParticleRefVector genParticles_
Definition: CaloParticle.h:232
int status() const
Status word.
Definition: CaloParticle.h:183
std::vector< float > fractions_
Definition: CaloParticle.h:226
math::XYZVectorD Vector
point in the space
Definition: CaloParticle.h:24
EncodedEventId event_
Definition: CaloParticle.h:221
int numberOfSimHits() const
Gives the total number of SimHits, in the cluster.
Definition: CaloParticle.h:193
const std::vector< SimTrack > & g4Tracks() const
Definition: CaloParticle.h:76
float eta() const
Momentum pseudorapidity. Note this is taken from the simtrack before the calorimeter.
Definition: CaloParticle.h:166
float massSqr() const
Mass squared. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:121
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:253
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:104
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: CaloParticle.h:22
SimClusterRefVector simClusters_
Definition: CaloParticle.h:234
std::vector< SimTrack >::const_iterator g4t_iterator
Definition: CaloParticle.h:28
genp_iterator genParticle_begin() const
iterators
Definition: CaloParticle.h:65
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:248
math::XYZTLorentzVectorF theMomentum_
Definition: CaloParticle.h:228
int numberOfRecHits() const
Gives the total number of SimHits, in the cluster.
Definition: CaloParticle.h:196
uint64_t particleId() const
Definition: CaloParticle.h:56
math::XYZPointD Point
point in the space
Definition: CaloParticle.h:23
reco::GenParticleRefVector::iterator genp_iterator
reference to reco::GenParticle
Definition: CaloParticle.h:27
sc_iterator simCluster_end() const
Definition: CaloParticle.h:70
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:86
static const unsigned int longLivedTag
long lived flag
Definition: CaloParticle.h:187
uint32_t particleId_
Definition: CaloParticle.h:223
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
sc_iterator simCluster_begin() const
Definition: CaloParticle.h:69
const reco::GenParticleRefVector & genParticles() const
Definition: CaloParticle.h:73
const SimClusterRefVector & simClusters() const
Definition: CaloParticle.h:74
void addGenParticle(const reco::GenParticleRef &ref)
Definition: CaloParticle.h:61
void addSimHit(const PCaloHit &hit)
add simhit&#39;s energy to cluster
Definition: CaloParticle.h:217
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:106
float theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:161
double f[11][100]
float simhit_energy_
Definition: CaloParticle.h:224
void addSimCluster(const SimClusterRef &ref)
Definition: CaloParticle.h:62
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:17
math::XYZVectorF boostToCM() const
Vector to boost to the particle centre of mass frame.
Definition: CaloParticle.h:96
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this CaloParticle.
Definition: CaloParticle.h:205
float mt() const
Transverse mass. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:126
float y() const
Same as rapidity().
Definition: CaloParticle.h:176
math::XYZVectorF momentum() const
spatial momentum vector
Definition: CaloParticle.h:91
float px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: CaloParticle.h:136
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:37
unsigned long long uint64_t
Definition: Time.h:15
void clear()
Clear the vector.
Definition: RefVector.h:147
float pt() const
Transverse momentum. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:151
int pdgId() const
PDG ID.
Definition: CaloParticle.h:43
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: CaloParticle.h:231
float phi() const
Momentum azimuthal angle. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:156
float rapidity() const
Rapidity. Note this is taken from the simtrack before the calorimeter.
Definition: CaloParticle.h:171
g4t_iterator g4Track_end() const
Definition: CaloParticle.h:68
bool longLived() const
is long lived?
Definition: CaloParticle.h:190
float mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:131
float mass() const
Mass. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:116
g4t_iterator g4Track_begin() const
Definition: CaloParticle.h:67
void clearSimClusters()
Definition: CaloParticle.h:78
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
int threeCharge() const
Gives charge in unit of quark charge (should be 3 times "charge()")
Definition: CaloParticle.h:83
float et() const
Transverse energy. Note this is taken from the first SimTrack only.
Definition: CaloParticle.h:111
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
float py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: CaloParticle.h:141
genp_iterator genParticle_end() const
Definition: CaloParticle.h:66
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40