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 
11 
19 
20 
22  SiPixelPhase1Base(iConfig)
23 {
24  pixelSrcToken_ = consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelSrc"));
25 
26  stripSrcToken_ = consumes<edmNew::DetSetVector<SiStripCluster>>(iConfig.getParameter<edm::InputTag>("stripSrc"));
27 }
28 
30  if( !checktrigger(iEvent,iSetup,DCS) ) return;
31 
33  iEvent.getByToken(pixelSrcToken_, inputPixel);
34  if (!inputPixel.isValid()) return;
35 
37  iEvent.getByToken(stripSrcToken_, inputStrip);
38  if (inputStrip.isValid())
39  {
40  if (!inputStrip.product()->data().empty())
41  {
42  histo[PIXEL_TO_STRIP_RATIO].fill((double)inputPixel.product()->data().size() / (double) inputStrip.product()->data().size(), DetId(0), &iEvent);
43  }
44  }
45 
46  bool hasClusters=false;
47 
49  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
50  assert(tracker.isValid());
51 
53  for (it = inputPixel->begin(); it != inputPixel->end(); ++it) {
54  auto id = DetId(it->detId());
55 
56  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*> ( tracker->idToDet(id) );
57  const PixelTopology& topol = theGeomDet->specificTopology();
58 
59  for(SiPixelCluster const& cluster : *it) {
60  int row = cluster.x()-0.5, col = cluster.y()-0.5;
61  histo[READOUT_CHARGE].fill(double(cluster.charge()), id, &iEvent, col, row);
62  histo[CHARGE].fill(double(cluster.charge()), id, &iEvent, col, row);
63  histo[SIZE ].fill(double(cluster.size() ), id, &iEvent, col, row);
64  histo[SIZEX ].fill(double(cluster.sizeX() ), id, &iEvent, col, row);
65  histo[SIZEY ].fill(double(cluster.sizeY() ), id, &iEvent, col, row);
66  histo[NCLUSTERS].fill(id, &iEvent, col, row);
67  histo[NCLUSTERSINCLUSIVE].fill(id, &iEvent);
68  hasClusters=true;
69  if (cluster.size()>1){
70  histo[READOUT_NCLUSTERS].fill(id, &iEvent);
71  }
72 
73  LocalPoint clustlp = topol.localPosition(MeasurementPoint(cluster.x(), cluster.y()));
74  GlobalPoint clustgp = theGeomDet->surface().toGlobal(clustlp);
75  histo[POSITION_B ].fill(clustgp.z(), clustgp.phi(), id, &iEvent);
76  histo[POSITION_F ].fill(clustgp.x(), clustgp.y(), id, &iEvent);
77  histo[POSITION_XZ].fill(clustgp.x(), clustgp.z(), id, &iEvent);
78  histo[POSITION_YZ].fill(clustgp.y(), clustgp.z(), id, &iEvent);
79  histo[SIZE_VS_ETA].fill(clustgp.eta(), cluster.sizeY(), id, &iEvent);
80 
81  }
82  }
83 
84 
85  if (hasClusters) histo[EVENTRATE].fill(DetId(0), &iEvent);
86 
87  histo[NCLUSTERS].executePerEventHarvesting(&iEvent);
88  histo[READOUT_NCLUSTERS].executePerEventHarvesting(&iEvent);
89  histo[NCLUSTERSINCLUSIVE].executePerEventHarvesting(&iEvent);
90 
91 }
92 
T getParameter(std::string const &) const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
SiPixelPhase1Clusters(const edm::ParameterSet &conf)
bool checktrigger(const edm::Event &iEvent, const edm::EventSetup &iSetup, const unsigned trgidx) const
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
int iEvent
Definition: GenABIO.cc:230
T z() const
Definition: PV3DBase.h:64
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripSrcToken_
bool isValid() const
Definition: HandleBase.h:74
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:58
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
std::vector< HistogramManager > histo
T eta() const
Definition: PV3DBase.h:76
Pixel cluster – collection of neighboring pixels above threshold.
col
Definition: cuy.py:1008
const TrackerGeomDet * idToDet(DetId) const override
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelSrcToken_
void analyze(const edm::Event &, const edm::EventSetup &) override
bool isValid() const
Definition: ESHandle.h:47
T x() const
Definition: PV3DBase.h:62