![]() |
![]() |
00001 /* \class SUSYHighPtPhotonSkim 00002 * 00003 * High Energy Photon SUSY Skim 00004 * one(two) photon(s) > xx GeV in barrel + isolation 00005 * 00006 * $Date: 2007/08/31 14:17:19 $ 00007 * $Revision: 1.5 $ 00008 * 00009 * \author Daniele del Re - Univ. La Sapienza & INFN 00010 * 00011 */ 00012 00013 #include <iostream> 00014 #include <string> 00015 #include <list> 00016 #include <cmath> 00017 #include <cstdio> 00018 #include <vector> 00019 #include <memory> 00020 00021 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00022 #include "DataFormats/Common/interface/Handle.h" 00023 00024 #include "DataFormats/EgammaCandidates/interface/Photon.h" 00025 #include "DataFormats/Candidate/interface/Candidate.h" 00026 #include "DataFormats/Common/interface/AssociationVector.h" 00027 00028 #include "SUSYBSMAnalysis/CSA07Skims/interface/SUSYHighPtPhotonSkim.h" 00029 00030 using namespace edm; 00031 using namespace std; 00032 using namespace reco; 00033 00034 SUSYHighPtPhotonSkim::SUSYHighPtPhotonSkim( const edm::ParameterSet& iConfig ) : 00035 nEvents_(0), nAccepted_(0) 00036 { 00037 Photonsrc_ = iConfig.getParameter<InputTag>( "Photonsrc" ); 00038 Photon1Ptmin_ = 00039 iConfig.getParameter<double>( "Photon1Ptmin"); 00040 Photon2Ptmin_ = 00041 iConfig.getParameter<double>( "Photon2Ptmin"); 00042 IsIsolated_ = iConfig.getParameter<bool>( "IsIsolated"); 00043 IsolationCut_ = iConfig.getParameter<double>( "IsolationCut"); 00044 } 00045 00046 /*------------------------------------------------------------------------*/ 00047 00048 SUSYHighPtPhotonSkim::~SUSYHighPtPhotonSkim() 00049 {} 00050 00051 /*------------------------------------------------------------------------*/ 00052 00053 bool SUSYHighPtPhotonSkim::filter( edm::Event& iEvent, 00054 const edm::EventSetup& iSetup ) 00055 { 00056 nEvents_++; 00057 00058 typedef AssociationVector<RefProd<CandidateCollection>, vector<double> > PhotonMapCollection; 00059 Handle<PhotonMapCollection> PhotonHandle; 00060 00061 iEvent.getByLabel( Photonsrc_, PhotonHandle ); 00062 00063 if ( PhotonHandle->empty() ) return false; 00064 00065 int nPhoton1 = 0; 00066 int nPhoton2 = 0; 00067 00068 for ( PhotonMapCollection::const_iterator it = PhotonHandle->begin(); 00069 it != PhotonHandle->end(); it++ ) { 00070 00071 bool iso = it->second < IsolationCut_; 00072 if(!IsIsolated_) iso = 1; 00073 00074 if (iso && fabs(it->first->eta()) < 1.479) { 00075 if (it->first->pt() > Photon1Ptmin_) nPhoton1++; 00076 if (it->first->pt() > Photon2Ptmin_) nPhoton2++; 00077 } 00078 } 00079 00080 if ( !nPhoton1 ) return false; 00081 if ( nPhoton2<2 ) return false; 00082 00083 nAccepted_++; 00084 00085 return true; 00086 } 00087 00088 /*------------------------------------------------------------------------*/ 00089 00090 void SUSYHighPtPhotonSkim::endJob() 00091 { 00092 edm::LogVerbatim( "SUSYHighPtPhotonSkim" ) 00093 << "Events read " << nEvents_ 00094 << " Events accepted " << nAccepted_ 00095 << "\nEfficiency " << (double)(nAccepted_)/(double)(nEvents_) 00096 << endl; 00097 }