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