Belle II Software development
DQMHistAnalysisSVDDose.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#include <dqm/analysis/modules/DQMHistAnalysisSVDDose.h>
10#include <framework/utilities/Utils.h>
11#include <TROOT.h>
12#include <TText.h>
13
14using namespace std;
15using namespace Belle2;
16
17// Utility function
18inline double getClockSeconds() { return Utils::getClock() / 1e9; }
19
20REG_MODULE(DQMHistAnalysisSVDDose);
21
22DQMHistAnalysisSVDDoseModule::DQMHistAnalysisSVDDoseModule()
23{
24 setDescription("Monitoring of SVD Dose with events from Poisson trigger w/o inj. veto. See also SVDDQMDoseModule.");
25 // THIS MODULE CAN NOT BE RUN IN PARALLEL
26 addParam("pvPrefix", m_pvPrefix, "Prefix for EPICS PVs.", std::string("DQM:SVD:"));
27 addParam("epicsUpdateSeconds", m_epicsUpdateSeconds,
28 "Minimum interval between two successive PV updates (in seconds).", 1000.0);
29 addParam("pvSuffix", m_pvSuffix, "Suffix for EPICS PVs.", std::string(":Occ:Pois:Avg"));
30 addParam("deltaTPVSuffix", m_deltaTPVSuffix, "Suffix for the PV that monitors the update interval of the PVs.",
31 std::string("Occ:Pois:UpdateInterval"));
32 addParam("statePVSuffix", m_statePVSuffix, "Suffix for the PV with the state of the monitoring.",
33 std::string("Occ:Pois:State"));
34}
35
36DQMHistAnalysisSVDDoseModule::~DQMHistAnalysisSVDDoseModule()
37{
38}
39
41{
42 B2DEBUG(18, "DQMHistAnalysisSVDDose: initialize");
43
44 gROOT->cd(); // Don't know why I need this, but DQMHistAnalysisSVDOnMiraBelle uses it
45
46 m_monObj = getMonitoringObject("svd"); // To write to MiraBelle
47
48 m_c_instOccu.reserve(c_sensorGroups.size());
49 m_c_occuLER.reserve(c_sensorGroups.size());
50 m_h_occuLER.reserve(c_sensorGroups.size());
51 m_c_occuHER.reserve(c_sensorGroups.size());
52 m_h_occuHER.reserve(c_sensorGroups.size());
53 m_c_occuLER1.reserve(c_sensorGroups.size());
54 m_h_occuLER1.reserve(c_sensorGroups.size());
55 m_c_occuHER1.reserve(c_sensorGroups.size());
56 m_h_occuHER1.reserve(c_sensorGroups.size());
57 m_c_instOccuAll.reserve(c_sensorGroups.size());
58 m_c_occuLERAll.reserve(c_sensorGroups.size());
59 m_h_occuLERAll.reserve(c_sensorGroups.size());
60 m_c_occuHERAll.reserve(c_sensorGroups.size());
61 m_h_occuHERAll.reserve(c_sensorGroups.size());
62 m_c_occuLER1All.reserve(c_sensorGroups.size());
63 m_h_occuLER1All.reserve(c_sensorGroups.size());
64 m_c_occuHER1All.reserve(c_sensorGroups.size());
65 m_h_occuHER1All.reserve(c_sensorGroups.size());
66 for (const auto& group : c_sensorGroups) {
67 TCanvas* c = new TCanvas("SVDDOSE/c_svd_instOccupancy_" + group.nameSuffix + "_pois",
68 "Instantaneous occupancy (Pois. trig.) " + group.titleSuffix);
69 m_c_instOccu.push_back(c);
70 c = new TCanvas("SVDDOSE/c_svd_occuLER_" + group.nameSuffix + "_pois",
71 "Occupancy vs time since LER inj. (Pois. trig.) " + group.titleSuffix);
72 m_c_occuLER.push_back(c);
73 m_h_occuLER.push_back(nullptr);
74 c = new TCanvas("SVDDOSE/c_svd_occuHER_" + group.nameSuffix + "_pois",
75 "Occupancy vs time since HER inj. (Pois. trig.) " + group.titleSuffix);
76 m_c_occuHER.push_back(c);
77 m_h_occuHER.push_back(nullptr);
78 c = new TCanvas("SVDDOSE/c_svd_1DoccuLER_" + group.nameSuffix + "_pois",
79 "Occupancy vs time since LER inj. (Pois. trig.) " + group.titleSuffix);
80 m_c_occuLER1.push_back(c);
81 m_h_occuLER1.push_back(nullptr);
82 c = new TCanvas("SVDDOSE/c_svd_1DoccuHER_" + group.nameSuffix + "_pois",
83 "Occupancy vs time since HER inj. (Pois. trig.) " + group.titleSuffix);
84 m_c_occuHER1.push_back(c);
85 m_h_occuHER1.push_back(nullptr);
86
87 c = new TCanvas("SVDDOSE/c_svd_instOccupancy_" + group.nameSuffix + "_all",
88 "Instantaneous occupancy (all events) " + group.titleSuffix);
89 m_c_instOccuAll.push_back(c);
90 c = new TCanvas("SVDDOSE/c_svd_occuLER_" + group.nameSuffix + "_all",
91 "Occupancy vs time since LER inj. (all events) " + group.titleSuffix);
92 m_c_occuLERAll.push_back(c);
93 m_h_occuLERAll.push_back(nullptr);
94 c = new TCanvas("SVDDOSE/c_svd_occuHER_" + group.nameSuffix + "_all",
95 "Occupancy vs time since HER inj. (all events) " + group.titleSuffix);
96 m_c_occuHERAll.push_back(c);
97 m_h_occuHERAll.push_back(nullptr);
98 c = new TCanvas("SVDDOSE/c_svd_1DoccuLER_" + group.nameSuffix + "_all",
99 "Occupancy vs time since LER inj. (all events) " + group.titleSuffix);
100 m_c_occuLER1All.push_back(c);
101 m_h_occuLER1All.push_back(nullptr);
102 c = new TCanvas("SVDDOSE/c_svd_1DoccuHER_" + group.nameSuffix + "_all",
103 "Occupancy vs time since HER inj. (all events) " + group.titleSuffix);
104 m_c_occuHER1All.push_back(c);
105 m_h_occuHER1All.push_back(nullptr);
106 }
107
108 m_legend = new TPaveText(0.53, 0.73, 0.68, 0.88, "brNDC");
109 m_legend->AddText("LER inj."); ((TText*)m_legend->GetListOfLines()->Last())->SetTextColor(kRed);
110 m_legend->AddText("HER inj."); ((TText*)m_legend->GetListOfLines()->Last())->SetTextColor(kAzure);
111 m_legend->AddText("No inj."); ((TText*)m_legend->GetListOfLines()->Last())->SetTextColor(kBlack);
112
113#ifdef _BELLE2_EPICS
114 if (getUseEpics()) {
115 // Channels for the occupancies
116 m_myPVs.resize(c_sensorGroups.size());
117 for (unsigned int g = 0; g < c_sensorGroups.size(); g++)
118 SEVCHK(ca_create_channel((m_pvPrefix + c_sensorGroups[g].pvMiddle + m_pvSuffix).data(),
119 NULL, NULL, 10, &m_myPVs[g].mychid), "ca_create_channel");
120 // Channels for update interval and state
121 SEVCHK(ca_create_channel((m_pvPrefix + m_deltaTPVSuffix).data(),
122 NULL, NULL, 10, &m_timeSinceLastPVUpdateChan), "ca_create_channel");
123 SEVCHK(ca_create_channel((m_pvPrefix + m_statePVSuffix).data(),
124 NULL, NULL, 10, &m_stateChan), "ca_create_channel");
125 // Actually do create the channels, communicating with the IOC
126 SEVCHK(ca_pend_io(2.0), "ca_pend_io");
127 // Get the state enum
128 if (m_stateChan) {
129 SEVCHK(ca_get(DBR_CTRL_ENUM, m_stateChan, &m_stateCtrl), "ca_get");
130 B2DEBUG(19, "State PV initialized (ca_get)" << LogVar("value", m_stateCtrl.value));
131 SEVCHK(ca_pend_io(2.0), "ca_pend_io");
132 } else {
133 B2DEBUG(18, "State PV failed to initialize, will retry in beginRun(), event() and endRun().");
134 }
135 // First update should happen m_epicsUpdateSeconds from now
136 m_lastPVUpdate = getClockSeconds();
137 }
138#endif
139}
140
142{
143 // Set status PV to running, reset last update time
144#ifdef _BELLE2_EPICS
145 if (getUseEpics()) {
146 B2DEBUG(19, "beginRun: setting state PV to RUNNING");
147 m_stateCtrl.value = 1;
148 if (m_stateChan) {
149 SEVCHK(ca_put(DBR_ENUM, m_stateChan, &m_stateCtrl.value), "ca_put");
150 } else {
151 SEVCHK(ca_create_channel((m_pvPrefix + m_statePVSuffix).data(),
152 NULL, NULL, 10, &m_stateChan), "ca_create_channel (reconnection)");
153 }
154 SEVCHK(ca_pend_io(2.0), "ca_pend_io");
155 // First update should happen m_epicsUpdateSeconds from now
156 m_lastPVUpdate = getClockSeconds();
157 }
158#endif
159}
160
162{
163 // Update PVs ("write" to EPICS)
164#ifdef _BELLE2_EPICS
165 double timeSinceLastPVUpdate = getClockSeconds() - m_lastPVUpdate;
166 if (getUseEpics() && timeSinceLastPVUpdate >= m_epicsUpdateSeconds) {
167 // Dummy ca_get to ensure reconnection to the IOC in case of past errors
168 if (m_stateChan) {
169 SEVCHK(ca_get(DBR_CTRL_ENUM, m_stateChan, &m_stateCtrl), "ca_get");
170 } else {
171 SEVCHK(ca_create_channel((m_pvPrefix + m_statePVSuffix).data(),
172 NULL, NULL, 10, &m_stateChan), "ca_create_channel (reconnection)");
173 }
174 SEVCHK(ca_pend_io(2.0), "ca_pend_io");
175
176 // Write updateInterval PV and state PV first
177 if (m_timeSinceLastPVUpdateChan) {
178 SEVCHK(ca_put(DBR_DOUBLE, m_timeSinceLastPVUpdateChan, (void*)&timeSinceLastPVUpdate), "ca_put");
179 } else {
180 SEVCHK(ca_create_channel((m_pvPrefix + m_deltaTPVSuffix).data(), NULL, NULL, 10, &m_timeSinceLastPVUpdateChan),
181 "ca_create_channel (reconnection)");
182 }
183 m_stateCtrl.value = 1; // If IOC is restarted this PV must be re-updated
184 if (m_stateChan) {
185 SEVCHK(ca_put(DBR_ENUM, m_stateChan, &m_stateCtrl.value), "ca_put");
186 } else {
187 SEVCHK(ca_create_channel((m_pvPrefix + m_statePVSuffix).data(),
188 NULL, NULL, 10, &m_stateChan), "ca_create_channel (reconnection)");
189 }
190 // Update occupancy PVs
191 for (unsigned int g = 0; g < c_sensorGroups.size(); g++) {
192 const auto& group = c_sensorGroups[g];
193 double nHits = 0.0, nEvts = 0.0;
194 for (TString dir : {"SVDDoseLERInjPois", "SVDDoseHERInjPois", "SVDDoseNoInjPois"}) {
195 auto hHits = findHistT<TH2F>(dir + "/SVDHitsVsTime_" + group.nameSuffix);
196 auto hEvts = findHistT<TH2F>(dir + "/SVDEvtsVsTime");
197 if (!hHits || !hEvts) {
198 B2WARNING("Histograms needed for Average Poisson Occupancy U-side not found.");
199 nEvts = 0.0;
200 break;
201 }
202 nHits += hHits->GetEntries();
203 nEvts += hEvts->GetEntries();
204 }
205
206 B2DEBUG(19, "DQMHistAnalysisSVDDose: PV write"
207 << LogVar("group", group.nameSuffix.Data())
208 << LogVar("nEvts", nEvts) << LogVar("nHits", nHits));
209
210 auto& pv = m_myPVs[g];
211 double delta_nHits = nHits - pv.lastNHits;
212 double delta_nEvts = nEvts - pv.lastNEvts;
213 double occ = delta_nEvts > 0.0 ? (delta_nHits / delta_nEvts * 100.0 / group.nStrips) : -1.0;
214 if (pv.mychid) {
215 SEVCHK(ca_put(DBR_DOUBLE, pv.mychid, (void*)&occ), "ca_put");
216 } else {
217 SEVCHK(ca_create_channel((m_pvPrefix + c_sensorGroups[g].pvMiddle + m_pvSuffix).data(),
218 NULL, NULL, 10, &m_myPVs[g].mychid), "ca_create_channel (reconnection)");
219 }
220 pv.lastNEvts = nEvts;
221 pv.lastNHits = nHits;
222 }
223
224 // Actually write all the PVs
225 SEVCHK(ca_pend_io(2.0), "ca_pend_io");
226 m_lastPVUpdate = getClockSeconds();
227 }
228#endif
229
231}
232
234{
235 B2DEBUG(18, "DQMHistAnalysisSVDDose: endRun");
236
237 // EPICS: reset the counters used for the delta computation, set state to NOT RUNNING
238#ifdef _BELLE2_EPICS
239 if (getUseEpics()) {
240 B2DEBUG(19, "endRun: setting state PV to NOT RUNNING");
241 m_stateCtrl.value = 0;
242 if (m_stateChan) {
243 SEVCHK(ca_put(DBR_ENUM, m_stateChan, &m_stateCtrl.value), "ca_put");
244 } else {
245 SEVCHK(ca_create_channel((m_pvPrefix + m_statePVSuffix).data(),
246 NULL, NULL, 10, &m_stateChan), "ca_create_channel (reconnection)");
247 }
248 SEVCHK(ca_pend_io(2.0), "ca_pend_io");
249 // Reset events and hits counters
250 for (auto& pv : m_myPVs)
251 pv.lastNEvts = pv.lastNHits = 0.0;
252 }
253#endif
254
255 // Write to MiraBelle
256 for (unsigned int g = 0; g < c_sensorGroups.size(); g++) {
257 const auto& group = c_sensorGroups[g];
258 double nHits = 0.0, nEvts = 0.0;
259 for (TString dir : {"SVDDoseLERInjPois", "SVDDoseHERInjPois", "SVDDoseNoInjPois"}) {
260 auto hHits = findHistT<TH2F>(dir + "/SVDHitsVsTime_" + group.nameSuffix);
261 auto hEvts = findHistT<TH2F>(dir + "/SVDEvtsVsTime");
262 if (!hHits || !hEvts) {
263 B2WARNING("Histograms needed for Average Poisson Occupancy U-side not found.");
264 nEvts = 0.0;
265 break;
266 }
267 nHits += hHits->GetEntries();
268 nEvts += hEvts->GetEntries();
269 }
270
271 B2DEBUG(19, "DQMHistAnalysisSVDDose: MonObj write"
272 << LogVar("group", group.nameSuffix.Data())
273 << LogVar("nEvts", nEvts) << LogVar("nHits", nHits));
274
275 if (nEvts > 0) {
276 const double occ = nHits / nEvts * 100.0 / group.nStrips;
277 TString vName = group.nameSuffix + "OccPoisAvg"; // e.g. L3XXUOccPoisAvg
278 m_monObj->setVariable(vName.Data(), occ);
279 }
280 }
281
283}
284
286{
287 B2DEBUG(18, "DQMHistAnalysisSVDDose: updating canvases");
288
289 for (unsigned int g = 0; g < c_sensorGroups.size(); g++) {
290 const auto& group = c_sensorGroups[g];
291
292 auto c = m_c_instOccu[g];
293 auto hLER = findHistT<TH1F>("SVDDoseLERInjPois/SVDInstOccu_" + group.nameSuffix);
294 auto hHER = findHistT<TH1F>("SVDDoseHERInjPois/SVDInstOccu_" + group.nameSuffix);
295 auto hNo = findHistT<TH1F>("SVDDoseNoInjPois/SVDInstOccu_" + group.nameSuffix);
296 if (hLER && hHER && hNo) {
297 hLER->SetLineColor(kRed);
298 hHER->SetLineColor(kAzure);
299 hNo->SetLineColor(kBlack);
300 carryOverflowOver(hLER);
301 carryOverflowOver(hHER);
303 c->Clear();
304 c->cd(0);
305 hNo->SetTitle("SVD instantaneous occu. " + group.titleSuffix + " U-side Pois. trig.");
306 hNo->Draw("hist"); // hNo usually has the larger maximum by far
307 hLER->Draw("hist same");
308 hHER->Draw("hist same");
309 c->SetLogy();
310 m_legend->Draw();
311 }
312
313 c = m_c_occuLER[g];
314 auto hHits = findHistT<TH2F>("SVDDoseLERInjPois/SVDHitsVsTime_" + group.nameSuffix);
315 auto hEvts = findHistT<TH2F>("SVDDoseLERInjPois/SVDEvtsVsTime");
316 if (hHits && hEvts) {
317 m_h_occuLER[g] = divide(hHits, hEvts, 100.0f / group.nStrips, m_h_occuLER[g]);
318 m_h_occuLER[g]->SetTitle("SVD Occupancy " + group.titleSuffix + " - LER inj. Pois. trig."
319 ";Time since last injection [#mus];Time in beam cycle [#mus]"
320 ";Occupancy [%]");
321 m_h_occuLER[g]->SetMinimum(1e-3);
322 m_h_occuLER[g]->SetMaximum(10);
323 c->Clear();
324 c->cd(0);
325 c->SetRightMargin(0.16); // For the colorbar
326 m_h_occuLER[g]->Draw("COLZ");
327 c->SetLogz();
328 }
329
330 c = m_c_occuLER1[g];
331 auto hpEvts = findHistT<TH1F>("SVDDoseLERInjPois/SVDEvtsVsTime1");
332 auto hpHits = findHistT<TH1F>("SVDDoseLERInjPois/SVDHitsVsTime1_" + group.nameSuffix);
333 if (hpHits && hpEvts) {
334 m_h_occuLER1[g] = divide(hpHits, hpEvts, 100.0f / group.nStrips, m_h_occuLER1[g]);
335 m_h_occuLER1[g]->SetTitle("SVD Occupancy " + group.titleSuffix + " - LER inj. Pois. trig."
336 ";Time since last injection [#mus];Occupancy [%]");
337 m_h_occuLER1[g]->SetMinimum(1e-3);
338 m_h_occuLER1[g]->SetMaximum(10);
339 c->Clear();
340 c->cd(0);
341 m_h_occuLER1[g]->SetMarkerStyle(7);
342 m_h_occuLER1[g]->Draw("hist P");
343 c->SetLogy();
344 }
345
346 c = m_c_occuHER[g];
347 hHits = findHistT<TH2F>("SVDDoseHERInjPois/SVDHitsVsTime_" + group.nameSuffix);
348 hEvts = findHistT<TH2F>("SVDDoseHERInjPois/SVDEvtsVsTime");
349 if (hHits && hEvts) {
350 m_h_occuHER[g] = divide(hHits, hEvts, 100.0f / group.nStrips, m_h_occuHER[g]);
351 m_h_occuHER[g]->SetTitle("SVD Occupancy " + group.titleSuffix + " - HER inj. Pois. trig."
352 ";Time since last injection [#mus];Time in beam cycle [#mus]"
353 ";Occupancy [%]");
354 m_h_occuHER[g]->SetMinimum(1e-3);
355 m_h_occuHER[g]->SetMaximum(10);
356 c->Clear();
357 c->cd(0);
358 c->SetRightMargin(0.16); // For the colorbar
359 m_h_occuHER[g]->Draw("COLZ");
360 c->SetLogz();
361 }
362
363 c = m_c_occuHER1[g];
364 hpEvts = findHistT<TH1F>("SVDDoseHERInjPois/SVDEvtsVsTime1");
365 hpHits = findHistT<TH1F>("SVDDoseHERInjPois/SVDHitsVsTime1_" + group.nameSuffix);
366 if (hpHits && hpEvts) {
367 m_h_occuHER1[g] = divide(hpHits, hpEvts, 100.0f / group.nStrips, m_h_occuHER1[g]);
368 m_h_occuHER1[g]->SetTitle("SVD Occupancy " + group.titleSuffix + " - HER inj. Pois. trig."
369 ";Time since last injection [#mus];Occupancy [%]");
370 m_h_occuHER1[g]->SetMinimum(1e-3);
371 m_h_occuHER1[g]->SetMaximum(10);
372 c->Clear();
373 c->cd(0);
374 m_h_occuHER1[g]->SetMarkerStyle(7);
375 m_h_occuHER1[g]->Draw("hist P");
376 c->SetLogy();
377 }
378
379 // ========== All events =============
380 c = m_c_instOccuAll[g];
381 hLER = findHistT<TH1F>("SVDDoseLERInjAll/SVDInstOccu_" + group.nameSuffix);
382 hHER = findHistT<TH1F>("SVDDoseHERInjAll/SVDInstOccu_" + group.nameSuffix);
383 hNo = findHistT<TH1F>("SVDDoseNoInjAll/SVDInstOccu_" + group.nameSuffix);
384 if (hLER && hHER && hNo) {
385 hLER->SetLineColor(kRed);
386 hHER->SetLineColor(kAzure);
387 hNo->SetLineColor(kBlack);
388 carryOverflowOver(hLER);
389 carryOverflowOver(hHER);
391 c->Clear();
392 c->cd(0);
393 hNo->SetTitle("SVD instantaneous occu. " + group.titleSuffix + " U-side all events");
394 hNo->Draw("hist"); // hNo usually has the larger maximum by far
395 hLER->Draw("hist same");
396 hHER->Draw("hist same");
397 c->SetLogy();
398 m_legend->Draw();
399 }
400
401 c = m_c_occuLERAll[g];
402 hHits = findHistT<TH2F>("SVDDoseLERInjAll/SVDHitsVsTime_" + group.nameSuffix);
403 hEvts = findHistT<TH2F>("SVDDoseLERInjAll/SVDEvtsVsTime");
404 if (hHits && hEvts) {
405 m_h_occuLERAll[g] = divide(hHits, hEvts, 100.0f / group.nStrips, m_h_occuLERAll[g]);
406 m_h_occuLERAll[g]->SetTitle("SVD Occupancy " + group.titleSuffix + " - LER inj. all events"
407 ";Time since last injection [#mus];Time in beam cycle [#mus]"
408 ";Occupancy [%]");
409 m_h_occuLERAll[g]->SetMinimum(1e-3);
410 m_h_occuLERAll[g]->SetMaximum(10);
411 c->Clear();
412 c->cd(0);
413 c->SetRightMargin(0.16); // For the colorbar
414 m_h_occuLERAll[g]->Draw("COLZ");
415 c->SetLogz();
416 }
417
418 c = m_c_occuLER1All[g];
419 hpEvts = findHistT<TH1F>("SVDDoseLERInjAll/SVDEvtsVsTime1");
420 hpHits = findHistT<TH1F>("SVDDoseLERInjAll/SVDHitsVsTime1_" + group.nameSuffix);
421 if (hpHits && hpEvts) {
422 m_h_occuLER1All[g] = divide(hpHits, hpEvts, 100.0f / group.nStrips, m_h_occuLER1All[g]);
423 m_h_occuLER1All[g]->SetTitle("SVD Occupancy " + group.titleSuffix + " - LER inj. all events"
424 ";Time since last injection [#mus];Occupancy [%]");
425 m_h_occuLER1All[g]->SetMinimum(1e-3);
426 m_h_occuLER1All[g]->SetMaximum(10);
427 c->Clear();
428 c->cd(0);
429 m_h_occuLER1All[g]->SetMarkerStyle(7);
430 m_h_occuLER1All[g]->Draw("hist P");
431 c->SetLogy();
432 }
433
434 c = m_c_occuHERAll[g];
435 hHits = findHistT<TH2F>("SVDDoseHERInjAll/SVDHitsVsTime_" + group.nameSuffix);
436 hEvts = findHistT<TH2F>("SVDDoseHERInjAll/SVDEvtsVsTime");
437 if (hHits && hEvts) {
438 m_h_occuHERAll[g] = divide(hHits, hEvts, 100.0f / group.nStrips, m_h_occuHERAll[g]);
439 m_h_occuHERAll[g]->SetTitle("SVD Occupancy " + group.titleSuffix + " - HER inj. all events"
440 ";Time since last injection [#mus];Time in beam cycle [#mus]"
441 ";Occupancy [%]");
442 m_h_occuHERAll[g]->SetMinimum(1e-3);
443 m_h_occuHERAll[g]->SetMaximum(10);
444 c->Clear();
445 c->cd(0);
446 c->SetRightMargin(0.16); // For the colorbar
447 m_h_occuHERAll[g]->Draw("COLZ");
448 c->SetLogz();
449 }
450
451 c = m_c_occuHER1All[g];
452 hpEvts = findHistT<TH1F>("SVDDoseHERInjAll/SVDEvtsVsTime1");
453 hpHits = findHistT<TH1F>("SVDDoseHERInjAll/SVDHitsVsTime1_" + group.nameSuffix);
454 if (hpHits && hpEvts) {
455 m_h_occuHER1All[g] = divide(hpHits, hpEvts, 100.0f / group.nStrips, m_h_occuHER1All[g]);
456 m_h_occuHER1All[g]->SetTitle("SVD Occupancy " + group.titleSuffix + " - HER inj. all events"
457 ";Time since last injection [#mus];Occupancy [%]");
458 m_h_occuHER1All[g]->SetMinimum(1e-3);
459 m_h_occuHER1All[g]->SetMaximum(10);
460 c->Clear();
461 c->cd(0);
462 m_h_occuHER1All[g]->SetMarkerStyle(7);
463 m_h_occuHER1All[g]->Draw("hist P");
464 c->SetLogy();
465 }
466 }
467}
468
470{
471 int i = h->GetNbinsX();
472 float t = h->GetBinContent(i) + h->GetBinContent(i + 1);
473 h->SetBinContent(i, t);
474 h->SetBinContent(i + 1, 0);
475}
476
477const vector<DQMHistAnalysisSVDDoseModule::SensorGroup> DQMHistAnalysisSVDDoseModule::c_sensorGroups = {
478 {"L31XU", "L3.1", "L3:1", 768 * 2},
479 {"L32XU", "L3.2", "L3:2", 768 * 2},
480 {"L3XXU", "L3 avg.", "L3", 768 * 14},
481 {"L4XXU", "L4 avg.", "L4", 768 * 30},
482 {"L5XXU", "L5 avg.", "L5", 768 * 48},
483 {"L6XXU", "L6 avg.", "L6", 768 * 80}
484};
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
bool getUseEpics(void)
Getter for EPICS usage.
void initialize() override final
Initializer.
std::vector< TCanvas * > m_c_occuLER
Canvases for the occu.
std::vector< TH2F * > m_h_occuHERAll
Histograms for the occu.
std::vector< TH1F * > m_h_occuLER1
Histograms for the 1D occu.
std::string m_deltaTPVSuffix
Suffix of the update-time monitoring PV.
std::vector< TCanvas * > m_c_occuHERAll
Canvases for the occu.
std::vector< TCanvas * > m_c_instOccu
Canvases for the instantaneous occupancy.
std::vector< TH2F * > m_h_occuLER
Histograms for the occ.
double m_epicsUpdateSeconds
Minimum interval between successive PV updates.
std::vector< TH2F * > m_h_occuLERAll
Histograms for the occu.
std::vector< TCanvas * > m_c_occuLER1
Canvases for the 1D occu.
std::string m_pvPrefix
Prefix for EPICS PVs.
std::vector< TH1F * > m_h_occuHER1All
Histograms for the 1D occu.
MonitoringObject * m_monObj
Monitoring object for MiraBelle.
std::vector< TCanvas * > m_c_occuLER1All
Canvases for the 1D occu.
void event() override final
This method is called for each event.
std::string m_statePVSuffix
Suffix of the state PV.
std::vector< TCanvas * > m_c_occuLERAll
Canvases for the occu.
std::vector< TCanvas * > m_c_occuHER
Canvases for the occu.
std::vector< TCanvas * > m_c_occuHER1All
Canvases for the 1D occu.
static void carryOverflowOver(TH1F *h)
Carries the content of the overflow bin into the last bin.
static const std::vector< SensorGroup > c_sensorGroups
List of sensors groups.
std::vector< TH1F * > m_h_occuLER1All
Histograms for the 1D occu.
void endRun() override final
This method is called if the current run ends.
void beginRun() override final
Called when entering a new run.
std::vector< TH1F * > m_h_occuHER1
Histograms for the 1D occu.
std::vector< TCanvas * > m_c_instOccuAll
Canvases for the instantaneous occupancy.
static T * divide(T *num, T *den, float scale=1.0f, T *res=nullptr)
Divide two histograms, ignoring errors on the second histogram.
TPaveText * m_legend
Legend of the inst.
std::vector< TCanvas * > m_c_occuHER1
Canvases for the 1D occu.
std::vector< TH2F * > m_h_occuHER
Histograms for the occu.
std::string m_pvSuffix
Suffix for EPICS PVs.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setVariable(const std::string &var, float val, float upErr=-1., float dwErr=-1)
set value to float variable (new variable is made if not yet existing)
Class to store variables with their name which were sent to the logging service.
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 getClock()
Return current value of the real-time clock.
Definition: Utils.cc:66
Abstract base class for different kinds of events.
STL namespace.