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 
18 
21 
23 
25 
26 #include <fstream>
27 #include <iostream>
28 #include <vector>
29 #include <algorithm>
30 
31 using namespace std;
32 
33 #include "TROOT.h"
34 #include "TH1F.h"
35 #include "TF1.h"
36 
37 class VertexHit
38 {
39  public:
40  float z;
41  float r;
42  int w;
43 };
44 
45 /*****************************************************************************/
47  (const edm::ParameterSet& ps) : theConfig(ps)
48 {
49  // Product
50  produces<reco::VertexCollection>();
51 }
52 
53 
54 /*****************************************************************************/
56 {
57 }
58 
59 /*****************************************************************************/
61  (edm::Run const & run, edm::EventSetup const & es)
62 {
63  // Get tracker geometry
64  edm::ESHandle<TrackerGeometry> trackerHandle;
65  es.get<TrackerDigiGeometryRecord>().get(trackerHandle);
66  theTracker = trackerHandle.product();
67 }
68 
69 /*****************************************************************************/
71  (vector<VertexHit> hits, float z0, float & chi)
72 {
73  int n = 0;
74  chi = 0.;
75 
76  for(vector<VertexHit>::const_iterator hit = hits.begin();
77  hit!= hits.end(); hit++)
78  {
79  // Predicted cluster width in y direction
80  float p = 2 * fabs(hit->z - z0)/hit->r + 0.5; // FIXME
81 
82  if(fabs(p - hit->w) <= 1.)
83  {
84  chi += fabs(p - hit->w);
85  n++;
86  }
87  }
88 
89  return n;
90 }
91 
92 /*****************************************************************************/
94  (edm::Event& ev, const edm::EventSetup& es)
95 {
96  // Get pixel hit collections
98  ev.getByLabel("siPixelRecHits", pixelColl);
99 
100  const SiPixelRecHitCollection* thePixelHits = pixelColl.product();
101 
102  std::auto_ptr<reco::VertexCollection> vertices(new reco::VertexCollection);
103 
104  if(thePixelHits->size() > 0)
105  {
106  vector<VertexHit> hits;
107 
109  recHit = thePixelHits->data().begin(),
110  recHitEnd = thePixelHits->data().end();
111  recHit != recHitEnd; ++recHit)
112  {
113  if(recHit->isValid())
114  {
115 // if(!(recHit->isOnEdge() || recHit->hasBadPixels()))
116  DetId id = recHit->geographicalId();
117  const PixelGeomDetUnit* pgdu =
118  dynamic_cast<const PixelGeomDetUnit*>(theTracker->idToDetUnit(id));
119  const PixelTopology* theTopol = ( &(pgdu->specificTopology()) );
120  vector<SiPixelCluster::Pixel> pixels = recHit->cluster()->pixels();
121 
122  bool pixelOnEdge = false;
123  for(vector<SiPixelCluster::Pixel>::const_iterator
124  pixel = pixels.begin(); pixel!= pixels.end(); pixel++)
125  {
126  int pos_x = (int)pixel->x;
127  int pos_y = (int)pixel->y;
128 
129  if(theTopol->isItEdgePixelInX(pos_x) ||
130  theTopol->isItEdgePixelInY(pos_y))
131  { pixelOnEdge = true; break; }
132  }
133 
134  if(!pixelOnEdge)
135  if(id.subdetId() == int(PixelSubdetector::PixelBarrel))
136  {
137  PXBDetId pid(id);
138 
139  LocalPoint lpos = LocalPoint(recHit->localPosition().x(),
140  recHit->localPosition().y(),
141  recHit->localPosition().z());
142 
143  GlobalPoint gpos = theTracker->idToDet(id)->toGlobal(lpos);
144 
145  VertexHit hit;
146  hit.z = gpos.z();
147  hit.r = gpos.perp();
148  hit.w = recHit->cluster()->sizeY();
149 
150  hits.push_back(hit);
151  }
152  }
153  }
154 
155  int nhits; int nhits_max = 0;
156  float chi; float chi_max = 1e+9;
157 
158  float zest = 0;
159 
160  for(float z0 = -15.9; z0 <= 15.95; z0 += 0.1)
161  {
162  nhits = getContainedHits(hits, z0, chi);
163 
164  if(nhits > 0)
165  {
166  if(nhits > nhits_max)
167  { chi_max = 1e+9; nhits_max = nhits; }
168 
169  if(nhits >= nhits_max)
170  if(chi < chi_max)
171  { chi_max = chi; zest = z0; }
172  }
173  }
174 
175  LogTrace("MinBiasTracking")
176  << " [vertex position] estimated = " << zest
177  << " | pixel barrel hits = " << thePixelHits->size();
178 
180  err(2,2) = 0.6 * 0.6;
181  reco::Vertex ver(reco::Vertex::Point(0,0,zest), err, 0, 1, 1);
182  vertices->push_back(ver);
183  }
184 
185  ev.put(vertices);
186 }
187 
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
void beginRun(edm::Run const &run, edm::EventSetup const &es)
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
virtual bool isItEdgePixelInX(int ixbin) const =0
int getContainedHits(std::vector< VertexHit > hits, float z0, float &chi)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
data_type const * data(size_t cell) const
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
#define LogTrace(id)
Definition: DetId.h:20
const T & get() const
Definition: EventSetup.h:55
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
T const * product() const
Definition: ESHandle.h:62
T const * product() const
Definition: Handle.h:74
size_type size() const
virtual bool isItEdgePixelInY(int iybin) const =0
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
virtual void produce(edm::Event &ev, const edm::EventSetup &es)
PixelVertexProducerClusters(const edm::ParameterSet &ps)
Definition: Run.h:33