CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoPixelVertexing/PixelLowPtUtilities/src/TrackCleaner.cc

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