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 
18 
25 
31 
32 //
33 // class declaration
34 //
35 
37 public:
40 
41  void produce(edm::Event &, const edm::EventSetup &) override;
42 
43 private:
45  edm::InputTag inputTag_; // input tag identifying product containing pixel clusters
46  double minZ_; // beginning z-vertex position
47  double maxZ_; // end z-vertex position
48  double zStep_; // size of steps in z-vertex test
49 
50  struct VertexHit {
51  float z;
52  float r;
53  float w;
54  };
55 
56  struct ContainedHits {
57  float z0;
58  int nHit;
59  float chi;
60  };
61 
62  ContainedHits getContainedHits(const std::vector<VertexHit> &hits, double z0) const;
63 };
64 
66  : inputTag_(config.getParameter<edm::InputTag>("inputTag")),
67  minZ_(config.getParameter<double>("minZ")),
68  maxZ_(config.getParameter<double>("maxZ")),
69  zStep_(config.getParameter<double>("zStep")) {
70  inputToken_ = consumes<SiPixelRecHitCollection>(inputTag_);
71  LogDebug("") << "Using the " << inputTag_ << " input collection";
72  produces<reco::ClusterCompatibility>();
73 }
74 
76 
78  auto creco = std::make_unique<reco::ClusterCompatibility>();
79 
80  // get hold of products from Event
82  iEvent.getByToken(inputToken_, hRecHits);
83 
84  // get tracker geometry
85  if (hRecHits.isValid()) {
86  edm::ESHandle<TrackerGeometry> trackerHandle;
87  iSetup.get<TrackerDigiGeometryRecord>().get(trackerHandle);
88  const TrackerGeometry *tgeo = trackerHandle.product();
89  const SiPixelRecHitCollection *hits = hRecHits.product();
90 
91  // loop over pixel rechits
92  int nPxlHits = 0;
93  std::vector<VertexHit> vhits;
94  for (SiPixelRecHitCollection::DataContainer::const_iterator hit = hits->data().begin(), end = hits->data().end();
95  hit != end;
96  ++hit) {
97  if (!hit->isValid())
98  continue;
99  ++nPxlHits;
100  DetId id(hit->geographicalId());
101  if (id.subdetId() != int(PixelSubdetector::PixelBarrel))
102  continue;
103  const PixelGeomDetUnit *pgdu = static_cast<const PixelGeomDetUnit *>(tgeo->idToDet(id));
104  const PixelTopology *pixTopo = &(pgdu->specificTopology());
105  std::vector<SiPixelCluster::Pixel> pixels(hit->cluster()->pixels());
106  bool pixelOnEdge = false;
107  for (std::vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin(); pixel != pixels.end(); ++pixel) {
108  int pixelX = pixel->x;
109  int pixelY = pixel->y;
110  if (pixTopo->isItEdgePixelInX(pixelX) || pixTopo->isItEdgePixelInY(pixelY)) {
111  pixelOnEdge = true;
112  break;
113  }
114  }
115  if (pixelOnEdge)
116  continue;
117 
118  LocalPoint lpos = LocalPoint(hit->localPosition().x(), hit->localPosition().y(), hit->localPosition().z());
119  GlobalPoint gpos = pgdu->toGlobal(lpos);
120  VertexHit vh;
121  vh.z = gpos.z();
122  vh.r = gpos.perp();
123  vh.w = hit->cluster()->sizeY();
124  vhits.push_back(vh);
125  }
126 
127  creco->setNValidPixelHits(nPxlHits);
128 
129  // append cluster compatibility for each z-position
130  for (double z0 = minZ_; z0 <= maxZ_; z0 += zStep_) {
132  creco->append(c.z0, c.nHit, c.chi);
133  }
134  }
135  iEvent.put(std::move(creco));
136 }
137 
139  const std::vector<VertexHit> &hits, double z0) const {
140  // Calculate number of hits contained in v-shaped window in cluster y-width vs. z-position.
141  int n = 0;
142  double chi = 0.;
143 
144  for (std::vector<VertexHit>::const_iterator hit = hits.begin(); hit != hits.end(); hit++) {
145  // the calculation of the predicted cluster width p was
146  // marked 'FIXME' in the HLTPixelClusterShapeFilter. It should
147  // be revisited but is retained as it was for compatibility with the
148  // older filter.
149  double p = 2 * fabs(hit->z - z0) / hit->r + 0.5;
150  if (fabs(p - hit->w) <= 1.) {
151  chi += fabs(p - hit->w);
152  n++;
153  }
154  }
156  output.z0 = z0;
157  output.nHit = n;
158  output.chi = chi;
159  return output;
160 }
161 
162 //define this as a plug-in
#define LogDebug(id)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
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
virtual bool isItEdgePixelInY(int iybin) const =0
Definition: config.py:1
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ClusterCompatibilityProducer(const edm::ParameterSet &)
ContainedHits getContainedHits(const std::vector< VertexHit > &hits, double z0) const
data_type const * data(size_t cell) const
#define end
Definition: vmac.h:39
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.
edm::EDGetTokenT< SiPixelRecHitCollection > inputToken_
HLT enums.
T get() const
Definition: EventSetup.h:73
const TrackerGeomDet * idToDet(DetId) const override
void produce(edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511