CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/HLTrigger/special/src/HLTPhysicsDeclared.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:   HLTPhysicsDeclared
00004 // Class:     HLTPhysicsDeclared
00005 //
00006 // Original Author:     Luca Malgeri
00007 // Adapted for HLT by:  Andrea Bocci
00008 
00009 // user include files
00010 #include "FWCore/Framework/interface/Frameworkfwd.h"
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "FWCore/Utilities/interface/InputTag.h"
00014 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
00015 
00016 //
00017 // class declaration
00018 //
00019 
00020 class HLTPhysicsDeclared : public edm::EDFilter {
00021 public:
00022   explicit HLTPhysicsDeclared( const edm::ParameterSet & );
00023   ~HLTPhysicsDeclared();
00024   
00025 private:
00026   virtual bool filter( edm::Event &, const edm::EventSetup & );
00027 
00028   bool          m_invert;
00029   edm::InputTag m_gtDigis;
00030   
00031 };
00032 
00033 // system include files
00034 #include <memory>
00035 
00036 // user include files
00037 #include "FWCore/Framework/interface/Frameworkfwd.h"
00038 #include "FWCore/Framework/interface/EDFilter.h"
00039 #include "FWCore/Framework/interface/Event.h"
00040 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00041 #include "FWCore/Utilities/interface/InputTag.h"
00042 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00043 #include "DataFormats/L1GlobalTrigger/interface/L1GtFdlWord.h"
00044 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00045 
00046 using namespace edm;
00047 
00048 HLTPhysicsDeclared::HLTPhysicsDeclared(const edm::ParameterSet & config) :
00049   m_invert(  config.getParameter<bool>("invert") ),
00050   m_gtDigis( config.getParameter<edm::InputTag>("L1GtReadoutRecordTag") )
00051 {
00052 }
00053 
00054 HLTPhysicsDeclared::~HLTPhysicsDeclared()
00055 {
00056 }
00057 
00058 bool HLTPhysicsDeclared::filter( edm::Event & event, const edm::EventSetup & setup)
00059 {
00060   bool accept = false;
00061 
00062   if (event.isRealData()) {
00063     // for real data, access the "physics enabled" bit in the L1 GT data
00064     edm::Handle<L1GlobalTriggerReadoutRecord> h_gtDigis;
00065     if (not event.getByLabel(m_gtDigis, h_gtDigis)) {
00066       edm::LogWarning(h_gtDigis.whyFailed()->category()) << h_gtDigis.whyFailed()->what();
00067       // not enough informations to make a decision - reject the event
00068       return false;
00069     } else {
00070       L1GtFdlWord fdlWord = h_gtDigis->gtFdlWord();
00071       if (fdlWord.physicsDeclared() == 1) 
00072         accept = true;
00073     }
00074   } else {
00075     // for MC, assume the "physics enabled" bit to be always set
00076     accept = true;
00077   }
00078 
00079   // if requested, invert the filter decision
00080   if (m_invert)
00081     accept = not accept;
00082 
00083   return accept;
00084 }
00085 
00086 // define this as a framework plug-in
00087 #include "FWCore/Framework/interface/MakerMacros.h"
00088 DEFINE_FWK_MODULE(HLTPhysicsDeclared);