Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:40

0001 #ifndef SimDataFormats_CaloParticle_h
0002 #define SimDataFormats_CaloParticle_h
0003 
0004 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0005 #include "DataFormats/Math/interface/LorentzVector.h"
0006 #include "DataFormats/Math/interface/Point3D.h"
0007 #include "DataFormats/Math/interface/Vector3D.h"
0008 #include "SimDataFormats/CaloAnalysis/interface/SimClusterFwd.h"
0009 #include "SimDataFormats/CaloHit/interface/PCaloHit.h"
0010 #include "SimDataFormats/Track/interface/SimTrack.h"
0011 #include <vector>
0012 
0013 class SimTrack;
0014 class EncodedEventId;
0015 
0016 class CaloParticle {
0017   friend std::ostream &operator<<(std::ostream &s, CaloParticle const &tp);
0018 
0019 public:
0020   typedef int Charge;                                       ///< electric charge type
0021   typedef math::XYZTLorentzVectorD LorentzVector;           ///< Lorentz vector
0022   typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;  ///< Lorentz vector
0023   typedef math::XYZPointD Point;                            ///< point in the space
0024   typedef math::XYZVectorD Vector;                          ///< point in the space
0025 
0026   /// reference to reco::GenParticle
0027   typedef reco::GenParticleRefVector::iterator genp_iterator;
0028   typedef std::vector<SimTrack>::const_iterator g4t_iterator;
0029   typedef SimClusterRefVector::iterator sc_iterator;
0030 
0031   CaloParticle();
0032 
0033   CaloParticle(const SimTrack &simtrk);
0034   CaloParticle(EncodedEventId eventID, uint32_t particleID);  // for PU
0035 
0036   // destructor
0037   ~CaloParticle();
0038 
0039   /** @brief PDG ID.
0040    *
0041    * Returns the PDG ID of the first associated gen particle. If there are no
0042    * gen particles associated then it returns type() from the first SimTrack. */
0043   int pdgId() const {
0044     if (genParticles_.empty())
0045       return g4Tracks_[0].type();
0046     else
0047       return (*genParticles_.begin())->pdgId();
0048   }
0049 
0050   /** @brief Signal source, crossing number.
0051    *
0052    * Note this is taken from the first SimTrack only, but there shouldn't be any
0053    * SimTracks from different crossings in the CaloParticle. */
0054   EncodedEventId eventId() const { return event_; }
0055 
0056   uint64_t particleId() const { return particleId_; }
0057 
0058   // Setters for G4 and reco::GenParticle
0059   void addGenParticle(const reco::GenParticleRef &ref) { genParticles_.push_back(ref); }
0060   void addSimCluster(const SimClusterRef &ref) { simClusters_.push_back(ref); }
0061   void addG4Track(const SimTrack &t) { g4Tracks_.push_back(t); }
0062   /// iterators
0063   genp_iterator genParticle_begin() const { return genParticles_.begin(); }
0064   genp_iterator genParticle_end() const { return genParticles_.end(); }
0065   g4t_iterator g4Track_begin() const { return g4Tracks_.begin(); }
0066   g4t_iterator g4Track_end() const { return g4Tracks_.end(); }
0067   sc_iterator simCluster_begin() const { return simClusters_.begin(); }
0068   sc_iterator simCluster_end() const { return simClusters_.end(); }
0069 
0070   // Getters for Embd and Sim Tracks
0071   const reco::GenParticleRefVector &genParticles() const { return genParticles_; }
0072   const SimClusterRefVector &simClusters() const { return simClusters_; }
0073   // Only for clusters from the signal vertex
0074   const std::vector<SimTrack> &g4Tracks() const { return g4Tracks_; }
0075 
0076   void clearSimClusters() { simClusters_.clear(); }
0077 
0078   /// @brief Electric charge. Note this is taken from the first SimTrack only.
0079   float charge() const { return g4Tracks_[0].charge(); }
0080   /// Gives charge in unit of quark charge (should be 3 times "charge()")
0081   int threeCharge() const { return lrintf(3.f * charge()); }
0082 
0083   /// @brief Four-momentum Lorentz vector. Note this is taken from the first
0084   /// SimTrack only.
0085   const math::XYZTLorentzVectorF &p4() const { return theMomentum_; }
0086 
0087   /// @brief spatial momentum vector
0088   math::XYZVectorF momentum() const { return p4().Vect(); }
0089 
0090   /// @brief Vector to boost to the particle centre of mass frame.
0091   math::XYZVectorF boostToCM() const { return p4().BoostToCM(); }
0092 
0093   /// @brief Magnitude of momentum vector. Note this is taken from the first
0094   /// SimTrack only.
0095   float p() const { return p4().P(); }
0096 
0097   /// @brief Energy. Note this is taken from the first SimTrack only.
0098   float energy() const { return p4().E(); }
0099 
0100   /// @brief Transverse energy. Note this is taken from the first SimTrack only.
0101   float et() const { return p4().Et(); }
0102 
0103   /// @brief Mass. Note this is taken from the first SimTrack only.
0104   float mass() const { return p4().M(); }
0105 
0106   /// @brief Mass squared. Note this is taken from the first SimTrack only.
0107   float massSqr() const { return pow(mass(), 2); }
0108 
0109   /// @brief Transverse mass. Note this is taken from the first SimTrack only.
0110   float mt() const { return p4().Mt(); }
0111 
0112   /// @brief Transverse mass squared. Note this is taken from the first SimTrack
0113   /// only.
0114   float mtSqr() const { return p4().Mt2(); }
0115 
0116   /// @brief x coordinate of momentum vector. Note this is taken from the first
0117   /// SimTrack only.
0118   float px() const { return p4().Px(); }
0119 
0120   /// @brief y coordinate of momentum vector. Note this is taken from the first
0121   /// SimTrack only.
0122   float py() const { return p4().Py(); }
0123 
0124   /// @brief z coordinate of momentum vector. Note this is taken from the first
0125   /// SimTrack only.
0126   float pz() const { return p4().Pz(); }
0127 
0128   /// @brief Transverse momentum. Note this is taken from the first SimTrack
0129   /// only.
0130   float pt() const { return p4().Pt(); }
0131 
0132   /// @brief Momentum azimuthal angle. Note this is taken from the first
0133   /// SimTrack only.
0134   float phi() const { return p4().Phi(); }
0135 
0136   /// @brief Momentum polar angle. Note this is taken from the first SimTrack
0137   /// only.
0138   float theta() const { return p4().Theta(); }
0139 
0140   /// @brief Momentum pseudorapidity. Note this is taken from the simtrack
0141   /// before the calorimeter
0142   float eta() const { return p4().Eta(); }
0143 
0144   /// @brief Rapidity. Note this is taken from the simtrack before the
0145   /// calorimeter
0146   float rapidity() const { return p4().Rapidity(); }
0147 
0148   /// @brief Same as rapidity().
0149   float y() const { return rapidity(); }
0150 
0151   /** @brief Status word.
0152    *
0153    * Returns status() from the first gen particle, or -99 if there are no gen
0154    * particles attached. */
0155   int status() const { return genParticles_.empty() ? -99 : (*genParticles_[0]).status(); }
0156 
0157   static const unsigned int longLivedTag;  ///< long lived flag
0158 
0159   /// is long lived?
0160   bool longLived() const { return status() & longLivedTag; }
0161 
0162   /** @brief Gives the total number of SimHits, in the cluster */
0163   int numberOfSimHits() const { return nsimhits_; }
0164 
0165   /** @brief Gives the total number of SimHits, in the cluster */
0166   int numberOfRecHits() const { return hits_.size(); }
0167 
0168   /** @brief add rechit with fraction */
0169   void addRecHitAndFraction(uint32_t hit, float fraction) {
0170     hits_.emplace_back(hit);
0171     fractions_.emplace_back(fraction);
0172   }
0173 
0174   /** @brief Returns list of rechit IDs and fractions for this CaloParticle */
0175   std::vector<std::pair<uint32_t, float>> hits_and_fractions() const {
0176     std::vector<std::pair<uint32_t, float>> result;
0177     for (size_t i = 0; i < hits_.size(); ++i) {
0178       result.emplace_back(hits_[i], fractions_[i]);
0179     }
0180     return result;
0181   }
0182 
0183   /** @brief returns the accumulated sim energy in the cluster */
0184   float simEnergy() const { return simhit_energy_; }
0185 
0186   /** @brief add simhit's energy to cluster */
0187   void addSimHit(const PCaloHit &hit) {
0188     simhit_energy_ += hit.energy();
0189     ++nsimhits_;
0190   }
0191 
0192 protected:
0193   uint64_t nsimhits_{0};
0194   EncodedEventId event_;
0195 
0196   uint32_t particleId_{0};
0197   float simhit_energy_{0.f};
0198   std::vector<uint32_t> hits_;
0199   std::vector<float> fractions_;
0200 
0201   math::XYZTLorentzVectorF theMomentum_;
0202 
0203   /// references to G4 and reco::GenParticle tracks
0204   std::vector<SimTrack> g4Tracks_;
0205   reco::GenParticleRefVector genParticles_;
0206 
0207   SimClusterRefVector simClusters_;
0208 };
0209 
0210 #endif  // SimDataFormats_CaloParticle_H