CMS 3D CMS Logo

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  public:
33  float z;
34  float r;
35  int w;
36  };
37 
38  /*****************************************************************************/
39  int getContainedHits(const vector<VertexHit>& hits, float z0, float& chi) {
40  int n = 0;
41  chi = 0.;
42 
43  for (vector<VertexHit>::const_iterator hit = hits.begin(); hit != hits.end(); hit++) {
44  // Predicted cluster width in y direction
45  float p = 2 * fabs(hit->z - z0) / hit->r + 0.5; // FIXME
46 
47  if (fabs(p - hit->w) <= 1.) {
48  chi += fabs(p - hit->w);
49  n++;
50  }
51  }
52 
53  return n;
54  }
55 } // namespace
56 
57 /*****************************************************************************/
59  : pixelToken_(consumes<SiPixelRecHitCollection>(edm::InputTag("siPixelRecHits"))) {
60  // Product
61  produces<reco::VertexCollection>();
62 }
63 
64 /*****************************************************************************/
66 
67 /*****************************************************************************/
69  //Retrieve tracker topology from geometry
71  es.get<TrackerTopologyRcd>().get(tTopo);
72 
73  // Get tracker geometry
74  edm::ESHandle<TrackerGeometry> trackerHandle;
75  es.get<TrackerDigiGeometryRecord>().get(trackerHandle);
76  const TrackerGeometry* theTracker = trackerHandle.product();
77 
78  // Get pixel hit collections
80  ev.getByToken(pixelToken_, pixelColl);
81 
82  const SiPixelRecHitCollection* thePixelHits = pixelColl.product();
83 
84  auto vertices = std::make_unique<reco::VertexCollection>();
85 
86  if (!thePixelHits->empty()) {
87  vector<VertexHit> hits;
88 
90  recHitEnd = thePixelHits->data().end();
91  recHit != recHitEnd;
92  ++recHit) {
93  if (recHit->isValid()) {
94  // if(!(recHit->isOnEdge() || recHit->hasBadPixels()))
95  DetId id = recHit->geographicalId();
96  const PixelGeomDetUnit* pgdu = dynamic_cast<const PixelGeomDetUnit*>(theTracker->idToDetUnit(id));
97  const PixelTopology* theTopol = (&(pgdu->specificTopology()));
98  vector<SiPixelCluster::Pixel> pixels = recHit->cluster()->pixels();
99 
100  bool pixelOnEdge = false;
101  for (vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin(); pixel != pixels.end(); pixel++) {
102  int pos_x = (int)pixel->x;
103  int pos_y = (int)pixel->y;
104 
105  if (theTopol->isItEdgePixelInX(pos_x) || theTopol->isItEdgePixelInY(pos_y)) {
106  pixelOnEdge = true;
107  break;
108  }
109  }
110 
111  if (!pixelOnEdge)
112  if (id.subdetId() == int(PixelSubdetector::PixelBarrel)) {
113  LocalPoint lpos =
114  LocalPoint(recHit->localPosition().x(), recHit->localPosition().y(), recHit->localPosition().z());
115 
116  GlobalPoint gpos = theTracker->idToDet(id)->toGlobal(lpos);
117 
118  VertexHit hit;
119  hit.z = gpos.z();
120  hit.r = gpos.perp();
121  hit.w = recHit->cluster()->sizeY();
122 
123  hits.push_back(hit);
124  }
125  }
126  }
127 
128  int nhits;
129  int nhits_max = 0;
130  float chi;
131  float chi_max = 1e+9;
132 
133  float zest = 0;
134 
135  for (float z0 = -15.9; z0 <= 15.95; z0 += 0.1) {
136  nhits = getContainedHits(hits, z0, chi);
137 
138  if (nhits > 0) {
139  if (nhits > nhits_max) {
140  chi_max = 1e+9;
141  nhits_max = nhits;
142  }
143 
144  if (nhits >= nhits_max)
145  if (chi < chi_max) {
146  chi_max = chi;
147  zest = z0;
148  }
149  }
150  }
151 
152  LogTrace("MinBiasTracking") << " [vertex position] estimated = " << zest
153  << " | pixel barrel hits = " << thePixelHits->size();
154 
156  err(2, 2) = 0.6 * 0.6;
157  reco::Vertex ver(reco::Vertex::Point(0, 0, zest), err, 0, 1, 1);
158  vertices->push_back(ver);
159  }
160 
161  ev.put(std::move(vertices));
162 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
const double w
Definition: UKUtility.cc:23
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::EDGetTokenT< SiPixelRecHitCollection > pixelToken_
virtual bool isItEdgePixelInX(int ixbin) const =0
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
virtual bool isItEdgePixelInY(int iybin) const =0
bool ev
float float float z
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
data_type const * data(size_t cell) const
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
#define LogTrace(id)
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
T const * product() const
Definition: Handle.h:69
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
size_type size() const
HLT enums.
T get() const
Definition: EventSetup.h:73
void produce(edm::StreamID, edm::Event &ev, const edm::EventSetup &es) const override
const TrackerGeomDet * idToDet(DetId) const override
PixelVertexProducerClusters(const edm::ParameterSet &ps)
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511