Go to the documentation of this file. 1 #ifndef SimDataFormats_SimCluster_h
2 #define SimDataFormats_SimCluster_h
101 float p()
const {
return p4().P(); }
107 float et()
const {
return p4().Et(); }
116 float mt()
const {
return p4().Mt(); }
124 float px()
const {
return p4().Px(); }
128 float py()
const {
return p4().Py(); }
132 float pz()
const {
return p4().Pz(); }
136 float pt()
const {
return p4().Pt(); }
140 float phi()
const {
return p4().Phi(); }
148 float eta()
const {
return p4().Eta(); }
185 std::vector<std::pair<uint32_t, float>>
result;
186 for (
size_t i = 0;
i <
hits_.size(); ++
i) {
195 std::vector<std::pair<uint32_t, float>>
result;
197 for (
size_t i = 0;
i <
hits_.size(); ++
i) {
205 std::vector<uint32_t>().swap(
hits_);
235 #endif // SimDataFormats_SimCluster_H
std::vector< float > fractions_
float px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only.
g4t_iterator g4Track_end() const
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
int threeCharge() const
Gives charge in unit of quark charge (should be 3 times "charge()")
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
genp_iterator genParticle_end() const
const_iterator begin() const
Initialize an iterator over the RefVector.
float energy() const
Energy. Note this is taken from the first SimTrack only.
Monte Carlo truth information used for tracking validation.
const reco::GenParticleRefVector & genParticles() const
genp_iterator genParticle_begin() const
iterators
bool longLived() const
is long lived?
void addSimHit(const PCaloHit &hit)
add simhit's energy to cluster
float rapidity() const
Rapidity. Note this is taken from the simtrack before the calorimeter.
const_iterator end() const
Termination of iteration.
bool empty() const
Is the RefVector empty.
int numberOfSimHits() const
Gives the total number of SimHits, in the cluster.
float mt() const
Transverse mass. Note this is taken from the first SimTrack only.
float p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
const std::vector< SimTrack > & g4Tracks() const
std::vector< std::pair< uint32_t, float > > hits_and_energies() const
Returns list of rechit IDs and energies for this SimCluster.
void addHitEnergy(float energy)
add rechit energy
std::vector< float > energies_
static const unsigned int longLivedTag
long lived flag
reco::GenParticleRefVector genParticles_
math::XYZTLorentzVectorF theMomentum_
EncodedEventId eventId() const
Signal source, crossing number.
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
float phi() const
Momentum azimuthal angle. Note this is taken from the first SimTrack only.
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
float et() const
Transverse energy. Note this is taken from the first SimTrack only.
float pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only.
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
math::XYZPointD Point
point in the space
void addG4Track(const SimTrack &t)
int status() const
Status word.
math::XYZVectorD Vector
point in the space
float eta() const
Momentum pseudorapidity. Note this is taken from the simtrack before the calorimeter.
const math::XYZTLorentzVectorF & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
g4t_iterator g4Track_begin() const
friend std::ostream & operator<<(std::ostream &s, SimCluster const &tp)
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
math::XYZVectorF boostToCM() const
Vector to boost to the particle centre of mass frame.
std::vector< SimTrack >::const_iterator g4t_iterator
float massSqr() const
Mass squared. Note this is taken from the first SimTrack only.
int Charge
electric charge type
math::XYZVectorF momentum() const
spatial momentum vector
void addGenParticle(const reco::GenParticleRef &ref)
uint64_t particleId() const
float mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
float simEnergy() const
returns the accumulated sim energy in the cluster
int numberOfRecHits() const
Gives the total number of SimHits, in the cluster.
void clearHitsEnergy()
clear the energies list
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
unsigned long long uint64_t
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
float y() const
Same as rapidity().
Power< A, B >::type pow(const A &a, const B &b)
float py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only.
void addRecHitAndFraction(uint32_t hit, float fraction)
add rechit with fraction
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
std::vector< uint32_t > hits_
float theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
float mass() const
Mass. Note this is taken from the first SimTrack only.
reco::GenParticleRefVector::iterator genp_iterator
reference to reco::GenParticle
void clearHitsAndFractions()
clear the hits and fractions list
float pt() const
Transverse momentum. Note this is taken from the first SimTrack only.