CMS 3D CMS Logo

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