test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTnoiseDBValidation Class Reference

#include <DTnoiseDBValidation.h>

Inheritance diagram for DTnoiseDBValidation:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &setup)
 
void beginRun (const edm::Run &run, const edm::EventSetup &setup)
 Operations. More...
 
 DTnoiseDBValidation (const edm::ParameterSet &pset)
 Constructor. More...
 
void endJob ()
 
void endRun (edm::Run const &, edm::EventSetup const &)
 
virtual ~DTnoiseDBValidation ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

void bookHisto (const DTChamberId &)
 

Private Attributes

DQMStoredbe_
 
MonitorElementdiffHisto_
 
std::string diffTestName_
 
edm::ESHandle< DTGeometrydtGeom_
 
std::string labelDB_
 
std::string labelDBRef_
 
MonitorElementlayerHisto_
 
std::string layerTestName_
 
std::map< DTChamberId,
MonitorElement * > 
noiseHistoMap_
 
const DTStatusFlagnoiseMap_
 
const DTStatusFlagnoiseRefMap_
 
int noisyCellsRef_
 
int noisyCellsValid_
 
std::string outputFileName_
 
bool outputMEsInRootFile_
 
MonitorElementsectorHisto_
 
std::string sectorTestName_
 
MonitorElementstationHisto_
 
std::string stationTestName_
 
MonitorElementwheelHisto_
 
std::string wheelTestName_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- 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 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

Plot the noise from the DB comparaison

Author
G. Mila - INFN Torino

Definition at line 29 of file DTnoiseDBValidation.h.

Constructor & Destructor Documentation

DTnoiseDBValidation::DTnoiseDBValidation ( const edm::ParameterSet pset)

Constructor.

Definition at line 40 of file DTnoiseDBValidation.cc.

References dbe_, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), cppFunctionSkipper::operator, DQMStore::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

40  {
41 
42  LogVerbatim("NoiseDBValidation") << "[DTnoiseDBValidation] Constructor called!";
43 
44  // Get the DQM needed services
46  dbe_->setCurrentFolder("DT/DtCalib/NoiseDBValidation");
47 
48  // Get dataBase label
49  labelDBRef_ = pset.getParameter<string>("labelDBRef");
50  labelDB_ = pset.getParameter<string>("labelDB");
51 
52  diffTestName_ = "noiseDifferenceInRange";
53  if( pset.exists("diffTestName") ) diffTestName_ = pset.getParameter<string>("diffTestName");
54 
55  wheelTestName_ = "noiseWheelOccInRange";
56  if( pset.exists("wheelTestName") ) wheelTestName_ = pset.getParameter<string>("wheelTestName");
57 
58  stationTestName_ = "noiseStationOccInRange";
59  if( pset.exists("stationTestName") ) stationTestName_ = pset.getParameter<string>("stationTestName");
60 
61  sectorTestName_ = "noiseSectorOccInRange";
62  if( pset.exists("sectorTestName") ) sectorTestName_ = pset.getParameter<string>("sectorTestName");
63 
64  layerTestName_ = "noiseLayerOccInRange";
65  if( pset.exists("layerTestName") ) layerTestName_ = pset.getParameter<string>("layerTestName");
66 
67  outputMEsInRootFile_ = false;
68  if( pset.exists("OutputFileName") ){
69  outputMEsInRootFile_ = true;
70  outputFileName_ = pset.getParameter<std::string>("OutputFileName");
71  }
72 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:667
DTnoiseDBValidation::~DTnoiseDBValidation ( )
virtual

Destructor.

Definition at line 75 of file DTnoiseDBValidation.cc.

75 {}

Member Function Documentation

void DTnoiseDBValidation::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
inlinevirtual

Implements edm::EDAnalyzer.

Definition at line 41 of file DTnoiseDBValidation.h.

41 {}
void DTnoiseDBValidation::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)
virtual

Operations.

Reimplemented from edm::EDAnalyzer.

Definition at line 77 of file DTnoiseDBValidation.cc.

References funct::abs(), plotTTrigFromHistos::binNumber(), DQMStore::book1D(), DTSuperLayerId::chamberId(), dbe_, edm::EventSetup::get(), DTTopology::lastChannel(), DTLayerId::layer(), DTWireId::layerId(), pileupReCalc_HLTpaths::scale, MonitorElement::setBinLabel(), relativeConstraints::station, DTSuperLayerId::superLayer(), DTLayerId::superlayerId(), and DTWireId::wire().

77  {
78  ESHandle<DTStatusFlag> noiseRef;
79  setup.get<DTStatusFlagRcd>().get(labelDBRef_, noiseRef);
80  noiseRefMap_ = &*noiseRef;
81 
82  ESHandle<DTStatusFlag> noiseValid;
83  setup.get<DTStatusFlagRcd>().get(labelDB_, noiseValid);
84  noiseMap_ = &*noiseValid;
85 
86  // Get the geometry
87  setup.get<MuonGeometryRecord>().get(dtGeom_);
88 
89  LogVerbatim("NoiseDBValidation")<<"[DTnoiseDBValidation] Parameters initialization";
90 
91  noisyCellsRef_ = 0;
92  noisyCellsValid_ = 0;
93 
94  // Histo booking
95  diffHisto_ = dbe_->book1D("noisyCellDiff", "percentual (wrt the previous db) total number of noisy cells",1, 0.5, 1.5);
96  diffHisto_->setBinLabel(1,"Diff");
97  wheelHisto_ = dbe_->book1D("wheelOccupancy", "percentual noisy cells occupancy per wheel",5, -2.5, 2.5);
98  wheelHisto_->setBinLabel(1,"Wh-2");
99  wheelHisto_->setBinLabel(2,"Wh-1");
100  wheelHisto_->setBinLabel(3,"Wh0");
101  wheelHisto_->setBinLabel(4,"Wh1");
102  wheelHisto_->setBinLabel(5,"Wh2");
103  stationHisto_ = dbe_->book1D("stationOccupancy", "percentual noisy cells occupancy per station",4, 0.5, 4.5);
104  stationHisto_->setBinLabel(1,"St1");
105  stationHisto_->setBinLabel(2,"St2");
106  stationHisto_->setBinLabel(3,"St3");
107  stationHisto_->setBinLabel(4,"St4");
108  sectorHisto_ = dbe_->book1D("sectorOccupancy", "percentual noisy cells occupancy per sector",12, 0.5, 12.5);
109  sectorHisto_->setBinLabel(1,"Sect1");
110  sectorHisto_->setBinLabel(2,"Sect2");
111  sectorHisto_->setBinLabel(3,"Sect3");
112  sectorHisto_->setBinLabel(4,"Sect4");
113  sectorHisto_->setBinLabel(5,"Sect5");
114  sectorHisto_->setBinLabel(6,"Sect6");
115  sectorHisto_->setBinLabel(7,"Sect7");
116  sectorHisto_->setBinLabel(8,"Sect8");
117  sectorHisto_->setBinLabel(9,"Sect9");
118  sectorHisto_->setBinLabel(10,"Sect10");
119  sectorHisto_->setBinLabel(11,"Sect11");
120  sectorHisto_->setBinLabel(12,"Sect12");
121  layerHisto_ = dbe_->book1D("layerOccupancy", "percentual noisy cells occupancy per layer",3, 0.5, 3.5);
122  layerHisto_->setBinLabel(1,"First 10 bins");
123  layerHisto_->setBinLabel(2,"Middle bins");
124  layerHisto_->setBinLabel(3,"Last 10 bins");
125 
126  // map initialization
127  map<int, int> whMap;
128  whMap.clear();
129  map<int, int> stMap;
130  stMap.clear();
131  map<int, int> sectMap;
132  sectMap.clear();
133  map<int, int> layerMap;
134  layerMap.clear();
135 
136  // Loop over reference DB entries
138  noise != noiseRefMap_->end(); noise++) {
139  DTWireId wireId((*noise).first.wheelId,
140  (*noise).first.stationId,
141  (*noise).first.sectorId,
142  (*noise).first.slId,
143  (*noise).first.layerId,
144  (*noise).first.cellId);
145  LogVerbatim("NoiseDBValidation") << "Ref. noisy wire: " << wireId;
146  ++noisyCellsRef_;
147  }
148 
149  // Loop over validation DB entries
151  noise != noiseMap_->end(); noise++) {
152  DTWireId wireId((*noise).first.wheelId,
153  (*noise).first.stationId,
154  (*noise).first.sectorId,
155  (*noise).first.slId,
156  (*noise).first.layerId,
157  (*noise).first.cellId);
158  LogVerbatim("NoiseDBValidation") << "Valid. noisy wire: " << wireId;
160 
161  whMap[(*noise).first.wheelId]++;
162  stMap[(*noise).first.stationId]++;
163  sectMap[(*noise).first.sectorId]++;
164 
165  const DTTopology& dtTopo = dtGeom_->layer(wireId.layerId())->specificTopology();
166  const int lastWire = dtTopo.lastChannel();
167  if((*noise).first.cellId<=10)
168  layerMap[1]++;
169  if((*noise).first.cellId>10 && (*noise).first.cellId<(lastWire-10))
170  layerMap[2]++;
171  if((*noise).first.cellId>=(lastWire-10))
172  layerMap[3]++;
173 
174  const DTChamberId chId = wireId.layerId().superlayerId().chamberId();
175  if( noiseHistoMap_.find(chId) == noiseHistoMap_.end() ) bookHisto(chId);
176  int binNumber = 4*(wireId.superLayer() - 1) + wireId.layer();
177  noiseHistoMap_[chId]->Fill(wireId.wire(),binNumber);
178  }
179 
180  //histo filling
181  double scale = 1/double(noisyCellsRef_);
183 
184  scale = 1/double(noisyCellsValid_);
185  for(map<int, int >::const_iterator wheel = whMap.begin();
186  wheel != whMap.end();
187  wheel++) {
188  wheelHisto_->Fill((*wheel).first, ((*wheel).second)*scale);
189  }
190  for(map<int, int >::const_iterator station = stMap.begin();
191  station != stMap.end();
192  station++) {
193  stationHisto_->Fill((*station).first, ((*station).second)*scale);
194  }
195  for(map<int, int >::const_iterator sector = sectMap.begin();
196  sector != sectMap.end();
197  sector++) {
198  sectorHisto_->Fill((*sector).first, ((*sector).second)*scale);
199  }
200  for(map<int, int >::const_iterator layer = layerMap.begin();
201  layer != layerMap.end();
202  layer++) {
203  layerHisto_->Fill((*layer).first, ((*layer).second)*scale);
204  }
205 
206 }
MonitorElement * wheelHisto_
const_iterator end() const
std::map< DTChamberId, MonitorElement * > noiseHistoMap_
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:954
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
void bookHisto(const DTChamberId &)
const DTStatusFlag * noiseRefMap_
const DTStatusFlag * noiseMap_
void Fill(long long x)
int lastChannel() const
Returns the wire number of the last wire.
Definition: DTTopology.h:80
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::pair< DTStatusFlagId, DTStatusFlagData > >::const_iterator const_iterator
Access methods to data.
Definition: DTStatusFlag.h:262
MonitorElement * diffHisto_
edm::ESHandle< DTGeometry > dtGeom_
const T & get() const
Definition: EventSetup.h:55
MonitorElement * stationHisto_
MonitorElement * layerHisto_
const_iterator begin() const
MonitorElement * sectorHisto_
void DTnoiseDBValidation::bookHisto ( const DTChamberId chId)
private

Definition at line 271 of file DTnoiseDBValidation.cc.

References plotTTrigFromHistos::binNumber(), DQMStore::book2D(), dbe_, diffTwoXMLs::label, LayerTriplets::layers(), DTChamberId::sector(), DTChamberId::station(), DTChamber::superLayers(), and DTChamberId::wheel().

271  {
272  stringstream histoName;
273  histoName << "NoiseOccupancy"
274  << "_W" << chId.wheel()
275  <<"_St" << chId.station()
276  << "_Sec" << chId.sector();
277 
278  if( noiseHistoMap_.find(chId) == noiseHistoMap_.end() ){ // Redundant check
279  // Get the chamber from the geometry
280  int nWiresMax = 0;
281  const DTChamber* dtchamber = dtGeom_->chamber(chId);
282  const vector<const DTSuperLayer*> superlayers = dtchamber->superLayers();
283 
284  // Loop over layers and find the max # of wires
285  for(vector<const DTSuperLayer*>::const_iterator sl = superlayers.begin();
286  sl != superlayers.end(); ++sl) { // loop over SLs
287  vector<const DTLayer*> layers = (*sl)->layers();
288  for(vector<const DTLayer*>::const_iterator lay = layers.begin();
289  lay != layers.end(); ++lay) { // loop over layers
290  int nWires = (*lay)->specificTopology().channels();
291  if(nWires > nWiresMax) nWiresMax = nWires;
292  }
293  }
294 
295  noiseHistoMap_[chId] = dbe_->book2D(histoName.str(),"Noise occupancy",nWiresMax,1,(nWiresMax + 1),12,1,13);
296  for(int i_sl = 1; i_sl <= 3; ++i_sl) {
297  for(int i_lay = 1; i_lay <= 4; ++i_lay) {
298  int binNumber = 4*(i_sl - 1) + i_lay;
299  stringstream label;
300  label << "SL" << i_sl << ": L" << i_lay;
301  noiseHistoMap_[chId]->setBinLabel(binNumber,label.str(),2);
302  }
303  }
304  }
305 
306 }
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
std::map< DTChamberId, MonitorElement * > noiseHistoMap_
const std::vector< const DTSuperLayer * > & superLayers() const
Return the superlayers in the chamber.
Definition: DTChamber.cc:60
edm::ESHandle< DTGeometry > dtGeom_
int sector() const
Definition: DTChamberId.h:61
int station() const
Return the station number.
Definition: DTChamberId.h:51
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:1082
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void DTnoiseDBValidation::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 266 of file DTnoiseDBValidation.cc.

References dbe_, and DQMStore::save().

266  {
267  // Write the histos in a ROOT file
269 }
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, const uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE", const bool resetMEsAfterWriting=false)
Definition: DQMStore.cc:2540
void DTnoiseDBValidation::endRun ( edm::Run const &  run,
edm::EventSetup const &  setup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 208 of file DTnoiseDBValidation.cc.

References QReport::getBadChannels().

208  {
209 
210  // test on difference histo
211  //string testCriterionName;
212  //testCriterionName = parameters.getUntrackedParameter<string>("diffTestName","noiseDifferenceInRange");
213  const QReport * theDiffQReport = diffHisto_->getQReport(diffTestName_);
214  if(theDiffQReport) {
215  vector<dqm::me_util::Channel> badChannels = theDiffQReport->getBadChannels();
216  for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
217  channel != badChannels.end(); channel++) {
218  LogWarning("NoiseDBValidation") << " Bad partial difference of noisy channels! Contents : " << (*channel).getContents();
219  }
220  }
221  //testCriterionName = parameters.getUntrackedParameter<string>("wheelTestName","noiseWheelOccInRange");
222  const QReport * theDiffQReport2 = wheelHisto_->getQReport(wheelTestName_);
223  if(theDiffQReport2) {
224  vector<dqm::me_util::Channel> badChannels = theDiffQReport2->getBadChannels();
225  for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
226  channel != badChannels.end(); channel++) {
227  int wheel = (*channel).getBin()-3;
228  LogWarning("NoiseDBValidation") << " Bad percentual occupancy for wheel : " << wheel << " Contents : " << (*channel).getContents();
229  }
230  }
231  //testCriterionName = parameters.getUntrackedParameter<string>("stationTestName","noiseStationOccInRange");
232  const QReport * theDiffQReport3 = stationHisto_->getQReport(stationTestName_);
233  if(theDiffQReport3) {
234  vector<dqm::me_util::Channel> badChannels = theDiffQReport3->getBadChannels();
235  for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
236  channel != badChannels.end(); channel++) {
237  LogWarning("NoiseDBValidation") << " Bad percentual occupancy for station : " << (*channel).getBin() << " Contents : " << (*channel).getContents();
238  }
239  }
240  //testCriterionName = parameters.getUntrackedParameter<string>("sectorTestName","noiseSectorOccInRange");
241  const QReport * theDiffQReport4 = sectorHisto_->getQReport(sectorTestName_);
242  if(theDiffQReport4) {
243  vector<dqm::me_util::Channel> badChannels = theDiffQReport4->getBadChannels();
244  for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
245  channel != badChannels.end(); channel++) {
246  LogWarning("NoiseDBValidation") << " Bad percentual occupancy for sector : " << (*channel).getBin() << " Contents : " << (*channel).getContents();
247  }
248  }
249  //testCriterionName = parameters.getUntrackedParameter<string>("layerTestName","noiseLayerOccInRange");
250  const QReport * theDiffQReport5 = layerHisto_->getQReport(layerTestName_);
251  if(theDiffQReport5) {
252  vector<dqm::me_util::Channel> badChannels = theDiffQReport5->getBadChannels();
253  for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
254  channel != badChannels.end(); channel++) {
255  if((*channel).getBin()==1)
256  LogWarning("NoiseDBValidation") << " Bad percentual occupancy for the first 10 wires! Contents : " << (*channel).getContents();
257  if((*channel).getBin()==2)
258  LogWarning("NoiseDBValidation") << " Bad percentual occupancy for the middle wires! Contents : " << (*channel).getContents();
259  if((*channel).getBin()==3)
260  LogWarning("NoiseDBValidation") << " Bad percentual occupancy for the last 10 wires! Contents : "<<(*channel).getContents();
261  }
262  }
263 
264 }
MonitorElement * wheelHisto_
const QReport * getQReport(const std::string &qtname) const
get QReport corresponding to &lt;qtname&gt; (null pointer if QReport does not exist)
MonitorElement * diffHisto_
const std::vector< DQMChannel > & getBadChannels(void) const
Definition: QReport.h:33
MonitorElement * stationHisto_
MonitorElement * layerHisto_
MonitorElement * sectorHisto_

Member Data Documentation

DQMStore* DTnoiseDBValidation::dbe_
private

Definition at line 48 of file DTnoiseDBValidation.h.

MonitorElement* DTnoiseDBValidation::diffHisto_
private

Definition at line 69 of file DTnoiseDBValidation.h.

std::string DTnoiseDBValidation::diffTestName_
private

Definition at line 52 of file DTnoiseDBValidation.h.

edm::ESHandle<DTGeometry> DTnoiseDBValidation::dtGeom_
private

Definition at line 59 of file DTnoiseDBValidation.h.

std::string DTnoiseDBValidation::labelDB_
private

Definition at line 51 of file DTnoiseDBValidation.h.

std::string DTnoiseDBValidation::labelDBRef_
private

Definition at line 50 of file DTnoiseDBValidation.h.

MonitorElement* DTnoiseDBValidation::layerHisto_
private

Definition at line 73 of file DTnoiseDBValidation.h.

std::string DTnoiseDBValidation::layerTestName_
private

Definition at line 52 of file DTnoiseDBValidation.h.

std::map<DTChamberId, MonitorElement*> DTnoiseDBValidation::noiseHistoMap_
private

Definition at line 74 of file DTnoiseDBValidation.h.

const DTStatusFlag* DTnoiseDBValidation::noiseMap_
private

Definition at line 62 of file DTnoiseDBValidation.h.

const DTStatusFlag* DTnoiseDBValidation::noiseRefMap_
private

Definition at line 63 of file DTnoiseDBValidation.h.

int DTnoiseDBValidation::noisyCellsRef_
private

Definition at line 66 of file DTnoiseDBValidation.h.

int DTnoiseDBValidation::noisyCellsValid_
private

Definition at line 67 of file DTnoiseDBValidation.h.

std::string DTnoiseDBValidation::outputFileName_
private

Definition at line 56 of file DTnoiseDBValidation.h.

bool DTnoiseDBValidation::outputMEsInRootFile_
private

Definition at line 55 of file DTnoiseDBValidation.h.

MonitorElement* DTnoiseDBValidation::sectorHisto_
private

Definition at line 72 of file DTnoiseDBValidation.h.

std::string DTnoiseDBValidation::sectorTestName_
private

Definition at line 52 of file DTnoiseDBValidation.h.

MonitorElement* DTnoiseDBValidation::stationHisto_
private

Definition at line 71 of file DTnoiseDBValidation.h.

std::string DTnoiseDBValidation::stationTestName_
private

Definition at line 52 of file DTnoiseDBValidation.h.

MonitorElement* DTnoiseDBValidation::wheelHisto_
private

Definition at line 70 of file DTnoiseDBValidation.h.

std::string DTnoiseDBValidation::wheelTestName_
private

Definition at line 52 of file DTnoiseDBValidation.h.