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