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 
11 #include <memory>
12 #include <string>
13 #include <iostream>
14 #include <cmath>
15 #include <vector>
16 
18 
24 
26 
28 
32 
33 //#include "DataFormats/TrackReco/src/classes.h"
34 
36 
38  // retrieve producer name of input TrackCollection(s)
39  trackProducer1 = consumes<reco::ConversionTrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer1"));
40  trackProducer2 = consumes<reco::ConversionTrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer2"));
41 
42  produces<reco::ConversionTrackCollection>();
43 }
44 
45 // Virtual destructor needed.
47 
48 // Functions that gets called by framework every event
50  double shareFrac = conf_.getParameter<double>("ShareFrac");
51  bool allowFirstHitShare = conf_.getParameter<bool>("allowFirstHitShare");
52  bool checkCharge = conf_.getParameter<bool>("checkCharge");
53  double minProb = conf_.getParameter<double>("minProb");
54 
55  int outputPreferCollection = conf_.getParameter<int>("outputPreferCollection");
56  int trackerOnlyPreferCollection = conf_.getParameter<int>("trackerOnlyPreferCollection");
57  int arbitratedEcalSeededPreferCollection = conf_.getParameter<int>("arbitratedEcalSeededPreferCollection");
58  int arbitratedMergedPreferCollection = conf_.getParameter<int>("arbitratedMergedPreferCollection");
60  conf_.getParameter<int>("arbitratedMergedEcalGeneralPreferCollection");
61 
62  // get Inputs
63  // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
64  // this allows ConversionTrackMerger to be used as a cleaner only if handed just one list
65  // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
66  //
67  const reco::ConversionTrackCollection* TC1 = nullptr;
68  static const reco::ConversionTrackCollection s_empty1, s_empty2;
70  e.getByToken(trackProducer1, trackCollection1);
71  if (trackCollection1.isValid()) {
72  TC1 = trackCollection1.product();
73  //std::cout << "1st collection " << trackProducer1 << " has "<< TC1->size() << " tracks" << std::endl ;
74  } else {
75  TC1 = &s_empty1;
76  edm::LogWarning("ConversionTrackMerger") << "1st TrackCollection not found;"
77  << " will only clean 2nd TrackCollection ";
78  }
80 
81  const reco::ConversionTrackCollection* TC2 = nullptr;
83  e.getByToken(trackProducer2, trackCollection2);
84  if (trackCollection2.isValid()) {
85  TC2 = trackCollection2.product();
86  //std::cout << "2nd collection " << trackProducer2 << " has "<< TC2->size() << " tracks" << std::endl ;
87  } else {
88  TC2 = &s_empty2;
89  edm::LogWarning("ConversionTrackMerger") << "2nd TrackCollection not found;"
90  << " will only clean 1st TrackCollection ";
91  }
93 
94  // Step B: create empty output collection
95  outputTrks = std::make_unique<reco::ConversionTrackCollection>();
96  int i;
97 
98  std::vector<int> selected1;
99  for (unsigned int i = 0; i < tC1.size(); ++i) {
100  selected1.push_back(1);
101  }
102  std::vector<int> selected2;
103  for (unsigned int i = 0; i < tC2.size(); ++i) {
104  selected2.push_back(1);
105  }
106 
107  std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh1;
108  std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh2;
109  for (reco::ConversionTrackCollection::const_iterator track = tC1.begin(); track != tC1.end(); ++track) {
110  trackingRecHit_iterator itB = track->track()->recHitsBegin();
111  trackingRecHit_iterator itE = track->track()->recHitsEnd();
112  for (trackingRecHit_iterator it = itB; it != itE; ++it) {
113  const TrackingRecHit* hit = &(**it);
114  rh1[track].push_back(hit);
115  }
116  }
117  for (reco::ConversionTrackCollection::const_iterator track = tC2.begin(); track != tC2.end(); ++track) {
118  trackingRecHit_iterator jtB = track->track()->recHitsBegin();
119  trackingRecHit_iterator jtE = track->track()->recHitsEnd();
120  for (trackingRecHit_iterator jt = jtB; jt != jtE; ++jt) {
121  const TrackingRecHit* hit = &(**jt);
122  rh2[track].push_back(hit);
123  }
124  }
125 
126  if ((!tC1.empty()) && (!tC2.empty())) {
127  i = -1;
128  for (reco::ConversionTrackCollection::iterator track = tC1.begin(); track != tC1.end(); ++track) {
129  i++;
130 
131  //clear flags if preferCollection was set to 0
132  selected1[i] = selected1[i] && outputPreferCollection != 0;
133  track->setIsTrackerOnly(track->isTrackerOnly() && trackerOnlyPreferCollection != 0);
134  track->setIsArbitratedEcalSeeded(track->isArbitratedEcalSeeded() && arbitratedEcalSeededPreferCollection != 0);
135  track->setIsArbitratedMerged(track->isArbitratedMerged() && arbitratedMergedPreferCollection != 0);
136  track->setIsArbitratedMergedEcalGeneral(track->isArbitratedMergedEcalGeneral() &&
138 
139  std::vector<const TrackingRecHit*>& iHits = rh1[track];
140  unsigned nh1 = iHits.size();
141  int j = -1;
142  for (reco::ConversionTrackCollection::iterator track2 = tC2.begin(); track2 != tC2.end(); ++track2) {
143  j++;
144 
145  //clear flags if preferCollection was set to 0
146  selected2[j] = selected2[j] && outputPreferCollection != 0;
147  track2->setIsTrackerOnly(track2->isTrackerOnly() && trackerOnlyPreferCollection != 0);
148  track2->setIsArbitratedEcalSeeded(track2->isArbitratedEcalSeeded() &&
150  track2->setIsArbitratedMerged(track2->isArbitratedMerged() && arbitratedMergedPreferCollection != 0);
151  track2->setIsArbitratedMergedEcalGeneral(track2->isArbitratedMergedEcalGeneral() &&
153 
154  std::vector<const TrackingRecHit*>& jHits = rh2[track2];
155  unsigned nh2 = jHits.size();
156  int noverlap = 0;
157  int firstoverlap = 0;
158  for (unsigned ih = 0; ih < nh1; ++ih) {
159  const TrackingRecHit* it = iHits[ih];
160  if (it->isValid()) {
161  int jj = -1;
162  for (unsigned jh = 0; jh < nh2; ++jh) {
163  const TrackingRecHit* jt = jHits[jh];
164  jj++;
165  if (jt->isValid()) {
166  if (it->sharesInput(jt, TrackingRecHit::some)) {
167  noverlap++;
168  if (allowFirstHitShare && (ih == 0) && (jh == 0))
169  firstoverlap = 1;
170  }
171  }
172  }
173  }
174  }
175  int nhit1 = track->track()->numberOfValidHits();
176  int nhit2 = track2->track()->numberOfValidHits();
177  //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());
178  //std::cout << " trk1 trk2 nhits1 nhits2 nover " << i << " " << j << " " << track->numberOfValidHits() << " " << track2->numberOfValidHits() << " " << noverlap << " " << fi << " " << fj <<std::endl;
179  if ((noverlap - firstoverlap) > (std::min(nhit1, nhit2) - firstoverlap) * shareFrac &&
180  (!checkCharge || track->track()->charge() * track2->track()->charge() > 0)) {
181  //printf("overlapping tracks\n");
182  //printf ("ndof1 = %5f, chisq1 = %5f, ndof2 = %5f, chisq2 = %5f\n",track->track()->ndof(),track->track()->chi2(),track2->track()->ndof(),track2->track()->chi2());
183 
184  double probFirst = ChiSquaredProbability(track->track()->chi2(), track->track()->ndof());
185  double probSecond = ChiSquaredProbability(track2->track()->chi2(), track2->track()->ndof());
186 
187  //arbitrate by number of hits and reduced chisq
188  bool keepFirst = (nhit1 > nhit2 ||
189  (nhit1 == nhit2 && track->track()->normalizedChi2() < track2->track()->normalizedChi2()));
190 
191  //override decision in case one track is radically worse quality than the other
192  keepFirst |= (probFirst > minProb && probSecond <= minProb);
193  keepFirst &= !(probFirst <= minProb && probSecond > minProb);
194 
195  bool keepSecond = !keepFirst;
196 
197  //set flags based on arbitration decision and precedence settings
198 
199  selected1[i] = selected1[i] && ((keepFirst && outputPreferCollection == 3) || outputPreferCollection == -1 ||
201  track->setIsTrackerOnly(track->isTrackerOnly() &&
202  ((keepFirst && trackerOnlyPreferCollection == 3) ||
204  track->setIsArbitratedEcalSeeded(track->isArbitratedEcalSeeded() &&
205  ((keepFirst && arbitratedEcalSeededPreferCollection == 3) ||
208  track->setIsArbitratedMerged(track->isArbitratedMerged() &&
209  ((keepFirst && arbitratedMergedPreferCollection == 3) ||
212  track->setIsArbitratedMergedEcalGeneral(track->isArbitratedMergedEcalGeneral() &&
213  ((keepFirst && arbitratedMergedEcalGeneralPreferCollection == 3) ||
216 
217  selected2[j] = selected2[j] && ((keepSecond && outputPreferCollection == 3) || outputPreferCollection == -1 ||
219  track2->setIsTrackerOnly(track2->isTrackerOnly() &&
220  ((keepSecond && trackerOnlyPreferCollection == 3) ||
222  track2->setIsArbitratedEcalSeeded(track2->isArbitratedEcalSeeded() &&
223  ((keepSecond && arbitratedEcalSeededPreferCollection == 3) ||
226  track2->setIsArbitratedMerged(track2->isArbitratedMerged() &&
227  ((keepSecond && arbitratedMergedPreferCollection == 3) ||
230  track2->setIsArbitratedMergedEcalGeneral(track2->isArbitratedMergedEcalGeneral() &&
231  ((keepSecond && arbitratedMergedEcalGeneralPreferCollection == 3) ||
234 
235  } //end got a duplicate
236  } //end track2 loop
237  } //end track loop
238  } //end more than 1 track
239 
240  //
241  // output selected tracks - if any
242  //
243 
244  if (!tC1.empty()) {
245  i = 0;
246  for (reco::ConversionTrackCollection::const_iterator track = tC1.begin(); track != tC1.end(); ++track, ++i) {
247  //don't store tracks rejected as duplicates
248  if (!selected1[i]) {
249  continue;
250  }
251  //fill the TrackCollection
252  outputTrks->push_back(*track);
253  } //end faux loop over tracks
254  } //end more than 0 track
255 
256  //Fill the trajectories, etc. for 1st collection
257 
258  if (!tC2.empty()) {
259  i = 0;
260  for (reco::ConversionTrackCollection::const_iterator track = tC2.begin(); track != tC2.end(); ++track, ++i) {
261  //don't store tracks rejected as duplicates
262  if (!selected2[i]) {
263  continue;
264  }
265  //fill the TrackCollection
266  outputTrks->push_back(*track);
267  } //end faux loop over tracks
268  } //end more than 0 track
269 
270  e.put(std::move(outputTrks));
271  return;
272 
273 } //end produce
PA_ZEESkim_cff.checkCharge
checkCharge
Definition: PA_ZEESkim_cff.py:27
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
TrackerGeometry.h
edm::Handle::product
T const * product() const
Definition: Handle.h:70
TrackCandidateCollection.h
ESHandle.h
SiStripMatchedRecHit2DCollection.h
min
T min(T a, T b)
Definition: MathUtil.h:58
ChiSquaredProbability
float ChiSquaredProbability(double chiSquared, double nrDOF)
Definition: ChiSquaredProbability.cc:13
reco::ConversionTrackCollection
std::vector< ConversionTrack > ConversionTrackCollection
collection of ConversionTracks
Definition: ConversionTrackFwd.h:6
ChiSquaredProbability.h
conversionTrackMerger_cfi.arbitratedEcalSeededPreferCollection
arbitratedEcalSeededPreferCollection
Definition: conversionTrackMerger_cfi.py:19
edm::Handle
Definition: AssociativeIterator.h:50
ConversionTrackMerger::outputTrks
std::unique_ptr< reco::ConversionTrackCollection > outputTrks
Definition: ConversionTrackMerger.h:43
ConversionTrackMerger::~ConversionTrackMerger
~ConversionTrackMerger() override
Definition: ConversionTrackMerger.cc:46
TrackingRecHit::some
Definition: TrackingRecHit.h:59
HLT_FULL_cff.shareFrac
shareFrac
Definition: HLT_FULL_cff.py:84042
edm::OwnVector::const_iterator
Definition: OwnVector.h:41
edm::LogWarning
Definition: MessageLogger.h:141
ConversionTrackMerger::produce
void produce(edm::Event &e, const edm::EventSetup &c) override
Definition: ConversionTrackMerger.cc:49
TrackerDigiGeometryRecord.h
edm::ParameterSet
Definition: ParameterSet.h:36
ConversionTrackMerger::trackProducer2
edm::EDGetTokenT< reco::ConversionTrackCollection > trackProducer2
Definition: ConversionTrackMerger.h:41
SiStripRecHit1DCollection.h
SiStripRecHit2DCollection.h
conversionTrackMerger_cfi.minProb
minProb
Definition: conversionTrackMerger_cfi.py:10
edm::EventSetup
Definition: EventSetup.h:57
TrajectorySeedCollection.h
HLT_2018_cff.allowFirstHitShare
allowFirstHitShare
Definition: HLT_2018_cff.py:8871
TrackingRecHit
Definition: TrackingRecHit.h:21
ConversionTrackMerger::ConversionTrackMerger
ConversionTrackMerger(const edm::ParameterSet &conf)
Definition: ConversionTrackMerger.cc:37
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
GeomDet.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
conversionTrackMerger_cfi.arbitratedMergedEcalGeneralPreferCollection
arbitratedMergedEcalGeneralPreferCollection
Definition: conversionTrackMerger_cfi.py:21
conversionTrackMerger_cfi.outputPreferCollection
outputPreferCollection
Definition: conversionTrackMerger_cfi.py:17
findQualityFiles.jj
string jj
Definition: findQualityFiles.py:188
TrackingRecHit::sharesInput
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
Definition: TrackingRecHit.cc:12
conversionTrackMerger_cfi.trackerOnlyPreferCollection
trackerOnlyPreferCollection
Definition: conversionTrackMerger_cfi.py:18
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
ConversionTrackMerger::conf_
edm::ParameterSet conf_
Definition: ConversionTrackMerger.h:38
edm::Event
Definition: Event.h:73
conversionTrackMerger_cfi.arbitratedMergedPreferCollection
arbitratedMergedPreferCollection
Definition: conversionTrackMerger_cfi.py:20
ConversionTrackMerger.h
TrackingRecHit::isValid
bool isValid() const
Definition: TrackingRecHit.h:141
edm::InputTag
Definition: InputTag.h:15
hit
Definition: SiStripHitEffFromCalibTree.cc:88
ConversionTrackMerger::trackProducer1
edm::EDGetTokenT< reco::ConversionTrackCollection > trackProducer1
Definition: ConversionTrackMerger.h:40
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37