CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
TrackCleaner Class Reference

#include <TrackCleaner.h>

Inheritance diagram for TrackCleaner:
PixelTrackCleaner

Public Member Functions

TracksWithRecHits cleanTracks (const TracksWithRecHits &tracksWithRecHits) const override
 
 TrackCleaner (const TrackerTopology *tTopo)
 
 ~TrackCleaner () override
 
- Public Member Functions inherited from PixelTrackCleaner
virtual void cleanTracks (TracksWithTTRHs &tracksWithRecHits) const
 
bool fast () const
 
virtual ~PixelTrackCleaner ()
 

Private Member Functions

bool areSame (const TrackingRecHit *a, const TrackingRecHit *b) const
 
bool canBeMerged (const std::vector< const TrackingRecHit * > &recHitsA, const std::vector< const TrackingRecHit * > &recHitsB) const
 
bool isCompatible (const DetId &i1, const DetId &i2) const
 
std::vector< const TrackingRecHit * > ttrhs (const SeedingHitSet &h) const
 

Private Attributes

const TrackerTopologytTopo_
 

Additional Inherited Members

- Public Types inherited from PixelTrackCleaner
using Record = TrajectoryFilter::Record
 
typedef pixeltrackfitting::TracksWithRecHits TracksWithRecHits
 
using TracksWithTTRHs = pixeltrackfitting::TracksWithTTRHs
 
- Protected Member Functions inherited from PixelTrackCleaner
 PixelTrackCleaner (bool fast=false)
 

Detailed Description

Definition at line 15 of file TrackCleaner.h.

Constructor & Destructor Documentation

TrackCleaner::TrackCleaner ( const TrackerTopology tTopo)
explicit

Definition at line 75 of file TrackCleaner.cc.

75  : tTopo_(tTopo)
76 {
77 }
const TrackerTopology * tTopo_
Definition: TrackCleaner.h:34
TrackCleaner::~TrackCleaner ( )
override

Definition at line 80 of file TrackCleaner.cc.

81 {
82 }

Member Function Documentation

bool TrackCleaner::areSame ( const TrackingRecHit a,
const TrackingRecHit b 
) const
private

Definition at line 85 of file TrackCleaner.cc.

References MillePedeFileConverter_cfg::e, TrackingRecHit::geographicalId(), TrackingRecHit::localPosition(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by canBeMerged(), and cleanTracks().

87 {
88  if(a->geographicalId() != b->geographicalId())
89  return false;
90 
91  if(fabs(a->localPosition().x() - b->localPosition().x()) < 1e-5 &&
92  fabs(a->localPosition().y() - b->localPosition().y()) < 1e-5)
93  return true;
94  else
95  return false;
96 }
T y() const
Definition: PV3DBase.h:63
virtual LocalPoint localPosition() const =0
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:62
bool TrackCleaner::canBeMerged ( const std::vector< const TrackingRecHit * > &  recHitsA,
const std::vector< const TrackingRecHit * > &  recHitsB 
) const
private

Definition at line 139 of file TrackCleaner.cc.

References areSame(), cleanTracks(), isCompatible(), and convertSQLiteXML::ok.

Referenced by cleanTracks(), and isCompatible().

141 {
142  bool ok = true;
143 
144  for(vector<const TrackingRecHit *>::const_iterator
145  recHitA = recHitsA.begin(); recHitA!= recHitsA.end(); recHitA++)
146  for(vector<const TrackingRecHit *>::const_iterator
147  recHitB = recHitsB.begin(); recHitB!= recHitsB.end(); recHitB++)
148  if(!areSame(*recHitA,*recHitB))
149  if(!isCompatible((*recHitA)->geographicalId(),
150  (*recHitB)->geographicalId()))
151  ok = false;
152 
153  return ok;
154 }
bool isCompatible(const DetId &i1, const DetId &i2) const
Definition: TrackCleaner.cc:99
bool areSame(const TrackingRecHit *a, const TrackingRecHit *b) const
Definition: TrackCleaner.cc:85
TracksWithRecHits TrackCleaner::cleanTracks ( const TracksWithRecHits tracksWithRecHits) const
overridevirtual

Reimplemented from PixelTrackCleaner.

Definition at line 158 of file TrackCleaner.cc.

References areSame(), canBeMerged(), plotBeamSpotDB::first, HitInfo::getInfo(), mps_fire::i, keep, LogTrace, min(), convertSQLiteXML::ok, rpcPointValidation_cfi::recHit, edm::second(), l1t::tracks, tracks_, and tTopo_.

Referenced by canBeMerged().

159 {
160  // Local copy
162 
163  typedef map<const TrackingRecHit*,vector<unsigned int>,HitComparator>
164  RecHitMap;
165 
166  vector<bool> keep(tracks.size(),true);
167 
168  int changes;
169 
170  LogTrace("MinBiasTracking")
171  << " [TrackCleaner] initial tracks : " << tracks.size();
172 
173  for(unsigned int i = 0; i < tracks.size(); i++)
174  LogTrace("TrackCleaner")
175  << " Track #" << i << " : " << HitInfo::getInfo(tracks[i].second, tTopo_);
176 
177  do
178  {
179  changes = 0;
180 
181  RecHitMap recHitMap;
182 
183 /*
184  LogTrace("MinBiasTracking")
185  << " [TrackCleaner] fill rechit map";
186 */
187 
188  // Fill the rechit map
189  for(unsigned int i = 0; i < tracks.size(); i++)
190  if(keep[i])
191  {
192  for(vector<const TrackingRecHit *>::const_iterator
193  recHit = tracks[i].second.begin();
194  recHit!= tracks[i].second.end(); recHit++)
195  recHitMap[*recHit].push_back(i);
196  }
197 
198  // Look at each track
199  typedef map<unsigned int,int,less<unsigned int> > TrackMap;
200 
201  for(unsigned int i = 0; i < tracks.size(); i++)
202  {
203  // Skip if 'i' already removed
204  if(!keep[i]) continue;
205 
206 /*
207 
208  bool addedNewHit = false;
209  do
210  {
211 */
212  TrackMap trackMap;
213 
214  // Go trough all rechits of this track
215  for(vector<const TrackingRecHit *>::const_iterator
216  recHit = tracks[i].second.begin();
217  recHit!= tracks[i].second.end(); recHit++)
218  {
219  // Get tracks sharing this rechit
220  vector<unsigned int> sharing = recHitMap[*recHit];
221 
222  for(vector<unsigned int>::iterator j = sharing.begin();
223  j!= sharing.end(); j++)
224  if(i < *j)
225  trackMap[*j]++;
226  }
227 
228  // Check for tracks with shared rechits
229  for(TrackMap::iterator sharing = trackMap.begin();
230  sharing!= trackMap.end(); sharing++)
231  {
232  unsigned int j = (*sharing).first;
233  if(!keep[i] || !keep[j]) continue;
234 
235  if(tracks[i].second.size() >=3)
236  { // triplet tracks
237  if((*sharing).second > min(int(tracks[i].second.size()),
238  int(tracks[j].second.size()))/2)
239  { // more than min(hits1,hits2)/2 rechits are shared
240  if(canBeMerged(tracks[i].second,tracks[j].second))
241  { // no common layer
242  // merge tracks, add separate hits of the second to the first one
243  for(vector<const TrackingRecHit *>::const_iterator
244  recHit = tracks[j].second.begin();
245  recHit!= tracks[j].second.end(); recHit++)
246  {
247  bool ok = true;
248  for(vector<const TrackingRecHit *>::const_iterator
249  recHitA = tracks[i].second.begin();
250  recHitA!= tracks[i].second.end(); recHitA++)
251  if(areSame(*recHit,*recHitA)) ok = false;
252 
253  if(ok)
254  {
255  tracks[i].second.push_back(*recHit);
256  recHitMap[*recHit].push_back(i);
257 
258  sort(tracks[i].second.begin(),
259  tracks[i].second.end(),
261 
262  //addedNewHit = true;
263  }
264  }
265 
266  LogTrace("TrackCleaner")
267  << " Merge #" << i << " #" << j
268  << ", first now has " << tracks[i].second.size();
269 
270  // Remove second track
271  keep[j] = false;
272 
273  changes++;
274  }
275  else
276  { // there is a common layer, keep smaller impact
277  if(fabs(tracks[i].first->d0())
278  < fabs(tracks[j].first->d0()))
279  keep[j] = false;
280  else
281  keep[i] = false;
282 
283  LogTrace("TrackCleaner")
284  << " Clash #" << i << " #" << j
285  << " keep lower d0 " << tracks[i].first->d0()
286  << " " << tracks[j].first->d0()
287  << ", keep #" << (keep[i] ? i : ( keep[j] ? j : 9999 ) );
288 
289  changes++;
290  }
291  }
292  else
293  { // note more than 50%, but at least two are shared
294  if((*sharing).second > 1)
295  {
296  if(tracks[i].second.size() != tracks[j].second.size())
297  { // keep longer
298  if(tracks[i].second.size() > tracks[j].second.size())
299  keep[j] = false; else keep[i] = false;
300  changes++;
301 
302  LogTrace("TrackCleaner")
303  << " Sharing " << (*sharing).second << " remove by size";
304  }
305  else
306  { // keep smaller impact
307  if(fabs(tracks[i].first->d0())
308  < fabs(tracks[j].first->d0()))
309  keep[j] = false; else keep[i] = false;
310  changes++;
311 
312  LogTrace("TrackCleaner")
313  << " Sharing " << (*sharing).second << " remove by d0";
314  }
315  }
316  }
317  }
318  else
319  { // pair tracks
320  if((*sharing).second > 0)
321  {
322  // Remove second track
323  keep[j] = false;
324 
325  changes++;
326  }
327  }
328  }
329 /*
330  }
331  while(addedNewHit);
332 */
333  }
334  }
335  while(changes > 0);
336 
337  // Final copy
338  TracksWithRecHits cleaned;
339 
340  for(unsigned int i = 0; i < tracks.size(); i++)
341  if(keep[i])
342  cleaned.push_back(tracks[i]);
343  else delete tracks_[i].first;
344 
345  LogTrace("MinBiasTracking")
346  << " [TrackCleaner] cleaned tracks : " << cleaned.size();
347 
348  for(unsigned int i = 0; i < cleaned.size(); i++)
349  LogTrace("TrackCleaner")
350  << " Track #" << i << " : " << HitInfo::getInfo(cleaned[i].second,tTopo_);
351 
352  return cleaned;
353 }
bool canBeMerged(const std::vector< const TrackingRecHit * > &recHitsA, const std::vector< const TrackingRecHit * > &recHitsB) const
const std::vector< reco::PFCandidatePtr > & tracks_
U second(std::pair< T, U > const &p)
const int keep
T min(T a, T b)
Definition: MathUtil.h:58
#define LogTrace(id)
const TrackerTopology * tTopo_
Definition: TrackCleaner.h:34
std::vector< TrackWithRecHits > TracksWithRecHits
bool areSame(const TrackingRecHit *a, const TrackingRecHit *b) const
Definition: TrackCleaner.cc:85
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:23
bool TrackCleaner::isCompatible ( const DetId i1,
const DetId i2 
) const
private

Definition at line 99 of file TrackCleaner.cc.

References funct::abs(), canBeMerged(), constexpr, runTauDisplay::dr, PVValHelper::dz, createfilelist::int, SiStripPI::max, PixelSubdetector::PixelBarrel, TrackerTopology::pxbLadder(), TrackerTopology::pxbLayer(), TrackerTopology::pxbModule(), TrackerTopology::pxfBlade(), TrackerTopology::pxfDisk(), TrackerTopology::pxfModule(), TrackerTopology::pxfPanel(), TrackerTopology::pxfSide(), DetId::subdetId(), and tTopo_.

Referenced by canBeMerged().

101 {
102  // different subdet
103  if(i1.subdetId() != i2.subdetId()) return true;
104 
106  { // barrel
107 
108  if(tTopo_->pxbLayer(i1) != tTopo_->pxbLayer(i2)) return true;
109 
110  int dphi = abs(int(tTopo_->pxbLadder(i1) - tTopo_->pxbLadder(i2)));
111  constexpr int max[3] = {20, 32, 44};
112  if(dphi > max[tTopo_->pxbLayer(i1)-1] / 2) dphi = max[tTopo_->pxbLayer(i1)-1] - dphi;
113 
114  int dz = abs(int(tTopo_->pxbModule(i1) - tTopo_->pxbModule(i2)));
115 
116  if(dphi == 1 && dz <= 1) return true;
117  }
118  else
119  { // endcap
120 
121  if(tTopo_->pxfSide(i1) != tTopo_->pxfSide(i2) ||
122  tTopo_->pxfDisk(i1) != tTopo_->pxfDisk(i2)) return true;
123 
124  int dphi = abs(int(tTopo_->pxfBlade(i1) - tTopo_->pxfBlade(i2)));
125  constexpr int max = 24;
126  if(dphi > max / 2) dphi = max - dphi;
127 
128  int dr = abs(int( ((tTopo_->pxfModule(i1)-1) * 2 + (tTopo_->pxfPanel(i1)-1)) -
129  ((tTopo_->pxfModule(i2)-1) * 2 + (tTopo_->pxfPanel(i2)-1)) ));
130 
131  if(dphi <= 1 && dr <= 1 && !(dphi == 0 && dr == 0)) return true;
132  }
133 
134  return false;
135 }
unsigned int pxfDisk(const DetId &id) const
unsigned int pxbLadder(const DetId &id) const
unsigned int pxbModule(const DetId &id) const
#define constexpr
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int pxfModule(const DetId &id) const
unsigned int pxbLayer(const DetId &id) const
const TrackerTopology * tTopo_
Definition: TrackCleaner.h:34
unsigned int pxfSide(const DetId &id) const
unsigned int pxfPanel(const DetId &id) const
unsigned int pxfBlade(const DetId &id) const
std::vector<const TrackingRecHit*> TrackCleaner::ttrhs ( const SeedingHitSet h) const
private

Member Data Documentation

const TrackerTopology* TrackCleaner::tTopo_
private

Definition at line 34 of file TrackCleaner.h.

Referenced by cleanTracks(), and isCompatible().