CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
OccupancyPlots.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: OccupancyPlots
4 // Class: OccupancyPlots
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Mon Oct 27 17:37:53 CET 2008
16 // $Id: OccupancyPlots.cc,v 1.3 2013/02/27 19:49:47 wmtan Exp $
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
24 
25 #include <vector>
26 #include <map>
27 #include <limits>
28 
29 #include "TProfile.h"
30 
41 
46 
51 
61 
62 //
63 // class decleration
64 //
65 
66 class OccupancyPlots : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
67 public:
68  explicit OccupancyPlots(const edm::ParameterSet&);
69  ~OccupancyPlots() override;
70 
71 private:
72  void beginJob() override;
73  void analyze(const edm::Event&, const edm::EventSetup&) override;
74  void beginRun(const edm::Run&, const edm::EventSetup&) override;
75  void endRun(const edm::Run&, const edm::EventSetup&) override;
76  void endJob() override;
77 
78  // ----------member data ---------------------------
79 
80  std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > > m_multiplicityMapTokens;
81  std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > > m_occupancyMapTokens;
86 
88  std::map<unsigned int, DetIdSelector> m_wantedsubdets;
89 
90  TProfile** m_avemultiplicity;
91  TProfile** m_aveoccupancy;
92 
95 
96  TProfile** m_averadius;
97  TProfile** m_avez;
98  TProfile** m_avex;
99  TProfile** m_avey;
100  TProfile** m_zavedr;
101  TProfile** m_zavedz;
102  TProfile** m_zavedrphi;
103  TProfile** m_yavedr;
104  TProfile** m_yavedz;
105  TProfile** m_yavedrphi;
106  TProfile** m_xavedr;
107  TProfile** m_xavedz;
108  TProfile** m_xavedrphi;
109 };
110 
111 //
112 // constants, enums and typedefs
113 //
114 
115 //
116 // static data member definitions
117 //
118 
119 //
120 // constructors and destructor
121 //
123  : m_multiplicityMapTokens(edm::vector_transform(
124  iConfig.getParameter<std::vector<edm::InputTag> >("multiplicityMaps"),
125  [this](edm::InputTag const& tag) { return consumes<std::map<unsigned int, int> >(tag); })),
127  iConfig.getParameter<std::vector<edm::InputTag> >("occupancyMaps"),
128  [this](edm::InputTag const& tag) { return consumes<std::map<unsigned int, int> >(tag); })),
129  m_fp(iConfig.getUntrackedParameter<edm::FileInPath>(
130  "file", edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt"))),
131  m_tkGeomToken(esConsumes<edm::Transition::EndRun>()),
132  m_stripQualityToken(esConsumes<edm::Transition::EndRun>()),
133  m_pixelQualityToken(esConsumes<edm::Transition::EndRun>()),
134  m_rhm(consumesCollector()),
136  //now do what ever initialization is needed
137 
138  m_avemultiplicity = m_rhm.makeTProfile("avemult", "Average Multiplicty", 6000, 0.5, 6000.5);
139  m_aveoccupancy = m_rhm.makeTProfile("aveoccu", "Average Occupancy", 6000, 0.5, 6000.5);
140 
141  m_nchannels_ideal = m_rhm.makeTH1F("nchannels_ideal", "Number of channels (ideal)", 6000, 0.5, 6000.5);
142  m_nchannels_real = m_rhm.makeTH1F("nchannels_real", "Number of channels (real)", 6000, 0.5, 6000.5);
143 
144  m_averadius = m_rhm.makeTProfile("averadius", "Average Module Radius", 6000, 0.5, 6000.5);
145  m_avez = m_rhm.makeTProfile("avez", "Average Module z coordinate", 6000, 0.5, 6000.5);
146  m_avex = m_rhm.makeTProfile("avex", "Average Module x coordinate", 6000, 0.5, 6000.5);
147  m_avey = m_rhm.makeTProfile("avey", "Average Module y coordinate", 6000, 0.5, 6000.5);
148 
149  m_zavedr = m_rhm.makeTProfile("zavedr", "Average z unit vector dr", 6000, 0.5, 6000.5);
150  m_zavedz = m_rhm.makeTProfile("zavedz", "Average z unit vector dz", 6000, 0.5, 6000.5);
151  m_zavedrphi = m_rhm.makeTProfile("zavedrphi", "Average z unit vector drphi", 6000, 0.5, 6000.5);
152  m_xavedr = m_rhm.makeTProfile("xavedr", "Average x unit vector dr", 6000, 0.5, 6000.5);
153  m_xavedz = m_rhm.makeTProfile("xavedz", "Average x unit vctor dz", 6000, 0.5, 6000.5);
154  m_xavedrphi = m_rhm.makeTProfile("xavedrphi", "Average Module x unit vector drphi", 6000, 0.5, 6000.5);
155  m_yavedr = m_rhm.makeTProfile("yavedr", "Average y unit vector dr", 6000, 0.5, 6000.5);
156  m_yavedz = m_rhm.makeTProfile("yavedz", "Average y unit vector dz", 6000, 0.5, 6000.5);
157  m_yavedrphi = m_rhm.makeTProfile("yavedrphi", "Average y unit vector drphi", 6000, 0.5, 6000.5);
158 
159  std::vector<edm::ParameterSet> wantedsubdets_ps =
160  iConfig.getParameter<std::vector<edm::ParameterSet> >("wantedSubDets");
161 
162  for (std::vector<edm::ParameterSet>::const_iterator wsdps = wantedsubdets_ps.begin(); wsdps != wantedsubdets_ps.end();
163  ++wsdps) {
164  unsigned int detsel = wsdps->getParameter<unsigned int>("detSelection");
165  std::vector<std::string> selstr = wsdps->getUntrackedParameter<std::vector<std::string> >("selection");
166  m_wantedsubdets[detsel] = DetIdSelector(selstr);
167  }
168 }
169 
171  // do anything here that needs to be done at desctruction time
172  // (e.g. close files, deallocate resources etc.)
173 }
174 
175 //
176 // member functions
177 //
178 
179 // ------------ method called to for each event ------------
181  using namespace edm;
182 
183  for (std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > >::const_iterator mapToken =
184  m_multiplicityMapTokens.begin();
185  mapToken != m_multiplicityMapTokens.end();
186  ++mapToken) {
188  iEvent.getByToken(*mapToken, mults);
189 
190  for (std::map<unsigned int, int>::const_iterator mult = mults->begin(); mult != mults->end(); mult++) {
192  (*m_avemultiplicity)->Fill(mult->first, mult->second);
193  }
194  }
195 
196  for (std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > >::const_iterator mapToken =
197  m_occupancyMapTokens.begin();
198  mapToken != m_occupancyMapTokens.end();
199  ++mapToken) {
201  iEvent.getByToken(*mapToken, occus);
202 
203  for (std::map<unsigned int, int>::const_iterator occu = occus->begin(); occu != occus->end(); occu++) {
205  (*m_aveoccupancy)->Fill(occu->first, occu->second);
206  }
207  }
208 }
209 
210 // ------------ method called once each job just before starting event loop ------------
212 
213 void OccupancyPlots::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) { m_rhm.beginRun(iRun); }
214 
215 void OccupancyPlots::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
216  const auto& trkgeo = iSetup.getData(m_tkGeomToken);
217 
218  // Test new TrackerGeometry features
219  LogDebug("IsThereTest") << "Test of TrackerGeometry::isThere";
220  LogTrace("IsThereTest") << " is there PixelBarrel: " << trkgeo.isThere(GeomDetEnumerators::PixelBarrel);
221  LogTrace("IsThereTest") << " is there PixelEndcap: " << trkgeo.isThere(GeomDetEnumerators::PixelEndcap);
222  LogTrace("IsThereTest") << " is there P1PXB: " << trkgeo.isThere(GeomDetEnumerators::P1PXB);
223  LogTrace("IsThereTest") << " is there P1PXEC: " << trkgeo.isThere(GeomDetEnumerators::P1PXEC);
224  LogTrace("IsThereTest") << " is there P2PXB: " << trkgeo.isThere(GeomDetEnumerators::P2PXB);
225  LogTrace("IsThereTest") << " is there P2PXEC: " << trkgeo.isThere(GeomDetEnumerators::P2PXEC);
226  LogTrace("IsThereTest") << " is there TIB: " << trkgeo.isThere(GeomDetEnumerators::TIB);
227  LogTrace("IsThereTest") << " is there TID: " << trkgeo.isThere(GeomDetEnumerators::TID);
228  LogTrace("IsThereTest") << " is there TOB: " << trkgeo.isThere(GeomDetEnumerators::TOB);
229  LogTrace("IsThereTest") << " is there TEC: " << trkgeo.isThere(GeomDetEnumerators::TEC);
230  LogTrace("IsThereTest") << " is there P2OTB: " << trkgeo.isThere(GeomDetEnumerators::P2OTB);
231  LogTrace("IsThereTest") << " is there P2OTEC: " << trkgeo.isThere(GeomDetEnumerators::P2OTEC);
232 
233  const Local2DPoint center(0., 0.);
234  const Local3DPoint locz(0., 0., 1.);
235  const Local3DPoint locx(1., 0., 0.);
236  const Local3DPoint locy(0., 1., 0.);
237  const GlobalPoint origin(0., 0., 0.);
238 
239  TrackingGeometry::DetIdContainer detunits = trkgeo.detUnitIds();
240 
241  for (TrackingGeometry::DetIdContainer::const_iterator det = detunits.begin(); det != detunits.end(); ++det) {
242  if (det->det() != DetId::Tracker)
243  continue;
244 
245  edm::LogInfo("DetIdFromGeometry") << det->rawId();
246 
247  GlobalPoint position = trkgeo.idToDet(*det)->toGlobal(center);
248  GlobalPoint zpos = trkgeo.idToDet(*det)->toGlobal(locz);
249  GlobalPoint xpos = trkgeo.idToDet(*det)->toGlobal(locx);
250  GlobalPoint ypos = trkgeo.idToDet(*det)->toGlobal(locy);
251  GlobalVector posvect = position - origin;
252  GlobalVector dz = zpos - position;
253  GlobalVector dx = xpos - position;
254  GlobalVector dy = ypos - position;
255 
256  double dzdr = posvect.perp() > 0 ? (dz.x() * posvect.x() + dz.y() * posvect.y()) / posvect.perp() : 0.;
257  double dxdr = posvect.perp() > 0 ? (dx.x() * posvect.x() + dx.y() * posvect.y()) / posvect.perp() : 0.;
258  double dydr = posvect.perp() > 0 ? (dy.x() * posvect.x() + dy.y() * posvect.y()) / posvect.perp() : 0.;
259 
260  double dzdrphi = posvect.perp() > 0 ? (dz.y() * posvect.x() - dz.x() * posvect.y()) / posvect.perp() : 0.;
261  double dxdrphi = posvect.perp() > 0 ? (dx.y() * posvect.x() - dx.x() * posvect.y()) / posvect.perp() : 0.;
262  double dydrphi = posvect.perp() > 0 ? (dy.y() * posvect.x() - dy.x() * posvect.y()) / posvect.perp() : 0.;
263 
264  for (std::map<unsigned int, DetIdSelector>::const_iterator sel = m_wantedsubdets.begin();
265  sel != m_wantedsubdets.end();
266  ++sel) {
267  if (sel->second.isSelected(*det)) {
268  edm::LogInfo("SelectedDetId") << sel->first;
269  // average positions
270  if (m_averadius && *m_averadius)
271  (*m_averadius)->Fill(sel->first, position.perp());
272  if (m_avez && *m_avez)
273  (*m_avez)->Fill(sel->first, position.z());
274  if (m_avex && *m_avex)
275  (*m_avex)->Fill(sel->first, position.x());
276  if (m_avey && *m_avey)
277  (*m_avey)->Fill(sel->first, position.y());
278  if (m_zavedr && *m_zavedr)
279  (*m_zavedr)->Fill(sel->first, dzdr);
280  if (m_zavedz && *m_zavedz)
281  (*m_zavedz)->Fill(sel->first, dz.z());
282  if (m_zavedrphi && *m_zavedrphi)
283  (*m_zavedrphi)->Fill(sel->first, dzdrphi);
284  if (m_xavedr && *m_xavedr)
285  (*m_xavedr)->Fill(sel->first, dxdr);
286  if (m_xavedz && *m_xavedz)
287  (*m_xavedz)->Fill(sel->first, dx.z());
288  if (m_xavedrphi && *m_xavedrphi)
289  (*m_xavedrphi)->Fill(sel->first, dxdrphi);
290  if (m_yavedr && *m_yavedr)
291  (*m_yavedr)->Fill(sel->first, dydr);
292  if (m_yavedz && *m_yavedz)
293  (*m_yavedz)->Fill(sel->first, dy.z());
294  if (m_yavedrphi && *m_yavedrphi)
295  (*m_yavedrphi)->Fill(sel->first, dydrphi);
296  }
297  }
298  }
299 
300  // counting the number of channels per module subset
301 
302  // the histograms have to be reset to avoid double counting if endRun is called more than once
303 
305  (*m_nchannels_ideal)->Reset();
307  (*m_nchannels_real)->Reset();
308 
309  const auto& stripQuality = iSetup.getData(m_stripQualityToken);
310 
311  for (const auto det : trkgeo.detUnits()) {
312  const StripGeomDetUnit* stripDet = dynamic_cast<const StripGeomDetUnit*>(det);
313  if (stripDet != nullptr) {
314  const DetId detid = stripDet->geographicalId();
315 
316  int nchannideal = stripDet->specificTopology().nstrips();
317  // int nchannreal = stripDet->specificTopology().nstrips();
318  int nchannreal = 0;
319  for (int strip = 0; strip < nchannideal; ++strip) {
320  if (!stripQuality.IsStripBad(detid, strip))
321  ++nchannreal;
322  }
323 
324  for (std::map<unsigned int, DetIdSelector>::const_iterator sel = m_wantedsubdets.begin();
325  sel != m_wantedsubdets.end();
326  ++sel) {
327  if (sel->second.isSelected(detid)) {
329  (*m_nchannels_ideal)->Fill(sel->first, nchannideal);
331  (*m_nchannels_real)->Fill(sel->first, nchannreal);
332  }
333  }
334  }
335  }
336 
337  const auto& pxlquality = iSetup.getData(m_pixelQualityToken);
338 
340 
341  const std::vector<uint32_t>& pxldetids = pxlreader.getAllDetIds();
342 
343  for (std::vector<uint32_t>::const_iterator detid = pxldetids.begin(); detid != pxldetids.end(); ++detid) {
344  int nchannideal = pxlreader.getDetUnitDimensions(*detid).first * pxlreader.getDetUnitDimensions(*detid).second;
345  int nchannreal = 0;
346  if (!pxlquality.IsModuleBad(*detid)) {
347  nchannreal = pxlreader.getDetUnitDimensions(*detid).first * pxlreader.getDetUnitDimensions(*detid).second;
348  }
349  /*
350  int nchannreal = 0;
351  for(int strip = 0; strip < nchannideal; ++strip) {
352  if(!stripquality.IsStripBad(*detid,strip)) ++nchannreal;
353  }
354  */
355 
356  for (std::map<unsigned int, DetIdSelector>::const_iterator sel = m_wantedsubdets.begin();
357  sel != m_wantedsubdets.end();
358  ++sel) {
359  if (sel->second.isSelected(*detid)) {
361  (*m_nchannels_ideal)->Fill(sel->first, nchannideal);
363  (*m_nchannels_real)->Fill(sel->first, nchannreal);
364  }
365  }
366  }
367 }
368 // ------------ method called once each job just after ending the event loop ------------
370 //define this as a plug-in
OccupancyPlots(const edm::ParameterSet &)
virtual int nstrips() const =0
TProfile ** m_zavedrphi
TProfile ** m_xavedz
T perp() const
Definition: PV3DBase.h:69
const std::vector< uint32_t > & getAllDetIds() const
TH1F ** m_nchannels_ideal
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
void beginRun(const edm::Run &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T y() const
Definition: PV3DBase.h:60
TProfile ** m_zavedr
std::map< unsigned int, DetIdSelector > m_wantedsubdets
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
void beginJob() override
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
RunHistogramManager m_rhm
TH1F ** m_nchannels_real
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_multiplicityMapTokens
#define LogTrace(id)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
m_fp(iConfig.getUntrackedParameter< edm::FileInPath >("file", edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt")))
int iEvent
Definition: GenABIO.cc:224
TProfile ** m_avex
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > m_tkGeomToken
TProfile ** m_aveoccupancy
TProfile ** m_averadius
m_rhm(consumesCollector())
T z() const
Definition: PV3DBase.h:61
if(conf_.getParameter< bool >("UseStripCablingDB"))
TProfile ** m_avez
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
TProfile ** m_xavedr
TProfile ** m_avey
TProfile ** m_yavedr
m_occupancyMapTokens(edm::vector_transform(iConfig.getParameter< std::vector< edm::InputTag > >("occupancyMaps"), [this](edm::InputTag const &tag){return consumes< std::map< unsigned int, int > >(tag);}))
Log< level::Info, false > LogInfo
Definition: DetId.h:17
TProfile ** m_yavedz
TProfile ** m_avemultiplicity
m_wantedsubdets()
edm::ESGetToken< SiStripQuality, SiStripQualityRcd > m_stripQualityToken
TProfile ** m_zavedz
void beginRun(const edm::Run &iRun)
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_occupancyMapTokens
m_stripQualityToken(esConsumes< edm::Transition::EndRun >())
TProfile ** m_yavedrphi
std::vector< DetId > DetIdContainer
~OccupancyPlots() override
edm::FileInPath m_fp
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::string fullPath() const
Definition: FileInPath.cc:161
T x() const
Definition: PV3DBase.h:59
m_tkGeomToken(esConsumes< edm::Transition::EndRun >())
edm::ESGetToken< SiPixelQuality, SiPixelQualityRcd > m_pixelQualityToken
void endJob() override
void endRun(const edm::Run &, const edm::EventSetup &) override
Definition: Run.h:45
TProfile ** m_xavedrphi
m_pixelQualityToken(esConsumes< edm::Transition::EndRun >())
#define LogDebug(id)
void analyze(const edm::Event &, const edm::EventSetup &) override