CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SimplePi0DiscAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoEgamma/Examples
4 // Class: SimplePi0DiscAnalyzer
5 //
13 //
14 // Original Author: Aristotelis Kyriakis NCSR "Demokritos" Athens
15 // D Maletic, "Vinca" Belgrade
16 // Created: May 26 13:22:06 CEST 2009
17 //
18 //
19 
31 
32 #include "CLHEP/Units/PhysicalConstants.h"
33 
34 #include "TFile.h"
35 #include "TH1F.h"
36 #include "TH1I.h"
37 #include "TH2F.h"
38 #include "TProfile.h"
39 #include "TTree.h"
40 
41 #include <iostream>
42 
44 public:
45  explicit SimplePi0DiscAnalyzer(const edm::ParameterSet& conf);
46 
47  ~SimplePi0DiscAnalyzer() override;
48 
49  void beginJob() override;
50  void endJob() override;
51  void analyze(const edm::Event& e, const edm::EventSetup& c) override;
52 
53 private:
54  // ----------member data ---------------------------
55 
58 
60  TFile* rootFile_;
61 
63 
76 };
77 
80 
81 using namespace reco;
82 
84  outputFile_ = conf.getParameter<std::string>("outputFile");
85  rootFile_ = new TFile(outputFile_.c_str(), "RECREATE");
86 
87  photonCollectionProducer_ = conf.getParameter<std::string>("phoProducer");
88  photonCollection_ = conf.getParameter<std::string>("photonCollection");
89 }
90 
92  // do anything here that needs to be done at desctruction time
93  // (e.g. close files, deallocate resources etc.)
94  rootFile_->Write();
95  rootFile_->Close();
96 }
97 
99  rootFile_->cd();
100  std::cout << "beginJob() -> Book the Histograms" << std::endl;
101 
102  hConv_ntracks_ = new TH1F("nConvTracks", "Number of tracks of converted Photons ", 10, 0., 10);
103  hAll_nnout_Assoc_ = new TH1F("All_nnout_Assoc", "NNout for All Photons(AssociationMap)", 100, 0., 1.);
104  hAll_nnout_NoConv_Assoc_ =
105  new TH1F("All_nnout_NoConv_Assoc", "NNout for Unconverted Photons(AssociationMap)", 100, 0., 1.);
106  hAll_nnout_NoConv_Assoc_R9_ =
107  new TH1F("All_nnout_NoConv_Assoc_R9", "NNout for Unconverted Photons with R9>0.93 (AssociationMap)", 100, 0., 1.);
108  hBarrel_nnout_Assoc_ = new TH1F("barrel_nnout_Assoc", "NNout for Barrel Photons(AssociationMap)", 100, 0., 1.);
109  hBarrel_nnout_NoConv_Assoc_ =
110  new TH1F("barrel_nnout_NoConv_Assoc", "NNout for Barrel Unconverted Photons(AssociationMap)", 100, 0., 1.);
111  hBarrel_nnout_NoConv_Assoc_R9_ = new TH1F(
112  "barrel_nnout_NoConv_Assoc_R9", "NNout for Barrel Unconverted Photons with R9>0.93 (AssociationMap)", 100, 0., 1.);
113  hEndcNoPresh_nnout_Assoc_ =
114  new TH1F("endcNoPresh_nnout_Assoc", "NNout for Endcap NoPresh Photons(AssociationMap)", 100, 0., 1.);
115  hEndcNoPresh_nnout_NoConv_Assoc_ = new TH1F(
116  "endcNoPresh_nnout_NoConv_Assoc", "NNout for Endcap Unconverted NoPresh Photons(AssociationMap)", 100, 0., 1.);
117  hEndcNoPresh_nnout_NoConv_Assoc_R9_ =
118  new TH1F("endcNoPresh_nnout_NoConv_Assoc_R9",
119  "NNout for Endcap Unconverted NoPresh Photons with R9>0.93 (AssociationMap)",
120  100,
121  0.,
122  1.);
123  hEndcWithPresh_nnout_Assoc_ =
124  new TH1F("endcWithPresh_nnout_Assoc", "NNout for Endcap WithPresh Photons(AssociationMap)", 100, 0., 1.);
125  hEndcWithPresh_nnout_NoConv_Assoc_ = new TH1F(
126  "endcWithPresh_nnout_NoConv_Assoc", "NNout for Endcap Unconverted WithPresh Photons(AssociationMap)", 100, 0., 1.);
127  hEndcWithPresh_nnout_NoConv_Assoc_R9_ =
128  new TH1F("endcWithPresh_nnout_NoConv_Assoc_R9",
129  "NNout for Endcap Unconverted WithPresh Photons with R9>0.93 (AssociationMap)",
130  100,
131  0.,
132  1.);
133 }
134 
136  rootFile_->cd();
137  std::cout << "endJob() -> Write the Histograms" << std::endl;
138  hConv_ntracks_->Write();
139 
140  hAll_nnout_Assoc_->Write();
141  hAll_nnout_NoConv_Assoc_->Write();
142  hAll_nnout_NoConv_Assoc_R9_->Write();
143  hBarrel_nnout_Assoc_->Write();
144  hBarrel_nnout_NoConv_Assoc_->Write();
145  hBarrel_nnout_NoConv_Assoc_R9_->Write();
146  hEndcNoPresh_nnout_Assoc_->Write();
147  hEndcNoPresh_nnout_NoConv_Assoc_->Write();
148  hEndcNoPresh_nnout_NoConv_Assoc_R9_->Write();
149  hEndcWithPresh_nnout_Assoc_->Write();
150  hEndcWithPresh_nnout_NoConv_Assoc_->Write();
151  hEndcWithPresh_nnout_NoConv_Assoc_R9_->Write();
152 }
153 
155  std::cout << std::endl;
156  std::cout << " -------------- NEW EVENT : Run, Event = " << iEvent.id() << std::endl;
157 
159  iEvent.getByLabel(photonCollectionProducer_, photonCollection_, PhotonHandle);
160  const reco::PhotonCollection photons = *(PhotonHandle.product());
161 
162  std::cout << "----> Photons size: " << photons.size() << std::endl;
163 
165  iEvent.getByLabel("piZeroDiscriminators", "PhotonPi0DiscriminatorAssociationMap", map);
167 
168  // int PhoInd = 0;
169 
170  for (reco::PhotonCollection::const_iterator iPho = photons.begin(); iPho != photons.end();
171  iPho++) { // Loop over Photons
172 
173  reco::Photon localPho(*iPho);
174 
175  float Photon_et = localPho.et();
176  float Photon_eta = localPho.eta();
177  float Photon_phi = localPho.phi();
178  float Photon_r9 = localPho.r9();
179  bool isPhotConv = localPho.hasConversionTracks();
180  // std::cout << "Photon Id = " << PhoInd
181  std::cout << "Photon Id = " << iPho - photons.begin() << " with Et = " << Photon_et << " Eta = " << Photon_eta
182  << " Phi = " << Photon_phi << " R9 = " << Photon_r9 << " and conv_id = " << isPhotConv << std::endl;
183 
184  auto it_super = localPho.superCluster(); // get the SC related to the Photon candidate
185 
186  // hConv_ntracks_->Fill(Ntrk_conv);
187 
188  float nn = -10;
189  // mapIter = map->find(edm::Ref<reco::PhotonCollection>(PhotonHandle,PhoInd));
190  mapIter = map->find(edm::Ref<reco::PhotonCollection>(PhotonHandle, iPho - photons.begin()));
191  if (mapIter != map->end()) {
192  nn = mapIter->val;
193  }
194  if (fabs(it_super->eta()) <= 1.442) {
195  hBarrel_nnout_Assoc_->Fill(nn);
196  hAll_nnout_Assoc_->Fill(nn);
197  std::cout << "AssociationMap Barrel NN = " << nn << std::endl;
198  if (!isPhotConv) {
199  hBarrel_nnout_NoConv_Assoc_->Fill(nn);
200  hAll_nnout_NoConv_Assoc_->Fill(nn);
201  }
202  if (Photon_r9 > 0.93) {
203  hBarrel_nnout_NoConv_Assoc_R9_->Fill(nn);
204  hAll_nnout_NoConv_Assoc_R9_->Fill(nn);
205  }
206  } else if ((fabs(it_super->eta()) >= 1.556 && fabs(it_super->eta()) < 1.65) || fabs(it_super->eta()) > 2.5) {
207  hEndcNoPresh_nnout_Assoc_->Fill(nn);
208  hAll_nnout_Assoc_->Fill(nn);
209  std::cout << "AssociationMap EndcNoPresh NN = " << nn << std::endl;
210  if (!isPhotConv) {
211  hEndcNoPresh_nnout_NoConv_Assoc_->Fill(nn);
212  hAll_nnout_NoConv_Assoc_->Fill(nn);
213  }
214  if (Photon_r9 > 0.93) {
215  hEndcNoPresh_nnout_NoConv_Assoc_R9_->Fill(nn);
216  hAll_nnout_NoConv_Assoc_R9_->Fill(nn);
217  }
218  } else if (fabs(it_super->eta()) >= 1.65 && fabs(it_super->eta()) <= 2.5) {
219  hEndcWithPresh_nnout_Assoc_->Fill(nn);
220  hAll_nnout_Assoc_->Fill(nn);
221  std::cout << "AssociationMap EndcWithPresh NN = " << nn << std::endl;
222  if (!isPhotConv) {
223  hEndcWithPresh_nnout_NoConv_Assoc_->Fill(nn);
224  hAll_nnout_NoConv_Assoc_->Fill(nn);
225  }
226  if (Photon_r9 > 0.93) {
227  hEndcWithPresh_nnout_NoConv_Assoc_R9_->Fill(nn);
228  hAll_nnout_NoConv_Assoc_R9_->Fill(nn);
229  }
230  }
231 
232  // PhoInd++;
233  } // End Loop over Photons
234 }
const edm::EventSetup & c
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void analyze(const edm::Event &e, const edm::EventSetup &c) override
SimplePi0DiscAnalyzer(const edm::ParameterSet &conf)
int iEvent
Definition: GenABIO.cc:224
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
bool hasConversionTracks() const
Bool flagging photons with a vector of refereces to conversions with size &gt;0.
Definition: Photon.h:70
T const * product() const
Definition: Handle.h:70
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EventID id() const
Definition: EventBase.h:59
double et() const final
transverse energy
float r9() const
Definition: Photon.h:276
int32_t *__restrict__ nn
tuple cout
Definition: gather_cfg.py:144
double phi() const final
momentum azimuthal angle
double eta() const final
momentum pseudorapidity