CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HIPixelClusterVtxProducer.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 /*****************************************************************************/
33  : srcPixels_(ps.getParameter<std::string>("pixelRecHits")),
34  minZ_(ps.getParameter<double>("minZ")),
35  maxZ_(ps.getParameter<double>("maxZ")),
36  zStep_(ps.getParameter<double>("zStep"))
37 
38 {
39  // Constructor
40  produces<reco::VertexCollection>();
41 }
42 
43 
44 /*****************************************************************************/
46 {
47  // Destructor
48 }
49 
50 
51 /*****************************************************************************/
53 {
54 
55  // new vertex collection
56  std::auto_ptr<reco::VertexCollection> vertices(new reco::VertexCollection);
57 
58  // get pixel rechits
60  try {
61  ev.getByLabel(edm::InputTag(srcPixels_),hRecHits);
62  } catch (...) {}
63 
64  // get tracker geometry
65  if (hRecHits.isValid()) {
66  edm::ESHandle<TrackerGeometry> trackerHandle;
67  es.get<TrackerDigiGeometryRecord>().get(trackerHandle);
68  const TrackerGeometry *tgeo = trackerHandle.product();
69  const SiPixelRecHitCollection *hits = hRecHits.product();
70 
71  // loop over pixel rechits
72  std::vector<VertexHit> vhits;
74  end = hits->data().end(); hit != end; ++hit) {
75  if (!hit->isValid())
76  continue;
77  DetId id(hit->geographicalId());
78  if(id.subdetId() != int(PixelSubdetector::PixelBarrel))
79  continue;
80  const PixelGeomDetUnit *pgdu = static_cast<const PixelGeomDetUnit*>(tgeo->idToDet(id));
81  if (1) {
82  const PixelTopology *pixTopo = &(pgdu->specificTopology());
83  std::vector<SiPixelCluster::Pixel> pixels(hit->cluster()->pixels());
84  bool pixelOnEdge = false;
85  for(std::vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin();
86  pixel != pixels.end(); ++pixel) {
87  int pixelX = pixel->x;
88  int pixelY = pixel->y;
89  if(pixTopo->isItEdgePixelInX(pixelX) || pixTopo->isItEdgePixelInY(pixelY)) {
90  pixelOnEdge = true;
91  break;
92  }
93  }
94  if (pixelOnEdge)
95  continue;
96  }
97 
98  LocalPoint lpos = LocalPoint(hit->localPosition().x(),
99  hit->localPosition().y(),
100  hit->localPosition().z());
101  GlobalPoint gpos = pgdu->toGlobal(lpos);
102  VertexHit vh;
103  vh.z = gpos.z();
104  vh.r = gpos.perp();
105  vh.w = hit->cluster()->sizeY();
106  vhits.push_back(vh);
107  }
108 
109  // estimate z-position from cluster lengths
110  double zest = 0.0;
111  int nhits = 0, nhits_max = 0;
112  double chi = 0, chi_max = 1e+9;
113  for(double z0 = minZ_; z0 <= maxZ_; z0 += zStep_) {
114  nhits = getContainedHits(vhits, z0, chi);
115  if(nhits == 0)
116  continue;
117  if(nhits > nhits_max) {
118  chi_max = 1e+9;
119  nhits_max = nhits;
120  }
121  if(nhits >= nhits_max && chi < chi_max) {
122  chi_max = chi;
123  zest = z0;
124  }
125  }
126 
127  LogTrace("MinBiasTracking")
128  << " [vertex position] estimated = " << zest
129  << " | pixel barrel hits = " << vhits.size();
130 
131  // put 1-d vertex and dummy errors into collection
133  err(2,2) = 0.6 * 0.6;
134  reco::Vertex ver(reco::Vertex::Point(0,0,zest), err, 0, 1, 1);
135  vertices->push_back(ver);
136  }
137 
138  ev.put(vertices);
139 }
140 
141 /*****************************************************************************/
142 int HIPixelClusterVtxProducer::getContainedHits(const std::vector<VertexHit> &hits, double z0, double &chi)
143 {
144  // Calculate number of hits contained in v-shaped window in cluster y-width vs. z-position.
145  int n = 0;
146  chi = 0.;
147 
148  for(std::vector<VertexHit>::const_iterator hit = hits.begin(); hit!= hits.end(); hit++) {
149  double p = 2 * fabs(hit->z - z0)/hit->r + 0.5; // FIXME
150  if(TMath::Abs(p - hit->w) <= 1.) {
151  chi += fabs(p - hit->w);
152  n++;
153  }
154  }
155  return n;
156 }
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
int getContainedHits(const std::vector< VertexHit > &hits, double z0, double &chi)
HIPixelClusterVtxProducer(const edm::ParameterSet &ps)
virtual void produce(edm::Event &ev, const edm::EventSetup &es)
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
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
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
#define end
Definition: vmac.h:38
bool isValid() const
Definition: HandleBase.h:76
virtual const GeomDet * idToDet(DetId) const
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
virtual bool isItEdgePixelInY(int iybin) const =0
Local3DPoint LocalPoint
Definition: LocalPoint.h:11