CMS 3D CMS Logo

PFMatchedCandidateRefExtractor.cc
Go to the documentation of this file.
6 
12 
14 
16 public:
17  explicit PFMatchedCandidateRefExtractor(const edm::ParameterSet& iConfig);
19 
20  void produce(edm::StreamID iID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
21 
22 private:
27 
29 };
30 
32  : col1Token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("col1"))),
33  col2Token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("col2"))),
34  moduleLabel_(iConfig.getParameter<std::string>("@module_label")) {
35  //register products
36 
37  extractPFCands_ = iConfig.getParameter<bool>("extractPFCandidates");
38 
39  produces<edm::PtrVector<reco::Candidate>>("col1");
40  produces<edm::PtrVector<reco::Candidate>>("col2");
41  if (extractPFCands_) {
42  pfCandToken_ = mayConsume<edm::View<reco::PFCandidate>>(iConfig.getParameter<edm::InputTag>("pfCandCollection"));
43  produces<edm::PtrVector<reco::PFCandidate>>("pfCandCol1");
44  produces<edm::PtrVector<reco::PFCandidate>>("pfCandCol2");
45  }
46 }
47 
49 
50 // ------------ method called to produce the data ------------
52  std::unique_ptr<edm::PtrVector<reco::Candidate>> outcol1(new edm::PtrVector<reco::Candidate>());
53  std::unique_ptr<edm::PtrVector<reco::Candidate>> outcol2(new edm::PtrVector<reco::Candidate>());
54 
55  std::unique_ptr<edm::PtrVector<reco::PFCandidate>> outPFCandCol1(new edm::PtrVector<reco::PFCandidate>());
56  std::unique_ptr<edm::PtrVector<reco::PFCandidate>> outPFCandCol2(new edm::PtrVector<reco::PFCandidate>());
57 
61 
62  iEvent.getByToken(col1Token_, col1Handle);
63  iEvent.getByToken(col2Token_, col2Handle);
64  if (extractPFCands_) {
65  iEvent.getByToken(pfCandToken_, pfCandHandle);
66  }
67 
68  for (size_t iC1 = 0; iC1 < col1Handle->size(); iC1++) {
69  for (size_t iC2 = 0; iC2 < col2Handle->size(); iC2++) {
70  bool matched = (col1Handle->ptrAt(iC1) == col2Handle->ptrAt(iC2));
71  if (!matched && deltaR2(col1Handle->ptrAt(iC1)->p4(), col2Handle->ptrAt(iC2)->p4()) < 0.000001) {
72  matched = true;
73  }
74 
75  if (matched) {
76  outcol1->push_back(col1Handle->ptrAt(iC1));
77  outcol2->push_back(col2Handle->ptrAt(iC2));
78 
79  if (!extractPFCands_)
80  continue;
81  std::set<reco::CandidatePtr> sc1s;
82  std::set<reco::CandidatePtr> sc2s;
83  for (size_t ics1 = 0; ics1 < col1Handle->ptrAt(iC1)->numberOfSourceCandidatePtrs(); ics1++) {
84  sc1s.insert(col1Handle->ptrAt(iC1)->sourceCandidatePtr(ics1));
85  }
86  for (size_t ics2 = 0; ics2 < col2Handle->ptrAt(iC2)->numberOfSourceCandidatePtrs(); ics2++) {
87  sc2s.insert(col2Handle->ptrAt(iC2)->sourceCandidatePtr(ics2));
88  }
89 
90  for (size_t ic = 0; ic < pfCandHandle->size(); ++ic) {
91  reco::PFCandidatePtr c = pfCandHandle->ptrAt(ic);
92 
93  bool match1 = (sc1s.find(c) != sc1s.end());
94  bool match2 = (sc2s.find(c) != sc2s.end());
95 
96  if (!match1) { //recovery when pfcandidate sources are not equivalent
97  for (size_t ics1 = 0; ics1 < sc1s.size(); ics1++) {
98  if (deltaR2(c->p4(),
99  col1Handle->ptrAt(iC1)->sourceCandidatePtr(ics1)->p4()) <
100  0.0000001) { //tight dR, pfcandidates should be the same
101  match1 = true;
102  break;
103  }
104  }
105  }
106  if (!match2) { //recovery when pfcandidate sources are not equivalent
107  for (size_t ics2 = 0; ics2 < sc2s.size(); ics2++) {
108  if (deltaR2(c->p4(),
109  col2Handle->ptrAt(iC2)->sourceCandidatePtr(ics2)->p4()) < 0.00001) { //tight dR
110  match2 = true;
111  break;
112  }
113  }
114  }
115 
116  if (match1) {
117  outPFCandCol1->push_back(c);
118  }
119  if (match2) {
120  outPFCandCol2->push_back(c);
121  }
122  } //pfcand loop
123 
124  } //matching
125  } //col2
126  } //col1
127 
128  iEvent.put(std::move(outcol1), "col1");
129  iEvent.put(std::move(outcol2), "col2");
130  if (extractPFCands_) {
131  iEvent.put(std::move(outPFCandCol1), "pfCandCol1");
132  iEvent.put(std::move(outPFCandCol2), "pfCandCol2");
133  }
134 }
135 
136 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
PFMatchedCandidateRefExtractor(const edm::ParameterSet &iConfig)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< edm::View< reco::Candidate > > col1Token_
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< edm::View< reco::Candidate > > col2Token_
edm::EDGetTokenT< edm::View< reco::PFCandidate > > pfCandToken_
fixed size matrix
HLT enums.
void produce(edm::StreamID iID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
def move(src, dest)
Definition: eostools.py:511