CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TripletFilter.cc
Go to the documentation of this file.
2 
5 
9 
12 
15 
18 
19 using namespace std;
20 
21 /*****************************************************************************/
23 {
24  // Get cluster shape hit filter
26  es.get<CkfComponentsRecord>().get("ClusterShapeHitFilter",shape);
27  theFilter = shape.product();
28 }
29 
30 /*****************************************************************************/
32 {
33 }
34 
35 /*****************************************************************************/
37 (const vector<const TrackingRecHit*>& recHits, const vector<LocalVector>& localDirs, const TrackerTopology *tTopo, const SiPixelClusterShapeCache& clusterShapeCache)
38 {
39  bool ok = true;
40 
41  vector<LocalVector>::const_iterator localDir = localDirs.begin();
42  for(vector<const TrackingRecHit*>::const_iterator recHit = recHits.begin();
43  recHit!= recHits.end();
44  recHit++)
45  {
46  const SiPixelRecHit* pixelRecHit =
47  dynamic_cast<const SiPixelRecHit *>(*recHit);
48 
49  if(! pixelRecHit->isValid())
50  { ok = false; break; }
51 
52  if(! theFilter->isCompatible(*pixelRecHit, *localDir, clusterShapeCache))
53  {
54  LogTrace("MinBiasTracking")
55  << " [TripletFilter] clusShape problem"
56  << HitInfo::getInfo(**recHit,tTopo);
57 
58  ok = false; break;
59  }
60 
61  localDir++;
62  }
63 
64  return ok;
65 }
66 
67 /*****************************************************************************/
69 (const vector<const TrackingRecHit*>& recHits, const vector<GlobalVector>& globalDirs, const TrackerTopology *tTopo, const SiPixelClusterShapeCache& clusterShapeCache)
70 {
71  bool ok = true;
72 
73  vector<GlobalVector>::const_iterator globalDir = globalDirs.begin();
74  for(vector<const TrackingRecHit*>::const_iterator recHit = recHits.begin();
75  recHit!= recHits.end();
76  recHit++)
77  {
78  const SiPixelRecHit* pixelRecHit =
79  dynamic_cast<const SiPixelRecHit *>(*recHit);
80 
81  if(! pixelRecHit->isValid())
82  { ok = false; break; }
83 
84  if(! theFilter->isCompatible(*pixelRecHit, *globalDir, clusterShapeCache))
85  {
86  LogTrace("MinBiasTracking")
87  << " [TripletFilter] clusShape problem"
88  << HitInfo::getInfo(**recHit,tTopo);
89 
90  ok = false; break;
91  }
92 
93  globalDir++;
94  }
95 
96  return ok;
97 }
98 
bool checkTrack(const std::vector< const TrackingRecHit * > &recHits, const std::vector< LocalVector > &localDirs, const TrackerTopology *tTopo, const SiPixelClusterShapeCache &clusterShapeCache)
#define LogTrace(id)
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
bool isValid() const
TripletFilter(const edm::EventSetup &es)
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:23
Our base class.
Definition: SiPixelRecHit.h:23