CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // $Id: SiStripElectron.h,v 1.15 2012/10/14 08:45:37 innocent Exp $
20 //
21 
22 // system include files
23 
24 #include <vector>
25 
26 // user include files
27 
31 
32 // forward declarations
33 
34 namespace reco {
35 
36  class SiStripElectron : public RecoCandidate {
37  public:
42  Charge q,
43  const std::vector<SiStripRecHit2D>& rphiRecHits,
44  const std::vector<SiStripRecHit2D>& stereoRecHits,
46  double phiVsRSlope,
47  double phiAtOrigin,
48  double chi2,
49  int ndof,
50  double pt,
51  double pz,
52  double zVsRSlope,
53  unsigned int numberOfStereoHits,
54  unsigned int numberOfBarrelRphiHits,
55  unsigned int numberOfEndcapZphiHits)
56  : RecoCandidate(q, PtEtaPhiMass(pt,etaFromRZ(pt,pz), phiAtOrigin, 0.000510f), Point(0,0,0), -11 * q )
57  , superCluster_(superCluster)
58  , rphiRecHits_(rphiRecHits)
59  , stereoRecHits_(stereoRecHits)
60  , superClusterPhiVsRSlope_(superClusterPhiVsRSlope)
61  , phiVsRSlope_(phiVsRSlope)
62  , phiAtOrigin_(phiAtOrigin)
63  , chi2_(chi2)
64  , ndof_(ndof)
65  , zVsRSlope_(zVsRSlope)
66  , numberOfStereoHits_(numberOfStereoHits)
67  , numberOfBarrelRphiHits_(numberOfBarrelRphiHits)
68  , numberOfEndcapZphiHits_(numberOfEndcapZphiHits) { }
69 
71  template<typename P4>
72  SiStripElectron( Charge q, const P4 & p4, const Point & vtx = Point( 0, 0, 0 ) ) :
73  RecoCandidate( q, p4, vtx, -11 * q ) { }
75  virtual ~SiStripElectron();
77  virtual SiStripElectron * clone() const;
79  virtual reco::SuperClusterRef superCluster() const;
80 
82  const std::vector<SiStripRecHit2D>& rphiRecHits() const { return rphiRecHits_; }
84  const std::vector<SiStripRecHit2D>& stereoRecHits() const { return stereoRecHits_; }
85 
89  double phiVsRSlope() const { return phiVsRSlope_; }
91  double phiAtOrigin() const { return phiAtOrigin_; }
93  double chi2() const { return chi2_; }
95  int ndof() const { return ndof_; }
96 
98  double zVsRSlope() const { return zVsRSlope_; }
99 
101  unsigned int numberOfStereoHits() const { return numberOfStereoHits_; }
103  unsigned int numberOfBarrelRphiHits() const { return numberOfBarrelRphiHits_; }
105  unsigned int numberOfEndcapZphiHits() const { return numberOfEndcapZphiHits_; }
106 
107  bool isElectron() const;
108  private:
110  virtual bool overlap( const Candidate & ) const;
113  std::vector<SiStripRecHit2D> rphiRecHits_;
114  std::vector<SiStripRecHit2D> stereoRecHits_;
115 
117  double phiVsRSlope_;
118  double phiAtOrigin_;
119  double chi2_;
120  int ndof_;
121 
122  double zVsRSlope_;
123 
124  unsigned int numberOfStereoHits_;
127  };
128 }
129 
130 #endif
int Charge
electric charge type
Definition: Candidate.h:39
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
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
SiStripElectron(Charge q, const P4 &p4, const Point &vtx=Point(0, 0, 0))
constructor from RecoCandidate
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
std::vector< SiStripRecHit2D > rphiRecHits_
virtual double pz() const GCC11_FINAL
z coordinate of momentum vector
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
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 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]
unsigned int numberOfStereoHits_
reco::SuperClusterRef superCluster_
reference to a SuperCluster
virtual bool overlap(const Candidate &) const
check overlap with another candidate
virtual reco::SuperClusterRef superCluster() const
reference to a SuperCluster
unsigned int numberOfEndcapZphiHits_
virtual ~SiStripElectron()
destructor
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
std::vector< SiStripRecHit2D > stereoRecHits_
virtual float pt() const GCC11_FINAL
transverse momentum
bool isElectron() const
virtual SiStripElectron * clone() const
returns a clone of the candidate
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:30
const std::vector< SiStripRecHit2D > & stereoRecHits() const
reference to the stereoRecHits identified as belonging to an electron