CMS 3D CMS Logo

SiPixelPhase1Clusters.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1Clusters
4 // Class: SiPixelPhase1Clusters
5 //
6 
7 // Original Author: Marcel Schneider
8 
14 
21 
22 namespace {
23 
24  class SiPixelPhase1Clusters final : public SiPixelPhase1Base {
25  enum {
26  CHARGE,
27  BIGPIXELCHARGE,
28  NOTBIGPIXELCHARGE,
29  SIZE,
30  SIZEX,
31  SIZEY,
32  NCLUSTERS,
33  NCLUSTERSINCLUSIVE,
34  EVENTRATE,
35  POSITION_B,
36  POSITION_F,
37  POSITION_XZ,
38  POSITION_YZ,
39  SIZE_VS_ETA,
40  READOUT_CHARGE,
41  READOUT_NCLUSTERS,
42  PIXEL_TO_STRIP_RATIO
43  };
44 
45  public:
46  explicit SiPixelPhase1Clusters(const edm::ParameterSet& conf);
47  void analyze(const edm::Event&, const edm::EventSetup&) override;
48 
49  private:
53  };
54 
55  SiPixelPhase1Clusters::SiPixelPhase1Clusters(const edm::ParameterSet& iConfig)
56  : SiPixelPhase1Base(iConfig), trackerGeometryToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()) {
57  pixelSrcToken_ = consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelSrc"));
58 
59  stripSrcToken_ = consumes<edmNew::DetSetVector<SiStripCluster>>(iConfig.getParameter<edm::InputTag>("stripSrc"));
60  }
61 
63  if (!checktrigger(iEvent, iSetup, DCS))
64  return;
65 
67  iEvent.getByToken(pixelSrcToken_, inputPixel);
68  if (!inputPixel.isValid())
69  return;
70 
72  iEvent.getByToken(stripSrcToken_, inputStrip);
73  if (inputStrip.isValid()) {
74  if (!inputStrip.product()->data().empty()) {
75  histo[PIXEL_TO_STRIP_RATIO].fill(
76  (double)inputPixel.product()->data().size() / (double)inputStrip.product()->data().size(),
77  DetId(0),
78  &iEvent);
79  }
80  }
81 
82  bool hasClusters = false;
83 
84  const TrackerGeometry& tracker = iSetup.getData(trackerGeometryToken_);
85 
87  for (it = inputPixel->begin(); it != inputPixel->end(); ++it) {
88  auto id = DetId(it->detId());
89 
90  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(tracker.idToDet(id));
91  const PixelTopology& topol = theGeomDet->specificTopology();
92 
93  for (SiPixelCluster const& cluster : *it) {
94  int row = cluster.x() - 0.5, col = cluster.y() - 0.5;
95  const std::vector<SiPixelCluster::Pixel> pixelsVec = cluster.pixels();
96 
97  for (unsigned int i = 0; i < pixelsVec.size(); ++i) {
98  float pixx = pixelsVec[i].x; // index as float=iteger, row index
99  float pixy = pixelsVec[i].y; // same, col index
100 
101  bool bigInX = topol.isItBigPixelInX(int(pixx));
102  bool bigInY = topol.isItBigPixelInY(int(pixy));
103  float pixel_charge = pixelsVec[i].adc;
104 
105  if (bigInX == true || bigInY == true) {
106  histo[BIGPIXELCHARGE].fill(pixel_charge, id, &iEvent, col, row);
107  }
108 
109  else {
110  histo[NOTBIGPIXELCHARGE].fill(pixel_charge, id, &iEvent, col, row);
111  }
112  }
113  histo[READOUT_CHARGE].fill(double(cluster.charge()), id, &iEvent, col, row);
114  histo[CHARGE].fill(double(cluster.charge()), id, &iEvent, col, row);
115  histo[SIZE].fill(double(cluster.size()), id, &iEvent, col, row);
116  histo[SIZEX].fill(double(cluster.sizeX()), id, &iEvent, col, row);
117  histo[SIZEY].fill(double(cluster.sizeY()), id, &iEvent, col, row);
118  histo[NCLUSTERS].fill(id, &iEvent, col, row);
119  histo[NCLUSTERSINCLUSIVE].fill(id, &iEvent);
120  hasClusters = true;
121  if (cluster.size() > 1) {
122  histo[READOUT_NCLUSTERS].fill(id, &iEvent);
123  }
124 
125  LocalPoint clustlp = topol.localPosition(MeasurementPoint(cluster.x(), cluster.y()));
126  GlobalPoint clustgp = theGeomDet->surface().toGlobal(clustlp);
127  histo[POSITION_B].fill(clustgp.z(), clustgp.phi(), id, &iEvent);
128  histo[POSITION_F].fill(clustgp.x(), clustgp.y(), id, &iEvent);
129  histo[POSITION_XZ].fill(clustgp.x(), clustgp.z(), id, &iEvent);
130  histo[POSITION_YZ].fill(clustgp.y(), clustgp.z(), id, &iEvent);
131  histo[SIZE_VS_ETA].fill(clustgp.eta(), cluster.sizeY(), id, &iEvent);
132  }
133  }
134 
135  if (hasClusters)
136  histo[EVENTRATE].fill(DetId(0), &iEvent);
137 
138  histo[NCLUSTERS].executePerEventHarvesting(&iEvent);
139  histo[READOUT_NCLUSTERS].executePerEventHarvesting(&iEvent);
140  histo[NCLUSTERSINCLUSIVE].executePerEventHarvesting(&iEvent);
141  }
142 
143 } // namespace
144 
145 DEFINE_FWK_MODULE(SiPixelPhase1Clusters);
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
T z() const
Definition: PV3DBase.h:61
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
T eta() const
Definition: PV3DBase.h:73
T const * product() const
Definition: Handle.h:70
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
int iEvent
Definition: GenABIO.cc:224
virtual bool isItBigPixelInX(int ixbin) const =0
virtual bool isItBigPixelInY(int iybin) const =0
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
void analyze(edm::Event const &e, edm::EventSetup const &) override=0
bool isValid() const
Definition: HandleBase.h:70
Pixel cluster – collection of neighboring pixels above threshold.
col
Definition: cuy.py:1009