Belle II Software  release-05-01-25
eclDQMExtended.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * ECL Data Quality Monitor (Second Module) *
6  * *
7  * This module provides histograms to check out ECL electronics logic *
8  * *
9  * Author: The Belle II Collaboration *
10  * Contributors: Dmitry Matvienko (d.v.matvienko@inp.nsk.su) *
11  * *
12  * This software is provided "as is" without any warranty. *
13  **************************************************************************/
14 
15 //THIS MODULE
16 #include <ecl/modules/eclDQMExtended/eclDQMExtended.h>
17 
18 //FRAMEWORK
19 #include <framework/core/HistoModule.h>
20 #include <framework/datastore/StoreArray.h>
21 #include <framework/logging/Logger.h>
22 
23 //ECL
24 #include <ecl/dataobjects/ECLDigit.h>
25 #include <ecl/dataobjects/ECLTrig.h>
26 #include <ecl/dataobjects/ECLDsp.h>
27 #include <ecl/utility/ECLDspUtilities.h>
28 #include <ecl/utility/ECLDspEmulator.h>
29 
30 //STL
31 #include <regex>
32 #include <map>
33 #include <vector>
34 #include <string>
35 #include <iostream>
36 
37 //ROOT
38 #include <TH1F.h>
39 #include <TH2F.h>
40 #include <TDirectory.h>
41 
42 //BOOST
43 #include <boost/filesystem.hpp>
44 #include <boost/format.hpp>
45 
46 
47 //NAMESPACE(S)
48 using namespace Belle2;
49 using namespace ECL;
50 
51 //-----------------------------------------------------------------
52 // Register the Module
53 //-----------------------------------------------------------------
54 REG_MODULE(ECLDQMEXTENDED)
55 
56 
57 //-----------------------------------------------------------------
58 // Implementation
59 //-----------------------------------------------------------------
60 
62  : HistoModule(),
63  m_ECLDspDataArray0("ECLDSPPars0"),
64  m_ECLDspDataArray1("ECLDSPPars1"),
65  m_ECLDspDataArray2("ECLDSPPars2"),
66  m_calibrationThrA0("ECL_FPGA_LowAmp"),
67  m_calibrationThrAhard("ECL_FPGA_HitThresh"),
68  m_calibrationThrAskip("ECL_FPGA_StoreDigit")
69 
70 {
71 
72  //Set module properties
73  setDescription("ECL Data Quality Monitor. Logic Test");
74  setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
75  addParam("histogramDirectoryName", m_histogramDirectoryName,
76  "histogram directory in ROOT file", std::string("ECL"));
77  addParam("InitKey", m_InitKey,
78  "How to initialize DSP coeffs: ''DB'' or ''File'' are acceptable ", std::string("DB"));
79  addParam("DSPDirectoryName", m_DSPDirectoryName,
80  "directory for DSP coeffs", std::string("/hsm/belle2/bdata/users/dmitry/dsp/"));
81  addParam("RunName", m_RunName,
82  "Name of run with DSP files", std::string("run0000/"));
83  addParam("SaveDetailedFitData", m_SaveDetailedFitData, "Save detailed data "
84  "(ampdiff_{cellid,shaper}, timediff_{cellid,shaper} histograms) for "
85  "failed fits.", false);
86  addParam("AdjustedTiming", m_adjusted_timing, "Use improved procedure for "
87  "time determination in ShaperDSP emulator", true);
88 }
89 
91 {
92 }
93 
94 //------------------------------------------------------------------
95 // Function to define histograms
96 //-----------------------------------------------------------------
97 
99 {
100  TDirectory* oldDir = gDirectory;
101 
102  // Create a separate histogram directory and cd into it.
103 
104  TDirectory* dirDAQ = dynamic_cast<TDirectory*>(oldDir->Get(m_histogramDirectoryName.c_str()));
105  if (!dirDAQ) dirDAQ = oldDir->mkdir(m_histogramDirectoryName.c_str());
106  dirDAQ->cd();
107 
108  //1D histograms creation.
109 
110  std::vector<TH1F*> h_amp_qualityfail_raw, h_time_qualityfail_raw;
111 
112  for (int i = 0; i < 4; i++) {
113  std::string h_name_a, h_title_a, h_name_t, h_title_t;
114  h_name_a = str(boost::format("amp_timefail_q%1%") % i);
115  h_title_a = str(boost::format("Amp for time mismatches w/ FPGA fit qual=%1%") % i);
116  h_name_t = str(boost::format("time_ampfail_q%1%") % i);
117  h_title_t = str(boost::format("Time for amp mismatches w/ FPGA fit qual=%1%") % i);
118  TH1F* h_a = new TH1F(h_name_a.c_str(), h_title_a.c_str(), 1200, 0, 262144);
119  TH1F* h_t = new TH1F(h_name_t.c_str(), h_title_t.c_str(), 240, -2050, 2050);
120  h_a->SetOption("LIVE");
121  h_t->SetOption("LIVE");
122  h_amp_timefail.push_back(h_a);
123  h_time_ampfail.push_back(h_t);
124  }
125 
126  for (int i = 0; i < 4; i++) {
127  for (int j = 0; j < 4; j++) {
128  if (j != i) {
129  std::string h_name_a, h_title_a, h_name_t, h_title_t;
130  h_name_a = str(boost::format("amp_qf%1%_qd%2%") % i % j);
131  h_title_a = str(boost::format("Amp for C++ fit qual=%1% and FPGA fit qual=%2%") % i % j);
132  h_name_t = str(boost::format("time_qf%1%_qd%2%") % i % j);
133  h_title_t = str(boost::format("Time for C++ fit qual=%1% and FPGA fit qual=%2%") % i % j);
134  TH1F* h_a = new TH1F(h_name_a.c_str(), h_title_a.c_str(), 1200, 0, 262144);
135  TH1F* h_t = new TH1F(h_name_t.c_str(), h_title_t.c_str(), 240, -2050, 2050);
136  h_a->SetOption("LIVE");
137  h_t->SetOption("LIVE");
138  h_amp_qualityfail_raw.push_back(h_a);
139  h_time_qualityfail_raw.push_back(h_t);
140  }
141  }
142  h_amp_qualityfail.push_back(h_amp_qualityfail_raw);
143  h_time_qualityfail.push_back(h_time_qualityfail_raw);
144  h_amp_qualityfail_raw.clear();
145  h_time_qualityfail_raw.clear();
146  }
147 
148  h_ampfail_quality = new TH1F("ampfail_quality", "Number of FPGA <-> C++ fitter #bf{amp} inconsistencies vs fit qual", 5, -1, 4);
149  h_ampfail_quality->GetXaxis()->SetTitle("FPGA fit qual. -1-all evts,0-good,1-int overflow,2-low amp,3-bad chi2");
150  h_ampfail_quality->SetFillColor(kPink - 4);
151  h_ampfail_quality->SetOption("LIVE");
152 
153  h_timefail_quality = new TH1F("timefail_quality", "Number of FPGA <-> C++ fitter #bf{time} inconsistencies vs fit qual", 5, -1, 4);
154  h_timefail_quality->SetFillColor(kPink - 4);
155  h_timefail_quality->GetXaxis()->SetTitle("FPGA fit qual. -1-all evts,0-good,1-int overflow,2-low amp,3-bad chi2");
156  h_timefail_quality->SetOption("LIVE");
157 
158  h_ampfail_cellid = new TH1F("ampfail_cellid", "Cell IDs w/ amp inconsistencies", 8736, 1, 8737);
159  h_ampfail_cellid->GetXaxis()->SetTitle("Cell ID");
160  h_ampfail_cellid->SetOption("LIVE");
161 
162  h_timefail_cellid = new TH1F("timefail_cellid", "Cell IDs w/ time inconsistencies", 8736, 1, 8737);
163  h_timefail_cellid->GetXaxis()->SetTitle("Cell ID");
164  h_timefail_cellid->SetOption("LIVE");
165 
166  h_amptimefail_cellid = new TH1F("amptimefail_cellid", "Cell IDs w/ time and amp inconsistencies", 8736, 1, 8737);
167  h_amptimefail_cellid->GetXaxis()->SetTitle("Cell ID");
168  h_amptimefail_cellid->SetOption("LIVE");
169 
170  h_ampfail_shaperid = new TH1F("ampfail_shaperid", "Shaper IDs w/ amp inconsistencies", 624, 1, 625);
171  h_ampfail_shaperid->GetXaxis()->SetTitle("Shaper ID");
172  h_ampfail_shaperid->SetOption("LIVE");
173 
174  h_timefail_shaperid = new TH1F("timefail_shaperid", "Shaper IDs w/ time inconsistencies", 624, 1, 625);
175  h_timefail_shaperid->GetXaxis()->SetTitle("Shaper ID");
176  h_timefail_shaperid->SetOption("LIVE");
177 
178  h_amptimefail_shaperid = new TH1F("amptimefail_shaperid", "Shaper IDs w/ time and amp inconsistencies", 624, 1, 625);
179  h_amptimefail_shaperid->GetXaxis()->SetTitle("Shaper ID");
180  h_amptimefail_shaperid->SetOption("LIVE");
181 
182  h_ampfail_crateid = new TH1F("ampfail_crateid", "Crate IDs w/ amp inconsistencies", 52, 1, 53);
183  h_ampfail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
184  h_ampfail_crateid->SetOption("LIVE");
185 
186  h_timefail_crateid = new TH1F("timefail_crateid", "Crate IDs w/ time inconsistencies", 52, 1, 53);
187  h_timefail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
188  h_timefail_crateid->SetOption("LIVE");
189 
190  h_amptimefail_crateid = new TH1F("amptimefail_crateid", "Crate IDs w/ time and amp inconsistencies", 52, 1, 53);
191  h_amptimefail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
192  h_amptimefail_crateid->SetOption("LIVE");
193 
194  h_qualityfail_cellid = new TH1F("qualityfail_cellid", "Cell IDs w/ fit qual inconsistencies", 8736, 1, 8737);
195  h_qualityfail_cellid->GetXaxis()->SetTitle("Cell ID");
196  h_qualityfail_cellid->SetOption("LIVE");
197 
198  h_qualityfail_shaperid = new TH1F("qualityfail_shaperid", "Shaper IDs w/ fit qual inconsistencies", 624, 1, 625);
199  h_qualityfail_shaperid->GetXaxis()->SetTitle("Shaper ID");
200  h_qualityfail_shaperid->SetOption("LIVE");
201 
202  h_qualityfail_crateid = new TH1F("qualityfail_crateid", "Crate IDs w/ fit qual inconsistencies", 52, 1, 53);
203  h_qualityfail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
204  h_qualityfail_crateid->SetOption("LIVE");
205 
206  h_fail_shaperid = new TH1F("fail_shaperid", "Shaper IDs w/ inconsistencies", 624, 1, 625);
207  h_fail_shaperid->GetXaxis()->SetTitle("Shaper ID");
208  h_fail_shaperid->SetOption("LIVE");
209 
210  h_fail_crateid = new TH1F("fail_crateid", "Crate IDs w/ inconsistencies", 52, 1, 53);
211  h_fail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
212  h_fail_crateid->SetOption("LIVE");
213 
214 
215  //2D histograms creation.
216 
217  if (m_SaveDetailedFitData) {
218  h_ampdiff_cellid = new TH2F("ampdiff_cellid", "Amp. diff. (Emulator-Data) for amp inconsistencies",
219  8736, 1, 8737, 239, -262143, 262143);
220  h_ampdiff_cellid->GetXaxis()->SetTitle("Cell ID");
221  h_ampdiff_cellid->GetYaxis()->SetTitle("Amplitude difference");
222  h_ampdiff_cellid->SetOption("LIVE");
223 
224  h_timediff_cellid = new TH2F("timediff_cellid", "Time diff. (Emulator-Data) for time inconsistencies",
225  8736, 1, 8737, 239, -4095, 4095);
226  h_timediff_cellid->GetXaxis()->SetTitle("Cell ID");
227  h_timediff_cellid->GetYaxis()->SetTitle("Time difference");
228  h_timediff_cellid->SetOption("LIVE");
229 
230  h_ampdiff_shaperid = new TH2F("ampdiff_shaper", "Amp. diff. (Emulator-Data) "
231  "for amp inconsistencies vs Shaper Id",
232  624, 1, 625, 239, -262143, 262143);
233  h_ampdiff_shaperid->GetXaxis()->SetTitle("Shaper Id");
234  h_ampdiff_shaperid->GetYaxis()->SetTitle("Amplitude difference");
235  h_ampdiff_shaperid->SetOption("LIVE");
236 
237  h_timediff_shaperid = new TH2F("timediff_shaper", "Time diff. (Emulator-Data) "
238  "for time inconsistencies vs Shaper Id",
239  624, 1, 625, 239, -4095, 4095);
240  h_timediff_shaperid->GetXaxis()->SetTitle("Shaper Id");
241  h_timediff_shaperid->GetYaxis()->SetTitle("Time difference");
242  h_timediff_shaperid->SetOption("LIVE");
243  }
244 
245  h_ampdiff_quality = new TH2F("ampdiff_quality", "Amp. diff. (Emulator-Data) for amp. inconsistencies", 4, 0, 4, 239,
246  -262143, 262143);
247  h_ampdiff_quality->GetXaxis()->SetTitle("FPGA fit quality. 0-good, 1-int overflow, 2-low amp, 3-bad chi2");
248  h_ampdiff_quality->GetYaxis()->SetTitle("Amplitude difference");
249  h_ampdiff_quality->SetOption("LIVE");
250 
251  h_timediff_quality = new TH2F("timediff_quality", "Time diff. (Emulator-Data) for time inconsistencies", 4, 0, 4, 239,
252  -4095, 4095);
253  h_timediff_quality->GetXaxis()->SetTitle("FPGA fit quality. 0-good, 1-int overflow, 2-low amp, 3-bad chi2");
254  h_timediff_quality->GetYaxis()->SetTitle("Time difference");
255  h_timediff_quality->SetOption("LIVE");
256 
257  h_quality_fit_data = new TH2F("quality_fit_data", "C++ fitter vs FPGA, fit quality inconsistencies", 4, 0, 4, 4, 0, 4);
258  h_quality_fit_data->GetXaxis()->SetTitle("C++ fit qual. 0-good,1-int overflow,2-low amp,3-bad chi2");
259  h_quality_fit_data->GetYaxis()->SetTitle("FPGA fit qual. 0-good,1-int overflow,2-low amp,3-bad chi2");
260  h_quality_fit_data->SetOption("LIVE");
261 
262  h_ampflag_qualityfail = new TH2F("ampflag_qualityfail", "Amp flag (0/1) for fit qual inconsistencies", 4, 0, 4, 4, -1,
263  3);
264  h_ampflag_qualityfail->GetXaxis()->SetTitle("FPGA fit quality. 0-good,1-int overflow,2-low amp,3-bad chi2");
265  h_ampflag_qualityfail->GetYaxis()->SetTitle("Amp flag (0-amp consistent)");
266  h_ampflag_qualityfail->SetOption("LIVE");
267 
268  h_timeflag_qualityfail = new TH2F("timeflag_qualityfail", "Time flag (0/1) for fit qual inconsistencies", 4, 0, 4, 4,
269  -1, 3);
270  h_timeflag_qualityfail->GetXaxis()->SetTitle("FPGA fit quality. 0-good,1-int overflow,2-low amp,3-bad chi2");
271  h_timeflag_qualityfail->GetYaxis()->SetTitle("Time flag (0-time consistent)");
272  h_timeflag_qualityfail->SetOption("LIVE");
273 
274  oldDir->cd();
275 }
276 
277 
279 {
280  REG_HISTOGRAM; // required to register histograms to HistoManager
281 
282  StoreArray<ECLDigit> ECLDigits;
283  ECLDigits.isRequired();
284 
285  StoreArray<ECLTrig> ECLTrigs;
286  ECLTrigs.isOptional();
287 
288  StoreArray<ECLDsp> ECLDsps;
289  ECLDsps.isOptional();
290 
291  if (!mapper.initFromDB()) B2FATAL("ECL DQM logic test FATAL:: Can't initialize eclChannelMapper");
292 
293  if (m_InitKey == "DB") initDspfromDB();
294  else if (m_InitKey == "File") initDspfromFile();
295  else B2FATAL("ECL DQM logic test FATAL: No way to initialize DSP coeffs!!! Please choose InitKey = DB or InitKey = File");
296 }
297 
298 void ECLDQMEXTENDEDModule::callbackCalibration(DBObjPtr<ECLCrystalCalib>& cal, std::vector<short int>& constants)
299 {
300  const std::vector<float> intermediate = cal->getCalibVector();
301  constants.resize(intermediate.size());
302  for (size_t i = 0; i < constants.size(); i++) constants[i] = (short int)intermediate[i];
303 }
304 
305 
306 void ECLDQMEXTENDEDModule::callbackCalibration(ECLDspData* dspdata, std::map<std::string, std::vector<short int>>& map1,
307  std::map<std::string, short int>& map2)
308 {
309 
310  dspdata->getF(map1["F"]);
311  dspdata->getF1(map1["F1"]);
312  dspdata->getF31(map1["F31"]);
313  dspdata->getF32(map1["F32"]);
314  dspdata->getF33(map1["F33"]);
315  dspdata->getF41(map1["F41"]);
316  dspdata->getF43(map1["F43"]);
317 
318  map2["k_a"] = (short int)dspdata->getka();
319  map2["k_b"] = (short int)dspdata->getkb();
320  map2["k_c"] = (short int)dspdata->getkc();
321  map2["k_16"] = (short int)dspdata->gety0Startr();
322  map2["k_1"] = (short int)dspdata->getk1();
323  map2["k_2"] = (short int)dspdata->getk2();
324  map2["chi_thres"] = dspdata->getchiThresh();
325 }
326 
327 
329 {
330  int iCrate = mapper.getCrateID(cellID);
331  int iShaperPosition = mapper.getShaperPosition(cellID);
332  return (iCrate - 1) * 12 + iShaperPosition;
333 
334 }
335 short int* ECLDQMEXTENDEDModule::vectorsplit(std::vector<short int>& vectorFrom, int iChannel)
336 {
337  size_t size = vectorFrom.size();
338  if (size % 16) B2ERROR("ECL DQM logic test error: Split is impossible!" << LogVar("Vector size", size));
339  return (vectorFrom.data() + (size / 16) * (iChannel - 1));
340 }
341 
342 
344 {
345  size_t iShaper = 0;
346 
347  callbackCalibration(m_calibrationThrA0, v_totalthrA0);
348  callbackCalibration(m_calibrationThrAhard, v_totalthrAhard);
349  callbackCalibration(m_calibrationThrAskip, v_totalthrAskip);
350 
351  for (const auto& dspdata : m_ECLDspDataArray0) { //iCrate = 1, ..., 18
352  iShaper++;
353  map_vec.clear();
354  map_coef.clear();
355  ECLDspData* dspointer = &(const_cast<ECLDspData&>(dspdata));
356  callbackCalibration(dspointer, map_vec, map_coef);
357  map_container_vec[iShaper] = map_vec;
358  map_container_coef[iShaper] = map_coef;
359  }
360 
361  for (const auto& dspdata : m_ECLDspDataArray1) { //iCrate = 19, ..., 36
362  iShaper++;
363  map_vec.clear();
364  map_coef.clear();
365  ECLDspData* dspointer = &(const_cast<ECLDspData&>(dspdata));
366  callbackCalibration(dspointer, map_vec, map_coef);
367  map_container_vec[iShaper] = map_vec;
368  map_container_coef[iShaper] = map_coef;
369  }
370 
371  for (const auto& dspdata : m_ECLDspDataArray2) { //iCrate = 37, ..., 52
372  iShaper++;
373  if (iShaper < 529) {
374  if (iShaper - (iShaper - 1) / 12 * 12 > 10) continue;
375  } else {
376  if (iShaper - (iShaper - 1) / 12 * 12 > 8) continue;
377  }
378  map_vec.clear();
379  map_coef.clear();
380  ECLDspData* dspointer = &(const_cast<ECLDspData&>(dspdata));
381  callbackCalibration(dspointer, map_vec, map_coef);
382  map_container_vec[iShaper] = map_vec;
383  map_container_coef[iShaper] = map_coef;
384  }
385 }
386 
388 {
389  const boost::filesystem::path MainDir(m_DSPDirectoryName);
390  const boost::filesystem::path RunSubDir(m_RunName);
391  const std::regex Filter(".*(crate)([0-9]{2})/.*(dsp)([0-9]{2})(.dat)");
392  if (!exists(MainDir / RunSubDir)) B2FATAL("ECL DQM logic test FATAL: Directory w/ DSP files don't exist" << LogVar("Directory",
393  MainDir / RunSubDir));
394  for (boost::filesystem::directory_entry& x : boost::filesystem::recursive_directory_iterator(MainDir / RunSubDir)) {
395  if (!std::regex_match(x.path().string(), Filter) || !boost::filesystem::is_regular_file(x.path())) continue;
396  int iCrate = atoi(std::regex_replace(x.path().string(), Filter, "$2").c_str());
397  int iShaperPosition = atoi(std::regex_replace(x.path().string(), Filter, "$4").c_str());
398  int iShaper = (iCrate - 1) * 12 + iShaperPosition;
399  if (iCrate > 36 && iCrate < 45) {
400  if (iShaperPosition > 10) continue;
401  } else if (iCrate > 44) {
402  if (iShaperPosition > 8) continue;
403  }
404  ECLDspData* dspdata = ECLDspUtilities::readEclDsp(x.path().string().c_str(), iShaperPosition - 1);
405  callbackCalibration(dspdata, map_vec, map_coef);
406  map_container_vec[iShaper] = map_vec;
407  map_container_coef[iShaper] = map_coef;
408  callbackCalibration(m_calibrationThrA0, v_totalthrA0);
409  callbackCalibration(m_calibrationThrAhard, v_totalthrAhard);
410  callbackCalibration(m_calibrationThrAskip, v_totalthrAskip);
411  }
412 }
413 
414 void ECLDQMEXTENDEDModule::emulator(int cellID, int trigger_time, std::vector<int> adc_data)
415 {
416 
417  int iShaper = conversion(cellID);
418  int iChannelPosition = mapper.getShaperChannel(cellID);
419  short int* f, *f1, *fg41, *fg43, *fg31, *fg32, *fg33;
420  int k_a, k_b, k_c, k_1, k_2, k_16, chi_thres;
421  int A0, Ahard, Askip;
422 
423  map_vec = map_container_vec[iShaper];
424  f = vectorsplit(map_vec["F"], iChannelPosition);
425  f1 = vectorsplit(map_vec["F1"], iChannelPosition);
426  fg31 = vectorsplit(map_vec["F31"], iChannelPosition);
427  fg32 = vectorsplit(map_vec["F32"], iChannelPosition);
428  fg33 = vectorsplit(map_vec["F33"], iChannelPosition);
429  fg41 = vectorsplit(map_vec["F41"], iChannelPosition);
430  fg43 = vectorsplit(map_vec["F43"], iChannelPosition);
431 
432  map_coef = map_container_coef[iShaper];
433  k_a = map_coef["k_a"];
434  k_b = map_coef["k_b"];
435  k_c = map_coef["k_c"];
436  k_1 = map_coef["k_1"];
437  k_2 = map_coef["k_2"];
438  k_16 = map_coef["k_16"];
439  chi_thres = map_coef["chi_thres"];
440 
441  A0 = (int)v_totalthrA0[cellID - 1];
442  Ahard = (int)v_totalthrAhard[cellID - 1];
443  Askip = (int)v_totalthrAskip[cellID - 1];
444 
445  int* y = adc_data.data();
446  int ttrig2 = trigger_time - 2 * (trigger_time / 8);
447 
448  auto result = lftda_(f, f1, fg41, fg43, fg31, fg32, fg33, y, ttrig2, A0,
449  Ahard, Askip, k_a, k_b, k_c, k_16, k_1, k_2, chi_thres,
450  m_adjusted_timing);
451  m_AmpFit = result.amp;
452  m_TimeFit = result.time;
453  m_QualityFit = result.quality;
454 
455  if (result.skip_thr || result.hit_thr) m_QualityFit += 4;
456 }
457 
459 {
460  for (int i = 0; i < 4; i++) {
461  h_amp_timefail[i]->Reset();
462  h_time_ampfail[i]->Reset();
463  }
464  for (int i = 0; i < 4; i++) {
465  for (int j = 0; j < 3; j++) {
466  h_amp_qualityfail[i][j]->Reset();
467  h_time_qualityfail[i][j]->Reset();
468  }
469  }
470  h_ampfail_quality->Reset();
471  h_timefail_quality->Reset();
472  h_ampfail_cellid->Reset();
473  h_timefail_cellid->Reset();
474  h_amptimefail_cellid->Reset();
475  h_ampfail_shaperid->Reset();
476  h_timefail_shaperid->Reset();
477  h_amptimefail_shaperid->Reset();
478  h_ampfail_crateid->Reset();
479  h_timefail_crateid->Reset();
480  h_amptimefail_crateid->Reset();
481  h_qualityfail_cellid->Reset();
482  h_qualityfail_shaperid->Reset();
483  h_qualityfail_crateid->Reset();
484  if (m_SaveDetailedFitData) {
485  h_ampdiff_cellid->Reset();
486  h_timediff_cellid->Reset();
487  h_ampdiff_shaperid->Reset();
488  h_timediff_shaperid->Reset();
489  }
490  h_ampdiff_quality->Reset();
491  h_timediff_quality->Reset();
492  h_quality_fit_data->Reset();
493  h_ampflag_qualityfail->Reset();
494  h_timeflag_qualityfail->Reset();
495 }
496 
498 {
499  StoreArray<ECLDsp> ECLDsps;
500  StoreArray<ECLTrig> ECLTrigs;
501  StoreArray<ECLDigit> ECLDigits;
502 
503  int iAmpflag_qualityfail = 0;
504  int iTimeflag_qualityfail = 0;
505 
506  for (auto& aECLDsp : ECLDsps) {
507  m_CellId = aECLDsp.getCellId();
508  std::vector<int> DspArray = aECLDsp.getDspA();
509  if (!ECLTrigs.isValid()) B2FATAL("ECL DQM logic test FATAL: Trigger time information is not available");
510  for (auto& aECLTrig : ECLTrigs) {
511  if (aECLTrig.getTrigId() == mapper.getCrateID(m_CellId)) {
512  m_TrigTime = aECLTrig.getTimeTrig();
513  break;
514  }
515  }
516 
517  emulator(m_CellId, m_TrigTime, DspArray);
518  ECLDigit* aECLDigit = aECLDsp.getRelated<ECLDigit>();
519 
520  if ((m_AmpFit >= (int)v_totalthrAskip[m_CellId - 1]) && m_QualityFit < 4 && !aECLDigit)
521  B2ERROR("ECL DQM logic test error: ECL Digit does not exist for A_emulator > Thr_skip"
522  << LogVar("Thr_skip", (int)v_totalthrAskip[m_CellId - 1])
523  << LogVar("A_emulator", m_AmpFit)
524  << LogVar("Quality_emulator", m_QualityFit));
525 
526  if ((m_AmpFit >= (int)v_totalthrAskip[m_CellId - 1]) && aECLDigit) {
527 
528  m_AmpData = aECLDigit->getAmp();
529  m_TimeData = aECLDigit->getTimeFit();
530  m_QualityData = aECLDigit->getQuality();
531 
532  if (m_AmpFit != m_AmpData) {
533  for (int i = 0; i < 4; i++) if (m_QualityData == i && m_TimeFit == m_TimeData) h_time_ampfail[i]->Fill(m_TimeData);
534  h_ampfail_quality->Fill(m_QualityData);
535  h_ampfail_cellid->Fill(m_CellId);
536  h_ampfail_shaperid->Fill(conversion(m_CellId));
537  h_ampfail_crateid->Fill(mapper.getCrateID(m_CellId));
538  if (m_SaveDetailedFitData) {
539  h_ampdiff_cellid->Fill(m_CellId, m_AmpFit - m_AmpData);
540  h_ampdiff_shaperid->Fill(conversion(m_CellId), m_AmpFit - m_AmpData);
541  }
542  h_ampdiff_quality->Fill(m_QualityData, m_AmpFit - m_AmpData);
543  }
544  if (m_TimeFit != m_TimeData) {
545  for (int i = 0; i < 4; i++) if (m_QualityData == i && m_AmpFit == m_AmpData) h_amp_timefail[i]->Fill(m_AmpData);
546  h_timefail_quality->Fill(m_QualityData);
547  h_timefail_cellid->Fill(m_CellId);
548  h_timefail_shaperid->Fill(conversion(m_CellId));
549  h_timefail_crateid->Fill(mapper.getCrateID(m_CellId));
550  if (m_SaveDetailedFitData) {
551  h_timediff_cellid->Fill(m_CellId, m_TimeFit - m_TimeData);
552  h_timediff_shaperid->Fill(conversion(m_CellId), m_TimeFit - m_TimeData);
553  }
554  h_timediff_quality->Fill(m_QualityData, m_TimeFit - m_TimeData);
555  }
556  if (m_AmpFit != m_AmpData && m_TimeFit != m_TimeData) {
557  h_amptimefail_cellid->Fill(m_CellId);
558  h_amptimefail_shaperid->Fill(conversion(m_CellId));
559  h_amptimefail_crateid->Fill(mapper.getCrateID(m_CellId));
560  }
561  if (m_QualityFit != m_QualityData) {
562  for (int i = 0; i < 4; i++) {
563  for (int j = 0; j < 3; j++) {
564  int k = 0;
565  if (i == 0) k = j + 1;
566  if (i == 1) { if (j == 0) k = j; else k = j + 1; }
567  if (i == 2) { if (j == 2) k = j + 1; else k = j; }
568  if (i == 3) k = j;
569  if (m_QualityFit == i && m_QualityData == k) { h_amp_qualityfail[i][j]->Fill(m_AmpData); h_time_qualityfail[i][j]->Fill(m_TimeData); }
570  }
571  }
572  if (m_AmpFit != m_AmpData) iAmpflag_qualityfail = 1;
573  if (m_TimeFit != m_TimeData) iTimeflag_qualityfail = 1;
574  h_qualityfail_cellid->Fill(m_CellId);
575  h_qualityfail_shaperid->Fill(conversion(m_CellId));
576  h_qualityfail_crateid->Fill(mapper.getCrateID(m_CellId));
577  h_ampflag_qualityfail->Fill(m_QualityData, iAmpflag_qualityfail);
578  h_timeflag_qualityfail->Fill(m_QualityData, iTimeflag_qualityfail);
579  h_quality_fit_data->Fill(m_QualityFit, m_QualityData);
580  }
581  if (m_AmpFit != m_AmpData || m_TimeFit != m_TimeData || m_QualityFit != m_QualityData) {
582  h_fail_shaperid->Fill(conversion(m_CellId));
583  h_fail_crateid->Fill(mapper.getCrateID(m_CellId));
584  }
585  h_ampfail_quality->Fill(-1);
586  h_timefail_quality->Fill(-1);
587  } //aECLDigit
588  } //aECLDsp
589 } //event
590 
591 
593 {
594 }
595 
596 
597 
599 {
600 }
Belle2::ECLDspData::getF31
void getF31(std::vector< short int > &dst) const
Array FG31, used to estimate signal amplitude.
Definition: ECLDspData.h:166
Belle2::ECLDQMEXTENDEDModule::initDspfromDB
void initDspfromDB()
Get DSP coeffs and auxiliary constants from DB.
Definition: eclDQMExtended.cc:343
Belle2::ECLDspData::getk1
unsigned char getk1() const
multipliers power of 2 for f, f1
Definition: ECLDspData.h:211
Belle2::ECLDQMEXTENDEDModule::conversion
int conversion(int)
Convert a CellID number to the global Shaper number.
Definition: eclDQMExtended.cc:328
Belle2::ECLDigit::getQuality
int getQuality() const
Get Fitting Quality.
Definition: ECLDigit.h:90
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::ECLDspData::getchiThresh
short int getchiThresh() const
chi2 threshold for fit quality flag
Definition: ECLDspData.h:209
Belle2::ECLDQMEXTENDEDModule
This module is created to monitor ECL electronics logic in frame of DQM system.
Definition: eclDQMExtended.h:47
Belle2::ECLDQMEXTENDEDModule::beginRun
virtual void beginRun() override
Call when a run begins.
Definition: eclDQMExtended.cc:458
Belle2::ECLDspData::getF41
void getF41(std::vector< short int > &dst) const
Alternative for FG31 for signals with small amplitude.
Definition: ECLDspData.h:183
Belle2::RelationsInterface::getRelated
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
Definition: RelationsObject.h:280
Belle2::ECLDspData::getkc
unsigned char getkc() const
Number of bits for FG33, FG43.
Definition: ECLDspData.h:219
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::ECLDspData
This object contains ECL DSP coefs – electromagnetic calorimeter digital signal processing coefficien...
Definition: ECLDspData.h:44
Belle2::ECLDQMEXTENDEDModule::defineHisto
virtual void defineHisto() override
Function to define histograms.
Definition: eclDQMExtended.cc:98
Belle2::ECLDspData::getF
void getF(std::vector< short int > &dst) const
Array with tabulated signal waveform.
Definition: ECLDspData.h:157
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::ECLDspData::getka
unsigned char getka() const
Number of bits for FG31, FG41.
Definition: ECLDspData.h:215
Belle2::ECLDspData::getF33
void getF33(std::vector< short int > &dst) const
Array FG33, used to estimate pedestal height in signal.
Definition: ECLDspData.h:179
Belle2::ECLDspData::getF43
void getF43(std::vector< short int > &dst) const
Alternative for FG33 for signals with small amplitude.
Definition: ECLDspData.h:187
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
Belle2::ECLDspData::getF1
void getF1(std::vector< short int > &dst) const
Array with tabulated derivative of signal waveform.
Definition: ECLDspData.h:161
Belle2::ECLDQMEXTENDEDModule::callbackCalibration
void callbackCalibration(DBObjPtr< ECLCrystalCalib > &, std::vector< short int > &)
Read calibration values for thresholds from DBObject.
Definition: eclDQMExtended.cc:298
Belle2::ECLDQMEXTENDEDModule::initialize
virtual void initialize() override
Initialize the module.
Definition: eclDQMExtended.cc:278
Belle2::StoreArray::isValid
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:298
Belle2::ECLDigit
Class to store ECL digitized hits (output of ECLDigi) relation to ECLHit filled in ecl/modules/eclDig...
Definition: ECLDigit.h:34
Belle2::ECLDQMEXTENDEDModule::vectorsplit
short int * vectorsplit(std::vector< short int > &, int)
Select from vector of DSP coeffs a subvector corresponding to accurate channel number.
Definition: eclDQMExtended.cc:335
Belle2::ECLDspData::gety0Startr
unsigned char gety0Startr() const
start point for pedestal calculation
Definition: ECLDspData.h:221
Belle2::ECLDigit::getTimeFit
int getTimeFit() const
Get Fitting Time.
Definition: ECLDigit.h:85
Belle2::ECLDQMEXTENDEDModule::~ECLDQMEXTENDEDModule
virtual ~ECLDQMEXTENDEDModule() override
Destructor.
Definition: eclDQMExtended.cc:90
Belle2::ECLDspData::getk2
unsigned char getk2() const
multipliers power of 2 for chi2 calculation
Definition: ECLDspData.h:213
Belle2::ECLDigit::getAmp
int getAmp() const
Get Fitting Amplitude.
Definition: ECLDigit.h:80
Belle2::ECLDspData::getF32
void getF32(std::vector< short int > &dst) const
Array FG32, used to estimate A * delta_t.
Definition: ECLDspData.h:173
Belle2::ECL::ECLDspUtilities::readEclDsp
static ECLDspData * readEclDsp(const char *raw_file, int boardNumber)
Convert ECLDspData from *.dat file to Root object.
Definition: ECLDspUtilities.cc:48
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::ECLDspData::getkb
unsigned char getkb() const
Number of bits for FG32.
Definition: ECLDspData.h:217
Belle2::ECLDQMEXTENDEDModule::initDspfromFile
void initDspfromFile()
Get DSP coeffs and auxiliary constants from Files.
Definition: eclDQMExtended.cc:387
Belle2::ECLDQMEXTENDEDModule::emulator
void emulator(int, int, std::vector< int >)
Call for DSP emulator.
Definition: eclDQMExtended.cc:414
Belle2::ECLDQMEXTENDEDModule::endRun
virtual void endRun() override
Call when a run ends.
Definition: eclDQMExtended.cc:592
Belle2::ECLDQMEXTENDEDModule::terminate
virtual void terminate() override
Terminate.
Definition: eclDQMExtended.cc:598
Belle2::ECLDQMEXTENDEDModule::event
virtual void event() override
Event processor.
Definition: eclDQMExtended.cc:497
Belle2::Filter
This class is used to select pairs, triplets...
Definition: Filter.h:44
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29