CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PackedCandidate.h
Go to the documentation of this file.
1 #ifndef __DataFormats_PatCandidates_PackedCandidate_h__
2 #define __DataFormats_PatCandidates_PackedCandidate_h__
3 
12 /* #include "DataFormats/Math/interface/PtEtaPhiMass.h" */
13 
14 namespace pat {
16  public:
27 
28  typedef unsigned int index;
29 
32  : p4_(0,0,0,0), p4c_(0,0,0,0), vertex_(0,0,0), dphi_(0), pdgId_(0), qualityFlags_(0), unpacked_(false),dxydxy_(0),dzdz_(0),dxydz_(0),dlambdadz_(0),dphidxy_(0),packedHits_(0) { }
33  explicit PackedCandidate( const reco::Candidate & c, const reco::VertexRef &pv)
34  : p4_(c.pt(), c.eta(), c.phi(), c.mass()), p4c_(p4_), vertex_(c.vertex()), dphi_(0), pdgId_(c.pdgId()), qualityFlags_(0), pvRef_(pv), unpacked_(true) , unpackedVtx_(true),dxydxy_(0),dzdz_(0),dxydz_(0),dlambdadz_(0),dphidxy_(0),packedHits_(0){ packBoth(); }
35 
36  explicit PackedCandidate( const PolarLorentzVector &p4, const Point &vtx, float phiAtVtx, int pdgId, const reco::VertexRef &pv)
37  : p4_(p4), p4c_(p4_), vertex_(vtx), dphi_(reco::deltaPhi(phiAtVtx,p4_.phi())), pdgId_(pdgId), qualityFlags_(0), pvRef_(pv), unpacked_(true), unpackedVtx_(true),dxydxy_(0),dzdz_(0),dxydz_(0),dlambdadz_(0),dphidxy_(0),packedHits_(0) { packBoth(); }
38 
39  explicit PackedCandidate( const LorentzVector &p4, const Point &vtx, float phiAtVtx, int pdgId, const reco::VertexRef &pv)
40  : p4_(p4.Pt(), p4.Eta(), p4.Phi(), p4.M()), p4c_(p4), vertex_(vtx), dphi_(reco::deltaPhi(phiAtVtx,p4_.phi())), pdgId_(pdgId), qualityFlags_(0), pvRef_(pv), unpacked_(true), unpackedVtx_(true),dxydxy_(0),dzdz_(0),dxydz_(0),dlambdadz_(0),dphidxy_(0),packedHits_(0) { packBoth(); }
41 
42 
43 
44 
46  virtual ~PackedCandidate();
48  virtual const_iterator begin() const;
50  virtual const_iterator end() const;
52  virtual iterator begin();
54  virtual iterator end();
56  virtual size_t numberOfDaughters() const;
58  virtual const reco::Candidate * daughter( size_type ) const;
60  virtual size_t numberOfMothers() const;
62  virtual const reco::Candidate * mother( size_type ) const;
64  virtual reco::Candidate * daughter( size_type );
66  virtual reco::Candidate * daughter(const std::string& s );
68  virtual const reco::Candidate * daughter(const std::string& s ) const;
71  virtual size_t numberOfSourceCandidatePtrs() const {return 0;}
75  return reco::CandidatePtr();
76  }
77 
79  virtual int charge() const {
80  switch (abs(pdgId_)) {
81  case 211: return (pdgId_>0)-(pdgId_<0);
82  case 11: return (-1)*(pdgId_>0)-(pdgId_<0); //e
83  case 13: return (-1)*(pdgId_>0)-(pdgId_<0); //mu
84  case 15: return (-1)*(pdgId_>0)-(pdgId_<0); //tau
85  case 24: return (-1)*(pdgId_>0)-(pdgId_<0); //W
86  default: return 0; //FIXME: charge is not defined
87  }
88  }
90  virtual void setCharge( int charge) {}
92  virtual int threeCharge() const {return charge()*3;}
94  virtual void setThreeCharge( int threecharge) {}
96  virtual const LorentzVector & p4() const { if (!unpacked_) unpack(); return p4c_; }
98  virtual const PolarLorentzVector & polarP4() const { if (!unpacked_) unpack(); return p4_; }
100  virtual Vector momentum() const { if (!unpacked_) unpack(); return p4c_.Vect(); }
103  virtual Vector boostToCM() const { if (!unpacked_) unpack(); return p4c_.BoostToCM(); }
105  virtual double p() const { if (!unpacked_) unpack(); return p4c_.P(); }
107  virtual double energy() const { if (!unpacked_) unpack(); return p4c_.E(); }
109  virtual double et() const { if (!unpacked_) unpack(); return p4_.Et(); }
111  virtual float mass() const { if (!unpacked_) unpack(); return p4_.M(); }
113  virtual float massSqr() const { if (!unpacked_) unpack(); return p4_.M()*p4_.M(); }
114 
116  virtual double mt() const { if (!unpacked_) unpack(); return p4_.Mt(); }
118  virtual double mtSqr() const { if (!unpacked_) unpack(); return p4_.Mt2(); }
120  virtual double px() const { if (!unpacked_) unpack(); return p4c_.Px(); }
122  virtual double py() const { if (!unpacked_) unpack(); return p4c_.Py(); }
124  virtual double pz() const { if (!unpacked_) unpack(); return p4c_.Pz(); }
126  virtual float pt() const { if (!unpacked_) unpack(); return p4_.Pt();}
128  virtual float phi() const { if (!unpacked_) unpack(); return p4_.Phi(); }
130  virtual float phiAtVtx() const {
131  maybeUnpackBoth();
132  float ret = p4_.Phi() + dphi_;
133  while (ret > float(M_PI)) ret -= 2*float(M_PI);
134  while (ret < -float(M_PI)) ret += 2*float(M_PI);
135  return ret;
136  }
138  virtual double theta() const { if (!unpacked_) unpack(); return p4_.Theta(); }
140  virtual float eta() const { if (!unpacked_) unpack(); return p4_.Eta(); }
142  virtual double rapidity() const { if (!unpacked_) unpack(); return p4_.Rapidity(); }
144  virtual double y() const { if (!unpacked_) unpack(); return p4_.Rapidity(); }
146  virtual void setP4( const LorentzVector & p4 ) {
147  maybeUnpackBoth(); // changing px,py,pz changes also mapping between dxy,dz and x,y,z
148  p4_ = PolarLorentzVector(p4.Pt(), p4.Eta(), p4.Phi(), p4.M());
149  packBoth();
150  }
152  virtual void setP4( const PolarLorentzVector & p4 ) {
153  maybeUnpackBoth(); // changing px,py,pz changes also mapping between dxy,dz and x,y,z
154  p4_ = p4;
155  packBoth();
156  }
158  virtual void setMass( double m ) {
159  if (!unpacked_) unpack();
160  p4_ = PolarLorentzVector(p4_.Pt(), p4_.Eta(), p4_.Phi(), m);
161  pack();
162  }
163  virtual void setPz( double pz ) {
164  maybeUnpackBoth(); // changing px,py,pz changes also mapping between dxy,dz and x,y,z
165  p4c_ = LorentzVector(p4c_.Px(), p4c_.Py(), pz, p4c_.E());
166  p4_ = PolarLorentzVector(p4c_.Pt(), p4c_.Eta(), p4c_.Phi(), p4c_.M());
167  packBoth();
168  }
170 
171  virtual void setTrackProperties( const reco::Track & tk, const reco::Track::CovarianceMatrix & covariance) {
172  dxydxy_ = covariance(3,3);
173  dxydz_ = covariance(3,4);
174  dzdz_ = covariance(4,4);
175  dphidxy_ = covariance(2,3);
176  dlambdadz_ = covariance(1,4);
177  dptdpt_ = covariance(0,0)*pt()*pt();
178  detadeta_ = covariance(1,1);
179  dphidphi_ = covariance(2,2)*pt()*pt();
180 
182  int numberOfPixelHits_ = tk.hitPattern().numberOfValidPixelHits();
183  if (numberOfPixelHits_ > 7) numberOfPixelHits_ = 7;
184  int numberOfStripHits_ = tk.hitPattern().numberOfValidHits() - numberOfPixelHits_;
185  if (numberOfStripHits_ > 31) numberOfStripHits_ = 31;
186  packedHits_ = (numberOfPixelHits_&0x7) | (numberOfStripHits_ << 3);
187  packBoth();
188  }
189 
190  virtual void setTrackProperties( const reco::Track & tk ) {
192  }
193 
194  int numberOfPixelHits() const { return packedHits_ & 0x7; }
195  int numberOfHits() const { return (packedHits_ >> 3) + numberOfPixelHits(); }
196 
198  virtual const Point & vertex() const { maybeUnpackBoth(); return vertex_; }//{ if (fromPV_) return Point(0,0,0); else return Point(0,0,100); }
200  virtual double vx() const { maybeUnpackBoth(); return vertex_.X(); }//{ return 0; }
202  virtual double vy() const { maybeUnpackBoth(); return vertex_.Y(); }//{ return 0; }
204  virtual double vz() const { maybeUnpackBoth(); return vertex_.Z(); }//{ if (fromPV_) return 0; else return 100; }
206  virtual void setVertex( const Point & vertex ) { maybeUnpackBoth(); vertex_ = vertex; packVtx(); }
207 
208  enum PVAssoc { NoPV=0, PVLoose=1, PVTight=2, PVUsedInFit=3 } ;
209  const PVAssoc fromPV() const { return PVAssoc((qualityFlags_ & fromPVMask)>>fromPVShift); }
211 
214  const reco::VertexRef vertexRef() const { return pvRef_; }
215 
217  virtual float dxy() const { maybeUnpackBoth(); return dxy_; }
219  virtual float dz() const { maybeUnpackBoth(); return dz_; }
221  virtual float dxy(const Point &p) const ;
223  virtual float dz(const Point &p) const ;
224 
226  virtual float dzError() const { maybeUnpackBoth(); return sqrt(dzdz_); }
228  virtual float dxyError() const { maybeUnpackBoth(); return sqrt(dxydxy_); }
229 
230 
232  virtual reco::Track pseudoTrack() const;
233 
238 
242  noLostInnerHits=0, // it could still not have a hit in the first layer, e.g. if it crosses an inactive sensor
245  };
248  }
250  int lost = hits; if (lost > 2) lost = 2; // protection against misuse
251  lost++; // shift so it's 0 .. 3 instead of (-1) .. 2
253  }
254 
255  void setMuonID(bool isStandAlone, bool isGlobal) {
256  int16_t muonFlags = isStandAlone | (2*isGlobal);
258  }
259 
261  virtual int pdgId() const { return pdgId_; }
262  // set PDG identifier
263  virtual void setPdgId( int pdgId ) { pdgId_ = pdgId; }
265  virtual int status() const { return qualityFlags_; } /*FIXME*/
267  virtual void setStatus( int status ) {} /*FIXME*/
269  static const unsigned int longLivedTag = 0; /*FIXME*/
271  virtual void setLongLived() {} /*FIXME*/
273  virtual bool longLived() const;
275  static const unsigned int massConstraintTag = 0; /*FIXME*/
277  virtual void setMassConstraint() {} /*FIXME*/
279  virtual bool massConstraint() const;
280 
282  virtual PackedCandidate * clone() const {
283  return new PackedCandidate( *this );
284  }
285 
287  virtual double vertexChi2() const;
294  virtual double vertexNdof() const;
296  virtual double vertexNormalizedChi2() const;
298  virtual double vertexCovariance(int i, int j) const;
302  virtual void fillVertexCovariance(CovarianceMatrix & v) const;
305  virtual bool hasMasterClone() const;
308  virtual const reco::CandidateBaseRef & masterClone() const;
311  virtual bool hasMasterClonePtr() const;
314 
315  virtual const reco::CandidatePtr & masterClonePtr() const;
316 
318  template<typename Ref>
319  Ref masterRef() const { return masterClone().template castTo<Ref>(); }
321 
322  /* template<typename T> T get() const { */
323  /* if ( hasMasterClone() ) return masterClone()->get<T>(); */
324  /* else return reco::get<T>( * this ); */
325  /* } */
326  /* /// get a component */
327  /* template<typename T, typename Tag> T get() const { */
328  /* if ( hasMasterClone() ) return masterClone()->get<T, Tag>(); */
329  /* else return reco::get<T, Tag>( * this ); */
330  /* } */
331  /* /// get a component */
332  /* template<typename T> T get( size_type i ) const { */
333  /* if ( hasMasterClone() ) return masterClone()->get<T>( i ); */
334  /* else return reco::get<T>( * this, i ); */
335  /* } */
336  /* /// get a component */
337  /* template<typename T, typename Tag> T get( size_type i ) const { */
338  /* if ( hasMasterClone() ) return masterClone()->get<T, Tag>( i ); */
339  /* else return reco::get<T, Tag>( * this, i ); */
340  /* } */
341  /* /// number of components */
342  /* template<typename T> size_type numberOf() const { */
343  /* if ( hasMasterClone() ) return masterClone()->numberOf<T>(); */
344  /* else return reco::numberOf<T>( * this ); */
345  /* } */
346  /* /// number of components */
347  /* template<typename T, typename Tag> size_type numberOf() const { */
348  /* if ( hasMasterClone() ) return masterClone()->numberOf<T, Tag>(); */
349  /* else return reco::numberOf<T, Tag>( * this ); */
350  /* } */
351 
352  /* template<typename S> */
353  /* struct daughter_iterator { */
354  /* typedef boost::filter_iterator<S, const_iterator> type; */
355  /* }; */
356 
357  /* template<typename S> */
358  /* typename daughter_iterator<S>::type beginFilter( const S & s ) const { */
359  /* return boost::make_filter_iterator(s, begin(), end()); */
360  /* } */
361  /* template<typename S> */
362  /* typename daughter_iterator<S>::type endFilter( const S & s ) const { */
363  /* return boost::make_filter_iterator(s, end(), end()); */
364  /* } */
365 
366 
367  virtual bool isElectron() const { return false; }
368  virtual bool isMuon() const { return false; }
369  virtual bool isStandAloneMuon() const { return ((qualityFlags_ & muonFlagsMask) >> muonFlagsShift) & 1; }
370  virtual bool isGlobalMuon() const { return ((qualityFlags_ & muonFlagsMask) >> muonFlagsShift) & 2; }
371  virtual bool isTrackerMuon() const { return false; }
372  virtual bool isCaloMuon() const { return false; }
373  virtual bool isPhoton() const { return false; }
374  virtual bool isConvertedPhoton() const { return false; }
375  virtual bool isJet() const { return false; }
376 
377  protected:
383  void pack(bool unpackAfterwards=true) ;
384  void unpack() const ;
385  void packVtx(bool unpackAfterwards=true) ;
386  void unpackVtx() const ;
387  void maybeUnpackBoth() const { if (!unpacked_) unpack(); if (!unpackedVtx_) unpackVtx(); }
388  void packBoth() { pack(false); packVtx(false); unpack(); unpackVtx(); } // do it this way, so that we don't loose precision on the angles before computing dxy,dz
389 
394  mutable Point vertex_;
395  mutable float dxy_, dz_, dphi_;
397  int pdgId_;
398  uint16_t qualityFlags_;
401  // is the momentum p4 unpacked
402  mutable bool unpacked_;
403  // are the dxy, dz and vertex unpacked
404  mutable bool unpackedVtx_;
407  uint8_t packedHits_;
409  uint8_t normalizedChi2_;
410 // uint8_t numberOfPixelHits_;
411  // uint8_t numberOfHits_;
412 
414  virtual bool overlap( const reco::Candidate & ) const;
415  template<typename, typename, typename> friend struct component;
416  friend class ::OverlapChecker;
417  friend class ShallowCloneCandidate;
419 
425  };
426  private:
427  // const iterator implementation
429  // iterator implementation
431  };
432 
433  typedef std::vector<pat::PackedCandidate> PackedCandidateCollection;
436 }
437 
438 #endif
virtual void setP4(const LorentzVector &p4)
set 4-momentum
void setMuonID(bool isStandAlone, bool isGlobal)
virtual void setVertex(const Point &vertex)
set vertex
int i
Definition: DBlmapReader.cc:9
virtual void setTrackProperties(const reco::Track &tk, const reco::Track::CovarianceMatrix &covariance)
set impact parameters covariance
virtual size_t numberOfMothers() const
number of mothers
virtual bool hasMasterClonePtr() const
virtual float phi() const
momentum azimuthal angle
virtual void setPdgId(int pdgId)
virtual float massSqr() const
mass squared
math::XYZVector Vector
point in the space
Definition: Candidate.h:47
virtual float dzError() const
uncertainty on dz
int numberOfHits() const
size_t size_type
Definition: Candidate.h:34
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:109
void setFromPV(PVAssoc fromPV)
PackedCandidate(const PolarLorentzVector &p4, const Point &vtx, float phiAtVtx, int pdgId, const reco::VertexRef &pv)
candidate::const_iterator const_iterator
Definition: Candidate.h:35
virtual const_iterator begin() const
first daughter const_iterator
virtual bool isCaloMuon() const
void maybeUnpackBoth() const
static const unsigned int longLivedTag
long lived flag
int numberOfValidHits() const
Definition: HitPattern.h:589
virtual double y() const
rapidity
virtual void fillVertexCovariance(CovarianceMatrix &v) const
fill SMatrix
static const unsigned int massConstraintTag
do mass constraint flag
virtual double mtSqr() const
transverse mass squared
std::vector< pat::PackedCandidate > PackedCandidateCollection
math::XYZPoint Point
point in the space
virtual double vertexNdof() const
virtual void setStatus(int status)
set status word
virtual double vx() const
x coordinate of vertex position
virtual bool overlap(const reco::Candidate &) const
check overlap with another Candidate
virtual bool isPhoton() const
virtual void setP4(const PolarLorentzVector &p4)
set 4-momentum
virtual float pt() const
transverse momentum
const reco::VertexRef vertexRef() const
void unpackVtx() const
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:180
Ref masterRef() const
cast master clone reference to a concrete type
PolarLorentzVector p4_
the four vector
void setVertexRef(const reco::VertexRef &vertexRef)
set reference to the primary vertex
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
reco::CandidateCollection daughters
collection of daughter candidates
virtual const reco::Candidate * daughter(size_type) const
return daughter at a given position (throws an exception)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
virtual void setTrackProperties(const reco::Track &tk)
virtual bool isTrackerMuon() const
virtual void setLongLived()
set long lived flag
virtual bool isConvertedPhoton() const
friend struct component
virtual void setMassConstraint()
set mass constraint flag
virtual bool isMuon() const
virtual PackedCandidate * clone() const
returns a clone of the Candidate object
math::XYZTLorentzVector LorentzVector
Lorentz vector.
T sqrt(T t)
Definition: SSEVec.h:48
virtual double py() const
y coordinate of momentum vector
bool trackHighPurity() const
true if the track had the highPurity quality bit
LostInnerHits lostInnerHits() const
virtual int status() const
status word
PackedCandidate()
default constructor
virtual double rapidity() const
rapidity
math::XYZPoint Point
virtual const reco::CandidateBaseRef & masterClone() const
edm::Ref< pat::PackedCandidateCollection > PackedCandidateRef
friend class ShallowCloneCandidate
virtual double energy() const
energy
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual bool isJet() const
int j
Definition: DBlmapReader.cc:9
virtual int charge() const
electric charge
virtual double theta() const
momentum polar angle
virtual int threeCharge() const
electric charge
virtual double vertexNormalizedChi2() const
chi-squared divided by n.d.o.f.
Point vertex_
vertex position
virtual void setThreeCharge(int threecharge)
set electric charge
virtual double et() const
transverse energy
virtual void setPz(double pz)
virtual double mt() const
transverse mass
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:221
int numberOfPixelHits() const
virtual const reco::CandidatePtr & masterClonePtr() const
#define M_PI
virtual double vy() const
y coordinate of vertex position
virtual const PolarLorentzVector & polarP4() const
four-momentum Lorentz vector
LostInnerHits
Enumerator specifying the.
uint8_t normalizedChi2_
track quality information
candidate::iterator iterator
Definition: Candidate.h:36
virtual size_t numberOfSourceCandidatePtrs() const
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
virtual Vector momentum() const
spatial momentum vector
virtual bool isGlobalMuon() const
PackedCandidate(const LorentzVector &p4, const Point &vtx, float phiAtVtx, int pdgId, const reco::VertexRef &pv)
virtual bool massConstraint() const
do mass constraint?
float dxydxy_
IP covariance.
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::RefVector< pat::PackedCandidateCollection > PackedCandidateRefVector
edm::Ref< reco::VertexCollection > pvRef_
Ref to primary vertex.
virtual Vector boostToCM() const
uint16_t packedCovarianceDxyDxy_
virtual double vz() const
z coordinate of vertex position
reco::candidate::const_iterator_imp_specific< daughters > const_iterator_imp_specific
virtual void setMass(double m)
set particle mass
virtual reco::CandidatePtr sourceCandidatePtr(size_type i) const
virtual float phiAtVtx() const
momentum azimuthal angle from the track (normally identical to phi())
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
int pdgId_
PDG identifier.
virtual const Point & vertex() const
vertex position
virtual size_t numberOfDaughters() const
number of daughters
CovarianceMatrix vertexCovariance() const
return SMatrix
virtual bool isStandAloneMuon() const
friend class ShallowClonePtrCandidate
virtual bool longLived() const
is long lived?
virtual const reco::Candidate * mother(size_type) const
return mother at a given position (throws an exception)
virtual void setCharge(int charge)
set electric charge
virtual float mass() const
mass
void packVtx(bool unpackAfterwards=true)
void setLostInnerHits(LostInnerHits hits)
int numberOfValidPixelHits() const
Definition: HitPattern.h:601
virtual const LorentzVector & p4() const
four-momentum Lorentz vecto r
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
void pack(bool unpackAfterwards=true)
virtual bool isElectron() const
get a component
virtual double px() const
x coordinate of momentum vector
volatile std::atomic< bool > shutdown_flag false
virtual float eta() const
momentum pseudorapidity
virtual float dxy() const
dxy with respect to the PV ref
virtual const_iterator end() const
last daughter const_iterator
virtual double vertexChi2() const
chi-squares
void setTrackHighPurity(bool highPurity)
set to true if the track had the highPurity quality bit
const PVAssoc fromPV() const
virtual ~PackedCandidate()
destructor
virtual int pdgId() const
PDG identifier.
virtual reco::Track pseudoTrack() const
Return by value (no caching heavy function) a pseudo track made with candidate kinematics, parameterized error for eta,phi,pt and full IP covariance.
virtual double pz() const
z coordinate of momentum vector
virtual float dxyError() const
uncertainty on dxy
reco::candidate::iterator_imp_specific< daughters > iterator_imp_specific
virtual float dz() const
dz with respect to the PV ref
math::XYZVector Vector
point in the space
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:70
virtual double p() const
magnitude of momentum vector
virtual bool hasMasterClone() const
PackedCandidate(const reco::Candidate &c, const reco::VertexRef &pv)
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:43