CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/EgammaAnalysis/CSA07Skims/src/EgammaProbeSelector.cc

Go to the documentation of this file.
00001 
00010 #include "EgammaAnalysis/CSA07Skims/interface/EgammaProbeSelector.h"
00011 #include "DataFormats/JetReco/interface/Jet.h"
00012 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00013 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00014 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 #include <iostream> 
00017 #include <TMath.h>
00018 
00019 using namespace edm;
00020 using namespace std;
00021 using namespace reco;
00022 
00023 
00024 EgammaProbeSelector::EgammaProbeSelector(const edm::ParameterSet& iConfig) {
00025   
00026   jetLabel        = iConfig.getParameter<InputTag>("JetCollection");
00027   minNumberOfjets = iConfig.getParameter<int>("MinNumberOfJets");
00028   jetEtMin        = iConfig.getParameter<double>("JetEtMin");
00029   jetEtaMin       = iConfig.getParameter<double>("JetEtaMin");
00030   jetEtaMax       = iConfig.getParameter<double>("JetEtaMax");
00031   
00032   nEvents         = 0;
00033   nSelectedEvents = 0;
00034   
00035   scLabel          = iConfig.getParameter<InputTag>("SuperClusterBarrelCollection");
00036   scEELabel        = iConfig.getParameter<InputTag>("SuperClusterEndCapCollection");
00037   minNumberOfSuperClusters = iConfig.getParameter<int>("MinNumberOfSuperClusters");
00038   scEtMin          = iConfig.getParameter<double>("ScEtMin");
00039   scEtaMin         = iConfig.getParameter<double>("ScEtaMin");
00040   scEtaMax         = iConfig.getParameter<double>("ScEtaMax");
00041 }
00042 
00043 
00044 EgammaProbeSelector::~EgammaProbeSelector(){
00045   edm::LogVerbatim("EgammaProbeSelector") 
00046     << " Number_events_read " << nEvents
00047     << " Number_events_kept " << nSelectedEvents
00048     << " Efficiency         " << ((double)nSelectedEvents)/((double) nEvents + 0.01) << std::endl;
00049   
00050 }
00051 
00052 
00053 bool EgammaProbeSelector::filter(edm::Event& iEvent, const edm::EventSetup& iSetup ){
00054   int nSC   = 0;
00055   int nJets = 0;  
00056   nEvents++;
00057   
00058   // Get reconstructed Jets, check their Et and eta, and count their number
00059   Handle<CaloJetCollection> jetHandle;
00060   iEvent.getByLabel(jetLabel,jetHandle);
00061   
00062   if(jetHandle.isValid()){
00063     const reco::CaloJetCollection & jets = *(jetHandle.product());
00064     CaloJetCollection::const_iterator ijet;
00065     for(ijet = jets.begin(); ijet!= jets.end(); ijet++){
00066       if(ijet->pt()  > jetEtMin  && 
00067          ijet->eta() > jetEtaMin &&
00068          ijet->eta() < jetEtaMax ) 
00069         nJets++;                
00070     }
00071   }
00072   
00073   // Get Super Clusters, check their Et and eta, and count their number
00074   // EB
00075   Handle<reco::SuperClusterCollection> scHandle;
00076   iEvent.getByLabel(scLabel,scHandle);
00077   
00078   if(scHandle.isValid()){
00079     const reco::SuperClusterCollection & SCs = *(scHandle.product());
00080     
00081     for(reco::SuperClusterCollection::const_iterator scIt = SCs.begin(); 
00082         scIt!= SCs.end(); scIt++){
00083 
00084       double scEt   = scIt->energy() * TMath::Sin(scIt->position().theta());
00085       double _eta = scIt->position().eta();
00086       if(scEt  > scEtMin && _eta > scEtaMin && _eta < scEtaMax) nSC++;          
00087     }
00088   }
00089 
00090   // EE
00091   Handle<reco::SuperClusterCollection> scEEHandle;
00092   iEvent.getByLabel(scEELabel,scEEHandle);
00093   
00094   if(scEEHandle.isValid()){
00095     const reco::SuperClusterCollection & SCs = *(scEEHandle.product());
00096     for(reco::SuperClusterCollection::const_iterator scIt = SCs.begin(); 
00097         scIt!= SCs.end(); scIt++){
00098       
00099       double scEt   = scIt->energy() * TMath::Sin(scIt->position().theta()); 
00100       double _eta = scIt->position().eta();
00101       if(scEt  > scEtMin && _eta > scEtaMin && _eta < scEtaMax) nSC++; 
00102     }
00103   } 
00104   
00105   // Make final desicion on passed events
00106   bool accepted = false;
00107   if(nJets >= minNumberOfjets || nSC >=minNumberOfSuperClusters) {
00108     accepted = true;
00109     nSelectedEvents++;
00110   }
00111   
00112   return accepted;
00113 }