Belle II Software development
eclDQM.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/eclDQM/eclDQM.h>
11
12/* ECL headers. */
13#include <ecl/dataobjects/ECLCalDigit.h>
14#include <ecl/dataobjects/ECLDigit.h>
15#include <ecl/dataobjects/ECLDsp.h>
16#include <ecl/dataobjects/ECLElementNumbers.h>
17#include <ecl/dataobjects/ECLTrig.h>
18#include <ecl/geometry/ECLGeometryPar.h>
19#include <ecl/mapper/ECLChannelMapper.h>
20
21/* Basf2 headers. */
22#include <framework/core/HistoModule.h>
23#include <framework/dataobjects/EventMetaData.h>
24#include <framework/gearbox/Unit.h>
25#include <framework/logging/Logger.h>
26#include <mdst/dataobjects/TRGSummary.h>
27
28/* Boost headers. */
29#include <boost/format.hpp>
30#include <boost/range/combine.hpp>
31
32/* ROOT headers. */
33#include <TDirectory.h>
34#include <TH1F.h>
35#include <TH2F.h>
36#include <TProfile.h>
37
38/* C++ headers. */
39#include <iostream>
40#include <iterator>
41#include <cmath>
42#include <stdexcept>
43
44//NAMESPACE(S)
45using namespace Belle2;
46using namespace ECL;
47
48REG_MODULE(ECLDQM);
49
51 : HistoModule(),
52 m_calibrationThrApsd("ECL_FPGA_StoreWaveform")
53{
54 //Set module properties.
56 "Primary module for ECL Data Quality Monitor.\n"
57 "This module provides a large set of low-level histograms: occupancy, time distribution, number of saved waveforms, etc.");
58 setPropertyFlags(c_ParallelProcessingCertified); // specify parallel processing.
59
60 m_WaveformOption = {"psd", "logic", "rand", "dphy", "other"};
61
62 addParam("histogramDirectoryName", m_histogramDirectoryName,
63 "histogram directory in ROOT file", std::string("ECL"));
64 addParam("EnergyUpperThr", m_EnergyUpperThr, "Upper threshold of energy deposition in event, [GeV]", 20.0 * Belle2::Unit::GeV);
65 addParam("HitThresholds", m_HitThresholds, "Thresholds to display hit occupancy, MeV", std::vector<double> {0, 5, 10, 50});
66 addParam("TotalEnergyThresholds", m_TotalEnergyThresholds, "Thresholds to display total energy, MeV", std::vector<double> {0, 5, 7});
67 addParam("TimingThresholds", m_TimingThresholds, "Thresholds (MeV) to display ECL timing", std::vector<double> {5, 10, 50});
68 addParam("HitNumberUpperlimits", m_HitNumberUpperLimits,
69 "Upper limit (# of hits) to display hit multiplicity", std::vector<double> {10000, 1000, 700, 200});
70 addParam("WaveformOption", m_WaveformOption, "Option (all,psd,logic,rand,dphy,other) to display waveform flow",
72 addParam("DPHYTTYP", m_DPHYTTYP,
73 "Flag to control trigger of delayed bhabha events; 0 - select events by 'bha_delay' trigger bit, 1 - select by TTYP_DPHY", false);
74}
75
77{
78}
79
80
82{
83 TDirectory* oldDir = gDirectory;
84
85 // Create a separate histogram directory and cd into it.
86
87 TDirectory* dirDAQ = dynamic_cast<TDirectory*>(oldDir->Get(m_histogramDirectoryName.c_str()));
88 if (!dirDAQ) dirDAQ = oldDir->mkdir(m_histogramDirectoryName.c_str());
89 dirDAQ->cd();
90
91 //1D histograms creation.
92 h_evtot = new TH1F("event", "Total event bank", 1, 0, 1);
93
94 h_quality = new TH1F("quality", "Fit quality flag. 0-good, 1-integer overflow, 2-low amplitude, 3-bad chi2", 4, 0, 4);
95 h_quality->GetXaxis()->SetTitle("Quality flag");
96 h_quality->GetYaxis()->SetTitle("ECL hits count");
97 h_quality->SetFillColor(kPink - 4);
98
99 h_quality_other = new TH1F("quality_other", "Fit quality flag for unexpectedly saved waveforms", 4, 0, 4);
100 h_quality_other->GetXaxis()->SetTitle("Quality flag. 0-good,1-int overflow,2-low amplitude,3-bad chi2");
101 h_quality_other->SetFillColor(kPink - 4);
102
103 h_bad_quality = new TH1F("bad_quality", "Fraction of hits with bad chi2 (qual=3) and E > 1 GeV vs Cell ID",
105 h_bad_quality->GetXaxis()->SetTitle("Cell ID");
106 h_bad_quality->GetYaxis()->SetTitle("ECL hits count");
107
108 h_trigtag1 = new TH1F("trigtag1", "Consistency b/w global event number and trigger tag. 0-good, 1-DQM error", 2, 0, 2);
109 h_trigtag1->GetXaxis()->SetTitle("Flag value");
110 h_trigtag1->GetYaxis()->SetTitle("Events count");
111 h_trigtag1->SetDrawOption("hist");
112 h_trigtag1->SetFillColor(kPink - 4);
113
114 h_adc_hits = new TH1F("adc_hits", "Fraction of high-energy hits (E > 50 MeV)", 1001, 0, 1.001);
115 h_adc_hits->GetXaxis()->SetTitle("Fraction");
116 h_adc_hits->GetYaxis()->SetTitle("Events count");
117
118 h_time_crate_Thr1GeV_large = new TH1F("time_crate_Thr1GeV_large",
119 "Number of hits with timing outside #pm 100 ns per Crate ID (E > 1 GeV)",
120 52, 1, 53);
121 h_time_crate_Thr1GeV_large->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
122 h_time_crate_Thr1GeV_large->GetYaxis()->SetTitle("ECL hits count");
123
124 for (const auto& id : m_HitThresholds) {
125 std::string h_name, h_title;
126 h_name = str(boost::format("cid_Thr%1%MeV") % id);
127 h_title = str(boost::format("Occupancy per Cell ID (E > %1% MeV)") % id);
128 TH1F* h = new TH1F(h_name.c_str(), h_title.c_str(),
130 h->GetXaxis()->SetTitle("Cell ID");
131 h->GetYaxis()->SetTitle("Occupancy (hits / events_count)");
132 h_cids.push_back(h);
133 }
134
135 for (const auto& id : m_TotalEnergyThresholds) {
136 std::string h_name, h_title;
137 h_name = str(boost::format("edep_Thr%1%MeV") % id);
138 h_title = str(boost::format("Total energy (thr = %1% MeV)") % id);
139 TH1F* h = new TH1F(h_name.c_str(), h_title.c_str(), (int)(100 * m_EnergyUpperThr), 0, m_EnergyUpperThr);
140 h->GetXaxis()->SetTitle("Energy, [GeV]");
141 h_edeps.push_back(h);
142 }
143
144 for (const auto& id : m_TimingThresholds) {
145 std::string h_bar_name, h_bar_title;
146 std::string h_end_name, h_end_title;
147 h_bar_name = str(boost::format("time_barrel_Thr%1%MeV") % id);
148 h_bar_title = str(boost::format("Reconstructed time for ECL barrel (E > %1% MeV)") % id);
149 h_end_name = str(boost::format("time_endcaps_Thr%1%MeV") % id);
150 h_end_title = str(boost::format("Reconstructed time for ECL endcaps (E > %1% MeV)") % id);
151 TH1F* h_time_barrel = new TH1F(h_bar_name.c_str(), h_bar_title.c_str(), 206, -1030, 1030);
152 TH1F* h_time_endcap = new TH1F(h_end_name.c_str(), h_end_title.c_str(), 206, -1030, 1030);
153 h_time_barrel->GetXaxis()->SetTitle("Time, [ns]");
154 h_time_endcap->GetXaxis()->SetTitle("Time, [ns]");
155 h_time_barrels.push_back(h_time_barrel);
156 h_time_endcaps.push_back(h_time_endcap);
157 }
158
159 for (const auto& id : boost::combine(m_HitThresholds, m_HitNumberUpperLimits)) {
160 double id1 = 0, id2 = 0;
161 boost::tie(id1, id2) = id;
162 std::string h_name, h_title;
163 h_name = str(boost::format("ncev_Thr%1%MeV") % id1);
164 h_title = str(boost::format("Number of hits in event (E > %1% MeV)") % id1);
165 TH1F* h = new TH1F(h_name.c_str(), h_title.c_str(), id2, 0, id2);
166 h->GetXaxis()->SetTitle("Number of hits");
167 h_ncevs.push_back(h);
168 }
169
170 for (int i = 0; i < ECL_CRATES; i++) {
171 int crate = i + 1;
172 std::string h_name, h_title;
173 h_name = str(boost::format("time_crate_%1%_Thr1GeV") % (crate));
174 h_title = str(boost::format("Reconstructed time for ECL crate #%1% with E > 1 GeV") % (crate));
175 TH1F* h = new TH1F(h_name.c_str(), h_title.c_str(), 400, -100, 100);
176 h->GetXaxis()->SetTitle("Time [ns]");
177 h_time_crate_Thr1GeV.push_back(h);
178 }
179
180 for (const auto& id : m_WaveformOption) {
181 if (id != "all" && id != "psd" && id != "logic" && id != "rand" && id != "dphy" && id != "other")
182 B2WARNING("Waveform Options are not correctly assigned. They must be 'all', 'psd', 'logic', 'rand', 'dphy', 'other'!");
183 std::string h_title;
184 std::string h_cell_name;
185 if (id == "other") h_title = "Unexpectedly saved waveforms";
186 if (id == "psd") h_title = "#frac{Saved}{Expected} waveforms for high-energy hits (E > 50 MeV)";
187 if (id == "logic") h_title = "#frac{Saved}{Expected} waveforms for every 1000th event";
188 if (id == "rand") h_title = "#frac{Saved}{Expected} waveforms for random trigger events";
189 if (id == "dphy") h_title = "#frac{Saved}{Expected} waveforms for delayed bhabha (DPHY) events";
190 if (id == "all") h_title = "#frac{Saved}{Expected} waveforms for all events";
191 h_cell_name = str(boost::format("wf_cid_%1%") % (id));
192 TH1F* h_cell = new TH1F(h_cell_name.c_str(), h_title.c_str(),
194 h_cell->GetXaxis()->SetTitle("Cell ID");
195 if (id == "psd") {
196 h_cell_psd_norm = new TH1F("psd_cid", "Normalization to psd hits for cid",
198 }
199 if (id == "logic") {
200 h_evtot_logic = new TH1F("event_logic", "Event bank for logic", 1, 0, 1);
201 }
202 if (id == "rand") {
203 h_evtot_rand = new TH1F("event_rand", "Event bank for rand", 1, 0, 1);
204 }
205 if (id == "dphy") {
206 h_evtot_dphy = new TH1F("event_dphy", "Event bank for dphy", 1, 0, 1);
207 }
208 h_cells.push_back(h_cell);
209 }
210
211 //2D histograms creation.
212
213 h_trigtag2_trigid = new TH2F("trigtag2_trigid", "Internal data consistency vs crate. 0-good, 1-data corruption",
214 52, 1, 53, 11, -1, 10);
215 h_trigtag2_trigid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
216 h_trigtag2_trigid->GetYaxis()->SetTitle("Data consistency flag");
217
218 h_pedmean_cellid = new TProfile("pedmean_cellid", "Pedestal vs Cell ID",
220 h_pedmean_cellid->GetXaxis()->SetTitle("Cell ID");
221 h_pedmean_cellid->GetYaxis()->SetTitle("Ped. average (ADC units, #approx 0.05 MeV)");
222
223 h_pedrms_cellid = new TProfile("pedrms_cellid", "Pedestal stddev vs Cell ID",
225 h_pedrms_cellid->GetXaxis()->SetTitle("Cell ID");
226 h_pedrms_cellid->GetYaxis()->SetTitle("Ped. stddev (ADC units, #approx 0.05 MeV)");
227
228 h_pedrms_thetaid = new TProfile("pedrms_thetaid", "Pedestal stddev vs #theta ID",
229 68, 0, 68);
230 h_pedrms_thetaid->GetXaxis()->SetTitle("#theta ID (0-12=FWD, 59-67=BWD endcap)");
231 h_pedrms_thetaid->GetYaxis()->SetTitle("Ped. stddev (ADC units, #approx 0.05 MeV)");
232
233 h_trigtime_trigid = new TH2F("trigtime_trigid", "Trigger time vs Crate ID", 52, 1, 53, 145, 0, 145);
234 h_trigtime_trigid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
235 h_trigtime_trigid->GetYaxis()->SetTitle("Trigger time (only even, 0-142)");
236
237 //cd into parent directory.
238
239 oldDir->cd();
240}
241
243{
244 REG_HISTOGRAM; // required to register histograms to HistoManager.
245 m_ECLDigits.isRequired();
246 m_ECLCalDigits.isOptional();
247 m_ECLTrigs.isOptional();
248 m_ECLDsps.isOptional();
249 m_l1Trigger.isOptional();
250
251 if (!mapper.initFromDB()) B2FATAL("ECL DQM: Can't initialize eclChannelMapper");
252
253 ecltot.resize(m_TotalEnergyThresholds.size());
254 nhits.resize(m_HitNumberUpperLimits.size());
255
257
258 v_totalthrApsd.resize((m_calibrationThrApsd->getCalibVector()).size());
259 for (size_t i = 0; i < v_totalthrApsd.size(); i++) v_totalthrApsd[i] = (int)(m_calibrationThrApsd->getCalibVector())[i];
260}
261
263{
264 h_evtot->Reset();
265 h_evtot_logic->Reset();
266 h_evtot_rand->Reset();
267 h_evtot_dphy->Reset();
268 h_quality->Reset();
269 h_quality_other->Reset();
270 h_bad_quality->Reset();
271 h_trigtag1->Reset();
272 h_adc_hits->Reset();
274 h_cell_psd_norm->Reset();
275 std::for_each(h_cids.begin(), h_cids.end(), [](auto & it) {it->Reset();});
276 std::for_each(h_edeps.begin(), h_edeps.end(), [](auto & it) {it->Reset();});
277 std::for_each(h_time_barrels.begin(), h_time_barrels.end(), [](auto & it) {it->Reset();});
278 std::for_each(h_time_endcaps.begin(), h_time_endcaps.end(), [](auto & it) {it->Reset();});
279 std::for_each(h_ncevs.begin(), h_ncevs.end(), [](auto & it) {it->Reset();});
280 std::for_each(h_cells.begin(), h_cells.end(), [](auto & it) {it->Reset();});
281 for (int i = 0; i < ECL_CRATES; i++) h_time_crate_Thr1GeV[i]->Reset();
282 h_trigtag2_trigid->Reset();
283 h_pedmean_cellid->Reset();
284 h_pedrms_cellid->Reset();
285 h_pedrms_thetaid->Reset();
286 h_trigtime_trigid->Reset();
287}
288
290{
291 int trigtag1 = 0;
292 int NDigits = 0;
293 for (auto& value : ecltot) value = 0;
294 for (auto& value : nhits) value = 0;
295 bool bhatrig = false;
296
297 if (m_l1Trigger.isValid() && m_DPHYTTYP) bhatrig = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_DPHY;
298 else if (m_l1Trigger.isValid() && !m_DPHYTTYP) {
299 try { bhatrig = m_l1Trigger->testInput("bha_delay"); }
300 catch (const std::exception&) { bhatrig = false; }
301 }
302
303 m_iEvent = -1;
304 if (m_eventmetadata.isValid()) {
305 if (m_eventmetadata->getErrorFlag() != 0x10) {
306 m_iEvent = m_eventmetadata->getEvent();
307 h_evtot->Fill(0);
308 for (const auto& id : m_WaveformOption) {
309 if (id == "logic" && m_iEvent % 1000 == 999) h_evtot_logic->Fill(0);
310 if (id == "rand" && isRandomTrigger()) h_evtot_rand->Fill(0);
311 if (id == "dphy" && bhatrig) h_evtot_dphy->Fill(0);
312 }
313 }
314 }
315
316 for (auto& aECLDigit : m_ECLDigits) {
317 int i = aECLDigit.getCellId() - 1;
318 h_quality->Fill(aECLDigit.getQuality()); //Fit quality histogram filling.
319 if (aECLDigit.getAmp() > 2.e04 && aECLDigit.getQuality() == 3) h_bad_quality->Fill(aECLDigit.getCellId());
320 if (aECLDigit.getAmp() >= (v_totalthrApsd[i] / 4 * 4)) NDigits ++;
321 for (const auto& id : m_WaveformOption) {
322 if (id != "psd") continue;
323 else if (id == "psd" && (m_iEvent % 1000 == 999 ||
324 isRandomTrigger() ||
325 bhatrig ||
326 aECLDigit.getAmp() < (v_totalthrApsd[i] / 4 * 4))) continue;
327 h_cell_psd_norm->Fill(aECLDigit.getCellId());
328 }
329 }
330
331 for (auto& aECLTrig : m_ECLTrigs) {
332 double itrg = aECLTrig.getTimeTrig();
333 //trigger time conversion to acceptable units in range (0, ..., 142).
334 //one trigger time clock corresponds to 0.567/144*1000 = 3.93 ns
335 int tg = (int)itrg - 2 * ((int)itrg / 8);
336 h_trigtime_trigid->Fill(aECLTrig.getTrigId(), tg); //Trigger time histogram filling.
337 trigtag1 += aECLTrig.getTrigTag();
338 h_trigtag2_trigid->Fill(aECLTrig.getTrigId(), aECLTrig.getTrigTagQualityFlag()); //Data consistency histogram filling.
339 }
340
341 if (m_ECLTrigs.getEntries() > 0) {
342 int flagtag = 1;
343 trigtag1 /= m_ECLTrigs.getEntries();
344 int compar = (65535 & m_iEvent);
345 if (compar == trigtag1) flagtag = 0;
346 h_trigtag1->Fill(flagtag); //Trigger tag flag #1 histogram filling.
347 }
348
349 for (auto& aECLCalDigit : m_ECLCalDigits) {
350 int cid = aECLCalDigit.getCellId();
351 double energy = aECLCalDigit.getEnergy(); //get calibrated energy.
352 double timing = aECLCalDigit.getTime(); //get calibrated time.
353
354 for (const auto& id : m_HitThresholds) {
355 auto scale = id / 1000.;
356 auto index = std::distance(m_HitThresholds.begin(), std::find(m_HitThresholds.begin(), m_HitThresholds.end(), id));
357 if (energy > scale) {
358 h_cids[index]->Fill(cid);
359 nhits[index] += 1;
360 }
361 }
362
363 for (const auto& id : m_TotalEnergyThresholds) {
364 auto scale = id / 1000.;
365 auto index = std::distance(m_TotalEnergyThresholds.begin(), std::find(m_TotalEnergyThresholds.begin(),
366 m_TotalEnergyThresholds.end(), id));
367 if (energy > scale) ecltot[index] += energy;
368 }
369
370 for (const auto& id : m_TimingThresholds) {
371 auto scale = id / 1000.;
372 auto index = std::distance(m_TimingThresholds.begin(), std::find(m_TimingThresholds.begin(), m_TimingThresholds.end(), id));
373 if (energy > scale) {
374 if (cid > ECL_FWD_CHANNELS && cid <= ECL_FWD_CHANNELS + ECL_BARREL_CHANNELS) h_time_barrels[index]->Fill(timing);
375 else h_time_endcaps[index]->Fill(timing);
376 }
377 }
378
379 if (energy > 1.000 && std::abs(timing) < 100.) h_time_crate_Thr1GeV[mapper.getCrateID(cid) - 1]->Fill(timing);
380 if (energy > 1.000 && std::abs(timing) > 100.) h_time_crate_Thr1GeV_large->Fill(mapper.getCrateID(cid));
381 }
382
383 for (auto& h : h_edeps) {
384 auto index = std::distance(h_edeps.begin(), std::find(h_edeps.begin(), h_edeps.end(), h));
385 h->Fill(ecltot[index]);
386 }
387
388 for (auto& h : h_ncevs) {
389 auto index = std::distance(h_ncevs.begin(), std::find(h_ncevs.begin(), h_ncevs.end(), h));
390 h->Fill(nhits[index]);
391 }
392
393 for (auto& aECLDsp : m_ECLDsps) {
394 int i = aECLDsp.getCellId() - 1; //get number of Cell ID in m_DspArray.
395 aECLDsp.getDspA(m_DspArray[i]);
396 m_PedestalMean[i] = 0;
397 m_PedestalRms[i] = 0;
398
399 for (int j = 0; j < 16; j++) m_PedestalMean[i] += m_DspArray[i][j];
400 m_PedestalMean[i] /= 16;
401 h_pedmean_cellid->Fill(aECLDsp.getCellId(), m_PedestalMean[i]); //Pedestal Avg histogram filling.
402
403 for (int j = 0; j < 16; j++) m_PedestalRms[i] += pow(m_DspArray[i][j] - m_PedestalMean[i], 2);
404 m_PedestalRms[i] = sqrt(m_PedestalRms[i] / 15.);
405 h_pedrms_cellid->Fill(aECLDsp.getCellId(), m_PedestalRms[i]); //Pedestal stddev histogram filling.
406 m_geom->Mapping(i);
408
409 ECLDigit* aECLDigit = ECLDigit::getByCellID(aECLDsp.getCellId());
410
411 for (const auto& id : m_WaveformOption) {
412 auto index = std::distance(m_WaveformOption.begin(), std::find(m_WaveformOption.begin(), m_WaveformOption.end(), id));
413 if (id != "all" && id != "psd" && id != "logic" && id != "rand" && id != "dphy" && id != "other") continue;
414 else if (id == "psd" && (m_iEvent % 1000 == 999 || isRandomTrigger() || bhatrig ||
415 !aECLDigit || aECLDigit->getAmp() < (v_totalthrApsd[i] / 4 * 4))) continue;
416 else if (id == "logic" && m_iEvent % 1000 != 999) continue;
417 else if (id == "rand" && (m_iEvent % 1000 == 999 || !isRandomTrigger())) continue;
418 else if (id == "dphy" && (m_iEvent % 1000 == 999 || !bhatrig)) continue;
419 else if (id == "other" && (m_iEvent % 1000 == 999 || isRandomTrigger() || bhatrig ||
420 (aECLDigit && aECLDigit->getAmp() >= (v_totalthrApsd[i] / 4 * 4)))) continue;
421 h_cells[index]->Fill(aECLDsp.getCellId());
422 if (id == "other" && aECLDigit) h_quality_other->Fill(aECLDigit->getQuality());
423 }
424 }
425 if (m_ECLDigits.getEntries() > 0)
426 h_adc_hits->Fill((double)NDigits / (double)m_ECLDigits.getEntries()); //Fraction of high-energy hits
427}
428
430{
431}
432
433
435{
436}
437
439{
440 if (!m_l1Trigger.isValid()) return false;
441 return m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_RAND ||
443}
444
double m_PedestalMean[ECLElementNumbers::c_NCrystals]
Pedestal average values.
Definition: eclDQM.h:123
TH1F * h_trigtag1
Histogram: Trigger tag flag #1.
Definition: eclDQM.h:142
std::vector< double > nhits
Container for channel multiplicity.
Definition: eclDQM.h:115
TH1F * h_bad_quality
Histogram: Cell IDs w/ bad fit quality flag.
Definition: eclDQM.h:140
double m_PedestalRms[ECLElementNumbers::c_NCrystals]
Pedestal rms error values.
Definition: eclDQM.h:125
bool m_DPHYTTYP
Flag to select events triggered by delayed bhabha.
Definition: eclDQM.h:100
std::vector< double > m_TimingThresholds
Parameters for timing histograms.
Definition: eclDQM.h:107
TH1F * h_adc_hits
Histogram: Fraction of digits above ADC threshold.
Definition: eclDQM.h:144
StoreArray< ECLDsp > m_ECLDsps
StoreArray ECLDsp.
Definition: eclDQM.h:85
virtual void initialize() override
Initialize the module.
Definition: eclDQM.cc:242
virtual ~ECLDQMModule()
Destructor.
Definition: eclDQM.cc:76
StoreArray< ECLCalDigit > m_ECLCalDigits
StoreArray ECLCalDigit.
Definition: eclDQM.h:89
StoreObjPtr< TRGSummary > m_l1Trigger
StoreObjPtr TRGSummary
Definition: eclDQM.h:79
TH1F * h_cell_psd_norm
Histogram: Normalize to psd hits for CellID.
Definition: eclDQM.h:161
virtual void event() override
Event processor.
Definition: eclDQM.cc:289
TH2F * h_trigtime_trigid
Histogram: Trigger time vs.
Definition: eclDQM.h:165
ECLDQMModule()
< derived from HistoModule class.
Definition: eclDQM.cc:50
std::vector< TH1F * > h_edeps
Histogram vector: Total energy.
Definition: eclDQM.h:151
virtual void endRun() override
Call when a run ends.
Definition: eclDQM.cc:429
std::vector< TH1F * > h_ncevs
Histogram vector: Channel multiplicity.
Definition: eclDQM.h:157
virtual void terminate() override
Terminate.
Definition: eclDQM.cc:434
TH1F * h_evtot
Histogram: Total event no (auxiliary) to normalize hit map .
Definition: eclDQM.h:128
int m_iEvent
Global event number.
Definition: eclDQM.h:94
ECL::ECLGeometryPar * m_geom
Geometry.
Definition: eclDQM.h:75
StoreArray< ECLDigit > m_ECLDigits
StoreArray ECLDigit.
Definition: eclDQM.h:83
std::vector< double > ecltot
Container for energy.
Definition: eclDQM.h:113
std::vector< TH1F * > h_time_barrels
Histogram vector: Reconstructed time for barrel.
Definition: eclDQM.h:153
double m_EnergyUpperThr
Upper threshold of energy deposition in event, [GeV].
Definition: eclDQM.h:98
TProfile * h_pedmean_cellid
Histogram: Pedestal Average vs.
Definition: eclDQM.h:169
TH1F * h_evtot_rand
Histogram: Event no for rand (auxiliary) to normalize rand waveform flow.
Definition: eclDQM.h:132
std::vector< TH1F * > h_time_crate_Thr1GeV
Histogram vector: Reconstructed signal time for all ECL crates above the threshold = 1 GeV.
Definition: eclDQM.h:163
TH1F * h_time_crate_Thr1GeV_large
Histogram: Entries with crate time offsets > 100 ns (E > 1 GeV).
Definition: eclDQM.h:146
virtual void beginRun() override
Call when a run begins.
Definition: eclDQM.cc:262
std::string m_histogramDirectoryName
Histogram directory in ROOT file.
Definition: eclDQM.h:96
ECL::ECLChannelMapper mapper
ECL channel mapper.
Definition: eclDQM.h:81
TH1F * h_evtot_logic
Histogram: Event no for logic (auxiliary) to normalize logic waveform flow.
Definition: eclDQM.h:130
std::vector< double > m_HitNumberUpperLimits
Parameters for number of hits histograms.
Definition: eclDQM.h:109
std::vector< double > m_TotalEnergyThresholds
Parameters for histograms w/ total energy.
Definition: eclDQM.h:105
TProfile * h_pedrms_thetaid
Histogram: Pedestal rms error vs.
Definition: eclDQM.h:173
std::vector< std::string > m_WaveformOption
Parameters for waveform histograms.
Definition: eclDQM.h:111
int m_DspArray[ECLElementNumbers::c_NCrystals][31]
WF sampling points for digit array.
Definition: eclDQM.h:121
TH1F * h_evtot_dphy
Histogram: Event no for dphy (auxiliary) to normalize dphy waveform flow.
Definition: eclDQM.h:134
std::vector< double > m_HitThresholds
Parameters for hit occ.
Definition: eclDQM.h:103
TH1F * h_quality
Histogram: Fit quality flag (0 - good, 1 - large amplitude, 3 - bad chi2).
Definition: eclDQM.h:136
StoreArray< ECLTrig > m_ECLTrigs
StoreArray ECLTrig.
Definition: eclDQM.h:87
StoreObjPtr< EventMetaData > m_eventmetadata
StoreObjPtr EventMetaData.
Definition: eclDQM.h:77
TH2F * h_trigtag2_trigid
Histogram: Trigger tag flag #2 vs.
Definition: eclDQM.h:167
std::vector< int > v_totalthrApsd
Vector to store psd wf amplitude threshold.
Definition: eclDQM.h:117
TH1F * h_quality_other
Histogram: Fit quality flag for waveform type 'other'.
Definition: eclDQM.h:138
std::vector< TH1F * > h_cids
Histogram vector: Hit map.
Definition: eclDQM.h:149
TProfile * h_pedrms_cellid
Histogram: Pedestal rms error vs.
Definition: eclDQM.h:171
DBObjPtr< ECLCrystalCalib > m_calibrationThrApsd
PSD waveform amplitude threshold.
Definition: eclDQM.h:91
std::vector< TH1F * > h_cells
Histogram vector: Waveforms vs CellID.
Definition: eclDQM.h:159
bool isRandomTrigger()
Definition: eclDQM.cc:438
virtual void defineHisto() override
Function to define histograms.
Definition: eclDQM.cc:81
std::vector< TH1F * > h_time_endcaps
Histogram vector: Reconstructed time for endcaps.
Definition: eclDQM.h:155
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
bool initFromDB()
Initialize channel mapper from the conditions database.
int getCrateID(int iCOPPERNode, int iFINESSE, bool pcie40=false)
Get crate number by given COPPER node number and FINESSE number.
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
void Mapping(int cid)
Mapping theta, phi Id.
int GetThetaID()
Get Theta Id.
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_POIS
poisson random trigger
Definition: TRGSummary.h:73
@ TTYP_RAND
random trigger events
Definition: TRGSummary.h:67
static const double GeV
Standard of [energy, momentum, mass].
Definition: Unit.h:51
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.