Belle II Software development
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/* Own header. */
10#include <ecl/modules/eclDQMInjection/eclDQMInjection.h>
11
12/* ECL headers. */
13#include <ecl/dataobjects/ECLDsp.h>
14#include <ecl/utility/ECLDspUtilities.h>
15
16/* Boost headers. */
17#include <boost/format.hpp>
18#include <boost/range/combine.hpp>
19
20/* ROOT headers. */
21#include <TDirectory.h>
22
23/* C++ headers. */
24#include <stdexcept>
25
26using namespace std;
27using namespace Belle2;
28using namespace Belle2::ECL;
29
30//-----------------------------------------------------------------
31// Register the Module
32//-----------------------------------------------------------------
33REG_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");
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 acquisition 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
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,
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,
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 hInjkickTimeShift[0] = new TH2F("ECLInjkickTimeShiftLER",
91 "LER Injection peak position in ECL data;"
92 "Time within beam cycle [ADC ticks];"
93 "Inj peak position [ADC ticks]",
94 18, 0, 18, 16, 0, 16);
95 hInjkickTimeShift[1] = new TH2F("ECLInjkickTimeShiftHER",
96 "HER Injection peak position in ECL data;"
97 "Time within beam cycle [ADC ticks];"
98 "Inj peak position [ADC ticks]",
99 18, 0, 18, 16, 0, 16);
100
101 //== Fill h_ped_peak vector
102
104 1.0, 1.5, 2.0, 4, 6, 8, 10
105 };
106 int ped_peak_range_count = m_ped_peak_range.size() - 1;
107
108 static const std::string part_names[] = {"fwd", "bar", "bwd"};
109 static const std::string title_suffix[] = {
110 "in fwd endcap", "in barrel", "in bwd endcap"
111 };
112
113 for (int ler_her = 0; ler_her < 2; ler_her++) {
114 std::string ring_name = (ler_her == 0) ? "LER" : "HER";
115 for (int part = 0; part < 3; part++) {
116 std::string suffix = title_suffix[part];
117 for (int i = 0; i < ped_peak_range_count; i++) {
118 float min_time = m_ped_peak_range[i];
119 float max_time = m_ped_peak_range[i + 1];
120 std::string name, title;
121 name = str(boost::format("ped_peak_%s_%s_%d") %
122 ring_name % part_names[part] % i);
123 title = str(boost::format("Peak height %.1f-%.1f ms after %s inj %s") %
124 min_time % max_time % ring_name % suffix);
125
126 auto h = new TH1F(name.c_str(), title.c_str(), 300, 0.0, 0.3);
127 h->GetXaxis()->SetTitle("Peak height in first 16 points [GeV]");
128
129 h_ped_peak.push_back(h);
130 }
131 }
132 }
133
134 // Initialize coefficients used by pedestalFit function
136
137 // cd back to root directory
138 oldDir->cd();
139}
140
142{
143 REG_HISTOGRAM
144 m_rawTTD.isOptional();
145 m_storeHits.isRequired(m_ECLDigitsName);
146 m_ECLTrigs.isOptional();
147 m_ECLDsps.isOptional();
148 m_l1Trigger.isOptional();
149
150 if (!mapper.initFromDB()) B2FATAL("ECL Display:: Can't initialize eclChannelMapper");
151
152 v_totalthrApsd.resize((m_calibrationThrApsd->getCalibVector()).size());
153 for (size_t i = 0; i < v_totalthrApsd.size(); i++) v_totalthrApsd[i] = (int)(m_calibrationThrApsd->getCalibVector())[i];
154}
155
157{
158 // Assume that everthing is non-yero ;-)
159 hHitsAfterInjLER->Reset();
160 hHitsAfterInjHER->Reset();
161 hEHitsAfterInjLER->Reset();
162 hEHitsAfterInjHER->Reset();
163 hBurstsAfterInjLER->Reset();
164 hBurstsAfterInjHER->Reset();
165 hEBurstsAfterInjLER->Reset();
166 hEBurstsAfterInjHER->Reset();
167 hVetoAfterInjLER->Reset();
168 hVetoAfterInjHER->Reset();
169 hOccAfterInjHER->Reset();
170 hOccAfterInjLER->Reset();
171}
172
174{
175 bool bhatrig = false;
176
177 if (m_l1Trigger.isValid() && m_DPHYTTYP) bhatrig = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_DPHY;
178 else if (m_l1Trigger.isValid() && !m_DPHYTTYP) {
179 try { bhatrig = m_l1Trigger->testInput("bha_delay"); }
180 catch (const std::exception&) { bhatrig = false; }
181 }
182
183 if (m_eventmetadata.isValid() && m_eventmetadata->getErrorFlag() != 0x10) {
184 m_iEvent = m_eventmetadata->getEvent();
185 } else m_iEvent = -1;
186
187 int amps[ECLElementNumbers::c_NCrystals] = {};
188 unsigned int ECLDigitsAboveThr = 0; // Threshold is set to 20 MeV
189 unsigned int ECLDigitsAboveThr1MeV = 0;
190 for (const auto& aECLDigit : m_storeHits) {
191 int amp = aECLDigit.getAmp();
192 amps[aECLDigit.getCellId() - 1] = amp;
193 if (amp > m_ECLThresholdforVetoTuning) ECLDigitsAboveThr += 1;
194 if (amp > 20) ECLDigitsAboveThr1MeV += 1;
195 }
196
197 int discarded_wfs = 0;
198 for (auto& aECLTrig : m_ECLTrigs) {
199 int crate = aECLTrig.getTrigId();
200 int suppress = aECLTrig.getBurstSuppressionMask();
201 int shaper_pos = 0;
202 while (suppress) {
203 shaper_pos ++;
204 bool shaper_bit = suppress & 1;
205 if (shaper_bit) {
206 if (m_iEvent % 1000 == 999 || (m_l1Trigger.isValid() && m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_RAND) ||
207 (m_l1Trigger.isValid() && bhatrig)) {
208 for (int channel_pos = 0; channel_pos < 16; channel_pos ++) {
209 if (mapper.getCellId(crate, shaper_pos, channel_pos) > 0) discarded_wfs += 1;
210 }
211 } else {
212 for (int channel_pos = 0; channel_pos < 16; channel_pos ++) {
213 int cid = mapper.getCellId(crate, shaper_pos, channel_pos);
214 if (cid > 0 && amps[cid - 1] >= (v_totalthrApsd[cid - 1] / 4 * 4)) discarded_wfs += 1;
215 }
216 }
217 }
218 suppress >>= 1;
219 }
220 }
221
222 for (auto& it : m_rawTTD) {
223 B2DEBUG(29, "TTD FTSW : " << hex << it.GetTTUtime(0) << " " << it.GetTTCtime(0) << " EvtNr " << it.GetEveNo(0) << " Type " <<
224 (it.GetTTCtimeTRGType(0) & 0xF) << " TimeSincePrev " << it.GetTimeSincePrevTrigger(0) << " TimeSinceInj " <<
225 it.GetTimeSinceLastInjection(0) << " IsHER " << it.GetIsHER(0) << " Bunch " << it.GetBunchNumber(0));
226
227 // get last injection time
228 auto difference = it.GetTimeSinceLastInjection(0);
229 // check time overflow, too long ago
230 if (difference == 0x7FFFFFFF) continue;
231
232 unsigned int all = m_storeHits.getEntries();
233 float diff2 = difference / 127.; // 127MHz clock ticks to us, inexact rounding
234
235 // Time within beam revolution (in 127 MHz ticks)
236 int time_within_cycle = difference % 1280;
237 // Time within beam revolution (in microseconds)
238 double time_in_cycle_us = time_within_cycle / 127.;
239 // Time within beam revolution (in ADC ticks)
240 // https://confluence.desy.de/pages/viewpage.action?spaceKey=BI&title=ECL+Technical+Notes
241 int time_within_cycle_adc_ticks = (1280 - time_within_cycle) / 72;
242
243 int is_her = it.GetIsHER(0);
244
245 if (is_her < 0 || is_her > 1) continue;
246
247 if (is_her) {
248 hHitsAfterInjHER->Fill(diff2, all);
249 hEHitsAfterInjHER->Fill(diff2);
250 hBurstsAfterInjHER->Fill(diff2, discarded_wfs);
251 hEBurstsAfterInjHER->Fill(diff2);
252 hVetoAfterInjHER->Fill(diff2, time_in_cycle_us, ECLDigitsAboveThr);
253 if (all > 0) hOccAfterInjHER->Fill(diff2, ECLDigitsAboveThr1MeV * 100.0 / ECLElementNumbers::c_NCrystals);
254 } else {
255 hHitsAfterInjLER->Fill(diff2, all);
256 hEHitsAfterInjLER->Fill(diff2);
257 hBurstsAfterInjLER->Fill(diff2, discarded_wfs);
258 hEBurstsAfterInjLER->Fill(diff2);
259 hVetoAfterInjLER->Fill(diff2, time_in_cycle_us, ECLDigitsAboveThr);
260 if (all > 0) hOccAfterInjLER->Fill(diff2, ECLDigitsAboveThr1MeV * 100.0 / ECLElementNumbers::c_NCrystals);
261 }
262
263 //== Filling h_ped_peak histograms
264 int range_count = m_ped_peak_range.size() - 1;
265 if (diff2 < m_ped_peak_range[range_count] * 1000) {
266 //== Identify which histogram to fill (according to inj time range)
267 int range_id;
268 for (range_id = 0; range_id < range_count; range_id++) {
269 // Converting from ms to us
270 float min_time = m_ped_peak_range[range_id ] * 1000;
271 float max_time = m_ped_peak_range[range_id + 1] * 1000;
272 if (diff2 > min_time && diff2 < max_time) break;
273 }
274 //== Find pedestal peaks in all available waveforms
275 if (range_id < range_count) {
276 for (auto& aECLDsp : m_ECLDsps) {
277 auto result = ECLDspUtilities::pedestalFit(aECLDsp.getDspA());
278
279 //== Identify which histogram to fill (HER/LER,FWD/BAR/BWD)
280 int cid = aECLDsp.getCellId();
281 int part_id = 0; // forward endcap
282 if (ECLElementNumbers::isBarrel(cid)) part_id = 1; // barrel
283 if (ECLElementNumbers::isBackward(cid)) part_id = 2; // backward endcap
284
285 int hist_id = is_her * 3 * range_count + part_id * range_count + range_id;
286 // NOTE: We are using the approximate conversion to energy here.
287 // (20'000 ADC counts ~= 1 GeV)
288 h_ped_peak[hist_id]->Fill(result.amp / 2e4);
289 }
290 }
291 }
292
293 //== Filling hInjkickTimeShift histograms
294
295 if (diff2 < 10e3) {
296 for (auto& aECLDsp : m_ECLDsps) {
297 int adc[31];
298 aECLDsp.getDspA(adc);
299 // Do a naive estimate of inj peak position by
300 // searching for the maximum ADC sample in the
301 // pedestal part of the waveform.
302 int* ped_max = std::max_element(adc, adc + 16);
303 int* ped_min = std::min_element(adc, adc + 16);
304 // The waveform should have at least ~10 MeV peak amplitude
305 if (*ped_max - *ped_min < 200) continue;
306 int max_ped_id = ped_max - adc;
307 hInjkickTimeShift[is_her]->Fill(time_within_cycle_adc_ticks, max_ped_id);
308 }
309 }
310
311
312 break;
313 }
314}
bool initFromDB()
Initialize channel mapper from the conditions database.
int getCellId(int iCrate, int iShaper, int iChannel)
Get CellId by given crate number, shaper position in the crate and DSP channel number in the shaper.
TH1F * hHitsAfterInjHER
Histogram Hits after HER injection.
TH2F * hVetoAfterInjLER
Histogram Veto tuning w/ ECL hits after LER injection.
TH2F * hInjkickTimeShift[2]
Histograms to determine injkick signal time offset for LER/HER injections.
TH1F * hEHitsAfterInjLER
Histogram for Nr Entries (=Triggrs) for normalization after LER injection.
void initialize() override final
initialize function
TH1F * hEHitsAfterInjHER
Histogram for Nr Entries (=Triggrs) for normalization after HER injection.
TH1F * hHitsAfterInjLER
Histogram Hits after LER injection.
bool m_DPHYTTYP
Flag to select events triggered by delayed bhabha.
TH2F * hVetoAfterInjHER
Histogram Veto tuning w/ ECL hits after HER injection.
TH1F * hBurstsAfterInjLER
Histogram Bursts suppression after LER injection.
TH1F * hEBurstsAfterInjLER
Histogram Bursts suppression for normalization after LER injection.
ECLDQMInjectionModule()
Constructor defining the parameters.
StoreArray< ECLDsp > m_ECLDsps
Input array for ECL waveform data.
StoreObjPtr< TRGSummary > m_l1Trigger
StoreObjPtr TRGSummary
double m_revolutionTime
The beam revolution cycle time in .
StoreArray< RawFTSW > m_rawTTD
Input array for DAQ Status.
void defineHisto() override final
defineHisto function
int m_iEvent
Global event number.
void event() override final
event function
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
StoreArray< ECLDigit > m_storeHits
Input array for ECL Raw Hits.
std::string m_ECLDigitsName
The name of the StoreArray of ECLRawHits to be generated.
std::vector< float > m_ped_peak_range
Injection time range (in ms) for h_ped_peak histograms.
double m_ECLThresholdforVetoTuning
ECL threshold for injection veto tuning, ADC channels.
TH1F * hBurstsAfterInjHER
Histogram Bursts suppression after HER injection.
TH1F * hEBurstsAfterInjHER
Histogram Bursts suppression for normalization after HER injection.
void beginRun() override final
beginRun function
StoreArray< ECLTrig > m_ECLTrigs
Input array for ECL burst suppresions.
StoreObjPtr< EventMetaData > m_eventmetadata
StoreObjPtr EventMetaData.
std::vector< int > v_totalthrApsd
Vector to store psd wf amplitude threshold.
std::vector< TH1F * > h_ped_peak
Distribution of pedestal peak (peak in first 16 waveform samples) after HER/LER injection,...
DBObjPtr< ECLCrystalCalib > m_calibrationThrApsd
PSD waveform amplitude threshold.
TH2F * hOccAfterInjLER
Histogram Occupancy after LER injection.
TH2F * hOccAfterInjHER
Histogram Occupancy after HER injection.
ECLChannelMapper mapper
ECL channel mapper.
static ECLPedestalFit pedestalFit(std::vector< int > adc)
Fit pedestal part of the signal waveform (first 16 samples) This method will fit the first 16 samples...
static void initPedestalFit()
Load DSP coefficients used in the pedestal fit function.
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
@ TTYP_DPHY
delayed physics events for background
Definition: TRGSummary.h:65
@ TTYP_RAND
random trigger events
Definition: TRGSummary.h:67
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
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
bool isBarrel(int cellId)
Check whether the crystal is in barrel ECL.
const int c_NCrystals
Number of crystals.
bool isBackward(int cellId)
Check whether the crystal is in backward ECL.
Abstract base class for different kinds of events.
STL namespace.