Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:03:48

0001 #ifndef Candidate_LeafCandidate_h
0002 #define Candidate_LeafCandidate_h
0003 /** \class reco::LeafCandidate
0004  *
0005  * particle candidate with no constituent nor daughters
0006  *
0007  * \author Luca Lista, INFN
0008  *
0009  *
0010  */
0011 #include "DataFormats/Candidate/interface/Candidate.h"
0012 #include "ParticleState.h"
0013 
0014 namespace reco {
0015 
0016   class LeafCandidate : public Candidate {
0017   public:
0018     /// collection of daughter candidates
0019     typedef CandidateCollection daughters;
0020     /// electric charge type
0021     typedef int Charge;
0022     /// Lorentz vector
0023     typedef math::XYZTLorentzVector LorentzVector;
0024     /// Lorentz vector
0025     typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;
0026     /// point in the space
0027     typedef math::XYZPoint Point;
0028     /// point in the space
0029     typedef math::XYZVector Vector;
0030 
0031     typedef unsigned int index;
0032 
0033     LeafCandidate() {}
0034 
0035     // constructor from candidate
0036     explicit LeafCandidate(const Candidate& c) : m_state(c.charge(), c.polarP4(), c.vertex(), c.pdgId(), c.status()) {}
0037 
0038 #if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
0039     template <typename... Args>
0040     explicit LeafCandidate(Args&&... args) : m_state(std::forward<Args>(args)...) {}
0041 
0042     LeafCandidate(LeafCandidate& rh) : m_state(rh.m_state) {}
0043 
0044     LeafCandidate(LeafCandidate&&) = default;
0045     LeafCandidate(LeafCandidate const&) = default;
0046     LeafCandidate& operator=(LeafCandidate&&) = default;
0047     LeafCandidate& operator=(LeafCandidate const&) = default;
0048 #else
0049     // for Reflex to parse...  (compilation will use the above)
0050     LeafCandidate(Charge q,
0051                   const PtEtaPhiMass& p4,
0052                   const Point& vtx = Point(0, 0, 0),
0053                   int pdgId = 0,
0054                   int status = 0,
0055                   bool integerCharge = true);
0056     LeafCandidate(Charge q,
0057                   const LorentzVector& p4,
0058                   const Point& vtx = Point(0, 0, 0),
0059                   int pdgId = 0,
0060                   int status = 0,
0061                   bool integerCharge = true);
0062     LeafCandidate(Charge q,
0063                   const PolarLorentzVector& p4,
0064                   const Point& vtx = Point(0, 0, 0),
0065                   int pdgId = 0,
0066                   int status = 0,
0067                   bool integerCharge = true);
0068     LeafCandidate(Charge q,
0069                   const GlobalVector& p3,
0070                   float iEnergy,
0071                   float imass,
0072                   const Point& vtx = Point(0, 0, 0),
0073                   int pdgId = 0,
0074                   int status = 0,
0075                   bool integerCharge = true);
0076 #endif
0077 
0078     void construct(int qx3, float pt, float eta, float phi, float mass, const Point& vtx, int pdgId, int status) {
0079       m_state = ParticleState(qx3, PolarLorentzVector(pt, eta, phi, mass), vtx, pdgId, status, false);
0080     }
0081 
0082     /// destructor
0083     ~LeafCandidate() override;
0084     /// number of daughters
0085     size_t numberOfDaughters() const override;
0086     /// return daughter at a given position (throws an exception)
0087     const Candidate* daughter(size_type) const override;
0088     /// number of mothers
0089     size_t numberOfMothers() const override;
0090     /// return mother at a given position (throws an exception)
0091     const Candidate* mother(size_type) const override;
0092     /// return daughter at a given position (throws an exception)
0093     Candidate* daughter(size_type) override;
0094     /// return daughter with a specified role name
0095     Candidate* daughter(const std::string& s) override;
0096     /// return daughter with a specified role name
0097     const Candidate* daughter(const std::string& s) const override;
0098     /// return the number of source Candidates
0099     /// ( the candidates used to construct this Candidate)
0100     size_t numberOfSourceCandidatePtrs() const override { return 0; }
0101     /// return a Ptr to one of the source Candidates
0102     /// ( the candidates used to construct this Candidate)
0103     CandidatePtr sourceCandidatePtr(size_type i) const override { return CandidatePtr(); }
0104 
0105     /// electric charge
0106     int charge() const final { return m_state.charge(); }
0107     /// set electric charge
0108     void setCharge(Charge q) final { m_state.setCharge(q); }
0109     /// electric charge
0110     int threeCharge() const final { return m_state.threeCharge(); }
0111     /// set electric charge
0112     void setThreeCharge(Charge qx3) final { m_state.setThreeCharge(qx3); }
0113     /// four-momentum Lorentz vector
0114     const LorentzVector& p4() const final { return m_state.p4(); }
0115     /// four-momentum Lorentz vector
0116     const PolarLorentzVector& polarP4() const final { return m_state.polarP4(); }
0117     /// spatial momentum vector
0118     Vector momentum() const final { return m_state.momentum(); }
0119     /// boost vector to boost a Lorentz vector
0120     /// to the particle center of mass system
0121     Vector boostToCM() const final { return m_state.boostToCM(); }
0122     /// magnitude of momentum vector
0123     double p() const final { return m_state.p(); }
0124     /// energy
0125     double energy() const final { return m_state.energy(); }
0126     /// transverse energy
0127     double et() const final { return m_state.et(); }
0128     /// transverse energy squared (use this for cut!)
0129     double et2() const final { return m_state.et2(); }
0130     /// mass
0131     double mass() const final { return m_state.mass(); }
0132     /// mass squared
0133     double massSqr() const final { return mass() * mass(); }
0134 
0135     /// transverse mass
0136     double mt() const final { return m_state.mt(); }
0137     /// transverse mass squared
0138     double mtSqr() const final { return m_state.mtSqr(); }
0139     /// x coordinate of momentum vector
0140     double px() const final { return m_state.px(); }
0141     /// y coordinate of momentum vector
0142     double py() const final { return m_state.py(); }
0143     /// z coordinate of momentum vector
0144     double pz() const final { return m_state.pz(); }
0145     /// transverse momentum
0146     double pt() const final { return m_state.pt(); }
0147     /// momentum azimuthal angle
0148     double phi() const final { return m_state.phi(); }
0149     /// momentum polar angle
0150     double theta() const final { return m_state.theta(); }
0151     /// momentum pseudorapidity
0152     double eta() const final { return m_state.eta(); }
0153     /// rapidity
0154     double rapidity() const final { return m_state.rapidity(); }
0155     /// rapidity
0156     double y() const final { return rapidity(); }
0157     /// set 4-momentum
0158     void setP4(const LorentzVector& p4) final { m_state.setP4(p4); }
0159     /// set 4-momentum
0160     void setP4(const PolarLorentzVector& p4) final { m_state.setP4(p4); }
0161     /// set particle mass
0162     void setMass(double m) final { m_state.setMass(m); }
0163     void setPz(double pz) final { m_state.setPz(pz); }
0164     /// vertex position                 (overwritten by PF...)
0165     const Point& vertex() const override { return m_state.vertex(); }
0166     /// x coordinate of vertex position
0167     double vx() const override { return m_state.vx(); }
0168     /// y coordinate of vertex position
0169     double vy() const override { return m_state.vy(); }
0170     /// z coordinate of vertex position
0171     double vz() const override { return m_state.vz(); }
0172     /// set vertex
0173     void setVertex(const Point& vertex) override { m_state.setVertex(vertex); }
0174 
0175     /// PDG identifier
0176     int pdgId() const final { return m_state.pdgId(); }
0177     // set PDG identifier
0178     void setPdgId(int pdgId) final { m_state.setPdgId(pdgId); }
0179     /// status word
0180     int status() const final { return m_state.status(); }
0181     /// set status word
0182     void setStatus(int status) final { m_state.setStatus(status); }
0183     /// long lived flag
0184     /// set long lived flag
0185     void setLongLived() final { m_state.setLongLived(); }
0186     /// is long lived?
0187     bool longLived() const final { return m_state.longLived(); }
0188     /// do mass constraint flag
0189     /// set mass constraint flag
0190     void setMassConstraint() final { m_state.setMassConstraint(); }
0191     /// do mass constraint?
0192     bool massConstraint() const final { return m_state.massConstraint(); }
0193 
0194     /// returns a clone of the Candidate object
0195     LeafCandidate* clone() const override { return new LeafCandidate(*this); }
0196 
0197     /// chi-squares
0198     double vertexChi2() const override;
0199     /** Number of degrees of freedom                                                                                   
0200      *  Meant to be Double32_t for soft-assignment fitters:                                                            
0201      *  tracks may contribute to the vertex with fractional weights.                                                   
0202      *  The ndof is then = to the sum of the track weights.                                                            
0203      *  see e.g. CMS NOTE-2006/032, CMS NOTE-2004/002                                                                  
0204      */
0205     double vertexNdof() const override;
0206     /// chi-squared divided by n.d.o.f.
0207     double vertexNormalizedChi2() const override;
0208     /// (i, j)-th element of error matrix, i, j = 0, ... 2
0209     double vertexCovariance(int i, int j) const override;
0210     /// return SMatrix
0211     CovarianceMatrix vertexCovariance() const final {
0212       CovarianceMatrix m;
0213       fillVertexCovariance(m);
0214       return m;
0215     }
0216     /// fill SMatrix
0217     void fillVertexCovariance(CovarianceMatrix& v) const override;
0218     /// returns true if this candidate has a reference to a master clone.
0219     /// This only happens if the concrete Candidate type is ShallowCloneCandidate
0220     bool hasMasterClone() const override;
0221     /// returns ptr to master clone, if existing.
0222     /// Throws an exception unless the concrete Candidate type is ShallowCloneCandidate
0223     const CandidateBaseRef& masterClone() const override;
0224     /// returns true if this candidate has a ptr to a master clone.
0225     /// This only happens if the concrete Candidate type is ShallowClonePtrCandidate
0226     bool hasMasterClonePtr() const override;
0227     /// returns ptr to master clone, if existing.
0228     /// Throws an exception unless the concrete Candidate type is ShallowClonePtrCandidate
0229     const CandidatePtr& masterClonePtr() const override;
0230 
0231     /// cast master clone reference to a concrete type
0232     template <typename Ref>
0233     Ref masterRef() const {
0234       return masterClone().template castTo<Ref>();
0235     }
0236     /// get a component
0237 
0238     template <typename T>
0239     T get() const {
0240       if (hasMasterClone())
0241         return masterClone()->get<T>();
0242       else
0243         return reco::get<T>(*this);
0244     }
0245     /// get a component
0246     template <typename T, typename Tag>
0247     T get() const {
0248       if (hasMasterClone())
0249         return masterClone()->get<T, Tag>();
0250       else
0251         return reco::get<T, Tag>(*this);
0252     }
0253     /// get a component
0254     template <typename T>
0255     T get(size_type i) const {
0256       if (hasMasterClone())
0257         return masterClone()->get<T>(i);
0258       else
0259         return reco::get<T>(*this, i);
0260     }
0261     /// get a component
0262     template <typename T, typename Tag>
0263     T get(size_type i) const {
0264       if (hasMasterClone())
0265         return masterClone()->get<T, Tag>(i);
0266       else
0267         return reco::get<T, Tag>(*this, i);
0268     }
0269     /// number of components
0270     template <typename T>
0271     size_type numberOf() const {
0272       if (hasMasterClone())
0273         return masterClone()->numberOf<T>();
0274       else
0275         return reco::numberOf<T>(*this);
0276     }
0277     /// number of components
0278     template <typename T, typename Tag>
0279     size_type numberOf() const {
0280       if (hasMasterClone())
0281         return masterClone()->numberOf<T, Tag>();
0282       else
0283         return reco::numberOf<T, Tag>(*this);
0284     }
0285 
0286     bool isElectron() const override;
0287     bool isMuon() const override;
0288     bool isStandAloneMuon() const override;
0289     bool isGlobalMuon() const override;
0290     bool isTrackerMuon() const override;
0291     bool isCaloMuon() const override;
0292     bool isPhoton() const override;
0293     bool isConvertedPhoton() const override;
0294     bool isJet() const override;
0295 
0296   private:
0297     ParticleState m_state;
0298 
0299   private:
0300     /// check overlap with another Candidate
0301     bool overlap(const Candidate&) const override;
0302     template <typename, typename, typename>
0303     friend struct component;
0304     friend class ::OverlapChecker;
0305     friend class ShallowCloneCandidate;
0306     friend class ShallowClonePtrCandidate;
0307   };
0308 
0309 }  // namespace reco
0310 
0311 #endif