CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ClusterShapeTrackFilter.cc
Go to the documentation of this file.
2 
6 
12 
16 
20 
23 
27 
28 inline float sqr(float x) { return x*x; }
29 
30 using namespace std;
31 
32 /*****************************************************************************/
34  theClusterShapeCacheToken(iC.consumes<SiPixelClusterShapeCache>(ps.getParameter<edm::InputTag>("clusterShapeCacheSrc"))),
35  theTracker(nullptr),
36  theFilter(nullptr),
37  theClusterShapeCache(nullptr)
38 {
39  // Get ptMin if available
40  ptMin = (ps.exists("ptMin") ? ps.getParameter<double>("ptMin") : 0.);
41  ptMax = (ps.exists("ptMax") ? ps.getParameter<double>("ptMax") : 999999.);
42 }
43 
44 /*****************************************************************************/
46 {
47 }
48 
49 /*****************************************************************************/
53  theClusterShapeCache = cache.product();
54 
55  // Get tracker geometry
57  es.get<TrackerDigiGeometryRecord>().get(tracker);
58  theTracker = tracker.product();
59 
60  // Get cluster shape hit filter
62  es.get<CkfComponentsRecord>().get("ClusterShapeHitFilter",shape);
63  theFilter = shape.product();
64 }
65 
66 /*****************************************************************************/
68  (const Global2DVector& a, const Global2DVector& b) const
69 {
70  return a.x() * b.y() - a.y() * b.x();
71 }
72 
73 /*****************************************************************************/
74 vector<GlobalVector> ClusterShapeTrackFilter::getGlobalDirs
75  (const vector<GlobalPoint> & g) const
76 {
77  // Get 2d points
78  vector<Global2DVector> p;
79  for(vector<GlobalPoint>::const_iterator ig = g.begin();
80  ig!= g.end(); ig++)
81  p.push_back( Global2DVector(ig->x(), ig->y()) );
82 
83  //
84  vector<GlobalVector> globalDirs;
85 
86  // Determine circle
87  CircleFromThreePoints circle(g[0],g[1],g[2]);
88 
89  if(circle.curvature() != 0.)
90  {
91  Global2DVector c (circle.center().x(), circle.center().y());
92 
93  float rad2 = (p[0] - c).mag2();
94  float a12 = asin(fabsf(areaParallelogram(p[0] - c, p[1] - c)) / rad2);
95 
96  float slope = (g[1].z() - g[0].z()) / a12;
97 
98  float cotTheta = slope * circle.curvature(); // == sinhEta
99  float coshEta = sqrt(1 + sqr(cotTheta)); // == 1/sinTheta
100 
101  // Calculate globalDirs
102  float sinTheta = 1. / coshEta;
103  float cosTheta = cotTheta * sinTheta;
104 
105  int dir;
106  if(areaParallelogram(p[0] - c, p[1] - c) > 0) dir = 1; else dir = -1;
107 
108  float curvature = circle.curvature();
109 
110  for(vector<Global2DVector>::const_iterator ip = p.begin();
111  ip!= p.end(); ip++)
112  {
113  Global2DVector v = (*ip - c)*curvature*dir;
114  globalDirs.push_back(GlobalVector(-v.y()*sinTheta,
115  v.x()*sinTheta,
116  cosTheta));
117  }
118  }
119 
120  return globalDirs;
121 }
122 
123 /*****************************************************************************/
124 vector<GlobalPoint> ClusterShapeTrackFilter::getGlobalPoss
125  (const vector<const TrackingRecHit *> & recHits) const
126 {
127  vector<GlobalPoint> globalPoss;
128 
129  for(vector<const TrackingRecHit *>::const_iterator recHit = recHits.begin();
130  recHit!= recHits.end();
131  recHit++)
132  {
133  DetId detId = (*recHit)->geographicalId();
134 
135  GlobalPoint gpos =
136  theTracker->idToDet(detId)->toGlobal((*recHit)->localPosition());
137 
138  globalPoss.push_back(gpos);
139  }
140 
141  return globalPoss;
142 }
143 
144 /*****************************************************************************/
146  (const reco::Track* track,
147  const vector<const TrackingRecHit *> & recHits,
148  const TrackerTopology *tTopo ) const
149 {
150  // Do not even look at pairs
151  if(recHits.size() <= 2) return true;
152 
153  // Check pt
154  if(track->pt() < ptMin ||
155  track->pt() > ptMax)
156  {
157  LogTrace("ClusterShapeTrackFilter")
158  << " [ClusterShapeTrackFilter] pt not in range: "
159  << ptMin << " " << track->pt() << " " << ptMax;
160  return false;
161  }
162 
163  // Get global positions
164  vector<GlobalPoint> globalPoss = getGlobalPoss(recHits);
165 
166  // Get global directions
167  vector<GlobalVector> globalDirs = getGlobalDirs(globalPoss);
168 
169  bool ok = true;
170 
171  // Check whether shape of pixel cluster is compatible
172  // with local track direction
173  for(unsigned int i = 0; i < recHits.size(); i++)
174  {
175  const SiPixelRecHit* pixelRecHit =
176  dynamic_cast<const SiPixelRecHit *>(recHits[i]);
177 
178  if(!pixelRecHit->isValid())
179  {
180  ok = false; break;
181  }
182 
183  if(! theFilter->isCompatible(*pixelRecHit, globalDirs[i], *theClusterShapeCache) )
184  {
185  LogTrace("ClusterShapeTrackFilter")
186  << " [ClusterShapeTrackFilter] clusShape problem"
187  << HitInfo::getInfo(*recHits[i],tTopo);
188 
189  ok = false; break;
190  }
191  }
192 
193  return ok;
194 }
195 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
T y() const
Definition: PV2DBase.h:46
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
static const double slope[3]
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool ev
const TrackerGeometry * theTracker
#define nullptr
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 x() const
Cartesian x coordinate.
T curvature(T InversePt, const edm::EventSetup &iSetup)
const ClusterShapeHitFilter * theFilter
T sqrt(T t)
Definition: SSEVec.h:48
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.
edm::EDGetTokenT< SiPixelClusterShapeCache > theClusterShapeCacheToken
#define LogTrace(id)
float areaParallelogram(const Global2DVector &a, const Global2DVector &b) const
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
bool isValid() const
double b
Definition: hdecay.h:120
std::vector< GlobalVector > getGlobalDirs(const std::vector< GlobalPoint > &globalPoss) const
ClusterShapeTrackFilter(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
double a
Definition: hdecay.h:121
Square< F >::type sqr(const F &f)
Definition: Square.h:13
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 x() const
Cartesian x coordinate.
Global3DVector GlobalVector
Definition: GlobalVector.h:10
void update(const edm::Event &ev, const edm::EventSetup &es) override
Our base class.
Definition: SiPixelRecHit.h:23