CMS 3D CMS Logo

PatPhotonSimpleAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PatPhotonSimpleAnalyzer
4 // Class: PatPhotonSimpleAnalyzer
5 //
14 //
15 // Original Author: M.B. Anderson
16 // based on simple photon analyzer by: J. Stilley, A. Askew
17 //
18 // Created: Wed Sep 23 12:00:01 CDT 2008
19 //
20 
26 
27 #include "TFile.h"
28 #include "TH1.h"
29 #include "TTree.h"
30 
31 #include <memory>
32 #include <string>
33 
35 public:
37  ~PatPhotonSimpleAnalyzer() override;
38 
39  void analyze(const edm::Event&, const edm::EventSetup&) override;
40  void beginJob() override;
41  void endJob() override;
42 
43 private:
44  std::string outputFile_; // output file
45  double minPhotonEt_; // minimum photon Et
46  double minPhotonAbsEta_; // min and
47  double maxPhotonAbsEta_; // max abs(eta)
48  double minPhotonR9_; // minimum R9 = E(3x3)/E(SuperCluster)
49  double maxPhotonHoverE_; // maximum HCAL / ECAL
50  bool createPhotonTTree_; // Create a TTree of photon variables
51 
52  // Will be used for creating TTree of photons.
53  // These names did not have to match those from a phtn->...
54  // but do match for clarity.
62  float isEBGap;
63  float isEEGap;
64  float isEBEEGap;
65  float r9;
66  float et;
67  float eta;
68  float phi;
70  float ecalIso;
71  float hcalIso;
72  float trackIso;
73  };
75 
76  // root file to store histograms
77  TFile* rootFile_;
78 
79  // data members for histograms to be filled
80 
81  // PhotonID Histograms
88  TH1F* h_ebgap_;
89  TH1F* h_eeGap_;
90  TH1F* h_ebeeGap_;
91  TH1F* h_r9_;
92 
93  // Photon Histograms
94  TH1F* h_photonEt_;
95  TH1F* h_photonEta_;
96  TH1F* h_photonPhi_;
97  TH1F* h_hadoverem_;
98 
99  // Photon's SuperCluster Histograms
104 
105  // Composite or Other Histograms
108  TH1F* h_nPho_;
109 
110  // TTree
112 };
113 
116 
117 using namespace std;
118 
120 // Constructor //
123  // Read Parameters from configuration file
124 
125  // output filename
126  outputFile_ = ps.getParameter<std::string>("outputFile");
127  // Read variables that must be passed to allow a
128  // supercluster to be placed in histograms as a photon.
129  minPhotonEt_ = ps.getParameter<double>("minPhotonEt");
130  minPhotonAbsEta_ = ps.getParameter<double>("minPhotonAbsEta");
131  maxPhotonAbsEta_ = ps.getParameter<double>("maxPhotonAbsEta");
132  minPhotonR9_ = ps.getParameter<double>("minPhotonR9");
133  maxPhotonHoverE_ = ps.getParameter<double>("maxPhotonHoverE");
134 
135  // Read variable to that decidedes whether
136  // a TTree of photons is created or not
137  createPhotonTTree_ = ps.getParameter<bool>("createPhotonTTree");
138 
139  // open output file to store histograms
140  rootFile_ = TFile::Open(outputFile_.c_str(), "RECREATE");
141 }
142 
144 // Destructor //
147  // do anything here that needs to be done at desctruction time
148  // (e.g. close files, deallocate resources etc.)
149 
150  delete rootFile_;
151 }
152 
154 // method called once each job just before starting event loop //
157  // go to *OUR* rootfile
158  rootFile_->cd();
159 
160  // Book Histograms
161  // PhotonID Histograms
162  h_isoEcalRecHit_ = new TH1F("photonEcalIso", "Ecal Rec Hit Isolation", 100, 0, 100);
163  h_isoHcalRecHit_ = new TH1F("photonHcalIso", "Hcal Rec Hit Isolation", 100, 0, 100);
164  h_trk_pt_solid_ = new TH1F("photonTrackSolidIso", "Sum of track pT in a cone of #DeltaR", 100, 0, 100);
165  h_trk_pt_hollow_ = new TH1F("photonTrackHollowIso", "Sum of track pT in a hollow cone", 100, 0, 100);
166  h_ntrk_solid_ = new TH1F("photonTrackCountSolid", "Number of tracks in a cone of #DeltaR", 100, 0, 100);
167  h_ntrk_hollow_ = new TH1F("photonTrackCountHollow", "Number of tracks in a hollow cone", 100, 0, 100);
168  h_ebgap_ = new TH1F("photonInEBgap", "Ecal Barrel gap flag", 2, -0.5, 1.5);
169  h_eeGap_ = new TH1F("photonInEEgap", "Ecal Endcap gap flag", 2, -0.5, 1.5);
170  h_ebeeGap_ = new TH1F("photonInEEgap", "Ecal Barrel/Endcap gap flag", 2, -0.5, 1.5);
171  h_r9_ = new TH1F("photonR9", "R9 = E(3x3) / E(SuperCluster)", 300, 0, 3);
172 
173  // Photon Histograms
174  h_photonEt_ = new TH1F("photonEt", "Photon E_{T}", 200, 0, 200);
175  h_photonEta_ = new TH1F("photonEta", "Photon #eta", 200, -4, 4);
176  h_photonPhi_ = new TH1F("photonPhi", "Photon #phi", 200, -1. * M_PI, M_PI);
177  h_hadoverem_ = new TH1F("photonHoverE", "Hadronic over EM", 200, 0, 1);
178 
179  // Photon's SuperCluster Histograms
180  h_photonScEt_ = new TH1F("photonScEt", "Photon SuperCluster E_{T}", 200, 0, 200);
181  h_photonScEta_ = new TH1F("photonScEta", "Photon #eta", 200, -4, 4);
182  h_photonScPhi_ = new TH1F("photonScPhi", "Photon #phi", 200, -1. * M_PI, M_PI);
183  h_photonScEtaWidth_ = new TH1F("photonScEtaWidth", "#eta-width", 100, 0, .1);
184 
185  // Composite or Other Histograms
186  h_photonInAnyGap_ = new TH1F("photonInAnyGap", "Photon in any gap flag", 2, -0.5, 1.5);
187  h_nPassingPho_ = new TH1F("photonPassingCount", "Total number photons (0=NotPassing, 1=Passing)", 2, -0.5, 1.5);
188  h_nPho_ = new TH1F("photonCount", "Number of photons passing cuts in event", 10, 0, 10);
189 
190  // Create a TTree of photons if set to 'True' in config file
191  if (createPhotonTTree_) {
192  tree_PhotonAll_ = new TTree("TreePhotonAll", "Reconstructed Photon");
193  tree_PhotonAll_->Branch(
194  "recPhoton",
195  &recPhoton.isolationEcalRecHit,
196  "isolationEcalRecHit/"
197  "F:isolationHcalRecHit:isolationSolidTrkCone:isolationHollowTrkCone:nTrkSolidCone:nTrkHollowCone:isEBGap:"
198  "isEEGap:isEBEEGap:r9:et:eta:phi:hadronicOverEm:ecalIso:hcalIso:trackIso");
199  }
200 }
201 
203 // method called to for each event //
206  using namespace std;
207  using namespace edm;
208 
209  // Grab pat::Photon
210  Handle<View<pat::Photon> > photonHandle;
211  evt.getByLabel("selectedLayer1Photons", photonHandle);
212  View<pat::Photon> photons = *photonHandle;
213 
214  int photonCounter = 0;
215 
216  for (int i = 0; i < int(photons.size()); i++) {
217  pat::Photon currentPhoton = photons.at(i);
218 
219  float photonEt = currentPhoton.et();
220  float superClusterEt =
221  (currentPhoton.superCluster()->energy()) / (cosh(currentPhoton.superCluster()->position().eta()));
222 
223  // Only store photon candidates (SuperClusters) that pass some simple cuts
224  bool passCuts = (photonEt > minPhotonEt_) && (fabs(currentPhoton.eta()) > minPhotonAbsEta_) &&
225  (fabs(currentPhoton.eta()) < maxPhotonAbsEta_) && (currentPhoton.r9() > minPhotonR9_) &&
226  (currentPhoton.hadronicOverEm() < maxPhotonHoverE_);
227 
228  if (passCuts) {
230  // fill histograms //
232  // PhotonID Variables
233  h_isoEcalRecHit_->Fill(currentPhoton.ecalRecHitSumEtConeDR04());
234  h_isoHcalRecHit_->Fill(currentPhoton.hcalTowerSumEtConeDR04());
235  h_trk_pt_solid_->Fill(currentPhoton.trkSumPtSolidConeDR04());
236  h_trk_pt_hollow_->Fill(currentPhoton.trkSumPtHollowConeDR04());
237  h_ntrk_solid_->Fill(currentPhoton.nTrkSolidConeDR04());
238  h_ntrk_hollow_->Fill(currentPhoton.nTrkHollowConeDR04());
239  h_ebgap_->Fill(currentPhoton.isEBGap());
240  h_eeGap_->Fill(currentPhoton.isEEGap());
241  h_ebeeGap_->Fill(currentPhoton.isEBEEGap());
242  h_r9_->Fill(currentPhoton.r9());
243 
244  // Photon Variables
245  h_photonEt_->Fill(photonEt);
246  h_photonEta_->Fill(currentPhoton.eta());
247  h_photonPhi_->Fill(currentPhoton.phi());
248  h_hadoverem_->Fill(currentPhoton.hadronicOverEm());
249 
250  // Photon's SuperCluster Variables
251  // eta is with respect to detector (not physics) vertex,
252  // thus Et and eta are different from photon.
253  h_photonScEt_->Fill(superClusterEt);
254  h_photonScEta_->Fill(currentPhoton.superCluster()->position().eta());
255  h_photonScPhi_->Fill(currentPhoton.superCluster()->position().phi());
256  h_photonScEtaWidth_->Fill(currentPhoton.superCluster()->etaWidth());
257 
258  // It passed photon cuts, mark it
259  h_nPassingPho_->Fill(1.0);
260 
262  // fill TTree (optional) //
264  if (createPhotonTTree_) {
265  recPhoton.isolationEcalRecHit = currentPhoton.ecalRecHitSumEtConeDR04();
266  recPhoton.isolationHcalRecHit = currentPhoton.hcalTowerSumEtConeDR04();
267  recPhoton.isolationSolidTrkCone = currentPhoton.trkSumPtSolidConeDR04();
268  recPhoton.isolationHollowTrkCone = currentPhoton.trkSumPtHollowConeDR04();
269  recPhoton.nTrkSolidCone = currentPhoton.nTrkSolidConeDR04();
270  recPhoton.nTrkHollowCone = currentPhoton.nTrkHollowConeDR04();
271  recPhoton.isEBGap = currentPhoton.isEBGap();
272  recPhoton.isEEGap = currentPhoton.isEEGap();
273  recPhoton.isEBEEGap = currentPhoton.isEBEEGap();
274  recPhoton.r9 = currentPhoton.r9();
275  recPhoton.et = currentPhoton.et();
276  recPhoton.eta = currentPhoton.eta();
277  recPhoton.phi = currentPhoton.phi();
278  recPhoton.hadronicOverEm = currentPhoton.hadronicOverEm();
279  recPhoton.ecalIso = currentPhoton.ecalIso();
280  recPhoton.hcalIso = currentPhoton.hcalIso();
281  recPhoton.trackIso = currentPhoton.trackIso();
282 
283  // Fill the tree (this records all the recPhoton.* since
284  // tree_PhotonAll_ was set to point at that.
285  tree_PhotonAll_->Fill();
286  }
287 
288  // Record whether it was near any module gap.
289  // Very convoluted at the moment.
290  bool inAnyGap = currentPhoton.isEBEEGap() || (currentPhoton.isEB() && currentPhoton.isEBGap()) ||
291  (currentPhoton.isEE() && currentPhoton.isEEGap());
292  if (inAnyGap) {
293  h_photonInAnyGap_->Fill(1.0);
294  } else {
295  h_photonInAnyGap_->Fill(0.0);
296  }
297 
298  photonCounter++;
299  } else {
300  // This didn't pass photon cuts, mark it
301  h_nPassingPho_->Fill(0.0);
302  }
303 
304  } // End Loop over photons
305  h_nPho_->Fill(photonCounter);
306 }
307 
309 // method called once each job just after ending the event loop //
312  // go to *OUR* root file and store histograms
313  rootFile_->cd();
314 
315  // PhotonID Histograms
316  h_isoEcalRecHit_->Write();
317  h_isoHcalRecHit_->Write();
318  h_trk_pt_solid_->Write();
319  h_trk_pt_hollow_->Write();
320  h_ntrk_solid_->Write();
321  h_ntrk_hollow_->Write();
322  h_ebgap_->Write();
323  h_eeGap_->Write();
324  h_ebeeGap_->Write();
325  h_r9_->Write();
326 
327  // Photon Histograms
328  h_photonEt_->Write();
329  h_photonEta_->Write();
330  h_photonPhi_->Write();
331  h_hadoverem_->Write();
332 
333  // Photon's SuperCluster Histograms
334  h_photonScEt_->Write();
335  h_photonScEta_->Write();
336  h_photonScPhi_->Write();
337  h_photonScEtaWidth_->Write();
338 
339  // Composite or Other Histograms
340  h_photonInAnyGap_->Write();
341  h_nPassingPho_->Write();
342  h_nPho_->Write();
343 
344  // Write the root file (really writes the TTree)
345  rootFile_->Write();
346  rootFile_->Close();
347 }
348 
349 //define this as a plug-in
350 // DEFINE_FWK_MODULE(PatPhotonSimpleAnalyzer);
Analysis-level Photon class.
Definition: Photon.h:46
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
float ecalRecHitSumEtConeDR04() const
Definition: Photon.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
PatPhotonSimpleAnalyzer(const edm::ParameterSet &)
bool isEB() const
Definition: Photon.h:123
float trackIso() const
Definition: Photon.h:111
float trkSumPtSolidConeDR04() const
Definition: Photon.h:495
float hcalTowerSumEtConeDR04(int depth=0) const
Definition: Photon.h:475
void analyze(const edm::Event &, const edm::EventSetup &) override
bool isEE() const
Definition: Photon.h:125
bool isEEGap() const
true if photon is in EE, and inside the boundaries in supercrystal/D
Definition: Photon.h:131
#define M_PI
float ecalIso() const
Definition: Photon.h:114
int nTrkSolidConeDR04() const
Definition: Photon.h:499
float trkSumPtHollowConeDR04() const
Definition: Photon.h:497
float hadronicOverEm(int depth=0) const
Definition: Photon.h:236
float hcalIso() const
Definition: Photon.h:117
reco::SuperClusterRef superCluster() const override
override the superCluster method from CaloJet, to access the internal storage of the supercluster ...
float r9() const
Definition: Photon.h:276
bool isEBGap() const
true if photon is in EB, and inside the boundaries in super crystals/modules
Definition: Photon.h:127
double et() const final
transverse energy
HLT enums.
int nTrkHollowConeDR04() const
Definition: Photon.h:501
double phi() const final
momentum azimuthal angle
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
bool isEBEEGap() const
true if photon is in boundary between EB and EE
Definition: Photon.h:135
double eta() const final
momentum pseudorapidity