CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackCleaner.cc
Go to the documentation of this file.
3 
5 
8 
12 
13 using namespace std;
14 using namespace pixeltrackfitting;
15 
16 /*****************************************************************************/
18 { // No access to geometry!
19  public:
20  bool operator() (const TrackingRecHit* a, const TrackingRecHit* b) const
21  {
22  DetId i1 = a->geographicalId();
23  DetId i2 = b->geographicalId();
24 
25  bool isBarrel = (i1.subdetId() == int(PixelSubdetector::PixelBarrel));
26 
27  if(i1.subdetId() != i2.subdetId())
28  {
29  return isBarrel;
30  }
31  else
32  {
33  if(isBarrel)
34  {
35  PXBDetId p1(i1);
36  PXBDetId p2(i2);
37 
38  int r1 = (p1.layer() - 1)*2 + (p1.ladder() - 1)%2;
39  int r2 = (p2.layer() - 1)*2 + (p2.ladder() - 1)%2;
40 
41  return (r1 < r2);
42  }
43  else
44  {
45  PXFDetId p1(i1);
46  PXFDetId p2(i2);
47 
48  int r1 = (p1.disk() - 1)*2 + (p1.panel() - 1)%2;
49  int r2 = (p2.disk() - 1)*2 + (p2.panel() - 1)%2;
50 
51  return (r1 < r2);
52  }
53  }
54  }
55 };
56 
57 /*****************************************************************************/
59 {
60 public:
61  bool operator() (const TrackingRecHit* a, const TrackingRecHit* b) const
62  {
63  if(a->geographicalId() < b->geographicalId()) return true;
64  if(b->geographicalId() < a->geographicalId()) return false;
65 
66  if(a->localPosition().x() < b->localPosition().x() - 1e-5) return true;
67  if(b->localPosition().x() < a->localPosition().x() - 1e-5) return false;
68 
69  if(a->localPosition().y() < b->localPosition().y() - 1e-5) return true;
70  if(b->localPosition().y() < a->localPosition().y() - 1e-5) return false;
71 
72  return false;
73  }
74 };
75 
76 /*****************************************************************************/
78  (const edm::ParameterSet& ps)
79 {
80 }
81 
82 /*****************************************************************************/
84 {
85 }
86 
87 /*****************************************************************************/
89  const TrackingRecHit * b)
90 {
91  if(a->geographicalId() != b->geographicalId())
92  return false;
93 
94  if(fabs(a->localPosition().x() - b->localPosition().x()) < 1e-5 &&
95  fabs(a->localPosition().y() - b->localPosition().y()) < 1e-5)
96  return true;
97  else
98  return false;
99 }
100 
101 /*****************************************************************************/
103  const DetId & i2)
104 {
105  // different subdet
106  if(i1.subdetId() != i2.subdetId()) return true;
107 
108  if(i1.subdetId() == int(PixelSubdetector::PixelBarrel))
109  { // barrel
110  PXBDetId p1(i1);
111  PXBDetId p2(i2);
112 
113  if(p1.layer() != p2.layer()) return true;
114 
115  int dphi = abs(int(p1.ladder() - p2.ladder()));
116  static int max[3] = {20, 32, 44};
117  if(dphi > max[p1.layer()-1] / 2) dphi = max[p1.layer()-1] - dphi;
118 
119  int dz = abs(int(p1.module() - p2.module()));
120 
121  if(dphi == 1 && dz <= 1) return true;
122  }
123  else
124  { // endcap
125  PXFDetId p1(i1);
126  PXFDetId p2(i2);
127 
128  if(p1.side() != p2.side() ||
129  p1.disk() != p2.disk()) return true;
130 
131  int dphi = abs(int(p1.blade() - p2.blade()));
132  static int max = 24;
133  if(dphi > max / 2) dphi = max - dphi;
134 
135  int dr = abs(int( ((p1.module()-1) * 2 + (p1.panel()-1)) -
136  ((p2.module()-1) * 2 + (p2.panel()-1)) ));
137 
138  if(dphi <= 1 && dr <= 1 && !(dphi == 0 && dr == 0)) return true;
139  }
140 
141  return false;
142 }
143 
144 /*****************************************************************************/
146  (vector<const TrackingRecHit *> recHitsA,
147  vector<const TrackingRecHit *> recHitsB)
148 {
149  bool ok = true;
150 
151  for(vector<const TrackingRecHit *>::const_iterator
152  recHitA = recHitsA.begin(); recHitA!= recHitsA.end(); recHitA++)
153  for(vector<const TrackingRecHit *>::const_iterator
154  recHitB = recHitsB.begin(); recHitB!= recHitsB.end(); recHitB++)
155  if(!areSame(*recHitA,*recHitB))
156  if(!isCompatible((*recHitA)->geographicalId(),
157  (*recHitB)->geographicalId()))
158  ok = false;
159 
160  return ok;
161 }
162 
163 /*****************************************************************************/
166 {
167  // Local copy
169 
170  typedef map<const TrackingRecHit*,vector<unsigned int>,HitComparator>
171  RecHitMap;
172 
173  vector<bool> keep(tracks.size(),true);
174 
175  int changes;
176 
177  LogTrace("MinBiasTracking")
178  << " [TrackCleaner] initial tracks : " << tracks.size();
179 
180  for(unsigned int i = 0; i < tracks.size(); i++)
181  LogTrace("TrackCleaner")
182  << " Track #" << i << " : " << HitInfo::getInfo(tracks[i].second);
183 
184  do
185  {
186  changes = 0;
187 
188  RecHitMap recHitMap;
189 
190 /*
191  LogTrace("MinBiasTracking")
192  << " [TrackCleaner] fill rechit map";
193 */
194 
195  // Fill the rechit map
196  for(unsigned int i = 0; i < tracks.size(); i++)
197  if(keep[i])
198  {
199  for(vector<const TrackingRecHit *>::const_iterator
200  recHit = tracks[i].second.begin();
201  recHit!= tracks[i].second.end(); recHit++)
202  recHitMap[*recHit].push_back(i);
203  }
204 
205  // Look at each track
206  typedef map<unsigned int,int,less<unsigned int> > TrackMap;
207 
208  for(unsigned int i = 0; i < tracks.size(); i++)
209  {
210  // Skip if 'i' already removed
211  if(!keep[i]) continue;
212 
213 /*
214 
215  bool addedNewHit = false;
216  do
217  {
218 */
219  TrackMap trackMap;
220 
221  // Go trough all rechits of this track
222  for(vector<const TrackingRecHit *>::const_iterator
223  recHit = tracks[i].second.begin();
224  recHit!= tracks[i].second.end(); recHit++)
225  {
226  // Get tracks sharing this rechit
227  vector<unsigned int> sharing = recHitMap[*recHit];
228 
229  for(vector<unsigned int>::iterator j = sharing.begin();
230  j!= sharing.end(); j++)
231  if(i < *j)
232  trackMap[*j]++;
233  }
234 
235  // Check for tracks with shared rechits
236  for(TrackMap::iterator sharing = trackMap.begin();
237  sharing!= trackMap.end(); sharing++)
238  {
239  unsigned int j = (*sharing).first;
240  if(!keep[i] || !keep[j]) continue;
241 
242  if(tracks[i].second.size() >=3)
243  { // triplet tracks
244  if((*sharing).second > min(int(tracks[i].second.size()),
245  int(tracks[j].second.size()))/2)
246  { // more than min(hits1,hits2)/2 rechits are shared
247  if(canBeMerged(tracks[i].second,tracks[j].second))
248  { // no common layer
249  // merge tracks, add separate hits of the second to the first one
250  for(vector<const TrackingRecHit *>::const_iterator
251  recHit = tracks[j].second.begin();
252  recHit!= tracks[j].second.end(); recHit++)
253  {
254  bool ok = true;
255  for(vector<const TrackingRecHit *>::const_iterator
256  recHitA = tracks[i].second.begin();
257  recHitA!= tracks[i].second.end(); recHitA++)
258  if(areSame(*recHit,*recHitA)) ok = false;
259 
260  if(ok)
261  {
262  tracks[i].second.push_back(*recHit);
263  recHitMap[*recHit].push_back(i);
264 
265  sort(tracks[i].second.begin(),
266  tracks[i].second.end(),
268 
269  //addedNewHit = true;
270  }
271  }
272 
273  LogTrace("TrackCleaner")
274  << " Merge #" << i << " #" << j
275  << ", first now has " << tracks[i].second.size();
276 
277  // Remove second track
278  keep[j] = false;
279 
280  changes++;
281  }
282  else
283  { // there is a common layer, keep smaller impact
284  if(fabs(tracks[i].first->d0())
285  < fabs(tracks[j].first->d0()))
286  keep[j] = false;
287  else
288  keep[i] = false;
289 
290  LogTrace("TrackCleaner")
291  << " Clash #" << i << " #" << j
292  << " keep lower d0 " << tracks[i].first->d0()
293  << " " << tracks[j].first->d0()
294  << ", keep #" << (keep[i] ? i : ( keep[j] ? j : 9999 ) );
295 
296  changes++;
297  }
298  }
299  else
300  { // note more than 50%, but at least two are shared
301  if((*sharing).second > 1)
302  {
303  if(tracks[i].second.size() != tracks[j].second.size())
304  { // keep longer
305  if(tracks[i].second.size() > tracks[j].second.size())
306  keep[j] = false; else keep[i] = false;
307  changes++;
308 
309  LogTrace("TrackCleaner")
310  << " Sharing " << (*sharing).second << " remove by size";
311  }
312  else
313  { // keep smaller impact
314  if(fabs(tracks[i].first->d0())
315  < fabs(tracks[j].first->d0()))
316  keep[j] = false; else keep[i] = false;
317  changes++;
318 
319  LogTrace("TrackCleaner")
320  << " Sharing " << (*sharing).second << " remove by d0";
321  }
322  }
323  }
324  }
325  else
326  { // pair tracks
327  if((*sharing).second > 0)
328  {
329  // Remove second track
330  keep[j] = false;
331 
332  changes++;
333  }
334  }
335  }
336 /*
337  }
338  while(addedNewHit);
339 */
340  }
341  }
342  while(changes > 0);
343 
344  // Final copy
345  TracksWithRecHits cleaned;
346 
347  for(unsigned int i = 0; i < tracks.size(); i++)
348  if(keep[i])
349  cleaned.push_back(tracks[i]);
350  else delete tracks_[i].first;
351 
352  LogTrace("MinBiasTracking")
353  << " [TrackCleaner] cleaned tracks : " << cleaned.size();
354 
355  for(unsigned int i = 0; i < cleaned.size(); i++)
356  LogTrace("TrackCleaner")
357  << " Track #" << i << " : " << HitInfo::getInfo(cleaned[i].second);
358 
359  return cleaned;
360 }
361 
static std::string getInfo(const DetId &id)
Definition: HitInfo.cc:24
TrackCleaner(const edm::ParameterSet &ps)
Definition: TrackCleaner.cc:78
int i
Definition: DBlmapReader.cc:9
unsigned int panel() const
panel id
Definition: PXFDetId.h:52
bool areSame(const TrackingRecHit *a, const TrackingRecHit *b)
Definition: TrackCleaner.cc:88
T y() const
Definition: PV3DBase.h:62
#define abs(x)
Definition: mlp_lapack.h:159
unsigned int ladder() const
ladder id
Definition: PXBDetId.h:39
#define min(a, b)
Definition: mlp_lapack.h:161
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
unsigned int blade() const
blade id
Definition: PXFDetId.h:48
U second(std::pair< T, U > const &p)
const int keep
bool isCompatible(const DetId &i1, const DetId &i2)
const T & max(const T &a, const T &b)
virtual ~TrackCleaner()
Definition: TrackCleaner.cc:83
bool canBeMerged(std::vector< const TrackingRecHit * > recHitsA, std::vector< const TrackingRecHit * > recHitsB)
unsigned int module() const
det id
Definition: PXBDetId.h:43
int j
Definition: DBlmapReader.cc:9
unsigned int module() const
det id
Definition: PXFDetId.h:56
bool first
Definition: L1TdeRCT.cc:94
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
unsigned int disk() const
disk id
Definition: PXFDetId.h:43
Definition: DetId.h:20
virtual TracksWithRecHits cleanTracks(const TracksWithRecHits &tracksWithRecHits)
tuple tracks
Definition: testEve_cfg.py:39
double b
Definition: hdecay.h:120
double p1[4]
Definition: TauolaWrapper.h:89
double a
Definition: hdecay.h:121
std::vector< TrackWithRecHits > TracksWithRecHits
unsigned int side() const
positive or negative id
Definition: PXFDetId.h:38
const reco::PFCandidateRefVector & tracks_
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:61
virtual LocalPoint localPosition() const =0
pixeltrackfitting::TracksWithRecHits TracksWithRecHits