CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LowPtClusterShapeSeedComparitor.cc
Go to the documentation of this file.
2 
6 
13 
18 
20 
22 
26 
27 #include<cmath>
28 
29 
30 
31 namespace {
32  typedef Basic2DVector<float> Vector2D;
33 
34  inline float sqr(float x) { return x*x; }
35 
36  /*****************************************************************************/
37  inline
38  float areaParallelogram
39  (const Vector2D& a, const Vector2D& b)
40  {
41  return a.x() * b.y() - a.y() * b.x();
42  }
43 
44  /*****************************************************************************/
45 
46  inline
47  bool getGlobalDirs(GlobalPoint const * g, GlobalVector * globalDirs)
48  {
49 
50 
51  // Determine circle
52  CircleFromThreePoints circle(g[0],g[1],g[2]);
53 
54  float curvature = circle.curvature();
55  if(0.f == curvature) {
56  LogDebug("LowPtClusterShapeSeedComparitor")<<"the curvature is null:"
57  <<"\n point1: "<<g[0]
58  <<"\n point2: "<<g[1]
59  <<"\n point3: "<<g[2];
60  return false;
61  }
62 
63  // Get 2d points
64  Vector2D p[3];
65  Vector2D c = circle.center();
66  for(int i=0; i!=3; i++)
67  p[i] = g[i].basicVector().xy() -c;
68 
69 
70  float area = std::abs(areaParallelogram(p[1] - p[0], p[1]));
71 
72  float a12 = std::asin(std::min(area*curvature*curvature,1.f));
73 
74  float slope = (g[1].z() - g[0].z()) / a12;
75 
76  // Calculate globalDirs
77 
78  float cotTheta = slope * curvature;
79  float sinTheta = 1.f/std::sqrt(1.f + sqr(cotTheta));
80  float cosTheta = cotTheta*sinTheta;
81 
82  if (areaParallelogram(p[0], p[1] ) < 0) sinTheta = - sinTheta;
83 
84  for(int i = 0; i!=3; i++) {
85  Vector2D vl = p[i]*(curvature*sinTheta);
86  globalDirs[i] = GlobalVector(-vl.y(),
87  vl.x(),
88  cosTheta
89  );
90  }
91  return true;
92  }
93 
94  /*****************************************************************************/
95 
96 
97  inline
98  void getGlobalPos(const SeedingHitSet &hits, GlobalPoint * globalPoss)
99  {
100 
101  for(unsigned int i=0; i!=hits.size(); ++i)
102  globalPoss[i] = hits[i]->globalPosition();
103  }
104 
105 } // namespace
106 
107 /*****************************************************************************/
109  thePixelClusterShapeCacheToken(iC.consumes<SiPixelClusterShapeCache>(ps.getParameter<edm::InputTag>("clusterShapeCacheSrc")))
110 {}
111 
112 /*****************************************************************************/
114  es.get<CkfComponentsRecord>().get("ClusterShapeHitFilter", theShapeFilter);
115  es.get<TrackerTopologyRcd>().get(theTTopo);
116 
118 }
119 
121 //(const reco::Track* track, const vector<const TrackingRecHit *> & recHits) const
122 {
123  assert(hits.size()==3);
124 
126  assert(filter != 0 && "LowPtClusterShapeSeedComparitor: init(EventSetup) method was not called");
127 
128  // Get global positions
129  GlobalPoint globalPoss[3];
130  getGlobalPos(hits, globalPoss);
131 
132  // Get global directions
133  GlobalVector globalDirs[3];
134 
135  bool ok = getGlobalDirs(globalPoss,globalDirs);
136 
137  // Check whether shape of pixel cluster is compatible
138  // with local track direction
139 
140  if (!ok)
141  {
142  LogDebug("LowPtClusterShapeSeedComparitor")<<"curvarture 0:"
143  <<"\nnHits: "<<hits.size()
144  <<" will say the seed is good anyway.";
145  return true;
146  }
147 
148  for(int i = 0; i < 3; i++)
149  {
150  const SiPixelRecHit* pixelRecHit =
151  dynamic_cast<const SiPixelRecHit *>(hits[i]->hit());
152 
153  if (!pixelRecHit){
154  edm::LogError("LowPtClusterShapeSeedComparitor")<<"this is not a pixel cluster";
155  ok = false; break;
156  }
157 
158  if(!pixelRecHit->isValid())
159  {
160  ok = false; break;
161  }
162 
163  LogDebug("LowPtClusterShapeSeedComparitor")<<"about to compute compatibility."
164  <<"hit ptr: "<<pixelRecHit
165  <<"global direction:"<< globalDirs[i];
166 
167 
168  if(! filter->isCompatible(*pixelRecHit, globalDirs[i], *thePixelClusterShapeCache) )
169  {
170  LogTrace("LowPtClusterShapeSeedComparitor")
171  << " clusShape is not compatible"
172  << HitInfo::getInfo(*hits[i]->hit(),theTTopo.product());
173 
174  ok = false; break;
175  }
176  }
177 
178  return ok;
179 }
180 
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
static const double slope[3]
assert(m_qm.get())
edm::ESHandle< ClusterShapeHitFilter > theShapeFilter
something
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
edm::ESHandle< TrackerTopology > theTTopo
T x() const
Cartesian x coordinate.
T curvature(T InversePt, const edm::EventSetup &iSetup)
virtual bool compatible(const SeedingHitSet &hits, const TrackingRegion &region) const
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
LowPtClusterShapeSeedComparitor(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
edm::Handle< SiPixelClusterShapeCache > thePixelClusterShapeCache
virtual void init(const edm::Event &e, const edm::EventSetup &es)
#define LogTrace(id)
edm::EDGetTokenT< SiPixelClusterShapeCache > thePixelClusterShapeCacheToken
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
string const
Definition: compareJSON.py:14
double a
Definition: hdecay.h:121
Square< F >::type sqr(const F &f)
Definition: Square.h:13
unsigned int size() const
Definition: SeedingHitSet.h:44
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:23
Global3DVector GlobalVector
Definition: GlobalVector.h:10
Our base class.
Definition: SiPixelRecHit.h:23