CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
OccupancyPlots Class Reference

#include <myTKAnalyses/DigiInvestigator/src/OccupancyPlots.cc>

Inheritance diagram for OccupancyPlots:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 OccupancyPlots (const edm::ParameterSet &)
 
 ~OccupancyPlots () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void beginRun (const edm::Run &, const edm::EventSetup &) override
 
void endJob () override
 
void endRun (const edm::Run &, const edm::EventSetup &) override
 

Private Attributes

TProfile ** m_avemultiplicity
 
TProfile ** m_aveoccupancy
 
TProfile ** m_averadius
 
TProfile ** m_avex
 
TProfile ** m_avey
 
TProfile ** m_avez
 
edm::FileInPath m_fp
 
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_multiplicityMapTokens
 
TH1F ** m_nchannels_ideal
 
TH1F ** m_nchannels_real
 
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_occupancyMapTokens
 
RunHistogramManager m_rhm
 
std::map< unsigned int, DetIdSelectorm_wantedsubdets
 
TProfile ** m_xavedr
 
TProfile ** m_xavedrphi
 
TProfile ** m_xavedz
 
TProfile ** m_yavedr
 
TProfile ** m_yavedrphi
 
TProfile ** m_yavedz
 
TProfile ** m_zavedr
 
TProfile ** m_zavedrphi
 
TProfile ** m_zavedz
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 67 of file OccupancyPlots.cc.

Constructor & Destructor Documentation

OccupancyPlots::OccupancyPlots ( const edm::ParameterSet iConfig)
explicit

Definition at line 123 of file OccupancyPlots.cc.

References edm::EDConsumerBase::consumesCollector(), m_avemultiplicity, m_aveoccupancy, m_averadius, m_avex, m_avey, m_avez, m_fp, m_nchannels_ideal, m_nchannels_real, m_occupancyMapTokens, m_rhm, m_wantedsubdets, m_xavedr, m_xavedrphi, m_xavedz, m_yavedr, m_yavedrphi, m_yavedz, m_zavedr, m_zavedrphi, m_zavedz, RunHistogramManager::makeTH1F(), RunHistogramManager::makeTProfile(), GlobalPosition_Frontier_DevDB_cff::tag, and edm::vector_transform().

123  :
124  m_multiplicityMapTokens(edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag> >("multiplicityMaps"), [this](edm::InputTag const & tag){return consumes<std::map<unsigned int, int> >(tag);})),
125  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);})),
126  m_fp(iConfig.getUntrackedParameter<edm::FileInPath>("file",edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt"))),
128 {
129  //now do what ever initialization is needed
130 
131  m_avemultiplicity = m_rhm.makeTProfile("avemult","Average Multiplicty",6000,0.5,6000.5);
132  m_aveoccupancy = m_rhm.makeTProfile("aveoccu","Average Occupancy",6000,0.5,6000.5);
133 
134  m_nchannels_ideal = m_rhm.makeTH1F("nchannels_ideal","Number of channels (ideal)",6000,0.5,6000.5);
135  m_nchannels_real = m_rhm.makeTH1F("nchannels_real","Number of channels (real)",6000,0.5,6000.5);
136 
137  m_averadius = m_rhm.makeTProfile("averadius","Average Module Radius",6000,0.5,6000.5);
138  m_avez = m_rhm.makeTProfile("avez","Average Module z coordinate",6000,0.5,6000.5);
139  m_avex = m_rhm.makeTProfile("avex","Average Module x coordinate",6000,0.5,6000.5);
140  m_avey = m_rhm.makeTProfile("avey","Average Module y coordinate",6000,0.5,6000.5);
141 
142  m_zavedr = m_rhm.makeTProfile("zavedr","Average z unit vector dr",6000,0.5,6000.5);
143  m_zavedz = m_rhm.makeTProfile("zavedz","Average z unit vector dz",6000,0.5,6000.5);
144  m_zavedrphi = m_rhm.makeTProfile("zavedrphi","Average z unit vector drphi",6000,0.5,6000.5);
145  m_xavedr = m_rhm.makeTProfile("xavedr","Average x unit vector dr",6000,0.5,6000.5);
146  m_xavedz = m_rhm.makeTProfile("xavedz","Average x unit vctor dz",6000,0.5,6000.5);
147  m_xavedrphi = m_rhm.makeTProfile("xavedrphi","Average Module x unit vector drphi",6000,0.5,6000.5);
148  m_yavedr = m_rhm.makeTProfile("yavedr","Average y unit vector dr",6000,0.5,6000.5);
149  m_yavedz = m_rhm.makeTProfile("yavedz","Average y unit vector dz",6000,0.5,6000.5);
150  m_yavedrphi = m_rhm.makeTProfile("yavedrphi","Average y unit vector drphi",6000,0.5,6000.5);
151 
152  std::vector<edm::ParameterSet> wantedsubdets_ps = iConfig.getParameter<std::vector<edm::ParameterSet> >("wantedSubDets");
153 
154  for(std::vector<edm::ParameterSet>::const_iterator wsdps = wantedsubdets_ps.begin();wsdps!=wantedsubdets_ps.end();++wsdps) {
155 
156  unsigned int detsel = wsdps->getParameter<unsigned int>("detSelection");
157  std::vector<std::string> selstr = wsdps->getUntrackedParameter<std::vector<std::string> >("selection");
158  m_wantedsubdets[detsel]=DetIdSelector(selstr);
159 
160  }
161 
162 
163 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
TProfile ** m_zavedrphi
TProfile ** m_xavedz
TH1F ** m_nchannels_ideal
TProfile ** m_zavedr
RunHistogramManager m_rhm
TH1F ** m_nchannels_real
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_multiplicityMapTokens
std::map< unsigned int, DetIdSelector > m_wantedsubdets
TH1F ** makeTH1F(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
TProfile ** m_avex
TProfile ** m_aveoccupancy
TProfile ** m_averadius
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
TProfile ** m_avez
TProfile ** makeTProfile(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
TProfile ** m_xavedr
TProfile ** m_avey
TProfile ** m_yavedr
TProfile ** m_yavedz
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
TProfile ** m_avemultiplicity
TProfile ** m_zavedz
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_occupancyMapTokens
TProfile ** m_yavedrphi
edm::FileInPath m_fp
TProfile ** m_xavedrphi
OccupancyPlots::~OccupancyPlots ( )
override

Definition at line 166 of file OccupancyPlots.cc.

167 {
168 
169  // do anything here that needs to be done at desctruction time
170  // (e.g. close files, deallocate resources etc.)
171 
172 }

Member Function Documentation

void OccupancyPlots::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 181 of file OccupancyPlots.cc.

References edm::Event::getByToken(), m_avemultiplicity, m_aveoccupancy, m_multiplicityMapTokens, and m_occupancyMapTokens.

182 {
183  using namespace edm;
184 
185 
186  for(std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > >::const_iterator mapToken = m_multiplicityMapTokens.begin();mapToken!=m_multiplicityMapTokens.end();++mapToken) {
187 
189  iEvent.getByToken(*mapToken,mults);
190 
191  for(std::map<unsigned int,int>::const_iterator mult=mults->begin();mult!=mults->end();mult++) {
192  if(m_avemultiplicity && *m_avemultiplicity) (*m_avemultiplicity)->Fill(mult->first,mult->second);
193  }
194  }
195 
196 
197 
198  for(std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > >::const_iterator mapToken = m_occupancyMapTokens.begin();mapToken!=m_occupancyMapTokens.end();++mapToken) {
199 
201  iEvent.getByToken(*mapToken,occus);
202 
203  for(std::map<unsigned int,int>::const_iterator occu=occus->begin();occu!=occus->end();occu++) {
204  if(m_aveoccupancy && *m_aveoccupancy) (*m_aveoccupancy)->Fill(occu->first,occu->second);
205  }
206  }
207 
208 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_multiplicityMapTokens
TProfile ** m_aveoccupancy
TProfile ** m_avemultiplicity
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_occupancyMapTokens
HLT enums.
void OccupancyPlots::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 213 of file OccupancyPlots.cc.

214 {
215 
216 }
void OccupancyPlots::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 219 of file OccupancyPlots.cc.

References RunHistogramManager::beginRun(), and m_rhm.

219  {
220 
221  m_rhm.beginRun(iRun);
222 
223 }
RunHistogramManager m_rhm
void beginRun(const edm::Run &iRun)
void OccupancyPlots::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 374 of file OccupancyPlots.cc.

References DEFINE_FWK_MODULE.

374  {
375 }
void OccupancyPlots::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 226 of file OccupancyPlots.cc.

References TrackerGeometry::detUnitIds(), PVValHelper::dx, PVValHelper::dy, PVValHelper::dz, edm::FileInPath::fullPath(), edm::EventSetup::get(), SiPixelDetInfoFileReader::getAllDetIds(), SiStripDetInfoFileReader::getAllDetIds(), SiStripDetInfoFileReader::getNumberOfApvsAndStripLength(), TrackerGeometry::idToDet(), SiPixelQuality::IsModuleBad(), SiStripQuality::IsStripBad(), TrackerGeometry::isThere(), LogDebug, LogTrace, m_averadius, m_avex, m_avey, m_avez, m_fp, m_nchannels_ideal, m_nchannels_real, m_wantedsubdets, m_xavedr, m_xavedrphi, m_xavedz, m_yavedr, m_yavedrphi, m_yavedz, m_zavedr, m_zavedrphi, m_zavedz, Utilities::operator, GeomDetEnumerators::P1PXB, GeomDetEnumerators::P1PXEC, GeomDetEnumerators::P2OTB, GeomDetEnumerators::P2OTEC, GeomDetEnumerators::P2PXB, GeomDetEnumerators::P2PXEC, PV3DBase< T, PVType, FrameType >::perp(), GeomDetEnumerators::PixelBarrel, GeomDetEnumerators::PixelEndcap, position, jets_cff::quality, matplotRender::reader, triggerObjects_cff::sel, digitizers_cfi::strip, GeomDetEnumerators::TEC, GeomDetEnumerators::TIB, GeomDetEnumerators::TID, GeomDetEnumerators::TOB, GeomDet::toGlobal(), DetId::Tracker, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

226  {
227 
228 
229  // edm::ESHandle<GlobalTrackingGeometry> trkgeo;
230  // iSetup.get<GlobalTrackingGeometryRecord>().get("",trkgeo);
232  iSetup.get<TrackerDigiGeometryRecord>().get("",trkgeo);
233 
234  // Test new TrackerGeometry features
235  LogDebug("IsThereTest") << "Test of TrackerGeometry::isThere";
236  LogTrace("IsThereTest") << " is there PixelBarrel: " << trkgeo->isThere(GeomDetEnumerators::PixelBarrel);
237  LogTrace("IsThereTest") << " is there PixelEndcap: " << trkgeo->isThere(GeomDetEnumerators::PixelEndcap);
238  LogTrace("IsThereTest") << " is there P1PXB: " << trkgeo->isThere(GeomDetEnumerators::P1PXB);
239  LogTrace("IsThereTest") << " is there P1PXEC: " << trkgeo->isThere(GeomDetEnumerators::P1PXEC);
240  LogTrace("IsThereTest") << " is there P2PXB: " << trkgeo->isThere(GeomDetEnumerators::P2PXB);
241  LogTrace("IsThereTest") << " is there P2PXEC: " << trkgeo->isThere(GeomDetEnumerators::P2PXEC);
242  LogTrace("IsThereTest") << " is there TIB: " << trkgeo->isThere(GeomDetEnumerators::TIB);
243  LogTrace("IsThereTest") << " is there TID: " << trkgeo->isThere(GeomDetEnumerators::TID);
244  LogTrace("IsThereTest") << " is there TOB: " << trkgeo->isThere(GeomDetEnumerators::TOB);
245  LogTrace("IsThereTest") << " is there TEC: " << trkgeo->isThere(GeomDetEnumerators::TEC);
246  LogTrace("IsThereTest") << " is there P2OTB: " << trkgeo->isThere(GeomDetEnumerators::P2OTB);
247  LogTrace("IsThereTest") << " is there P2OTEC: " << trkgeo->isThere(GeomDetEnumerators::P2OTEC);
248 
249 
250  const Local2DPoint center(0.,0.);
251  const Local3DPoint locz(0.,0.,1.);
252  const Local3DPoint locx(1.,0.,0.);
253  const Local3DPoint locy(0.,1.,0.);
254  const GlobalPoint origin(0.,0.,0.);
255 
256  TrackingGeometry::DetIdContainer detunits = trkgeo->detUnitIds();
257 
258  for(TrackingGeometry::DetIdContainer::const_iterator det = detunits.begin(); det!=detunits.end(); ++det) {
259 
260  if(det->det()!=DetId::Tracker) continue;
261 
262  edm::LogInfo("DetIdFromGeometry") << det->rawId();
263 
264  GlobalPoint position = trkgeo->idToDet(*det)->toGlobal(center);
265  GlobalPoint zpos = trkgeo->idToDet(*det)->toGlobal(locz);
266  GlobalPoint xpos = trkgeo->idToDet(*det)->toGlobal(locx);
267  GlobalPoint ypos = trkgeo->idToDet(*det)->toGlobal(locy);
268  GlobalVector posvect = position - origin;
269  GlobalVector dz = zpos - position;
270  GlobalVector dx = xpos - position;
271  GlobalVector dy = ypos - position;
272 
273  double dzdr = posvect.perp()>0 ? (dz.x()*posvect.x()+dz.y()*posvect.y())/posvect.perp() : 0. ;
274  double dxdr = posvect.perp()>0 ? (dx.x()*posvect.x()+dx.y()*posvect.y())/posvect.perp() : 0. ;
275  double dydr = posvect.perp()>0 ? (dy.x()*posvect.x()+dy.y()*posvect.y())/posvect.perp() : 0. ;
276 
277  double dzdrphi = posvect.perp()>0 ? (dz.y()*posvect.x()-dz.x()*posvect.y())/posvect.perp() : 0. ;
278  double dxdrphi = posvect.perp()>0 ? (dx.y()*posvect.x()-dx.x()*posvect.y())/posvect.perp() : 0. ;
279  double dydrphi = posvect.perp()>0 ? (dy.y()*posvect.x()-dy.x()*posvect.y())/posvect.perp() : 0. ;
280 
281  for(std::map<unsigned int,DetIdSelector>::const_iterator sel=m_wantedsubdets.begin();sel!=m_wantedsubdets.end();++sel) {
282 
283  if(sel->second.isSelected(*det)) {
284  edm::LogInfo("SelectedDetId") << sel->first;
285  // average positions
286  if(m_averadius && *m_averadius) (*m_averadius)->Fill(sel->first,position.perp());
287  if(m_avez && *m_avez) (*m_avez)->Fill(sel->first,position.z());
288  if(m_avex && *m_avex) (*m_avex)->Fill(sel->first,position.x());
289  if(m_avey && *m_avey) (*m_avey)->Fill(sel->first,position.y());
290  if(m_zavedr && *m_zavedr) (*m_zavedr)->Fill(sel->first,dzdr);
291  if(m_zavedz && *m_zavedz) (*m_zavedz)->Fill(sel->first,dz.z());
292  if(m_zavedrphi && *m_zavedrphi) (*m_zavedrphi)->Fill(sel->first,dzdrphi);
293  if(m_xavedr && *m_xavedr) (*m_xavedr)->Fill(sel->first,dxdr);
294  if(m_xavedz && *m_xavedz) (*m_xavedz)->Fill(sel->first,dx.z());
295  if(m_xavedrphi && *m_xavedrphi) (*m_xavedrphi)->Fill(sel->first,dxdrphi);
296  if(m_yavedr && *m_yavedr) (*m_yavedr)->Fill(sel->first,dydr);
297  if(m_yavedz && *m_yavedz) (*m_yavedz)->Fill(sel->first,dy.z());
298  if(m_yavedrphi && *m_yavedrphi) (*m_yavedrphi)->Fill(sel->first,dydrphi);
299  }
300  }
301  }
302 
303  // counting the number of channels per module subset
304 
305  // the histograms have to be reset to avoid double counting if endRun is called more than once
306 
307  if(m_nchannels_ideal && *m_nchannels_ideal) (*m_nchannels_ideal)->Reset();
308  if(m_nchannels_real && *m_nchannels_real) (*m_nchannels_real)->Reset();
309 
311  iSetup.get<SiStripQualityRcd>().get("",quality);
312 
313 
315 
316  const std::vector<uint32_t>& detids = reader->getAllDetIds();
317 
318  for(std::vector<uint32_t>::const_iterator detid=detids.begin();detid!=detids.end();++detid) {
319 
320  int nchannideal = reader->getNumberOfApvsAndStripLength(*detid).first*128;
321  // int nchannreal = reader->getNumberOfApvsAndStripLength(*detid).first*128;
322  int nchannreal = 0;
323  for(int strip = 0; strip < nchannideal; ++strip) {
324  if(!quality->IsStripBad(*detid,strip)) ++nchannreal;
325  }
326 
327 
328  for(std::map<unsigned int,DetIdSelector>::const_iterator sel=m_wantedsubdets.begin();sel!=m_wantedsubdets.end();++sel) {
329 
330  if(sel->second.isSelected(*detid)) {
331  if(m_nchannels_ideal && *m_nchannels_ideal) (*m_nchannels_ideal)->Fill(sel->first,nchannideal);
332  if(m_nchannels_real && *m_nchannels_real) (*m_nchannels_real)->Fill(sel->first,nchannreal);
333  }
334  }
335 
336  }
337 
338 
340  iSetup.get<SiPixelQualityRcd>().get("",pxlquality);
341 
342 
344 
345  const std::vector<uint32_t>& pxldetids = pxlreader.getAllDetIds();
346 
347  for(std::vector<uint32_t>::const_iterator detid=pxldetids.begin();detid!=pxldetids.end();++detid) {
348 
349  int nchannideal = pxlreader.getDetUnitDimensions(*detid).first*pxlreader.getDetUnitDimensions(*detid).second;
350  int nchannreal = 0;
351  if(!pxlquality->IsModuleBad(*detid)) {
352  nchannreal = pxlreader.getDetUnitDimensions(*detid).first*pxlreader.getDetUnitDimensions(*detid).second;
353  }
354  /*
355  int nchannreal = 0;
356  for(int strip = 0; strip < nchannideal; ++strip) {
357  if(!quality->IsStripBad(*detid,strip)) ++nchannreal;
358  }
359  */
360 
361  for(std::map<unsigned int,DetIdSelector>::const_iterator sel=m_wantedsubdets.begin();sel!=m_wantedsubdets.end();++sel) {
362 
363  if(sel->second.isSelected(*detid)) {
364  if(m_nchannels_ideal && *m_nchannels_ideal) (*m_nchannels_ideal)->Fill(sel->first,nchannideal);
365  if(m_nchannels_real && *m_nchannels_real) (*m_nchannels_real)->Fill(sel->first,nchannreal);
366  }
367  }
368 
369  }
370 
371 }
#define LogDebug(id)
const DetIdContainer & detUnitIds() const override
Returm a vector of all GeomDetUnit DetIds.
TProfile ** m_zavedrphi
TProfile ** m_xavedz
T perp() const
Definition: PV3DBase.h:72
const std::vector< uint32_t > & getAllDetIds() const
TH1F ** m_nchannels_ideal
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
T y() const
Definition: PV3DBase.h:63
TProfile ** m_zavedr
TH1F ** m_nchannels_real
std::map< unsigned int, DetIdSelector > m_wantedsubdets
bool isThere(GeomDetEnumerators::SubDetector subdet) const
TProfile ** m_avex
TProfile ** m_averadius
T z() const
Definition: PV3DBase.h:64
TProfile ** m_avez
TProfile ** m_xavedr
TProfile ** m_avey
bool IsModuleBad(const uint32_t &detid) const
#define LogTrace(id)
TProfile ** m_yavedr
TProfile ** m_yavedz
TProfile ** m_zavedz
std::vector< DetId > DetIdContainer
TProfile ** m_yavedrphi
edm::FileInPath m_fp
static int position[264][3]
Definition: ReadPGInfo.cc:509
T get() const
Definition: EventSetup.h:71
const TrackerGeomDet * idToDet(DetId) const override
std::string fullPath() const
Definition: FileInPath.cc:163
T x() const
Definition: PV3DBase.h:62
TProfile ** m_xavedrphi

Member Data Documentation

TProfile** OccupancyPlots::m_avemultiplicity
private

Definition at line 89 of file OccupancyPlots.cc.

Referenced by analyze(), and OccupancyPlots().

TProfile** OccupancyPlots::m_aveoccupancy
private

Definition at line 90 of file OccupancyPlots.cc.

Referenced by analyze(), and OccupancyPlots().

TProfile** OccupancyPlots::m_averadius
private

Definition at line 95 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_avex
private

Definition at line 97 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_avey
private

Definition at line 98 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_avez
private

Definition at line 96 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

edm::FileInPath OccupancyPlots::m_fp
private

Definition at line 84 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > > OccupancyPlots::m_multiplicityMapTokens
private

Definition at line 82 of file OccupancyPlots.cc.

Referenced by analyze().

TH1F** OccupancyPlots::m_nchannels_ideal
private

Definition at line 92 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TH1F** OccupancyPlots::m_nchannels_real
private

Definition at line 93 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

std::vector<edm::EDGetTokenT<std::map<unsigned int, int> > > OccupancyPlots::m_occupancyMapTokens
private

Definition at line 83 of file OccupancyPlots.cc.

Referenced by analyze(), and OccupancyPlots().

RunHistogramManager OccupancyPlots::m_rhm
private

Definition at line 86 of file OccupancyPlots.cc.

Referenced by beginRun(), and OccupancyPlots().

std::map<unsigned int,DetIdSelector> OccupancyPlots::m_wantedsubdets
private

Definition at line 87 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_xavedr
private

Definition at line 105 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_xavedrphi
private

Definition at line 107 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_xavedz
private

Definition at line 106 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_yavedr
private

Definition at line 102 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_yavedrphi
private

Definition at line 104 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_yavedz
private

Definition at line 103 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_zavedr
private

Definition at line 99 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_zavedrphi
private

Definition at line 101 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().

TProfile** OccupancyPlots::m_zavedz
private

Definition at line 100 of file OccupancyPlots.cc.

Referenced by endRun(), and OccupancyPlots().