CMS 3D CMS Logo

HadSUSYQCDSkim Class Reference

all hadronic SUSY Skim >= 2 barrel jets 100 GeV, dphi, R1, R2 cuts QCD control, test L1 acoplanar path (QCD trigger path) More...

#include <SUSYBSMAnalysis/CSA07Skims/interface/HadSUSYQCDSkim.h>

Inheritance diagram for HadSUSYQCDSkim:

edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual void endJob ()
virtual bool filter (edm::Event &, const edm::EventSetup &)
 HadSUSYQCDSkim (const edm::ParameterSet &)
 ~HadSUSYQCDSkim ()

Private Attributes

double CaloJetPtmin_
edm::InputTag CaloJetsrc_
edm::InputTag CaloMETsrc_
unsigned int nAccepted_
unsigned int nEvents_
int NminCaloJet_


Detailed Description

all hadronic SUSY Skim >= 2 barrel jets 100 GeV, dphi, R1, R2 cuts QCD control, test L1 acoplanar path (QCD trigger path)

Date
2007/07/12 09:18:45
Revision
1.1

Author:
Michael Tytgat, Maria Spiropulu - CERN
Date
2007/09/25 17:54:51
Revision
1.7

Author:
Michael Tytgat, Maria Spiropulu - CERN

Definition at line 24 of file HadSUSYQCDSkim.h.


Constructor & Destructor Documentation

HadSUSYQCDSkim::HadSUSYQCDSkim ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 50 of file HadSUSYQCDSkim.cc.

References CaloJetPtmin_, CaloJetsrc_, CaloMETsrc_, edm::ParameterSet::getParameter(), and NminCaloJet_.

00050                                                                :
00051   nEvents_(0), nAccepted_(0)
00052 {
00053   CaloJetsrc_ = iConfig.getParameter<InputTag>( "CaloJetsrc" );
00054   CaloJetPtmin_ = 
00055     iConfig.getParameter<double>( "CaloJetPtmin");
00056   NminCaloJet_ = iConfig.getParameter<int>( "NminCaloJet");
00057   if ( NminCaloJet_ < 2 ) {
00058     edm::LogError( "HadSUSYQCDSkim" ) 
00059       << "Setting NminCaloJet to 2 !!";
00060     NminCaloJet_ = 2;
00061   }
00062   CaloMETsrc_ = iConfig.getParameter<InputTag>( "CaloMETsrc" );
00063 }

HadSUSYQCDSkim::~HadSUSYQCDSkim (  ) 

Definition at line 67 of file HadSUSYQCDSkim.cc.

00068 {}


Member Function Documentation

void HadSUSYQCDSkim::endJob ( void   )  [virtual]

Reimplemented from edm::EDFilter.

Definition at line 136 of file HadSUSYQCDSkim.cc.

References lat::endl(), nAccepted_, and nEvents_.

00137 {
00138   edm::LogVerbatim( "HadSUSYQCDSkim" ) 
00139     << "Events read " << nEvents_
00140     << " Events accepted " << nAccepted_
00141     << "\nEfficiency " << (double)(nAccepted_)/(double)(nEvents_) 
00142     << endl;
00143 }

bool HadSUSYQCDSkim::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDFilter.

Definition at line 72 of file HadSUSYQCDSkim.cc.

References CaloJetPtmin_, CaloJetsrc_, CaloMETsrc_, DeltaPhi(), edm::Event::getByLabel(), it, nAccepted_, nEvents_, NminCaloJet_, and funct::sqrt().

00074 {
00075   nEvents_++;
00076 
00077   Handle<CaloJetCollection> CaloJetsHandle;
00078 //  try {
00079     iEvent.getByLabel( CaloJetsrc_, CaloJetsHandle );
00080 //  } 
00081 //  catch ( cms::Exception& ex ) {
00082 //    edm::LogError( "HadSUSYQCDSkim" ) 
00083 //      << "Unable to get CaloJet collection "
00084 //      << CaloJetsrc_.label();
00085 //    return false;
00086 //  }
00087   if ( CaloJetsHandle->empty() ) return false;
00088   CaloJetCollection TheCaloJets = *CaloJetsHandle;
00089   std::stable_sort( TheCaloJets.begin(), TheCaloJets.end(), PtSorter() );
00090 
00091   Handle<CaloMETCollection> METHandle;
00092 //  try {
00093     iEvent.getByLabel( CaloMETsrc_, METHandle );
00094 //  }
00095 //  catch ( cms::Exception& ex ) {
00096 //    edm::LogError( "HadSUSYQCDSkim" )
00097 //      << "Unable to get CaloMET collection "
00098 //      << CaloMETsrc_.label();
00099 //    return false;
00100 //  }
00101   if ( METHandle->empty() ) return false;
00102   double METphi = METHandle->begin()->phi();
00103 
00104   // Jet and phi cuts
00105   int nJet = 0;
00106   double phiJet[2];
00107   for ( CaloJetCollection::const_iterator it = TheCaloJets.begin(); 
00108         it != TheCaloJets.end() ; it++ ) {
00109     if ( (fabs(it->eta()) < 3.0) &&
00110          (it->pt() > CaloJetPtmin_) ) {
00111      if(nJet<2)  phiJet[nJet] = it->phi();
00112        nJet++;
00113       if ( DeltaPhi( it->phi(), METphi ) < 0.3 ) return false;
00114     }
00115   }
00116   if ( nJet < NminCaloJet_ ) return false;
00117 
00118   double dphi1 = DeltaPhi( phiJet[0], METphi );
00119   double dphi2 = DeltaPhi( phiJet[1], METphi );
00120   double R1 = sqrt( dphi2*dphi2 
00121                     + (2.*acos(-1.)-dphi1)*(2.*acos(-1.)-dphi1) );
00122   if ( R1 < 0.5 ) return false;
00123   double R2 = sqrt( dphi1*dphi1 
00124                     + (2.*acos(-1.)-dphi2)*(2.*acos(-1.)-dphi2) ); 
00125   if ( R2 < 0.5 ) return false;
00126   if ( DeltaPhi( METphi, phiJet[1] ) < (20./180.*acos(-1.)) )
00127     return false;
00128   
00129   nAccepted_++;
00130 
00131   return true;
00132 }


Member Data Documentation

double HadSUSYQCDSkim::CaloJetPtmin_ [private]

Definition at line 34 of file HadSUSYQCDSkim.h.

Referenced by filter(), and HadSUSYQCDSkim().

edm::InputTag HadSUSYQCDSkim::CaloJetsrc_ [private]

Definition at line 33 of file HadSUSYQCDSkim.h.

Referenced by filter(), and HadSUSYQCDSkim().

edm::InputTag HadSUSYQCDSkim::CaloMETsrc_ [private]

Definition at line 36 of file HadSUSYQCDSkim.h.

Referenced by filter(), and HadSUSYQCDSkim().

unsigned int HadSUSYQCDSkim::nAccepted_ [private]

Definition at line 38 of file HadSUSYQCDSkim.h.

Referenced by endJob(), and filter().

unsigned int HadSUSYQCDSkim::nEvents_ [private]

Definition at line 37 of file HadSUSYQCDSkim.h.

Referenced by endJob(), and filter().

int HadSUSYQCDSkim::NminCaloJet_ [private]

Definition at line 35 of file HadSUSYQCDSkim.h.

Referenced by filter(), and HadSUSYQCDSkim().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:23:26 2009 for CMSSW by  doxygen 1.5.4