Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:13:40

0001 // -*- C++ -*-
0002 //
0003 // Package:    MCProcessRangeFilter
0004 // Class:      MCProcessRangeFilter
0005 //
0006 /*
0007 
0008  Description: filter events based on the Pythia ProcessID and the Pt_hat
0009 
0010  Implementation: inherits from generic EDFilter
0011 
0012 */
0013 //
0014 // Original Author:  Filip Moortgat
0015 //         Created:  Mon Sept 11 10:57:54 CET 2006
0016 //
0017 
0018 #include "DataFormats/Common/interface/Handle.h"
0019 #include "FWCore/Framework/interface/Event.h"
0020 #include "FWCore/Framework/interface/Frameworkfwd.h"
0021 #include "FWCore/Framework/interface/global/EDFilter.h"
0022 #include "FWCore/Framework/interface/MakerMacros.h"
0023 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0024 #include "FWCore/Utilities/interface/EDGetToken.h"
0025 #include "FWCore/Utilities/interface/InputTag.h"
0026 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0027 
0028 #include <string>
0029 
0030 class MCProcessRangeFilter : public edm::global::EDFilter<> {
0031 public:
0032   explicit MCProcessRangeFilter(const edm::ParameterSet&);
0033 
0034   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0035 
0036 private:
0037   const edm::EDGetTokenT<edm::HepMCProduct> token_;
0038   const int minProcessID;
0039   const int maxProcessID;
0040   const double pthatMin;
0041   const double pthatMax;
0042 };
0043 
0044 MCProcessRangeFilter::MCProcessRangeFilter(const edm::ParameterSet& iConfig)
0045     : token_(consumes<edm::HepMCProduct>(
0046           edm::InputTag(iConfig.getUntrackedParameter("moduleLabel", std::string("generator")), "unsmeared"))),
0047       minProcessID(iConfig.getUntrackedParameter("MinProcessID", 0)),
0048       maxProcessID(iConfig.getUntrackedParameter("MaxProcessID", 500)),
0049       pthatMin(iConfig.getUntrackedParameter("MinPthat", 0)),
0050       pthatMax(iConfig.getUntrackedParameter("MaxPthat", 14000)) {}
0051 
0052 bool MCProcessRangeFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0053   bool accepted = false;
0054   edm::Handle<edm::HepMCProduct> evt;
0055   iEvent.getByToken(token_, evt);
0056 
0057   const HepMC::GenEvent* myGenEvent = evt->GetEvent();
0058 
0059   // do the selection -- processID 0 is always accepted
0060 
0061   if (myGenEvent->signal_process_id() > minProcessID && myGenEvent->signal_process_id() < maxProcessID) {
0062     if (myGenEvent->event_scale() > pthatMin && myGenEvent->event_scale() < pthatMax) {
0063       accepted = true;
0064     }
0065   }
0066   return accepted;
0067 }
0068 
0069 DEFINE_FWK_MODULE(MCProcessRangeFilter);