CMS 3D CMS Logo

MCPdgIndexFilter.cc
Go to the documentation of this file.
4 
6  token_(consumes<edm::HepMCProduct>(edm::InputTag(cfg.getUntrackedParameter("moduleLabel",std::string("generator")),"unsmeared"))),
7  pdgID(cfg.getParameter<std::vector<int> >("PdgId")),
8  index(cfg.getParameter<std::vector<unsigned> >("Index")),
9  maxIndex(*std::max_element(index.begin(),index.end())),
10  taggingMode(cfg.getUntrackedParameter<bool>("TagMode",false))
11 {
12  if (pdgID.size() != index.size())
13  edm::LogWarning("MCPdgIndexFilter")
14  << "Configuration Error :"
15  << "Sizes of array parameters 'PdgId' and 'Index' differ.";
16 
17  if (taggingMode) {
18  auto tag = cfg.getUntrackedParameter<std::string>("Tag","");
19  putToken_ = produces<bool>(tag);
20  edm::LogInfo("TagMode") << "Filter result in '" << tag << "', filtering disabled.";
21  }
22 }
23 
24 
26  bool result = pass(evt);
27  LogDebug("FilterResult") << (result?"Pass":"Fail");
28  if (!taggingMode) return result;
29  evt.emplace(putToken_,result);
30  return true;
31 }
32 
33 
34 bool MCPdgIndexFilter::pass(const edm::Event& evt) const {
36  evt.getByToken(token_, hepmc);
37 
38  const HepMC::GenEvent * genEvent = hepmc->GetEvent();
39 
40  HepMC::GenEvent::particle_const_iterator
41  p(genEvent->particles_begin()),
42  p_end(genEvent->particles_end());
43 
44  for ( unsigned i=0; p!=p_end && i<=maxIndex; ++p, i++ ) {
45  LogDebug("Particle") << "index: " << i << " pdgID: " << (*p)->pdg_id();
46  for (unsigned j = 0; j < pdgID.size(); j++) {
47  if (i==index[j] && pdgID[j] != (*p)->pdg_id())
48  return false;
49  }
50  }
51  return true;
52 }
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
const bool taggingMode
bool pass(const edm::Event &) const
MCPdgIndexFilter(const edm::ParameterSet &)
const std::vector< unsigned > index
#define end
Definition: vmac.h:39
const unsigned maxIndex
const edm::EDGetTokenT< edm::HepMCProduct > token_
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
OrphanHandle< PROD > emplace(EDPutTokenT< PROD > token, Args &&...args)
puts a new product
Definition: Event.h:413
#define begin
Definition: vmac.h:32
const std::vector< int > pdgID
HLT enums.
edm::EDPutTokenT< bool > putToken_