Belle II Software development
DQMHistAnalysisPXDCM.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// File : DQMHistAnalysisPXDCM.cc
10// Description : Analysis of PXD Common Modes
11//-
12
13
14#include <dqm/analysis/modules/DQMHistAnalysisPXDCM.h>
15#include <TROOT.h>
16#include <TLatex.h>
17#include <TPaveText.h>
18#include <vxd/geometry/GeoCache.h>
19#include <framework/core/ModuleParam.templateDetails.h>
20
21using namespace std;
22using namespace Belle2;
23
24//-----------------------------------------------------------------
25// Register the Module
26//-----------------------------------------------------------------
27REG_MODULE(DQMHistAnalysisPXDCM);
28
29//-----------------------------------------------------------------
30// Implementation
31//-----------------------------------------------------------------
32
35{
36 // This module CAN NOT be run in parallel!
37 setDescription("DQM Analysis for PXD Common Mode");
38
39 // Parameter definition
40 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("PXDDAQ"));
41 addParam("minEntries", m_minEntries, "minimum number of new entries for last time slot", 10000);
42
43 addParam("warnMean", m_warnMean, "warn level for peak position", 2.0);
44 addParam("errorMean", m_errorMean, "error level for peak position", 3.0);
45 addParam("warnOutside", m_warnOutside, "warn level for outside fraction", 1e-5);
46 addParam("errorOutside", m_errorOutside, "error level for outside fraction", 1e-4);
47 addParam("upperLine", m_upperLine, "upper threshold and line for outside fraction", 17);
48
49 addParam("gateMaskModuleList", m_parModuleList, "Module List for Gate Masking");
50 addParam("gateMaskGateList", m_parGateList, "Gate List for Gate Masking");
51 addParam("excluded", m_excluded, "excluded module (indizes starting from 0 to 39)");
52
53 B2DEBUG(99, "DQMHistAnalysisPXDCM: Constructor done.");
54}
55
57{
60
61 // collect the list of all PXD Modules in the geometry here
62 std::vector<VxdID> sensors = geo.getListOfSensors();
63 for (VxdID& aVxdID : sensors) {
64 VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
65 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
66 m_PXDModules.push_back(aVxdID); // reorder, sort would be better
67 }
68 std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
69
70 gROOT->cd(); // this seems to be important, or strange things happen
71
72 if (m_PXDModules.size() == 0) {
73 // Backup if no geometry is present (testing...)
74 B2WARNING("No PXDModules in Geometry found! Use hard-coded setup.");
75 std::vector <string> mod = {
76 "1.1.1", "1.1.2", "1.2.1", "1.2.2", "1.3.1", "1.3.2", "1.4.1", "1.4.2",
77 "1.5.1", "1.5.2", "1.6.1", "1.6.2", "1.7.1", "1.7.2", "1.8.1", "1.8.2",
78 "2.1.1", "2.1.2", "2.2.1", "2.2.2", "2.3.1", "2.3.2", "2.4.1", "2.4.2",
79 "2.5.1", "2.5.2", "2.6.1", "2.6.2", "2.7.1", "2.7.2", "2.8.1", "2.8.2",
80 "2.9.1", "2.9.2", "2.10.1", "2.10.2", "2.11.1", "2.11.2", "2.12.1", "2.12.2"
81 };
82 for (auto& it : mod) m_PXDModules.push_back(VxdID(it));
83 }
84 m_cCommonModeDelta = new TCanvas((m_histogramDirectoryName + "/c_CommonModeDelta").data());
85
86 m_hCommonModeDelta = new TH2D("hPXDCommonMode", "PXD CommonMode ; Module; CommonMode", m_PXDModules.size(), 0,
87 m_PXDModules.size(), 63, 0, 63);
88 m_hCommonModeDelta->SetDirectory(0);// dont mess with it, this is MY histogram
89 m_hCommonModeDelta->SetStats(false);
90
91 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
92 TString ModuleName = (std::string)m_PXDModules[i];
93 m_hCommonModeDelta->GetXaxis()->SetBinLabel(i + 1, ModuleName);
94 addDeltaPar(m_histogramDirectoryName, "PXDDAQCM_" + (std::string)m_PXDModules[i], HistDelta::c_Underflow, m_minEntries,
95 1); // register delta
96 }
97 //Unfortunately this only changes the labels, but can't fill the bins by the VxdIDs
98 m_hCommonModeDelta->Draw("colz");
99
101
102 if (m_parModuleList.size() != m_parGateList.size()) {
103 B2FATAL("Parameter list need same length");
104 return;
105 }
106 for (size_t i = 0; i < m_parModuleList.size(); i++) {
107 for (auto n : m_parGateList[i]) {
108 m_maskedGates[VxdID(m_parModuleList[i])].push_back(n);
109 }
110 }
111
113 m_line10 = new TLine(0, 10, m_PXDModules.size(), 10);
114 m_lineOutside = new TLine(0, m_upperLine, m_PXDModules.size(), m_upperLine);
115 m_line10->SetHorizontal(true);
116 m_line10->SetLineColor(3);// Green
117 m_line10->SetLineWidth(3);
118 m_lineOutside->SetHorizontal(true);
119 m_lineOutside->SetLineColor(1);// Black
120 m_lineOutside->SetLineWidth(3);
121
122
123 registerEpicsPV("PXD:CommonMode:Status", "Status");
124 registerEpicsPV("PXD:CommonMode:Outside", "Outside");
125 registerEpicsPV("PXD:CommonMode:CM63", "CM63");
126 //registerEpicsPV("PXD:CommonMode:CM62", "CM62");
127
128 for (VxdID& aPXDModule : m_PXDModules) {
129 auto buff = (std::string)aPXDModule;
130 replace(buff.begin(), buff.end(), '.', '_');
131 registerEpicsPV("PXD:CommonMode:Mean:" + buff, (std::string)aPXDModule);
132 }
133 B2DEBUG(99, "DQMHistAnalysisPXDCM: initialized.");
134}
135
137{
138 B2DEBUG(99, "DQMHistAnalysisPXDCM: beginRun called.");
139
140 m_cCommonModeDelta->Clear();
141 m_cCommonModeDelta->SetLogz();
142
143 // this is needed at least for the "Old" and "Delta" one or update doesn't work
144 m_hCommonModeDelta->Reset();
145}
146
148{
149 double all_outside = 0.0, all = 0.0;
150 double all_cm = 0.0;
151 bool error_flag = false;
152 bool warn_flag = false;
153 bool anyupdate = false;
154
155 static TPaveText* leg = nullptr;
156
157 if (leg == nullptr) {
158 leg = new TPaveText(0.1, 0.6, 0.90, 0.95, "NDC");
159 leg->SetFillStyle(0);
160 leg->SetBorderSize(0);
161 } else {
162 leg->Clear();
163 }
164
165 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
166 auto modname = (std::string)m_PXDModules[i];
167 std::string name = "PXDDAQCM_" + modname;
168 bool excluded = find(m_excluded.begin(), m_excluded.end(), i) != m_excluded.end();
169
170 auto hh1 = getDelta(m_histogramDirectoryName, name); // default, only updated
171 if (hh1) {
172 auto scale = hh1->GetBinContent(0); // misuse underflow as event counter
173 bool update = scale >= m_minEntries ; // filter initial sampling
174 anyupdate |= update;
175 if (scale > 0) scale = 1.0 / scale;
176 else scale = 1.; // worst case, no events at run start
177
178 auto& gm = m_maskedGates[m_PXDModules[i]];
179 // We loop over a 2d histogram!
180 // loop CM values
181 for (int bin = 1; bin <= 64; bin++) { // including CM63!!!
182 // loop gates*asics
183 double v = 0;
184 for (int gate = 0; gate < 192; gate++) {
185 // attention, gate is not bin nr!
186 if (std::find(gm.begin(), gm.end(), gate) == gm.end()) {
187 v += hh1->GetBinContent(hh1->GetBin(gate + 1 + 192 * 0, bin)) +
188 hh1->GetBinContent(hh1->GetBin(gate + 1 + 192 * 1, bin)) +
189 hh1->GetBinContent(hh1->GetBin(gate + 1 + 192 * 2, bin)) +
190 hh1->GetBinContent(hh1->GetBin(gate + 1 + 192 * 3, bin));
191 }
192 }
193 // integration intervals depend on CM default value, this seems to be agreed =10
194 // FIXME currently we have to much noise below the line ... thus excluding this to avoid false alarms
195 // outside_full += hh1->Integral(1 /*0*/, 5); /// FIXME we exclude bin 0 as we use it for debugging/timing pixels
196 // attention, n bins!
197 // we integrate up including value 62 (cm overflow), but not 63 (fifo full)
198 if (!excluded) {
199 if (bin == 63 + 1) { // CM63
200 all_cm += v;
201 } else { // excluding CM63
202 all += v;
203 if (bin > m_upperLine + 1) all_outside += v;
204 }
205 }
206 m_hCommonModeDelta->SetBinContent(i + 1, bin, v * scale); // attention, mixing bin nr and index
207 }
208
209 if (update) {
210 Double_t mean = 0.;
211 Double_t entries = 0.;
212 Double_t outside = 0.;
213
214 // Attention, Bins
215 // we do not need to re-scale it as the scale is the same for all bins
216 for (int cm_y = 0; cm_y < m_upperLine; cm_y++) {
217 auto v = m_hCommonModeDelta->GetBinContent(m_hCommonModeDelta->GetBin(i + 1, cm_y + 1));
218 entries += v;
219 mean += v * (cm_y + 1);
220 }
221 // Attention, Bins
222 // We ignore CM63 in outside and overall count
223 for (int cm_y = m_upperLine; cm_y < 63; cm_y++) {
224 auto v = m_hCommonModeDelta->GetBinContent(m_hCommonModeDelta->GetBin(i + 1, cm_y + 1));
225 entries += v;
226 outside += v;
227 }
228 if (entries > 0 and scale < 1e-3) { // ignore modules with minimum events
229 // scale <1e-3 == >1000 events
230 mean /= entries; // calculate mean
231 auto warn_tmp_m = fabs(10.0 - mean) > m_warnMean;
232 auto err_tmp_m = fabs(10.0 - mean) > m_errorMean;
233 auto warn_tmp_os = outside / entries > m_warnOutside;
234 auto err_tmp_os = outside / entries > m_errorOutside;
235 if (not excluded) {
236 warn_flag |= warn_tmp_m or warn_tmp_os;
237 error_flag |= err_tmp_m or err_tmp_os;
238
239 if (warn_tmp_m or err_tmp_m) {
240 TString tmp;
241 tmp.Form("%s: Mean %f", modname.c_str(), mean);
242 leg->AddText(tmp);
243 B2INFO(name << " Mean " << mean << " " << warn_tmp_m << err_tmp_m);
244 }
245 if (warn_tmp_os or err_tmp_os) {
246 TString tmp;
247 tmp.Form("%s: Outside %f %%", modname.c_str(), 100. * outside / entries);
248 leg->AddText(tmp);
249 B2INFO(name << " Outside " << outside / entries << " (" << outside << "/" << entries << ") " << warn_tmp_os
250 << err_tmp_os);
251 }
252 }
253 m_monObj->setVariable(("cm_" + modname).c_str(), mean);
254
255 setEpicsPV((std::string)m_PXDModules[i], mean);
256 }
257 }
258 }
259 }
260
261 auto status = makeStatus(all >= 10000, warn_flag, error_flag);
262
263 m_cCommonModeDelta->cd();
265
266 if (anyupdate) {
267 double dataoutside = all > 0 ? (all_outside / all) : 0;
268 double datacm = all > 0 ? (all_cm / all) : 0;
269 setEpicsPV("Status", status);
270 setEpicsPV("Outside", dataoutside);
271 setEpicsPV("CM63", datacm);
272 }
273 if (m_hCommonModeDelta) {
274 m_hCommonModeDelta->Draw("colz");
275 leg->Draw();
276 m_line10->Draw();
277 m_lineOutside->Draw();
278 for (auto& it : m_excluded) {
279 static std::map <int, TLatex*> ltmap;
280 auto tt = ltmap[it];
281 if (!tt) {
282 tt = new TLatex(it + 0.5, 0, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
283 tt->SetTextSize(0.035);
284 tt->SetTextAngle(90);// Rotated
285 tt->SetTextAlign(12);// Centered
286 ltmap[it] = tt;
287 }
288 tt->Draw();
289 }
290 }
291
293 m_cCommonModeDelta->Modified();
294 m_cCommonModeDelta->Update();
295}
296
298{
299 B2DEBUG(99, "DQMHistAnalysisPXDCM: terminate called");
300
303 if (m_line10) delete m_line10;
304 if (m_lineOutside) delete m_lineOutside;
305}
306
The base class for the histogram analysis module.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
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.
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 makeStatus(bool enough, bool warn_flag, bool error_flag)
Helper function to judge the status for coloring and EPICS.
int registerEpicsPV(std::string pvname, std::string keyname="")
EPICS related Functions.
void UpdateCanvas(std::string name, bool updated=true)
Mark canvas as updated (or not)
void terminate(void) override final
This method is called at the end of the event processing.
int m_minEntries
Update entry interval.
TLine * m_line10
Line in the Canvas to guide the eye, target CM.
double m_warnOutside
warn level for outside fraction
TCanvas * m_cCommonModeDelta
Final Canvas.
int m_upperLine
threshold level/line for outside fraction
void initialize(void) override final
Initializer.
TLine * m_lineOutside
Line in the Canvas to guide the eye, outside boundary.
MonitoringObject * m_monObj
Monitoring Object.
TH2D * m_hCommonModeDelta
histogram covering all modules
std::vector< VxdID > m_PXDModules
IDs of all PXD Modules to iterate over.
std::string m_histogramDirectoryName
name of histogram directory
double m_errorOutside
error level for outside fraction
double m_errorMean
error level for mean
std::map< VxdID, std::vector< int > > m_maskedGates
Module wise gate masking in CM plot and alarm.
double m_warnMean
warn level for mean
std::vector< std::vector< int > > m_parGateList
Gate list for masking.
std::vector< int > m_excluded
Indizes of excluded PXD Modules.
void beginRun(void) override final
Called when entering a new run.
void event(void) override final
This method is called for each event.
std::vector< std::string > m_parModuleList
Module list for masking.
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)
void addCanvas(TCanvas *canv)
Add Canvas to monitoring object.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:59
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
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
Abstract base class for different kinds of events.
STL namespace.