00001
00020 #include <iostream>
00021 #include <string>
00022 #include <list>
00023 #include <cmath>
00024 #include <cstdio>
00025 #include <vector>
00026 #include <memory>
00027
00028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00029 #include "DataFormats/Common/interface/Handle.h"
00030
00031 #include "SUSYBSMAnalysis/CSA07Skims/interface/LepSUSYSkim.h"
00032 #include "DataFormats/MuonReco/interface/Muon.h"
00033 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00034 #include "DataFormats/JetReco/interface/CaloJet.h"
00035 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00036 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00037 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00038 #include "DataFormats/METReco/interface/CaloMET.h"
00039
00040 using namespace edm;
00041 using namespace std;
00042 using namespace reco;
00043
00044 class PtSorter {
00045 public:
00046 template <class T> bool operator() ( const T& a, const T& b ) {
00047 return ( a.pt() > b.pt() );
00048 }
00049 };
00050
00051 LepSUSYSkim::LepSUSYSkim( const edm::ParameterSet& iConfig ) :
00052 nEvents_(0), nAccepted_(0)
00053 {
00054 Muonsrc_ = iConfig.getParameter<InputTag>( "Muonsrc" );
00055 Elecsrc_ = iConfig.getParameter<InputTag>( "Elecsrc" );
00056 CaloJetsrc_ = iConfig.getParameter<InputTag>( "CaloJetsrc" );
00057 CaloMETsrc_ = iConfig.getParameter<InputTag>( "CaloMETsrc" );
00058 NminMuon_ = iConfig.getParameter<int>( "NminMuon");
00059 MuonPtmin_ = iConfig.getParameter<double>( "MuonPtmin");
00060 NminElec_ = iConfig.getParameter<int>( "NminElec");
00061 ElecPtmin_ = iConfig.getParameter<double>( "ElecPtmin");
00062 NminCaloJet_ = iConfig.getParameter<int>( "NminCaloJet");
00063 CaloJet1Ptmin_ = iConfig.getParameter<double>( "CaloJet1Ptmin");
00064 CaloJet2Ptmin_ = iConfig.getParameter<double>( "CaloJet2Ptmin");
00065 CaloMetmin_ = iConfig.getParameter<double>( "CaloMetmin" );
00066
00067 }
00068
00069
00070
00071 LepSUSYSkim::~LepSUSYSkim()
00072 {}
00073
00074
00075
00076 bool LepSUSYSkim::filter( edm::Event& iEvent,
00077 const edm::EventSetup& iSetup )
00078 {
00079 nEvents_++;
00080
00081
00082
00083 Handle<MuonCollection> MuonHandle;
00084
00085 iEvent.getByLabel( Muonsrc_, MuonHandle );
00086
00087
00088
00089
00090
00091
00092
00093 MuonCollection TheMuons = *MuonHandle;
00094 std::stable_sort( TheMuons.begin(), TheMuons.end(), PtSorter() );
00095
00096
00097 Handle<GsfElectronCollection> ElecHandle;
00098
00099 iEvent.getByLabel( Elecsrc_, ElecHandle );
00100
00101
00102
00103
00104
00105
00106
00107 GsfElectronCollection TheElecs = *ElecHandle;
00108 std::stable_sort( TheElecs.begin(), TheElecs.end(), PtSorter() );
00109
00110
00111
00112 Handle<CaloJetCollection> CaloJetsHandle;
00113
00114 iEvent.getByLabel( CaloJetsrc_, CaloJetsHandle );
00115
00116
00117
00118
00119
00120
00121
00122 CaloJetCollection TheCaloJets = *CaloJetsHandle;
00123 std::stable_sort( TheCaloJets.begin(), TheCaloJets.end(), PtSorter() );
00124
00125
00126 Handle<CaloMETCollection> METHandle;
00127
00128 iEvent.getByLabel( CaloMETsrc_, METHandle );
00129
00130
00131
00132
00133
00134
00135
00136 if ( METHandle->empty() ) return false;
00137
00138
00139
00140
00141
00142 int nMuon = 0;
00143 for ( MuonCollection::const_iterator it = TheMuons.begin();
00144 it != TheMuons.end(); it++ ) {
00145
00146
00147 if ( (it->pt() > MuonPtmin_) &&
00148 (fabs(it->eta()) < 3.0) ) {
00149 nMuon++;
00150 }
00151 }
00152
00153 if ( nMuon < NminMuon_ ) return false;
00154
00155
00156 int nElec = 0;
00157 for ( GsfElectronCollection::const_iterator it = TheElecs.begin();
00158 it != TheElecs.end(); it++ ) {
00159
00160
00161 if ( (it->pt() > ElecPtmin_) &&
00162 (fabs(it->eta()) < 3.0) ) {
00163 nElec++;
00164 }
00165 }
00166
00167 if ( nElec < NminElec_ ) return false;
00168
00169
00170
00171 int nJet = 0;
00172 for ( CaloJetCollection::const_iterator it = TheCaloJets.begin();
00173 it != TheCaloJets.end(); it++ ) {
00174
00175
00176 if (fabs(it->eta()) < 3.0) nJet++;
00177 }
00178
00179 if ( nJet < NminCaloJet_ ) return false;
00180
00181
00182 if(NminCaloJet_ > 0) {if ( TheCaloJets[0].pt() < CaloJet1Ptmin_ ) return false;}
00183 if(NminCaloJet_ > 1) {if ( TheCaloJets[1].pt() < CaloJet2Ptmin_ ) return false;}
00184
00185
00186
00187
00188
00189
00190 double MetValue = METHandle->begin()->pt();
00191
00192 if(MetValue < CaloMetmin_ ) return false;
00193
00194 nAccepted_++;
00195
00196 return true;
00197 }
00198
00199
00200
00201 void LepSUSYSkim::endJob()
00202 {
00203 edm::LogVerbatim( "LepSUSYSkim" )
00204 << "Events read " << nEvents_
00205 << "\nEvents accepted " << nAccepted_
00206 << "\nEfficiency " << (double)(nAccepted_)/(double)(nEvents_)
00207 << endl;
00208 }