Belle II Software development
DQMHistAnalysisPXDReduction.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 : DQMHistAnalysisPXDReduction.cc
10// Description : Analysis of PXD Reduction
11//-
12
13
14#include <dqm/analysis/modules/DQMHistAnalysisPXDReduction.h>
15#include <TROOT.h>
16#include <TLatex.h>
17#include <vxd/geometry/GeoCache.h>
18
19using namespace std;
20using namespace Belle2;
21
22//-----------------------------------------------------------------
23// Register the Module
24//-----------------------------------------------------------------
25REG_MODULE(DQMHistAnalysisPXDReduction);
26
27//-----------------------------------------------------------------
28// Implementation
29//-----------------------------------------------------------------
30
33{
34 // This module CAN NOT be run in parallel!
35 setDescription("PXD DQM analysis module for ONSEN Data Reduction Monitoring");
36
37 // Parameter definition
38 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("PXDDAQ"));
39 addParam("LowerWarnLimit", m_meanLowerWarn, "Mean Reduction Low limit for warning", double(NAN)); // default is NAN =disable
40 addParam("LowerErrorLimit", m_meanLowerAlarm, "Mean Reduction Low limit for alarms", double(NAN)); // default is NAN =disable
41 addParam("UpperWarnLimit", m_meanUpperWarn, "Mean Reduction High limit for warning", double(NAN)); // default is NAN =disable
42 addParam("UpperErrorLimit", m_meanUpperAlarm, "Mean Reduction High limit for alarms", double(NAN)); // default is NAN =disable
43 addParam("minEntries", m_minEntries, "minimum number of new entries for last time slot", 1000);
44 addParam("excluded", m_excluded, "excluded module (indizes starting from 0 to 39)");
45 B2DEBUG(1, "DQMHistAnalysisPXDReduction: Constructor done.");
46}
47
49{
50}
51
53{
54 B2DEBUG(1, "DQMHistAnalysisPXDReduction: initialized.");
55
58
59 //collect the list of all PXD Modules in the geometry here
60 std::vector<VxdID> sensors = geo.getListOfSensors();
61 for (VxdID& aVxdID : sensors) {
62 VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
63 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
64 m_PXDModules.push_back(aVxdID); // reorder, sort would be better
65 }
66 std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
67
68 gROOT->cd(); // this seems to be important, or strange things happen
69
70 if (m_PXDModules.size() == 0) {
71 // Backup if no geometry is present (testing...)
72 B2WARNING("No PXDModules in Geometry found! Use hard-coded setup.");
73 std::vector <string> mod = {
74 "1.1.1", "1.1.2", "1.2.1", "1.2.2", "1.3.1", "1.3.2", "1.4.1", "1.4.2",
75 "1.5.1", "1.5.2", "1.6.1", "1.6.2", "1.7.1", "1.7.2", "1.8.1", "1.8.2",
76 "2.1.1", "2.1.2", "2.2.1", "2.2.2", "2.3.1", "2.3.2", "2.4.1", "2.4.2",
77 "2.5.1", "2.5.2", "2.6.1", "2.6.2", "2.7.1", "2.7.2", "2.8.1", "2.8.2",
78 "2.9.1", "2.9.2", "2.10.1", "2.10.2", "2.11.1", "2.11.2", "2.12.1", "2.12.2"
79 };
80 for (auto& it : mod) m_PXDModules.push_back(VxdID(it));
81 }
82
83 m_cReduction = new TCanvas((m_histogramDirectoryName + "/c_Reduction").data());
84 m_hReduction = new TH1F("hPXDReduction", "PXD Reduction; Module; Reduction", m_PXDModules.size(), 0, m_PXDModules.size());
85 m_hReduction->SetDirectory(0);// dont mess with it, this is MY histogram
86 m_hReduction->SetStats(false);
87 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
88 TString ModuleName = (std::string)m_PXDModules[i];
89 m_hReduction->GetXaxis()->SetBinLabel(i + 1, ModuleName);
91 "PXDDAQDHEDataReduction_" + (std::string)m_PXDModules[i])) addDeltaPar(m_histogramDirectoryName,
92 "PXDDAQDHEDataReduction_" + (std::string)m_PXDModules[i], HistDelta::c_Entries, m_minEntries,
93 1); // register delta
94 }
95 //Unfortunately this only changes the labels, but can't fill the bins by the VxdIDs
96 m_hReduction->Draw("");
98
99 m_meanLine = new TLine(0, 10, m_PXDModules.size(), 10);
100 m_meanUpperWarnLine = new TLine(0, 16, m_PXDModules.size(), 16);
101 m_meanLowerWarnLine = new TLine(0, 0.9, m_PXDModules.size(), 0.9);
102 m_meanUpperAlarmLine = new TLine(0, 20, m_PXDModules.size(), 20);
103 m_meanLowerAlarmLine = new TLine(0, 0.5, m_PXDModules.size(), 0.5);
104 m_meanLine->SetHorizontal(true);
105 m_meanLine->SetLineColor(kBlue);
106 m_meanLine->SetLineWidth(3);
107 m_meanUpperWarnLine->SetHorizontal(true);
108 m_meanUpperWarnLine->SetLineColor(c_ColorWarning + 2);
109 m_meanUpperWarnLine->SetLineWidth(3);
110 m_meanLowerWarnLine->SetHorizontal(true);
111 m_meanLowerWarnLine->SetLineColor(c_ColorWarning + 2);
112 m_meanLowerWarnLine->SetLineWidth(3);
113 m_meanUpperAlarmLine->SetHorizontal(true);
114 m_meanUpperAlarmLine->SetLineColor(c_ColorError + 2);
115 m_meanUpperAlarmLine->SetLineWidth(3);
116 m_meanLowerAlarmLine->SetHorizontal(true);
117 m_meanLowerAlarmLine->SetLineColor(c_ColorError + 2);
118 m_meanLowerAlarmLine->SetLineWidth(3);
119
120 registerEpicsPV("PXD:Red:Status", "Status");
121 registerEpicsPV("PXD:Red:Value", "Value");
122}
123
124
126{
127 B2DEBUG(1, "DQMHistAnalysisPXDReduction: beginRun called.");
128
129 m_cReduction->Clear();
130 m_hReduction->Reset(); // dont sum up!!!
132
133 // override with limits from EPICS. if they are set
135
136 if (!std::isnan(m_meanLowerAlarm)) {
139 }
140 if (!std::isnan(m_meanLowerWarn)) {
143 }
144 if (!std::isnan(m_meanUpperWarn)) {
147 }
148 if (!std::isnan(m_meanUpperAlarm)) {
151 }
152}
153
155{
156 if (!m_cReduction) return;
157
158 double ireduction = 0.0;
159 int ireductioncnt = 0;
160
161 bool anyupdate = false;
162 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
163 std::string name = "PXDDAQDHEDataReduction_" + (std::string)m_PXDModules[i ];
164 // std::replace( name.begin(), name.end(), '.', '_');
165
166 TH1* hh1 = getDelta(m_histogramDirectoryName, name);
167 // no initial sampling, we should get plenty of statistics
168 if (hh1) {
169 auto mean = hh1->GetMean();
170 m_hReduction->SetBinContent(i + 1, mean);
171 anyupdate = true;
172 }
173 }
174
175 if (!anyupdate) return; // nothing new -> no update
176
177 // calculate the mean of the mean
178 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
179 // ignore modules in exclude list
180 if (std::find(m_excluded.begin(), m_excluded.end(), i) != m_excluded.end()) continue;
181 auto mean = m_hReduction->GetBinContent(i + 1);
182 if (mean > 0) { // only for valid values
183 ireduction += mean; // well fit would be better
184 ireductioncnt++;
185 }
186 }
187
188 m_cReduction->cd();
189
190 double value = ireductioncnt > 0 ? ireduction / ireductioncnt : 0;
191
192 // if any if NaN, the comparison is false
193 auto stat_data = makeStatus(ireductioncnt >= 15,
194 value > m_meanUpperWarn || value < m_meanLowerWarn,
195 value > m_meanUpperAlarm || value < m_meanLowerAlarm);
196
197 if (m_hReduction) {
198 m_hReduction->Draw("");
199 if (stat_data != c_StatusTooFew) {
200 m_meanLine->SetY1(value);
201 m_meanLine->SetY2(value); // aka SetHorizontal
202 m_meanLine->Draw();
203 }
204 if (!std::isnan(m_meanLowerAlarm)) m_meanLowerAlarmLine->Draw();
205 if (!std::isnan(m_meanLowerWarn)) m_meanLowerWarnLine->Draw();
206 if (!std::isnan(m_meanUpperWarn)) m_meanUpperWarnLine->Draw();
207 if (!std::isnan(m_meanUpperAlarm)) m_meanUpperAlarmLine->Draw();
208 for (auto& it : m_excluded) {
209 static std::map <int, TLatex*> ltmap;
210 auto tt = ltmap[it];
211 if (!tt) {
212 tt = new TLatex(it + 0.5, 0, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
213 tt->SetTextSize(0.035);
214 tt->SetTextAngle(90);// Rotated
215 tt->SetTextAlign(12);// Centered
216 ltmap[it] = tt;
217 }
218 tt->Draw();
219 }
220 }
221
222 m_monObj->setVariable("reduction", value);
223
224 colorizeCanvas(m_cReduction, stat_data);
226 m_cReduction->Modified();
227 m_cReduction->Update();
228
229 // better only update if statistics is reasonable, we dont want "0" drops between runs!
230 setEpicsPV("Status", stat_data);
231 setEpicsPV("Value", value);
232}
233
235{
236 B2DEBUG(1, "DQMHistAnalysisPXDReduction: terminate called");
237
238 if (m_cReduction) delete m_cReduction;
239 if (m_hReduction) delete m_hReduction;
240
241 if (m_meanLine) delete m_meanLine;
246}
247
The base class for the histogram analysis module.
bool hasDeltaPar(const std::string &dirname, const std::string &histname)
Check if Delta histogram parameters exist for histogram.
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.
@ c_ColorWarning
Analysis result: Warning, there may be minor issues.
@ c_ColorError
Analysis result: Severe issue found.
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.
@ c_StatusTooFew
Not enough entries/event to judge.
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)
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
void terminate(void) override final
This method is called at the end of the event processing.
TH1F * m_hReduction
Histogram covering all modules.
TLine * m_meanUpperAlarmLine
Line in the Canvas to indicate limits.
TLine * m_meanLowerAlarmLine
Line in the Canvas to indicate limits.
double m_meanUpperAlarm
high error limit for alarm
void initialize(void) override final
Initializer.
MonitoringObject * m_monObj
Monitoring Object.
std::vector< VxdID > m_PXDModules
IDs of all PXD Modules to iterate over.
std::string m_histogramDirectoryName
name of histogram directory
TLine * m_meanLowerWarnLine
Line in the Canvas to indicate limits.
TLine * m_meanLine
Line in the Canvas to guide the eye.
TLine * m_meanUpperWarnLine
Line in the Canvas to indicate limits.
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.
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.