CMS 3D CMS Logo

DTOccupancyTestML.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author G. Cerminara - University and INFN Torino
5  *
6  * threadsafe version (//-) oct/nov 2014 - WATWanAbdullah -ncpp-um-my
7  *
8  */
9 
10 
12 
23 
24 #include "TMath.h"
25 
27 #include <vector>
28 
29 using namespace edm;
30 using namespace std;
31 
32 
33 
34 
36  LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTestML") << "[DTOccupancyTestML]: Constructor";
37 
38  // Get the DQM service
39 
40  lsCounter = 0;
41 
42  writeRootFile = ps.getUntrackedParameter<bool>("writeRootFile", false);
43  if(writeRootFile) {
44  rootFile = new TFile("MLDTOccupancyTest.root","RECREATE");
45  ntuple = new TNtuple("OccupancyNtuple", "OccupancyNtuple", "ls:wh:st:se:lay1MeanCell:lay1RMS:lay2MeanCell:lay2RMS:lay3MeanCell:lay3RMS:lay4MeanCell:lay4RMS:lay5MeanCell:lay5RMS:lay6MeanCell:lay6RMS:lay7MeanCell:lay7RMS:lay8MeanCell:lay8RMS:lay9MeanCell:lay9RMS:lay10MeanCell:lay10RMS:lay11MeanCell:lay11RMS:lay12MeanCell:lay12RMS");
46  }
47 
48  // switch on the mode for running on test pulses (different top folder)
49  tpMode = ps.getUntrackedParameter<bool>("testPulseMode", false);
50 
51  runOnAllHitsOccupancies = ps.getUntrackedParameter<bool>("runOnAllHitsOccupancies", true);
52  runOnNoiseOccupancies = ps.getUntrackedParameter<bool>("runOnNoiseOccupancies", false);
53  runOnInTimeOccupancies = ps.getUntrackedParameter<bool>("runOnInTimeOccupancies", false);
54  nMinEvts = ps.getUntrackedParameter<int>("nEventsCert", 5000);
55  nMinEvtsPC = ps.getUntrackedParameter<int>("nEventsMinPC", 2200);
56  nZeroEvtsPC = ps.getUntrackedParameter<int>("nEventsZeroPC", 30);
57 
58  bookingdone = false;
59 
60  // Event counter
61  nevents = 0;
62 
63 }
64 
66  LogVerbatim ("DTDQM|DTMonitorClient|MLDTOccupancyTest") << " destructor called" << endl;
67 }
68 
69 
70 void DTOccupancyTestML::beginRun(const edm::Run& run, const EventSetup& context){
71 
72  LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTestML") << "[DTOccupancyTestML]: BeginRun";
73 
74  // Get the geometry
75  context.get<MuonGeometryRecord>().get(muonGeom);
76 
77 }
78 
80  edm::LuminosityBlock const & lumiSeg, edm::EventSetup const & context) {
81  if (!bookingdone) {
82 
83  // Book the summary histos
84  // - one summary per wheel
85  for(int wh = -2; wh <= 2; ++wh) { // loop over wheels
86  bookHistos(ibooker,wh, string("MLDTOccupancies"), "MLOccupancySummary");
87  }
88 
89  ibooker.setCurrentFolder(topFolder(true));
90  string title = "Occupancy Summary";
91  if(tpMode) {
92  title = "Test Pulse Occupancy Summary";
93  }
94  // - global summary with alarms
95  summaryHisto = ibooker.book2D("MLOccupancySummary",title.c_str(),12,1,13,5,-2,3);
96  summaryHisto->setAxisTitle("sector",1);
97  summaryHisto->setAxisTitle("wheel",2);
98 
99  // - global summary with percentages
100  glbSummaryHisto = ibooker.book2D("MLOccupancyGlbSummary",title.c_str(),12,1,13,5,-2,3);
101  glbSummaryHisto->setAxisTitle("sector",1);
102  glbSummaryHisto->setAxisTitle("wheel",2);
103 
104 
105  // assign the name of the input histogram
106  if(runOnAllHitsOccupancies) {
107  nameMonitoredHisto = "OccupancyAllHits_perCh";
108  } else if(runOnNoiseOccupancies) {
109  nameMonitoredHisto = "OccupancyNoise_perCh";
110  } else if(runOnInTimeOccupancies) {
111  nameMonitoredHisto = "OccupancyInTimeHits_perCh";
112  } else { // default is AllHits histo
113  nameMonitoredHisto = "OccupancyAllHits_perCh";
114  }
115 
116  }
117  bookingdone = true;
118 
119 
120  LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTestML")
121  <<"[DTOccupancyTestML]: End of LS transition, performing the DQM client operation";
122  lsCounter++;
123 
124 
125 
126  // Reset the global summary
127  summaryHisto->Reset();
128  glbSummaryHisto->Reset();
129 
130  nChannelTotal = 0;
131  nChannelDead = 0;
132 
133  // Get all the DT chambers
134  vector<const DTChamber*> chambers = muonGeom->chambers();
135 
136  // Load graph
138  edm::FileInPath modelFilePath("DQM/DTMonitorClient/data/occupancy_cnn_v1.pb");
139  tensorflow::GraphDef* graphDef = tensorflow::loadGraphDef(modelFilePath.fullPath());
140 
141  // Create session
142  tensorflow::Session* session = tensorflow::createSession(graphDef);
143 
144  for(vector<const DTChamber*>::const_iterator chamber = chambers.begin();
145  chamber != chambers.end(); ++chamber) { // Loop over all chambers
146  DTChamberId chId = (*chamber)->id();
147 
148  MonitorElement * chamberOccupancyHisto = igetter.get(getMEName(nameMonitoredHisto, chId));
149 
150  // Run the tests on the plot for the various granularities
151  if(chamberOccupancyHisto != nullptr) {
152  // Get the 2D histo
153  TH2F *histo = chamberOccupancyHisto->getTH2F();
154 
155  float chamberPercentage = 1.;
156  int result = runOccupancyTest(histo, chId, chamberPercentage, graphDef, session);
157  int sector = chId.sector();
158 
159  if(sector == 13) {
160  sector = 4;
161  float resultSect4 = wheelHistos[chId.wheel()]->getBinContent(sector, chId.station());
162  if(resultSect4 > result) {
163  result = (int)resultSect4;
164  }
165  } else if(sector == 14) {
166  sector = 10;
167  float resultSect10 = wheelHistos[chId.wheel()]->getBinContent(sector, chId.station());
168  if(resultSect10 > result) {
169  result = (int)resultSect10;
170  }
171  }
172 
173  // the 2 MB4 of Sect 4 and 10 count as half a chamber
174  if((sector == 4 || sector == 10) && chId.station() == 4)
175  chamberPercentage = chamberPercentage/2.;
176 
177  wheelHistos[chId.wheel()]->setBinContent(sector, chId.station(),result);
178  if(result > summaryHisto->getBinContent(sector, chId.wheel()+3)) {
179  summaryHisto->setBinContent(sector, chId.wheel()+3, result);
180  }
181  glbSummaryHisto->Fill(sector, chId.wheel(), chamberPercentage*1./4.);
182  } else {
183  LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTestML") << "[DTOccupancyTestML] ME: "
184  << getMEName(nameMonitoredHisto, chId) << " not found!" << endl;
185  }
186 
187  }
188 
189  // Clean up neural network graph
190  tensorflow::closeSession(session);
191  delete graphDef;
192 
193  string nEvtsName = "DT/EventInfo/Counters/nProcessedEventsDigi";
194 
195  MonitorElement * meProcEvts = igetter.get(nEvtsName);
196 
197  if (meProcEvts) {
198  int nProcEvts = meProcEvts->getFloatValue();
199  glbSummaryHisto->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
200  summaryHisto->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
201  } else {
202  glbSummaryHisto->setEntries(nMinEvts +1);
203  summaryHisto->setEntries(nMinEvts + 1);
204  LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTestML") << "[DTOccupancyTestML] ME: "
205  << nEvtsName << " not found!" << endl;
206  }
207 
208  // Fill the global summary
209  // Check for entire sectors off and report them on the global summary
210  //FIXME: TODO
211 
212  if(writeRootFile) ntuple->AutoSave("SaveSelf");
213 
214 }
215 
217 
218  LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTestML") << "[DTOccupancyTestML] endjob called!";
219  if(writeRootFile) {
220  rootFile->cd();
221  ntuple->Write();
222  rootFile->Close();
223  }
224 }
225 
226 
227 
228 // --------------------------------------------------
229 
230 void DTOccupancyTestML::bookHistos(DQMStore::IBooker & ibooker, const int wheelId,
231  string folder, string histoTag) {
232  // Set the current folder
233  stringstream wheel; wheel << wheelId;
234 
235  ibooker.setCurrentFolder(topFolder(true));
236 
237  // build the histo name
238  string histoName = histoTag + "_W" + wheel.str();
239 
240 
241  LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTestML") <<"[DTOccupancyTestML]: booking wheel histo:"
242  << histoName
243  << " (tag "
244  << histoTag
245  << ") in: "
246  << topFolder(true) + "Wheel"+ wheel.str() + "/" + folder << endl;
247 
248  string histoTitle = "Occupancy summary WHEEL: "+wheel.str();
249  if(tpMode) {
250  histoTitle = "TP Occupancy summary WHEEL: "+wheel.str();
251  }
252 
253  wheelHistos[wheelId] = ibooker.book2D(histoName,histoTitle,12,1,13,4,1,5);
254  wheelHistos[wheelId]->setBinLabel(1,"MB1",2);
255  wheelHistos[wheelId]->setBinLabel(2,"MB2",2);
256  wheelHistos[wheelId]->setBinLabel(3,"MB3",2);
257  wheelHistos[wheelId]->setBinLabel(4,"MB4",2);
258  wheelHistos[wheelId]->setAxisTitle("sector",1);
259 }
260 
261 
262 
263 string DTOccupancyTestML::getMEName(string histoTag, const DTChamberId& chId) {
264 
265  stringstream wheel; wheel << chId.wheel();
266  stringstream station; station << chId.station();
267  stringstream sector; sector << chId.sector();
268 
269 
270  string folderRoot = topFolder(false) + "Wheel" + wheel.str() +
271  "/Sector" + sector.str() +
272  "/Station" + station.str() + "/";
273 
274  // build the histo name
275  string histoName = histoTag
276  + "_W" + wheel.str()
277  + "_St" + station.str()
278  + "_Sec" + sector.str();
279 
280  string histoname = folderRoot + histoName;
281 
282  return histoname;
283 }
284 
285 
286 int DTOccupancyTestML::getIntegral(TH2F *histo, int firstBinX, int lastBinX, int firstBinY, int lastBinY, bool doall) {
287 
288  int sum = 0;
289  for (Int_t i = firstBinX; i < lastBinX+1; i++) {
290  for (Int_t j = firstBinY; j < lastBinY+1; j++) {
291 
292  if (histo->GetBinContent(i,j) >0){
293  if (!doall) return 1;
294  sum += histo->GetBinContent(i,j);
295  }
296  }
297  }
298 
299  return sum;
300 }
301 
303  float& chamberPercentage,
304  tensorflow::GraphDef *graphDef,
305  tensorflow::Session *session) {
306 
307  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTestML")
308  << "--- Occupancy test for chamber: " << chId << endl;
309 
310  // Initialize counters
311  int totalLayers = 0;
312  int badLayers = 0;
313 
314  // Loop over the super layers
315  for (int superLayer = 1; superLayer <= 3; superLayer++) {
316  int binYlow = ((superLayer-1)*4)+1;
317 
318  // Skip for non-existent super layers
319  if (chId.station() == 4 && superLayer == 2) continue;
320 
321  // Loop over layers
322  for (int layer = 1; layer <= 4; layer++) {
323  DTLayerId layID(chId, superLayer, layer);
324  int firstWire = muonGeom->layer(layID)->
325  specificTopology().firstChannel();
326  int nWires = muonGeom->layer(layID)->specificTopology().channels();
327  int binY = binYlow+(layer-1);
328  std::vector<float> layerOccupancy(nWires);
329  int channelId = 0;
330 
331  // Loop over cells within a layer
332  for (int cell = firstWire; cell != (nWires + firstWire); cell++) {
333  double cellOccupancy = histo->GetBinContent(cell, binY);
334  layerOccupancy.at(channelId) = cellOccupancy;
335  channelId++;
336  }
337 
338  int targetSize = 47;
339  std::vector<float> reshapedLayerOccupancy = interpolateLayers(layerOccupancy,
340  nWires,
341  targetSize);
342 
343  // Scale occupancy
344  float maxOccupancyInLayer =
345  *std::max_element(reshapedLayerOccupancy.begin(),
346  reshapedLayerOccupancy.end());
347 
348  if (maxOccupancyInLayer != 0) {
349  for (auto &val : reshapedLayerOccupancy)
350  val /= maxOccupancyInLayer;
351  }
352 
353  // Define input
354  tensorflow::Tensor input(tensorflow::DT_FLOAT, {1, targetSize});
355  for (int i = 0; i < targetSize; i++)
356  input.matrix<float>()(0, i) = reshapedLayerOccupancy[i];
357 
358  std::vector<tensorflow::Tensor> outputs;
359  tensorflow::run(session, { { "input_cnn_input", input } },
360  { "output_cnn/Softmax" }, &outputs);
361 
362  totalLayers++;
363  bool isBad = outputs[0].matrix<float>()(0, 1) > 0.95;
364  if (isBad) badLayers++;
365  }
366  }
367 
368  // Calculate a fraction of good layers
369  chamberPercentage = 1.0 - static_cast<float>(badLayers)/totalLayers;
370 
371  if (badLayers > 8) return 3; // 3 SL
372  if (badLayers > 4) return 2; // 2 SL
373  if (badLayers > 0) return 1; // 1 SL
374  return 0;
375 }
376 
377 std::vector<float> DTOccupancyTestML::interpolateLayers(std::vector<float> const& inputs, int size, int targetSize) {
378  // Reshape layers with linear interpolation
379  int interpolationFloor = 0;
380  float interpolationPoint = 0.;
381  float step = static_cast<float>(size) / targetSize;
382  std::vector<float> reshapedInput(targetSize);
383 
384  for (int i = 0; i < targetSize; i++) {
385  interpolationFloor = static_cast<int>(std::floor(interpolationPoint));
386  // Interpolating here
387  reshapedInput.at(i) = (interpolationPoint - interpolationFloor) *
388  (inputs[interpolationFloor + 1] - inputs[interpolationFloor]) +
389  inputs[interpolationFloor];
390  interpolationPoint = step + interpolationPoint;
391  }
392  return reshapedInput;
393 }
394 
395 string DTOccupancyTestML::topFolder(bool isBooking) const {
396  if(tpMode) return string("DT/10-TestPulses/");
397  if(isBooking) return string("DT/01-Digi/ML");
398  return string("DT/01-Digi/");
399 }
size
Write out results.
Session * createSession(SessionOptions &sessionOptions)
Definition: TensorFlow.cc:87
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
Endjob.
T getUntrackedParameter(std::string const &, T const &) const
int getIntegral(TH2F *histo, int, int, int, int, bool)
double getFloatValue() const
MonitorElement * get(const std::string &path)
Definition: DQMStore.cc:302
GraphDef * loadGraphDef(const std::string &pbFile)
Definition: TensorFlow.cc:68
void bookHistos(DQMStore::IBooker &, const int wheelId, std::string folder, std::string histoTag)
book the summary histograms
std::vector< float > interpolateLayers(std::vector< float > const &inputs, int size, int targetSize)
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)
Definition: ntuple.py:1
static std::string const input
Definition: EdmProvDump.cc:44
void bookHistos()
Definition: Histogram.h:33
void dqmEndLuminosityBlock(DQMStore::IBooker &, DQMStore::IGetter &, edm::LuminosityBlock const &, edm::EventSetup const &) override
DQM Client Diagnostic.
std::string getMEName(std::string histoTag, const DTChamberId &chId)
Get the ME name.
bool closeSession(Session *&session)
Definition: TensorFlow.cc:193
void beginRun(edm::Run const &run, edm::EventSetup const &context) override
BeginRun.
std::string topFolder(bool isBooking) const
int runOccupancyTest(TH2F *histo, const DTChamberId &chId, float &chamberPercentage, tensorflow::GraphDef *graphDef, tensorflow::Session *session)
#define LogTrace(id)
TH2F * getTH2F() const
void setLogging(const std::string &level="3")
Definition: TensorFlow.cc:14
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
~DTOccupancyTestML() override
Destructor.
HLT enums.
int sector() const
Definition: DTChamberId.h:61
T get() const
Definition: EventSetup.h:63
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, const std::vector< std::string > &targetNodes, std::vector< Tensor > *outputs)
Definition: TensorFlow.cc:210
std::string fullPath() const
Definition: FileInPath.cc:197
step
int station() const
Return the station number.
Definition: DTChamberId.h:51
static char chambers[264][20]
Definition: ReadPGInfo.cc:243
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
DTOccupancyTestML(const edm::ParameterSet &ps)
Constructor.
Definition: Run.h:44