CMS 3D CMS Logo

SUSYControlHighPtPhotonSkim.cc

Go to the documentation of this file.
00001 /* \class SUSYControlHighPtPhotonSkim
00002  *
00003  * High Energy Photon SUSY Skim (control sample)
00004  * one photon and one electron > xx GeV in barrel + isolation 
00005  *
00006  * $Date: 2007/12/09 10:54:57 $
00007  * $Revision: 1.6 $
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/EgammaCandidates/interface/GsfElectron.h"
00026 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00027 
00028 #include "DataFormats/Candidate/interface/Candidate.h"
00029 #include "DataFormats/Common/interface/AssociationVector.h"
00030 
00031 #include "SUSYBSMAnalysis/CSA07Skims/interface/SUSYControlHighPtPhotonSkim.h"
00032 
00033 using namespace edm;
00034 using namespace std;
00035 using namespace reco;
00036 
00037 SUSYControlHighPtPhotonSkim::SUSYControlHighPtPhotonSkim( const edm::ParameterSet& iConfig ) :
00038   nEvents_(0), nAccepted_(0)
00039 {
00040   Photonsrc_ = iConfig.getParameter<InputTag>( "Photonsrc" );
00041   Electronsrc_ = iConfig.getParameter<InputTag>( "Electronsrc" );
00042   PhotonPtmin_ = 
00043     iConfig.getParameter<double>( "PhotonPtmin");
00044   ElectronPtmin_ = 
00045     iConfig.getParameter<double>( "ElectronPtmin");
00046   IsIsolated_ = iConfig.getParameter<bool>( "IsIsolated");
00047   IsolationCut_ = iConfig.getParameter<double>( "IsolationCut");
00048 }
00049 
00050 /*------------------------------------------------------------------------*/
00051 
00052 SUSYControlHighPtPhotonSkim::~SUSYControlHighPtPhotonSkim() 
00053 {}
00054 
00055 /*------------------------------------------------------------------------*/
00056 
00057 bool SUSYControlHighPtPhotonSkim::filter( edm::Event& iEvent, 
00058                                        const edm::EventSetup& iSetup )
00059 {
00060   nEvents_++;
00061 
00062   typedef AssociationVector<RefProd<CandidateCollection>, vector<double> > PhotonMapCollection;
00063   Handle<PhotonMapCollection>  PhotonHandle;
00064 
00065   iEvent.getByLabel( Photonsrc_, PhotonHandle );
00066 
00067   if ( PhotonHandle->empty() ) return false;
00068 
00069   Handle<GsfElectronCollection>  ElectronHandle;
00070 
00071   iEvent.getByLabel( Electronsrc_, ElectronHandle );
00072 
00073   if ( ElectronHandle->empty() ) return false;
00074 
00075   int nPhoton = 0;
00076   int nElectron = 0;
00077 
00078   for ( PhotonMapCollection::const_iterator it = PhotonHandle->begin(); 
00079         it != PhotonHandle->end(); it++ ) {
00080 
00081     bool iso = it->second < IsolationCut_;
00082     if(!IsIsolated_) iso = 1; 
00083 
00084     if (iso && fabs(it->first->eta()) < 1.479 && it->first->pt() > PhotonPtmin_) { 
00085 
00086       int overlap(0);
00087       for ( GsfElectronCollection::const_iterator itel = ElectronHandle->begin(); 
00088             itel != ElectronHandle->end(); itel++ ) {
00089         
00090         if( itel->pt() > ElectronPtmin_ ){
00091           double dr = sqrt( pow(it->first->eta() - itel->eta(),2) +
00092                             pow(it->first->phi() - itel->phi(),2) );
00093           if ( dr<0.1 ) overlap = 1;
00094         }
00095       }
00096       if (!overlap) nPhoton++;
00097 
00098     }
00099 
00100   }
00101 
00102   for ( GsfElectronCollection::const_iterator itel = ElectronHandle->begin(); 
00103         itel != ElectronHandle->end(); itel++ ) {
00104 
00105     if (fabs(itel->eta()) < 1.479 && itel->pt() > ElectronPtmin_) 
00106       nElectron++;
00107  }
00108 
00109   if (!nPhoton) return false;
00110   if (!nElectron) return false;
00111 
00112   nAccepted_++;
00113 
00114   return true;
00115 }
00116 
00117 /*------------------------------------------------------------------------*/
00118 
00119 void SUSYControlHighPtPhotonSkim::endJob()
00120 {
00121   edm::LogVerbatim( "SUSYControlHighPtPhotonSkim" ) 
00122     << "Events read " << nEvents_
00123     << " Events accepted " << nAccepted_
00124     << "\nEfficiency " << (double)(nAccepted_)/(double)(nEvents_) 
00125     << endl;
00126 }

Generated on Tue Jun 9 17:48:04 2009 for CMSSW by  doxygen 1.5.4