CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTDeadChannelTest.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * \author G. Mila - INFN Torino
6  */
7 
8 
10 
11 // Framework
13 
14 
15 // Geometry
20 
23 
27 
28 #include <stdio.h>
29 #include <sstream>
30 #include <math.h>
31 
32 
33 using namespace edm;
34 using namespace std;
35 
36 
38 
39  edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: Constructor";
40 
41  parameters = ps;
42 
43  dbe = edm::Service<DQMStore>().operator->();
44 
45  prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1);
46 
47 }
48 
50 
51  edm::LogVerbatim ("deadChannel") << "DTDeadChannelTest: analyzed " << nevents << " events";
52 
53 }
54 
55 
57 
58  edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: BeginJob";
59 
60  nevents = 0;
61 
62 }
63 
64 void DTDeadChannelTest::beginRun(Run const& run, EventSetup const& context) {
65 
66  edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: BeginRun";
67 
68  // Get the geometry
69  context.get<MuonGeometryRecord>().get(muonGeom);
70 
71 }
72 
73 
75 
76  edm::LogVerbatim ("deadChannel") <<"[DTDeadChannelTest]: Begin of LS transition";
77 
78  // Get the run number
79  run = lumiSeg.run();
80 
81 }
82 
83 
84 
86 
87  nevents++;
88  edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: "<<nevents<<" events";
89 
90 }
91 
92 
93 
95 
96  // counts number of updats (online mode) or number of events (standalone mode)
97  //nevents++;
98  // if running in standalone perform diagnostic only after a reasonalbe amount of events
99  //if ( parameters.getUntrackedParameter<bool>("runningStandalone", false) &&
100  // nevents%parameters.getUntrackedParameter<int>("diagnosticPrescale", 1000) != 0 ) return;
101  //edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: "<<nevents<<" updates";
102 
103 
104  edm::LogVerbatim ("deadChannel") <<"[DTDeadChannelTest]: End of LS transition, performing the DQM client operation";
105 
106  // counts number of lumiSegs
107  nLumiSegs = lumiSeg.id().luminosityBlock();
108 
109  // prescale factor
110  if ( nLumiSegs%prescaleFactor != 0 ) return;
111 
112  edm::LogVerbatim ("deadChannel") <<"[DTDeadChannelTest]: "<<nLumiSegs<<" updates";
113 
114 
115  vector<const DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
116  vector<const DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
117 
118  edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: Occupancy tests results";
119 
120  // Loop over the chambers
121  for (; ch_it != ch_end; ++ch_it) {
122  DTChamberId chID = (*ch_it)->id();
123  vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
124  vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
125 
126  stringstream wheel; wheel << chID.wheel();
127  stringstream station; station << chID.station();
128  stringstream sector; sector << chID.sector();
129 
130  context.get<DTTtrigRcd>().get(tTrigMap);
131 
132  string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
133 
134  // Get the ME produced by DigiTask Source
135  MonitorElement * noise_histo = dbe->get(getMEName("OccupancyNoise_perCh", chID));
136  MonitorElement * hitInTime_histo = dbe->get(getMEName("OccupancyInTimeHits_perCh", chID));
137 
138  // ME -> TH2F
139  if(noise_histo && hitInTime_histo) {
140  TH2F * noise_histo_root = noise_histo->getTH2F();
141  TH2F * hitInTime_histo_root = hitInTime_histo->getTH2F();
142 
143  // Loop over the SuperLayers
144  for(; sl_it != sl_end; ++sl_it) {
145  DTSuperLayerId slID = (*sl_it)->id();
146  vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin();
147  vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
148 
149  // ttrig and rms are counts
150  float tTrig, tTrigRMS, kFactor;
151  tTrigMap->get(slID, tTrig, tTrigRMS, kFactor, DTTimeUnits::counts);
152 
153  // Loop over the layers
154  for(; l_it != l_end; ++l_it) {
155  DTLayerId lID = (*l_it)->id();
156 
157  //Parameters to fill histos
158  stringstream superLayer; superLayer << slID.superlayer();
159  stringstream layer; layer << lID.layer();
160  string HistoNameTest = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() + "_SL" + superLayer.str() + "_L" + layer.str();
161 
162  const int firstWire = muonGeom->layer(lID)->specificTopology().firstChannel();
163  const int lastWire = muonGeom->layer(lID)->specificTopology().lastChannel();
164 
165  int entry=-1;
166  if(slID.superlayer() == 1) entry=0;
167  if(slID.superlayer() == 2) entry=4;
168  if(slID.superlayer() == 3) entry=8;
169  int YBinNumber = entry+lID.layer();
170 
171 
172  // Loop over the TH2F bin and fill the ME to be used for the Quality Test
173  for(int bin=firstWire; bin <= lastWire; bin++) {
174  if (OccupancyDiffHistos.find(HistoNameTest) == OccupancyDiffHistos.end()) bookHistos(lID, firstWire, lastWire);
175  // tMax default value
176  float tMax = 450.0;
177 
178  float difference = (hitInTime_histo_root->GetBinContent(bin, YBinNumber) / tMax)
179  - (noise_histo_root->GetBinContent(bin, YBinNumber) / tTrig);
180  OccupancyDiffHistos.find(HistoNameTest)->second->setBinContent(bin, difference);
181  }
182  } // loop on layers
183  } // loop on superlayers
184  }
185  } // loop on chambers
186 
187  // Occupancy Difference test
188  string OccupancyDiffCriterionName = parameters.getUntrackedParameter<string>("OccupancyDiffTestName","OccupancyDiffInRange");
189  for(map<string, MonitorElement*>::const_iterator hOccDiff = OccupancyDiffHistos.begin();
190  hOccDiff != OccupancyDiffHistos.end();
191  hOccDiff++) {
192  const QReport * theOccupancyDiffQReport = (*hOccDiff).second->getQReport(OccupancyDiffCriterionName);
193  if(theOccupancyDiffQReport) {
194  vector<dqm::me_util::Channel> badChannels = theOccupancyDiffQReport->getBadChannels();
195  for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
196  channel != badChannels.end(); channel++) {
197  edm::LogError ("deadChannel") << "Layer : "<<(*hOccDiff).first<<" Bad occupancy difference channels: "<<(*channel).getBin()<<" Contents : "<<(*channel).getContents();
198  }
199  // FIXME: getMessage() sometimes returns and invalid string (null pointer inside QReport data member)
200  // edm::LogWarning("deadChannel")<< "-------- Layer : "<<(*hOccDiff).first<<" "<<theOccupancyDiffQReport->getMessage()<<" ------- "<<theOccupancyDiffQReport->getStatus();
201  }
202  }
203 
204 }
205 
206 
208 
209  edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest] endjob called!";
210 
211  dbe->rmdir("DT/Tests/DTDeadChannel");
212 
213 }
214 
215 string DTDeadChannelTest::getMEName(string histoTag, const DTChamberId & chId) {
216 
217  stringstream wheel; wheel << chId.wheel();
218  stringstream station; station << chId.station();
219  stringstream sector; sector << chId.sector();
220 
221  string folderRoot = parameters.getUntrackedParameter<string>("folderRoot", "Collector/FU0/");
222  string folderName =
223  folderRoot + "DT/DTDigiTask/Wheel" + wheel.str() +
224  "/Station" + station.str() +
225  "/Sector" + sector.str() +
226  "/Occupancies" + "/";
227 
228  string histoname = folderName + histoTag
229  + "_W" + wheel.str()
230  + "_St" + station.str()
231  + "_Sec" + sector.str();
232 
233  return histoname;
234 
235 }
236 
237 
238 void DTDeadChannelTest::bookHistos(const DTLayerId & lId, int firstWire, int lastWire) {
239 
240  stringstream wheel; wheel << lId.superlayerId().wheel();
241  stringstream station; station << lId.superlayerId().station();
242  stringstream sector; sector << lId.superlayerId().sector();
243  stringstream superLayer; superLayer << lId.superlayerId().superlayer();
244  stringstream layer; layer << lId.layer();
245 
246  string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() + "_SL" + superLayer.str() + "_L" + layer.str();
247  string OccupancyDiffHistoName = "OccupancyDiff_" + HistoName;
248 
249  dbe->setCurrentFolder("DT/Tests/DTDeadChannel/Wheel" + wheel.str() +
250  "/Station" + station.str() +
251  "/Sector" + sector.str());
252 
253  OccupancyDiffHistos[HistoName] = dbe->book1D(OccupancyDiffHistoName.c_str(),OccupancyDiffHistoName.c_str(),lastWire-firstWire+1, firstWire-0.5, lastWire+0.5);
254 
255 }
LuminosityBlockID id() const
dictionary parameters
Definition: Parameters.py:2
std::string getMEName(std::string histoTag, const DTChamberId &chId)
Get the ME name.
void endJob()
Endjob.
int layer() const
Return the layer number.
Definition: DTLayerId.h:53
DTDeadChannelTest(const edm::ParameterSet &ps)
Constructor.
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:59
void bookHistos()
Definition: Histogram.h:33
void analyze(const edm::Event &e, const edm::EventSetup &c)
Analyze.
RunNumber_t run() const
int nevents
void beginJob()
BeginJob.
const std::vector< DQMChannel > & getBadChannels(void) const
Definition: QReport.h:33
void beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context)
int superlayer() const
Return the superlayer number (deprecated method name)
const T & get() const
Definition: EventSetup.h:55
LuminosityBlockNumber_t luminosityBlock() const
string histoname
Definition: dt-layouts.py:48
std::string HistoName
void bookHistos(const DTLayerId &ch, int firstWire, int lastWire)
book the new ME
int sector() const
Definition: DTChamberId.h:61
int station() const
Return the station number.
Definition: DTChamberId.h:51
TH2F * getTH2F(void) const
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
virtual ~DTDeadChannelTest()
Destructor.
Definition: Run.h:41
void endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &c)
DQM Client Diagnostic.
void beginRun(edm::Run const &run, edm::EventSetup const &context)
BeginRun.