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, "Thresholds to display hit occupancy, MeV", std::vector<double> {0, 5, 10, 50});
72 addParam("TotalEnergyThresholds", m_TotalEnergyThresholds, "Thresholds to display total energy, MeV", std::vector<double> {0, 5, 7});
73 addParam("TimingThresholds", m_TimingThresholds, "Thresholds (MeV) to display ECL timing", std::vector<double> {5, 10, 50});
74 addParam("HitNumberUpperlimits", m_HitNumberUpperLimits,
75 "Upper limit (# of hits) to display hit multiplicity", std::vector<double> {10000, 1000, 700, 200});
76 addParam("WaveformOption", m_WaveformOption, "Option (all,psd,logic,rand,dphy,other) to display waveform flow",
78 addParam("DPHYTTYP", m_DPHYTTYP,
79 "Flag to control trigger of delayed bhabha events; 0 - select events by 'bha_delay' trigger bit, 1 - select by TTYP_DPHY", false);
80 addParam("PI0PListName", m_pi0PListName, "Name of the pi0 particle list", std::string("pi0:eclDQM"));
81}
82
84{
85}
86
87
89{
90 TDirectory* oldDir = gDirectory;
91
92 // Create a separate histogram directory and cd into it.
93
94 TDirectory* dirDAQ = dynamic_cast<TDirectory*>(oldDir->Get(m_histogramDirectoryName.c_str()));
95 if (!dirDAQ) dirDAQ = oldDir->mkdir(m_histogramDirectoryName.c_str());
96 dirDAQ->cd();
97
98 //1D histograms creation.
99 h_evtot = new TH1F("event", "Total event bank", 1, 0, 1);
100
101 h_quality = new TH1F("quality", "Fit quality flag. 0-good, 1-integer overflow, 2-low amplitude, 3-bad chi2", 4, 0, 4);
102 h_quality->GetXaxis()->SetTitle("Quality flag");
103 h_quality->GetYaxis()->SetTitle("ECL hits count");
104 h_quality->SetFillColor(kPink - 4);
105
106 h_quality_other = new TH1F("quality_other", "Fit quality flag for unexpectedly saved waveforms", 4, 0, 4);
107 h_quality_other->GetXaxis()->SetTitle("Quality flag. 0-good,1-int overflow,2-low amplitude,3-bad chi2");
108 h_quality_other->SetFillColor(kPink - 4);
109
110 h_bad_quality = new TH1F("bad_quality", "Fraction of hits with bad chi2 (qual=3) and E > 1 GeV vs Cell ID",
112 h_bad_quality->GetXaxis()->SetTitle("Cell ID");
113 h_bad_quality->GetYaxis()->SetTitle("ECL hits count");
114
115 h_trigtag1 = new TH1F("trigtag1", "Consistency b/w global event number and trigger tag. 0-good, 1-DQM error", 2, 0, 2);
116 h_trigtag1->GetXaxis()->SetTitle("Flag value");
117 h_trigtag1->GetYaxis()->SetTitle("Events count");
118 h_trigtag1->SetDrawOption("hist");
119 h_trigtag1->SetFillColor(kPink - 4);
120
121 h_adc_hits = new TH1F("adc_hits", "Fraction of high-energy hits (E > 50 MeV)", 1001, 0, 1.001);
122 h_adc_hits->GetXaxis()->SetTitle("Fraction");
123 h_adc_hits->GetYaxis()->SetTitle("Events count");
124
125 h_time_crate_Thr1GeV_large = new TH1F("time_crate_Thr1GeV_large",
126 "Number of hits with timing outside #pm 100 ns per Crate ID (E > 1 GeV)",
127 52, 1, 53);
128 h_time_crate_Thr1GeV_large->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
129 h_time_crate_Thr1GeV_large->GetYaxis()->SetTitle("ECL hits count");
130
131 for (const auto& id : m_HitThresholds) {
132 std::string h_name, h_title;
133 h_name = str(boost::format("cid_Thr%1%MeV") % id);
134 h_title = str(boost::format("Occupancy per Cell ID (E > %1% MeV)") % id);
135 TH1F* h = new TH1F(h_name.c_str(), h_title.c_str(),
137 h->GetXaxis()->SetTitle("Cell ID");
138 h->GetYaxis()->SetTitle("Occupancy (hits / events_count)");
139 h_cids.push_back(h);
140 }
141
142 for (const auto& id : m_TotalEnergyThresholds) {
143 std::string h_name, h_title;
144 h_name = str(boost::format("edep_Thr%1%MeV") % id);
145 h_title = str(boost::format("Total energy (thr = %1% MeV)") % id);
146 TH1F* h = new TH1F(h_name.c_str(), h_title.c_str(), (int)(100 * m_EnergyUpperThr), 0, m_EnergyUpperThr);
147 h->GetXaxis()->SetTitle("Energy, [GeV]");
148 h_edeps.push_back(h);
149 }
150
151 for (const auto& id : m_TimingThresholds) {
152 std::string h_bar_name, h_bar_title;
153 std::string h_end_name, h_end_title;
154 h_bar_name = str(boost::format("time_barrel_Thr%1%MeV") % id);
155 h_bar_title = str(boost::format("Reconstructed time for ECL barrel (E > %1% MeV)") % id);
156 h_end_name = str(boost::format("time_endcaps_Thr%1%MeV") % id);
157 h_end_title = str(boost::format("Reconstructed time for ECL endcaps (E > %1% MeV)") % id);
158 TH1F* h_time_barrel = new TH1F(h_bar_name.c_str(), h_bar_title.c_str(), 206, -1030, 1030);
159 TH1F* h_time_endcap = new TH1F(h_end_name.c_str(), h_end_title.c_str(), 206, -1030, 1030);
160 h_time_barrel->GetXaxis()->SetTitle("Time, [ns]");
161 h_time_endcap->GetXaxis()->SetTitle("Time, [ns]");
162 h_time_barrels.push_back(h_time_barrel);
163 h_time_endcaps.push_back(h_time_endcap);
164 }
165
166 for (const auto& id : boost::combine(m_HitThresholds, m_HitNumberUpperLimits)) {
167 double id1 = 0, id2 = 0;
168 boost::tie(id1, id2) = id;
169 std::string h_name, h_title;
170 h_name = str(boost::format("ncev_Thr%1%MeV") % id1);
171 h_title = str(boost::format("Number of hits in event (E > %1% MeV)") % id1);
172 TH1F* h = new TH1F(h_name.c_str(), h_title.c_str(), id2, 0, id2);
173 h->GetXaxis()->SetTitle("Number of hits");
174 h_ncevs.push_back(h);
175 }
176
177 for (int i = 0; i < ECL_CRATES; i++) {
178 int crate = i + 1;
179 std::string h_name, h_title;
180 h_name = str(boost::format("time_crate_%1%_Thr1GeV") % (crate));
181 h_title = str(boost::format("Reconstructed time for ECL crate #%1% with E > 1 GeV") % (crate));
182 TH1F* h = new TH1F(h_name.c_str(), h_title.c_str(), 400, -100, 100);
183 h->GetXaxis()->SetTitle("Time [ns]");
184 h_time_crate_Thr1GeV.push_back(h);
185 }
186
187 for (const auto& id : m_WaveformOption) {
188 if (id != "all" && id != "psd" && id != "logic" && id != "rand" && id != "dphy" && id != "other")
189 B2WARNING("Waveform Options are not correctly assigned. They must be 'all', 'psd', 'logic', 'rand', 'dphy', 'other'!");
190 std::string h_title;
191 std::string h_cell_name;
192 if (id == "other") h_title = "Unexpectedly saved waveforms";
193 if (id == "psd") h_title = "#frac{Saved}{Expected} waveforms for high-energy hits (E > 50 MeV)";
194 if (id == "logic") h_title = "#frac{Saved}{Expected} waveforms for every 1000th event";
195 if (id == "rand") h_title = "#frac{Saved}{Expected} waveforms for random trigger events";
196 if (id == "dphy") h_title = "#frac{Saved}{Expected} waveforms for delayed bhabha (DPHY) events";
197 if (id == "all") h_title = "#frac{Saved}{Expected} waveforms for all events";
198 h_cell_name = str(boost::format("wf_cid_%1%") % (id));
199 TH1F* h_cell = new TH1F(h_cell_name.c_str(), h_title.c_str(),
201 h_cell->GetXaxis()->SetTitle("Cell ID");
202 if (id == "psd") {
203 h_cell_psd_norm = new TH1F("psd_cid", "Normalization to psd hits for cid",
205 }
206 if (id == "logic") {
207 h_evtot_logic = new TH1F("event_logic", "Event bank for logic", 1, 0, 1);
208 }
209 if (id == "rand") {
210 h_evtot_rand = new TH1F("event_rand", "Event bank for rand", 1, 0, 1);
211 }
212 if (id == "dphy") {
213 h_evtot_dphy = new TH1F("event_dphy", "Event bank for dphy", 1, 0, 1);
214 }
215 h_cells.push_back(h_cell);
216 }
217
218 //2D histograms creation.
219
220 h_trigtag2_trigid = new TH2F("trigtag2_trigid", "Internal data consistency vs crate. 0-good, 1-data corruption",
221 52, 1, 53, 11, -1, 10);
222 h_trigtag2_trigid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
223 h_trigtag2_trigid->GetYaxis()->SetTitle("Data consistency flag");
224
225 h_pedmean_cellid = new TProfile("pedmean_cellid", "Pedestal vs Cell ID",
227 h_pedmean_cellid->GetXaxis()->SetTitle("Cell ID");
228 h_pedmean_cellid->GetYaxis()->SetTitle("Ped. average (ADC units, #approx 0.05 MeV)");
229
230 h_pedrms_cellid = new TProfile("pedrms_cellid", "Pedestal stddev vs Cell ID",
232 h_pedrms_cellid->GetXaxis()->SetTitle("Cell ID");
233 h_pedrms_cellid->GetYaxis()->SetTitle("Ped. stddev (ADC units, #approx 0.05 MeV)");
234
235 h_pedrms_thetaid = new TProfile("pedrms_thetaid", "Pedestal stddev vs #theta ID",
236 68, 0, 68);
237 h_pedrms_thetaid->GetXaxis()->SetTitle("#theta ID (0-12=FWD, 59-67=BWD endcap)");
238 h_pedrms_thetaid->GetYaxis()->SetTitle("Ped. stddev (ADC units, #approx 0.05 MeV)");
239
240 h_trigtime_trigid = new TH2F("trigtime_trigid", "Trigger time vs Crate ID", 52, 1, 53, 145, 0, 145);
241 h_trigtime_trigid->GetXaxis()->SetTitle("Crate ID (same as ECLCollector ID)");
242 h_trigtime_trigid->GetYaxis()->SetTitle("Trigger time (only even, 0-142)");
243
244 h_pi0_mass = new TH1F("ecl_pi0_mass", "ecl_pi0_mass", 120, 0.08, 0.20);
245
246 //cd into parent directory.
247
248 oldDir->cd();
249}
250
252{
253 REG_HISTOGRAM; // required to register histograms to HistoManager.
254 m_ECLDigits.isRequired();
255 m_ECLCalDigits.isOptional();
256 m_ECLTrigs.isOptional();
257 m_ECLDsps.isOptional();
258 m_l1Trigger.isOptional();
259
260 if (!mapper.initFromDB()) B2FATAL("ECL DQM: Can't initialize eclChannelMapper");
261
262 ecltot.resize(m_TotalEnergyThresholds.size());
263 nhits.resize(m_HitNumberUpperLimits.size());
264
266
267 v_totalthrApsd.resize((m_calibrationThrApsd->getCalibVector()).size());
268 for (size_t i = 0; i < v_totalthrApsd.size(); i++) v_totalthrApsd[i] = (int)(m_calibrationThrApsd->getCalibVector())[i];
269}
270
272{
273 h_evtot->Reset();
274 h_evtot_logic->Reset();
275 h_evtot_rand->Reset();
276 h_evtot_dphy->Reset();
277 h_quality->Reset();
278 h_quality_other->Reset();
279 h_bad_quality->Reset();
280 h_trigtag1->Reset();
281 h_adc_hits->Reset();
283 h_cell_psd_norm->Reset();
284 for (TH1F* histogram : h_cids)
285 histogram->Reset();
286 for (TH1F* histogram : h_edeps)
287 histogram->Reset();
288 for (TH1F* histogram : h_time_barrels)
289 histogram->Reset();
290 for (TH1F* histogram : h_time_endcaps)
291 histogram->Reset();
292 for (TH1F* histogram : h_ncevs)
293 histogram->Reset();
294 for (TH1F* histogram : h_cells)
295 histogram->Reset();
296 for (TH1F* histogram : h_time_crate_Thr1GeV)
297 histogram->Reset();
298 h_trigtag2_trigid->Reset();
299 h_pedmean_cellid->Reset();
300 h_pedrms_cellid->Reset();
301 h_pedrms_thetaid->Reset();
302 h_trigtime_trigid->Reset();
303 h_pi0_mass->Reset();
304}
305
307{
308 int trigtag1 = 0;
309 int NDigits = 0;
310 for (auto& value : ecltot) value = 0;
311 for (auto& value : nhits) value = 0;
312 bool bhatrig = false;
313
314 if (m_l1Trigger.isValid() && m_DPHYTTYP) bhatrig = m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_DPHY;
315 else if (m_l1Trigger.isValid() && !m_DPHYTTYP) {
316 try { bhatrig = m_l1Trigger->testInput("bha_delay"); }
317 catch (const std::exception&) { bhatrig = false; }
318 }
319
320 m_iEvent = -1;
321 if (m_eventmetadata.isValid()) {
322 if (m_eventmetadata->getErrorFlag() != 0x10) {
323 m_iEvent = m_eventmetadata->getEvent();
324 h_evtot->Fill(0);
325 for (const auto& id : m_WaveformOption) {
326 if (id == "logic" && m_iEvent % 1000 == 999) h_evtot_logic->Fill(0);
327 if (id == "rand" && isRandomTrigger()) h_evtot_rand->Fill(0);
328 if (id == "dphy" && bhatrig) h_evtot_dphy->Fill(0);
329 }
330 }
331 }
332
333 for (auto& aECLDigit : m_ECLDigits) {
334 int i = aECLDigit.getCellId() - 1;
335 h_quality->Fill(aECLDigit.getQuality()); //Fit quality histogram filling.
336 if (aECLDigit.getAmp() > 2.e04 && aECLDigit.getQuality() == 3) h_bad_quality->Fill(aECLDigit.getCellId());
337 if (aECLDigit.getAmp() >= (v_totalthrApsd[i] / 4 * 4)) NDigits ++;
338 for (const auto& id : m_WaveformOption) {
339 if (id != "psd") continue;
340 else if (id == "psd" && (m_iEvent % 1000 == 999 ||
341 isRandomTrigger() ||
342 bhatrig ||
343 aECLDigit.getAmp() < (v_totalthrApsd[i] / 4 * 4))) continue;
344 h_cell_psd_norm->Fill(aECLDigit.getCellId());
345 }
346 }
347
348 for (auto& aECLTrig : m_ECLTrigs) {
349 double itrg = aECLTrig.getTimeTrig();
350 //trigger time conversion to acceptable units in range (0, ..., 142).
351 //one trigger time clock corresponds to 0.567/144*1000 = 3.93 ns
352 int tg = (int)itrg - 2 * ((int)itrg / 8);
353 h_trigtime_trigid->Fill(aECLTrig.getTrigId(), tg); //Trigger time histogram filling.
354 trigtag1 += aECLTrig.getTrigTag();
355 h_trigtag2_trigid->Fill(aECLTrig.getTrigId(), aECLTrig.getTrigTagQualityFlag()); //Data consistency histogram filling.
356 }
357
358 if (m_ECLTrigs.getEntries() > 0) {
359 int flagtag = 1;
360 trigtag1 /= m_ECLTrigs.getEntries();
361 int compar = (65535 & m_iEvent);
362 if (compar == trigtag1) flagtag = 0;
363 h_trigtag1->Fill(flagtag); //Trigger tag flag #1 histogram filling.
364 }
365
366 for (auto& aECLCalDigit : m_ECLCalDigits) {
367 int cid = aECLCalDigit.getCellId();
368 double energy = aECLCalDigit.getEnergy(); //get calibrated energy.
369 double timing = aECLCalDigit.getTime(); //get calibrated time.
370
371 for (size_t i = 0; i < m_HitThresholds.size(); i++) {
372 auto thrGeV = m_HitThresholds[i] / 1000.;
373 if (energy > thrGeV) {
374 h_cids[i]->Fill(cid);
375 nhits[i] += 1;
376 }
377 }
378
379 for (const auto& thr : m_TotalEnergyThresholds | indexed(0)) {
380 auto thrGeV = thr.value() / 1000.;
381 if (energy > thrGeV) ecltot[thr.index()] += energy;
382 }
383
384 for (const auto& thr : m_TimingThresholds | indexed(0)) {
385 auto thrGeV = thr.value() / 1000.;
386 if (energy > thrGeV) {
387 if (cid > ECL_FWD_CHANNELS && cid <= ECL_FWD_CHANNELS + ECL_BARREL_CHANNELS) h_time_barrels[thr.index()]->Fill(timing);
388 else h_time_endcaps[thr.index()]->Fill(timing);
389 }
390 }
391
392 if (energy > 1.000 && std::abs(timing) < 100.) h_time_crate_Thr1GeV[mapper.getCrateID(cid) - 1]->Fill(timing);
393 if (energy > 1.000 && std::abs(timing) > 100.) h_time_crate_Thr1GeV_large->Fill(mapper.getCrateID(cid));
394 }
395
396 for (const auto& h : h_edeps | indexed(0)) {
397 h.value()->Fill(ecltot[h.index()]);
398 }
399
400 for (const auto& h : h_ncevs | indexed(0)) {
401 h.value()->Fill(nhits[h.index()]);
402 }
403
404 for (auto& aECLDsp : m_ECLDsps) {
405 int i = aECLDsp.getCellId() - 1; //get number of Cell ID in m_DspArray.
406 aECLDsp.getDspA(m_DspArray[i]);
407 m_PedestalMean[i] = 0;
408 m_PedestalRms[i] = 0;
409
410 for (int j = 0; j < 16; j++) m_PedestalMean[i] += m_DspArray[i][j];
411 m_PedestalMean[i] /= 16;
412 h_pedmean_cellid->Fill(aECLDsp.getCellId(), m_PedestalMean[i]); //Pedestal Avg histogram filling.
413
414 for (int j = 0; j < 16; j++) m_PedestalRms[i] += pow(m_DspArray[i][j] - m_PedestalMean[i], 2);
415 m_PedestalRms[i] = sqrt(m_PedestalRms[i] / 15.);
416 h_pedrms_cellid->Fill(aECLDsp.getCellId(), m_PedestalRms[i]); //Pedestal stddev histogram filling.
417 m_geom->Mapping(i);
419
420 ECLDigit* aECLDigit = ECLDigit::getByCellID(aECLDsp.getCellId());
421
422 for (const auto& iter : m_WaveformOption | indexed(0)) {
423 const auto& index = iter.index();
424 const auto& wf_opt = iter.value();
425 if (wf_opt != "all" && wf_opt != "psd" && wf_opt != "logic" && wf_opt != "rand" && wf_opt != "dphy" && wf_opt != "other") continue;
426 else if (wf_opt == "psd" && (m_iEvent % 1000 == 999 || isRandomTrigger() || bhatrig ||
427 !aECLDigit || aECLDigit->getAmp() < (v_totalthrApsd[i] / 4 * 4))) continue;
428 else if (wf_opt == "logic" && m_iEvent % 1000 != 999) continue;
429 else if (wf_opt == "rand" && (m_iEvent % 1000 == 999 || !isRandomTrigger())) continue;
430 else if (wf_opt == "dphy" && (m_iEvent % 1000 == 999 || !bhatrig)) continue;
431 else if (wf_opt == "other" && (m_iEvent % 1000 == 999 || isRandomTrigger() || bhatrig ||
432 (aECLDigit && aECLDigit->getAmp() >= (v_totalthrApsd[i] / 4 * 4)))) continue;
433 h_cells[index]->Fill(aECLDsp.getCellId());
434 if (wf_opt == "other" && aECLDigit) h_quality_other->Fill(aECLDigit->getQuality());
435 }
436 }
437 if (m_ECLDigits.getEntries() > 0)
438 h_adc_hits->Fill((double)NDigits / (double)m_ECLDigits.getEntries()); //Fraction of high-energy hits
439
441}
442
444{
445}
446
447
449{
450}
451
453{
454 if (!m_l1Trigger.isValid()) return false;
455 return m_l1Trigger->getTimType() == TRGSummary::ETimingType::TTYP_RAND ||
457}
458
460{
461 const std::string trg_identifier = "software_trigger_cut&skim&accept_hadron";
463 if (!result.isValid()) return false;
464
465 // check if trigger identifier is available
466 const std::map<std::string, int>& results = result->getResults();
467 if (results.find(trg_identifier) == results.end()) return false;
468 // apply software trigger
469 const bool accepted = (result->getResult(trg_identifier) == SoftwareTriggerCutResult::c_accept);
470 if (!accepted) return false;
471
473 if (!particles.isValid()) return false;
474
475 for (unsigned int i = 0; i < particles->getListSize(); i++) {
476 const Particle* part = particles->getParticle(i);
477 auto inv_mass = Variable::eclClusterOnlyInvariantMass(part);
478 h_pi0_mass->Fill(inv_mass);
479 }
480
481 return true;
482}
483
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:459
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:251
virtual ~ECLDQMModule()
Destructor.
Definition: eclDQM.cc:83
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:306
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:443
std::vector< TH1F * > h_ncevs
Histogram vector: Channel multiplicity.
Definition: eclDQM.h:168
virtual void terminate() override
Terminate.
Definition: eclDQM.cc:448
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:271
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
bool isRandomTrigger()
Definition: eclDQM.cc:452
virtual void defineHisto() override
Function to define histograms.
Definition: eclDQM.cc:88
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: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
Class to store reconstructed particles.
Definition: Particle.h:75
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: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
@ c_accept
Accept this event.
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.