Belle II Software  release-08-01-10
eclDQMExtended.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 /* Own header. */
10 #include <ecl/modules/eclDQMExtended/eclDQMExtended.h>
11 
12 /* ECL headers. */
13 #include <ecl/dataobjects/ECLElementNumbers.h>
14 #include <ecl/utility/ECLDspEmulator.h>
15 #include <ecl/utility/ECLDspUtilities.h>
16 
17 /* Basf2 headers. */
18 #include <framework/core/HistoModule.h>
19 #include <framework/logging/Logger.h>
20 
21 /* ROOT headers. */
22 #include <TDirectory.h>
23 #include <TH1F.h>
24 #include <TH2F.h>
25 
26 /* Boost headers. */
27 #include <boost/format.hpp>
28 
29 /* C++ headers. */
30 #include <filesystem>
31 #include <iostream>
32 #include <map>
33 #include <regex>
34 #include <string>
35 #include <vector>
36 
37 
38 //NAMESPACE(S)
39 using namespace Belle2;
40 using namespace ECL;
41 
42 //-----------------------------------------------------------------
43 // Register the Module
44 //-----------------------------------------------------------------
45 REG_MODULE(ECLDQMEXTENDED);
46 
47 
48 //-----------------------------------------------------------------
49 // Implementation
50 //-----------------------------------------------------------------
51 
53  : HistoModule(),
54  m_ECLDspDataArray0("ECLDSPPars0"),
55  m_ECLDspDataArray1("ECLDSPPars1"),
56  m_ECLDspDataArray2("ECLDSPPars2"),
57  m_calibrationThrA0("ECL_FPGA_LowAmp"),
58  m_calibrationThrAhard("ECL_FPGA_HitThresh"),
59  m_calibrationThrAskip("ECL_FPGA_StoreDigit")
60 
61 {
62 
63  //Set module properties
64  setDescription("ECL Data Quality Monitor. Logic Test");
65  setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
66  addParam("histogramDirectoryName", m_histogramDirectoryName,
67  "histogram directory in ROOT file", std::string("ECL"));
68  addParam("InitKey", m_InitKey,
69  "How to initialize DSP coeffs: ''DB'' or ''File'' are acceptable ", std::string("DB"));
70  addParam("DSPDirectoryName", m_DSPDirectoryName,
71  "directory for DSP coeffs", std::string("/hsm/belle2/bdata/users/dmitry/dsp/"));
72  addParam("RunName", m_RunName,
73  "Name of run with DSP files", std::string("run0000/"));
74  addParam("SaveDetailedFitData", m_SaveDetailedFitData, "Save detailed data "
75  "(ampdiff_{cellid,shaper}, timediff_{cellid,shaper} histograms) for "
76  "failed fits.", false);
77  addParam("AdjustedTiming", m_adjusted_timing, "Use improved procedure for "
78  "time determination in ShaperDSP emulator", true);
79 
80  std::vector<int> default_skip = {
81  // By default, skip delayed bhabha and random trigger events
82  TRGSummary::ETimingType::TTYP_DPHY,
83  TRGSummary::ETimingType::TTYP_RAND,
84  TRGSummary::ETimingType::TTYP_POIS,
85  };
86  addParam("skipEvents", m_skipEvents, "Skip events that have listed timing "
87  "source (from TRGSummary)", default_skip);
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_amp_timefail.push_back(h_a);
121  h_time_ampfail.push_back(h_t);
122  }
123 
124  for (int i = 0; i < 4; i++) {
125  for (int j = 0; j < 4; j++) {
126  if (j != i) {
127  std::string h_name_a, h_title_a, h_name_t, h_title_t;
128  h_name_a = str(boost::format("amp_qf%1%_qd%2%") % i % j);
129  h_title_a = str(boost::format("Amp for C++ fit qual=%1% and FPGA fit qual=%2%") % i % j);
130  h_name_t = str(boost::format("time_qf%1%_qd%2%") % i % j);
131  h_title_t = str(boost::format("Time for C++ fit qual=%1% and FPGA fit qual=%2%") % i % j);
132  TH1F* h_a = new TH1F(h_name_a.c_str(), h_title_a.c_str(), 1200, 0, 262144);
133  TH1F* h_t = new TH1F(h_name_t.c_str(), h_title_t.c_str(), 240, -2050, 2050);
134  h_amp_qualityfail_raw.push_back(h_a);
135  h_time_qualityfail_raw.push_back(h_t);
136  }
137  }
138  h_amp_qualityfail.push_back(h_amp_qualityfail_raw);
139  h_time_qualityfail.push_back(h_time_qualityfail_raw);
140  h_amp_qualityfail_raw.clear();
141  h_time_qualityfail_raw.clear();
142  }
143 
144  h_ampfail_quality = new TH1F("ampfail_quality", "Number of FPGA <-> C++ fitter #bf{amp} inconsistencies vs fit qual", 5, -1, 4);
145  h_ampfail_quality->SetFillColor(kPink - 4);
146  h_ampfail_quality->GetXaxis()->SetTitle("FPGA fit qual. -1-all evts,0-good,1-int overflow,2-low amp,3-bad chi2");
147  h_ampfail_quality->SetDrawOption("hist");
148 
149  h_timefail_quality = new TH1F("timefail_quality", "Number of FPGA <-> C++ fitter #bf{time} inconsistencies vs fit qual", 5, -1, 4);
150  h_timefail_quality->SetFillColor(kPink - 4);
151  h_timefail_quality->GetXaxis()->SetTitle("FPGA fit qual. -1-all evts,0-good,1-int overflow,2-low amp,3-bad chi2");
152 
153  h_ampfail_cellid = new TH1F("ampfail_cellid", "Cell IDs w/ amp inconsistencies",
155  h_ampfail_cellid->GetXaxis()->SetTitle("Cell ID");
156 
157  h_timefail_cellid = new TH1F("timefail_cellid", "Cell IDs w/ time inconsistencies",
159  h_timefail_cellid->GetXaxis()->SetTitle("Cell ID");
160 
161  h_amptimefail_cellid = new TH1F("amptimefail_cellid", "Cell IDs w/ time and amp inconsistencies",
163  h_amptimefail_cellid->GetXaxis()->SetTitle("Cell ID");
164 
165  h_ampfail_shaperid = new TH1F("ampfail_shaperid", "Shaper IDs w/ amp inconsistencies", 624, 1, 625);
166  h_ampfail_shaperid->GetXaxis()->SetTitle("Shaper ID");
167 
168  h_timefail_shaperid = new TH1F("timefail_shaperid", "Shaper IDs w/ time inconsistencies", 624, 1, 625);
169  h_timefail_shaperid->GetXaxis()->SetTitle("Shaper ID");
170 
171  h_amptimefail_shaperid = new TH1F("amptimefail_shaperid", "Shaper IDs w/ time and amp inconsistencies", 624, 1, 625);
172  h_amptimefail_shaperid->GetXaxis()->SetTitle("Shaper ID");
173 
174  h_ampfail_crateid = new TH1F("ampfail_crateid", "Crate IDs w/ amp inconsistencies", 52, 1, 53);
175  h_ampfail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
176 
177  h_timefail_crateid = new TH1F("timefail_crateid", "Crate IDs w/ time inconsistencies", 52, 1, 53);
178  h_timefail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
179 
180  h_amptimefail_crateid = new TH1F("amptimefail_crateid", "Crate IDs w/ time and amp inconsistencies", 52, 1, 53);
181  h_amptimefail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
182 
183  h_qualityfail_cellid = new TH1F("qualityfail_cellid", "Cell IDs w/ fit qual inconsistencies",
185  h_qualityfail_cellid->GetXaxis()->SetTitle("Cell ID");
186 
187  h_qualityfail_shaperid = new TH1F("qualityfail_shaperid", "Shaper IDs w/ fit qual inconsistencies", 624, 1, 625);
188  h_qualityfail_shaperid->GetXaxis()->SetTitle("Shaper ID");
189 
190  h_qualityfail_crateid = new TH1F("qualityfail_crateid", "Crate IDs w/ fit qual inconsistencies", 52, 1, 53);
191  h_qualityfail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
192 
193  h_fail_shaperid = new TH1F("fail_shaperid", "Shaper IDs w/ inconsistencies", 624, 1, 625);
194  h_fail_shaperid->GetXaxis()->SetTitle("Shaper ID");
195 
196  h_fail_crateid = new TH1F("fail_crateid", "Crate IDs w/ inconsistencies", 52, 1, 53);
197  h_fail_crateid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
198 
199 
200  //2D histograms creation.
201 
202  if (m_SaveDetailedFitData) {
203  h_ampdiff_cellid = new TH2F("ampdiff_cellid", "Amp. diff. (Emulator-Data) for amp inconsistencies",
204  ECLElementNumbers::c_NCrystals, 1, ECLElementNumbers::c_NCrystals + 1, 239, -262143, 262143);
205  h_ampdiff_cellid->GetXaxis()->SetTitle("Cell ID");
206  h_ampdiff_cellid->GetYaxis()->SetTitle("Amplitude difference");
207 
208  h_timediff_cellid = new TH2F("timediff_cellid", "Time diff. (Emulator-Data) for time inconsistencies",
210  h_timediff_cellid->GetXaxis()->SetTitle("Cell ID");
211  h_timediff_cellid->GetYaxis()->SetTitle("Time difference");
212 
213  h_ampdiff_shaperid = new TH2F("ampdiff_shaper", "Amp. diff. (Emulator-Data) "
214  "for amp inconsistencies vs Shaper Id",
215  624, 1, 625, 239, -262143, 262143);
216  h_ampdiff_shaperid->GetXaxis()->SetTitle("Shaper Id");
217  h_ampdiff_shaperid->GetYaxis()->SetTitle("Amplitude difference");
218 
219  h_timediff_shaperid = new TH2F("timediff_shaper", "Time diff. (Emulator-Data) "
220  "for time inconsistencies vs Shaper Id",
221  624, 1, 625, 239, -4095, 4095);
222  h_timediff_shaperid->GetXaxis()->SetTitle("Shaper Id");
223  h_timediff_shaperid->GetYaxis()->SetTitle("Time difference");
224  }
225 
226  h_ampdiff_quality = new TH2F("ampdiff_quality", "Amp. diff. (Emulator-Data) for amp. inconsistencies", 4, 0, 4, 239,
227  -262143, 262143);
228  h_ampdiff_quality->GetXaxis()->SetTitle("FPGA fit quality. 0-good, 1-int overflow, 2-low amp, 3-bad chi2");
229  h_ampdiff_quality->GetYaxis()->SetTitle("Amplitude difference");
230 
231  h_timediff_quality = new TH2F("timediff_quality", "Time diff. (Emulator-Data) for time inconsistencies", 4, 0, 4, 239,
232  -4095, 4095);
233  h_timediff_quality->GetXaxis()->SetTitle("FPGA fit quality. 0-good, 1-int overflow, 2-low amp, 3-bad chi2");
234  h_timediff_quality->GetYaxis()->SetTitle("Time difference");
235 
236  h_quality_fit_data = new TH2F("quality_fit_data", "C++ fitter vs FPGA, fit quality inconsistencies", 4, 0, 4, 4, 0, 4);
237  h_quality_fit_data->GetXaxis()->SetTitle("C++ fit qual. 0-good,1-int overflow,2-low amp,3-bad chi2");
238  h_quality_fit_data->GetYaxis()->SetTitle("FPGA fit qual. 0-good,1-int overflow,2-low amp,3-bad chi2");
239 
240  h_ampflag_qualityfail = new TH2F("ampflag_qualityfail", "Amp flag (0/1) for fit qual inconsistencies", 4, 0, 4, 4, -1,
241  3);
242  h_ampflag_qualityfail->GetXaxis()->SetTitle("FPGA fit quality. 0-good,1-int overflow,2-low amp,3-bad chi2");
243  h_ampflag_qualityfail->GetYaxis()->SetTitle("Amp flag (0-amp consistent)");
244 
245  h_timeflag_qualityfail = new TH2F("timeflag_qualityfail", "Time flag (0/1) for fit qual inconsistencies", 4, 0, 4, 4,
246  -1, 3);
247  h_timeflag_qualityfail->GetXaxis()->SetTitle("FPGA fit quality. 0-good,1-int overflow,2-low amp,3-bad chi2");
248  h_timeflag_qualityfail->GetYaxis()->SetTitle("Time flag (0-time consistent)");
249 
250  h_missing_ecldigits = new TH1F("missing_ecldigits", "", 4, 0, 4);
251  h_missing_ecldigits->SetTitle("Fit quality flag for missing ECLDigits (missing fit results)");
252  h_missing_ecldigits->GetXaxis()->SetTitle("C++ fit quality. 0-good,1-int overflow,2-low amplitude,3-bad chi2");
253 
254  oldDir->cd();
255 }
256 
257 
259 {
260  REG_HISTOGRAM; // required to register histograms to HistoManager
261 
262  m_ECLDigits.isRequired();
263  m_ECLTrigs.isOptional();
264  m_ECLDsps.isOptional();
265 
266  if (!mapper.initFromDB()) B2FATAL("ECL DQM logic test FATAL:: Can't initialize eclChannelMapper");
267 
268  if (m_InitKey == "DB") initDspfromDB();
269  else if (m_InitKey == "File") initDspfromFile();
270  else B2FATAL("ECL DQM logic test FATAL: No way to initialize DSP coeffs!!! Please choose InitKey = DB or InitKey = File");
271 }
272 
273 void ECLDQMEXTENDEDModule::callbackCalibration(DBObjPtr<ECLCrystalCalib>& cal, std::vector<short int>& constants)
274 {
275  const std::vector<float> intermediate = cal->getCalibVector();
276  constants.resize(intermediate.size());
277  for (size_t i = 0; i < constants.size(); i++) constants[i] = (short int)intermediate[i];
278 }
279 
280 
282  std::map<std::string, std::vector<short int>>& map1,
283  std::map<std::string, short int>& map2)
284 {
285 
286  dspdata->getF(map1["F"]);
287  dspdata->getF1(map1["F1"]);
288  dspdata->getF31(map1["F31"]);
289  dspdata->getF32(map1["F32"]);
290  dspdata->getF33(map1["F33"]);
291  dspdata->getF41(map1["F41"]);
292  dspdata->getF43(map1["F43"]);
293 
294  map2["k_a"] = (short int)dspdata->getka();
295  map2["k_b"] = (short int)dspdata->getkb();
296  map2["k_c"] = (short int)dspdata->getkc();
297  map2["k_16"] = (short int)dspdata->gety0Startr();
298  map2["k_1"] = (short int)dspdata->getk1();
299  map2["k_2"] = (short int)dspdata->getk2();
300  map2["chi_thres"] = dspdata->getchiThresh();
301 }
302 
303 
305 {
306  int iCrate = mapper.getCrateID(cellID);
307  int iShaperPosition = mapper.getShaperPosition(cellID);
308  return (iCrate - 1) * 12 + iShaperPosition;
309 
310 }
311 const short int* ECLDQMEXTENDEDModule::vectorsplit(const std::vector<short int>& vectorFrom, int iChannel)
312 {
313  size_t size = vectorFrom.size();
314  if (size % 16) B2ERROR("ECL DQM logic test error: Split is impossible!" << LogVar("Vector size", size));
315  return (vectorFrom.data() + (size / 16) * (iChannel - 1));
316 }
317 
318 
320 {
321  size_t iShaper = 0;
322 
326 
327  for (const auto& dspdata : m_ECLDspDataArray0) { //iCrate = 1, ..., 18
328  iShaper++;
329  callbackCalibration(&dspdata, map_container_vec[iShaper],
330  map_container_coef[iShaper]);
331  }
332 
333  for (const auto& dspdata : m_ECLDspDataArray1) { //iCrate = 19, ..., 36
334  iShaper++;
335  callbackCalibration(&dspdata, map_container_vec[iShaper],
336  map_container_coef[iShaper]);
337  }
338 
339  for (const auto& dspdata : m_ECLDspDataArray2) { //iCrate = 37, ..., 52
340  iShaper++;
341  if (iShaper < 529) {
342  if (iShaper - (iShaper - 1) / 12 * 12 > 10) continue;
343  } else {
344  if (iShaper - (iShaper - 1) / 12 * 12 > 8) continue;
345  }
346  callbackCalibration(&dspdata, map_container_vec[iShaper],
347  map_container_coef[iShaper]);
348  }
349 }
350 
352 {
353  const std::filesystem::path MainDir(m_DSPDirectoryName);
354  const std::filesystem::path RunSubDir(m_RunName);
355  const std::regex Filter(".*(crate)([0-9]{2})/.*(dsp)([0-9]{2})(.dat)");
356  if (!exists(MainDir / RunSubDir)) B2FATAL("ECL DQM logic test FATAL: Directory w/ DSP files don't exist" << LogVar("Directory",
357  MainDir / RunSubDir));
358  for (const std::filesystem::directory_entry& x : std::filesystem::recursive_directory_iterator(MainDir / RunSubDir)) {
359  if (!std::regex_match(x.path().string(), Filter) || !std::filesystem::is_regular_file(x.path())) continue;
360  int iCrate = atoi(std::regex_replace(x.path().string(), Filter, "$2").c_str());
361  int iShaperPosition = atoi(std::regex_replace(x.path().string(), Filter, "$4").c_str());
362  int iShaper = (iCrate - 1) * 12 + iShaperPosition;
363  if (iCrate > 36 && iCrate < 45) {
364  if (iShaperPosition > 10) continue;
365  } else if (iCrate > 44) {
366  if (iShaperPosition > 8) continue;
367  }
368  ECLDspData* dspdata = ECLDspUtilities::readEclDsp(x.path().string().c_str(), iShaperPosition - 1);
369  callbackCalibration(dspdata, map_container_vec[iShaper],
370  map_container_coef[iShaper]);
374  }
375 }
376 
377 void ECLDQMEXTENDEDModule::emulator(int cellID, int trigger_time, std::vector<int> adc_data)
378 {
379  int iShaper = conversion(cellID);
380  int iChannelPosition = mapper.getShaperChannel(cellID);
381 
382  const auto& map_vec = map_container_vec[iShaper];
383  const short int* f = vectorsplit(map_vec.at("F"), iChannelPosition);
384  const short int* f1 = vectorsplit(map_vec.at("F1"), iChannelPosition);
385  const short int* fg31 = vectorsplit(map_vec.at("F31"), iChannelPosition);
386  const short int* fg32 = vectorsplit(map_vec.at("F32"), iChannelPosition);
387  const short int* fg33 = vectorsplit(map_vec.at("F33"), iChannelPosition);
388  const short int* fg41 = vectorsplit(map_vec.at("F41"), iChannelPosition);
389  const short int* fg43 = vectorsplit(map_vec.at("F43"), iChannelPosition);
390 
391  const auto& map_coef = map_container_coef[iShaper];
392  int k_a = map_coef.at("k_a");
393  int k_b = map_coef.at("k_b");
394  int k_c = map_coef.at("k_c");
395  int k_1 = map_coef.at("k_1");
396  int k_2 = map_coef.at("k_2");
397  int k_16 = map_coef.at("k_16");
398  int chi_thres = map_coef.at("chi_thres");
399 
400  int A0 = (int)v_totalthrA0[cellID - 1];
401  int Ahard = (int)v_totalthrAhard[cellID - 1];
402  int Askip = (int)v_totalthrAskip[cellID - 1];
403 
404  int* y = adc_data.data();
405  int ttrig2 = trigger_time - 2 * (trigger_time / 8);
406 
407  auto result = lftda_(f, f1, fg41, fg43, fg31, fg32, fg33, y, ttrig2, A0,
408  Ahard, Askip, k_a, k_b, k_c, k_16, k_1, k_2, chi_thres,
410  m_AmpFit = result.amp;
411  m_TimeFit = result.time;
412  m_QualityFit = result.quality;
413 
414  if (result.skip_thr || result.hit_thr) m_QualityFit += 4;
415 }
416 
418 {
419  for (int i = 0; i < 4; i++) {
420  h_amp_timefail[i]->Reset();
421  h_time_ampfail[i]->Reset();
422  }
423  for (int i = 0; i < 4; i++) {
424  for (int j = 0; j < 3; j++) {
425  h_amp_qualityfail[i][j]->Reset();
426  h_time_qualityfail[i][j]->Reset();
427  }
428  }
429  h_ampfail_quality->Reset();
430  h_timefail_quality->Reset();
431  h_ampfail_cellid->Reset();
432  h_timefail_cellid->Reset();
433  h_amptimefail_cellid->Reset();
434  h_ampfail_shaperid->Reset();
435  h_timefail_shaperid->Reset();
436  h_amptimefail_shaperid->Reset();
437  h_ampfail_crateid->Reset();
438  h_timefail_crateid->Reset();
439  h_amptimefail_crateid->Reset();
440  h_qualityfail_cellid->Reset();
441  h_qualityfail_shaperid->Reset();
442  h_qualityfail_crateid->Reset();
443  if (m_SaveDetailedFitData) {
444  h_ampdiff_cellid->Reset();
445  h_timediff_cellid->Reset();
446  h_ampdiff_shaperid->Reset();
447  h_timediff_shaperid->Reset();
448  }
449  h_ampdiff_quality->Reset();
450  h_timediff_quality->Reset();
451  h_quality_fit_data->Reset();
452  h_ampflag_qualityfail->Reset();
453  h_timeflag_qualityfail->Reset();
454  h_missing_ecldigits->Reset();
455 }
456 
458 {
459  if (m_TRGSummary.isValid()) {
460  // Skip events for several types of timing sources
461  // to reduce CPU time usage on HLT (random trigger, delayed bhahba).
462  int timing_type = m_TRGSummary->getTimType();
463  for (auto& skipped_timing_type : m_skipEvents) {
464  if (timing_type == skipped_timing_type) return;
465  }
466  }
467 
468  int iAmpflag_qualityfail = 0;
469  int iTimeflag_qualityfail = 0;
470 
471  if (!m_ECLTrigs.isValid()) B2FATAL("ECL DQM logic test FATAL: Trigger time information is not available");
472 
473  for (auto& aECLDsp : m_ECLDsps) {
474  m_CellId = aECLDsp.getCellId();
475  std::vector<int> DspArray = aECLDsp.getDspA();
477 
478  emulator(m_CellId, m_TrigTime, DspArray);
480 
481  if ((m_AmpFit >= (int)v_totalthrAskip[m_CellId - 1]) && m_QualityFit < 4 && !aECLDigit) {
482  if (h_missing_ecldigits->GetEntries() == 0) {
483  B2ERROR("ECL DQM logic test error: ECL Digit does not exist for A_emulator > Thr_skip"
484  << LogVar("Thr_skip", (int)v_totalthrAskip[m_CellId - 1])
485  << LogVar("A_emulator", m_AmpFit)
486  << LogVar("Quality_emulator", m_QualityFit));
487  }
489  }
490 
491  if ((m_AmpFit >= (int)v_totalthrAskip[m_CellId - 1]) && aECLDigit) {
492 
493  m_AmpData = aECLDigit->getAmp();
494  m_TimeData = aECLDigit->getTimeFit();
495  m_QualityData = aECLDigit->getQuality();
496 
497  if (m_AmpFit != m_AmpData) {
498  for (int i = 0; i < 4; i++) if (m_QualityData == i && m_TimeFit == m_TimeData) h_time_ampfail[i]->Fill(m_TimeData);
500  h_ampfail_cellid->Fill(m_CellId);
503  if (m_SaveDetailedFitData) {
506  }
508  }
509  if (m_TimeFit != m_TimeData) {
510  for (int i = 0; i < 4; i++) if (m_QualityData == i && m_AmpFit == m_AmpData) h_amp_timefail[i]->Fill(m_AmpData);
515  if (m_SaveDetailedFitData) {
518  }
520  }
521  if (m_AmpFit != m_AmpData && m_TimeFit != m_TimeData) {
525  }
526  if (m_QualityFit != m_QualityData) {
527  for (int i = 0; i < 4; i++) {
528  for (int j = 0; j < 3; j++) {
529  int k = 0;
530  if (i == 0) k = j + 1;
531  if (i == 1) { if (j == 0) k = j; else k = j + 1; }
532  if (i == 2) { if (j == 2) k = j + 1; else k = j; }
533  if (i == 3) k = j;
534  if (m_QualityFit == i && m_QualityData == k) { h_amp_qualityfail[i][j]->Fill(m_AmpData); h_time_qualityfail[i][j]->Fill(m_TimeData); }
535  }
536  }
537  if (m_AmpFit != m_AmpData) iAmpflag_qualityfail = 1;
538  if (m_TimeFit != m_TimeData) iTimeflag_qualityfail = 1;
542  h_ampflag_qualityfail->Fill(m_QualityData, iAmpflag_qualityfail);
543  h_timeflag_qualityfail->Fill(m_QualityData, iTimeflag_qualityfail);
545  }
549  }
550  h_ampfail_quality->Fill(-1);
551  h_timefail_quality->Fill(-1);
552  } //aECLDigit
553  } //aECLDsp
554 } //event
555 
556 
558 {
559 }
560 
561 
562 
564 {
565 }
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
TH2F * h_timediff_cellid
Histogram: Time diff.
virtual ~ECLDQMEXTENDEDModule() override
Destructor.
TH1F * h_amptimefail_shaperid
Histogram: ShaperIDs w/ failed amplitudes and times.
int m_TimeData
Signal time from ECL data.
void initDspfromDB()
Get DSP coeffs and auxiliary constants from DB.
TH2F * h_timediff_shaperid
Histogram: Time diff.
TH1F * h_ampfail_crateid
Histogram: CrateIDs w/ failed amplitudes.
TH1F * h_fail_crateid
Histogram: CrateIDs w/ failed logic.
bool m_SaveDetailedFitData
Save detailed fit data for failed fits.
std::vector< short int > v_totalthrA0
Vector to store Low amplitude thresholds.
TH1F * h_timefail_quality
Histogram: Time control flags in bins of QualityData.
bool m_adjusted_timing
Use modified time determination algorithm in emulator, same as in ShaperDSP version >= 1....
int m_AmpFit
Signal amplitude obtaining from DSP emulator.
int m_TrigTime
Trigger time value.
int m_CellId
Cell ID number.
StoreArray< ECLDsp > m_ECLDsps
ECL DSP data.
virtual void initialize() override
Initialize the module.
TH1F * h_timefail_shaperid
Histogram: ShaperIDs w/ failed times.
std::vector< short int > v_totalthrAhard
Vector to store hit thresholds.
int m_QualityData
Quality flag from ECL data.
TH1F * h_qualityfail_shaperid
Histogram: ShaperIDs w/ failed qualities.
virtual void event() override
Event processor.
TH1F * h_fail_shaperid
Histogram: ShaperIDs w/ failed logic.
void emulator(int, int, std::vector< int >)
Call for DSP emulator.
TH2F * h_timediff_quality
Histogram: Time diff.
virtual void endRun() override
Call when a run ends.
TH1F * h_qualityfail_crateid
Histogram: CrateIDs w/ failed qualities.
TH1F * h_amptimefail_crateid
Histogram: CrateIDs w/ failed amplitudes and times.
virtual void terminate() override
Terminate.
int m_AmpData
Signal amplitude from ECL data.
void callbackCalibration(DBObjPtr< ECLCrystalCalib > &, std::vector< short int > &)
Read calibration values for thresholds from DBObject.
std::vector< short int > v_totalthrAskip
Vector to store skip amplitude threshold.
std::map< int, std::map< std::string, std::vector< short int > > > map_container_vec
Map to store DSP coeffs.
std::vector< int > m_skipEvents
Skip events with specified type of timing source (see TRGSummary class, ETimingType enum)
StoreArray< ECLDigit > m_ECLDigits
ECL digits.
DBArray< ECLDspData > m_ECLDspDataArray2
DBArray for payload 'ECLDSPPars2'.
std::vector< TH1F * > h_amp_timefail
Histogram vector: Amplitude for time mismacthes (AmplitudeFit == AmplitudeData) w/ various QualityDat...
TH1F * h_amptimefail_cellid
Histogram: CellIDs w/ failed amplitudes and times.
TH2F * h_ampflag_qualityfail
Histogram: Amp flag (0/1) w/ failed qualities in bins of QualityData.
DBObjPtr< ECLCrystalCalib > m_calibrationThrAhard
Hit threshold.
TH1F * h_timefail_crateid
Histogram: CrateIDs w/ failed times.
const short int * vectorsplit(const std::vector< short int > &, int)
Select from vector of DSP coeffs a subvector corresponding to accurate channel number.
std::string m_InitKey
Key to initialize DSP coeffs.
TH2F * h_quality_fit_data
Histogram: QualityFit vs QualityData for quality fails.
DBObjPtr< ECLCrystalCalib > m_calibrationThrAskip
Skip amplitude threshold.
int conversion(int)
Convert a CellID number to the global Shaper number.
virtual void beginRun() override
Call when a run begins.
std::string m_histogramDirectoryName
Histogram directory in ROOT file.
ECL::ECLChannelMapper mapper
ECL channel mapper.
ECLDQMEXTENDEDModule()
< derived from HistoModule class.
TH2F * h_ampdiff_quality
Histogram: Amp.
DBArray< ECLDspData > m_ECLDspDataArray1
DBArray for payload 'ECLDSPPars1'.
TH1F * h_ampfail_cellid
Histogram: CellIDs w/ failed amplitudes.
TH2F * h_ampdiff_shaperid
Histogram: Amp.
TH1F * h_ampfail_quality
Histogram: Amp.
int m_QualityFit
Quality flag obtaining from DSP emulator.
std::string m_RunName
Run with valid DSP coeffs.
DBArray< ECLDspData > m_ECLDspDataArray0
DBArray for payload 'ECLDSPPars0'.
TH2F * h_ampdiff_cellid
Histogram: Amplitude diff.
TH1F * h_ampfail_shaperid
Histogram: ShaperIDs w/ failed amplitudes.
TH1F * h_timefail_cellid
Histogram: CellIDs w/ failed times.
StoreArray< ECLTrig > m_ECLTrigs
ECL trigger data.
DBObjPtr< ECLCrystalCalib > m_calibrationThrA0
Low amplitude threshold.
TH1F * h_qualityfail_cellid
Histogram: CellIDs w/ failed qualities.
std::map< int, std::map< std::string, short int > > map_container_coef
Map to store auxiliary constants for all shapers.
std::vector< std::vector< TH1F * > > h_amp_qualityfail
Histogram vector: AmplitudeData for quality mismathes w/ various QualityFit (raw) and QualityData (co...
StoreObjPtr< TRGSummary > m_TRGSummary
StoreObjPtr TRGSummary
std::vector< TH1F * > h_time_ampfail
Histogram vector: Time for amplitude mismathces (TimeFit == TimeData) w/ various QualityData values.
std::string m_DSPDirectoryName
Directory name consisting of DSP coeffs.
int m_TimeFit
Signal time obtaining from DSP emulator.
void initDspfromFile()
Get DSP coeffs and auxiliary constants from Files.
TH2F * h_timeflag_qualityfail
Histogram: Time flag (0/1) w/ failed qualities in bins of Quality Data.
TH1F * h_missing_ecldigits
Histogram: Quality flag (from emulator) for the waveforms where fit result (ECLDigit) should have bee...
std::vector< std::vector< TH1F * > > h_time_qualityfail
Histogram vector: TimeData for quality mismacthes w/ various QualitFit (raw) and QualityData (column)...
virtual void defineHisto() override
Function to define histograms.
Class to store ECL digitized hits (output of ECLDigi) relation to ECLHit filled in ecl/modules/eclDig...
Definition: ECLDigit.h:24
int getAmp() const
Get Fitting Amplitude.
Definition: ECLDigit.h:70
int getQuality() const
Get Fitting Quality.
Definition: ECLDigit.h:80
static ECLDigit * getByCellID(int cid)
Find ECLDigit by Cell ID using linear search.
Definition: ECLDigit.cc:14
int getTimeFit() const
Get Fitting Time.
Definition: ECLDigit.h:75
This object contains ECL DSP coefs – electromagnetic calorimeter digital signal processing coefficien...
Definition: ECLDspData.h:34
unsigned char getkb() const
Number of bits for FG32.
Definition: ECLDspData.h:200
short int getchiThresh() const
chi2 threshold for fit quality flag
Definition: ECLDspData.h:192
void getF43(std::vector< short int > &dst) const
Alternative for FG33 for signals with small amplitude.
Definition: ECLDspData.h:170
unsigned char gety0Startr() const
start point for pedestal calculation
Definition: ECLDspData.h:204
void getF(std::vector< short int > &dst) const
Array with tabulated signal waveform.
Definition: ECLDspData.h:140
void getF31(std::vector< short int > &dst) const
Array FG31, used to estimate signal amplitude.
Definition: ECLDspData.h:149
unsigned char getk2() const
multipliers power of 2 for chi2 calculation
Definition: ECLDspData.h:196
void getF33(std::vector< short int > &dst) const
Array FG33, used to estimate pedestal height in signal.
Definition: ECLDspData.h:162
unsigned char getka() const
Number of bits for FG31, FG41.
Definition: ECLDspData.h:198
unsigned char getk1() const
multipliers power of 2 for f, f1
Definition: ECLDspData.h:194
void getF32(std::vector< short int > &dst) const
Array FG32, used to estimate A * delta_t.
Definition: ECLDspData.h:156
void getF41(std::vector< short int > &dst) const
Alternative for FG31 for signals with small amplitude.
Definition: ECLDspData.h:166
void getF1(std::vector< short int > &dst) const
Array with tabulated derivative of signal waveform.
Definition: ECLDspData.h:144
unsigned char getkc() const
Number of bits for FG33, FG43.
Definition: ECLDspData.h:202
double getTimeTrig() const
Get Trig Time.
Definition: ECLTrig.h:111
static ECLTrig * getByCellID(int cid)
Find ECLTrig by Cell ID using linear search.
Definition: ECLTrig.cc:15
bool initFromDB()
Initialize channel mapper from the conditions database.
int getShaperChannel(int cellID)
Get number of DSP channel in the shaper by given number of CellId.
int getShaperPosition(int cellID)
Get position of the shaper in the crate by given CellId.
int getCrateID(int iCOPPERNode, int iFINESSE, bool pcie40=false)
Get crate number by given COPPER node number and FINESSE number.
static ECLDspData * readEclDsp(const char *raw_file, int boardNumber)
Convert ECLDspData from *.dat file to Root object.
This class is used to select pairs, triplets...
Definition: Filter.h:34
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
Class to store variables with their name which were sent to the logging service.
REG_MODULE(arichBtest)
Register the Module.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.