CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelVertexProducerClusters.cc
Go to the documentation of this file.
2 
7 
10 
16 
19 
22 
24 
25 #include <vector>
26 #include <algorithm>
27 
28 using namespace std;
29 
30 namespace {
31 class VertexHit
32 {
33  public:
34  float z;
35  float r;
36  int w;
37 };
38 
39 /*****************************************************************************/
40 int getContainedHits
41  (const vector<VertexHit>& hits, float z0, float & chi)
42 {
43  int n = 0;
44  chi = 0.;
45 
46  for(vector<VertexHit>::const_iterator hit = hits.begin();
47  hit!= hits.end(); hit++)
48  {
49  // Predicted cluster width in y direction
50  float p = 2 * fabs(hit->z - z0)/hit->r + 0.5; // FIXME
51 
52  if(fabs(p - hit->w) <= 1.)
53  {
54  chi += fabs(p - hit->w);
55  n++;
56  }
57  }
58 
59  return n;
60 }
61 }
62 
63 /*****************************************************************************/
65  (const edm::ParameterSet& ps) : pixelToken_(consumes<SiPixelRecHitCollection>(edm::InputTag("siPixelRecHits")))
66 {
67  // Product
68  produces<reco::VertexCollection>();
69 }
70 
71 /*****************************************************************************/
73 {
74 }
75 
76 
77 /*****************************************************************************/
80 {
81  //Retrieve tracker topology from geometry
83  es.get<TrackerTopologyRcd>().get(tTopo);
84 
85  // Get tracker geometry
86  edm::ESHandle<TrackerGeometry> trackerHandle;
87  es.get<TrackerDigiGeometryRecord>().get(trackerHandle);
88  const TrackerGeometry* theTracker = trackerHandle.product();
89 
90 
91  // Get pixel hit collections
93  ev.getByToken(pixelToken_, pixelColl);
94 
95  const SiPixelRecHitCollection* thePixelHits = pixelColl.product();
96 
97  std::auto_ptr<reco::VertexCollection> vertices(new reco::VertexCollection);
98 
99  if(thePixelHits->size() > 0)
100  {
101  vector<VertexHit> hits;
102 
104  recHit = thePixelHits->data().begin(),
105  recHitEnd = thePixelHits->data().end();
106  recHit != recHitEnd; ++recHit)
107  {
108  if(recHit->isValid())
109  {
110 // if(!(recHit->isOnEdge() || recHit->hasBadPixels()))
111  DetId id = recHit->geographicalId();
112  const PixelGeomDetUnit* pgdu =
113  dynamic_cast<const PixelGeomDetUnit*>(theTracker->idToDetUnit(id));
114  const PixelTopology* theTopol = ( &(pgdu->specificTopology()) );
115  vector<SiPixelCluster::Pixel> pixels = recHit->cluster()->pixels();
116 
117  bool pixelOnEdge = false;
118  for(vector<SiPixelCluster::Pixel>::const_iterator
119  pixel = pixels.begin(); pixel!= pixels.end(); pixel++)
120  {
121  int pos_x = (int)pixel->x;
122  int pos_y = (int)pixel->y;
123 
124  if(theTopol->isItEdgePixelInX(pos_x) ||
125  theTopol->isItEdgePixelInY(pos_y))
126  { pixelOnEdge = true; break; }
127  }
128 
129  if(!pixelOnEdge)
130  if(id.subdetId() == int(PixelSubdetector::PixelBarrel))
131  {
132 
133 
134  LocalPoint lpos = LocalPoint(recHit->localPosition().x(),
135  recHit->localPosition().y(),
136  recHit->localPosition().z());
137 
138  GlobalPoint gpos = theTracker->idToDet(id)->toGlobal(lpos);
139 
140  VertexHit hit;
141  hit.z = gpos.z();
142  hit.r = gpos.perp();
143  hit.w = recHit->cluster()->sizeY();
144 
145  hits.push_back(hit);
146  }
147  }
148  }
149 
150  int nhits; int nhits_max = 0;
151  float chi; float chi_max = 1e+9;
152 
153  float zest = 0;
154 
155  for(float z0 = -15.9; z0 <= 15.95; z0 += 0.1)
156  {
157  nhits = getContainedHits(hits, z0, chi);
158 
159  if(nhits > 0)
160  {
161  if(nhits > nhits_max)
162  { chi_max = 1e+9; nhits_max = nhits; }
163 
164  if(nhits >= nhits_max)
165  if(chi < chi_max)
166  { chi_max = chi; zest = z0; }
167  }
168  }
169 
170  LogTrace("MinBiasTracking")
171  << " [vertex position] estimated = " << zest
172  << " | pixel barrel hits = " << thePixelHits->size();
173 
175  err(2,2) = 0.6 * 0.6;
176  reco::Vertex ver(reco::Vertex::Point(0,0,zest), err, 0, 1, 1);
177  vertices->push_back(ver);
178  }
179 
180  ev.put(vertices);
181 }
182 
const TrackerGeomDet * idToDetUnit(DetId) const
Return the pointer to the GeomDetUnit corresponding to a given DetId.
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const double w
Definition: UKUtility.cc:23
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
virtual void produce(edm::StreamID, edm::Event &ev, const edm::EventSetup &es) const override
float float float z
virtual bool isItEdgePixelInX(int ixbin) const =0
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
data_type const * data(size_t cell) const
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
#define LogTrace(id)
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:56
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
T const * product() const
Definition: ESHandle.h:86
size_type size() const
virtual bool isItEdgePixelInY(int iybin) const =0
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
PixelVertexProducerClusters(const edm::ParameterSet &ps)
const TrackerGeomDet * idToDet(DetId) const