ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/Mods/src/PhotonIDMod.cc
Revision: 1.25
Committed: Wed Aug 3 17:15:43 2011 UTC (13 years, 9 months ago) by bendavid
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_025, Mit_025pre2, Mit_024b, Mit_025pre1, Mit_024a, Mit_024
Changes since 1.24: +16 -112 lines
Log Message:
Reorganize photon code and add new PhotonTreeWriter class

File Contents

# Content
1 // $Id: PhotonIDMod.cc,v 1.24 2011/07/27 15:17:37 bendavid Exp $
2
3 #include "TDataMember.h"
4 #include "TTree.h"
5 #include "MitPhysics/Mods/interface/PhotonIDMod.h"
6 #include "MitAna/DataTree/interface/PhotonCol.h"
7 #include "MitPhysics/Init/interface/ModNames.h"
8 #include "MitPhysics/Utils/interface/IsolationTools.h"
9 #include "MitPhysics/Utils/interface/PhotonTools.h"
10
11 using namespace mithep;
12
13 ClassImp(mithep::PhotonIDMod)
14
15 //--------------------------------------------------------------------------------------------------
16 PhotonIDMod::PhotonIDMod(const char *name, const char *title) :
17 BaseMod(name,title),
18 fPhotonBranchName (Names::gkPhotonBrn),
19 fGoodPhotonsName (ModNames::gkGoodPhotonsName),
20 fTrackBranchName (Names::gkTrackBrn),
21 fBeamspotBranchName(Names::gkBeamSpotBrn),
22 fPileUpDenName (Names::gkPileupEnergyDensityBrn),
23 fConversionName ("MergedConversions"),
24 fElectronName ("Electrons"),
25 fGoodElectronName ("Electrons"),
26 fPVName (Names::gkPVBeamSpotBrn),
27 // MC specific stuff...
28 fMCParticleName (Names::gkMCPartBrn),
29 fPileUpName (Names::gkPileupInfoBrn),
30
31 fPhotonIDType ("Custom"),
32 fPhotonIsoType ("Custom"),
33 fPhotonPtMin (15.0),
34 fHadOverEmMax (0.02),
35 fApplySpikeRemoval (kFALSE),
36 fApplyPixelSeed (kTRUE),
37 fApplyElectronVeto (kFALSE),
38 fInvertElectronVeto(kFALSE),
39 fApplyElectronVetoConvRecovery(kFALSE),
40 fApplyConversionId (kFALSE),
41 fApplyTriggerMatching(kFALSE),
42 fPhotonR9Min (0.5),
43 fPhIdType (kIdUndef),
44 fPhIsoType (kIsoUndef),
45 fFiduciality (kTRUE),
46 fEtaWidthEB (0.01),
47 fEtaWidthEE (0.028),
48 fAbsEtaMax (999.99),
49 fApplyR9Min (kFALSE),
50
51 fEffAreaEcalEE (0.089),
52 fEffAreaHcalEE (0.156),
53 fEffAreaTrackEE (0.261),
54
55 fEffAreaEcalEB (0.183),
56 fEffAreaHcalEB (0.062),
57 fEffAreaTrackEB (0.306),
58
59 fPhotons(0),
60 fTracks(0),
61 fBeamspots(0),
62 fPileUpDen(0),
63 fConversions(0),
64 fElectrons(0),
65 fGoodElectrons(0),
66 fPV(0),
67 fMCParticles (0),
68 fPileUp (0),
69
70 fPVFromBranch (true),
71 fGoodElectronsFromBranch (kTRUE),
72 fIsData(false)
73
74 {
75 // Constructor.
76 }
77
78 //--------------------------------------------------------------------------------------------------
79 void PhotonIDMod::Process()
80 {
81 // Process entries of the tree.
82
83 LoadEventObject(fPhotonBranchName, fPhotons);
84
85 const BaseVertex *bsp = NULL;
86 Double_t _tRho = -1.;
87 const TriggerObjectCol *trigObjs = 0;
88 if (fPhotons->GetEntries()>0) {
89 LoadEventObject(fTrackBranchName, fTracks);
90 LoadEventObject(fBeamspotBranchName, fBeamspots);
91 LoadEventObject(fPileUpDenName, fPileUpDen);
92 LoadEventObject(fConversionName, fConversions);
93 LoadEventObject(fElectronName, fElectrons);
94 LoadEventObject(fGoodElectronName, fGoodElectrons);
95 LoadEventObject(fPVName, fPV);
96
97 if(fBeamspots->GetEntries() > 0)
98 bsp = fBeamspots->At(0);
99
100 if(fPileUpDen->GetEntries() > 0)
101 _tRho = (Double_t) fPileUpDen->At(0)->RhoRandomLowEta();
102
103 //get trigger object collection if trigger matching is enabled
104 if (fApplyTriggerMatching) {
105 trigObjs = GetHLTObjects(fTrigObjectsName);
106 }
107
108 }
109
110 PhotonOArr *GoodPhotons = new PhotonOArr;
111 GoodPhotons->SetName(fGoodPhotonsName);
112
113 for (UInt_t i=0; i<fPhotons->GetEntries(); ++i) {
114 const Photon *ph = fPhotons->At(i);
115
116
117 if (fFiduciality == kTRUE &&
118 (ph->SCluster()->AbsEta()>=2.5 || (ph->SCluster()->AbsEta()>=1.4442 && ph->SCluster()->AbsEta()<=1.566) ) )
119 continue;
120
121 if (fInvertElectronVeto && PhotonTools::PassElectronVeto(ph,fGoodElectrons)) {
122 continue;
123 }
124
125 // ---------------------------------------------------------------------
126 // check if we use the CiC Selection. If yes, bypass all the below...
127 if(fPhIdType == kBaseLineCiC) {
128 if( PhotonTools::PassCiCSelection(ph, fPV->At(0), fTracks, fElectrons, fPV, _tRho, fPhotonPtMin, fApplyElectronVeto) )
129 GoodPhotons->Add(fPhotons->At(i));
130 continue; // go to next Photons
131 }
132 // ---------------------------------------------------------------------
133
134 if (ph->Pt() <= fPhotonPtMin)
135 continue; // add good electron
136
137 Bool_t isbarrel = ph->SCluster()->AbsEta()<1.5;
138
139 Bool_t passSpikeRemovalFilter = kTRUE;
140
141 if (ph->SCluster() && ph->SCluster()->Seed()) {
142 if(ph->SCluster()->Seed()->Energy() > 5.0 &&
143 ph->SCluster()->Seed()->EMax() / ph->SCluster()->Seed()->E3x3() > 0.95
144 ) {
145 passSpikeRemovalFilter = kFALSE;
146 }
147 }
148
149 // For Now Only use the EMax/E3x3 prescription.
150 //if(ph->SCluster()->Seed()->Energy() > 5.0 &&
151 // (1 - (ph->SCluster()->Seed()->E1x3() + ph->SCluster()->Seed()->E3x1() - 2*ph->SCluster()->Seed()->EMax())) > 0.95
152 // ) {
153 // passSpikeRemovalFilter = kFALSE;
154 //}
155 if (fApplySpikeRemoval && !passSpikeRemovalFilter) continue;
156
157 if (ph->HadOverEm() >= fHadOverEmMax)
158 continue;
159
160 if (fApplyPixelSeed == kTRUE &&
161 ph->HasPixelSeed() == kTRUE)
162 continue;
163
164 if (fApplyElectronVeto && !PhotonTools::PassElectronVeto(ph,fElectrons) ) continue;
165
166 if (fApplyElectronVetoConvRecovery && !PhotonTools::PassElectronVetoConvRecovery(ph,fElectrons,fConversions,bsp) ) continue;
167
168 if (fApplyConversionId && !PhotonTools::PassConversionId(ph,PhotonTools::MatchedConversion(ph,fConversions,bsp))) continue;
169
170 if (fApplyTriggerMatching && !PhotonTools::PassTriggerMatching(ph,trigObjs)) continue;
171
172 Bool_t idcut = kFALSE;
173 switch (fPhIdType) {
174 case kTight:
175 idcut = ph->IsTightPhoton();
176 break;
177 case kLoose:
178 idcut = ph->IsLoosePhoton();
179 break;
180 case kLooseEM:
181 idcut = ph->IsLooseEM();
182 break;
183 case kCustomId:
184 idcut = kTRUE;
185 default:
186 break;
187 }
188
189 if (!idcut)
190 continue;
191
192 Bool_t isocut = kFALSE;
193 switch (fPhIsoType) {
194 case kNoIso:
195 isocut = kTRUE;
196 break;
197 case kCombinedIso:
198 {
199 Double_t totalIso = ph->HollowConeTrkIsoDr04()+
200 ph->EcalRecHitIsoDr04() +
201 ph->HcalTowerSumEtDr04();
202 if (totalIso/ph->Pt() < 0.25)
203 isocut = kTRUE;
204 }
205 break;
206 case kCustomIso:
207 {
208 if ( ph->HollowConeTrkIsoDr04() < (1.5 + 0.001*ph->Pt()) && ph->EcalRecHitIsoDr04()<(2.0+0.006*ph->Pt()) && ph->HcalTowerSumEtDr04()<(2.0+0.0025*ph->Pt()) )
209 isocut = kTRUE;
210 }
211 break;
212
213 case kMITPUCorrected:
214 {
215 // compute the PU corrections only if Rho is available
216 // ... otherwise (_tRho = -1.0) it's the std isolation
217 isocut = kTRUE;
218 Double_t fEffAreaEcal = fEffAreaEcalEB;
219 Double_t fEffAreaHcal = fEffAreaHcalEB;
220 Double_t fEffAreaTrack = fEffAreaTrackEB;
221
222 if( !isbarrel ) {
223 fEffAreaEcal = fEffAreaEcalEE;
224 fEffAreaHcal = fEffAreaHcalEE;
225 fEffAreaTrack = fEffAreaTrackEE;
226 }
227
228 Double_t EcalCorrISO = ph->EcalRecHitIsoDr04();
229 if(_tRho > -0.5 ) EcalCorrISO -= _tRho * fEffAreaEcal;
230 if ( EcalCorrISO > (2.0+0.006*ph->Pt()) ) isocut = kFALSE;
231 if ( isocut ) {
232 Double_t HcalCorrISO = ph->HcalTowerSumEtDr04();
233 if(_tRho > -0.5 ) HcalCorrISO -= _tRho * fEffAreaHcal;
234 if ( HcalCorrISO > (2.0+0.0025*ph->Pt()) ) isocut = kFALSE;
235 }
236 if ( isocut ) {
237 Double_t TrackCorrISO = IsolationTools::TrackIsolationNoPV(ph, bsp, 0.4, 0.04, 0.0, 0.015, 0.1, TrackQuality::highPurity, fTracks);
238 if(_tRho > -0.5 )
239 TrackCorrISO -= _tRho * fEffAreaTrack;
240 if ( TrackCorrISO > (1.5 + 0.001*ph->Pt()) ) isocut = kFALSE;
241 }
242 break;
243 }
244 default:
245 break;
246 }
247
248 if (!isocut)
249 continue;
250
251 if ( fApplyR9Min && ph->R9() <= fPhotonR9Min)
252 continue;
253
254 if ((isbarrel && ph->CoviEtaiEta() >= fEtaWidthEB) ||
255 (!isbarrel && ph->CoviEtaiEta() >= fEtaWidthEE))
256 continue;
257
258 if (ph->AbsEta() >= fAbsEtaMax)
259 continue;
260
261
262 // add good electron
263 GoodPhotons->Add(fPhotons->At(i));
264 }
265
266 // sort according to pt
267 GoodPhotons->Sort();
268
269 // add to event for other modules to use
270 AddObjThisEvt(GoodPhotons);
271
272 return;
273
274 }
275
276 //--------------------------------------------------------------------------------------------------
277 void PhotonIDMod::SlaveBegin()
278 {
279 // Run startup code on the computer (slave) doing the actual analysis. Here,
280 // we just request the photon collection branch.
281
282 ReqEventObject(fPhotonBranchName, fPhotons, kTRUE);
283 ReqEventObject(fTrackBranchName, fTracks, kTRUE);
284 ReqEventObject(fBeamspotBranchName, fBeamspots, kTRUE);
285 ReqEventObject(fConversionName, fConversions, kTRUE);
286 ReqEventObject(fElectronName, fElectrons, kTRUE);
287 ReqEventObject(fGoodElectronName, fGoodElectrons, fGoodElectronsFromBranch);
288 ReqEventObject(fPVName, fPV, fPVFromBranch);
289 ReqEventObject(fPileUpDenName, fPileUpDen, kTRUE);
290 if (!fIsData) {
291 ReqBranch(fPileUpName, fPileUp);
292 ReqBranch(fMCParticleName, fMCParticles);
293 }
294
295 if (fPhotonIDType.CompareTo("Tight") == 0)
296 fPhIdType = kTight;
297 else if (fPhotonIDType.CompareTo("Loose") == 0)
298 fPhIdType = kLoose;
299 else if (fPhotonIDType.CompareTo("LooseEM") == 0)
300 fPhIdType = kLooseEM;
301 else if (fPhotonIDType.CompareTo("Custom") == 0)
302 fPhIdType = kCustomId;
303 else if (fPhotonIDType.CompareTo("BaseLineCiC") == 0) {
304 fPhIdType = kBaseLineCiC;
305 fPhotonIsoType = "NoIso";
306 } else {
307 SendError(kAbortAnalysis, "SlaveBegin",
308 "The specified photon identification %s is not defined.",
309 fPhotonIDType.Data());
310 return;
311 }
312
313 if (fPhotonIsoType.CompareTo("NoIso") == 0 )
314 fPhIsoType = kNoIso;
315 else if (fPhotonIsoType.CompareTo("CombinedIso") == 0 )
316 fPhIsoType = kCombinedIso;
317 else if (fPhotonIsoType.CompareTo("Custom") == 0 )
318 fPhIsoType = kCustomIso;
319 else if (fPhotonIsoType.CompareTo("MITPUCorrected") == 0 )
320 fPhIsoType = kMITPUCorrected;
321 else {
322 SendError(kAbortAnalysis, "SlaveBegin",
323 "The specified photon isolation %s is not defined.",
324 fPhotonIsoType.Data());
325 return;
326 }
327
328
329 }