CMS 3D CMS Logo

HadSUSYdiElecSkim.cc

Go to the documentation of this file.
00001 
00014 #include <iostream>
00015 #include <string>
00016 #include <list>
00017 #include <cmath>
00018 #include <cstdio>
00019 #include <vector>
00020 #include <memory>
00021 
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 #include "DataFormats/Common/interface/Handle.h"    
00024 
00025 #include "SUSYBSMAnalysis/CSA07Skims/interface/HadSUSYdiElecSkim.h"
00026 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00027 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00028 
00029 using namespace edm;
00030 using namespace std;
00031 using namespace reco;
00032 
00033 class PtSorter {
00034 public:
00035   template <class T> bool operator() ( const T& a, const T& b ) {
00036     return ( a.pt() > b.pt() );
00037   }
00038 };
00039 
00040 HadSUSYdiElecSkim::HadSUSYdiElecSkim( const edm::ParameterSet& iConfig ) :
00041   nEvents_(0), nAccepted_(0)
00042 {
00043   Elecsrc_ = iConfig.getParameter<InputTag>( "Elecsrc" );
00044   NminElec_ = iConfig.getParameter<int>( "NminElec");
00045   ElecPtmin_ = iConfig.getParameter<double>( "ElecPtmin");
00046   PtmindiElec_ = iConfig.getParameter<double>( "PtmindiElec");
00047 }
00048 
00049 /*------------------------------------------------------------------------*/
00050 
00051 HadSUSYdiElecSkim::~HadSUSYdiElecSkim() 
00052 {}
00053 
00054 /*------------------------------------------------------------------------*/
00055 
00056 bool HadSUSYdiElecSkim::filter( edm::Event& iEvent, 
00057                                 const edm::EventSetup& iSetup )
00058 {
00059   nEvents_++;
00060 
00061   Handle<GsfElectronCollection> ElecHandle;
00062 //  try {
00063     iEvent.getByLabel( Elecsrc_, ElecHandle );
00064 //  }
00065 //  catch ( cms::Exception& ex ) {
00066 //    edm::LogError( "HadSUSYdiElecSkim" ) 
00067 //      << "Unable to get Elec collection "
00068 //      << Elecsrc_.label();
00069 //    return false;
00070 //  }
00071   if ( ElecHandle->empty() ) return false;
00072   GsfElectronCollection TheElecs = *ElecHandle;
00073   std::stable_sort( TheElecs.begin(), TheElecs.end(), PtSorter() );
00074   
00075   int nElec = 0;
00076   double Pxdielec = 0., Pydielec = 0.;
00077   for ( GsfElectronCollection::const_iterator it = TheElecs.begin();
00078         it != TheElecs.end(); it++ ) {
00079     if ( (it->pt() > ElecPtmin_) 
00080          && (fabs(it->eta()) < 3.0) ) {
00081       if ( nElec < 2 ) {
00082         Pxdielec += it->p()*sin(it->theta())*cos(it->phi());
00083         Pydielec += it->p()*sin(it->theta())*sin(it->phi());
00084       }
00085       nElec++;
00086     }
00087   }
00088   
00089   if ( nElec < NminElec_ ) return false;
00090 
00091   double PtdiElec = sqrt( Pxdielec*Pxdielec + Pydielec*Pydielec );
00092   if ( PtdiElec < PtmindiElec_ ) return false;
00093     
00094   nAccepted_++;
00095 
00096   return true;
00097 }
00098 
00099 /*------------------------------------------------------------------------*/
00100 
00101 void HadSUSYdiElecSkim::endJob()
00102 {
00103   edm::LogVerbatim( "HadSUSYdiElecSkim" ) 
00104     << "Events read " << nEvents_
00105     << " Events accepted " << nAccepted_
00106     << "\nEfficiency " << (double)(nAccepted_)/(double)(nEvents_) 
00107     << endl;
00108 }

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