Belle II Software development
DQMHistAnalysisARICH.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 <dqm/analysis/modules/DQMHistAnalysisARICH.h>
11
12#include <TH1F.h>
13#include <TH2F.h>
14#include <TCanvas.h>
15#include <TLine.h>
16#include <TClass.h>
17#include <TROOT.h>
18
19#include <fstream>
20#include <vector>
21#include <algorithm>
22
23using namespace std;
24using namespace Belle2;
25
26//-----------------------------------------------------------------
27// Register module
28//-----------------------------------------------------------------
29
30REG_MODULE(DQMHistAnalysisARICH);
31
34{
35 // set module description (e.g. insert text)
36 setDescription("Modify and analyze the data quality histograms of ARICH");
38 addParam("debug", m_debug, "debug mode", false);
39 addParam("alert", m_enableAlert, "Enable color alert", true);
40 addParam("hotLimit", m_hotLimit, "Occupancy limit for hot channels", 0.005);
41 addParam("badApdOccLimit", m_badApdOccLimit, "Occupancy limit for bad APDs (in fraction of average occupancy)", 0.3);
42
43}
44
46{
47}
48
50{
51 gROOT->cd();
52
53 //definition of new TObjects for modification and analysis
54 for (int i = 0; i < 5; i++) {
55 m_LineForMB[i] = new TLine();
56 m_LineForMB[i]->SetLineStyle(3);
57 m_LineForMB[i]->SetLineWidth(1);
58 m_LineForMB[i]->SetLineColor(kBlack);
59 }
60
61 m_apdHist = new ARICHChannelHist("tmpChHist", "tmpChHist", 2);
62 m_apdPoly = new TH2Poly();
63 m_apdPoly->SetName("ARICH/apdHitMap");
64 m_apdPoly->SetTitle("# of hits/APD/event");
65 m_apdPoly->SetOption("colz");
66 m_c_apdHist = new TCanvas("ARICH/c_apdHist");
67
68 addDeltaPar("ARICH", "bitsPerChannel", HistDelta::c_Events,
69 1000000, 1); // update each 100k events (from daq histogram)
70
71 // addDeltaPar(m_histogramDirectoryName, m_histogramName, HistDelta::c_Entries, 10000, 1); // update each 10000 entries
72
73 registerEpicsPV("ARI:badAPDs", "badAPDs");
74 registerEpicsPV("ARI:deadHAPDs", "deadHAPDs");
75 registerEpicsPV("ARI:hotChannels", "hotChannels");
76
77 registerEpicsPV("ARI:hotLimit", "hotLimit");
78 registerEpicsPV("ARI:badApdOccLimit", "badApdOccLimit");
79 registerEpicsPV("ARI:sigBitFracLimit", "sigBitFracLimit");
80
81 B2DEBUG(20, "DQMHistAnalysisARICH: initialized.");
82}
83
85{
86
87 // read the alarm limits
88 double unused = NAN;
91 requestLimitsFromEpicsPVs("sigBitFracLimit", m_sigBitFracLowAlarm, m_sigBitFracLowWarn, unused, unused);
92
93 double hot = getEpicsPV("hotLimit");
94 if (hot != NAN) m_hotLimit = hot;
95
96 double bad = getEpicsPV("badApdOccLimit");
97 if (bad != NAN) m_badApdOccLimit = bad;
98
99}
100
102{
103
104 //Show alert by empty bins = red and strange entries = yellow
105 //Draw lines on mergerHits histogram for shifters to divide sectors
106 TH1* m_h_mergerHit = findHist("ARICH/mergerHit");
107 m_c_mergerHit = findCanvas("ARICH/c_mergerHit");
108 if (m_h_mergerHit != NULL && m_c_mergerHit != NULL) {
109 m_c_mergerHit->Clear();
110 m_c_mergerHit->cd();
111 m_h_mergerHit->SetMinimum(0);
112 m_h_mergerHit->Draw("hist");
113 gPad->Update();
114
115 int alertMerger = 0;
116 double mean = m_h_mergerHit->Integral() / 72;
117 for (int i = 0; i < 72; i++) {
118 int hit = m_h_mergerHit->GetBinContent(i + 1);
119 if ((bool)hit ^ (bool)m_h_mergerHit->GetEntries()) {
120 //only if the empty bin is not a masked merger, show alert.
121 auto itr = std::find(maskedMergers.begin(), maskedMergers.end(), i + 1);
122 if (itr == maskedMergers.end()) {
123 alertMerger = 2;
124 break;
125 }
126 }
127 if (hit > mean * 100 && alertMerger < 1) alertMerger = 1;
128 }
129 if (m_enableAlert && m_minStats < m_h_mergerHit->GetEntries()) m_c_mergerHit->SetFillColor(alertColor[alertMerger]);
130
131 //Draw lines divide the sectors
132 for (int i = 0; i < 5; i++) {
133 m_LineForMB[i]->DrawLine(12 * (i + 1) + 0.5, 0, 12 * (i + 1) + 0.5, gPad->GetUymax());
134 }
135
136 m_c_mergerHit->Modified();
137 } else {
138 B2INFO("Histogram/canvas named mergerHit is not found.");
139 }
140
141
142 //Show alert by the ratio of center 2 bins to side 2bins. <1.5 = red, <2 = yellow
143 TH1* m_h_bits = findHist("ARICH/bits");
144 m_c_bits = findCanvas("ARICH/c_bits");
145 if (m_h_bits != NULL && m_c_bits != NULL) {
146 m_c_bits->Clear();
147 m_c_bits->cd();
148 m_h_bits->SetMinimum(0);
149 m_h_bits->Draw("hist");
150 gPad->Update();
151
152 double side = m_h_bits->GetBinContent(2) + m_h_bits->GetBinContent(5);
153 double center = m_h_bits->GetBinContent(3) + m_h_bits->GetBinContent(4);
154 EStatus status = c_StatusGood;
155 if (center / side < m_sigBitFracLowWarn) status = c_StatusWarning;
156 if (center / side < m_sigBitFracLowAlarm) status = c_StatusError;
157 if (m_enableAlert && m_minStats < m_h_bits->GetEntries()) colorizeCanvas(m_c_bits, status);
158
159 m_c_bits->Modified();
160 } else {
161 B2INFO("Histogram/canvas named bits is not found.");
162 }
163
164 //Show alert by no entry = red and 0 peak = yellow
165 TH1* m_h_hitsPerEvent = findHist("ARICH/hitsPerEvent");
166 m_c_hitsPerEvent = findCanvas("ARICH/c_hitsPerEvent");
167 int nEvents = 0;
168 if (m_h_hitsPerEvent != NULL && m_c_hitsPerEvent != NULL) {
169 m_c_hitsPerEvent->Clear();
170 m_c_hitsPerEvent->cd();
171 m_h_hitsPerEvent->SetMinimum(0);
172 m_h_hitsPerEvent->Draw("hist");
173 gPad->Update();
174 nEvents = m_h_hitsPerEvent->GetEntries();
175 int alertHitsPerEvent = 0;
176 double mean = m_h_hitsPerEvent->GetMean();
177 if (mean < 1) alertHitsPerEvent = 1;
178 double entry = m_h_hitsPerEvent->GetEntries();
179 if (entry == 0) alertHitsPerEvent = 2;
180 if (m_enableAlert) m_c_hitsPerEvent->SetFillColor(alertColor[alertHitsPerEvent]);
181
182 m_c_hitsPerEvent->Modified();
183 } else {
184 B2INFO("Histogram/canvas named hitsPerEvent is not found.");
185 }
186
187 //Draw 2D hit map of channels and APDs
188 TH1* m_h_chHit = findHist("ARICH/chipHit");
189 if (m_h_chHit != NULL) {
190 int nevt = 0;
191 TH1* htmp = findHist("ARICH/hitsPerEvent");
192 if (htmp) nevt = htmp->GetEntries();
193 m_apdHist->fillFromTH1(m_h_chHit);
194 if (nevt) m_apdHist->Scale(1. / float(nevt));
195 m_apdPoly->SetMaximum(0.1);
197 m_apdPoly->SetMinimum(0.0001);
198 m_c_apdHist->Clear();
199 m_c_apdHist->cd();
200 m_apdPoly->Draw("colz");
201 m_apdPoly->GetXaxis()->SetTickLength(0);
202 m_apdPoly->GetYaxis()->SetTickLength(0);
203 gPad->SetLogz();
204 m_c_apdHist->Update();
205 } else {
206 B2INFO("Histogram named chipHit is not found.");
207 }
208
209 TH1F* chDigit = (TH1F*)findHist("ARICH/chDigit");
210 int nhot = 0;
211 double avgOcc = 0;
212 if (chDigit != NULL && nEvents != 0) {
213 avgOcc = chDigit->GetEntries() / 60480.;
214 if (avgOcc > 100.) {
215 for (int i = 0; i < chDigit->GetNbinsX(); i++) {
216 int nhit = chDigit->GetBinContent(i + 1);
217 if ((nhit - 3. * sqrt(nhit)) / float(nEvents) > m_hotLimit) nhot++;
218 }
219 }
220 }
221 setEpicsPV("hotChannels", nhot);
222
223 int ndeadHapd = 0;
224 TH1F* hapdDigit = (TH1F*)findHist("ARICH/hapdDigit");
225 if (hapdDigit != NULL && avgOcc * 144. > 100.) {
226 for (int i = 0; i < hapdDigit->GetNbinsX(); i++) {
227 if (hapdDigit->GetBinContent(i + 1) == 0) ndeadHapd++;
228 }
229 }
230 setEpicsPV("deadHAPDs", ndeadHapd);
231
232 if (avgOcc * 36. < 100.) {
233 setEpicsPV("badAPDs", 0);
234 return;
235 }
236
237 auto h_bitsPerChannel = getDelta("ARICH", "bitsPerChannel", 0, true);
238 if (h_bitsPerChannel == NULL) return;
239 TH1F* apdHits = new TH1F("apdHits", "nSigHits/nevt for all apds", 1680, -0.5, 1679.5);
240 double avg = 0.;
241 for (int i = 1; i < 60481; i++) {
242 int nnoise = h_bitsPerChannel->GetBinContent(1, i) + h_bitsPerChannel->GetBinContent(4, i);
243 int nsig = h_bitsPerChannel->GetBinContent(2, i) + h_bitsPerChannel->GetBinContent(3, i) - nnoise;
244 double sig2noise = float(nsig) / float(nsig + nnoise);
245 if (sig2noise > 0.2) avg += nsig;
246 if (i % 36 == 0 && i > 1) {
247 apdHits->SetBinContent(i / 36 + 1, avg);
248 avg = 0;
249 }
250 }
251
252 // for convenience I make this independent of the DB payloads... (it will never change)
253 double ringApdAvg[7] = {0.};
254
255 const int hapdInRing[7] = {42, 48, 54, 60, 66, 72, 78};
256
257 // average in hapd ring
258 for (int i = 1; i < 1681; i++) {
259 int hapd = (i - 1) / 4 + 1;
260 ringApdAvg[getRing(hapd)] += apdHits->GetBinContent(i);
261 }
262
263 for (int i = 0; i < 7; i++) {
264 ringApdAvg[i] /= float(4.*hapdInRing[i]);
265 }
266 int nbadApd = 0;
267 for (int i = 0; i < 1680; i++) {
268 int ring = getRing(i / 4 + 1);
269 int nhits = apdHits->GetBinContent(i + 1);
270 if (nhits + 3. * sqrt(nhits) < ringApdAvg[ring] * m_badApdOccLimit) nbadApd++;
271 }
272 setEpicsPV("badAPDs", nbadApd);
273
274 delete apdHits;
275
276}
277
278
280{
281 if (modID <= 42) return 0;
282 if (modID <= 90) return 1;
283 if (modID <= 144) return 2;
284 if (modID <= 204) return 3;
285 if (modID <= 270) return 4;
286 if (modID <= 342) return 5;
287 if (modID <= 420) return 6;
288 return -1; // -1 if invalid input
289}
290
292{
293 B2DEBUG(20, "DQMHistAnalysisARICH : endRun called");
294}
295
297{
298
299 B2DEBUG(20, "terminate called");
300}
ARICH histogram with HAPD plane 3 options for bin segmentation are available type 0 - one bin per HAP...
void fillFromTH1(TH1 *hist)
Fill the channelHist from the histogram Type 0 channelHist has to be filled with 420*144bin TH1 (each...
void setPoly(TH2Poly *poly)
Fill pure TH2Poly from ARICHChannelHist, makes bins and fills content.
double m_hotLimit
Occupancy limit for hot channels.
double m_badApdOccLimit
Occupancy threshold for bad APDs, in units of average APD occupancy.
void initialize() override final
Initialize the Module.
TCanvas * m_c_apdHist
Canvas for 2D hit map of APDs.
TH2Poly * m_apdPoly
hit map for each APD
TLine * m_LineForMB[5]
Lines to divide the sectors on mergerHit histogram.
TCanvas * m_c_mergerHit
Canvas for modified mergerHit histogram.
int getRing(int modID)
Returns ring number of HAPD with given moduleID.
bool m_enableAlert
Enable alert by base color of canvases.
Belle2::ARICHChannelHist * m_apdHist
ARICH TObject to draw hit map for each APD.
void terminate() override final
Termination action.
void event() override final
Event processor.
double m_sigBitFracLowAlarm
Alarm limit for overall signal/background fraction.
void endRun() override final
End-of-run action.
TCanvas * m_c_bits
Canvas for modified bits histogram.
std::vector< int > maskedMergers
The id numbers of masked merger boards to avoid unnecessary alert.
void beginRun() override final
begin-of-run action.
TCanvas * m_c_hitsPerEvent
Canvas for modified hitsPerTrack histogram.
double m_sigBitFracLowWarn
Warning limit for overall signal/background fraction.
int alertColor[3]
Alert color of canvases.
The base class for the histogram analysis module.
TCanvas * findCanvas(TString cname)
Find canvas by name.
void addDeltaPar(const std::string &dirname, const std::string &histname, HistDelta::EDeltaType t, int p, unsigned int a=1)
Add Delta histogram parameters.
void colorizeCanvas(TCanvas *canvas, EStatus status)
Helper function for Canvas colorization.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
double getEpicsPV(std::string keyname)
Read value from a EPICS PV.
TH1 * getDelta(const std::string &fullname, int n=0, bool onlyIfUpdated=true)
Get Delta histogram.
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
EStatus
Status flag of histogram/canvas.
@ c_StatusError
Analysis result: Severe issue found.
@ c_StatusWarning
Analysis result: Warning, there may be minor issues.
@ c_StatusGood
Analysis result: Good.
int registerEpicsPV(std::string pvname, std::string keyname="")
EPICS related Functions.
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
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
TH1F * hapd[6]
histogram of hits for each hapd
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
Abstract base class for different kinds of events.
STL namespace.