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 
34 
36  // retrieve producer name of input TrackCollection(s)
37  trackProducer1 = consumes<reco::ConversionTrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer1"));
38  trackProducer2 = consumes<reco::ConversionTrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer2"));
39 
40  produces<reco::ConversionTrackCollection>();
41 }
42 
43 // Virtual destructor needed.
45 
46 // Functions that gets called by framework every event
48  double shareFrac = conf_.getParameter<double>("ShareFrac");
49  bool allowFirstHitShare = conf_.getParameter<bool>("allowFirstHitShare");
50  bool checkCharge = conf_.getParameter<bool>("checkCharge");
51  double minProb = conf_.getParameter<double>("minProb");
52 
53  int outputPreferCollection = conf_.getParameter<int>("outputPreferCollection");
54  int trackerOnlyPreferCollection = conf_.getParameter<int>("trackerOnlyPreferCollection");
55  int arbitratedEcalSeededPreferCollection = conf_.getParameter<int>("arbitratedEcalSeededPreferCollection");
56  int arbitratedMergedPreferCollection = conf_.getParameter<int>("arbitratedMergedPreferCollection");
58  conf_.getParameter<int>("arbitratedMergedEcalGeneralPreferCollection");
59 
60  // get Inputs
61  // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
62  // this allows ConversionTrackMerger to be used as a cleaner only if handed just one list
63  // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
64  //
65  const reco::ConversionTrackCollection* TC1 = nullptr;
66  static const reco::ConversionTrackCollection s_empty1, s_empty2;
68  e.getByToken(trackProducer1, trackCollection1);
69  if (trackCollection1.isValid()) {
70  TC1 = trackCollection1.product();
71  //std::cout << "1st collection " << trackProducer1 << " has "<< TC1->size() << " tracks" << std::endl ;
72  } else {
73  TC1 = &s_empty1;
74  edm::LogWarning("ConversionTrackMerger") << "1st TrackCollection not found;"
75  << " will only clean 2nd TrackCollection ";
76  }
78 
79  const reco::ConversionTrackCollection* TC2 = nullptr;
81  e.getByToken(trackProducer2, trackCollection2);
82  if (trackCollection2.isValid()) {
83  TC2 = trackCollection2.product();
84  //std::cout << "2nd collection " << trackProducer2 << " has "<< TC2->size() << " tracks" << std::endl ;
85  } else {
86  TC2 = &s_empty2;
87  edm::LogWarning("ConversionTrackMerger") << "2nd TrackCollection not found;"
88  << " will only clean 1st TrackCollection ";
89  }
91 
92  // Step B: create empty output collection
93  outputTrks = std::make_unique<reco::ConversionTrackCollection>();
94  int i;
95 
96  std::vector<int> selected1;
97  for (unsigned int i = 0; i < tC1.size(); ++i) {
98  selected1.push_back(1);
99  }
100  std::vector<int> selected2;
101  for (unsigned int i = 0; i < tC2.size(); ++i) {
102  selected2.push_back(1);
103  }
104 
105  std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh1;
106  std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh2;
107  for (reco::ConversionTrackCollection::const_iterator track = tC1.begin(); track != tC1.end(); ++track) {
108  trackingRecHit_iterator itB = track->track()->recHitsBegin();
109  trackingRecHit_iterator itE = track->track()->recHitsEnd();
110  for (trackingRecHit_iterator it = itB; it != itE; ++it) {
111  const TrackingRecHit* hit = &(**it);
112  rh1[track].push_back(hit);
113  }
114  }
115  for (reco::ConversionTrackCollection::const_iterator track = tC2.begin(); track != tC2.end(); ++track) {
116  trackingRecHit_iterator jtB = track->track()->recHitsBegin();
117  trackingRecHit_iterator jtE = track->track()->recHitsEnd();
118  for (trackingRecHit_iterator jt = jtB; jt != jtE; ++jt) {
119  const TrackingRecHit* hit = &(**jt);
120  rh2[track].push_back(hit);
121  }
122  }
123 
124  if ((!tC1.empty()) && (!tC2.empty())) {
125  i = -1;
126  for (reco::ConversionTrackCollection::iterator track = tC1.begin(); track != tC1.end(); ++track) {
127  i++;
128 
129  //clear flags if preferCollection was set to 0
130  selected1[i] = selected1[i] && outputPreferCollection != 0;
131  track->setIsTrackerOnly(track->isTrackerOnly() && trackerOnlyPreferCollection != 0);
132  track->setIsArbitratedEcalSeeded(track->isArbitratedEcalSeeded() && arbitratedEcalSeededPreferCollection != 0);
133  track->setIsArbitratedMerged(track->isArbitratedMerged() && arbitratedMergedPreferCollection != 0);
134  track->setIsArbitratedMergedEcalGeneral(track->isArbitratedMergedEcalGeneral() &&
136 
137  std::vector<const TrackingRecHit*>& iHits = rh1[track];
138  unsigned nh1 = iHits.size();
139  int j = -1;
140  for (reco::ConversionTrackCollection::iterator track2 = tC2.begin(); track2 != tC2.end(); ++track2) {
141  j++;
142 
143  //clear flags if preferCollection was set to 0
144  selected2[j] = selected2[j] && outputPreferCollection != 0;
145  track2->setIsTrackerOnly(track2->isTrackerOnly() && trackerOnlyPreferCollection != 0);
146  track2->setIsArbitratedEcalSeeded(track2->isArbitratedEcalSeeded() &&
148  track2->setIsArbitratedMerged(track2->isArbitratedMerged() && arbitratedMergedPreferCollection != 0);
149  track2->setIsArbitratedMergedEcalGeneral(track2->isArbitratedMergedEcalGeneral() &&
151 
152  std::vector<const TrackingRecHit*>& jHits = rh2[track2];
153  unsigned nh2 = jHits.size();
154  int noverlap = 0;
155  int firstoverlap = 0;
156  for (unsigned ih = 0; ih < nh1; ++ih) {
157  const TrackingRecHit* it = iHits[ih];
158  if (it->isValid()) {
159  int jj = -1;
160  for (unsigned jh = 0; jh < nh2; ++jh) {
161  const TrackingRecHit* jt = jHits[jh];
162  jj++;
163  if (jt->isValid()) {
164  if (it->sharesInput(jt, TrackingRecHit::some)) {
165  noverlap++;
166  if (allowFirstHitShare && (ih == 0) && (jh == 0))
167  firstoverlap = 1;
168  }
169  }
170  }
171  }
172  }
173  int nhit1 = track->track()->numberOfValidHits();
174  int nhit2 = track2->track()->numberOfValidHits();
175  //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());
176  //std::cout << " trk1 trk2 nhits1 nhits2 nover " << i << " " << j << " " << track->numberOfValidHits() << " " << track2->numberOfValidHits() << " " << noverlap << " " << fi << " " << fj <<std::endl;
177  if ((noverlap - firstoverlap) > (std::min(nhit1, nhit2) - firstoverlap) * shareFrac &&
178  (!checkCharge || track->track()->charge() * track2->track()->charge() > 0)) {
179  //printf("overlapping tracks\n");
180  //printf ("ndof1 = %5f, chisq1 = %5f, ndof2 = %5f, chisq2 = %5f\n",track->track()->ndof(),track->track()->chi2(),track2->track()->ndof(),track2->track()->chi2());
181 
182  double probFirst = ChiSquaredProbability(track->track()->chi2(), track->track()->ndof());
183  double probSecond = ChiSquaredProbability(track2->track()->chi2(), track2->track()->ndof());
184 
185  //arbitrate by number of hits and reduced chisq
186  bool keepFirst = (nhit1 > nhit2 ||
187  (nhit1 == nhit2 && track->track()->normalizedChi2() < track2->track()->normalizedChi2()));
188 
189  //override decision in case one track is radically worse quality than the other
190  keepFirst |= (probFirst > minProb && probSecond <= minProb);
191  keepFirst &= !(probFirst <= minProb && probSecond > minProb);
192 
193  bool keepSecond = !keepFirst;
194 
195  //set flags based on arbitration decision and precedence settings
196 
197  selected1[i] = selected1[i] && ((keepFirst && outputPreferCollection == 3) || outputPreferCollection == -1 ||
199  track->setIsTrackerOnly(track->isTrackerOnly() &&
200  ((keepFirst && trackerOnlyPreferCollection == 3) ||
202  track->setIsArbitratedEcalSeeded(track->isArbitratedEcalSeeded() &&
203  ((keepFirst && arbitratedEcalSeededPreferCollection == 3) ||
206  track->setIsArbitratedMerged(track->isArbitratedMerged() &&
207  ((keepFirst && arbitratedMergedPreferCollection == 3) ||
210  track->setIsArbitratedMergedEcalGeneral(track->isArbitratedMergedEcalGeneral() &&
211  ((keepFirst && arbitratedMergedEcalGeneralPreferCollection == 3) ||
214 
215  selected2[j] = selected2[j] && ((keepSecond && outputPreferCollection == 3) || outputPreferCollection == -1 ||
217  track2->setIsTrackerOnly(track2->isTrackerOnly() &&
218  ((keepSecond && trackerOnlyPreferCollection == 3) ||
220  track2->setIsArbitratedEcalSeeded(track2->isArbitratedEcalSeeded() &&
221  ((keepSecond && arbitratedEcalSeededPreferCollection == 3) ||
224  track2->setIsArbitratedMerged(track2->isArbitratedMerged() &&
225  ((keepSecond && arbitratedMergedPreferCollection == 3) ||
228  track2->setIsArbitratedMergedEcalGeneral(track2->isArbitratedMergedEcalGeneral() &&
229  ((keepSecond && arbitratedMergedEcalGeneralPreferCollection == 3) ||
232 
233  } //end got a duplicate
234  } //end track2 loop
235  } //end track loop
236  } //end more than 1 track
237 
238  //
239  // output selected tracks - if any
240  //
241 
242  if (!tC1.empty()) {
243  i = 0;
244  for (reco::ConversionTrackCollection::const_iterator track = tC1.begin(); track != tC1.end(); ++track, ++i) {
245  //don't store tracks rejected as duplicates
246  if (!selected1[i]) {
247  continue;
248  }
249  //fill the TrackCollection
250  outputTrks->push_back(*track);
251  } //end faux loop over tracks
252  } //end more than 0 track
253 
254  //Fill the trajectories, etc. for 1st collection
255 
256  if (!tC2.empty()) {
257  i = 0;
258  for (reco::ConversionTrackCollection::const_iterator track = tC2.begin(); track != tC2.end(); ++track, ++i) {
259  //don't store tracks rejected as duplicates
260  if (!selected2[i]) {
261  continue;
262  }
263  //fill the TrackCollection
264  outputTrks->push_back(*track);
265  } //end faux loop over tracks
266  } //end more than 0 track
267 
268  e.put(std::move(outputTrks));
269  return;
270 
271 } //end produce
PA_ZEESkim_cff.checkCharge
checkCharge
Definition: PA_ZEESkim_cff.py:27
mps_fire.i
i
Definition: mps_fire.py:428
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11776
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
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
ConversionTrackMerger::outputTrks
std::unique_ptr< reco::ConversionTrackCollection > outputTrks
Definition: ConversionTrackMerger.h:43
ConversionTrackMerger::~ConversionTrackMerger
~ConversionTrackMerger() override
Definition: ConversionTrackMerger.cc:44
TrackingRecHit::some
Definition: TrackingRecHit.h:59
HLT_FULL_cff.shareFrac
shareFrac
Definition: HLT_FULL_cff.py:84539
HLT_FULL_cff.allowFirstHitShare
allowFirstHitShare
Definition: HLT_FULL_cff.py:10275
edm::OwnVector::const_iterator
Definition: OwnVector.h:41
ConversionTrackMerger::produce
void produce(edm::Event &e, const edm::EventSetup &c) override
Definition: ConversionTrackMerger.cc:47
TrackerDigiGeometryRecord.h
edm::ParameterSet
Definition: ParameterSet.h:47
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
TrackingRecHit
Definition: TrackingRecHit.h:21
ConversionTrackMerger::ConversionTrackMerger
ConversionTrackMerger(const edm::ParameterSet &conf)
Definition: ConversionTrackMerger.cc:35
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
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
conversionTrackMerger_cfi.trackerOnlyPreferCollection
trackerOnlyPreferCollection
Definition: conversionTrackMerger_cfi.py:18
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