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