CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Attributes | Private Attributes | Friends
SimCluster Class Reference

Monte Carlo truth information used for tracking validation. More...

#include <SimCluster.h>

Public Types

typedef int Charge
 electric charge type More...
 
typedef std::vector< SimTrack >::const_iterator g4t_iterator
 
typedef reco::GenParticleRefVector::iterator genp_iterator
 reference to reco::GenParticle More...
 
typedef math::XYZTLorentzVectorD LorentzVector
 Lorentz vector. More...
 
typedef math::XYZPointD Point
 point in the space More...
 
typedef math::PtEtaPhiMLorentzVector PolarLorentzVector
 Lorentz vector. More...
 
typedef math::XYZVectorD Vector
 point in the space More...
 

Public Member Functions

void addG4Track (const SimTrack &t)
 
void addGenParticle (const reco::GenParticleRef &ref)
 
void addHitEnergy (float energy)
 add rechit energy More...
 
void addRecHitAndFraction (uint32_t hit, float fraction)
 add rechit with fraction More...
 
void addSimHit (const PCaloHit &hit)
 add simhit's energy to cluster More...
 
math::XYZVectorF boostToCM () const
 Vector to boost to the particle centre of mass frame. More...
 
float charge () const
 Electric charge. Note this is taken from the first SimTrack only. More...
 
void clearHitsAndFractions ()
 clear the hits and fractions list More...
 
void clearHitsEnergy ()
 clear the energies list More...
 
float energy () const
 Energy. Note this is taken from the first SimTrack only. More...
 
float et () const
 Transverse energy. Note this is taken from the first SimTrack only. More...
 
float eta () const
 Momentum pseudorapidity. Note this is taken from the simtrack before the calorimeter. More...
 
EncodedEventId eventId () const
 Signal source, crossing number. More...
 
g4t_iterator g4Track_begin () const
 
g4t_iterator g4Track_end () const
 
const std::vector< SimTrack > & g4Tracks () const
 
genp_iterator genParticle_begin () const
 iterators More...
 
genp_iterator genParticle_end () const
 
const reco::GenParticleRefVectorgenParticles () const
 
std::vector< std::pair< uint32_t, float > > hits_and_energies () const
 Returns list of rechit IDs and energies for this SimCluster. More...
 
std::vector< std::pair< uint32_t, float > > hits_and_fractions () const
 Returns list of rechit IDs and fractions for this SimCluster. More...
 
bool longLived () const
 is long lived? More...
 
float mass () const
 Mass. Note this is taken from the first SimTrack only. More...
 
float massSqr () const
 Mass squared. Note this is taken from the first SimTrack only. More...
 
math::XYZVectorF momentum () const
 spatial momentum vector More...
 
float mt () const
 Transverse mass. Note this is taken from the first SimTrack only. More...
 
float mtSqr () const
 Transverse mass squared. Note this is taken from the first SimTrack only. More...
 
int numberOfRecHits () const
 Gives the total number of SimHits, in the cluster. More...
 
int numberOfSimHits () const
 Gives the total number of SimHits, in the cluster. More...
 
float p () const
 Magnitude of momentum vector. Note this is taken from the first SimTrack only. More...
 
const math::XYZTLorentzVectorFp4 () const
 Four-momentum Lorentz vector. Note this is taken from the first SimTrack only. More...
 
uint64_t particleId () const
 
int pdgId () const
 PDG ID. More...
 
float phi () const
 Momentum azimuthal angle. Note this is taken from the first SimTrack only. More...
 
float pt () const
 Transverse momentum. Note this is taken from the first SimTrack only. More...
 
float px () const
 x coordinate of momentum vector. Note this is taken from the first SimTrack only. More...
 
float py () const
 y coordinate of momentum vector. Note this is taken from the first SimTrack only. More...
 
float pz () const
 z coordinate of momentum vector. Note this is taken from the first SimTrack only. More...
 
float rapidity () const
 Rapidity. Note this is taken from the simtrack before the calorimeter. More...
 
 SimCluster ()
 
 SimCluster (const SimTrack &simtrk)
 
 SimCluster (EncodedEventId eventID, uint32_t particleID)
 
float simEnergy () const
 returns the accumulated sim energy in the cluster More...
 
int status () const
 Status word. More...
 
float theta () const
 Momentum polar angle. Note this is taken from the first SimTrack only. More...
 
int threeCharge () const
 Gives charge in unit of quark charge (should be 3 times "charge()") More...
 
float y () const
 Same as rapidity(). More...
 
 ~SimCluster ()
 

Static Public Attributes

static const unsigned int longLivedTag = 65536
 long lived flag More...
 

Private Attributes

std::vector< float > energies_
 
EncodedEventId event_
 
std::vector< float > fractions_
 
std::vector< SimTrackg4Tracks_
 references to G4 and reco::GenParticle tracks More...
 
reco::GenParticleRefVector genParticles_
 
std::vector< uint32_t > hits_
 
uint64_t nsimhits_ {0}
 
uint32_t particleId_ {0}
 
float simhit_energy_ {0.f}
 
math::XYZTLorentzVectorF theMomentum_
 

Friends

std::ostream & operator<< (std::ostream &s, SimCluster const &tp)
 

Detailed Description

Monte Carlo truth information used for tracking validation.

Object with references to the original SimTrack and parent and daughter TrackingVertices. Simulation with high (~100) pileup was taking too much memory so the class was slimmed down and copies of the SimHits were removed.

Author
original author unknown, re-engineering and slimming by Subir Sarkar (subir.nosp@m..sar.nosp@m.kar@c.nosp@m.ern..nosp@m.ch), some tweaking and documentation by Mark Grimes (mark..nosp@m.grim.nosp@m.es@br.nosp@m.isto.nosp@m.l.ac..nosp@m.uk).
Date
original date unknown, re-engineering Jan-May 2013

Definition at line 29 of file SimCluster.h.

Member Typedef Documentation

◆ Charge

typedef int SimCluster::Charge

electric charge type

Definition at line 33 of file SimCluster.h.

◆ g4t_iterator

typedef std::vector<SimTrack>::const_iterator SimCluster::g4t_iterator

Definition at line 41 of file SimCluster.h.

◆ genp_iterator

reference to reco::GenParticle

Definition at line 40 of file SimCluster.h.

◆ LorentzVector

Lorentz vector.

Definition at line 34 of file SimCluster.h.

◆ Point

point in the space

Definition at line 36 of file SimCluster.h.

◆ PolarLorentzVector

Lorentz vector.

Definition at line 35 of file SimCluster.h.

◆ Vector

point in the space

Definition at line 37 of file SimCluster.h.

Constructor & Destructor Documentation

◆ SimCluster() [1/3]

SimCluster::SimCluster ( )

Definition at line 11 of file SimCluster.cc.

11  {
12  // No operation
13 }

◆ SimCluster() [2/3]

SimCluster::SimCluster ( const SimTrack simtrk)

Definition at line 15 of file SimCluster.cc.

References addG4Track(), event_, CoreSimTrack::eventId(), CoreSimTrack::momentum(), particleId_, theMomentum_, and CoreSimTrack::trackId().

15  {
16  addG4Track(simtrk);
17  event_ = simtrk.eventId();
18  particleId_ = simtrk.trackId();
19 
20  theMomentum_.SetPxPyPzE(
21  simtrk.momentum().px(), simtrk.momentum().py(), simtrk.momentum().pz(), simtrk.momentum().E());
22 }
void addG4Track(const SimTrack &t)
Definition: SimCluster.h:72
EncodedEventId event_
Definition: SimCluster.h:223
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
EncodedEventId eventId() const
Definition: CoreSimTrack.h:28
unsigned int trackId() const
Definition: CoreSimTrack.h:31
uint32_t particleId_
Definition: SimCluster.h:225
math::XYZTLorentzVectorF theMomentum_
Definition: SimCluster.h:231

◆ SimCluster() [3/3]

SimCluster::SimCluster ( EncodedEventId  eventID,
uint32_t  particleID 
)

Definition at line 24 of file SimCluster.cc.

References event_, EgammaObjectsElectrons_cfi::particleID, and particleId_.

24  {
25  event_ = eventID;
27 }
EncodedEventId event_
Definition: SimCluster.h:223
uint32_t particleId_
Definition: SimCluster.h:225

◆ ~SimCluster()

SimCluster::~SimCluster ( )

Definition at line 29 of file SimCluster.cc.

29 {}

Member Function Documentation

◆ addG4Track()

void SimCluster::addG4Track ( const SimTrack t)
inline

Definition at line 72 of file SimCluster.h.

References g4Tracks_, and submitPVValidationJobs::t.

Referenced by SimCluster().

72 { g4Tracks_.push_back(t); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: SimCluster.h:234

◆ addGenParticle()

void SimCluster::addGenParticle ( const reco::GenParticleRef ref)
inline

Definition at line 71 of file SimCluster.h.

References genParticles_, and edm::RefVector< C, T, F >::push_back().

71 { genParticles_.push_back(ref); }
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
reco::GenParticleRefVector genParticles_
Definition: SimCluster.h:235

◆ addHitEnergy()

void SimCluster::addHitEnergy ( float  energy)
inline

add rechit energy

Definition at line 181 of file SimCluster.h.

References energies_, and energy().

181 { energies_.emplace_back(energy); }
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
std::vector< float > energies_
Definition: SimCluster.h:229

◆ addRecHitAndFraction()

void SimCluster::addRecHitAndFraction ( uint32_t  hit,
float  fraction 
)
inline

add rechit with fraction

Definition at line 175 of file SimCluster.h.

References HLT_2022v12_cff::fraction, fractions_, and hits_.

175  {
176  hits_.emplace_back(hit);
177  fractions_.emplace_back(fraction);
178  }
std::vector< uint32_t > hits_
Definition: SimCluster.h:227
std::vector< float > fractions_
Definition: SimCluster.h:228

◆ addSimHit()

void SimCluster::addSimHit ( const PCaloHit hit)
inline

add simhit's energy to cluster

Definition at line 216 of file SimCluster.h.

References nsimhits_, and simhit_energy_.

216  {
217  simhit_energy_ += hit.energy();
218  ++nsimhits_;
219  }
float simhit_energy_
Definition: SimCluster.h:226
uint64_t nsimhits_
Definition: SimCluster.h:222

◆ boostToCM()

math::XYZVectorF SimCluster::boostToCM ( ) const
inline

Vector to boost to the particle centre of mass frame.

Definition at line 97 of file SimCluster.h.

References p4().

97 { return p4().BoostToCM(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ charge()

float SimCluster::charge ( void  ) const
inline

Electric charge. Note this is taken from the first SimTrack only.

Definition at line 85 of file SimCluster.h.

References g4Tracks_.

Referenced by threeCharge().

85 { return g4Tracks_[0].charge(); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: SimCluster.h:234

◆ clearHitsAndFractions()

void SimCluster::clearHitsAndFractions ( )
inline

clear the hits and fractions list

Definition at line 204 of file SimCluster.h.

References fractions_, hits_, and edm::swap().

204  {
205  std::vector<uint32_t>().swap(hits_);
206  std::vector<float>().swap(fractions_);
207  }
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
std::vector< uint32_t > hits_
Definition: SimCluster.h:227
std::vector< float > fractions_
Definition: SimCluster.h:228

◆ clearHitsEnergy()

void SimCluster::clearHitsEnergy ( )
inline

clear the energies list

Definition at line 210 of file SimCluster.h.

References energies_, and edm::swap().

210 { std::vector<float>().swap(energies_); }
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
std::vector< float > energies_
Definition: SimCluster.h:229

◆ energy()

float SimCluster::energy ( void  ) const
inline

Energy. Note this is taken from the first SimTrack only.

Definition at line 104 of file SimCluster.h.

References p4().

Referenced by addHitEnergy(), HGCalShowerSeparation::analyze(), HGCalHitCalibration::analyze(), and Jet.Jet::rawEnergy().

104 { return p4().E(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ et()

float SimCluster::et ( ) const
inline

Transverse energy. Note this is taken from the first SimTrack only.

Definition at line 107 of file SimCluster.h.

References p4().

107 { return p4().Et(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ eta()

float SimCluster::eta ( void  ) const
inline

Momentum pseudorapidity. Note this is taken from the simtrack before the calorimeter.

Definition at line 148 of file SimCluster.h.

References p4().

Referenced by Particle.Particle::__str__(), Jet.Jet::jetID(), and Jet.Jet::puJetId().

148 { return p4().Eta(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ eventId()

EncodedEventId SimCluster::eventId ( ) const
inline

Signal source, crossing number.

Note this is taken from the first SimTrack only, but there shouldn't be any SimTracks from different crossings in the SimCluster.

Definition at line 66 of file SimCluster.h.

References event_.

Referenced by ntupleDataFormat.Event::eventIdStr().

66 { return event_; }
EncodedEventId event_
Definition: SimCluster.h:223

◆ g4Track_begin()

g4t_iterator SimCluster::g4Track_begin ( ) const
inline

Definition at line 76 of file SimCluster.h.

References g4Tracks_.

76 { return g4Tracks_.begin(); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: SimCluster.h:234

◆ g4Track_end()

g4t_iterator SimCluster::g4Track_end ( ) const
inline

Definition at line 77 of file SimCluster.h.

References g4Tracks_.

77 { return g4Tracks_.end(); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: SimCluster.h:234

◆ g4Tracks()

const std::vector<SimTrack>& SimCluster::g4Tracks ( ) const
inline

Definition at line 82 of file SimCluster.h.

References g4Tracks_.

82 { return g4Tracks_; }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: SimCluster.h:234

◆ genParticle_begin()

genp_iterator SimCluster::genParticle_begin ( ) const
inline

iterators

Definition at line 74 of file SimCluster.h.

References edm::RefVector< C, T, F >::begin(), and genParticles_.

74 { return genParticles_.begin(); }
reco::GenParticleRefVector genParticles_
Definition: SimCluster.h:235
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223

◆ genParticle_end()

genp_iterator SimCluster::genParticle_end ( ) const
inline

Definition at line 75 of file SimCluster.h.

References edm::RefVector< C, T, F >::end(), and genParticles_.

75 { return genParticles_.end(); }
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
reco::GenParticleRefVector genParticles_
Definition: SimCluster.h:235

◆ genParticles()

const reco::GenParticleRefVector& SimCluster::genParticles ( ) const
inline

Definition at line 80 of file SimCluster.h.

References genParticles_.

80 { return genParticles_; }
reco::GenParticleRefVector genParticles_
Definition: SimCluster.h:235

◆ hits_and_energies()

std::vector<std::pair<uint32_t, float> > SimCluster::hits_and_energies ( ) const
inline

Returns list of rechit IDs and energies for this SimCluster.

Definition at line 193 of file SimCluster.h.

References cms::cuda::assert(), energies_, hits_, mps_fire::i, and mps_fire::result.

193  {
194  assert(hits_.size() == energies_.size());
195  std::vector<std::pair<uint32_t, float>> result;
196  result.reserve(hits_.size());
197  for (size_t i = 0; i < hits_.size(); ++i) {
198  result.emplace_back(hits_[i], energies_[i]);
199  }
200  return result;
201  }
assert(be >=bs)
std::vector< uint32_t > hits_
Definition: SimCluster.h:227
std::vector< float > energies_
Definition: SimCluster.h:229

◆ hits_and_fractions()

std::vector<std::pair<uint32_t, float> > SimCluster::hits_and_fractions ( ) const
inline

Returns list of rechit IDs and fractions for this SimCluster.

Definition at line 184 of file SimCluster.h.

References fractions_, hits_, mps_fire::i, and mps_fire::result.

Referenced by HGCalShowerSeparation::analyze(), HGCalHitCalibration::analyze(), MultiClusterAssociatorByEnergyScoreImpl::makeConnections(), and LCToCPAssociatorByEnergyScoreImpl::makeConnections().

184  {
185  std::vector<std::pair<uint32_t, float>> result;
186  for (size_t i = 0; i < hits_.size(); ++i) {
187  result.emplace_back(hits_[i], fractions_[i]);
188  }
189  return result;
190  }
std::vector< uint32_t > hits_
Definition: SimCluster.h:227
std::vector< float > fractions_
Definition: SimCluster.h:228

◆ longLived()

bool SimCluster::longLived ( ) const
inline

is long lived?

Definition at line 166 of file SimCluster.h.

References longLivedTag, and status().

166 { return status() & longLivedTag; }
int status() const
Status word.
Definition: SimCluster.h:161
static const unsigned int longLivedTag
long lived flag
Definition: SimCluster.h:163

◆ mass()

float SimCluster::mass ( ) const
inline

Mass. Note this is taken from the first SimTrack only.

Definition at line 110 of file SimCluster.h.

References p4().

Referenced by Particle.Particle::__str__(), DiObject.DiMuon::__str__(), and massSqr().

110 { return p4().M(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ massSqr()

float SimCluster::massSqr ( ) const
inline

Mass squared. Note this is taken from the first SimTrack only.

Definition at line 113 of file SimCluster.h.

References mass(), and funct::pow().

113 { return pow(mass(), 2); }
float mass() const
Mass. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:110
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ momentum()

math::XYZVectorF SimCluster::momentum ( ) const
inline

spatial momentum vector

Definition at line 94 of file SimCluster.h.

References p4().

94 { return p4().Vect(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ mt()

float SimCluster::mt ( ) const
inline

Transverse mass. Note this is taken from the first SimTrack only.

Definition at line 116 of file SimCluster.h.

References p4().

116 { return p4().Mt(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ mtSqr()

float SimCluster::mtSqr ( ) const
inline

Transverse mass squared. Note this is taken from the first SimTrack only.

Definition at line 120 of file SimCluster.h.

References p4().

120 { return p4().Mt2(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ numberOfRecHits()

int SimCluster::numberOfRecHits ( void  ) const
inline

Gives the total number of SimHits, in the cluster.

Definition at line 172 of file SimCluster.h.

References hits_.

172 { return hits_.size(); }
std::vector< uint32_t > hits_
Definition: SimCluster.h:227

◆ numberOfSimHits()

int SimCluster::numberOfSimHits ( ) const
inline

Gives the total number of SimHits, in the cluster.

Definition at line 169 of file SimCluster.h.

References nsimhits_.

169 { return nsimhits_; }
uint64_t nsimhits_
Definition: SimCluster.h:222

◆ p()

float SimCluster::p ( ) const
inline

Magnitude of momentum vector. Note this is taken from the first SimTrack only.

Definition at line 101 of file SimCluster.h.

References p4().

Referenced by Electron.Electron::ptErr().

101 { return p4().P(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ p4()

const math::XYZTLorentzVectorF& SimCluster::p4 ( ) const
inline

Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

Definition at line 91 of file SimCluster.h.

References theMomentum_.

Referenced by boostToCM(), Tau.Tau::dxy_approx(), Tau.Tau::dz(), energy(), et(), eta(), mass(), momentum(), mt(), mtSqr(), p(), Lepton.Lepton::p4WithFSR(), phi(), pt(), px(), py(), pz(), rapidity(), and theta().

91 { return theMomentum_; }
math::XYZTLorentzVectorF theMomentum_
Definition: SimCluster.h:231

◆ particleId()

uint64_t SimCluster::particleId ( ) const
inline

Definition at line 68 of file SimCluster.h.

References particleId_.

68 { return particleId_; }
uint32_t particleId_
Definition: SimCluster.h:225

◆ pdgId()

int SimCluster::pdgId ( ) const
inline

PDG ID.

Returns the PDG ID of the first associated gen particle. If there are no gen particles associated then it returns type() from the first SimTrack.

Definition at line 55 of file SimCluster.h.

References edm::RefVector< C, T, F >::begin(), edm::RefVector< C, T, F >::empty(), g4Tracks_, and genParticles_.

Referenced by Particle.Particle::__str__().

55  {
56  if (genParticles_.empty())
57  return g4Tracks_[0].type();
58  else
59  return (*genParticles_.begin())->pdgId();
60  }
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: SimCluster.h:234
int pdgId() const
PDG ID.
Definition: SimCluster.h:55
reco::GenParticleRefVector genParticles_
Definition: SimCluster.h:235
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223

◆ phi()

float SimCluster::phi ( void  ) const
inline

Momentum azimuthal angle. Note this is taken from the first SimTrack only.

Definition at line 140 of file SimCluster.h.

References p4().

Referenced by Particle.Particle::__str__(), and ntupleDataFormat.Track::phiPull().

140 { return p4().Phi(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ pt()

float SimCluster::pt ( ) const
inline

◆ px()

float SimCluster::px ( ) const
inline

x coordinate of momentum vector. Note this is taken from the first SimTrack only.

Definition at line 124 of file SimCluster.h.

References p4().

124 { return p4().Px(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ py()

float SimCluster::py ( ) const
inline

y coordinate of momentum vector. Note this is taken from the first SimTrack only.

Definition at line 128 of file SimCluster.h.

References p4().

128 { return p4().Py(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ pz()

float SimCluster::pz ( ) const
inline

z coordinate of momentum vector. Note this is taken from the first SimTrack only.

Definition at line 132 of file SimCluster.h.

References p4().

132 { return p4().Pz(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ rapidity()

float SimCluster::rapidity ( ) const
inline

Rapidity. Note this is taken from the simtrack before the calorimeter.

Definition at line 152 of file SimCluster.h.

References p4().

Referenced by y().

152 { return p4().Rapidity(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ simEnergy()

float SimCluster::simEnergy ( ) const
inline

returns the accumulated sim energy in the cluster

Definition at line 213 of file SimCluster.h.

References simhit_energy_.

Referenced by HGCalShowerSeparation::analyze(), and HGCalHitCalibration::analyze().

213 { return simhit_energy_; }
float simhit_energy_
Definition: SimCluster.h:226

◆ status()

int SimCluster::status ( void  ) const
inline

Status word.

Returns status() from the first gen particle, or -99 if there are no gen particles attached.

Definition at line 161 of file SimCluster.h.

References edm::RefVector< C, T, F >::empty(), genParticles_, and status().

Referenced by longLived(), and status().

161 { return genParticles_.empty() ? -99 : (*genParticles_[0]).status(); }
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
int status() const
Status word.
Definition: SimCluster.h:161
reco::GenParticleRefVector genParticles_
Definition: SimCluster.h:235

◆ theta()

float SimCluster::theta ( void  ) const
inline

Momentum polar angle. Note this is taken from the first SimTrack only.

Definition at line 144 of file SimCluster.h.

References p4().

Referenced by Tau.Tau::zImpact().

144 { return p4().Theta(); }
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:91

◆ threeCharge()

int SimCluster::threeCharge ( ) const
inline

Gives charge in unit of quark charge (should be 3 times "charge()")

Definition at line 87 of file SimCluster.h.

References charge(), and f.

87 { return lrintf(3.f * charge()); }
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:85
double f[11][100]

◆ y()

float SimCluster::y ( ) const
inline

Same as rapidity().

Definition at line 155 of file SimCluster.h.

References rapidity().

Referenced by svgfig.Ellipse::__repr__(), geometryXMLparser.Alignable::pos(), ntupleDataFormat._HitObject::r(), and ntupleDataFormat._HitObject::r3D().

155 { return rapidity(); }
float rapidity() const
Rapidity. Note this is taken from the simtrack before the calorimeter.
Definition: SimCluster.h:152

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  s,
SimCluster const &  tp 
)
friend

Definition at line 31 of file SimCluster.cc.

31  {
32  s << "CP momentum, q, ID, & Event #: " << tp.p4() << " " << tp.charge() << " " << tp.pdgId() << " "
33  << tp.eventId().bunchCrossing() << "." << tp.eventId().event() << std::endl;
34 
35  for (SimCluster::genp_iterator hepT = tp.genParticle_begin(); hepT != tp.genParticle_end(); ++hepT) {
36  s << " HepMC Track Momentum " << (*hepT)->momentum().rho() << std::endl;
37  }
38 
39  for (SimCluster::g4t_iterator g4T = tp.g4Track_begin(); g4T != tp.g4Track_end(); ++g4T) {
40  s << " Geant Track Momentum " << g4T->momentum() << std::endl;
41  s << " Geant Track ID & type " << g4T->trackId() << " " << g4T->type() << std::endl;
42  if (g4T->type() != tp.pdgId()) {
43  s << " Mismatch b/t SimCluster and Geant types" << std::endl;
44  }
45  }
46  s << " # of cells = " << tp.hits_.size()
47  << ", effective cells = " << std::accumulate(tp.fractions_.begin(), tp.fractions_.end(), 0.f) << std::endl;
48  return s;
49 }
std::vector< SimTrack >::const_iterator g4t_iterator
Definition: SimCluster.h:41

Member Data Documentation

◆ energies_

std::vector<float> SimCluster::energies_
private

Definition at line 229 of file SimCluster.h.

Referenced by addHitEnergy(), clearHitsEnergy(), and hits_and_energies().

◆ event_

EncodedEventId SimCluster::event_
private

Definition at line 223 of file SimCluster.h.

Referenced by eventId(), and SimCluster().

◆ fractions_

std::vector<float> SimCluster::fractions_
private

Definition at line 228 of file SimCluster.h.

Referenced by addRecHitAndFraction(), clearHitsAndFractions(), and hits_and_fractions().

◆ g4Tracks_

std::vector<SimTrack> SimCluster::g4Tracks_
private

references to G4 and reco::GenParticle tracks

Definition at line 234 of file SimCluster.h.

Referenced by addG4Track(), charge(), g4Track_begin(), g4Track_end(), g4Tracks(), and pdgId().

◆ genParticles_

reco::GenParticleRefVector SimCluster::genParticles_
private

◆ hits_

std::vector<uint32_t> SimCluster::hits_
private

◆ longLivedTag

const unsigned int SimCluster::longLivedTag = 65536
static

long lived flag

Definition at line 163 of file SimCluster.h.

Referenced by longLived().

◆ nsimhits_

uint64_t SimCluster::nsimhits_ {0}
private

Definition at line 222 of file SimCluster.h.

Referenced by addSimHit(), and numberOfSimHits().

◆ particleId_

uint32_t SimCluster::particleId_ {0}
private

Definition at line 225 of file SimCluster.h.

Referenced by particleId(), and SimCluster().

◆ simhit_energy_

float SimCluster::simhit_energy_ {0.f}
private

Definition at line 226 of file SimCluster.h.

Referenced by addSimHit(), and simEnergy().

◆ theMomentum_

math::XYZTLorentzVectorF SimCluster::theMomentum_
private

Definition at line 231 of file SimCluster.h.

Referenced by p4(), and SimCluster().