CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex >
const & 
esGetTokenRecordIndicesVector (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::array< std::vector< ModuleDescription const * > *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, 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
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
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
 
edm::ESGetToken
< SiPixelQuality,
SiPixelQualityRcd
m_pixelQualityToken
 
RunHistogramManager m_rhm
 
edm::ESGetToken
< SiStripQuality,
SiStripQualityRcd
m_stripQualityToken
 
edm::ESGetToken
< TrackerGeometry,
TrackerDigiGeometryRecord
m_tkGeomToken
 
std::map< unsigned int,
DetIdSelector
m_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 wantsInputProcessBlocks ()
 
static bool wantsProcessBlocks ()
 
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)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< B > consumes (edm::InputTag tag) noexcept
 
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<Transition Tr = Transition::Event>
constexpr auto esConsumes () noexcept
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag) noexcept
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

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

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

Definition at line 66 of file OccupancyPlots.cc.

Constructor & Destructor Documentation

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

Definition at line 122 of file OccupancyPlots.cc.

References GlobalPosition_Frontier_DevDB_cff::tag.

124  iConfig.getParameter<std::vector<edm::InputTag> >("multiplicityMaps"),
125  [this](edm::InputTag const& tag) { return consumes<std::map<unsigned int, int> >(tag); })),
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
std::vector< edm::EDGetTokenT< std::map< unsigned int, int > > > m_multiplicityMapTokens
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
OccupancyPlots::~OccupancyPlots ( )
override

Definition at line 170 of file OccupancyPlots.cc.

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

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 180 of file OccupancyPlots.cc.

References edm::Event::getByToken(), m_avemultiplicity, m_aveoccupancy, m_multiplicityMapTokens, m_occupancyMapTokens, VarParsing::mult, and trackerHitRTTI::vector.

180  {
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 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
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
void OccupancyPlots::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 211 of file OccupancyPlots.cc.

211 {}
void OccupancyPlots::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 213 of file OccupancyPlots.cc.

References RunHistogramManager::beginRun(), and m_rhm.

213 { m_rhm.beginRun(iRun); }
RunHistogramManager m_rhm
void beginRun(const edm::Run &iRun)
void OccupancyPlots::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 369 of file OccupancyPlots.cc.

369 {}
void OccupancyPlots::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 215 of file OccupancyPlots.cc.

References PVValHelper::dx, PVValHelper::dy, PVValHelper::dz, edm::FileInPath::fullPath(), GeomDet::geographicalId(), SiPixelDetInfoFileReader::getAllDetIds(), edm::EventSetup::getData(), if(), LogDebug, LogTrace, m_averadius, m_avex, m_avey, m_avez, m_fp, m_nchannels_ideal, m_nchannels_real, m_pixelQualityToken, m_stripQualityToken, m_tkGeomToken, m_wantedsubdets, m_xavedr, m_xavedrphi, m_xavedz, m_yavedr, m_yavedrphi, m_yavedz, m_zavedr, m_zavedrphi, m_zavedz, StripTopology::nstrips(), GeomDetEnumerators::P1PXB, GeomDetEnumerators::P1PXEC, GeomDetEnumerators::P2OTB, GeomDetEnumerators::P2OTEC, GeomDetEnumerators::P2PXB, GeomDetEnumerators::P2PXEC, PV3DBase< T, PVType, FrameType >::perp(), GeomDetEnumerators::PixelBarrel, GeomDetEnumerators::PixelEndcap, position, EgammaValidation_Wenu_cff::sel, StripGeomDetUnit::specificTopology(), digitizers_cfi::strip, GeomDetEnumerators::TEC, GeomDetEnumerators::TIB, GeomDetEnumerators::TID, GeomDetEnumerators::TOB, DetId::Tracker, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

215  {
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 }
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
T y() const
Definition: PV3DBase.h:60
TProfile ** m_zavedr
std::map< unsigned int, DetIdSelector > m_wantedsubdets
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
TH1F ** m_nchannels_real
#define LogTrace(id)
bool getData(T &iHolder) const
Definition: EventSetup.h:128
TProfile ** m_avex
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > m_tkGeomToken
TProfile ** m_averadius
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
Log< level::Info, false > LogInfo
Definition: DetId.h:17
TProfile ** m_yavedz
edm::ESGetToken< SiStripQuality, SiStripQualityRcd > m_stripQualityToken
TProfile ** m_zavedz
TProfile ** m_yavedrphi
std::vector< DetId > DetIdContainer
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
edm::ESGetToken< SiPixelQuality, SiPixelQualityRcd > m_pixelQualityToken
TProfile ** m_xavedrphi
#define LogDebug(id)

Member Data Documentation

TProfile** OccupancyPlots::m_avemultiplicity
private

Definition at line 90 of file OccupancyPlots.cc.

Referenced by analyze().

TProfile** OccupancyPlots::m_aveoccupancy
private

Definition at line 91 of file OccupancyPlots.cc.

Referenced by analyze().

TProfile** OccupancyPlots::m_averadius
private

Definition at line 96 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_avex
private

Definition at line 98 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_avey
private

Definition at line 99 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_avez
private

Definition at line 97 of file OccupancyPlots.cc.

Referenced by endRun().

edm::FileInPath OccupancyPlots::m_fp
private

Definition at line 82 of file OccupancyPlots.cc.

Referenced by endRun().

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

Definition at line 80 of file OccupancyPlots.cc.

Referenced by analyze().

TH1F** OccupancyPlots::m_nchannels_ideal
private

Definition at line 93 of file OccupancyPlots.cc.

Referenced by endRun().

TH1F** OccupancyPlots::m_nchannels_real
private

Definition at line 94 of file OccupancyPlots.cc.

Referenced by endRun().

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

Definition at line 81 of file OccupancyPlots.cc.

Referenced by analyze().

edm::ESGetToken<SiPixelQuality, SiPixelQualityRcd> OccupancyPlots::m_pixelQualityToken
private

Definition at line 85 of file OccupancyPlots.cc.

Referenced by endRun().

RunHistogramManager OccupancyPlots::m_rhm
private

Definition at line 87 of file OccupancyPlots.cc.

Referenced by beginRun().

edm::ESGetToken<SiStripQuality, SiStripQualityRcd> OccupancyPlots::m_stripQualityToken
private

Definition at line 84 of file OccupancyPlots.cc.

Referenced by endRun().

edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> OccupancyPlots::m_tkGeomToken
private

Definition at line 83 of file OccupancyPlots.cc.

Referenced by endRun().

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

Definition at line 88 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_xavedr
private

Definition at line 106 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_xavedrphi
private

Definition at line 108 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_xavedz
private

Definition at line 107 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_yavedr
private

Definition at line 103 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_yavedrphi
private

Definition at line 105 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_yavedz
private

Definition at line 104 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_zavedr
private

Definition at line 100 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_zavedrphi
private

Definition at line 102 of file OccupancyPlots.cc.

Referenced by endRun().

TProfile** OccupancyPlots::m_zavedz
private

Definition at line 101 of file OccupancyPlots.cc.

Referenced by endRun().