CMS 3D CMS Logo

DuplicateListMerger.cc
Go to the documentation of this file.
1 
9 
15 
19 
22 
24 
27 
28 #include <vector>
29 #include <algorithm>
30 #include <string>
31 #include <iostream>
32 #include <memory>
33 
34 // #include "TMVA/Reader.h"
35 
36 using namespace reco;
37 
38 namespace {
39  class DuplicateListMerger final : public edm::global::EDProducer<> {
40  public:
42  explicit DuplicateListMerger(const edm::ParameterSet& iPara);
44  ~DuplicateListMerger() override;
45 
47  using CandidateToDuplicate = std::vector<std::pair<int, int>>;
48 
49  using RecHitContainer = edm::OwnVector<TrackingRecHit>;
50 
51  using MVACollection = std::vector<float>;
52  using QualityMaskCollection = std::vector<unsigned char>;
53 
54  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
55 
56  private:
58  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
59 
60  private:
61  TrackCollectionCloner collectionCloner;
62  TrackCollectionCloner::Tokens mergedTrackSource_;
63  TrackCollectionCloner::Tokens originalTrackSource_;
64 
65  edm::EDGetTokenT<CandidateToDuplicate> candidateComponents_;
67 
68  edm::EDGetTokenT<MVACollection> originalMVAValsToken_;
69  edm::EDGetTokenT<MVACollection> mergedMVAValsToken_;
71 
72  std::string priorityName_;
73 
74  int diffHitsCut_;
75  };
76 } // namespace
77 
82 
85 
86 #include <tuple>
87 #include <array>
89 
90 namespace {
93  desc.add<edm::InputTag>("mergedSource", edm::InputTag());
94  desc.add<edm::InputTag>("originalSource", edm::InputTag());
95  desc.add<edm::InputTag>("mergedMVAVals", edm::InputTag());
96  desc.add<edm::InputTag>("originalMVAVals", edm::InputTag());
97  desc.add<edm::InputTag>("candidateSource", edm::InputTag());
98  desc.add<edm::InputTag>("candidateComponents", edm::InputTag());
99  desc.add<std::string>("trackAlgoPriorityOrder", "trackAlgoPriorityOrder");
100  desc.add<int>("diffHitsCut", 5);
102  descriptions.add("DuplicateListMerger", desc);
103  }
104 
105  DuplicateListMerger::DuplicateListMerger(const edm::ParameterSet& iPara)
106  : collectionCloner(producesCollector(), iPara, true),
107  mergedTrackSource_(iPara.getParameter<edm::InputTag>("mergedSource"), consumesCollector()),
108  originalTrackSource_(iPara.getParameter<edm::InputTag>("originalSource"), consumesCollector()),
109  priorityName_(iPara.getParameter<std::string>("trackAlgoPriorityOrder")) {
110  diffHitsCut_ = iPara.getParameter<int>("diffHitsCut");
111  candidateSource_ = consumes<std::vector<TrackCandidate>>(iPara.getParameter<edm::InputTag>("candidateSource"));
112  candidateComponents_ = consumes<CandidateToDuplicate>(iPara.getParameter<edm::InputTag>("candidateComponents"));
113 
114  mergedMVAValsToken_ = consumes<MVACollection>(iPara.getParameter<edm::InputTag>("mergedMVAVals"));
115  originalMVAValsToken_ = consumes<MVACollection>(iPara.getParameter<edm::InputTag>("originalMVAVals"));
116  priorityOrderToken_ = esConsumes<TrackAlgoPriorityOrder, CkfComponentsRecord>(edm::ESInputTag("", priorityName_));
117 
118  produces<MVACollection>("MVAValues");
119  produces<QualityMaskCollection>("QualityMasks");
120  }
121 
122  DuplicateListMerger::~DuplicateListMerger() { /* no op */
123  }
124 
125  void DuplicateListMerger::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
127 
128  auto const& originals = originalTrackSource_.tracks(iEvent);
129  auto const& merged = mergedTrackSource_.tracks(iEvent);
130  auto const& candIndices = mergedTrackSource_.indicesInput(iEvent);
131 
133  iEvent.getByToken(candidateSource_, candidateH);
134  auto const& candidates = *candidateH;
135 
136  edm::Handle<CandidateToDuplicate> candidateComponentsH;
137  iEvent.getByToken(candidateComponents_, candidateComponentsH);
138  auto const& candidateComponents = *candidateComponentsH;
139 
140  edm::Handle<MVACollection> originalMVAStore;
141  edm::Handle<MVACollection> mergedMVAStore;
142 
143  iEvent.getByToken(originalMVAValsToken_, originalMVAStore);
144  iEvent.getByToken(mergedMVAValsToken_, mergedMVAStore);
145 
146  edm::ESHandle<TrackAlgoPriorityOrder> priorityH = iSetup.getHandle(priorityOrderToken_);
147  auto const& trackAlgoPriorityOrder = *priorityH;
148 
149  MVACollection mvaVec;
150 
151  auto mergedMVA = *mergedMVAStore;
152 
153  //match new tracks to their candidates
154  std::vector<std::array<int, 3>> matches;
155  for (int i = 0; i < (int)merged.size(); ++i) {
156  auto cInd = candIndices[i];
157  auto const& cand = candidates[cInd];
158  const reco::Track& matchedTrack = merged[i];
159 
160  if (mergedMVA[i] < -0.7f)
161  continue; // at least "loose" ( FIXME: take cut value from CutSelector)
162 
163  // if( ChiSquaredProbability(matchedTrack.chi2(),matchedTrack.ndof()) < minTrkProbCut_)continue;
164  int dHits = cand.nRecHits() - matchedTrack.recHitsSize();
165  if (dHits > diffHitsCut_)
166  continue;
167  matches.push_back(std::array<int, 3>{{i, candidateComponents[cInd].first, candidateComponents[cInd].second}});
168  }
169 
170  //check for candidates/tracks that share merged tracks, select minimum chi2, remove the rest
171  if (matches.size() > 1)
172  for (auto matchIter0 = matches.begin(); matchIter0 != matches.end() - 1; ++matchIter0) {
173  if ((*matchIter0)[0] < 0)
174  continue;
175  auto nchi2 = merged[(*matchIter0)[0]].normalizedChi2();
176  for (auto matchIter1 = matchIter0 + 1; matchIter1 != matches.end(); ++matchIter1) {
177  if ((*matchIter1)[0] < 0)
178  continue;
179  if ((*matchIter0)[1] == (*matchIter1)[1] || (*matchIter0)[1] == (*matchIter1)[2] ||
180  (*matchIter0)[2] == (*matchIter1)[1] || (*matchIter0)[2] == (*matchIter1)[2]) {
181  auto nchi2_1 = merged[(*matchIter1)[0]].normalizedChi2();
182  if (nchi2_1 < nchi2) {
183  (*matchIter0)[0] = -1;
184  break;
185  } else {
186  (*matchIter1)[0] = -1;
187  }
188  }
189  }
190  }
191 
192  // products
193  auto pmvas = std::make_unique<MVACollection>();
194  auto pquals = std::make_unique<QualityMaskCollection>();
195 
196  //add the good merged tracks to the output list, remove input tracks
197  std::vector<int> inputTracks;
198 
199  std::vector<unsigned int> selId;
200  auto ntotTk = matches.size();
201  // declareDynArray(reco::TrackBase::TrackAlgorithm, ntotTk, algo);
203  declareDynArray(reco::TrackBase::AlgoMask, ntotTk, algoMask);
204 
205  auto nsel = 0U;
206  for (auto matchIter0 = matches.begin(); matchIter0 != matches.end(); matchIter0++) {
207  if ((*matchIter0)[0] < 0)
208  continue;
209  selId.push_back((*matchIter0)[0]);
210 
211  pmvas->push_back(mergedMVA[(*matchIter0)[0]]);
212 
213  const reco::Track& inTrk1 = originals[(*matchIter0)[1]];
214  const reco::Track& inTrk2 = originals[(*matchIter0)[2]];
215  oriAlgo[nsel] = std::min(
217  return trackAlgoPriorityOrder.priority(a) < trackAlgoPriorityOrder.priority(b);
218  });
219 
220  algoMask[nsel] = inTrk1.algoMask() | inTrk2.algoMask();
221 
222  pquals->push_back((inTrk1.qualityMask() | inTrk2.qualityMask()));
223  pquals->back() |= (1 << reco::TrackBase::confirmed);
224 
225  inputTracks.push_back((*matchIter0)[1]);
226  inputTracks.push_back((*matchIter0)[2]);
227 
228  ++nsel;
229  }
230 
231  producer(mergedTrackSource_, selId);
232  assert(producer.selTracks_->size() == pquals->size());
233 
234  for (auto isel = 0U; isel < nsel; ++isel) {
235  algoMask[isel].set(reco::TrackBase::duplicateMerge);
236  auto& otk = (*producer.selTracks_)[isel];
237  otk.setQualityMask((*pquals)[isel]);
238  otk.setAlgorithm(reco::TrackBase::duplicateMerge);
239  otk.setOriginalAlgorithm(oriAlgo[isel]);
240  otk.setAlgoMask(algoMask[isel]);
241  }
242 
243  selId.clear();
244  for (int i = 0; i < (int)originals.size(); i++) {
245  const reco::Track& origTrack = originals[i];
246  if (std::find(inputTracks.begin(), inputTracks.end(), i) != inputTracks.end())
247  continue;
248  selId.push_back(i);
249  pmvas->push_back((*originalMVAStore)[i]);
250  pquals->push_back(origTrack.qualityMask());
251  }
252 
253  producer(originalTrackSource_, selId);
254  assert(producer.selTracks_->size() == pquals->size());
255 
256  iEvent.put(std::move(pmvas), "MVAValues");
257  iEvent.put(std::move(pquals), "QualityMasks");
258  }
259 
260 } // namespace
261 
264 
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:97
AlgoMask algoMask() const
Definition: TrackBase.h:415
static void fill(edm::ParameterSetDescription &desc)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
TrackAlgorithm
track algorithm
Definition: TrackBase.h:89
int iEvent
Definition: GenABIO.cc:224
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
int qualityMask() const
Definition: TrackBase.h:843
std::bitset< algoSize > AlgoMask
algo mask
Definition: TrackBase.h:145
TrackAlgorithm algo() const
Definition: TrackBase.h:547
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
constexpr std::array< unsigned int, reco::TrackBase::algoSize > trackAlgoPriorityOrder
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
def move(src, dest)
Definition: eostools.py:511