Belle II Software  release-08-01-10
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 
17 using namespace std;
18 using namespace Belle2;
19 
20 //-----------------------------------------------------------------
21 // Register the Module
22 //-----------------------------------------------------------------
23 REG_MODULE(DQMHistAnalysisDeltaEpicsMonObjExample);
24 
25 //-----------------------------------------------------------------
26 // Implementation
27 //-----------------------------------------------------------------
28 
29 DQMHistAnalysisDeltaEpicsMonObjExampleModule::DQMHistAnalysisDeltaEpicsMonObjExampleModule()
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  // destructor not needed
46  // EPICS singleton deletion not urgent -> can be done by framework
47 }
48 
50 {
51  B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: initialized.");
52 
53  gROOT->cd(); // this seems to be important before creating new Canvas, or strange things happen
54  m_canvas = new TCanvas((m_histogramDirectoryName + "/c_Test").data());
55 
56  m_monObj = getMonitoringObject("test");
58 
59  // delta parameters, two examples:
61  100000, 1); // update each 100k events (from daq histogram)
62  // addDeltaPar(m_histogramDirectoryName, m_histogramName, HistDelta::c_Entries, 10000, 1); // update each 10000 entries
63 
64  registerEpicsPV(m_pvPrefix + "mean", "mean");
65  registerEpicsPV(m_pvPrefix + "width", "width");
66  registerEpicsPV(m_pvPrefix + "alarm", "alarm");
68  5.0); // -> now trigger update. this may be optional, framework can take care unless we want to now the result immediately
69 
70  m_meanLowerAlarmLine = new TLine();
71  m_meanLowerAlarmLine->SetLineColor(kRed + 3);
72  m_meanLowerAlarmLine->SetLineWidth(3);
73  m_meanLowerAlarmLine->SetVertical(true);
74  m_meanLowerWarnLine = new TLine();
75  m_meanLowerWarnLine->SetLineColor(kOrange - 3);
76  m_meanLowerWarnLine->SetLineWidth(3);
77  m_meanLowerWarnLine->SetVertical(true);
78  m_meanUpperWarnLine = new TLine();
79  m_meanUpperWarnLine->SetLineColor(kOrange - 3);
80  m_meanUpperWarnLine->SetLineWidth(3);
81  m_meanUpperWarnLine->SetVertical(true);
82  m_meanUpperAlarmLine = new TLine();
83  m_meanUpperAlarmLine->SetLineColor(kRed + 3);
84  m_meanUpperAlarmLine->SetLineWidth(3);
85  m_meanUpperAlarmLine->SetVertical(true);
86 }
87 
89 {
90  B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: beginRun called.");
91 
92  double unused = NAN;
93  m_meanLowerAlarm = -1.0; // this is some way to set default values, which are used
94  m_meanLowerWarn = -0.9; // when none are set in epics PV directly.
95  m_meanUpperWarn = 0.9;
96  m_meanUpperAlarm = 1.0;
98 
99  if (!std::isnan(m_meanLowerAlarm)) {
102  }
103  if (!std::isnan(m_meanLowerWarn)) {
106  }
107  if (!std::isnan(m_meanUpperWarn)) {
110  }
111  if (!std::isnan(m_meanUpperAlarm)) {
114  }
115 
117 }
118 
120 {
121  B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: endRun called.");
122  // MiraBelle export code should run at end of Run
123  // but it still "remembers" the state from last event call.
124  doHistAnalysis(true);
125 }
126 
128 {
129  B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: event called.");
130  doHistAnalysis(false);
131 }
132 
134 {
135  double data_mean = 0.0;
136  double data_width = 0.0;
137  int data_alarm = 0;
138 
139  m_canvas->Clear();
140  m_canvas->cd(0);
141 
142  // more handy to have it as a full name, but find/get functions
143  // can work with one or two parameters
144  // std::string fullname = m_histogramDirectoryName + "/" + m_histogramName;
145 
146  // In the following, enable one
147 
148  // get most recent delta, but only if updated since last call
149  auto hist = getDelta(m_histogramDirectoryName, m_histogramName, 0, true);// only if updated
150  // get basic histogram (run integrated up) but only if updated since last call
151  // auto hist = findHist(m_histogramDirectoryName, m_histogramName, true);// only if updated
152 
153  // the following cases do not make sense, unless you want to archieve something special.
154  // get most recent delta even if not updated
155  // auto hist = getDelta(m_histogramDirectoryName, m_histogramName, 0, false);// even if no update
156  // get basic histogram (run integrated up) even if not updated
157  // auto hist = findHist(m_histogramDirectoryName, m_histogramName, false);// even if no update
158 
159  if (hist) {
160  data_mean = hist->GetMean();
161  data_width = hist->GetRMS();
162 
163  data_alarm = 2; // good, green
164 
165  // your conditions for warn and alarm
166  // here we would have to check agains NaN, too ... just too lazy to write the code in this example
167  if (data_mean > m_meanUpperWarn || data_mean < m_meanLowerWarn || data_width > m_widthUpperWarn) data_alarm = 3; // warning, yellow
168  if (data_mean > m_meanUpperAlarm || data_mean < m_meanLowerAlarm || data_width > m_widthUpperAlarm) data_alarm = 4; // error, red
169  // maybe other condition, e.g. too low to judge -> 0
170 
171  hist->Draw("hist");
172 
173  // line position need to be adjusted to histogram content
174  double y1 = hist->GetMaximum();
175  double y2 = hist->GetMinimum();
176  auto ym = y1 + (y1 - y2) * 0.05;
177  m_meanLowerAlarmLine->SetY1(ym);
178  m_meanLowerWarnLine->SetY1(ym);
179  m_meanUpperWarnLine->SetY1(ym);
180  m_meanUpperAlarmLine->SetY1(ym);
181  auto yn = y2 + (y1 - y2) * 0.05;
182  m_meanLowerAlarmLine->SetY2(yn);
183  m_meanLowerWarnLine->SetY2(yn);
184  m_meanUpperWarnLine->SetY2(yn);
185  m_meanUpperAlarmLine->SetY2(yn);
186 
187  if (!std::isnan(m_meanLowerAlarm)) m_meanLowerAlarmLine->Draw();
188  if (!std::isnan(m_meanLowerWarn)) m_meanLowerWarnLine->Draw();
189  if (!std::isnan(m_meanUpperWarn)) m_meanUpperWarnLine->Draw();
190  if (!std::isnan(m_meanUpperAlarm)) m_meanUpperAlarmLine->Draw();
191 
192  // the following will be replaced by a base class function
193  auto color = kWhite;
194  switch (data_alarm) {
195  case 2:
196  color = kGreen;
197  break;
198  case 3:
199  color = kYellow;
200  break;
201  case 4:
202  color = kRed;
203  break;
204  default:
205  color = kGray;
206  break;
207  }
208  m_canvas->Pad()->SetFillColor(color);
209  }
210 
211  // Tag canvas as updated ONLY if things have changed.
212  UpdateCanvas(m_canvas->GetName(), hist != nullptr);
213  // Remark: if you do not tag, you may have an empty Canvas (as it is cleared above)
214 
215  // this if left over from jsroot, may not be needed anymore (to check)
216  m_canvas->Update();
217 
218  // actually, we would prefer to only update the epics variable
219  // if the main histogram or delta~ has changed.
220 
221  if (hist) {
222  if (forMiraBelle) {
223  // currently, monObj are only used at end of run, thus this is not necessary
224  // but, if in the future we want to use finer granularity, this code needs change.
225  m_monObj->setVariable("mean", data_mean);
226  m_monObj->setVariable("width", data_width);
227  } else {
228  // we do not want to fill epics again at end of run, as this would be identical to the last event called before.
229  setEpicsPV("mean", data_mean);
230  setEpicsPV("width", data_width);
231  setEpicsPV("alarm", data_alarm);
232  // until now, chnages are buffered only locally, not written out to network for performance reason
234  5.0); // -> now trigger update. this may be optional, framework can take care unless we want to now the result immediately
235  }
236  }
237 }
238 
240 {
241  B2DEBUG(1, "DQMHistAnalysisDeltaEpicsMonObjExample: terminate called");
242 }
243 
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.
The base class for the histogram analysis module.
int registerEpicsPV(std::string pvname, std::string keyname="", bool update_pvs=true)
EPICS related Functions.
void addDeltaPar(const std::string &dirname, const std::string &histname, HistDelta::EDeltaType t, int p, unsigned int a=1)
Add Delta histogram parameters.
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.
void UpdateCanvas(std::string name, bool updated=true)
Mark canvas as updated (or not)
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
int updateEpicsPVs(float timeout)
Update all EPICS PV (flush to network)
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.
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.