Go to the documentation of this file.00001 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/TrackCleaner.h"
00002 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/HitInfo.h"
00003
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005
00006 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00007 #include "DataFormats/TrackReco/interface/Track.h"
00008
00009 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00010 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00011 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00012
00013 using namespace std;
00014 using namespace pixeltrackfitting;
00015
00016
00017 class HitComparatorByRadius
00018 {
00019 public:
00020 HitComparatorByRadius(const TrackerTopology *tTopo) {tTopo_=tTopo;}
00021 private:
00022 const TrackerTopology *tTopo_;
00023 public:
00024 bool operator() (const TrackingRecHit* a, const TrackingRecHit* b) const
00025 {
00026 DetId i1 = a->geographicalId();
00027 DetId i2 = b->geographicalId();
00028
00029 bool isBarrel = (i1.subdetId() == int(PixelSubdetector::PixelBarrel));
00030
00031 if(i1.subdetId() != i2.subdetId())
00032 {
00033 return isBarrel;
00034 }
00035 else
00036 {
00037 if(isBarrel)
00038 {
00039 int r1 = (tTopo_->pxbLayer(i1) - 1)*2 + (tTopo_->pxbLadder(i1) - 1)%2;
00040 int r2 = (tTopo_->pxbLayer(i2) - 1)*2 + (tTopo_->pxbLadder(i2) - 1)%2;
00041
00042 return (r1 < r2);
00043 }
00044 else
00045 {
00046 int r1 = (tTopo_->pxfDisk(i1) - 1)*2 + (tTopo_->pxfPanel(i1) - 1)%2;
00047 int r2 = (tTopo_->pxfDisk(i2) - 1)*2 + (tTopo_->pxfPanel(i2) - 1)%2;
00048
00049 return (r1 < r2);
00050 }
00051 }
00052 }
00053 };
00054
00055
00056 class HitComparator
00057 {
00058 public:
00059 bool operator() (const TrackingRecHit* a, const TrackingRecHit* b) const
00060 {
00061 if(a->geographicalId() < b->geographicalId()) return true;
00062 if(b->geographicalId() < a->geographicalId()) return false;
00063
00064 if(a->localPosition().x() < b->localPosition().x() - 1e-5) return true;
00065 if(b->localPosition().x() < a->localPosition().x() - 1e-5) return false;
00066
00067 if(a->localPosition().y() < b->localPosition().y() - 1e-5) return true;
00068 if(b->localPosition().y() < a->localPosition().y() - 1e-5) return false;
00069
00070 return false;
00071 }
00072 };
00073
00074
00075 TrackCleaner::TrackCleaner
00076 (const edm::ParameterSet& ps)
00077 {
00078 }
00079
00080
00081 TrackCleaner::~TrackCleaner()
00082 {
00083 }
00084
00085
00086 bool TrackCleaner::areSame(const TrackingRecHit * a,
00087 const TrackingRecHit * b)
00088 {
00089 if(a->geographicalId() != b->geographicalId())
00090 return false;
00091
00092 if(fabs(a->localPosition().x() - b->localPosition().x()) < 1e-5 &&
00093 fabs(a->localPosition().y() - b->localPosition().y()) < 1e-5)
00094 return true;
00095 else
00096 return false;
00097 }
00098
00099
00100 bool TrackCleaner::isCompatible(const DetId & i1,
00101 const DetId & i2,
00102 const TrackerTopology *tTopo)
00103 {
00104
00105 if(i1.subdetId() != i2.subdetId()) return true;
00106
00107 if(i1.subdetId() == int(PixelSubdetector::PixelBarrel))
00108 {
00109
00110 if(tTopo->pxbLayer(i1) != tTopo->pxbLayer(i2)) return true;
00111
00112 int dphi = abs(int(tTopo->pxbLadder(i1) - tTopo->pxbLadder(i2)));
00113 static int max[3] = {20, 32, 44};
00114 if(dphi > max[tTopo->pxbLayer(i1)-1] / 2) dphi = max[tTopo->pxbLayer(i1)-1] - dphi;
00115
00116 int dz = abs(int(tTopo->pxbModule(i1) - tTopo->pxbModule(i2)));
00117
00118 if(dphi == 1 && dz <= 1) return true;
00119 }
00120 else
00121 {
00122
00123 if(tTopo->pxfSide(i1) != tTopo->pxfSide(i2) ||
00124 tTopo->pxfDisk(i1) != tTopo->pxfDisk(i2)) return true;
00125
00126 int dphi = abs(int(tTopo->pxfBlade(i1) - tTopo->pxfBlade(i2)));
00127 static int max = 24;
00128 if(dphi > max / 2) dphi = max - dphi;
00129
00130 int dr = abs(int( ((tTopo->pxfModule(i1)-1) * 2 + (tTopo->pxfPanel(i1)-1)) -
00131 ((tTopo->pxfModule(i2)-1) * 2 + (tTopo->pxfPanel(i2)-1)) ));
00132
00133 if(dphi <= 1 && dr <= 1 && !(dphi == 0 && dr == 0)) return true;
00134 }
00135
00136 return false;
00137 }
00138
00139
00140 bool TrackCleaner::canBeMerged
00141 (vector<const TrackingRecHit *> recHitsA,
00142 vector<const TrackingRecHit *> recHitsB,
00143 const TrackerTopology *tTopo)
00144 {
00145 bool ok = true;
00146
00147 for(vector<const TrackingRecHit *>::const_iterator
00148 recHitA = recHitsA.begin(); recHitA!= recHitsA.end(); recHitA++)
00149 for(vector<const TrackingRecHit *>::const_iterator
00150 recHitB = recHitsB.begin(); recHitB!= recHitsB.end(); recHitB++)
00151 if(!areSame(*recHitA,*recHitB))
00152 if(!isCompatible((*recHitA)->geographicalId(),
00153 (*recHitB)->geographicalId(),
00154 tTopo))
00155 ok = false;
00156
00157 return ok;
00158 }
00159
00160
00161 TracksWithRecHits TrackCleaner::cleanTracks
00162 (const TracksWithRecHits & tracks_, const TrackerTopology *tTopo)
00163 {
00164
00165 TracksWithRecHits tracks = tracks_;
00166
00167 typedef map<const TrackingRecHit*,vector<unsigned int>,HitComparator>
00168 RecHitMap;
00169
00170 vector<bool> keep(tracks.size(),true);
00171
00172 int changes;
00173
00174 LogTrace("MinBiasTracking")
00175 << " [TrackCleaner] initial tracks : " << tracks.size();
00176
00177 for(unsigned int i = 0; i < tracks.size(); i++)
00178 LogTrace("TrackCleaner")
00179 << " Track #" << i << " : " << HitInfo::getInfo(tracks[i].second,tTopo);
00180
00181 do
00182 {
00183 changes = 0;
00184
00185 RecHitMap recHitMap;
00186
00187
00188
00189
00190
00191
00192
00193 for(unsigned int i = 0; i < tracks.size(); i++)
00194 if(keep[i])
00195 {
00196 for(vector<const TrackingRecHit *>::const_iterator
00197 recHit = tracks[i].second.begin();
00198 recHit!= tracks[i].second.end(); recHit++)
00199 recHitMap[*recHit].push_back(i);
00200 }
00201
00202
00203 typedef map<unsigned int,int,less<unsigned int> > TrackMap;
00204
00205 for(unsigned int i = 0; i < tracks.size(); i++)
00206 {
00207
00208 if(!keep[i]) continue;
00209
00210
00211
00212
00213
00214
00215
00216 TrackMap trackMap;
00217
00218
00219 for(vector<const TrackingRecHit *>::const_iterator
00220 recHit = tracks[i].second.begin();
00221 recHit!= tracks[i].second.end(); recHit++)
00222 {
00223
00224 vector<unsigned int> sharing = recHitMap[*recHit];
00225
00226 for(vector<unsigned int>::iterator j = sharing.begin();
00227 j!= sharing.end(); j++)
00228 if(i < *j)
00229 trackMap[*j]++;
00230 }
00231
00232
00233 for(TrackMap::iterator sharing = trackMap.begin();
00234 sharing!= trackMap.end(); sharing++)
00235 {
00236 unsigned int j = (*sharing).first;
00237 if(!keep[i] || !keep[j]) continue;
00238
00239 if(tracks[i].second.size() >=3)
00240 {
00241 if((*sharing).second > min(int(tracks[i].second.size()),
00242 int(tracks[j].second.size()))/2)
00243 {
00244 if(canBeMerged(tracks[i].second,tracks[j].second,tTopo))
00245 {
00246
00247 for(vector<const TrackingRecHit *>::const_iterator
00248 recHit = tracks[j].second.begin();
00249 recHit!= tracks[j].second.end(); recHit++)
00250 {
00251 bool ok = true;
00252 for(vector<const TrackingRecHit *>::const_iterator
00253 recHitA = tracks[i].second.begin();
00254 recHitA!= tracks[i].second.end(); recHitA++)
00255 if(areSame(*recHit,*recHitA)) ok = false;
00256
00257 if(ok)
00258 {
00259 tracks[i].second.push_back(*recHit);
00260 recHitMap[*recHit].push_back(i);
00261
00262 sort(tracks[i].second.begin(),
00263 tracks[i].second.end(),
00264 HitComparatorByRadius(tTopo));
00265
00266
00267 }
00268 }
00269
00270 LogTrace("TrackCleaner")
00271 << " Merge #" << i << " #" << j
00272 << ", first now has " << tracks[i].second.size();
00273
00274
00275 keep[j] = false;
00276
00277 changes++;
00278 }
00279 else
00280 {
00281 if(fabs(tracks[i].first->d0())
00282 < fabs(tracks[j].first->d0()))
00283 keep[j] = false;
00284 else
00285 keep[i] = false;
00286
00287 LogTrace("TrackCleaner")
00288 << " Clash #" << i << " #" << j
00289 << " keep lower d0 " << tracks[i].first->d0()
00290 << " " << tracks[j].first->d0()
00291 << ", keep #" << (keep[i] ? i : ( keep[j] ? j : 9999 ) );
00292
00293 changes++;
00294 }
00295 }
00296 else
00297 {
00298 if((*sharing).second > 1)
00299 {
00300 if(tracks[i].second.size() != tracks[j].second.size())
00301 {
00302 if(tracks[i].second.size() > tracks[j].second.size())
00303 keep[j] = false; else keep[i] = false;
00304 changes++;
00305
00306 LogTrace("TrackCleaner")
00307 << " Sharing " << (*sharing).second << " remove by size";
00308 }
00309 else
00310 {
00311 if(fabs(tracks[i].first->d0())
00312 < fabs(tracks[j].first->d0()))
00313 keep[j] = false; else keep[i] = false;
00314 changes++;
00315
00316 LogTrace("TrackCleaner")
00317 << " Sharing " << (*sharing).second << " remove by d0";
00318 }
00319 }
00320 }
00321 }
00322 else
00323 {
00324 if((*sharing).second > 0)
00325 {
00326
00327 keep[j] = false;
00328
00329 changes++;
00330 }
00331 }
00332 }
00333
00334
00335
00336
00337 }
00338 }
00339 while(changes > 0);
00340
00341
00342 TracksWithRecHits cleaned;
00343
00344 for(unsigned int i = 0; i < tracks.size(); i++)
00345 if(keep[i])
00346 cleaned.push_back(tracks[i]);
00347 else delete tracks_[i].first;
00348
00349 LogTrace("MinBiasTracking")
00350 << " [TrackCleaner] cleaned tracks : " << cleaned.size();
00351
00352 for(unsigned int i = 0; i < cleaned.size(); i++)
00353 LogTrace("TrackCleaner")
00354 << " Track #" << i << " : " << HitInfo::getInfo(cleaned[i].second,tTopo);
00355
00356 return cleaned;
00357 }
00358