CMS 3D CMS Logo

SimCluster.h
Go to the documentation of this file.
1 #ifndef SimDataFormats_SimCluster_h
2 #define SimDataFormats_SimCluster_h
3 
10 #include <vector>
11 
12 //
13 // Forward declarations
14 //
15 class SimTrack;
16 class EncodedEventId;
17 
29 class SimCluster {
30  friend std::ostream &operator<<(std::ostream &s, SimCluster const &tp);
31 
32 public:
33  typedef int Charge;
38 
41  typedef std::vector<SimTrack>::const_iterator g4t_iterator;
42 
43  SimCluster();
44 
45  SimCluster(const SimTrack &simtrk);
46  SimCluster(EncodedEventId eventID, uint32_t particleID); // for PU
47 
48  // destructor
49  ~SimCluster();
50 
55  int pdgId() const {
56  if (genParticles_.empty())
57  return g4Tracks_[0].type();
58  else
59  return (*genParticles_.begin())->pdgId();
60  }
61 
66  EncodedEventId eventId() const { return event_; }
67 
68  uint64_t particleId() const { return particleId_; }
69 
70  // Setters for G4 and reco::GenParticle
72  void addG4Track(const SimTrack &t) { g4Tracks_.push_back(t); }
74  genp_iterator genParticle_begin() const { return genParticles_.begin(); }
75  genp_iterator genParticle_end() const { return genParticles_.end(); }
76  g4t_iterator g4Track_begin() const { return g4Tracks_.begin(); }
77  g4t_iterator g4Track_end() const { return g4Tracks_.end(); }
78 
79  // Getters for Embd and Sim Tracks
81  // Only for clusters from the signal vertex
82  const std::vector<SimTrack> &g4Tracks() const { return g4Tracks_; }
83 
85  float charge() const { return g4Tracks_[0].charge(); }
87  int threeCharge() const { return lrintf(3.f * charge()); }
88 
91  const math::XYZTLorentzVectorF &p4() const { return theMomentum_; }
92 
94  math::XYZVectorF momentum() const { return p4().Vect(); }
95 
97  math::XYZVectorF boostToCM() const { return p4().BoostToCM(); }
98 
101  float p() const { return p4().P(); }
102 
104  float energy() const { return p4().E(); }
105 
107  float et() const { return p4().Et(); }
108 
110  float mass() const { return p4().M(); }
111 
113  float massSqr() const { return pow(mass(), 2); }
114 
116  float mt() const { return p4().Mt(); }
117 
120  float mtSqr() const { return p4().Mt2(); }
121 
124  float px() const { return p4().Px(); }
125 
128  float py() const { return p4().Py(); }
129 
132  float pz() const { return p4().Pz(); }
133 
136  float pt() const { return p4().Pt(); }
137 
140  float phi() const { return p4().Phi(); }
141 
144  float theta() const { return p4().Theta(); }
145 
148  float eta() const { return p4().Eta(); }
149 
152  float rapidity() const { return p4().Rapidity(); }
153 
155  float y() const { return rapidity(); }
156 
161  int status() const { return genParticles_.empty() ? -99 : (*genParticles_[0]).status(); }
162 
163  static const unsigned int longLivedTag;
164 
166  bool longLived() const { return status() & longLivedTag; }
167 
169  int numberOfSimHits() const { return nsimhits_; }
170 
172  int numberOfRecHits() const { return hits_.size(); }
173 
175  void addRecHitAndFraction(uint32_t hit, float fraction) {
176  hits_.emplace_back(hit);
177  fractions_.emplace_back(fraction);
178  }
179 
181  std::vector<std::pair<uint32_t, float>> hits_and_fractions() const {
182  std::vector<std::pair<uint32_t, float>> result;
183  for (size_t i = 0; i < hits_.size(); ++i) {
184  result.emplace_back(hits_[i], fractions_[i]);
185  }
186  return result;
187  }
188 
191  std::vector<uint32_t>().swap(hits_);
192  std::vector<float>().swap(fractions_);
193  }
194 
196  float simEnergy() const { return simhit_energy_; }
197 
199  void addSimHit(const PCaloHit &hit) { simhit_energy_ += hit.energy(); }
200 
201 private:
204 
205  uint32_t particleId_;
207  std::vector<uint32_t> hits_;
208  std::vector<float> fractions_;
209 
211 
213  std::vector<SimTrack> g4Tracks_;
215 };
216 
217 #endif // SimDataFormats_SimCluster_H
type
Definition: HCALResponse.h:21
float y() const
Same as rapidity().
Definition: SimCluster.h:155
void addG4Track(const SimTrack &t)
Definition: SimCluster.h:72
float mt() const
Transverse mass. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:116
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
Definition: SimCluster.h:34
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:140
float massSqr() const
Mass squared. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:113
float simhit_energy_
Definition: SimCluster.h:206
float et() const
Transverse energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:107
int threeCharge() const
Gives charge in unit of quark charge (should be 3 times "charge()")
Definition: SimCluster.h:87
math::XYZPointD Point
point in the space
Definition: SimCluster.h:36
double energy() const
Definition: PCaloHit.h:24
float rapidity() const
Rapidity. Note this is taken from the simtrack before the calorimeter.
Definition: SimCluster.h:152
int pdgId() const
PDG ID.
Definition: SimCluster.h:55
math::XYZVectorD Vector
point in the space
Definition: SimCluster.h:37
friend std::ostream & operator<<(std::ostream &s, SimCluster const &tp)
Definition: SimCluster.cc:31
EncodedEventId event_
Definition: SimCluster.h:203
std::vector< SimTrack >::const_iterator g4t_iterator
Definition: SimCluster.h:41
uint64_t particleId() const
Definition: SimCluster.h:68
int Charge
electric charge type
Definition: SimCluster.h:33
void addGenParticle(const reco::GenParticleRef &ref)
Definition: SimCluster.h:71
EncodedEventId eventId() const
Signal source, crossing number.
Definition: SimCluster.h:66
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:175
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:202
genp_iterator genParticle_end() const
Definition: SimCluster.h:75
std::vector< uint32_t > hits_
Definition: SimCluster.h:207
float pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: SimCluster.h:132
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:136
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:29
float py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: SimCluster.h:128
math::XYZVectorF momentum() const
spatial momentum vector
Definition: SimCluster.h:94
float p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:101
bool longLived() const
is long lived?
Definition: SimCluster.h:166
void clearHitsAndFractions()
clear the hits and fractions list
Definition: SimCluster.h:190
reco::GenParticleRefVector::iterator genp_iterator
reference to reco::GenParticle
Definition: SimCluster.h:40
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
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:208
float theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:144
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: SimCluster.h:196
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:85
float px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
Definition: SimCluster.h:124
int numberOfSimHits() const
Gives the total number of SimHits, in the cluster.
Definition: SimCluster.h:169
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: SimCluster.h:213
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:161
float mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:120
const std::vector< SimTrack > & g4Tracks() const
Definition: SimCluster.h:82
void addSimHit(const PCaloHit &hit)
add simhit&#39;s energy to cluster
Definition: SimCluster.h:199
uint32_t particleId_
Definition: SimCluster.h:205
g4t_iterator g4Track_begin() const
Definition: SimCluster.h:76
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:110
const reco::GenParticleRefVector & genParticles() const
Definition: SimCluster.h:80
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:214
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:97
int numberOfRecHits() const
Gives the total number of SimHits, in the cluster.
Definition: SimCluster.h:172
genp_iterator genParticle_begin() const
iterators
Definition: SimCluster.h:74
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
g4t_iterator g4Track_end() const
Definition: SimCluster.h:77
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:181
static const unsigned int longLivedTag
long lived flag
Definition: SimCluster.h:163
float eta() const
Momentum pseudorapidity. Note this is taken from the simtrack before the calorimeter.
Definition: SimCluster.h:148
math::XYZTLorentzVectorF theMomentum_
Definition: SimCluster.h:210
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: SimCluster.h:35