CMS 3D CMS Logo

ClusterCompatibilityProducer.cc
Go to the documentation of this file.
1 //
2 // Derived from HLTrigger/special/src/HLTPixelClusterShapeFilter.cc
3 // at version 7_5_0_pre3
4 //
5 // Original Author (of Derivative Producer): Eric Appelt
6 // Created: Mon Apr 27, 2015
7 
8 #include <iostream>
9 
19 
26 
32 
33 //
34 // class declaration
35 //
36 
38 public:
41 
42  void produce(edm::Event &, const edm::EventSetup &) override;
43 
44 private:
47  edm::InputTag inputTag_; // input tag identifying product containing pixel clusters
48  double minZ_; // beginning z-vertex position
49  double maxZ_; // end z-vertex position
50  double zStep_; // size of steps in z-vertex test
51 
52  struct VertexHit {
53  float z;
54  float r;
55  float w;
56  };
57 
58  struct ContainedHits {
59  float z0;
60  int nHit;
61  float chi;
62  };
63 
64  ContainedHits getContainedHits(const std::vector<VertexHit> &hits, double z0) const;
65 };
66 
68  : inputTag_(config.getParameter<edm::InputTag>("inputTag")),
69  minZ_(config.getParameter<double>("minZ")),
70  maxZ_(config.getParameter<double>("maxZ")),
71  zStep_(config.getParameter<double>("zStep")) {
72  inputToken_ = consumes<SiPixelRecHitCollection>(inputTag_);
73  trackerToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
74  LogDebug("") << "Using the " << inputTag_ << " input collection";
75  produces<reco::ClusterCompatibility>();
76 }
77 
79 
81  auto creco = std::make_unique<reco::ClusterCompatibility>();
82 
83  // get hold of products from Event
85  iEvent.getByToken(inputToken_, hRecHits);
86 
87  // get tracker geometry
88  if (hRecHits.isValid()) {
90  const TrackerGeometry *tgeo = trackerHandle.product();
91  const SiPixelRecHitCollection *hits = hRecHits.product();
92 
93  // loop over pixel rechits
94  int nPxlHits = 0;
95  std::vector<VertexHit> vhits;
96  for (SiPixelRecHitCollection::DataContainer::const_iterator hit = hits->data().begin(), end = hits->data().end();
97  hit != end;
98  ++hit) {
99  if (!hit->isValid())
100  continue;
101  ++nPxlHits;
102  DetId id(hit->geographicalId());
103  if (id.subdetId() != int(PixelSubdetector::PixelBarrel))
104  continue;
105  const PixelGeomDetUnit *pgdu = static_cast<const PixelGeomDetUnit *>(tgeo->idToDet(id));
106  const PixelTopology *pixTopo = &(pgdu->specificTopology());
107  std::vector<SiPixelCluster::Pixel> pixels(hit->cluster()->pixels());
108  bool pixelOnEdge = false;
109  for (std::vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin(); pixel != pixels.end(); ++pixel) {
110  int pixelX = pixel->x;
111  int pixelY = pixel->y;
112  if (pixTopo->isItEdgePixelInX(pixelX) || pixTopo->isItEdgePixelInY(pixelY)) {
113  pixelOnEdge = true;
114  break;
115  }
116  }
117  if (pixelOnEdge)
118  continue;
119 
120  LocalPoint lpos = LocalPoint(hit->localPosition().x(), hit->localPosition().y(), hit->localPosition().z());
121  GlobalPoint gpos = pgdu->toGlobal(lpos);
122  VertexHit vh;
123  vh.z = gpos.z();
124  vh.r = gpos.perp();
125  vh.w = hit->cluster()->sizeY();
126  vhits.push_back(vh);
127  }
128 
129  creco->setNValidPixelHits(nPxlHits);
130 
131  // append cluster compatibility for each z-position
132  for (double z0 = minZ_; z0 <= maxZ_; z0 += zStep_) {
134  creco->append(c.z0, c.nHit, c.chi);
135  }
136  }
137  iEvent.put(std::move(creco));
138 }
139 
141  const std::vector<VertexHit> &hits, double z0) const {
142  // Calculate number of hits contained in v-shaped window in cluster y-width vs. z-position.
143  int n = 0;
144  double chi = 0.;
145 
146  for (std::vector<VertexHit>::const_iterator hit = hits.begin(); hit != hits.end(); hit++) {
147  // the calculation of the predicted cluster width p was
148  // marked 'FIXME' in the HLTPixelClusterShapeFilter. It should
149  // be revisited but is retained as it was for compatibility with the
150  // older filter.
151  double p = 2 * fabs(hit->z - z0) / hit->r + 0.5;
152  if (fabs(p - hit->w) <= 1.) {
153  chi += fabs(p - hit->w);
154  n++;
155  }
156  }
158  output.z0 = z0;
159  output.nHit = n;
160  output.chi = chi;
161  return output;
162 }
163 
164 //define this as a plug-in
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
Definition: config.py:1
virtual bool isItEdgePixelInX(int ixbin) const =0
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
ClusterCompatibilityProducer(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ContainedHits getContainedHits(const std::vector< VertexHit > &hits, double z0) const
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
const TrackerGeomDet * idToDet(DetId) const override
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerToken_
edm::EDGetTokenT< SiPixelRecHitCollection > inputToken_
virtual bool isItEdgePixelInY(int iybin) const =0
HLT enums.
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
Definition: output.py:1
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)