CMS 3D CMS Logo

ProbeMulteplicityProducer.cc
Go to the documentation of this file.
1 //
2 
14 
18 
22 
24 public:
25  explicit ProbeMulteplicityProducer(const edm::ParameterSet& iConfig);
26  ~ProbeMulteplicityProducer() override;
27 
28  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
29 
30 private:
33  pairCut_; // lazy parsing, to allow cutting on variables not in reco::Candidate class
35  probeCut_; // lazy parsing, to allow cutting on variables not in reco::Candidate class
36 };
37 
39  : pairs_(consumes<reco::CandidateView>(iConfig.getParameter<edm::InputTag>("pairs"))),
40  pairCut_(iConfig.existsAs<std::string>("pairSelection") ? iConfig.getParameter<std::string>("pairSelection") : "",
41  true),
42  probeCut_(
43  iConfig.existsAs<std::string>("probeSelection") ? iConfig.getParameter<std::string>("probeSelection") : "",
44  true) {
45  produces<edm::ValueMap<float>>();
46 }
47 
49 
51  using namespace edm;
52 
53  // read input
55  iEvent.getByToken(pairs_, pairs);
56 
57  // fill
58  unsigned int i = 0;
59  std::vector<unsigned int> tagKeys;
60  std::vector<float> values;
61  View<reco::Candidate>::const_iterator pair, endpairs = pairs->end();
62  for (pair = pairs->begin(); pair != endpairs; ++pair, ++i) {
63  reco::CandidateBaseRef probeRef = pair->daughter(1)->masterClone();
64  unsigned int tagKey = pair->daughter(0)->masterClone().key();
65  unsigned int copies = 1;
66  if (pairCut_(*pair) && probeCut_(*probeRef)) {
67  for (unsigned int j = 0; j < i; ++j)
68  if (tagKeys[j] == tagKey)
69  copies++;
70  for (unsigned int j = 0; j < i; ++j)
71  if (tagKeys[j] == tagKey)
72  values[j] = copies;
73  } else {
75  copies = 0;
76  }
77  tagKeys.push_back(tagKey);
78  values.push_back(copies);
79  }
80 
81  // convert into ValueMap and store
82  auto valMap = std::make_unique<ValueMap<float>>();
84  filler.insert(pairs, values.begin(), values.end());
85  filler.fill();
86  iEvent.put(std::move(valMap));
87 }
88 
edm::EDGetTokenT< reco::CandidateView > pairs_
StringCutObjectSelector< reco::Candidate, true > probeCut_
StringCutObjectSelector< reco::Candidate, true > pairCut_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Matcher of number of reconstructed objects in the event to probe.
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
ProbeMulteplicityProducer(const edm::ParameterSet &iConfig)
def move(src, dest)
Definition: eostools.py:511
edm::View< Candidate > CandidateView
view of a collection containing candidates
Definition: CandidateFwd.h:23