File indexing completed on 2024-04-06 12:13:41
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
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 <cmath>
0029 #include <cstdlib>
0030 #include <string>
0031
0032 class PythiaFilterHT : public edm::global::EDFilter<> {
0033 public:
0034 explicit PythiaFilterHT(const edm::ParameterSet&);
0035
0036 bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0037
0038 private:
0039 const edm::EDGetTokenT<edm::HepMCProduct> label_;
0040 const double minhtcut;
0041 const int motherID;
0042 };
0043
0044 using namespace std;
0045
0046 PythiaFilterHT::PythiaFilterHT(const edm::ParameterSet& iConfig)
0047 : label_(consumes<edm::HepMCProduct>(
0048 edm::InputTag(iConfig.getUntrackedParameter("moduleLabel", std::string("generator")), "unsmeared"))),
0049 minhtcut(iConfig.getUntrackedParameter("MinHT", 0.)),
0050 motherID(iConfig.getUntrackedParameter("MotherID", 0)) {}
0051
0052 bool PythiaFilterHT::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0053 bool accepted = false;
0054 edm::Handle<edm::HepMCProduct> evt;
0055 iEvent.getByToken(label_, evt);
0056
0057 const HepMC::GenEvent* myGenEvent = evt->GetEvent();
0058
0059 double HT = 0;
0060
0061 for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
0062 ++p) {
0063 if (((*p)->status() == 23) && ((abs((*p)->pdg_id()) < 6) || ((*p)->pdg_id() == 21))) {
0064 if (motherID == 0) {
0065 HT += (*p)->momentum().perp();
0066 } else {
0067 HepMC::GenParticle* mother = (*((*p)->production_vertex()->particles_in_const_begin()));
0068 if (abs(mother->pdg_id()) == abs(motherID)) {
0069 HT += (*p)->momentum().perp();
0070 }
0071 }
0072 }
0073 }
0074 if (HT > minhtcut)
0075 accepted = true;
0076 return accepted;
0077 }
0078
0079 DEFINE_FWK_MODULE(PythiaFilterHT);