MCProcessRangeFilter

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69
// -*- C++ -*-
//
// Package:    MCProcessRangeFilter
// Class:      MCProcessRangeFilter
//
/*

 Description: filter events based on the Pythia ProcessID and the Pt_hat

 Implementation: inherits from generic EDFilter

*/
//
// Original Author:  Filip Moortgat
//         Created:  Mon Sept 11 10:57:54 CET 2006
//

#include "DataFormats/Common/interface/Handle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/global/EDFilter.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/EDGetToken.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"

#include <string>

class MCProcessRangeFilter : public edm::global::EDFilter<> {
public:
  explicit MCProcessRangeFilter(const edm::ParameterSet&);

  bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;

private:
  const edm::EDGetTokenT<edm::HepMCProduct> token_;
  const int minProcessID;
  const int maxProcessID;
  const double pthatMin;
  const double pthatMax;
};

MCProcessRangeFilter::MCProcessRangeFilter(const edm::ParameterSet& iConfig)
    : token_(consumes<edm::HepMCProduct>(
          edm::InputTag(iConfig.getUntrackedParameter("moduleLabel", std::string("generator")), "unsmeared"))),
      minProcessID(iConfig.getUntrackedParameter("MinProcessID", 0)),
      maxProcessID(iConfig.getUntrackedParameter("MaxProcessID", 500)),
      pthatMin(iConfig.getUntrackedParameter("MinPthat", 0)),
      pthatMax(iConfig.getUntrackedParameter("MaxPthat", 14000)) {}

bool MCProcessRangeFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
  bool accepted = false;
  edm::Handle<edm::HepMCProduct> evt;
  iEvent.getByToken(token_, evt);

  const HepMC::GenEvent* myGenEvent = evt->GetEvent();

  // do the selection -- processID 0 is always accepted

  if (myGenEvent->signal_process_id() > minProcessID && myGenEvent->signal_process_id() < maxProcessID) {
    if (myGenEvent->event_scale() > pthatMin && myGenEvent->event_scale() < pthatMax) {
      accepted = true;
    }
  }
  return accepted;
}

DEFINE_FWK_MODULE(MCProcessRangeFilter);