CMS 3D CMS Logo

MLClient.cc
Go to the documentation of this file.
2 
4 
6 
8 
10 
12 
13 #include <fstream>
14 #include <string>
15 
16 using namespace cms::Ort;
17 
18 namespace ecaldqm {
19  MLClient::MLClient() : DQWorkerClient() { qualitySummaries_.insert("MLQualitySummary"); }
20 
22  EBThreshold_ = _params.getUntrackedParameter<double>("EBThreshold");
23  EEpThreshold_ = _params.getUntrackedParameter<double>("EEpThreshold");
24  EEmThreshold_ = _params.getUntrackedParameter<double>("EEmThreshold");
25  EB_PUcorr_slope_ = _params.getUntrackedParameter<double>("EB_PUcorr_slope");
26  EB_PUcorr_intercept_ = _params.getUntrackedParameter<double>("EB_PUcorr_intercept");
27  EEp_PUcorr_slope_ = _params.getUntrackedParameter<double>("EEp_PUcorr_slope");
28  EEp_PUcorr_intercept_ = _params.getUntrackedParameter<double>("EEp_PUcorr_intercept");
29  EEm_PUcorr_slope_ = _params.getUntrackedParameter<double>("EEm_PUcorr_slope");
30  EEm_PUcorr_intercept_ = _params.getUntrackedParameter<double>("EEm_PUcorr_intercept");
31 
32  if (!onlineMode_) {
33  MEs_.erase(std::string("MLQualitySummary"));
34  MEs_.erase(std::string("EventsperMLImage"));
35  sources_.erase(std::string("PU"));
36  sources_.erase(std::string("NumEvents"));
37  sources_.erase(std::string("DigiAllByLumi"));
38  sources_.erase(std::string("AELoss"));
39  sources_.erase(std::string("BadTowerCount"));
40  sources_.erase(std::string("BadTowerCountNorm"));
41  }
42  }
43 
45  if (!onlineMode_)
46  return;
47  nbadtowerEB = 0;
48  nbadtowerEE = 0;
49 
50  using namespace std;
51  MESet& meMLQualitySummary(MEs_.at("MLQualitySummary"));
52  MESet& meEventsperMLImage(MEs_.at("EventsperMLImage"));
53 
54  MESetNonObject const& sPU(static_cast<MESetNonObject&>(sources_.at("PU")));
55  MESetNonObject const& sNumEvents(static_cast<MESetNonObject&>(sources_.at("NumEvents")));
56 
57  //Get the no.of events and the PU per LS calculated in OccupancyTask
58  int nEv = sNumEvents.getFloatValue();
59  double pu = sPU.getFloatValue();
60 
61  //Do not compute ML quality if PU is non existent.
62  if (pu <= 0.) {
63  return;
64  }
68 
70  //Inorder to feed the data into the ML model we apply some preprocessing.
71  //We use the Digi Occupancy per Lumisection as the input source.
72  //The model was trained on each occupancy plot having 500 events.
73  //In apprehension of the low luminosity in the beginning of Run3, where in online DQM
74  //the no.of events per LS could be lower than 500, we sum the occupancies over a fixed no.of lumisections as a running sum,
75  //and require that the total no.of events on this summed occupancy to be atleast 200.
76  //(This no.of LS and the no.of events are parameters which would require tuning later)
77  //This summed occupancy is now the input image, which is then corrected for PileUp(PU) dependence and
78  //change in no.of events, which are derived from training.
79  //The input image is also padded by replicating the top and bottom rows so as to prevent the "edge effect"
80  //wherein the ML model's learning degrades near the edge of the data set it sees.
81  //This padding is then removed during inference on the model output.
82 
83  //Get the histogram of the input digi occupancy per lumisection.
84  TH2F* hEEmDigiMap((sources_.at("DigiAllByLumi")).getME(0)->getTH2F());
85  TH2F* hEbDigiMap((sources_.at("DigiAllByLumi")).getME(1)->getTH2F());
86  TH2F* hEEpDigiMap((sources_.at("DigiAllByLumi")).getME(2)->getTH2F());
87 
88  size_t nEBTowers = nEBEtaTowers * nEBPhiTowers; //Each EB occupancy map is of size 34x72 towers
89  size_t nEETowers = nEEEtaTowers * nEEPhiTowers; //Each EE occupancy map is of size 20x20 towers
90 
91  //Vectors to feed into the ML network
92  std::vector<float> ebOccMap1dCumulPad;
93  std::vector<float> eemOccMap1dCumulPad;
94  std::vector<float> eepOccMap1dCumulPad;
95 
96  //Array to store occupancy maps
97  std::valarray<float> ebOccMap1d(nEBTowers);
98  std::valarray<float> eemOccMap1d(nEETowers);
99  std::valarray<float> eepOccMap1d(nEETowers);
100 
101  //Store the values from the input histogram into the array
102  //to do preprocessing
103  for (int i = 0; i < hEbDigiMap->GetNbinsY(); i++) { //NbinsY = 34, NbinsX = 72
104  for (int j = 0; j < hEbDigiMap->GetNbinsX(); j++) {
105  int bin = hEbDigiMap->GetBin(j + 1, i + 1);
106  int k = (i * nEBPhiTowers) + j;
107  ebOccMap1d[k] = hEbDigiMap->GetBinContent(bin);
108  }
109  }
110  ebOccMap1dQ.push_back(ebOccMap1d); //Queue which stores input occupancy maps for nLS lumis
111 
112  for (int i = 0; i < hEEpDigiMap->GetNbinsY(); i++) { //NbinsY = 20, NbinsX = 20
113  for (int j = 0; j < hEEpDigiMap->GetNbinsX(); j++) {
114  int bin = hEEpDigiMap->GetBin(j + 1, i + 1);
115  int k = (i * nEEPhiTowers) + j;
116  eemOccMap1d[k] = hEEmDigiMap->GetBinContent(bin);
117  eepOccMap1d[k] = hEEpDigiMap->GetBinContent(bin);
118  }
119  }
120 
121  //Queue which stores input occupancy maps for nLS lumis
122  eemOccMap1dQ.push_back(eemOccMap1d);
123  eepOccMap1dQ.push_back(eepOccMap1d);
124 
125  NEventQ.push_back(nEv); //Queue which stores the no.of events per LS for nLS lumis
126 
127  if (NEventQ.size() < nLS) {
128  return; //Should have nLS lumis to add the occupancy over.
129  }
130  if (NEventQ.size() > nLS) {
131  NEventQ.pop_front(); //Keep only nLS consecutive LS. Pop the first one if size greater than nLS
132  }
133  if (ebOccMap1dQ.size() > nLS) {
134  ebOccMap1dQ.pop_front(); //Same conditon for the input occupancy maps.
135  eemOccMap1dQ.pop_front();
136  eepOccMap1dQ.pop_front();
137  }
138 
139  int TNum = 0;
140  for (size_t i = 0; i < nLS; i++) {
141  TNum += NEventQ[i]; //Total no.of events over nLS lumis
142  }
143 
144  if (TNum < 400) {
145  return; //The total no.of events should be atleast 400 over nLS for meaningful statistics
146  }
147  //Fill the ME to monitor the trend of the total no.of events in each input image to the ML model
148  meEventsperMLImage.fill(getEcalDQMSetupObjects(), EcalBarrel, double(timestamp_.iLumi), double(TNum));
149 
150  //Array to hold the sum of inputs, which make atleast 400 events.
151  std::valarray<float> ebOccMap1dCumul(0., nEBTowers);
152  std::valarray<float> eemOccMap1dCumul(0., nEETowers);
153  std::valarray<float> eepOccMap1dCumul(0., nEETowers);
154 
155  //Sum the input arrays of nLS.
156  for (size_t i = 0; i < ebOccMap1dQ.size(); i++) {
157  ebOccMap1dCumul += ebOccMap1dQ[i];
158  eemOccMap1dCumul += eemOccMap1dQ[i];
159  eepOccMap1dCumul += eepOccMap1dQ[i];
160  }
161 
162  //Applying PU correction derived from training
163  ebOccMap1dCumul = ebOccMap1dCumul / (EB_PUcorr_slope_ * pu + EB_PUcorr_intercept_);
164  eemOccMap1dCumul = eemOccMap1dCumul / (EEm_PUcorr_slope_ * pu + EEm_PUcorr_intercept_);
165  eepOccMap1dCumul = eepOccMap1dCumul / (EEp_PUcorr_slope_ * pu + EEp_PUcorr_intercept_);
166 
167  //Scaling up to match input dimensions.
168  ebOccMap1dCumul = ebOccMap1dCumul * (nEBEtaTowers * nEBPhiTowers);
169  eemOccMap1dCumul = eemOccMap1dCumul * nEEEtaTowers * nEEPhiTowers; //(nEETowersPad * nEETowersPad);
170  eepOccMap1dCumul = eepOccMap1dCumul * nEEEtaTowers * nEEPhiTowers; //(nEETowersPad * nEETowersPad);
171 
172  //Correction for no.of events in each input image as originally model trained with 500 events per image
173  ebOccMap1dCumul = ebOccMap1dCumul * (500. / TNum);
174  eemOccMap1dCumul = eemOccMap1dCumul * (500. / TNum);
175  eepOccMap1dCumul = eepOccMap1dCumul * (500. / TNum);
176 
177  std::vector<std::vector<float>> ebOccMap2dCumul(nEBEtaTowers, std::vector<float>(nEBPhiTowers, 0.));
178  //Convert 1dCumul to 2d
179  for (size_t i = 0; i < nEBEtaTowers; i++) {
180  for (size_t j = 0; j < nEBPhiTowers; j++) {
181  int k = (i * nEBPhiTowers) + j;
182  ebOccMap2dCumul[i][j] = ebOccMap1dCumul[k];
183  }
184  }
185 
186  std::vector<float> pad_top;
187  std::vector<float> pad_bottom;
188  std::vector<float> pad_left;
189  std::vector<float> pad_right;
190 
191  pad_top = ebOccMap2dCumul[0];
192  pad_bottom = ebOccMap2dCumul[ebOccMap2dCumul.size() - 1];
193 
194  ebOccMap2dCumul.insert(ebOccMap2dCumul.begin(), pad_top);
195  ebOccMap2dCumul.push_back(pad_bottom);
196 
198  std::vector<std::vector<float>> eemOccMap2dCumul(nEEEtaTowers, std::vector<float>(nEEPhiTowers, 0.));
199  std::vector<std::vector<float>> eepOccMap2dCumul(nEEEtaTowers, std::vector<float>(nEEPhiTowers, 0.));
200 
201  for (size_t i = 0; i < nEEEtaTowers; i++) {
202  for (size_t j = 0; j < nEEPhiTowers; j++) {
203  int k = (i * nEEPhiTowers) + j;
204  eemOccMap2dCumul[i][j] = eemOccMap1dCumul[k];
205  eepOccMap2dCumul[i][j] = eepOccMap1dCumul[k];
206  }
207  }
208 
209  // EE - //
210  pad_top.clear();
211  pad_bottom.clear();
212  pad_left.clear();
213  pad_right.clear();
214 
215  pad_top = eemOccMap2dCumul[0];
216  pad_bottom = eemOccMap2dCumul[eemOccMap2dCumul.size() - 1];
217 
218  eemOccMap2dCumul.insert(eemOccMap2dCumul.begin(), pad_top);
219  eemOccMap2dCumul.push_back(pad_bottom);
220 
221  for (auto& row : eemOccMap2dCumul) {
222  pad_left.push_back(row[0]);
223  pad_right.push_back(row[row.size() - 1]);
224  }
225 
226  std::size_t Lindex = 0;
227  std::size_t Rindex = 0;
228 
229  for (auto& row : eemOccMap2dCumul) {
230  row.insert(row.begin(), pad_left[Lindex++]);
231  row.insert(row.end(), pad_right[Rindex++]);
232  }
233 
234  // EE + //
235  pad_top.clear();
236  pad_bottom.clear();
237 
238  pad_top = eepOccMap2dCumul[0];
239  pad_bottom = eepOccMap2dCumul[eepOccMap2dCumul.size() - 1];
240 
241  eepOccMap2dCumul.insert(eepOccMap2dCumul.begin(), pad_top);
242  eepOccMap2dCumul.push_back(pad_bottom);
243 
244  for (auto& row : eepOccMap2dCumul) {
245  pad_left.push_back(row[0]);
246  pad_right.push_back(row[row.size() - 1]);
247  }
248 
249  Lindex = 0;
250  Rindex = 0;
251 
252  for (auto& row : eepOccMap2dCumul) {
253  row.insert(row.begin(), pad_left[Lindex++]);
254  row.insert(row.end(), pad_right[Rindex++]);
255  }
256 
257  //The pre-processed input is now fed into the 1D input tensor vector which will go into the ML model
258  for (auto& row : ebOccMap2dCumul) {
259  ebOccMap1dCumulPad.insert(ebOccMap1dCumulPad.end(), row.begin(), row.end());
260  }
261 
262  for (auto& row : eemOccMap2dCumul) {
263  eemOccMap1dCumulPad.insert(eemOccMap1dCumulPad.end(), row.begin(), row.end());
264  }
265 
266  for (auto& row : eepOccMap2dCumul) {
267  eepOccMap1dCumulPad.insert(eepOccMap1dCumulPad.end(), row.begin(), row.end());
268  }
269 
271  //An Autoencoder (AE) network with resnet architecture is used here which is trained on
272  //certified good data (EB digi occupancy) from Run 2018 data.
273  //On giving an input occupancy map, the encoder part of the AE compresses and reduces the input data, learning its features,
274  //and the decoder reconstructs the data from the encoded form into a representation as close to the original input as possible.
275  //We then compute the Mean squared error (MSE) between the input and output image, also called the Reconstruction Loss,
276  //calculated at a tower by tower basis.
277  //Thus, given an anomalous tower the loss should be significantly higher than the loss with respect to good towers, which the model
278  //has already seen --> anomaly detection.
279  //When calculating the loss we also apply a response correction by dividing each input and output image with the average occupancy from
280  //all 2018 data (also to be tuned),to accommodate the difference in response of crystals in different regions of the Ecal Barrel
281  //Further each loss map from each input image is then multiplied by the last N loss maps,
283  //A quality threshold is then applied on this time multiplied loss map, to mark them as GOOD or BAD,
284  //after which it is stored as a quality summary ME.
285 
287  std::string instanceName{"AE-DQM-inference"};
288  std::string modelFilepath = edm::FileInPath("DQM/EcalMonitorClient/data/onnxModels/resnet.onnx").fullPath();
289 
290  Ort::SessionOptions sessionOptions;
291  sessionOptions.SetIntraOpNumThreads(1);
292  Ort::Env env(OrtLoggingLevel::ORT_LOGGING_LEVEL_WARNING, instanceName.c_str());
293  Ort::Session session(env, modelFilepath.c_str(), sessionOptions);
294 
295  Ort::AllocatorWithDefaultOptions allocator;
296 
297  // Strings returned by session.GetInputNameAllocated are temporary, need to copy them before they are deallocated
298  std::string inputName{session.GetInputNameAllocated(0, allocator).get()};
299 
300  Ort::TypeInfo inputTypeInfo = session.GetInputTypeInfo(0);
301  auto inputTensorInfo = inputTypeInfo.GetTensorTypeAndShapeInfo();
302 
303  std::vector<int64_t> inputDims = inputTensorInfo.GetShape();
304 
305  std::string outputName{session.GetOutputNameAllocated(0, allocator).get()};
306 
307  Ort::TypeInfo outputTypeInfo = session.GetOutputTypeInfo(0);
308  auto outputTensorInfo = outputTypeInfo.GetTensorTypeAndShapeInfo();
309 
310  std::vector<int64_t> outputDims = outputTensorInfo.GetShape();
311 
312  size_t TensorSize = nEBEtaTowersPad * nEBPhiTowers;
313  std::vector<float> ebRecoOccMap1dPad(TensorSize); //To store the output reconstructed occupancy
314 
315  std::vector<const char*> inputNames{inputName.c_str()};
316  std::vector<const char*> outputNames{outputName.c_str()};
317  std::vector<Ort::Value> inputTensors;
318  std::vector<Ort::Value> outputTensors;
319 
320  Ort::MemoryInfo memoryInfo =
321  Ort::MemoryInfo::CreateCpu(OrtAllocatorType::OrtArenaAllocator, OrtMemType::OrtMemTypeDefault);
322  inputTensors.push_back(Ort::Value::CreateTensor<float>(
323  memoryInfo, ebOccMap1dCumulPad.data(), TensorSize, inputDims.data(), inputDims.size()));
324 
325  outputTensors.push_back(Ort::Value::CreateTensor<float>(
326  memoryInfo, ebRecoOccMap1dPad.data(), TensorSize, outputDims.data(), outputDims.size()));
327 
328  session.Run(Ort::RunOptions{nullptr},
329  inputNames.data(),
330  inputTensors.data(),
331  1,
332  outputNames.data(),
333  outputTensors.data(),
334  1);
335 
336  //Endcaps
337  // EE- //
338 
339  inputDims.clear();
340  outputDims.clear();
341  inputNames.clear();
342  outputNames.clear();
343  inputTensors.clear();
344  outputTensors.clear();
345 
346  modelFilepath = edm::FileInPath("DQM/EcalMonitorClient/data/onnxModels/EEm_resnet2018.onnx").fullPath();
347 
348  Ort::Session EEm_session(env, modelFilepath.c_str(), sessionOptions);
349 
350  inputName = EEm_session.GetInputNameAllocated(0, allocator).get();
351 
352  inputTypeInfo = EEm_session.GetInputTypeInfo(0);
353  auto EEm_inputTensorInfo = inputTypeInfo.GetTensorTypeAndShapeInfo();
354 
355  inputDims = EEm_inputTensorInfo.GetShape();
356 
357  outputName = EEm_session.GetOutputNameAllocated(0, allocator).get();
358 
359  //Ort::TypeInfo
360  outputTypeInfo = EEm_session.GetOutputTypeInfo(0);
361  auto EEm_outputTensorInfo = outputTypeInfo.GetTensorTypeAndShapeInfo();
362 
363  outputDims = EEm_outputTensorInfo.GetShape();
364 
365  size_t EE_TensorSize = nEETowersPad * nEETowersPad;
366  std::vector<float> eemRecoOccMap1dPad(EE_TensorSize); //To store the output reconstructed occupancy
367 
368  inputNames.push_back(inputName.c_str());
369  outputNames.push_back(outputName.c_str());
370 
371  //Ort::MemoryInfo
372  memoryInfo = Ort::MemoryInfo::CreateCpu(OrtAllocatorType::OrtArenaAllocator, OrtMemType::OrtMemTypeDefault);
373  inputTensors.push_back(Ort::Value::CreateTensor<float>(
374  memoryInfo, eemOccMap1dCumulPad.data(), EE_TensorSize, inputDims.data(), inputDims.size()));
375 
376  outputTensors.push_back(Ort::Value::CreateTensor<float>(
377  memoryInfo, eemRecoOccMap1dPad.data(), EE_TensorSize, outputDims.data(), outputDims.size()));
378 
379  EEm_session.Run(Ort::RunOptions{nullptr},
380  inputNames.data(),
381  inputTensors.data(),
382  1,
383  outputNames.data(),
384  outputTensors.data(),
385  1);
386 
387  // EE+ //
388  inputDims.clear();
389  outputDims.clear();
390  inputNames.clear();
391  outputNames.clear();
392  inputTensors.clear();
393  outputTensors.clear();
394 
395  modelFilepath = edm::FileInPath("DQM/EcalMonitorClient/data/onnxModels/EEp_resnet2018.onnx").fullPath();
396 
397  Ort::Session EEp_session(env, modelFilepath.c_str(), sessionOptions);
398 
399  inputName = EEp_session.GetInputNameAllocated(0, allocator).get();
400 
401  inputTypeInfo = EEp_session.GetInputTypeInfo(0);
402  auto EEp_inputTensorInfo = inputTypeInfo.GetTensorTypeAndShapeInfo();
403 
404  inputDims = EEp_inputTensorInfo.GetShape();
405 
406  outputName = EEp_session.GetOutputNameAllocated(0, allocator).get();
407 
408  outputTypeInfo = EEp_session.GetOutputTypeInfo(0);
409 
410  auto EEp_outputTensorInfo = outputTypeInfo.GetTensorTypeAndShapeInfo();
411 
412  outputDims = EEp_outputTensorInfo.GetShape();
413 
414  std::vector<float> eepRecoOccMap1dPad(EE_TensorSize); //To store the output reconstructed occupancy
415 
416  inputNames.push_back(inputName.c_str());
417  outputNames.push_back(outputName.c_str());
418 
419  //Ort::MemoryInfo
420  memoryInfo = Ort::MemoryInfo::CreateCpu(OrtAllocatorType::OrtArenaAllocator, OrtMemType::OrtMemTypeDefault);
421  inputTensors.push_back(Ort::Value::CreateTensor<float>(
422  memoryInfo, eepOccMap1dCumulPad.data(), EE_TensorSize, inputDims.data(), inputDims.size()));
423 
424  outputTensors.push_back(Ort::Value::CreateTensor<float>(
425  memoryInfo, eepRecoOccMap1dPad.data(), EE_TensorSize, outputDims.data(), outputDims.size()));
426 
427  EEp_session.Run(Ort::RunOptions{nullptr},
428  inputNames.data(),
429  inputTensors.data(),
430  1,
431  outputNames.data(),
432  outputTensors.data(),
433  1);
434 
436  //2D Loss map to store tower by tower loss between the output (reconstructed) and input occupancies,
437  //Have same dimensions as the occupancy plot
438  std::valarray<std::valarray<float>> EBlossMap2d(std::valarray<float>(nEBPhiTowers), nEBEtaTowers);
439  std::valarray<std::valarray<float>> EEmlossMap2d(std::valarray<float>(nEEPhiTowers), nEEEtaTowers);
440  std::valarray<std::valarray<float>> EEplossMap2d(std::valarray<float>(nEEPhiTowers), nEEEtaTowers);
441 
442  //1D val arrays to store row wise information corresponding to the reconstructed, input and average occupancies, and loss.
443  //and to do element wise (tower wise) operations on them to calculate the MSE loss between the reco and input occupancy.
444  std::valarray<float> EBrecoOcc1d(0., nEBPhiTowers);
445  std::valarray<float> EBinputOcc1d(0., nEBPhiTowers);
446  std::valarray<float> EBavgOcc1d(0., nEBPhiTowers);
447  std::valarray<float> EBloss_;
448 
449  std::valarray<float> EEmrecoOcc1d(0., nEEPhiTowers);
450  std::valarray<float> EEminputOcc1d(0., nEEPhiTowers);
451  std::valarray<float> EEmavgOcc1d(0., nEEPhiTowers);
452  std::valarray<float> EEmloss_;
453 
454  std::valarray<float> EEprecoOcc1d(0., nEEPhiTowers);
455  std::valarray<float> EEpinputOcc1d(0., nEEPhiTowers);
456  std::valarray<float> EEpavgOcc1d(0., nEEPhiTowers);
457  std::valarray<float> EEploss_;
458 
459  std::string EBOccpath =
460  edm::FileInPath("DQM/EcalMonitorClient/data/MLAvgOccupancy/EB_avgocc_Run2022_500ev.dat").fullPath();
461  std::ifstream inFile;
462  double val;
463  inFile.open((EBOccpath).c_str());
464  while (inFile) {
465  inFile >> val;
466  if (inFile.eof())
467  break;
468  EBavgOcc.push_back(val);
469  }
470  inFile.close();
471 
472  std::string EEmOccpath =
473  edm::FileInPath("DQM/EcalMonitorClient/data/MLAvgOccupancy/EEm_avgocc_Run2022_500ev.dat").fullPath();
474  inFile.open((EEmOccpath).c_str());
475  while (inFile) {
476  inFile >> val;
477  if (inFile.eof())
478  break;
479  EEmavgOcc.push_back(val);
480  }
481  inFile.close();
482 
483  std::string EEpOccpath =
484  edm::FileInPath("DQM/EcalMonitorClient/data/MLAvgOccupancy/EEp_avgocc_Run2022_500ev.dat").fullPath();
485  inFile.open((EEpOccpath).c_str());
486  while (inFile) {
487  inFile >> val;
488  if (inFile.eof())
489  break;
490  EEpavgOcc.push_back(val);
491  }
492  inFile.close();
493 
494  //Loss calculation
495  //Ignore the top and bottom replicated padded rows when doing inference
496  //by making index i run over (1,35) instead of (0,36) for EB, and over (1,21) for EE
497 
498  MESet const& sAEReco(sources_.at("AEReco"));
499  TH2F* hEBRecoMap2d(sAEReco.getME(1)->getTH2F());
500 
501  for (int i = 1; i < nEBEtaTowersPad - 1; i++) {
502  for (int j = 0; j < nEBPhiTowers; j++) {
503  int k = (i * nEBPhiTowers) + j;
504  int bin_ = hEBRecoMap2d->GetBin(j + 1, i);
505  EBrecoOcc1d[j] = ebRecoOccMap1dPad[k];
506  EBinputOcc1d[j] = ebOccMap1dCumulPad[k];
507  EBavgOcc1d[j] = EBavgOcc[k];
508  double content = ebRecoOccMap1dPad[k];
509  hEBRecoMap2d->SetBinContent(bin_, content);
510  }
511  //Calculate the MSE loss = (output-input)^2, with avg response correction
512  EBloss_ = std::pow((EBrecoOcc1d / EBavgOcc1d - EBinputOcc1d / EBavgOcc1d), 2);
513  EBlossMap2d[i - 1] = (EBloss_);
514  }
515 
516  TH2F* hEEmRecoMap2d(sAEReco.getME(0)->getTH2F());
517  TH2F* hEEpRecoMap2d(sAEReco.getME(2)->getTH2F());
518 
519  for (int i = 1; i < nEETowersPad - 1; i++) {
520  for (int j = 0; j < nEEPhiTowers; j++) {
521  int k = (i * nEETowersPad) + j + 1;
522  int bin_ = hEEmRecoMap2d->GetBin(j + 1, i);
523 
524  EEmrecoOcc1d[j] = eemRecoOccMap1dPad[k];
525  EEminputOcc1d[j] = eemOccMap1dCumulPad[k];
526  EEmavgOcc1d[j] = EEmavgOcc[k];
527  double EEmcontent = eemRecoOccMap1dPad[k];
528  hEEmRecoMap2d->SetBinContent(bin_, EEmcontent);
529 
530  EEprecoOcc1d[j] = eepRecoOccMap1dPad[k];
531  EEpinputOcc1d[j] = eepOccMap1dCumulPad[k];
532  EEpavgOcc1d[j] = EEpavgOcc[k];
533  double EEpcontent = eepRecoOccMap1dPad[k];
534  hEEpRecoMap2d->SetBinContent(bin_, EEpcontent);
535  }
536  //Calculate the MSE loss = (output-input)^2, with avg response correction
537  EEmloss_ = std::pow((EEmrecoOcc1d / EEmavgOcc1d - EEminputOcc1d / EEmavgOcc1d), 2);
538  EEmlossMap2d[i - 1] = (EEmloss_);
539 
540  EEploss_ = std::pow((EEprecoOcc1d / EEpavgOcc1d - EEpinputOcc1d / EEpavgOcc1d), 2);
541  EEplossMap2d[i - 1] = (EEploss_);
542  }
543 
544  //Store each loss map from the output in the queue
545  EBlossMap2dQ.push_back(EBlossMap2d);
546  EEmlossMap2dQ.push_back(EEmlossMap2d);
547  EEplossMap2dQ.push_back(EEplossMap2d);
548 
549  //Keep exactly nLSloss loss maps to multiply
550  if (EBlossMap2dQ.size() > nLSloss) {
551  EBlossMap2dQ.pop_front();
552  EEmlossMap2dQ.pop_front();
553  EEplossMap2dQ.pop_front();
554  }
555  if (EBlossMap2dQ.size() < nLSloss) { //Exit if there are not nLSloss loss maps
556  return;
557  }
558 
559  //To hold the final multiplied loss
560  std::valarray<std::valarray<float>> EBlossMap2dMult(std::valarray<float>(1., nEBPhiTowers), nEBEtaTowers);
561  std::valarray<std::valarray<float>> EEmlossMap2dMult(std::valarray<float>(1., nEEPhiTowers), nEEEtaTowers);
562  std::valarray<std::valarray<float>> EEplossMap2dMult(std::valarray<float>(1., nEEPhiTowers), nEEEtaTowers);
563 
564  //Multiply together the last nLSloss loss maps
565  //So that real anomalies which persist with time are enhanced and fluctuations are suppressed.
566  for (size_t i = 0; i < EBlossMap2dQ.size(); i++) {
567  EBlossMap2dMult *= EBlossMap2dQ[i];
568  EEmlossMap2dMult *= EEmlossMap2dQ[i];
569  EEplossMap2dMult *= EEplossMap2dQ[i];
570  }
571 
572  //Fill the AELoss ME with the values of this time multiplied loss map
573  //MESet const& sAELoss(sources_.at("AELoss"));
574  MESet& sAELoss(sources_.at("AELoss"));
575 
576  TH2F* hEBLossMap2dMult(sAELoss.getME(1)->getTH2F());
577 
578  for (int i = 0; i < hEBLossMap2dMult->GetNbinsY(); i++) {
579  for (int j = 0; j < hEBLossMap2dMult->GetNbinsX(); j++) {
580  int bin_ = hEBLossMap2dMult->GetBin(j + 1, i + 1);
581  double content = EBlossMap2dMult[i][j];
582  hEBLossMap2dMult->SetBinContent(bin_, content);
583  }
584  }
585 
586  TH2F* hEEmLossMap2dMult(sAELoss.getME(0)->getTH2F());
587  TH2F* hEEpLossMap2dMult(sAELoss.getME(2)->getTH2F());
588 
589  for (int i = 0; i < hEEmLossMap2dMult->GetNbinsY(); i++) {
590  for (int j = 0; j < hEEmLossMap2dMult->GetNbinsX(); j++) {
591  int bin_ = hEEmLossMap2dMult->GetBin(j + 1, i + 1);
592 
593  double EEmcontent = EEmlossMap2dMult[i][j];
594  hEEmLossMap2dMult->SetBinContent(bin_, EEmcontent);
595 
596  double EEpcontent = EEplossMap2dMult[i][j];
597  hEEpLossMap2dMult->SetBinContent(bin_, EEpcontent);
598  }
599  }
600 
602  //Apply the quality threshold on the time multiplied loss map stored in the ME AELoss
603  //If anomalous, the tower entry will have a large loss value. If good, the value will be close to zero.
604 
605  MESet& meBadTowerCount(sources_.at("BadTowerCount"));
606  MESet& meBadTowerCountNorm(sources_.at("BadTowerCountNorm"));
607  MESet& meTrendMLBadTower(MEs_.at("TrendMLBadTower"));
608 
609  LScount++;
610 
611  MESet::const_iterator dAEnd(sAELoss.end(GetElectronicsMap()));
612  for (MESet::const_iterator dItr(sAELoss.beginChannel(GetElectronicsMap())); dItr != dAEnd;
614  DetId id(dItr->getId());
615 
616  bool doMaskML(meMLQualitySummary.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
617 
618  float entries(dItr->getBinContent());
619 
620  int quality(doMaskML ? kMGood : kGood);
621  float MLThreshold;
622 
623  if (id.subdetId() == EcalEndcap) {
624  EEDetId eeid(id);
625  if (eeid.zside() > 0)
626  MLThreshold = EEpThreshold_;
627  else
628  MLThreshold = EEmThreshold_;
629  } else {
630  MLThreshold = EBThreshold_;
631  }
632 
633  //If a trigger tower entry is greater than the ML threshold, set it to Bad quality, otherwise Good.
634  if (entries > MLThreshold) {
635  quality = doMaskML ? kMBad : kBad;
636  meBadTowerCount.fill(getEcalDQMSetupObjects(), id);
637  if (id.subdetId() == EcalEndcap)
638  nbadtowerEE++;
639  else
640  nbadtowerEB++;
641  }
642  //Fill the quality summary with the quality of the given tower id.
643  meMLQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, double(quality));
644 
645  double badtowcount(meBadTowerCount.getBinContent(getEcalDQMSetupObjects(), id));
646  meBadTowerCountNorm.setBinContent(getEcalDQMSetupObjects(), id, double(badtowcount / LScount));
647  } // ML Quality Summary
648 
649  meTrendMLBadTower.fill(getEcalDQMSetupObjects(), EcalBarrel, double(timestamp_.iLumi), double(nbadtowerEB));
650  meTrendMLBadTower.fill(getEcalDQMSetupObjects(), EcalEndcap, double(timestamp_.iLumi), double(nbadtowerEE));
651 
652  } // producePlots()
653 
655 } // namespace ecaldqm
std::deque< std::valarray< std::valarray< float > > > EBlossMap2dQ
Definition: MLClient.h:60
static const int nEEPhiTowers
Definition: MLClient.h:27
float EEm_PUcorr_slope_
Definition: MLClient.h:38
void setParams(edm::ParameterSet const &) override
Definition: MLClient.cc:21
#define DEFINE_ECALDQM_WORKER(TYPE)
Definition: DQWorker.h:168
static const int nEBEtaTowers
Definition: MLClient.h:23
edm::LuminosityBlockNumber_t iLumi
Definition: DQWorker.h:48
static const int PHYSICS_BAD_CHANNEL_WARNING
std::string fullPath() const
Definition: FileInPath.cc:161
MESet & at(const std::string &key)
Definition: MESet.h:399
std::deque< std::valarray< std::valarray< float > > > EEplossMap2dQ
Definition: MLClient.h:61
static const int PEDESTAL_ONLINE_HIGH_GAIN_RMS_ERROR
static const int PHYSICS_BAD_CHANNEL_ERROR
void producePlots(ProcessType) override
Definition: MLClient.cc:44
static const int nEBEtaTowersPad
Definition: MLClient.h:31
string quality
double getFloatValue() const
std::set< std::string > qualitySummaries_
std::vector< double > EEpavgOcc
Definition: MLClient.h:57
float EEpThreshold_
Definition: MLClient.h:34
std::deque< std::valarray< float > > eemOccMap1dQ
Definition: MLClient.h:54
StatusManager const * statusManager_
float EEp_PUcorr_intercept_
Definition: MLClient.h:40
size_t nLSloss
Definition: MLClient.h:44
std::deque< std::valarray< float > > eepOccMap1dQ
Definition: MLClient.h:53
std::deque< std::valarray< std::valarray< float > > > EEmlossMap2dQ
Definition: MLClient.h:62
static const int nEETowersPad
Definition: MLClient.h:32
MESetCollection sources_
float EB_PUcorr_slope_
Definition: MLClient.h:36
float EEm_PUcorr_intercept_
Definition: MLClient.h:41
EcalElectronicsMapping const * GetElectronicsMap()
Definition: DQWorker.cc:150
EcalDQMSetupObjects const getEcalDQMSetupObjects()
Definition: DQWorker.cc:170
Definition: DetId.h:17
Timestamp timestamp_
Definition: DQWorker.h:134
int zside() const
Definition: EEDetId.h:71
std::vector< double > EEmavgOcc
Definition: MLClient.h:58
float EEmThreshold_
Definition: MLClient.h:35
static const int nEBPhiTowers
Definition: MLClient.h:24
MESetCollection MEs_
Definition: DQWorker.h:131
std::deque< int > NEventQ
Definition: MLClient.h:49
float EB_PUcorr_intercept_
Definition: MLClient.h:39
const_iterator & toNextChannel(EcalElectronicsMapping const *)
Definition: MESet.cc:413
std::vector< double > EBavgOcc
Definition: MLClient.h:56
EcalTrigTowerConstituentsMap const * GetTrigTowerMap()
Definition: DQWorker.cc:155
float EBThreshold_
Definition: MLClient.h:33
float EEp_PUcorr_slope_
Definition: MLClient.h:37
static const int nEEEtaTowers
Definition: MLClient.h:26
void erase(const std::string &key)
Definition: MESet.h:390
std::deque< std::valarray< float > > ebOccMap1dQ
Definition: MLClient.h:52
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29