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 | Protected Member Functions | Private Member Functions | Private Attributes
DTSegmentAnalysisTask Class Reference

#include <DTSegmentAnalysisTask.h>

Inheritance diagram for DTSegmentAnalysisTask:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &setup) override
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 BeginRun. More...
 
 DTSegmentAnalysisTask (const edm::ParameterSet &pset)
 Constructor. More...
 
 ~DTSegmentAnalysisTask () override
 Destructor. More...
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Protected Member Functions

void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 

Private Member Functions

void bookHistos (DQMStore::IBooker &ibooker, DTChamberId chamberId)
 
void fillHistos (DTChamberId chamberId, int nHits, float chi2)
 

Private Attributes

bool checkNoisyChannels
 
bool detailedAnalysis
 
const DTGeometrydtGeom
 
std::map< DTChamberId,
std::vector< MonitorElement * > > 
histosPerCh
 
bool hltDQMMode
 
edm::ESGetToken< DTGeometry,
MuonGeometryRecord
muonGeomToken_
 
MonitorElementnEventMonitor
 
int nevents
 
int nhitsCut
 
double phiSegmCut
 
edm::EDGetTokenT
< DTRecSegment4DCollection
recHits4DToken_
 
const DTStatusFlagstatusMap
 
edm::ESGetToken< DTStatusFlag,
DTStatusFlagRcd
statusMapToken_
 
std::map< int, MonitorElement * > summaryHistos
 
std::string topHistoFolder
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T...>
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T...>
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr
< DQMEDAnalyzerGlobalCache
initializeGlobalCache (edm::ParameterSet const &)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

DQM Analysis of 4D DT segments, it produces plots about:

Author
G. Cerminara - INFN Torino

Definition at line 42 of file DTSegmentAnalysisTask.h.

Constructor & Destructor Documentation

DTSegmentAnalysisTask::DTSegmentAnalysisTask ( const edm::ParameterSet pset)

Constructor.

Definition at line 33 of file DTSegmentAnalysisTask.cc.

References checkNoisyChannels, detailedAnalysis, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hltDQMMode, HLT_FULL_cff::InputTag, nhitsCut, phiSegmCut, recHits4DToken_, and topHistoFolder.

34  : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()), statusMapToken_(esConsumes()), nevents(0) {
35  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Constructor called!";
36 
37  // switch for detailed analysis
38  detailedAnalysis = pset.getUntrackedParameter<bool>("detailedAnalysis", false);
39  // the name of the 4D rec hits collection
40  recHits4DToken_ = consumes<DTRecSegment4DCollection>(edm::InputTag(pset.getParameter<string>("recHits4DLabel")));
41  // Get the map of noisy channels
42  checkNoisyChannels = pset.getUntrackedParameter<bool>("checkNoisyChannels", false);
43  phiSegmCut = pset.getUntrackedParameter<double>("phiSegmCut", 30.);
44  nhitsCut = pset.getUntrackedParameter<int>("nhitsCut", 12);
45 
46  // top folder for the histograms in DQMStore
47  topHistoFolder = pset.getUntrackedParameter<string>("topHistoFolder", "DT/02-Segments");
48  // hlt DQM mode
49  hltDQMMode = pset.getUntrackedParameter<bool>("hltDQMMode", false);
50 }
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
edm::ESGetToken< DTStatusFlag, DTStatusFlagRcd > statusMapToken_
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< DTRecSegment4DCollection > recHits4DToken_
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
DTSegmentAnalysisTask::~DTSegmentAnalysisTask ( )
override

Destructor.

Definition at line 52 of file DTSegmentAnalysisTask.cc.

52  {
53  //FR moved fron endjob
54  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Destructor called!";
55 }
Log< level::Info, true > LogVerbatim

Member Function Documentation

void DTSegmentAnalysisTask::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
overridevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 105 of file DTSegmentAnalysisTask.cc.

References DTStatusFlag::cellStatus(), checkNoisyChannels, HLT_FULL_cff::distance, edm::EventID::event(), dqm::impl::MonitorElement::Fill(), fillHistos(), edm::EventSetup::getData(), edm::EventBase::id(), edm::HandleBase::isValid(), nEventMonitor, nevents, nHits, nhitsCut, phiSegmCut, Pi, sistrip::SpyUtilities::range(), recHits4DToken_, findQualityFiles::size, DTRecSegment2D::specificRecHits(), statusMap, statusMapToken_, and xdir.

105  {
106  nevents++;
108 
109  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")
110  << "[DTSegmentAnalysisTask] Analyze #Run: " << event.id().run() << " #Event: " << event.id().event();
111  if (!(event.id().event() % 1000))
112  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")
113  << "[DTSegmentAnalysisTask] Analyze #Run: " << event.id().run() << " #Event: " << event.id().event();
114 
115  if (checkNoisyChannels) {
117  }
118 
119  // -- 4D segment analysis -----------------------------------------------------
120 
121  // Get the 4D segment collection from the event
123  event.getByToken(recHits4DToken_, all4DSegments);
124 
125  if (!all4DSegments.isValid())
126  return;
127 
128  // Loop over all chambers containing a segment
130  for (chamberId = all4DSegments->id_begin(); chamberId != all4DSegments->id_end(); ++chamberId) {
131  // Get the range for the corresponding ChamerId
132  DTRecSegment4DCollection::range range = all4DSegments->get(*chamberId);
133 
134  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")
135  << " Chamber: " << *chamberId << " has " << distance(range.first, range.second) << " 4D segments";
136 
137  // Loop over the rechits of this ChamerId
138  for (DTRecSegment4DCollection::const_iterator segment4D = range.first; segment4D != range.second; ++segment4D) {
139  //FOR NOISY CHANNELS////////////////////////////////
140  bool segmNoisy = false;
141  if (checkNoisyChannels) {
142  if ((*segment4D).hasPhi()) {
143  const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
144  vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
145  map<DTSuperLayerId, vector<DTRecHit1D> > hitsBySLMap;
146  for (vector<DTRecHit1D>::const_iterator hit = phiHits.begin(); hit != phiHits.end(); ++hit) {
147  DTWireId wireId = (*hit).wireId();
148 
149  // Check for noisy channels to skip them
150  bool isNoisy = false;
151  bool isFEMasked = false;
152  bool isTDCMasked = false;
153  bool isTrigMask = false;
154  bool isDead = false;
155  bool isNohv = false;
156  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
157  if (isNoisy) {
158  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")
159  << "Wire: " << wireId << " is noisy, skipping!";
160  segmNoisy = true;
161  }
162  }
163  }
164 
165  if ((*segment4D).hasZed()) {
166  const DTSLRecSegment2D* zSeg = (*segment4D).zSegment(); // zSeg lives in the SL RF
167  // Check for noisy channels to skip them
168  vector<DTRecHit1D> zHits = zSeg->specificRecHits();
169  for (vector<DTRecHit1D>::const_iterator hit = zHits.begin(); hit != zHits.end(); ++hit) {
170  DTWireId wireId = (*hit).wireId();
171  bool isNoisy = false;
172  bool isFEMasked = false;
173  bool isTDCMasked = false;
174  bool isTrigMask = false;
175  bool isDead = false;
176  bool isNohv = false;
177  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
178  if (isNoisy) {
179  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")
180  << "Wire: " << wireId << " is noisy, skipping!";
181  segmNoisy = true;
182  }
183  }
184  }
185 
186  } // end of switch on noisy channels
187  if (segmNoisy) {
188  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")
189  << "skipping the segment: it contains noisy cells";
190  continue;
191  }
192  //END FOR NOISY CHANNELS////////////////////////////////
193 
194  int nHits = 0;
195  if ((*segment4D).hasPhi())
196  nHits = (((*segment4D).phiSegment())->specificRecHits()).size();
197  if ((*segment4D).hasZed())
198  nHits = nHits + ((((*segment4D).zSegment())->specificRecHits()).size());
199 
200  double anglePhiSegm(0.);
201  if ((*segment4D).hasPhi()) {
202  double xdir = (*segment4D).phiSegment()->localDirection().x();
203  double zdir = (*segment4D).phiSegment()->localDirection().z();
204 
205  anglePhiSegm = atan(xdir / zdir) * 180. / TMath::Pi();
206  }
207  if (fabs(anglePhiSegm) > phiSegmCut)
208  continue;
209  // If the segment is in Wh+-2/SecX/MB1, get the DT chambers just above and check if there is a segment
210  // to validate the segment present in MB1
211  if (fabs((*chamberId).wheel()) == 2 && (*chamberId).station() == 1) {
212  bool segmOk = false;
213  int mb(2);
214  while (mb < 4) {
215  DTChamberId checkMB((*chamberId).wheel(), mb, (*chamberId).sector());
216  DTRecSegment4DCollection::range ckrange = all4DSegments->get(checkMB);
217 
218  for (DTRecSegment4DCollection::const_iterator cksegment4D = ckrange.first; cksegment4D != ckrange.second;
219  ++cksegment4D) {
220  int nHits = 0;
221  if ((*cksegment4D).hasPhi())
222  nHits = (((*cksegment4D).phiSegment())->specificRecHits()).size();
223  if ((*cksegment4D).hasZed())
224  nHits = nHits + ((((*cksegment4D).zSegment())->specificRecHits()).size());
225 
226  if (nHits >= nhitsCut)
227  segmOk = true;
228  }
229  mb++;
230  }
231 
232  if (!segmOk)
233  continue;
234  }
235  fillHistos(*chamberId, nHits, (*segment4D).chi2() / (*segment4D).degreesOfFreedom());
236  }
237  }
238 
239  // -----------------------------------------------------------------------------
240 }
const double Pi
Log< level::Info, true > LogVerbatim
EventNumber_t event() const
Definition: EventID.h:40
void fillHistos(DTChamberId chamberId, int nHits, float chi2)
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
identifier iterator
Definition: RangeMap.h:130
edm::ESGetToken< DTStatusFlag, DTStatusFlagRcd > statusMapToken_
const uint16_t range(const Frame &aFrame)
void Fill(long long x)
bool getData(T &iHolder) const
Definition: EventSetup.h:128
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
bool isValid() const
Definition: HandleBase.h:70
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
const DTStatusFlag * statusMap
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
int cellStatus(int wheelId, int stationId, int sectorId, int slId, int layerId, int cellId, bool &noiseFlag, bool &feMask, bool &tdcMask, bool &trigMask, bool &deadFlag, bool &nohvFlag) const
get content
Definition: DTStatusFlag.h:90
MonitorElement * nEventMonitor
edm::EventID id() const
Definition: EventBase.h:59
edm::EDGetTokenT< DTRecSegment4DCollection > recHits4DToken_
tuple size
Write out results.
void DTSegmentAnalysisTask::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  context 
)
overrideprotectedvirtual

Implements DQMEDAnalyzer.

Definition at line 62 of file DTSegmentAnalysisTask.cc.

References dqm::implementation::IBooker::book2D(), dqm::implementation::IBooker::bookFloat(), bookHistos(), DTGeometry::chambers(), chambers, dtGeom, hltDQMMode, nEventMonitor, dqm::implementation::NavigatorBase::setCurrentFolder(), summaryHistos, and topHistoFolder.

64  {
65  if (!hltDQMMode) {
66  ibooker.setCurrentFolder("DT/EventInfo/Counters");
67  nEventMonitor = ibooker.bookFloat("nProcessedEventsSegment");
68  }
69 
70  for (int wh = -2; wh <= 2; wh++) {
71  stringstream wheel;
72  wheel << wh;
73  ibooker.setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str());
74  string histoName = "numberOfSegments_W" + wheel.str();
75 
76  summaryHistos[wh] = ibooker.book2D(histoName.c_str(), histoName.c_str(), 12, 1, 13, 4, 1, 5);
77  summaryHistos[wh]->setAxisTitle("Sector", 1);
78  summaryHistos[wh]->setBinLabel(1, "1", 1);
79  summaryHistos[wh]->setBinLabel(2, "2", 1);
80  summaryHistos[wh]->setBinLabel(3, "3", 1);
81  summaryHistos[wh]->setBinLabel(4, "4", 1);
82  summaryHistos[wh]->setBinLabel(5, "5", 1);
83  summaryHistos[wh]->setBinLabel(6, "6", 1);
84  summaryHistos[wh]->setBinLabel(7, "7", 1);
85  summaryHistos[wh]->setBinLabel(8, "8", 1);
86  summaryHistos[wh]->setBinLabel(9, "9", 1);
87  summaryHistos[wh]->setBinLabel(10, "10", 1);
88  summaryHistos[wh]->setBinLabel(11, "11", 1);
89  summaryHistos[wh]->setBinLabel(12, "12", 1);
90  summaryHistos[wh]->setBinLabel(1, "MB1", 2);
91  summaryHistos[wh]->setBinLabel(2, "MB2", 2);
92  summaryHistos[wh]->setBinLabel(3, "MB3", 2);
93  summaryHistos[wh]->setBinLabel(4, "MB4", 2);
94  }
95 
96  // loop over all the DT chambers & book the histos
97  const vector<const DTChamber*>& chambers = dtGeom->chambers();
98  vector<const DTChamber*>::const_iterator ch_it = chambers.begin();
99  vector<const DTChamber*>::const_iterator ch_end = chambers.end();
100  for (; ch_it != ch_end; ++ch_it) {
101  bookHistos(ibooker, (*ch_it)->id());
102  }
103 }
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
MonitorElement * bookFloat(TString const &name, FUNC onbooking=NOOP())
Definition: DQMStore.h:80
const DTGeometry * dtGeom
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
void bookHistos(DQMStore::IBooker &ibooker, DTChamberId chamberId)
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
MonitorElement * nEventMonitor
std::map< int, MonitorElement * > summaryHistos
static char chambers[264][20]
Definition: ReadPGInfo.cc:243
void DTSegmentAnalysisTask::bookHistos ( DQMStore::IBooker ibooker,
DTChamberId  chamberId 
)
private

Definition at line 243 of file DTSegmentAnalysisTask.cc.

References dqm::implementation::IBooker::book1D(), detailedAnalysis, mergeVDriftHistosByStation::histos, histosPerCh, DTChamberId::sector(), dqm::implementation::NavigatorBase::setCurrentFolder(), DTChamberId::station(), relativeConstraints::station, topHistoFolder, and DTChamberId::wheel().

Referenced by bookHistograms().

243  {
244  edm::LogVerbatim("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << " Booking histos for chamber: " << chamberId;
245 
246  // Compose the chamber name
247  stringstream wheel;
248  wheel << chamberId.wheel();
249  stringstream station;
250  station << chamberId.station();
251  stringstream sector;
252  sector << chamberId.sector();
253 
254  string chamberHistoName = "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
255 
256  ibooker.setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" +
257  station.str());
258 
259  // Create the monitor elements
260  vector<MonitorElement*> histos;
261  histos.push_back(ibooker.book1D("h4DSegmNHits" + chamberHistoName, "# of hits per segment", 16, 0.5, 16.5));
262  if (detailedAnalysis) {
263  histos.push_back(ibooker.book1D("h4DChi2" + chamberHistoName, "4D Segment reduced Chi2", 20, 0, 20));
264  }
265  histosPerCh[chamberId] = histos;
266 }
Log< level::Info, true > LogVerbatim
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
int sector() const
Definition: DTChamberId.h:49
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
void DTSegmentAnalysisTask::dqmBeginRun ( const edm::Run run,
const edm::EventSetup context 
)
overridevirtual

BeginRun.

Reimplemented from DQMEDAnalyzer.

Definition at line 57 of file DTSegmentAnalysisTask.cc.

References dtGeom, edm::EventSetup::getData(), and muonGeomToken_.

57  {
58  // Get the DT Geometry
59  dtGeom = &context.getData(muonGeomToken_);
60 }
const DTGeometry * dtGeom
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_
bool getData(T &iHolder) const
Definition: EventSetup.h:128
void DTSegmentAnalysisTask::fillHistos ( DTChamberId  chamberId,
int  nHits,
float  chi2 
)
private

Definition at line 269 of file DTSegmentAnalysisTask.cc.

References detailedAnalysis, mergeVDriftHistosByStation::histos, histosPerCh, DTChamberId::sector(), DTChamberId::station(), summaryHistos, and DTChamberId::wheel().

Referenced by analyze().

269  {
270  int sector = chamberId.sector();
271  if (chamberId.sector() == 13) {
272  sector = 4;
273  } else if (chamberId.sector() == 14) {
274  sector = 10;
275  }
276 
277  summaryHistos[chamberId.wheel()]->Fill(sector, chamberId.station());
278 
279  vector<MonitorElement*> histos = histosPerCh[chamberId];
280  histos[0]->Fill(nHits);
281  if (detailedAnalysis) {
282  histos[1]->Fill(chi2);
283  }
284 }
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
int sector() const
Definition: DTChamberId.h:49
std::map< int, MonitorElement * > summaryHistos
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39

Member Data Documentation

bool DTSegmentAnalysisTask::checkNoisyChannels
private

Definition at line 76 of file DTSegmentAnalysisTask.h.

Referenced by analyze(), and DTSegmentAnalysisTask().

bool DTSegmentAnalysisTask::detailedAnalysis
private

Definition at line 62 of file DTSegmentAnalysisTask.h.

Referenced by bookHistos(), DTSegmentAnalysisTask(), and fillHistos().

const DTGeometry* DTSegmentAnalysisTask::dtGeom
private

Definition at line 66 of file DTSegmentAnalysisTask.h.

Referenced by bookHistograms(), and dqmBeginRun().

std::map<DTChamberId, std::vector<MonitorElement*> > DTSegmentAnalysisTask::histosPerCh
private

Definition at line 84 of file DTSegmentAnalysisTask.h.

Referenced by bookHistos(), and fillHistos().

bool DTSegmentAnalysisTask::hltDQMMode
private

Definition at line 91 of file DTSegmentAnalysisTask.h.

Referenced by bookHistograms(), and DTSegmentAnalysisTask().

edm::ESGetToken<DTGeometry, MuonGeometryRecord> DTSegmentAnalysisTask::muonGeomToken_
private

Definition at line 65 of file DTSegmentAnalysisTask.h.

Referenced by dqmBeginRun().

MonitorElement* DTSegmentAnalysisTask::nEventMonitor
private

Definition at line 97 of file DTSegmentAnalysisTask.h.

Referenced by analyze(), and bookHistograms().

int DTSegmentAnalysisTask::nevents
private

Definition at line 87 of file DTSegmentAnalysisTask.h.

Referenced by analyze().

int DTSegmentAnalysisTask::nhitsCut
private

Definition at line 95 of file DTSegmentAnalysisTask.h.

Referenced by analyze(), and DTSegmentAnalysisTask().

double DTSegmentAnalysisTask::phiSegmCut
private

Definition at line 93 of file DTSegmentAnalysisTask.h.

Referenced by analyze(), and DTSegmentAnalysisTask().

edm::EDGetTokenT<DTRecSegment4DCollection> DTSegmentAnalysisTask::recHits4DToken_
private

Definition at line 73 of file DTSegmentAnalysisTask.h.

Referenced by analyze(), and DTSegmentAnalysisTask().

const DTStatusFlag* DTSegmentAnalysisTask::statusMap
private

Definition at line 70 of file DTSegmentAnalysisTask.h.

Referenced by analyze().

edm::ESGetToken<DTStatusFlag, DTStatusFlagRcd> DTSegmentAnalysisTask::statusMapToken_
private

Definition at line 69 of file DTSegmentAnalysisTask.h.

Referenced by analyze().

std::map<int, MonitorElement*> DTSegmentAnalysisTask::summaryHistos
private

Definition at line 85 of file DTSegmentAnalysisTask.h.

Referenced by bookHistograms(), and fillHistos().

std::string DTSegmentAnalysisTask::topHistoFolder
private

Definition at line 89 of file DTSegmentAnalysisTask.h.

Referenced by bookHistograms(), bookHistos(), and DTSegmentAnalysisTask().