CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoEgamma/Examples/plugins/PatPhotonSimpleAnalyzer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    PatPhotonSimpleAnalyzer
00004 // Class:      PatPhotonSimpleAnalyzer
00005 //
00014 //
00015 // Original Author:  M.B. Anderson
00016 //   based on simple photon analyzer by:  J. Stilley, A. Askew
00017 //
00018 //         Created:  Wed Sep 23 12:00:01 CDT 2008
00019 //
00021 //                    header file for this analyzer                  //
00023 #include "RecoEgamma/Examples/plugins/PatPhotonSimpleAnalyzer.h"
00024 
00026 //                        CMSSW includes                             //
00028 #include "DataFormats/Common/interface/View.h"
00029 #include "DataFormats/PatCandidates/interface/Photon.h"
00030 
00032 //                      Root include files                           //
00034 #include "TH1.h"
00035 #include "TFile.h"
00036 #include "TMath.h"
00037 #include "TTree.h"
00038 
00039 using namespace std;
00040 
00042 //                           Constructor                             //
00044 PatPhotonSimpleAnalyzer::PatPhotonSimpleAnalyzer(const edm::ParameterSet& ps)
00045 {
00046   // Read Parameters from configuration file
00047 
00048   // output filename
00049   outputFile_   = ps.getParameter<std::string>("outputFile");
00050   // Read variables that must be passed to allow a
00051   //  supercluster to be placed in histograms as a photon.
00052   minPhotonEt_     = ps.getParameter<double>("minPhotonEt");
00053   minPhotonAbsEta_ = ps.getParameter<double>("minPhotonAbsEta");
00054   maxPhotonAbsEta_ = ps.getParameter<double>("maxPhotonAbsEta");
00055   minPhotonR9_     = ps.getParameter<double>("minPhotonR9");
00056   maxPhotonHoverE_ = ps.getParameter<double>("maxPhotonHoverE");
00057 
00058   // Read variable to that decidedes whether
00059   // a TTree of photons is created or not
00060   createPhotonTTree_ = ps.getParameter<bool>("createPhotonTTree");
00061 
00062   // open output file to store histograms
00063   rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
00064 }
00065 
00067 //                            Destructor                             //
00069 PatPhotonSimpleAnalyzer::~PatPhotonSimpleAnalyzer()
00070 {
00071 
00072 // do anything here that needs to be done at desctruction time
00073 // (e.g. close files, deallocate resources etc.)
00074 
00075   delete rootFile_;
00076 
00077 }
00078 
00080 //    method called once each job just before starting event loop    //
00082 void
00083 PatPhotonSimpleAnalyzer::beginJob()
00084 {
00085 
00086   // go to *OUR* rootfile
00087   rootFile_->cd();
00088 
00089   // Book Histograms
00090   // PhotonID Histograms
00091   h_isoEcalRecHit_ = new TH1F("photonEcalIso",          "Ecal Rec Hit Isolation", 100, 0, 100);
00092   h_isoHcalRecHit_ = new TH1F("photonHcalIso",          "Hcal Rec Hit Isolation", 100, 0, 100);
00093   h_trk_pt_solid_  = new TH1F("photonTrackSolidIso",    "Sum of track pT in a cone of #DeltaR" , 100, 0, 100);
00094   h_trk_pt_hollow_ = new TH1F("photonTrackHollowIso",   "Sum of track pT in a hollow cone" ,     100, 0, 100);
00095   h_ntrk_solid_    = new TH1F("photonTrackCountSolid",  "Number of tracks in a cone of #DeltaR", 100, 0, 100);
00096   h_ntrk_hollow_   = new TH1F("photonTrackCountHollow", "Number of tracks in a hollow cone",     100, 0, 100);
00097   h_ebgap_         = new TH1F("photonInEBgap",          "Ecal Barrel gap flag",  2, -0.5, 1.5);
00098   h_eeGap_         = new TH1F("photonInEEgap",          "Ecal Endcap gap flag",  2, -0.5, 1.5);
00099   h_ebeeGap_       = new TH1F("photonInEEgap",          "Ecal Barrel/Endcap gap flag",  2, -0.5, 1.5);
00100   h_r9_            = new TH1F("photonR9",               "R9 = E(3x3) / E(SuperCluster)", 300, 0, 3);
00101 
00102   // Photon Histograms
00103   h_photonEt_      = new TH1F("photonEt",     "Photon E_{T}",  200,  0, 200);
00104   h_photonEta_     = new TH1F("photonEta",    "Photon #eta",   200, -4,   4);
00105   h_photonPhi_     = new TH1F("photonPhi",    "Photon #phi",   200, -1.*TMath::Pi(), TMath::Pi());
00106   h_hadoverem_     = new TH1F("photonHoverE", "Hadronic over EM", 200, 0, 1);
00107 
00108   // Photon's SuperCluster Histograms
00109   h_photonScEt_       = new TH1F("photonScEt",  "Photon SuperCluster E_{T}", 200,  0, 200);
00110   h_photonScEta_      = new TH1F("photonScEta", "Photon #eta",               200, -4,   4);
00111   h_photonScPhi_      = new TH1F("photonScPhi", "Photon #phi", 200, -1.*TMath::Pi(), TMath::Pi());
00112   h_photonScEtaWidth_ = new TH1F("photonScEtaWidth","#eta-width",            100,  0,  .1);
00113 
00114   // Composite or Other Histograms
00115   h_photonInAnyGap_   = new TH1F("photonInAnyGap",     "Photon in any gap flag",  2, -0.5, 1.5);
00116   h_nPassingPho_      = new TH1F("photonPassingCount", "Total number photons (0=NotPassing, 1=Passing)", 2, -0.5, 1.5);
00117   h_nPho_             = new TH1F("photonCount",        "Number of photons passing cuts in event",  10,  0,  10);
00118 
00119   // Create a TTree of photons if set to 'True' in config file
00120   if ( createPhotonTTree_ ) {
00121     tree_PhotonAll_     = new TTree("TreePhotonAll", "Reconstructed Photon");
00122     tree_PhotonAll_->Branch("recPhoton", &recPhoton.isolationEcalRecHit, "isolationEcalRecHit/F:isolationHcalRecHit:isolationSolidTrkCone:isolationHollowTrkCone:nTrkSolidCone:nTrkHollowCone:isEBGap:isEEGap:isEBEEGap:r9:et:eta:phi:hadronicOverEm:ecalIso:hcalIso:trackIso");
00123   }
00124 }
00125 
00127 //                method called to for each event                    //
00129 void
00130 PatPhotonSimpleAnalyzer::analyze(const edm::Event& evt, const edm::EventSetup& es)
00131 {
00132 
00133   using namespace std;
00134   using namespace edm;
00135 
00136   // Grab pat::Photon
00137   Handle< View<pat::Photon> >  photonHandle;
00138   evt.getByLabel("selectedLayer1Photons", photonHandle);
00139   View<pat::Photon> photons = *photonHandle;
00140 
00141   int photonCounter = 0;
00142 
00143   for (int i=0; i<int(photons.size()); i++)
00144   {
00145 
00146     pat::Photon currentPhoton = photons.at(i);
00147 
00148     float photonEt       = currentPhoton.et();
00149     float superClusterEt = (currentPhoton.superCluster()->energy())/(cosh(currentPhoton.superCluster()->position().eta()));
00150 
00151     // Only store photon candidates (SuperClusters) that pass some simple cuts
00152     bool passCuts = (              photonEt > minPhotonEt_     ) &&
00153                     (      fabs(currentPhoton.eta()) > minPhotonAbsEta_ ) &&
00154                     (      fabs(currentPhoton.eta()) < maxPhotonAbsEta_ ) &&
00155                     (          currentPhoton.r9() > minPhotonR9_     ) &&
00156                     ( currentPhoton.hadronicOverEm() < maxPhotonHoverE_ ) ;
00157 
00158     if ( passCuts )
00159     {
00161       //                fill histograms                    //
00163       // PhotonID Variables
00164       h_isoEcalRecHit_->Fill(currentPhoton.ecalRecHitSumEtConeDR04());
00165       h_isoHcalRecHit_->Fill(currentPhoton.hcalTowerSumEtConeDR04());
00166       h_trk_pt_solid_ ->Fill(currentPhoton.trkSumPtSolidConeDR04());
00167       h_trk_pt_hollow_->Fill(currentPhoton.trkSumPtHollowConeDR04());
00168       h_ntrk_solid_->   Fill(currentPhoton.nTrkSolidConeDR04());
00169       h_ntrk_hollow_->  Fill(currentPhoton.nTrkHollowConeDR04());
00170       h_ebgap_->        Fill(currentPhoton.isEBGap());
00171       h_eeGap_->        Fill(currentPhoton.isEEGap());
00172       h_ebeeGap_->      Fill(currentPhoton.isEBEEGap());
00173       h_r9_->           Fill(currentPhoton.r9());
00174 
00175       // Photon Variables
00176       h_photonEt_->  Fill(photonEt);
00177       h_photonEta_-> Fill(currentPhoton.eta());
00178       h_photonPhi_-> Fill(currentPhoton.phi());
00179       h_hadoverem_-> Fill(currentPhoton.hadronicOverEm());
00180 
00181       // Photon's SuperCluster Variables
00182       // eta is with respect to detector (not physics) vertex,
00183       // thus Et and eta are different from photon.
00184       h_photonScEt_->      Fill(superClusterEt);
00185       h_photonScEta_->     Fill(currentPhoton.superCluster()->position().eta());
00186       h_photonScPhi_->     Fill(currentPhoton.superCluster()->position().phi());
00187       h_photonScEtaWidth_->Fill(currentPhoton.superCluster()->etaWidth());
00188 
00189       // It passed photon cuts, mark it
00190       h_nPassingPho_->Fill(1.0);
00191 
00193       //                fill TTree (optional)              //
00195       if ( createPhotonTTree_ ) {
00196         recPhoton.isolationEcalRecHit    = currentPhoton.ecalRecHitSumEtConeDR04();
00197         recPhoton.isolationHcalRecHit    = currentPhoton.hcalTowerSumEtConeDR04();
00198         recPhoton.isolationSolidTrkCone  = currentPhoton.trkSumPtSolidConeDR04();
00199         recPhoton.isolationHollowTrkCone = currentPhoton.trkSumPtHollowConeDR04();
00200         recPhoton.nTrkSolidCone          = currentPhoton.nTrkSolidConeDR04();
00201         recPhoton.nTrkHollowCone         = currentPhoton.nTrkHollowConeDR04();
00202         recPhoton.isEBGap                = currentPhoton.isEBGap();
00203         recPhoton.isEEGap                = currentPhoton.isEEGap();
00204         recPhoton.isEBEEGap              = currentPhoton.isEBEEGap();
00205         recPhoton.r9                     = currentPhoton.r9();
00206         recPhoton.et                     = currentPhoton.et();
00207         recPhoton.eta                    = currentPhoton.eta();
00208         recPhoton.phi                    = currentPhoton.phi();
00209         recPhoton.hadronicOverEm         = currentPhoton.hadronicOverEm();
00210         recPhoton.ecalIso                = currentPhoton.ecalIso();
00211         recPhoton.hcalIso                = currentPhoton.hcalIso();
00212         recPhoton.trackIso               = currentPhoton.trackIso();
00213 
00214         // Fill the tree (this records all the recPhoton.* since
00215         // tree_PhotonAll_ was set to point at that.
00216         tree_PhotonAll_->Fill();
00217       }
00218 
00219       // Record whether it was near any module gap.
00220       // Very convoluted at the moment.
00221       bool inAnyGap = currentPhoton.isEBEEGap() || (currentPhoton.isEB()&&currentPhoton.isEBGap()) || (currentPhoton.isEE()&&currentPhoton.isEEGap());
00222       if (inAnyGap) {
00223         h_photonInAnyGap_->Fill(1.0);
00224       } else {
00225         h_photonInAnyGap_->Fill(0.0);
00226       }
00227 
00228       photonCounter++;
00229     }
00230     else
00231     {
00232       // This didn't pass photon cuts, mark it
00233       h_nPassingPho_->Fill(0.0);
00234     }
00235 
00236   } // End Loop over photons
00237   h_nPho_->Fill(photonCounter);
00238 
00239 }
00240 
00242 //    method called once each job just after ending the event loop   //
00244 void
00245 PatPhotonSimpleAnalyzer::endJob()
00246 {
00247 
00248   // go to *OUR* root file and store histograms
00249   rootFile_->cd();
00250 
00251   // PhotonID Histograms
00252   h_isoEcalRecHit_->Write();
00253   h_isoHcalRecHit_->Write();
00254   h_trk_pt_solid_-> Write();
00255   h_trk_pt_hollow_->Write();
00256   h_ntrk_solid_->   Write();
00257   h_ntrk_hollow_->  Write();
00258   h_ebgap_->     Write();
00259   h_eeGap_->     Write();
00260   h_ebeeGap_->   Write();
00261   h_r9_->        Write();
00262 
00263   // Photon Histograms
00264   h_photonEt_->  Write();
00265   h_photonEta_-> Write();
00266   h_photonPhi_-> Write();
00267   h_hadoverem_-> Write();
00268 
00269   // Photon's SuperCluster Histograms
00270   h_photonScEt_->      Write();
00271   h_photonScEta_->     Write();
00272   h_photonScPhi_->     Write();
00273   h_photonScEtaWidth_->Write();
00274 
00275   // Composite or Other Histograms
00276   h_photonInAnyGap_->Write();
00277   h_nPassingPho_->   Write();
00278   h_nPho_->          Write();
00279 
00280   // Write the root file (really writes the TTree)
00281   rootFile_->Write();
00282   rootFile_->Close();
00283 
00284 }
00285 
00286 //define this as a plug-in
00287 // DEFINE_FWK_MODULE(PatPhotonSimpleAnalyzer);