CMS 3D CMS Logo

ConversionTrackMerger.cc
Go to the documentation of this file.
1 //
2 // Package: RecoTracker/FinalTrackSelectors
3 // Class: ConversionTrackMerger
4 //
5 // Description: Merger for ConversionTracks, adapted from SimpleTrackListMerger
6 //
7 // Original Author: J.Bendavid
8 //
9 //
10 
29 
30 #include <vector>
31 
33 public:
34  explicit ConversionTrackMerger(const edm::ParameterSet& conf);
35 
36  ~ConversionTrackMerger() override;
37 
38  void produce(edm::Event& e, const edm::EventSetup& c) override;
39 
40 private:
42 
45 
46  std::unique_ptr<reco::ConversionTrackCollection> outputTrks;
47 };
48 
51 
53  // retrieve producer name of input TrackCollection(s)
54  trackProducer1 = consumes<reco::ConversionTrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer1"));
55  trackProducer2 = consumes<reco::ConversionTrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer2"));
56 
57  produces<reco::ConversionTrackCollection>();
58 }
59 
60 // Virtual destructor needed.
62 
63 // Functions that gets called by framework every event
65  double shareFrac = conf_.getParameter<double>("ShareFrac");
66  bool allowFirstHitShare = conf_.getParameter<bool>("allowFirstHitShare");
67  bool checkCharge = conf_.getParameter<bool>("checkCharge");
68  double minProb = conf_.getParameter<double>("minProb");
69 
70  int outputPreferCollection = conf_.getParameter<int>("outputPreferCollection");
71  int trackerOnlyPreferCollection = conf_.getParameter<int>("trackerOnlyPreferCollection");
72  int arbitratedEcalSeededPreferCollection = conf_.getParameter<int>("arbitratedEcalSeededPreferCollection");
73  int arbitratedMergedPreferCollection = conf_.getParameter<int>("arbitratedMergedPreferCollection");
75  conf_.getParameter<int>("arbitratedMergedEcalGeneralPreferCollection");
76 
77  // get Inputs
78  // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
79  // this allows ConversionTrackMerger to be used as a cleaner only if handed just one list
80  // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
81  //
82  const reco::ConversionTrackCollection* TC1 = nullptr;
83  static const reco::ConversionTrackCollection s_empty1, s_empty2;
85  e.getByToken(trackProducer1, trackCollection1);
86  if (trackCollection1.isValid()) {
87  TC1 = trackCollection1.product();
88  //std::cout << "1st collection " << trackProducer1 << " has "<< TC1->size() << " tracks" << std::endl ;
89  } else {
90  TC1 = &s_empty1;
91  edm::LogWarning("ConversionTrackMerger") << "1st TrackCollection not found;"
92  << " will only clean 2nd TrackCollection ";
93  }
95 
96  const reco::ConversionTrackCollection* TC2 = nullptr;
98  e.getByToken(trackProducer2, trackCollection2);
99  if (trackCollection2.isValid()) {
100  TC2 = trackCollection2.product();
101  //std::cout << "2nd collection " << trackProducer2 << " has "<< TC2->size() << " tracks" << std::endl ;
102  } else {
103  TC2 = &s_empty2;
104  edm::LogWarning("ConversionTrackMerger") << "2nd TrackCollection not found;"
105  << " will only clean 1st TrackCollection ";
106  }
108 
109  // Step B: create empty output collection
110  outputTrks = std::make_unique<reco::ConversionTrackCollection>();
111  int i;
112 
113  std::vector<int> selected1;
114  for (unsigned int i = 0; i < tC1.size(); ++i) {
115  selected1.push_back(1);
116  }
117  std::vector<int> selected2;
118  for (unsigned int i = 0; i < tC2.size(); ++i) {
119  selected2.push_back(1);
120  }
121 
122  std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh1;
123  std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh2;
124  for (reco::ConversionTrackCollection::const_iterator track = tC1.begin(); track != tC1.end(); ++track) {
125  trackingRecHit_iterator itB = track->track()->recHitsBegin();
126  trackingRecHit_iterator itE = track->track()->recHitsEnd();
127  for (trackingRecHit_iterator it = itB; it != itE; ++it) {
128  const TrackingRecHit* hit = &(**it);
129  rh1[track].push_back(hit);
130  }
131  }
132  for (reco::ConversionTrackCollection::const_iterator track = tC2.begin(); track != tC2.end(); ++track) {
133  trackingRecHit_iterator jtB = track->track()->recHitsBegin();
134  trackingRecHit_iterator jtE = track->track()->recHitsEnd();
135  for (trackingRecHit_iterator jt = jtB; jt != jtE; ++jt) {
136  const TrackingRecHit* hit = &(**jt);
137  rh2[track].push_back(hit);
138  }
139  }
140 
141  if ((!tC1.empty()) && (!tC2.empty())) {
142  i = -1;
143  for (reco::ConversionTrackCollection::iterator track = tC1.begin(); track != tC1.end(); ++track) {
144  i++;
145 
146  //clear flags if preferCollection was set to 0
147  selected1[i] = selected1[i] && outputPreferCollection != 0;
148  track->setIsTrackerOnly(track->isTrackerOnly() && trackerOnlyPreferCollection != 0);
149  track->setIsArbitratedEcalSeeded(track->isArbitratedEcalSeeded() && arbitratedEcalSeededPreferCollection != 0);
150  track->setIsArbitratedMerged(track->isArbitratedMerged() && arbitratedMergedPreferCollection != 0);
151  track->setIsArbitratedMergedEcalGeneral(track->isArbitratedMergedEcalGeneral() &&
153 
154  std::vector<const TrackingRecHit*>& iHits = rh1[track];
155  unsigned nh1 = iHits.size();
156  int j = -1;
157  for (reco::ConversionTrackCollection::iterator track2 = tC2.begin(); track2 != tC2.end(); ++track2) {
158  j++;
159 
160  //clear flags if preferCollection was set to 0
161  selected2[j] = selected2[j] && outputPreferCollection != 0;
162  track2->setIsTrackerOnly(track2->isTrackerOnly() && trackerOnlyPreferCollection != 0);
163  track2->setIsArbitratedEcalSeeded(track2->isArbitratedEcalSeeded() &&
165  track2->setIsArbitratedMerged(track2->isArbitratedMerged() && arbitratedMergedPreferCollection != 0);
166  track2->setIsArbitratedMergedEcalGeneral(track2->isArbitratedMergedEcalGeneral() &&
168 
169  std::vector<const TrackingRecHit*>& jHits = rh2[track2];
170  unsigned nh2 = jHits.size();
171  int noverlap = 0;
172  int firstoverlap = 0;
173  for (unsigned ih = 0; ih < nh1; ++ih) {
174  const TrackingRecHit* it = iHits[ih];
175  if (it->isValid()) {
176  for (unsigned jh = 0; jh < nh2; ++jh) {
177  const TrackingRecHit* jt = jHits[jh];
178  if (jt->isValid()) {
179  if (it->sharesInput(jt, TrackingRecHit::some)) {
180  noverlap++;
181  if (allowFirstHitShare && (ih == 0) && (jh == 0))
182  firstoverlap = 1;
183  }
184  }
185  }
186  }
187  }
188  int nhit1 = track->track()->numberOfValidHits();
189  int nhit2 = track2->track()->numberOfValidHits();
190  //if (noverlap>0) printf("noverlap = %i, firstoverlap = %i, nhit1 = %i, nhit2 = %i, algo1 = %i, algo2 = %i, q1 = %i, q2 = %i\n",noverlap,firstoverlap,nhit1,nhit2,track->track()->algo(),track2->track()->algo(),track->track()->charge(),track2->track()->charge());
191  //std::cout << " trk1 trk2 nhits1 nhits2 nover " << i << " " << j << " " << track->numberOfValidHits() << " " << track2->numberOfValidHits() << " " << noverlap << " " << fi << " " << fj <<std::endl;
192  if ((noverlap - firstoverlap) > (std::min(nhit1, nhit2) - firstoverlap) * shareFrac &&
193  (!checkCharge || track->track()->charge() * track2->track()->charge() > 0)) {
194  //printf("overlapping tracks\n");
195  //printf ("ndof1 = %5f, chisq1 = %5f, ndof2 = %5f, chisq2 = %5f\n",track->track()->ndof(),track->track()->chi2(),track2->track()->ndof(),track2->track()->chi2());
196 
197  double probFirst = ChiSquaredProbability(track->track()->chi2(), track->track()->ndof());
198  double probSecond = ChiSquaredProbability(track2->track()->chi2(), track2->track()->ndof());
199 
200  //arbitrate by number of hits and reduced chisq
201  bool keepFirst = (nhit1 > nhit2 ||
202  (nhit1 == nhit2 && track->track()->normalizedChi2() < track2->track()->normalizedChi2()));
203 
204  //override decision in case one track is radically worse quality than the other
205  keepFirst |= (probFirst > minProb && probSecond <= minProb);
206  keepFirst &= !(probFirst <= minProb && probSecond > minProb);
207 
208  bool keepSecond = !keepFirst;
209 
210  //set flags based on arbitration decision and precedence settings
211 
212  selected1[i] = selected1[i] && ((keepFirst && outputPreferCollection == 3) || outputPreferCollection == -1 ||
214  track->setIsTrackerOnly(track->isTrackerOnly() &&
215  ((keepFirst && trackerOnlyPreferCollection == 3) ||
217  track->setIsArbitratedEcalSeeded(track->isArbitratedEcalSeeded() &&
218  ((keepFirst && arbitratedEcalSeededPreferCollection == 3) ||
221  track->setIsArbitratedMerged(track->isArbitratedMerged() &&
222  ((keepFirst && arbitratedMergedPreferCollection == 3) ||
225  track->setIsArbitratedMergedEcalGeneral(track->isArbitratedMergedEcalGeneral() &&
226  ((keepFirst && arbitratedMergedEcalGeneralPreferCollection == 3) ||
229 
230  selected2[j] = selected2[j] && ((keepSecond && outputPreferCollection == 3) || outputPreferCollection == -1 ||
232  track2->setIsTrackerOnly(track2->isTrackerOnly() &&
233  ((keepSecond && trackerOnlyPreferCollection == 3) ||
235  track2->setIsArbitratedEcalSeeded(track2->isArbitratedEcalSeeded() &&
236  ((keepSecond && arbitratedEcalSeededPreferCollection == 3) ||
239  track2->setIsArbitratedMerged(track2->isArbitratedMerged() &&
240  ((keepSecond && arbitratedMergedPreferCollection == 3) ||
243  track2->setIsArbitratedMergedEcalGeneral(track2->isArbitratedMergedEcalGeneral() &&
244  ((keepSecond && arbitratedMergedEcalGeneralPreferCollection == 3) ||
247 
248  } //end got a duplicate
249  } //end track2 loop
250  } //end track loop
251  } //end more than 1 track
252 
253  //
254  // output selected tracks - if any
255  //
256 
257  if (!tC1.empty()) {
258  i = 0;
259  for (reco::ConversionTrackCollection::const_iterator track = tC1.begin(); track != tC1.end(); ++track, ++i) {
260  //don't store tracks rejected as duplicates
261  if (!selected1[i]) {
262  continue;
263  }
264  //fill the TrackCollection
265  outputTrks->push_back(*track);
266  } //end faux loop over tracks
267  } //end more than 0 track
268 
269  //Fill the trajectories, etc. for 1st collection
270 
271  if (!tC2.empty()) {
272  i = 0;
273  for (reco::ConversionTrackCollection::const_iterator track = tC2.begin(); track != tC2.end(); ++track, ++i) {
274  //don't store tracks rejected as duplicates
275  if (!selected2[i]) {
276  continue;
277  }
278  //fill the TrackCollection
279  outputTrks->push_back(*track);
280  } //end faux loop over tracks
281  } //end more than 0 track
282 
283  e.put(std::move(outputTrks));
284  return;
285 
286 } //end produce
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T const * product() const
Definition: Handle.h:70
edm::EDGetTokenT< reco::ConversionTrackCollection > trackProducer1
std::vector< ConversionTrack > ConversionTrackCollection
collection of ConversionTracks
bool isValid() const
std::unique_ptr< reco::ConversionTrackCollection > outputTrks
float ChiSquaredProbability(double chiSquared, double nrDOF)
ConversionTrackMerger(const edm::ParameterSet &conf)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool isValid() const
Definition: HandleBase.h:70
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511
void produce(edm::Event &e, const edm::EventSetup &c) override
edm::EDGetTokenT< reco::ConversionTrackCollection > trackProducer2