CMS 3D CMS Logo

SiStripElectron.h
Go to the documentation of this file.
1 #ifndef EgammaCandidates_SiStripElectron_h
2 #define EgammaCandidates_SiStripElectron_h
3 // -*- C++ -*-
4 //
5 // Package: EgammaCandidates
6 // Class : SiStripElectron
7 //
16 //
17 // Original Author: Jim Pivarski
18 // Created: Fri May 26 15:43:14 EDT 2006
19 //
20 
21 // system include files
22 
23 #include <vector>
24 
25 // user include files
26 
30 
31 // forward declarations
32 
33 namespace reco {
34 
35  class SiStripElectron : public RecoCandidate {
36  public:
41  Charge q,
42  const std::vector<SiStripRecHit2D>& rphiRecHits,
43  const std::vector<SiStripRecHit2D>& stereoRecHits,
45  double phiVsRSlope,
46  double phiAtOrigin,
47  double chi2,
48  int ndof,
49  double pt,
50  double pz,
51  double zVsRSlope,
52  unsigned int numberOfStereoHits,
53  unsigned int numberOfBarrelRphiHits,
54  unsigned int numberOfEndcapZphiHits)
55  : RecoCandidate(q, PtEtaPhiMass(pt,etaFromRZ(pt,pz), phiAtOrigin, 0.000510f), Point(0,0,0), -11 * q )
56  , superCluster_(superCluster)
57  , rphiRecHits_(rphiRecHits)
58  , stereoRecHits_(stereoRecHits)
59  , superClusterPhiVsRSlope_(superClusterPhiVsRSlope)
60  , phiVsRSlope_(phiVsRSlope)
61  , phiAtOrigin_(phiAtOrigin)
62  , chi2_(chi2)
63  , ndof_(ndof)
64  , zVsRSlope_(zVsRSlope)
65  , numberOfStereoHits_(numberOfStereoHits)
66  , numberOfBarrelRphiHits_(numberOfBarrelRphiHits)
67  , numberOfEndcapZphiHits_(numberOfEndcapZphiHits) { }
68 
70  template<typename P4>
71  SiStripElectron( Charge q, const P4 & p4, const Point & vtx = Point( 0, 0, 0 ) ) :
72  RecoCandidate( q, p4, vtx, -11 * q ) { }
74  ~SiStripElectron() override;
76  SiStripElectron * clone() const override;
78  reco::SuperClusterRef superCluster() const override;
79 
81  const std::vector<SiStripRecHit2D>& rphiRecHits() const { return rphiRecHits_; }
83  const std::vector<SiStripRecHit2D>& stereoRecHits() const { return stereoRecHits_; }
84 
88  double phiVsRSlope() const { return phiVsRSlope_; }
90  double phiAtOrigin() const { return phiAtOrigin_; }
92  double chi2() const { return chi2_; }
94  int ndof() const { return ndof_; }
95 
97  double zVsRSlope() const { return zVsRSlope_; }
98 
100  unsigned int numberOfStereoHits() const { return numberOfStereoHits_; }
102  unsigned int numberOfBarrelRphiHits() const { return numberOfBarrelRphiHits_; }
104  unsigned int numberOfEndcapZphiHits() const { return numberOfEndcapZphiHits_; }
105 
106  bool isElectron() const override;
107  private:
109  bool overlap( const Candidate & ) const override;
112  std::vector<SiStripRecHit2D> rphiRecHits_;
113  std::vector<SiStripRecHit2D> stereoRecHits_;
114 
116  double phiVsRSlope_;
117  double phiAtOrigin_;
118  double chi2_;
119  int ndof_;
120 
121  double zVsRSlope_;
122 
123  unsigned int numberOfStereoHits_;
126  };
127 }
128 
129 #endif
int Charge
electric charge type
Definition: Candidate.h:35
double chi2() const
returns chi^2 of fit to tracker hits
SiStripElectron()
default constructor
int ndof() const
returns number of degrees of freedom of fit to tracker hits
double phiAtOrigin() const
returns phi(r=0) intercept from fit to tracker hits
SiStripElectron(Charge q, const P4 &p4, const Point &vtx=Point(0, 0, 0))
constructor from RecoCandidate
bool isElectron() const override
unsigned int numberOfBarrelRphiHits_
SiStripElectron(const reco::SuperClusterRef &superCluster, Charge q, const std::vector< SiStripRecHit2D > &rphiRecHits, const std::vector< SiStripRecHit2D > &stereoRecHits, double superClusterPhiVsRSlope, double phiVsRSlope, double phiAtOrigin, double chi2, int ndof, double pt, double pz, double zVsRSlope, unsigned int numberOfStereoHits, unsigned int numberOfBarrelRphiHits, unsigned int numberOfEndcapZphiHits)
constructor from band algorithm
double pt() const final
transverse momentum
std::vector< SiStripRecHit2D > rphiRecHits_
unsigned int numberOfBarrelRphiHits() const
returns number of barrel rphi hits in phi band
const std::vector< SiStripRecHit2D > & rphiRecHits() const
reference to the rphiRecHits identified as belonging to an electron
~SiStripElectron() override
destructor
unsigned int numberOfStereoHits() const
returns number of stereo hits in phi band (barrel + endcap)
double zVsRSlope() const
returns z(r) slope fit from stereo tracker hits (constrained to pass through supercluster) ...
double pz() const final
z coordinate of momentum vector
double phiVsRSlope() const
returns phi(r) slope from fit to tracker hits
double superClusterPhiVsRSlope() const
returns phi(r) projection from supercluster
unsigned int numberOfEndcapZphiHits() const
returns number of endcap zphi hits in phi band
double f[11][100]
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
reco::SuperClusterRef superCluster() const override
reference to a SuperCluster
unsigned int numberOfStereoHits_
reco::SuperClusterRef superCluster_
reference to a SuperCluster
fixed size matrix
unsigned int numberOfEndcapZphiHits_
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
std::vector< SiStripRecHit2D > stereoRecHits_
bool overlap(const Candidate &) const override
check overlap with another candidate
SiStripElectron * clone() const override
returns a clone of the candidate
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:27
const std::vector< SiStripRecHit2D > & stereoRecHits() const
reference to the stereoRecHits identified as belonging to an electron