Belle II Software  release-06-00-14
eclDQMInjection.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 //This module
10 #include <ecl/modules/eclDQMInjection/eclDQMInjection.h>
11 
12 //Boost
13 #include <boost/format.hpp>
14 #include <boost/range/combine.hpp>
15 
16 //ECL
17 #include <ecl/dataobjects/ECLDsp.h>
18 #include <ecl/utility/ECLDspUtilities.h>
19 
20 //ROOT
21 #include "TDirectory.h"
22 
23 //STL
24 #include <stdexcept>
25 
26 using namespace std;
27 using namespace Belle2;
28 using namespace Belle2::ECL;
29 
30 //-----------------------------------------------------------------
31 // Register the Module
32 //-----------------------------------------------------------------
33 REG_MODULE(ECLDQMInjection)
34 
35 //-----------------------------------------------------------------
36 // Implementation
37 //-----------------------------------------------------------------
38 
40  : HistoModule(),
41  m_calibrationThrApsd("ECL_FPGA_StoreWaveform")
42 {
43  //Set module properties
44  setDescription("Monitor Occupancy after Injection");
45  setPropertyFlags(c_ParallelProcessingCertified);
46  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
47  std::string("ECLINJ"));
48  addParam("ECLDigitsName", m_ECLDigitsName, "Name of ECL hits", std::string(""));
49  // BeamRevolutionCycle is set based on 'Timing distribution for the Belle II
50  // data acquistion system'. RF clock of 508 MHz is synchronized to
51  // beam-revolution cycle (5120 RF bunches in one cycle).
52  addParam("BeamRevolutionCycle", m_revolutionTime, "Beam revolution cycle in musec", 5120 / 508.);
53  addParam("ECLThresholdforVetoTuning", m_ECLThresholdforVetoTuning, "ECL Threshold for injection veto tuning, ADC channels", 400.);
54  addParam("DPHYTTYP", m_DPHYTTYP,
55  "Flag to control trigger of delayed bhabha events; 0 - select events by 'bha_delay' trigger bit, 1 - select by TTYP_DPHY", false);
56 
57 }
58 
59 void ECLDQMInjectionModule::defineHisto()
60 {
61  TDirectory* oldDir = gDirectory;
62  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not rely on return value, might be ZERO
63  oldDir->cd(m_histogramDirectoryName.c_str());//changing to the right directory
64 
65  hHitsAfterInjLER = new TH1F("ECLHitsInjLER", "ECLHitsInjLER/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
66  hHitsAfterInjHER = new TH1F("ECLHitsInjHER", "ECLHitsInjHER/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
67  hEHitsAfterInjLER = new TH1F("ECLEHitsInjLER", "ECLEHitsInjLER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0, 20000);
68  hEHitsAfterInjHER = new TH1F("ECLEHitsInjHER", "ECLEHitsInjHER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0, 20000);
69  hBurstsAfterInjLER = new TH1F("ECLBurstsInjLER", "ECLBurstsInjLER/Time;Time in #mus;Count/Time (1 #mus bins)", 20000, 0, 20000);
70  hBurstsAfterInjHER = new TH1F("ECLBurstsInjHER", "ECLBurstsInjHER/Time;Time in #mus;Count/Time (1 #mus bins)", 20000, 0, 20000);
71  hEBurstsAfterInjLER = new TH1F("ECLEBurstsInjLER", "ECLEBurstsInjLER/Time;Time in #mus;Triggers/Time (1 #mus bins)", 20000, 0,
72  20000);
73  hEBurstsAfterInjHER = new TH1F("ECLEBurstsInjHER", "ECLEBurstsInjHER/Time;Time in #mus;Triggers/Time (1 #mus bins)", 20000, 0,
74  20000);
75  hVetoAfterInjLER = new TH2F("ECLVetoAfterInjLER",
76  "ECL Hits for LER veto tuning (E > 20 MeV);"
77  "Time since last injection in #mus;Time within beam cycle in #mus", 500, 0, 30000, 100, 0,
78  m_revolutionTime);
79  hVetoAfterInjHER = new TH2F("ECLVetoAfterInjHER",
80  "ECL Hits for HER veto tuning (E > 20 MeV);"
81  "Time since last injection in #mus;Time within beam cycle in #mus", 500, 0, 30000, 100, 0,
82  m_revolutionTime);
83  hOccAfterInjLER = new TH2F("ECLOccAfterInjLER",
84  "ECL Occupancy after LER injection (E > 1 MeV); Time since last injection in #mus;Occupancy (Nhits/8736) [%]",
85  100, 0, 20000, 98, 2, 100);
86  hOccAfterInjHER = new TH2F("ECLOccAfterInjHER",
87  "ECL Occupancy after HER injection (E > 1 MeV); Time since last injection in #mus;Occupancy (Nhits/8736) [%]",
88  100, 0, 20000, 98, 2, 100);
89 
90 
91  //== Fill h_ped_peak vector
92 
93  m_ped_peak_range = {
94  1.0, 1.5, 2.0, 4, 6, 8, 10
95  };
96  int ped_peak_range_count = m_ped_peak_range.size() - 1;
97 
98  static const std::string part_names[] = {"fwd", "bar", "bwd"};
99  static const std::string title_suffix[] = {
100  "in fwd endcap", "in barrel", "in bwd endcap"
101  };
102 
103  for (int ler_her = 0; ler_her < 2; ler_her++) {
104  std::string ring_name = (ler_her == 0) ? "LER" : "HER";
105  for (int part = 0; part < 3; part++) {
106  std::string suffix = title_suffix[part];
107  for (int i = 0; i < ped_peak_range_count; i++) {
108  float min_time = m_ped_peak_range[i];
109  float max_time = m_ped_peak_range[i + 1];
110  std::string name, title;
111  name = str(boost::format("ped_peak_%s_%s_%d") %
112  ring_name % part_names[part] % i);
113  title = str(boost::format("Peak height %.1f-%.1f ms after %s inj %s") %
114  min_time % max_time % ring_name % suffix);
115 
116  auto h = new TH1F(name.c_str(), title.c_str(), 300, 0.0, 0.3);
117  h->GetXaxis()->SetTitle("Peak height in first 16 points [GeV]");
118 
119  h_ped_peak.push_back(h);
120  }
121  }
122  }
123 
124  // Initialize coefficients used by pedestalFit function
125  ECLDspUtilities::initPedestalFit();
126 
127  // cd back to root directory
128  oldDir->cd();
129 }
130 
131 void ECLDQMInjectionModule::initialize()
132 {
133  REG_HISTOGRAM
134  m_rawTTD.isOptional();
135  m_storeHits.isRequired(m_ECLDigitsName);
136  m_ECLTrigs.isOptional();
137  m_ECLDsps.isOptional();
138  m_l1Trigger.isOptional();
139 
140  if (!mapper.initFromDB()) B2FATAL("ECL Display:: Can't initialize eclChannelMapper");
141 
142  v_totalthrApsd.resize((m_calibrationThrApsd->getCalibVector()).size());
143  for (size_t i = 0; i < v_totalthrApsd.size(); i++) v_totalthrApsd[i] = (int)(m_calibrationThrApsd->getCalibVector())[i];
144 }
145 
146 void ECLDQMInjectionModule::beginRun()
147 {
148  // Assume that everthing is non-yero ;-)
149  hHitsAfterInjLER->Reset();
150  hHitsAfterInjHER->Reset();
151  hEHitsAfterInjLER->Reset();
152  hEHitsAfterInjHER->Reset();
153  hBurstsAfterInjLER->Reset();
154  hBurstsAfterInjHER->Reset();
155  hEBurstsAfterInjLER->Reset();
156  hEBurstsAfterInjHER->Reset();
157  hVetoAfterInjLER->Reset();
158  hVetoAfterInjHER->Reset();
159  hOccAfterInjHER->Reset();
160  hOccAfterInjLER->Reset();
161 }
162 
163 void ECLDQMInjectionModule::event()
164 {
165  bool bhatrig = false;
166 
167  if (m_l1Trigger.isValid() && m_DPHYTTYP) bhatrig = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_DPHY;
168  else if (m_l1Trigger.isValid() && !m_DPHYTTYP) {
169  try { bhatrig = m_l1Trigger->testInput("bha_delay"); }
170  catch (const std::exception&) { bhatrig = false; }
171  }
172 
173  if (m_eventmetadata.isValid() && m_eventmetadata->getErrorFlag() != 0x10) {
174  m_iEvent = m_eventmetadata->getEvent();
175  } else m_iEvent = -1;
176  int discarded_wfs = 0;
177  for (auto& aECLTrig : m_ECLTrigs) {
178  int crate = aECLTrig.getTrigId();
179  int suppress = aECLTrig.getBurstSuppressionMask();
180  int shaper_pos = 0;
181  while (suppress) {
182  shaper_pos ++;
183  bool shaper_bit = suppress & 1;
184  if (shaper_bit) {
185  if (m_iEvent % 1000 == 999 || (m_l1Trigger.isValid() && m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_RAND) ||
186  (m_l1Trigger.isValid() && bhatrig)) {
187  for (int channel_pos = 0; channel_pos < 16; channel_pos ++) {
188  if (mapper.getCellId(crate, shaper_pos, channel_pos) > 0) discarded_wfs += 1;
189  }
190  } else {
191  for (auto& aECLDigit : m_storeHits) {
192  if (crate == mapper.getCrateID(aECLDigit.getCellId()) && shaper_pos == mapper.getShaperPosition(aECLDigit.getCellId()) &&
193  aECLDigit.getAmp() >= (v_totalthrApsd[aECLDigit.getCellId() - 1] / 4 * 4)) discarded_wfs += 1;
194  }
195  }
196  }
197  suppress >>= 1;
198  }
199  }
200 
201  unsigned int ECLDigitsAboveThr = 0; // Threshold is set to 20 MeV
202  unsigned int ECLDigitsAboveThr1MeV = 0;
203  for (auto& aECLDigit : m_storeHits) {
204  if (aECLDigit.getAmp() > m_ECLThresholdforVetoTuning) ECLDigitsAboveThr += 1;
205  if (aECLDigit.getAmp() > 20) ECLDigitsAboveThr1MeV += 1;
206  }
207 
208  for (auto& it : m_rawTTD) {
209  B2DEBUG(29, "TTD FTSW : " << hex << it.GetTTUtime(0) << " " << it.GetTTCtime(0) << " EvtNr " << it.GetEveNo(0) << " Type " <<
210  (it.GetTTCtimeTRGType(0) & 0xF) << " TimeSincePrev " << it.GetTimeSincePrevTrigger(0) << " TimeSinceInj " <<
211  it.GetTimeSinceLastInjection(0) << " IsHER " << it.GetIsHER(0) << " Bunch " << it.GetBunchNumber(0));
212 
213  // get last injection time
214  auto difference = it.GetTimeSinceLastInjection(0);
215  // check time overflow, too long ago
216  if (difference != 0x7FFFFFFF) {
217  unsigned int all = m_storeHits.getEntries();
218  float diff2 = difference / 127.; // 127MHz clock ticks to us, inexact rounding
219  int is_her = it.GetIsHER(0);
220  if (is_her) {
221  hHitsAfterInjHER->Fill(diff2, all);
222  hEHitsAfterInjHER->Fill(diff2);
223  hBurstsAfterInjHER->Fill(diff2, discarded_wfs);
224  hEBurstsAfterInjHER->Fill(diff2);
225  hVetoAfterInjHER->Fill(diff2, diff2 - int(diff2 / m_revolutionTime)*m_revolutionTime, ECLDigitsAboveThr);
226  if (all > 0) hOccAfterInjHER->Fill(diff2, ECLDigitsAboveThr1MeV / 8736.*100.);
227  } else {
228  hHitsAfterInjLER->Fill(diff2, all);
229  hEHitsAfterInjLER->Fill(diff2);
230  hBurstsAfterInjLER->Fill(diff2, discarded_wfs);
231  hEBurstsAfterInjLER->Fill(diff2);
232  hVetoAfterInjLER->Fill(diff2, diff2 - int(diff2 / m_revolutionTime)*m_revolutionTime, ECLDigitsAboveThr);
233  if (all > 0) hOccAfterInjLER->Fill(diff2, ECLDigitsAboveThr1MeV / 8736.*100.);
234  }
235 
236  //== Filling h_ped_peak histograms
237  int range_count = m_ped_peak_range.size() - 1;
238  if (diff2 < m_ped_peak_range[range_count] * 1000) {
239  //== Identify which histogram to fill (according to inj time range)
240  int range_id;
241  for (range_id = 0; range_id < range_count; range_id++) {
242  // Converting from ms to us
243  float min_time = m_ped_peak_range[range_id ] * 1000;
244  float max_time = m_ped_peak_range[range_id + 1] * 1000;
245  if (diff2 > min_time && diff2 < max_time) break;
246  }
247  //== Find pedestal peaks in all available waveforms
248  if (range_id < range_count) {
249  for (auto& aECLDsp : m_ECLDsps) {
250  auto result = ECLDspUtilities::pedestalFit(aECLDsp.getDspA());
251 
252  //== Identify which histogram to fill (HER/LER,FWD/BAR/BWD)
253  int cid = aECLDsp.getCellId();
254  int part_id = 0; // forward endcap
255  if (cid >= 1153) part_id = 1; // barrel
256  if (cid >= 7777) part_id = 2; // backward endcap
257 
258  int hist_id = is_her * 3 * range_count + part_id * range_count + range_id;
259  // NOTE: We are using the approximate conversion to energy here.
260  // (20'000 ADC counts ~= 1 GeV)
261  h_ped_peak[hist_id]->Fill(result.amp / 2e4);
262  }
263  }
264  }
265  }
266 
267  break;
268  }
269 }
The ECL Occ after Injection DQM module.
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.