CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgHLTOffEle.h
Go to the documentation of this file.
1 #ifndef DQMOFFLINE_TRIGGER_EGHLTOFFELE
2 #define DQMOFFLINE_TRIGGER_EGHLTOFFELE
3 
4 //class: EgHLTOffEle
5 //
6 //author: Sam Harper (July 2008)
7 //
8 //
9 //aim: to allow easy access to electron ID variables
10 // currently the CMSSW electron classes are a mess with key electron selection variables not being accessable from GsfElectron
11 // this a stop gap to produce a simple electron class with all variables easily accessable via methods
12 // note as this is meant for HLT Offline DQM, I do not want the overhead of converting to pat
13 //
14 //implimentation: aims to be a wrapper for GsfElectron methods, it is hoped that in time these methods will be directly added to GsfElectron and so
15 // make this class obsolute
16 // unfortunately can not be a pure wrapper as needs to store isol and cluster shape
17 //
18 
19 
25 
28 
29 namespace egHLT {
30  class OffEle {
31 
32  public:
33  //helper struct to store the isolations
34  struct IsolData {
35  float em;
36  float hadDepth1;
37  float hadDepth2;
38  float ptTrks;
39  int nrTrks;
40  //possibly going to move these to hlt data
41  float hltHad;
42  float hltTrksEle;
43  float hltTrksPho;
44  float hltEm;
45  };
46 
47  public:
48  //helper struct to store the cluster shapes
49  struct ClusShapeData {
50  float sigmaEtaEta;
51  float sigmaIEtaIEta;
52  float sigmaPhiPhi;
53  float sigmaIPhiIPhi;
54  float e1x5Over5x5;
56  float r9;
57 
58  };
59 
60  public:
61  //helper struct to store reco approximations of variables made by HLT
62  struct HLTData {
63  float dEtaIn;
64  float dPhiIn;
65  float invEInvP;
66  };
67 
68  public:
69  //helper struct to store event-wide variables
70  struct EventData {
71  int NVertex;
72  };
73 
74  private:
75  const reco::GsfElectron* gsfEle_; //pointers to the underlying electron (we do not own this)
76 
81 
82  //these are bit-packed words telling me which cuts the electron fail (ie 0x0 is passed all cuts)
83  int cutCode_;
85  //the idea is that these are user definable cuts meant to be idenital to the specified trigger
86  //it is probably clear to the reader that I havent decided on the most efficient way to do this
87  std::vector<std::pair<TrigCodes::TrigBitSet,int> > trigCutsCutCodes_; //unsorted vector (may sort if have performance issues)
88 
89  //and these are the trigger bits stored
90  //note that the trigger bits are defined at the begining of each job
91  //and do not necessaryly map between jobs
93 
94  public:
95 
96  OffEle(const reco::GsfElectron& ele,const ClusShapeData& shapeData,const IsolData& isolData,const HLTData& hltData,const EventData& eventData):
97  gsfEle_(&ele),clusShapeData_(shapeData),isolData_(isolData),hltData_(hltData),eventData_(eventData),
99  ~OffEle(){}
100 
101 
102  //modifiers
103  int NVertex()const{return eventData_.NVertex;}
104  void setCutCode(int code){cutCode_=code;}
105  void setLooseCutCode(int code){looseCutCode_=code;}
106  //slightly inefficient way, think I can afford it and its a lot easier to just make the sorted vector outside the class
107  void setTrigCutsCutCodes(const std::vector<std::pair<TrigCodes::TrigBitSet,int> > trigCutsCutCodes){trigCutsCutCodes_=trigCutsCutCodes;}
109 
110  const reco::GsfElectron* gsfEle()const{return gsfEle_;}
111 
112  //kinematic and geometric methods
113  float et()const{return gsfEle_->et();}
114  // float et()const{return etSC();}
115  float energy()const{return gsfEle_->energy();}
116  float eta()const{return gsfEle_->eta();}
117  float phi()const{return gsfEle_->phi();}
118  float etSC()const{return gsfEle_->superCluster()->position().rho()/gsfEle_->superCluster()->position().r()*caloEnergy();}
119  float caloEnergy()const{return gsfEle_->caloEnergy();}
120  float etaSC()const{return gsfEle_->superCluster()->eta();}
121  float detEta()const{return etaSC();}
122  float phiSC()const{return gsfEle_->superCluster()->phi();}
123  float zVtx()const{return gsfEle_->TrackPositionAtVtx().z();}
124  const math::XYZTLorentzVector& p4()const{return gsfEle_->p4();}
125 
126  //classification (couldnt they have just named it 'type')
127  int classification()const{return gsfEle_->classification();}
128  bool isGap()const{return gsfEle_->isEBGap() || gsfEle_->isEEGap() || gsfEle_->isEBEEGap();}
129 
130  //track methods
131  int charge()const{return gsfEle_->charge();}
132  float pVtx()const{return gsfEle_->trackMomentumAtVtx().R();}
133  float pCalo()const{return gsfEle_->trackMomentumAtCalo().R();}
134  float ptVtx()const{return gsfEle_->trackMomentumAtVtx().rho();}
135  float ptCalo()const{return gsfEle_->trackMomentumAtCalo().rho();}
136 
137 
138  //abreviations of overly long GsfElectron methods, I'm sorry but if you cant figure out what hOverE() means, you shouldnt be using this class
139  float hOverE()const{return gsfEle_->hadronicOverEm();}
144  float epIn()const{return gsfEle_->eSuperClusterOverP();}
145  float epOut()const{return gsfEle_->eSeedClusterOverPout();}
146 
147  //variables with no direct method
148  float sigmaEtaEta()const;
151  float sigmaPhiPhi()const{return clusShapeData_.sigmaPhiPhi;}
152  //float sigmaIPhiIPhi()const{return clusShapeData_.sigmaIPhiIPhi;}
154  float e1x5Over5x5()const{return clusShapeData_.e1x5Over5x5;}
155 
156  float r9()const{return clusShapeData_.r9;}
157  //float sigmaPhiPhi()const{return clusShape_!=NULL ? sqrt(clusShape_->covPhiPhi()) : 999;}
158  float bremFrac()const{return (pVtx()-pCalo())/pVtx();}
159  float invEInvP()const{return gsfEle_->caloEnergy()!=0 && gsfEle_->trackMomentumAtVtx().R()!=0. ? 1./gsfEle_->caloEnergy() - 1./gsfEle_->trackMomentumAtVtx().R() : -999;}
160  //float e9OverE25()const{return clusShape_!=NULL ? clusShape_->e3x3()/clusShape_->e5x5() : -999;}
161 
162  //isolation
163  float isolEm()const{return isolData_.em;}
164  float isolHad()const{return isolHadDepth1()+isolHadDepth2();}
165  float isolHadDepth1()const{return isolData_.hadDepth1;}
166  float isolHadDepth2()const{return isolData_.hadDepth2;}
167  float isolPtTrks()const{return isolData_.ptTrks;}
168  int isolNrTrks()const{return isolData_.nrTrks;}
169  float hltIsolTrksEle()const{return isolData_.hltTrksEle;}
170  float hltIsolTrksPho()const{return isolData_.hltTrksPho;}
171  float hltIsolHad()const{return isolData_.hltHad;}
172  float hltIsolEm()const{return isolData_.hltEm;}
173 
174  //some hlt id variables (note these are reco approximations)
175  float hltDEtaIn()const{return hltData_.dEtaIn;}
176  float hltDPhiIn()const{return hltData_.dPhiIn;}
177  float hltInvEInvP()const{return hltData_.invEInvP;}
178 
179  //ctf track accessor and validatity checker
180  reco::TrackRef ctfTrack()const{return gsfEle_->closestCtfTrackRef();} //in theory lightweight (if they follow good design),return by value
181  //track is only valid if it exists and track extra exists (track extra is only stored in reco)
183 
184 
185  //ctf track varibles, used as hlt uses this algo
186  float ctfTrkP()const{return validCTFTrack() ? ctfTrack()->p() : -999.;}
187  float ctfTrkPt()const{return validCTFTrack() ? ctfTrack()->pt() : -999.;}
188  float ctfTrkEta()const{return validCTFTrack() ? ctfTrack()->eta() : -999.;}
189  float ctfTrkChi2()const{return validCTFTrack() ? ctfTrack()->chi2() : 999.;}
190  float ctfTrkNDof()const{return validCTFTrack() ? ctfTrack()->ndof() : 999.;} //this will give chi2/ndof a valid value, perhaps rethink
191  float ctfTrkPtOuter()const{return validCTFTrack() ? ctfTrack()->outerMomentum().Perp2() : -999.;}
192  float ctfTrkPtInner()const{return validCTFTrack() ? ctfTrack()->innerMomentum().Perp2() : -999.;}
193  float ctfTrkInnerRadius()const{return validCTFTrack() ? ctfTrack()->innerPosition().Rho() : 999.;}
194  float ctfTrkOuterRadius()const{return validCTFTrack() ? ctfTrack()->outerPosition().Rho() : -999.;}
195  int ctfTrkHitsFound()const{return validCTFTrack() ? static_cast<int>(ctfTrack()->found()) : -999;}
196  int ctfTrkHitsLost()const{return validCTFTrack() ? static_cast<int>(ctfTrack()->lost()) : -999;}
197  int ctfTrkNrHits()const{return validCTFTrack() ? static_cast<int>(ctfTrack()->recHitsSize()) : -999;}
198 
199  //selection cuts
200  int cutCode()const{return cutCode_;}
201  int looseCutCode()const{return looseCutCode_;}
202 
203 
204  //trigger codes are just used as a unique identifier of the trigger, it is an error to specify more than a single bit
205  //the idea here is to allow an arbitary number of electron triggers
206  int trigCutsCutCode(const TrigCodes::TrigBitSet& trigger)const;
207  //trigger
209 
210 
211  };
212 
213 }
214 
215 
216 
217 #endif
bool isEEGap() const
Definition: GsfElectron.h:347
float e2x5MaxOver5x5() const
Definition: EgHLTOffEle.h:153
OffEle(const reco::GsfElectron &ele, const ClusShapeData &shapeData, const IsolData &isolData, const HLTData &hltData, const EventData &eventData)
Definition: EgHLTOffEle.h:96
float eta() const
Definition: EgHLTOffEle.h:116
int ctfTrkHitsFound() const
Definition: EgHLTOffEle.h:195
reco::TrackRef ctfTrack() const
Definition: EgHLTOffEle.h:180
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:196
float eSuperClusterOverP() const
Definition: GsfElectron.h:215
virtual double et() const
transverse energy
EventData eventData_
Definition: EgHLTOffEle.h:80
int looseCutCode() const
Definition: EgHLTOffEle.h:201
bool isEBEEGap() const
Definition: GsfElectron.h:343
float detEta() const
Definition: EgHLTOffEle.h:121
TrigCodes::TrigBitSet trigBits() const
Definition: EgHLTOffEle.h:208
void setTrigCutsCutCodes(const std::vector< std::pair< TrigCodes::TrigBitSet, int > > trigCutsCutCodes)
Definition: EgHLTOffEle.h:107
void setLooseCutCode(int code)
Definition: EgHLTOffEle.h:105
float caloEnergy() const
Definition: EgHLTOffEle.h:119
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision bits
float pVtx() const
Definition: EgHLTOffEle.h:132
math::XYZVectorF trackMomentumAtVtx() const
Definition: GsfElectron.h:258
float ctfTrkEta() const
Definition: EgHLTOffEle.h:188
float epIn() const
Definition: EgHLTOffEle.h:144
int ctfTrkNrHits() const
Definition: EgHLTOffEle.h:197
float r9() const
Definition: EgHLTOffEle.h:156
float invEInvP() const
Definition: EgHLTOffEle.h:159
float et() const
Definition: EgHLTOffEle.h:113
float ctfTrkInnerRadius() const
Definition: EgHLTOffEle.h:193
const reco::GsfElectron * gsfEle_
Definition: EgHLTOffEle.h:75
float sigmaEtaEtaUnCorr() const
Definition: EgHLTOffEle.h:149
bool validCTFTrack() const
Definition: EgHLTOffEle.h:182
float dPhiIn() const
Definition: EgHLTOffEle.h:141
SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:168
float zVtx() const
Definition: EgHLTOffEle.h:123
float isolHadDepth1() const
Definition: EgHLTOffEle.h:165
math::XYZPointF TrackPositionAtVtx() const
Definition: GsfElectron.h:266
virtual double eta() const
momentum pseudorapidity
float sigmaPhiPhi() const
Definition: EgHLTOffEle.h:151
TrackRef closestCtfTrackRef() const
Definition: GsfElectron.h:309
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
int charge() const
Definition: EgHLTOffEle.h:131
float hltDPhiIn() const
Definition: EgHLTOffEle.h:176
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:249
float isolEm() const
Definition: EgHLTOffEle.h:163
math::XYZVectorF trackMomentumAtCalo() const
Definition: GsfElectron.h:259
float ctfTrkP() const
Definition: EgHLTOffEle.h:186
float hltInvEInvP() const
Definition: EgHLTOffEle.h:177
virtual double energy() const
energy
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:219
float hadronicOverEm() const
Definition: GsfElectron.h:399
float bremFrac() const
Definition: EgHLTOffEle.h:158
float dEtaOut() const
Definition: EgHLTOffEle.h:143
float dEtaIn() const
Definition: EgHLTOffEle.h:140
float deltaPhiSeedClusterTrackAtCalo() const
Definition: GsfElectron.h:223
HLTData hltData_
Definition: EgHLTOffEle.h:79
float phi() const
Definition: EgHLTOffEle.h:117
float ctfTrkPtOuter() const
Definition: EgHLTOffEle.h:191
float etaSC() const
Definition: EgHLTOffEle.h:120
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:222
#define INVALID
Definition: myFilter.cc:40
virtual int charge() const
electric charge
bool isGap() const
Definition: EgHLTOffEle.h:128
float sigmaIEtaIEta() const
Definition: EgHLTOffEle.h:150
std::vector< std::pair< TrigCodes::TrigBitSet, int > > trigCutsCutCodes_
Definition: EgHLTOffEle.h:87
float eSeedClusterOverPout() const
Definition: GsfElectron.h:217
float sigmaEtaEta() const
Definition: EgHLTOffEle.cc:5
void setCutCode(int code)
Definition: EgHLTOffEle.h:104
float hltIsolHad() const
Definition: EgHLTOffEle.h:171
int classification() const
Definition: EgHLTOffEle.h:127
TrigCodes::TrigBitSet trigBits_
Definition: EgHLTOffEle.h:92
bool isEBGap() const
Definition: GsfElectron.h:344
int cutCode() const
Definition: EgHLTOffEle.h:200
int ctfTrkHitsLost() const
Definition: EgHLTOffEle.h:196
float epOut() const
Definition: EgHLTOffEle.h:145
float energy() const
Definition: EgHLTOffEle.h:115
ClusShapeData clusShapeData_
Definition: EgHLTOffEle.h:77
Classification classification() const
Definition: GsfElectron.h:602
float ptVtx() const
Definition: EgHLTOffEle.h:134
float isolHadDepth2() const
Definition: EgHLTOffEle.h:166
float hltIsolTrksPho() const
Definition: EgHLTOffEle.h:170
float etSC() const
Definition: EgHLTOffEle.h:118
float ctfTrkPt() const
Definition: EgHLTOffEle.h:187
float isolPtTrks() const
Definition: EgHLTOffEle.h:167
float ctfTrkChi2() const
Definition: EgHLTOffEle.h:189
float ctfTrkOuterRadius() const
Definition: EgHLTOffEle.h:194
float hltIsolTrksEle() const
Definition: EgHLTOffEle.h:169
float hltDEtaIn() const
Definition: EgHLTOffEle.h:175
IsolData isolData_
Definition: EgHLTOffEle.h:78
float deltaEtaSeedClusterTrackAtCalo() const
Definition: GsfElectron.h:220
float ctfTrkNDof() const
Definition: EgHLTOffEle.h:190
float ptCalo() const
Definition: EgHLTOffEle.h:135
int trigCutsCutCode(const TrigCodes::TrigBitSet &trigger) const
Definition: EgHLTOffEle.cc:18
float ctfTrkPtInner() const
Definition: EgHLTOffEle.h:192
const reco::GsfElectron * gsfEle() const
Definition: EgHLTOffEle.h:110
const math::XYZTLorentzVector & p4() const
Definition: EgHLTOffEle.h:124
float e1x5Over5x5() const
Definition: EgHLTOffEle.h:154
int NVertex() const
Definition: EgHLTOffEle.h:103
float dPhiOut() const
Definition: EgHLTOffEle.h:142
float hOverE() const
Definition: EgHLTOffEle.h:139
float phiSC() const
Definition: EgHLTOffEle.h:122
float caloEnergy() const
Definition: GsfElectron.h:682
virtual double phi() const
momentum azimuthal angle
float isolHad() const
Definition: EgHLTOffEle.h:164
int isolNrTrks() const
Definition: EgHLTOffEle.h:168
float pCalo() const
Definition: EgHLTOffEle.h:133
float hltIsolEm() const
Definition: EgHLTOffEle.h:172
void setTrigBits(TrigCodes::TrigBitSet bits)
Definition: EgHLTOffEle.h:108
std::bitset< maxNrBits_ > TrigBitSet