CMS 3D CMS Logo

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