Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "GeneratorInterface/GenFilters/plugins/MCPdgIndexFilter.h"
0002 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 
0005 MCPdgIndexFilter::MCPdgIndexFilter(const edm::ParameterSet& cfg)
0006     : token_(consumes<edm::HepMCProduct>(
0007           edm::InputTag(cfg.getUntrackedParameter("moduleLabel", std::string("generator")), "unsmeared"))),
0008       pdgID(cfg.getParameter<std::vector<int> >("PdgId")),
0009       index(cfg.getParameter<std::vector<unsigned> >("Index")),
0010       maxIndex(*std::max_element(index.begin(), index.end())),
0011       taggingMode(cfg.getUntrackedParameter<bool>("TagMode", false)) {
0012   if (pdgID.size() != index.size())
0013     edm::LogWarning("MCPdgIndexFilter") << "Configuration Error :"
0014                                         << "Sizes of array parameters 'PdgId' and 'Index' differ.";
0015 
0016   if (taggingMode) {
0017     auto tag = cfg.getUntrackedParameter<std::string>("Tag", "");
0018     putToken_ = produces<bool>(tag);
0019     edm::LogInfo("TagMode") << "Filter result in '" << tag << "', filtering disabled.";
0020   }
0021 }
0022 
0023 bool MCPdgIndexFilter::filter(edm::StreamID, edm::Event& evt, const edm::EventSetup&) const {
0024   bool result = pass(evt);
0025   LogDebug("FilterResult") << (result ? "Pass" : "Fail");
0026   if (!taggingMode)
0027     return result;
0028   evt.emplace(putToken_, result);
0029   return true;
0030 }
0031 
0032 bool MCPdgIndexFilter::pass(const edm::Event& evt) const {
0033   edm::Handle<edm::HepMCProduct> hepmc;
0034   evt.getByToken(token_, hepmc);
0035 
0036   const HepMC::GenEvent* genEvent = hepmc->GetEvent();
0037 
0038   HepMC::GenEvent::particle_const_iterator p(genEvent->particles_begin()), p_end(genEvent->particles_end());
0039 
0040   for (unsigned i = 0; p != p_end && i <= maxIndex; ++p, i++) {
0041     LogDebug("Particle") << "index: " << i << "   pdgID: " << (*p)->pdg_id();
0042     for (unsigned j = 0; j < pdgID.size(); j++) {
0043       if (i == index[j] && pdgID[j] != (*p)->pdg_id())
0044         return false;
0045     }
0046   }
0047   return true;
0048 }