CMS 3D CMS Logo

DuplicateListMerger.cc
Go to the documentation of this file.
1 
10 
11 
16 
20 
21 
24 
25 
27 
30 
31 #include <vector>
32 #include <algorithm>
33 #include <string>
34 #include <iostream>
35 #include<memory>
36 
37 // #include "TMVA/Reader.h"
38 
39 using namespace reco;
40 
41 namespace {
42 class DuplicateListMerger final : public edm::global::EDProducer<> {
43  public:
45  explicit DuplicateListMerger(const edm::ParameterSet& iPara);
47  ~DuplicateListMerger() override;
48 
50  using CandidateToDuplicate = std::vector<std::pair<int, int>>;
51 
52  using RecHitContainer = edm::OwnVector<TrackingRecHit>;
53 
54  using MVACollection = std::vector<float>;
55  using QualityMaskCollection = std::vector<unsigned char>;
56 
57 
58  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
59 
60 
61  private:
63  void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
64 
65  private:
66 
67 
68  TrackCollectionCloner collectionCloner;
69  TrackCollectionCloner::Tokens mergedTrackSource_;
70  TrackCollectionCloner::Tokens originalTrackSource_;
71 
72 
73  edm::EDGetTokenT<CandidateToDuplicate> candidateComponents_;
75 
76  edm::EDGetTokenT<MVACollection> originalMVAValsToken_;
77  edm::EDGetTokenT<MVACollection> mergedMVAValsToken_;
78 
79  std::string priorityName_;
80 
81  int diffHitsCut_;
82 
83 
84 };
85 }
86 
91 
94 
95 #include<tuple>
96 #include<array>
98 
99 
100 namespace {
103  desc.add<edm::InputTag>("mergedSource",edm::InputTag());
104  desc.add<edm::InputTag>("originalSource",edm::InputTag());
105  desc.add<edm::InputTag>("mergedMVAVals",edm::InputTag());
106  desc.add<edm::InputTag>("originalMVAVals",edm::InputTag());
107  desc.add<edm::InputTag>("candidateSource",edm::InputTag());
108  desc.add<edm::InputTag>("candidateComponents",edm::InputTag());
109  desc.add<std::string>("trackAlgoPriorityOrder", "trackAlgoPriorityOrder");
110  desc.add<int>("diffHitsCut",5);
112  descriptions.add("DuplicateListMerger", desc);
113 }
114 
115 DuplicateListMerger::DuplicateListMerger(const edm::ParameterSet& iPara) :
116  collectionCloner(*this, iPara, true),
117  mergedTrackSource_(iPara.getParameter<edm::InputTag>("mergedSource"),consumesCollector()),
118  originalTrackSource_(iPara.getParameter<edm::InputTag>("originalSource"),consumesCollector()),
119  priorityName_(iPara.getParameter<std::string>("trackAlgoPriorityOrder"))
120 {
121 
122  diffHitsCut_ = iPara.getParameter<int>("diffHitsCut");
123  candidateSource_ = consumes<std::vector<TrackCandidate> >(iPara.getParameter<edm::InputTag>("candidateSource"));
124  candidateComponents_ = consumes<CandidateToDuplicate>(iPara.getParameter<edm::InputTag>("candidateComponents"));
125 
126  mergedMVAValsToken_ = consumes<MVACollection>(iPara.getParameter<edm::InputTag>("mergedMVAVals"));
127  originalMVAValsToken_ = consumes<MVACollection>(iPara.getParameter<edm::InputTag>("originalMVAVals"));
128 
129 
130 
131  produces<MVACollection>("MVAValues");
132  produces<QualityMaskCollection>("QualityMasks");
133 
134 
135 }
136 
137 DuplicateListMerger::~DuplicateListMerger()
138 {
139 
140  /* no op */
141 
142 }
143 
144 void DuplicateListMerger::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const
145 {
146 
147  TrackCollectionCloner::Producer producer(iEvent, collectionCloner);
148 
149 
150 
151  auto const & originals = originalTrackSource_.tracks(iEvent);
152  auto const & merged = mergedTrackSource_.tracks(iEvent);
153  auto const & candIndices = mergedTrackSource_.indicesInput(iEvent);
154 
156  iEvent.getByToken(candidateSource_,candidateH);
157  auto const & candidates = *candidateH;
158 
159  edm::Handle<CandidateToDuplicate> candidateComponentsH;
160  iEvent.getByToken(candidateComponents_,candidateComponentsH);
161  auto const & candidateComponents = *candidateComponentsH;
162 
163 
164 
165  edm::Handle<MVACollection> originalMVAStore;
166  edm::Handle<MVACollection> mergedMVAStore;
167 
168  iEvent.getByToken(originalMVAValsToken_,originalMVAStore);
169  iEvent.getByToken(mergedMVAValsToken_,mergedMVAStore);
170 
172  iSetup.get<CkfComponentsRecord>().get(priorityName_, priorityH);
173  auto const & trackAlgoPriorityOrder = *priorityH;
174 
175  MVACollection mvaVec;
176 
177 
178  auto mergedMVA = *mergedMVAStore;
179 
180  //match new tracks to their candidates
181  std::vector<std::array<int,3>> matches;
182  for(int i = 0; i < (int)merged.size(); ++i) {
183  auto cInd = candIndices[i];
184  auto const & cand = candidates[cInd];
185  const reco::Track& matchedTrack = merged[i];
186 
187  if (mergedMVA[i]< -0.7f) continue; // at least "loose" ( FIXME: take cut value from CutSelector)
188 
189  // if( ChiSquaredProbability(matchedTrack.chi2(),matchedTrack.ndof()) < minTrkProbCut_)continue;
190  int dHits = (cand.recHits().second - cand.recHits().first) - matchedTrack.recHitsSize();
191  if(dHits > diffHitsCut_)continue;
192  matches.push_back(std::array<int,3>{{i,candidateComponents[cInd].first,candidateComponents[cInd].second}});
193  }
194 
195  //check for candidates/tracks that share merged tracks, select minimum chi2, remove the rest
196  if (matches.size()>1)
197  for ( auto matchIter0 = matches.begin(); matchIter0 != matches.end()-1; ++matchIter0) {
198  if ( (*matchIter0)[0]<0) continue;
199  auto nchi2 = merged[(*matchIter0)[0]].normalizedChi2();
200  for( auto matchIter1 = matchIter0+1; matchIter1 != matches.end(); ++matchIter1){
201  if ( (*matchIter1)[0]<0) continue;
202  if ( (*matchIter0)[1]==(*matchIter1)[1]
203  || (*matchIter0)[1]==(*matchIter1)[2]
204  || (*matchIter0)[2]==(*matchIter1)[1]
205  || (*matchIter0)[2]==(*matchIter1)[2]
206  ) {
207  auto nchi2_1 = merged[(*matchIter1)[0]].normalizedChi2();
208  if(nchi2_1 < nchi2){
209  (*matchIter0)[0] = -1;
210  break;
211  }else{
212  (*matchIter1)[0] = -1;
213  }
214  }
215  }
216  }
217 
218 
219  // products
220  auto pmvas = std::make_unique<MVACollection>();
221  auto pquals = std::make_unique<QualityMaskCollection>();
222 
223 
224  //add the good merged tracks to the output list, remove input tracks
225  std::vector<int> inputTracks;
226 
227 
228  std::vector<unsigned int> selId;
229  auto ntotTk = matches.size();
230  // declareDynArray(reco::TrackBase::TrackAlgorithm, ntotTk, algo);
232  declareDynArray(reco::TrackBase::AlgoMask, ntotTk, algoMask);
233 
234  auto nsel=0U;
235  for(auto matchIter0 = matches.begin(); matchIter0 != matches.end(); matchIter0++){
236  if ( (*matchIter0)[0]<0) continue;
237  selId.push_back((*matchIter0)[0]);
238 
239  pmvas->push_back(mergedMVA[(*matchIter0)[0]]);
240 
241  const reco::Track& inTrk1 = originals[(*matchIter0)[1]];
242  const reco::Track& inTrk2 = originals[(*matchIter0)[2]];
243  oriAlgo[nsel] = std::min(inTrk1.algo(),inTrk2.algo(),
245  return trackAlgoPriorityOrder.priority(a) < trackAlgoPriorityOrder.priority(b);
246  });
247 
248 
249  algoMask[nsel] = inTrk1.algoMask() | inTrk2.algoMask();
250 
251  pquals->push_back((inTrk1.qualityMask() | inTrk2.qualityMask()));
252  pquals->back() |= (1<<reco::TrackBase::confirmed);
253 
254  inputTracks.push_back((*matchIter0)[1]);
255  inputTracks.push_back((*matchIter0)[2]);
256 
257 
258  ++nsel;
259 
260  }
261 
262 
263  producer(mergedTrackSource_,selId);
264  assert(producer.selTracks_->size()==pquals->size());
265 
266  for (auto isel=0U;isel<nsel;++isel) {
267  algoMask[isel].set(reco::TrackBase::duplicateMerge);
268  auto & otk = (*producer.selTracks_)[isel];
269  otk.setQualityMask((*pquals)[isel]);
270  otk.setAlgorithm(reco::TrackBase::duplicateMerge);
271  otk.setOriginalAlgorithm(oriAlgo[isel]);
272  otk.setAlgoMask(algoMask[isel]);
273  }
274 
275 
276  selId.clear();
277  for(int i = 0; i < (int)originals.size(); i++){
278  const reco::Track& origTrack = originals[i];
279  if (std::find(inputTracks.begin(),inputTracks.end(),i) != inputTracks.end()) continue;
280  selId.push_back(i);
281  pmvas->push_back((*originalMVAStore)[i]);
282  pquals->push_back(origTrack.qualityMask());
283  }
284 
285  producer(originalTrackSource_,selId);
286  assert(producer.selTracks_->size()==pquals->size());
287 
288  iEvent.put(std::move(pmvas),"MVAValues");
289  iEvent.put(std::move(pquals),"QualityMasks");
290 
291 
292 }
293 
294 }
295 
298 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:119
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:20
TrackAlgorithm
track algorithm
Definition: TrackBase.h:99
TrackAlgorithm algo() const
Definition: TrackBase.h:497
int iEvent
Definition: GenABIO.cc:230
int qualityMask() const
Definition: TrackBase.h:862
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(ConfigurationDescriptions &descriptions)
std::bitset< algoSize > AlgoMask
algo mask
Definition: TrackBase.h:146
AlgoMask algoMask() const
Definition: TrackBase.h:364
const T & get() const
Definition: EventSetup.h:55
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
std::array< unsigned int, reco::TrackBase::algoSize > trackAlgoPriorityOrder
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
def move(src, dest)
Definition: eostools.py:510