Go to the documentation of this file.00001
00015 #include <iostream>
00016 #include <string>
00017 #include <list>
00018 #include <cmath>
00019 #include <cstdio>
00020 #include <vector>
00021 #include <memory>
00022
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024 #include "DataFormats/Common/interface/Handle.h"
00025
00026 #include "TopQuarkAnalysis/TopSkimming/plugins/TopLeptonTauFilter.h"
00027
00028
00029 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00030 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00031
00032
00033 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00034 #include "DataFormats/MuonReco/interface/Muon.h"
00035
00036 #include "DataFormats/TauReco/interface/BaseTau.h"
00037
00038 #include "DataFormats/JetReco/interface/CaloJet.h"
00039 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00040
00041 #include "DataFormats/TrackReco/interface/Track.h"
00042
00043 class PtSorter {
00044 public:
00045 template <class T> bool operator() ( const T& a, const T& b ) {
00046 return ( a.pt() > b.pt() );
00047 }
00048 };
00049
00050 TopLeptonTauFilter::TopLeptonTauFilter( const edm::ParameterSet& iConfig ) :
00051 nEvents_(0), nAccepted_(0)
00052 {
00053
00054 ElecFilter_ = iConfig.getParameter<bool>( "ElecFilter" );
00055 MuonFilter_ = iConfig.getParameter<bool>( "MuonFilter" );
00056 TauFilter_ = iConfig.getParameter<bool>( "TauFilter" );
00057 JetFilter_ = iConfig.getParameter<bool>( "JetFilter" );
00058
00059 Elecsrc_ = iConfig.getParameter<edm::InputTag>( "Elecsrc" );
00060 Muonsrc_ = iConfig.getParameter<edm::InputTag>( "Muonsrc" );
00061 Tausrc_ = iConfig.getParameter<edm::InputTag>( "Tausrc" );
00062 CaloJetsrc_ = iConfig.getParameter<edm::InputTag>( "CaloJetsrc" );
00063
00064 NminMuon_ = iConfig.getParameter<int>( "NminMuon" );
00065 NminTau_ = iConfig.getParameter<int>( "NminTau" );
00066 NminElec_ = iConfig.getParameter<int>( "NminElec" );
00067 NminCaloJet_ = iConfig.getParameter<int>( "NminCaloJet" );
00068
00069 ElecPtmin_ = iConfig.getParameter<double>( "ElecPtmin" );
00070 MuonPtmin_ = iConfig.getParameter<double>( "MuonPtmin" );
00071 TauPtmin_ = iConfig.getParameter<double>( "TauLeadTkPtmin" );
00072 CaloJetPtmin_ = iConfig.getParameter<double>( "CaloJetPtmin" );
00073
00074 }
00075
00076
00077
00078 TopLeptonTauFilter::~TopLeptonTauFilter()
00079 {}
00080
00081
00082
00083 bool TopLeptonTauFilter::filter( edm::Event& iEvent,
00084 const edm::EventSetup& iSetup )
00085 {
00086
00087 bool filterResult(true);
00088 nEvents_++;
00089
00090 if(ElecFilter_){ filterResult = electronFilter(iEvent,iSetup); }
00091 if(MuonFilter_){ filterResult = muonFilter(iEvent,iSetup); }
00092 if(TauFilter_) { filterResult = tauFilter(iEvent,iSetup); }
00093 if(JetFilter_) { filterResult = jetFilter(iEvent,iSetup); }
00094
00095 if (filterResult) nAccepted_++;
00096
00097 return filterResult;
00098 }
00099
00100
00101 bool TopLeptonTauFilter::electronFilter( edm::Event& iEvent, const edm::EventSetup
00102 & iSetup ){
00103
00104
00105 edm::Handle<reco::GsfElectronCollection> ElecHandle;
00106 iEvent.getByLabel( Elecsrc_, ElecHandle );
00107 if ( ElecHandle->empty() && NminElec_!=0 ) return false;
00108 reco::GsfElectronCollection TheElecs = *ElecHandle;
00109 std::stable_sort( TheElecs.begin(), TheElecs.end(), PtSorter() );
00110
00111 int nElec = 0;
00112 for ( reco::GsfElectronCollection::const_iterator it = TheElecs.begin();
00113 it != TheElecs.end(); it++ ) {
00114 if ( (it->pt() > ElecPtmin_)
00115 && (fabs(it->eta()) < 3.0) ) {
00116 nElec++;
00117 }
00118 }
00119 if ( nElec < NminElec_ ) return false;
00120
00121
00122 return true;
00123 }
00124
00125 bool TopLeptonTauFilter::muonFilter( edm::Event& iEvent, const edm::EventSetup
00126 & iSetup ){
00127
00128 edm::Handle<reco::MuonCollection> MuonHandle;
00129 iEvent.getByLabel( Muonsrc_, MuonHandle );
00130 if ( MuonHandle->empty() && NminMuon_!=0 ) return false;
00131 reco::MuonCollection TheMuons = *MuonHandle;
00132 std::stable_sort( TheMuons.begin(), TheMuons.end(), PtSorter() );
00133
00134 int nMuon = 0;
00135
00136
00137 for ( reco::MuonCollection::const_iterator it = TheMuons.begin();
00138 it != TheMuons.end(); it++ ) {
00139 if ( (it->pt() > MuonPtmin_)
00140 && (fabs(it->eta()) < 3.0) ) {
00141 nMuon++;
00142 }
00143 }
00144
00145 if ( nMuon < NminMuon_ ) return false;
00146
00147 nAccepted_++;
00148
00149 return true;
00150 }
00151
00152
00153
00154 bool TopLeptonTauFilter::tauFilter( edm::Event& iEvent, const edm::EventSetup
00155 & iSetup ){
00156
00157 edm::Handle<reco::BaseTauCollection> TauHandle;
00158 iEvent.getByLabel(Tausrc_,TauHandle);
00159 const reco::BaseTauCollection& myTauCollection=*(TauHandle.product());
00160 if ( myTauCollection.empty() && NminTau_!=0 ) return false;
00161
00162 int nTau = 0;
00163 for(reco::BaseTauCollection::const_iterator it =myTauCollection.begin();it !=myTauCollection.end();it++)
00164 {
00165 reco::TrackRef theLeadTk = it->leadTrack();
00166 if(!theLeadTk) {}
00167 else{
00168 double leadTkPt = (*theLeadTk).pt();
00169 double leadTkEta = (*theLeadTk).eta();
00170
00171 if ( (leadTkPt> TauPtmin_)
00172 && (fabs(leadTkEta) < 2.4) ) {
00173 nTau++;
00174 }
00175 }
00176 }
00177 if ( nTau < NminTau_ ) return false;
00178
00179 return true;
00180 }
00181
00182
00183
00184
00185 bool TopLeptonTauFilter::jetFilter( edm::Event& iEvent, const edm::EventSetup
00186 & iSetup ){
00187
00188 edm::Handle<reco::CaloJetCollection> CaloJetsHandle;
00189 iEvent.getByLabel( CaloJetsrc_, CaloJetsHandle );
00190 if ( CaloJetsHandle->empty() && NminCaloJet_!=0 ) return false;
00191
00192 int nJet = 0;
00193 for ( reco::CaloJetCollection::const_iterator it = CaloJetsHandle->begin();
00194 it != CaloJetsHandle->end(); it++ ) {
00195 if ( (fabs(it->eta()) < 3.0) &&
00196 (it->pt() > CaloJetPtmin_) ) nJet++;
00197 }
00198
00199 if ( nJet < NminCaloJet_ ) return false;
00200
00201 return true;
00202
00203 }
00204
00205
00206
00207
00208
00209
00210 void TopLeptonTauFilter::endJob()
00211 {
00212
00213
00214 edm::LogInfo( "TopLeptonTauFilter" )
00215 << "\n Events read " << nEvents_
00216 << " Events accepted " << nAccepted_
00217 << "\nEfficiency " << (double)(nAccepted_)/(double)(nEvents_)
00218 << std::endl;
00219
00220 }
00221
00222 #include "FWCore/Framework/interface/MakerMacros.h"
00223
00224
00225 DEFINE_FWK_MODULE(TopLeptonTauFilter);