CMS 3D CMS Logo

HIPixelClusterVtxProducer.cc
Go to the documentation of this file.
2 
7 
10 
16 
19 
21 
22 #include <fstream>
23 #include <iostream>
24 #include <vector>
25 #include <algorithm>
26 
27 /*****************************************************************************/
29  : srcPixelsString_(ps.getParameter<std::string>("pixelRecHits")),
30  minZ_(ps.getParameter<double>("minZ")),
31  maxZ_(ps.getParameter<double>("maxZ")),
32  zStep_(ps.getParameter<double>("zStep"))
33 
34 {
35  // Constructor
36  produces<reco::VertexCollection>();
37  srcPixels_ = (consumes<SiPixelRecHitCollection>(srcPixelsString_));
38 }
39 
40 /*****************************************************************************/
42  // Destructor
43 }
44 
45 /*****************************************************************************/
47  // new vertex collection
48  auto vertices = std::make_unique<reco::VertexCollection>();
49 
50  // get pixel rechits
52  ev.getByToken(srcPixels_, hRecHits);
53 
54  // get tracker geometry
55  if (hRecHits.isValid()) {
56  edm::ESHandle<TrackerGeometry> trackerHandle;
57  es.get<TrackerDigiGeometryRecord>().get(trackerHandle);
58  const TrackerGeometry *tgeo = trackerHandle.product();
59  const SiPixelRecHitCollection *hits = hRecHits.product();
60 
61  // loop over pixel rechits
62  std::vector<VertexHit> vhits;
63  for (SiPixelRecHitCollection::DataContainer::const_iterator hit = hits->data().begin(), end = hits->data().end();
64  hit != end;
65  ++hit) {
66  if (!hit->isValid())
67  continue;
68  DetId id(hit->geographicalId());
69  if (id.subdetId() != int(PixelSubdetector::PixelBarrel))
70  continue;
71  const PixelGeomDetUnit *pgdu = static_cast<const PixelGeomDetUnit *>(tgeo->idToDet(id));
72  if (true) {
73  const PixelTopology *pixTopo = &(pgdu->specificTopology());
74  std::vector<SiPixelCluster::Pixel> pixels(hit->cluster()->pixels());
75  bool pixelOnEdge = false;
76  for (std::vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin(); pixel != pixels.end();
77  ++pixel) {
78  int pixelX = pixel->x;
79  int pixelY = pixel->y;
80  if (pixTopo->isItEdgePixelInX(pixelX) || pixTopo->isItEdgePixelInY(pixelY)) {
81  pixelOnEdge = true;
82  break;
83  }
84  }
85  if (pixelOnEdge)
86  continue;
87  }
88 
89  LocalPoint lpos = LocalPoint(hit->localPosition().x(), hit->localPosition().y(), hit->localPosition().z());
90  GlobalPoint gpos = pgdu->toGlobal(lpos);
91  VertexHit vh;
92  vh.z = gpos.z();
93  vh.r = gpos.perp();
94  vh.w = hit->cluster()->sizeY();
95  vhits.push_back(vh);
96  }
97 
98  // estimate z-position from cluster lengths
99  double zest = 0.0;
100  int nhits = 0, nhits_max = 0;
101  double chi = 0, chi_max = 1e+9;
102  for (double z0 = minZ_; z0 <= maxZ_; z0 += zStep_) {
103  nhits = getContainedHits(vhits, z0, chi);
104  if (nhits == 0)
105  continue;
106  if (nhits > nhits_max) {
107  chi_max = 1e+9;
108  nhits_max = nhits;
109  }
110  if (nhits >= nhits_max && chi < chi_max) {
111  chi_max = chi;
112  zest = z0;
113  }
114  }
115 
116  LogTrace("MinBiasTracking") << " [vertex position] estimated = " << zest
117  << " | pixel barrel hits = " << vhits.size();
118 
119  // put 1-d vertex and dummy errors into collection
121  err(2, 2) = 0.6 * 0.6;
122  reco::Vertex ver(reco::Vertex::Point(0, 0, zest), err, 0, 1, 1);
123  vertices->push_back(ver);
124  }
125 
126  ev.put(std::move(vertices));
127 }
128 
129 /*****************************************************************************/
130 int HIPixelClusterVtxProducer::getContainedHits(const std::vector<VertexHit> &hits, double z0, double &chi) {
131  // Calculate number of hits contained in v-shaped window in cluster y-width vs. z-position.
132  int n = 0;
133  chi = 0.;
134 
135  for (std::vector<VertexHit>::const_iterator hit = hits.begin(); hit != hits.end(); hit++) {
136  double p = 2 * fabs(hit->z - z0) / hit->r + 0.5; // FIXME
137  if (TMath::Abs(p - hit->w) <= 1.) {
138  chi += fabs(p - hit->w);
139  n++;
140  }
141  }
142  return n;
143 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
int getContainedHits(const std::vector< VertexHit > &hits, double z0, double &chi)
HIPixelClusterVtxProducer(const edm::ParameterSet &ps)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
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
void produce(edm::Event &ev, const edm::EventSetup &es) override
T Abs(T a)
Definition: MathUtil.h:49
data_type const * data(size_t cell) const
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
#define end
Definition: vmac.h:39
edm::EDGetTokenT< SiPixelRecHitCollection > srcPixels_
#define LogTrace(id)
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
T get() const
Definition: EventSetup.h:73
const TrackerGeomDet * idToDet(DetId) const override
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511