Belle II Software prerelease-11-00-00a
DQMHistAnalysisDeltaEpicsMonObjExample.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 : DQMHistAnalysisDeltaEpicsMonObjExample.cc
10// Description : Test Module for Delta Histogram Access
11//-
12
13
14#include <dqm/analysis/modules/DQMHistAnalysisDeltaEpicsMonObjExample.h>
15#include <TROOT.h>
16
17using namespace std;
18using namespace Belle2;
19
20//-----------------------------------------------------------------
21// Register the Module
22//-----------------------------------------------------------------
23REG_MODULE(DQMHistAnalysisDeltaEpicsMonObjExample);
24
25//-----------------------------------------------------------------
26// Implementation
27//-----------------------------------------------------------------
28
31{
32 // This module CAN NOT be run in parallel!
33 setDescription("Example module for using delta histogramming and sending values to EPICS");
34
35 // Parameter definition
36 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("test"));
37 addParam("histogramName", m_histogramName, "Name of Histogram", std::string("testHist"));
38 addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("DQM:TEST:"));
39 B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: Constructor done.");
40
41}
42
44{
45 B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: initialized.");
46
47 gROOT->cd(); // this seems to be important before creating new Canvas, or strange things happen
48 m_canvas = new TCanvas((m_histogramDirectoryName + "/c_Test").data());
49
51 m_monObj->addCanvas(m_canvas);
52
53 // delta parameters, two examples:
55 100000, 1); // update each 100k events (from daq histogram)
56 // addDeltaPar(m_histogramDirectoryName, m_histogramName, HistDelta::c_Entries, 10000, 1); // update each 10000 entries
57
58 registerEpicsPV(m_pvPrefix + "mean", "mean");
59 registerEpicsPV(m_pvPrefix + "width", "width");
60 registerEpicsPV(m_pvPrefix + "alarm", "alarm");
61
62 m_meanLowerAlarmLine = new TLine();
63 m_meanLowerAlarmLine->SetLineColor(kRed + 3);
64 m_meanLowerAlarmLine->SetLineWidth(3);
65 m_meanLowerAlarmLine->SetVertical(true);
66 m_meanLowerWarnLine = new TLine();
67 m_meanLowerWarnLine->SetLineColor(kOrange - 3);
68 m_meanLowerWarnLine->SetLineWidth(3);
69 m_meanLowerWarnLine->SetVertical(true);
70 m_meanUpperWarnLine = new TLine();
71 m_meanUpperWarnLine->SetLineColor(kOrange - 3);
72 m_meanUpperWarnLine->SetLineWidth(3);
73 m_meanUpperWarnLine->SetVertical(true);
74 m_meanUpperAlarmLine = new TLine();
75 m_meanUpperAlarmLine->SetLineColor(kRed + 3);
76 m_meanUpperAlarmLine->SetLineWidth(3);
77 m_meanUpperAlarmLine->SetVertical(true);
78}
79
81{
82 B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: beginRun called.");
83
84 double unused = NAN;
85 m_meanLowerAlarm = -1.0; // this is some way to set default values, which are used
86 m_meanLowerWarn = -0.9; // when none are set in epics PV directly.
87 m_meanUpperWarn = 0.9;
88 m_meanUpperAlarm = 1.0;
90
91 if (!std::isnan(m_meanLowerAlarm)) {
94 }
95 if (!std::isnan(m_meanLowerWarn)) {
98 }
99 if (!std::isnan(m_meanUpperWarn)) {
102 }
103 if (!std::isnan(m_meanUpperAlarm)) {
106 }
107
109}
110
112{
113 B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: endRun called.");
114 // MiraBelle export code should run at end of Run
115 // but it still "remembers" the state from last event call.
116 doHistAnalysis(true);
117}
118
120{
121 B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: event called.");
122 doHistAnalysis(false);
123}
124
126{
127 double data_mean = 0.0;
128 double data_width = 0.0;
129 EStatus status = c_StatusDefault;
130
131 m_canvas->Clear();
132 m_canvas->cd(0);
133
134 // more handy to have it as a full name, but find/get functions
135 // can work with one or two parameters
136 // std::string fullname = m_histogramDirectoryName + "/" + m_histogramName;
137
138 // In the following, enable one
139
140 // get most recent delta, but only if updated since last call
141 auto hist = getDelta(m_histogramDirectoryName, m_histogramName, 0, true);// only if updated
142 // get basic histogram (run integrated up) but only if updated since last call
143 // auto hist = findHist(m_histogramDirectoryName, m_histogramName, true);// only if updated
144
145 // the following cases do not make sense, unless you want to achieve something special.
146 // get most recent delta even if not updated
147 // auto hist = getDelta(m_histogramDirectoryName, m_histogramName, 0, false);// even if no update
148 // get basic histogram (run integrated up) even if not updated
149 // auto hist = findHist(m_histogramDirectoryName, m_histogramName, false);// even if no update
150
151 if (hist) {
152 data_mean = hist->GetMean();
153 data_width = hist->GetRMS();
154
155 // your conditions for warn and alarm
156 // here we would have to check against NaN, too ... just too lazy to write the code in this example
157 // maybe other condition, e.g. too low to judge -> 0
158 status = makeStatus(true,
159 data_mean > m_meanUpperWarn || data_mean < m_meanLowerWarn || data_width > m_widthUpperWarn,
160 data_mean > m_meanUpperAlarm || data_mean < m_meanLowerAlarm || data_width > m_widthUpperAlarm);
161
162 hist->Draw("hist");
163
164 // line position need to be adjusted to histogram content
165 double y1 = hist->GetMaximum();
166 double y2 = hist->GetMinimum();
167 auto ym = y1 + (y1 - y2) * 0.05;
168 m_meanLowerAlarmLine->SetY1(ym);
169 m_meanLowerWarnLine->SetY1(ym);
170 m_meanUpperWarnLine->SetY1(ym);
171 m_meanUpperAlarmLine->SetY1(ym);
172 auto yn = y2 + (y1 - y2) * 0.05;
173 m_meanLowerAlarmLine->SetY2(yn);
174 m_meanLowerWarnLine->SetY2(yn);
175 m_meanUpperWarnLine->SetY2(yn);
176 m_meanUpperAlarmLine->SetY2(yn);
177
178 if (!std::isnan(m_meanLowerAlarm)) m_meanLowerAlarmLine->Draw();
179 if (!std::isnan(m_meanLowerWarn)) m_meanLowerWarnLine->Draw();
180 if (!std::isnan(m_meanUpperWarn)) m_meanUpperWarnLine->Draw();
181 if (!std::isnan(m_meanUpperAlarm)) m_meanUpperAlarmLine->Draw();
182
183 colorizeCanvas(m_canvas, status);
184 }
185
186 // Tag canvas as updated ONLY if things have changed.
187 UpdateCanvas(m_canvas->GetName(), hist != nullptr);
188 // Remark: if you do not tag, you may have an empty Canvas (as it is cleared above)
189
190 // this if left over from jsroot, may not be needed anymore (to check)
191 m_canvas->Update();
192
193 // actually, we would prefer to only update the epics variable
194 // if the main histogram or delta~ has changed.
195
196 if (hist) {
197 if (forMiraBelle) {
198 // currently, monObj are only used at end of run, thus this is not necessary
199 // but, if in the future we want to use finer granularity, this code needs change.
200 m_monObj->setVariable("mean", data_mean);
201 m_monObj->setVariable("width", data_width);
202 } else {
203 // we do not want to fill epics again at end of run, as this would be identical to the last event called before.
204 setEpicsPV("mean", data_mean);
205 setEpicsPV("width", data_width);
206 setEpicsPV("alarm", int(status));
207 // until now, changes are buffered only locally, not waiting for them being written out to performance reason
208 // framework will take care after all analysis modules are done
209 }
210 }
211}
212
214{
215 B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: terminate called");
216}
217
void terminate(void) override final
This method is called at the end of the event processing.
void beginRun(void) override final
Called when entering a new run.
void event(void) override final
This method is called for each event.
void doHistAnalysis(bool forMiraBelle=false)
Do the actual processing.
int registerEpicsPV(const std::string &pvname, const std::string &keyname="")
EPICS related Functions.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
static void addDeltaPar(const std::string &dirname, const std::string &histname, HistDelta::EDeltaType t, int p, unsigned int a=1)
Add Delta histogram parameters.
static void colorizeCanvas(TCanvas *canvas, EStatus status)
Helper function for Canvas colorization.
static void UpdateCanvas(const std::string &name, bool updated=true)
Mark canvas as updated (or not)
TH1 * getDelta(const std::string &fullname, int n=0, bool onlyIfUpdated=true)
Get Delta histogram.
DQMHistAnalysisModule()
Constructor / Destructor.
EStatus
Status flag of histogram/canvas.
@ c_StatusDefault
default for non-coloring
static EStatus makeStatus(bool enough, bool warn_flag, bool error_flag)
Helper function to judge the status for coloring and EPICS.
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
void setEpicsPV(const std::string &keyname, double value)
Write value to a EPICS PV.
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
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
Abstract base class for different kinds of events.
STL namespace.