CMS 3D CMS Logo

ClusterShapeTrackFilter.cc
Go to the documentation of this file.
17 
18 inline float sqr(float x) { return x*x; }
19 
20 using namespace std;
21 
22 /*****************************************************************************/
24  theClusterShapeCache(cache),
25  ptMin(ptmin), ptMax(ptmax)
26 {
27  // Get tracker geometry
29  es.get<TrackerDigiGeometryRecord>().get(tracker);
30  theTracker = tracker.product();
31 
32  // Get cluster shape hit filter
34  es.get<CkfComponentsRecord>().get("ClusterShapeHitFilter",shape);
35  theFilter = shape.product();
36 
38  es.get<TrackerTopologyRcd>().get(tTopoHand);
39  tTopo = tTopoHand.product();
40 }
41 
42 /*****************************************************************************/
44 {
45 }
46 
47 /*****************************************************************************/
49  (const Global2DVector& a, const Global2DVector& b) const
50 {
51  return a.x() * b.y() - a.y() * b.x();
52 }
53 
54 /*****************************************************************************/
55 vector<GlobalVector> ClusterShapeTrackFilter::getGlobalDirs
56  (const vector<GlobalPoint> & g) const
57 {
58  // Get 2d points
59  vector<Global2DVector> p;
60  for(vector<GlobalPoint>::const_iterator ig = g.begin();
61  ig!= g.end(); ig++)
62  p.push_back( Global2DVector(ig->x(), ig->y()) );
63 
64  //
65  vector<GlobalVector> globalDirs;
66 
67  // Determine circle
68  CircleFromThreePoints circle(g[0],g[1],g[2]);
69 
70  if(circle.curvature() != 0.)
71  {
72  Global2DVector c (circle.center().x(), circle.center().y());
73 
74  float rad2 = (p[0] - c).mag2();
75  float a12 = asin(fabsf(areaParallelogram(p[0] - c, p[1] - c)) / rad2);
76 
77  float slope = (g[1].z() - g[0].z()) / a12;
78 
79  float cotTheta = slope * circle.curvature(); // == sinhEta
80  float coshEta = sqrt(1 + sqr(cotTheta)); // == 1/sinTheta
81 
82  // Calculate globalDirs
83  float sinTheta = 1. / coshEta;
84  float cosTheta = cotTheta * sinTheta;
85 
86  int dir;
87  if(areaParallelogram(p[0] - c, p[1] - c) > 0) dir = 1; else dir = -1;
88 
89  float curvature = circle.curvature();
90 
91  for(vector<Global2DVector>::const_iterator ip = p.begin();
92  ip!= p.end(); ip++)
93  {
94  Global2DVector v = (*ip - c)*curvature*dir;
95  globalDirs.push_back(GlobalVector(-v.y()*sinTheta,
96  v.x()*sinTheta,
97  cosTheta));
98  }
99  }
100 
101  return globalDirs;
102 }
103 
104 /*****************************************************************************/
105 vector<GlobalPoint> ClusterShapeTrackFilter::getGlobalPoss
106  (const vector<const TrackingRecHit *> & recHits) const
107 {
108  vector<GlobalPoint> globalPoss;
109 
110  for(vector<const TrackingRecHit *>::const_iterator recHit = recHits.begin();
111  recHit!= recHits.end();
112  recHit++)
113  {
114  DetId detId = (*recHit)->geographicalId();
115 
116  GlobalPoint gpos =
117  theTracker->idToDet(detId)->toGlobal((*recHit)->localPosition());
118 
119  globalPoss.push_back(gpos);
120  }
121 
122  return globalPoss;
123 }
124 
125 /*****************************************************************************/
128  const vector<const TrackingRecHit *> & recHits) const
129 {
130  // Do not even look at pairs
131  if(recHits.size() <= 2) return true;
132 
133  // Check pt
134  if(track->pt() < ptMin ||
135  track->pt() > ptMax)
136  {
137  LogTrace("ClusterShapeTrackFilter")
138  << " [ClusterShapeTrackFilter] pt not in range: "
139  << ptMin << " " << track->pt() << " " << ptMax;
140  return false;
141  }
142 
143  // Get global positions
144  vector<GlobalPoint> globalPoss = getGlobalPoss(recHits);
145 
146  // Get global directions
147  vector<GlobalVector> globalDirs = getGlobalDirs(globalPoss);
148  if ( globalDirs.empty() ) return false;
149 
150  bool ok = true;
151 
152  // Check whether shape of pixel cluster is compatible
153  // with local track direction
154  for(unsigned int i = 0; i < recHits.size(); i++)
155  {
156  const SiPixelRecHit* pixelRecHit =
157  dynamic_cast<const SiPixelRecHit *>(recHits[i]);
158 
159  if(!pixelRecHit->isValid())
160  {
161  ok = false; break;
162  }
163 
164  if(! theFilter->isCompatible(*pixelRecHit, globalDirs[i], *theClusterShapeCache) )
165  {
166  LogTrace("ClusterShapeTrackFilter")
167  << " [ClusterShapeTrackFilter] clusShape problem"
168  << HitInfo::getInfo(*recHits[i],tTopo);
169 
170  ok = false; break;
171  }
172  }
173 
174  return ok;
175 }
176 
T y() const
Definition: PV2DBase.h:46
static const double slope[3]
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
const TrackerGeometry * theTracker
ClusterShapeTrackFilter(const SiPixelClusterShapeCache *cache, double ptmin, double ptmax, const edm::EventSetup &es)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
T curvature(T InversePt, const edm::EventSetup &iSetup)
const ClusterShapeHitFilter * theFilter
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:660
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
std::vector< GlobalPoint > getGlobalPoss(const std::vector< const TrackingRecHit * > &recHits) const
T y() const
Cartesian y coordinate.
#define LogTrace(id)
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, const SiPixelClusterShapeCache &clusterShapeCache, PixelData const *pd=0) const
float areaParallelogram(const Global2DVector &a, const Global2DVector &b) const
Definition: DetId.h:18
const TrackerTopology * tTopo
bool isValid() const
double b
Definition: hdecay.h:120
std::vector< GlobalVector > getGlobalDirs(const std::vector< GlobalPoint > &globalPoss) const
double ptmin
Definition: HydjetWrapper.h:90
double a
Definition: hdecay.h:121
def cache(function)
Definition: utilities.py:3
T get() const
Definition: EventSetup.h:71
const TrackerGeomDet * idToDet(DetId) const override
dbl *** dir
Definition: mlp_gen.cc:35
const SiPixelClusterShapeCache * theClusterShapeCache
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:23
T x() const
Definition: PV2DBase.h:45
T const * product() const
Definition: ESHandle.h:86
float sqr(float x)
T x() const
Cartesian x coordinate.
Global3DVector GlobalVector
Definition: GlobalVector.h:10
Our base class.
Definition: SiPixelRecHit.h:23