Belle II Software development
DQMHistAnalysis.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 : DQMHistAnalysisModule.cc
10// Description : Baseclass for DQM histogram analysis module
11//-
12
13#include <dqm/core/DQMHistAnalysis.h>
14#include <boost/algorithm/string.hpp>
15#include <TROOT.h>
16#include <TClass.h>
17
18using namespace std;
19using namespace Belle2;
20
21//-----------------------------------------------------------------
22// Register the Module
23//-----------------------------------------------------------------
24REG_MODULE(DQMHistAnalysis);
25
26//-----------------------------------------------------------------
27// Implementation
28//-----------------------------------------------------------------
29
35#ifdef _BELLE2_EPICS
36std::vector <chid> DQMHistAnalysisModule::m_epicsChID;
37#endif
38
39bool DQMHistAnalysisModule::m_useEpics = false; // default to false, to enable EPICS, add special EPICS Module class into chain
41 false; // special for second "online" use (reading limits). default to false, to enable EPICS, add special EPICS Module parameter
42std::string DQMHistAnalysisModule::m_PVPrefix = "TEST:"; // default to "TEST:", for production, set in EPICS enabler to e.g. "DQM:"
43
45{
46 //Set module properties
47 setDescription("Histogram Analysis module base class");
48}
49
51{
52 s_histList.clear();
53 s_refList.clear();
54// s_monObjList;
55 s_deltaList.clear();
56 s_canvasUpdatedList.clear();
57}
58
59bool DQMHistAnalysisModule::addHist(const std::string& dirname, const std::string& histname, TH1* h)
60{
61 std::string fullname;
62 if (dirname.size() > 0) {
63 fullname = dirname + "/" + histname;
64 } else {
65 fullname = histname;
66 }
67
68 if (s_histList[fullname].update(h)) {
69 // only if histogram changed, check if delta histogram update needed
70 auto it = s_deltaList.find(fullname);
71 if (it != s_deltaList.end()) {
72 B2DEBUG(20, "Found Delta" << fullname);
73 it->second.update(h); // update
74 }
75 return true; // histogram changed
76 }
77
78 return false; // histogram didn't change
79}
80
81// void DQMHistAnalysisModule::addRef(const std::string& dirname, const std::string& histname, TH1* ref)
82// {
83// std::string fullname;
84// if (dirname.size() > 0) {
85// fullname = dirname + "/" + histname;
86// } else {
87// fullname = histname;
88// }
89// auto it = s_refList.find(fullname);
90// if (it == s_refList.end()) {
91// B2DEBUG(1, "Did not find histogram " << fullname << "in s_refList, so inserting now.");
92// s_refList.insert({fullname, ref});
93// }
94// }
95
96void DQMHistAnalysisModule::addDeltaPar(const std::string& dirname, const std::string& histname, HistDelta::EDeltaType t, int p,
97 unsigned int a)
98{
99 std::string fullname;
100 if (dirname.size() > 0) {
101 fullname = dirname + "/" + histname;
102 } else {
103 fullname = histname;
104 }
105 s_deltaList[fullname].set(t, p, a);
106}
107
108bool DQMHistAnalysisModule::hasDeltaPar(const std::string& dirname, const std::string& histname)
109{
110 std::string fullname;
111 if (dirname.size() > 0) {
112 fullname = dirname + "/" + histname;
113 } else {
114 fullname = histname;
115 }
116 return s_deltaList.find(fullname) != s_deltaList.end(); // contains() if we switch to C++20
117}
118
119TH1* DQMHistAnalysisModule::getDelta(const std::string& dirname, const std::string& histname, int n, bool onlyIfUpdated)
120{
121 std::string fullname;
122 if (dirname.size() > 0) {
123 fullname = dirname + "/" + histname;
124 } else {
125 fullname = histname;
126 }
127 return getDelta(fullname, n, onlyIfUpdated);
128}
129
130TH1* DQMHistAnalysisModule::getDelta(const std::string& fullname, int n, bool onlyIfUpdated)
131{
132 auto it = s_deltaList.find(fullname);
133 if (it != s_deltaList.end()) {
134 return it->second.getDelta(n, onlyIfUpdated);
135 }
136 B2WARNING("Delta hist " << fullname << " not found");
137 return nullptr;
138}
139
141{
142 auto obj = &s_monObjList[objName];
143 obj->SetName(objName.c_str());
144 return obj;
145}
146
147TCanvas* DQMHistAnalysisModule::findCanvas(TString canvas_name)
148{
149 TIter nextkey(gROOT->GetListOfCanvases());
150 TObject* obj{};
151
152 while ((obj = dynamic_cast<TObject*>(nextkey()))) {
153 if (obj->IsA()->InheritsFrom("TCanvas")) {
154 if (obj->GetName() == canvas_name)
155 return dynamic_cast<TCanvas*>(obj);
156 }
157 }
158 return nullptr;
159}
160
161TH1* DQMHistAnalysisModule::findHist(const std::string& histname, bool was_updated)
162{
163 if (s_histList.find(histname) != s_histList.end()) {
164 if (was_updated && !s_histList[histname].isUpdated()) return nullptr;
165 if (s_histList[histname].getHist()) {
166 return s_histList[histname].getHist();
167 } else {
168 B2ERROR("Histogram " << histname << " in histogram list but nullptr.");
169 }
170 }
171 B2INFO("Histogram " << histname << " not in list.");
172 return nullptr;
173}
174
175TH1* DQMHistAnalysisModule::findHist(const std::string& dirname, const std::string& histname, bool updated)
176{
177 if (dirname.size() > 0) {
178 return findHist(dirname + "/" + histname, updated);
179 }
180 return findHist(histname, updated);
181}
182
183TH1* DQMHistAnalysisModule::scaleReference(ERefScaling scaling, const TH1* hist, TH1* ref)
184{
185 // if hist/ref is nullptr, nothing to do
186 if (!hist || !ref)
187 return ref;
188
189 switch (scaling) {
190 // default: do nothing
191 case ERefScaling::c_RefScaleNone: //do nothing
192 break;
193 case ERefScaling::c_RefScaleEntries: // Integral
194 // only if we have entries in reference
195 if (hist->Integral() != 0 and ref->Integral() != 0) {
196 ref->Scale(hist->Integral() / ref->Integral());
197 }
198 break;
199 case ERefScaling::c_RefScaleMax: // Maximum
200 // only if we have entries in reference
201 if (hist->GetMaximum() != 0 and ref->GetMaximum() != 0) {
202 ref->Scale(hist->GetMaximum() / ref->GetMaximum());
203 }
204 break;
205 }
206 return ref;
207}
208
209TH1* DQMHistAnalysisModule::findRefHist(const std::string& histname, ERefScaling scaling, const TH1* hist)
210{
211 if (s_refList.find(histname) != s_refList.end()) {
212 // get a copy of the reference which we can modify
213 // (it is still owned and managed by the framework)
214 // then do the scaling
215 return scaleReference(scaling, hist, s_refList[histname].getReference());
216 }
217 return nullptr;
218}
219
220TH1* DQMHistAnalysisModule::findRefHist(const std::string& dirname, const std::string& histname, ERefScaling scaling,
221 const TH1* hist)
222{
223 if (dirname.size() > 0) {
224 return findRefHist(dirname + "/" + histname, scaling, hist);
225 }
226 return findRefHist(histname, scaling, hist);
227}
228
229TH1* DQMHistAnalysisModule::findHistInCanvas(const std::string& histo_name, TCanvas** cobj)
230{
231 TCanvas* cnv = nullptr;
232 // try to get canvas from outside
233 if (cobj) cnv = *cobj;
234 // if no canvas search for it
235 if (cnv == nullptr) {
236 // parse the dir+histo name and create the corresponding canvas name
237 auto s = StringSplit(histo_name, '/');
238 if (s.size() != 2) {
239 B2ERROR("findHistInCanvas: histoname not valid (missing dir?), should be 'dirname/histname': " << histo_name);
240 return nullptr;
241 }
242 auto dirname = s.at(0);
243 auto hname = s.at(1);
244 std::string canvas_name = dirname + "/c_" + hname;
245 cnv = findCanvas(canvas_name);
246 // set canvas pointer for outside
247 if (cnv && cobj) *cobj = cnv;
248 }
249
250 // get histogram pointer
251 if (cnv != nullptr) {
252 TIter nextkey(cnv->GetListOfPrimitives());
253 TObject* obj{};
254 while ((obj = dynamic_cast<TObject*>(nextkey()))) {
255 if (obj->IsA()->InheritsFrom("TH1")) {
256 if (obj->GetName() == histo_name)
257 return dynamic_cast<TH1*>(obj);
258 }
259 }
260 }
261 return nullptr;
262}
263
264TH1* DQMHistAnalysisModule::findHistInFile(TFile* file, const std::string& histname)
265{
266 // find histogram by name in file, histname CAN contain directory!
267 // will return nullptr if file is zeroptr, not found or not correct type
268 if (file && file->IsOpen()) {
269 auto obj = file->Get(histname.data());
270 if (obj != nullptr) {
271 // check class type
272 if (obj->IsA()->InheritsFrom("TH1")) {
273 B2DEBUG(20, "Histogram " << histname << " found in file");
274 return dynamic_cast<TH1*>(obj);
275 } else {
276 B2INFO("Found Object " << histname << " in file is not a histogram");
277 }
278 } else {
279 B2INFO("Histogram " << histname << " not found in file");
280 }
281 }
282 return nullptr;
283}
284
286{
287 if (s_monObjList.find(objName) != s_monObjList.end()) {
288 return &s_monObjList[objName];
289 }
290 B2INFO("MonitoringObject " << objName << " not in memfile.");
291 return nullptr;
292}
293
295{
296 double probs[2] = {0.16, 1 - 0.16};
297 double quant[2] = {0, 0};
298 h->GetQuantiles(2, quant, probs);
299 const double sigma68 = (-quant[0] + quant[1]) / 2;
300 return sigma68;
301}
302
303std::vector <std::string> DQMHistAnalysisModule::StringSplit(const std::string& in, const char delim)
304{
305 std::vector <std::string> out;
306 boost::split(out, in, [delim](char c) {return c == delim;});
307 return out;
308}
309
311{
312 TIter nextckey(gROOT->GetListOfCanvases());
313 TObject* cobj = nullptr;
314
315 while ((cobj = dynamic_cast<TObject*>(nextckey()))) {
316 if (cobj->IsA()->InheritsFrom("TCanvas")) {
317 TCanvas* cnv = dynamic_cast<TCanvas*>(cobj);
318 cnv->Clear();
320 }
321 }
322}
323
325{
326 for (auto& it : s_histList) {
327 // attention, we must use reference, otherwise we work on a copy
328 it.second.resetBeforeEvent();
329 }
330 for (auto& it : s_deltaList) {
331 // attention, we must use reference, otherwise we work on a copy
332 it.second.setNotUpdated();
333 }
334
335 s_canvasUpdatedList.clear();
336}
337
339{
340 s_histList.clear();
341}
342
344{
345 s_refList.clear();
346}
347
349{
350 for (auto& d : s_deltaList) {
351 d.second.reset();
352 }
353}
354
355void DQMHistAnalysisModule::UpdateCanvas(std::string name, bool updated)
356{
357 s_canvasUpdatedList[name] = updated;
358}
359
360void DQMHistAnalysisModule::UpdateCanvas(TCanvas* c, bool updated)
361{
362 if (c) UpdateCanvas(c->GetName(), updated);
363}
364
365void DQMHistAnalysisModule::ExtractRunType(std::vector <TH1*>& hs)
366{
367 s_runType = "";
368 for (size_t i = 0; i < hs.size(); i++) {
369 if (hs[i]->GetName() == std::string("DQMInfo/rtype")) {
370 s_runType = hs[i]->GetTitle();
371 return;
372 }
373 }
374 B2ERROR("ExtractRunType: Histogram \"DQMInfo/rtype\" missing");
375}
376
377void DQMHistAnalysisModule::ExtractEvent(std::vector <TH1*>& hs)
378{
380 for (size_t i = 0; i < hs.size(); i++) {
381 if (hs[i]->GetName() == std::string("DAQ/Nevent")) {
382 s_eventProcessed = hs[i]->GetEntries();
383 return;
384 }
385 }
386 B2ERROR("ExtractEvent: Histogram \"DAQ/Nevent\" missing");
387}
388
389int DQMHistAnalysisModule::registerEpicsPV(std::string pvname, std::string keyname)
390{
391 return registerEpicsPVwithPrefix(m_PVPrefix, pvname, keyname);
392}
393
394int DQMHistAnalysisModule::registerExternalEpicsPV(std::string pvname, std::string keyname)
395{
396 return registerEpicsPVwithPrefix(std::string(""), pvname, keyname);
397}
398
399int DQMHistAnalysisModule::registerEpicsPVwithPrefix(std::string prefix, std::string pvname, std::string keyname)
400{
401 if (!m_useEpics) return -1;
402#ifdef _BELLE2_EPICS
403 if (m_epicsNameToChID[pvname] != nullptr) {
404 B2ERROR("Epics PV " << pvname << " already registered!");
405 return -1;
406 }
407 if (keyname != "" && m_epicsNameToChID[keyname] != nullptr) {
408 B2ERROR("Epics PV with key " << keyname << " already registered!");
409 return -1;
410 }
411
412 m_epicsChID.emplace_back();
413 auto ptr = &m_epicsChID.back();
414 if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
415 // the subscribed name includes the prefix, the map below does *not*
416 CheckEpicsError(ca_create_channel((prefix + pvname).data(), NULL, NULL, 10, ptr), "ca_create_channel failure", pvname);
417
418 m_epicsNameToChID[pvname] = *ptr;
419 if (keyname != "") m_epicsNameToChID[keyname] = *ptr;
420 return m_epicsChID.size() - 1; // return index to last added item
421#else
422 return -1;
423#endif
424}
425
426void DQMHistAnalysisModule::setEpicsPV(std::string keyname, double value)
427{
428 if (!m_useEpics || m_epicsReadOnly) return;
429#ifdef _BELLE2_EPICS
430 if (m_epicsNameToChID[keyname] == nullptr) {
431 B2ERROR("Epics PV " << keyname << " not registered!");
432 return;
433 }
434 CheckEpicsError(ca_put(DBR_DOUBLE, m_epicsNameToChID[keyname], (void*)&value), "ca_set failure", keyname);
435#endif
436}
437
438void DQMHistAnalysisModule::setEpicsPV(std::string keyname, int value)
439{
440 if (!m_useEpics || m_epicsReadOnly) return;
441#ifdef _BELLE2_EPICS
442 if (m_epicsNameToChID[keyname] == nullptr) {
443 B2ERROR("Epics PV " << keyname << " not registered!");
444 return;
445 }
446 CheckEpicsError(ca_put(DBR_SHORT, m_epicsNameToChID[keyname], (void*)&value), "ca_set failure", keyname);
447#endif
448}
449
450void DQMHistAnalysisModule::setEpicsStringPV(std::string keyname, std::string value)
451{
452 if (!m_useEpics || m_epicsReadOnly) return;
453#ifdef _BELLE2_EPICS
454 if (m_epicsNameToChID[keyname] == nullptr) {
455 B2ERROR("Epics PV " << keyname << " not registered!");
456 return;
457 }
458 if (value.length() > 40) {
459 B2ERROR("Epics string PV " << keyname << " too long (>40 characters)!");
460 return;
461 }
462 char text[40];
463 strcpy(text, value.c_str());
464 CheckEpicsError(ca_put(DBR_STRING, m_epicsNameToChID[keyname], text), "ca_set failure", keyname);
465#endif
466}
467
468void DQMHistAnalysisModule::setEpicsPV(int index, double value)
469{
470 if (!m_useEpics || m_epicsReadOnly) return;
471#ifdef _BELLE2_EPICS
472 if (index < 0 || index >= (int)m_epicsChID.size()) {
473 B2ERROR("Epics PV with " << index << " not registered!");
474 return;
475 }
476 CheckEpicsError(ca_put(DBR_DOUBLE, m_epicsChID[index], (void*)&value), "ca_set failure", m_epicsChID[index]);
477#endif
478}
479
480void DQMHistAnalysisModule::setEpicsPV(int index, int value)
481{
482 if (!m_useEpics || m_epicsReadOnly) return;
483#ifdef _BELLE2_EPICS
484 if (index < 0 || index >= (int)m_epicsChID.size()) {
485 B2ERROR("Epics PV with " << index << " not registered!");
486 return;
487 }
488 CheckEpicsError(ca_put(DBR_SHORT, m_epicsChID[index], (void*)&value), "ca_set failure", m_epicsChID[index]);
489#endif
490}
491
492void DQMHistAnalysisModule::setEpicsStringPV(int index, std::string value)
493{
494 if (!m_useEpics || m_epicsReadOnly) return;
495#ifdef _BELLE2_EPICS
496 if (index < 0 || index >= (int)m_epicsChID.size()) {
497 B2ERROR("Epics PV with " << index << " not registered!");
498 return;
499 }
500 char text[41];
501 strncpy(text, value.c_str(), 40);
502 text[40] = 0;
503 CheckEpicsError(ca_put(DBR_STRING, m_epicsChID[index], text), "ca_set failure", m_epicsChID[index]);
504#endif
505}
506
507double DQMHistAnalysisModule::getEpicsPV(std::string keyname)
508{
509 double value{NAN};
510 if (!m_useEpics) return value;
511#ifdef _BELLE2_EPICS
512 if (m_epicsNameToChID[keyname] == nullptr) {
513 B2ERROR("Epics PV " << keyname << " not registered!");
514 return value;
515 }
516 // From EPICS doc. When ca_get or ca_array_get are invoked the returned channel value can't be assumed to be stable
517 // in the application supplied buffer until after ECA_NORMAL is returned from ca_pend_io. If a connection is lost
518 // outstanding get requests are not automatically reissued following reconnect.
519 auto r = ca_get(DBR_DOUBLE, m_epicsNameToChID[keyname], (void*)&value);
520 if (r == ECA_NORMAL) r = ca_pend_io(5.0); // this is needed!
521 if (r == ECA_NORMAL) {
522 return value;
523 } else {
524 CheckEpicsError(r, "Read PV failed in ca_get or ca_pend_io failure", keyname);
525 }
526#endif
527 return NAN;
528}
529
531{
532 double value{NAN};
533 if (!m_useEpics) return value;
534#ifdef _BELLE2_EPICS
535 if (index < 0 || index >= (int)m_epicsChID.size()) {
536 B2ERROR("Epics PV with " << index << " not registered!");
537 return value;
538 }
539 // From EPICS doc. When ca_get or ca_array_get are invoked the returned channel value can't be assumed to be stable
540 // in the application supplied buffer until after ECA_NORMAL is returned from ca_pend_io. If a connection is lost
541 // outstanding get requests are not automatically reissued following reconnect.
542 auto r = ca_get(DBR_DOUBLE, m_epicsChID[index], (void*)&value);
543 if (r == ECA_NORMAL) r = ca_pend_io(5.0); // this is needed!
544 if (r == ECA_NORMAL) {
545 return value;
546 } else {
547 CheckEpicsError(r, "Read PV failed in ca_get or ca_pend_io failure", m_epicsChID[index]);
548 }
549#endif
550 return NAN;
551}
552
553std::string DQMHistAnalysisModule::getEpicsStringPV(std::string keyname, bool& status)
554{
555 status = false;
556 char value[40] = "";
557 if (!m_useEpics) return std::string(value);
558#ifdef _BELLE2_EPICS
559 if (m_epicsNameToChID[keyname] == nullptr) {
560 B2ERROR("Epics PV " << keyname << " not registered!");
561 return std::string(value);
562 }
563 // From EPICS doc. When ca_get or ca_array_get are invoked the returned channel value can't be assumed to be stable
564 // in the application supplied buffer until after ECA_NORMAL is returned from ca_pend_io. If a connection is lost
565 // outstanding get requests are not automatically reissued following reconnect.
566 auto r = ca_get(DBR_STRING, m_epicsNameToChID[keyname], value);
567 if (r == ECA_NORMAL) r = ca_pend_io(5.0); // this is needed!
568 if (r == ECA_NORMAL) {
569 status = true;
570 return std::string(value);
571 } else {
572 CheckEpicsError(r, "Read PV (string) failed in ca_get or ca_pend_io failure", keyname);
573 }
574#endif
575 return std::string(value);
576}
577
578std::string DQMHistAnalysisModule::getEpicsStringPV(int index, bool& status)
579{
580 status = false;
581 char value[40] = "";
582 if (!m_useEpics) return std::string(value);
583#ifdef _BELLE2_EPICS
584 if (index < 0 || index >= (int)m_epicsChID.size()) {
585 B2ERROR("Epics PV with " << index << " not registered!");
586 return std::string(value);
587 }
588 // From EPICS doc. When ca_get or ca_array_get are invoked the returned channel value can't be assumed to be stable
589 // in the application supplied buffer until after ECA_NORMAL is returned from ca_pend_io. If a connection is lost
590 // outstanding get requests are not automatically reissued following reconnect.
591 auto r = ca_get(DBR_DOUBLE, m_epicsChID[index], value);
592 if (r == ECA_NORMAL) r = ca_pend_io(5.0); // this is needed!
593 if (r == ECA_NORMAL) {
594 status = true;
595 return std::string(value);
596 } else {
597 CheckEpicsError(r, "Read PV (string) failed in ca_get or ca_pend_io failure", m_epicsChID[index]);
598 }
599#endif
600 return std::string(value);
601}
602
604{
605#ifdef _BELLE2_EPICS
606 if (m_useEpics) {
607 if (m_epicsNameToChID[keyname] != nullptr) {
608 return m_epicsNameToChID[keyname];
609 } else {
610 B2ERROR("Epics PV " << keyname << " not registered!");
611 }
612 }
613#endif
614 return nullptr;
615}
616
618{
619#ifdef _BELLE2_EPICS
620 if (m_useEpics) {
621 if (index >= 0 && index < (int)m_epicsChID.size()) {
622 return m_epicsChID[index];
623 } else {
624 B2ERROR("Epics PV with " << index << " not registered!");
625 }
626 }
627#endif
628 return nullptr;
629}
630
632{
633 int state = ECA_NORMAL;
634 if (!m_useEpics) return state;
635#ifdef _BELLE2_EPICS
636 if (wait > 0.) {
637 state = ca_pend_io(wait);
638 SEVCHK(state, "ca_pend_io failure");
639 }
640#endif
641 return state;
642}
643
645{
646 // this should be called in terminate function of analysis modules
647#ifdef _BELLE2_EPICS
648 if (getUseEpics()) {
649 for (auto& it : m_epicsChID) CheckEpicsError(ca_clear_channel(it), "ca_clear_channel failure", it);
650 updateEpicsPVs(5.0);
651 // Make sure we clean up both afterwards!
652 m_epicsChID.clear();
653 m_epicsNameToChID.clear();
654 }
655#endif
656}
657
658bool DQMHistAnalysisModule::requestLimitsFromEpicsPVs(std::string name, double& lowerAlarm, double& lowerWarn, double& upperWarn,
659 double& upperAlarm)
660{
661 return requestLimitsFromEpicsPVs(getEpicsPVChID(name), lowerAlarm, lowerWarn, upperWarn, upperAlarm);
662}
663
664bool DQMHistAnalysisModule::requestLimitsFromEpicsPVs(int index, double& lowerAlarm, double& lowerWarn, double& upperWarn,
665 double& upperAlarm)
666{
667 return requestLimitsFromEpicsPVs(getEpicsPVChID(index), lowerAlarm, lowerWarn, upperWarn, upperAlarm);
668}
669
670bool DQMHistAnalysisModule::requestLimitsFromEpicsPVs(chid pv, double& lowerAlarm, double& lowerWarn, double& upperWarn,
671 double& upperAlarm)
672{
673 // get warn and error limit only if pv exists
674 // overwrite only if limit is defined (not NaN)
675 // user should initialize with NaN before calling, unless
676 // some "default" values should be set otherwise
677 if (pv != nullptr) {
678 struct dbr_ctrl_double tPvData;
679 // From EPICS doc. When ca_get or ca_array_get are invoked the returned channel value can't be assumed to be stable
680 // in the application supplied buffer until after ECA_NORMAL is returned from ca_pend_io. If a connection is lost
681 // outstanding get requests are not automatically reissued following reconnect.
682 auto r = ca_get(DBR_CTRL_DOUBLE, pv, &tPvData);
683 if (r == ECA_NORMAL) r = ca_pend_io(5.0); // this is needed!
684 if (r == ECA_NORMAL) {
685 if (!std::isnan(tPvData.lower_alarm_limit)) {
686 lowerAlarm = tPvData.lower_alarm_limit;
687 }
688 if (!std::isnan(tPvData.lower_warning_limit)) {
689 lowerWarn = tPvData.lower_warning_limit;
690 }
691 if (!std::isnan(tPvData.upper_warning_limit)) {
692 upperWarn = tPvData.upper_warning_limit;
693 }
694 if (!std::isnan(tPvData.upper_alarm_limit)) {
695 upperAlarm = tPvData.upper_alarm_limit;
696 }
697 return true;
698 } else {
699 CheckEpicsError(r, "Reading PV Limits failed in ca_get or ca_pend_io failure", pv);
700 }
701 }
702 return false;
703}
704
705DQMHistAnalysisModule::EStatus DQMHistAnalysisModule::makeStatus(bool enough, bool warn_flag, bool error_flag)
706{
707 // white color is the default, if no colorize
708 if (!enough) {
709 return (c_StatusTooFew);
710 } else {
711 if (error_flag) {
712 return (c_StatusError);
713 } else if (warn_flag) {
714 return (c_StatusWarning);
715 } else {
716 return (c_StatusGood);
717 }
718 }
719
720 return (c_StatusDefault); // default, but should not be reached
721}
722
724{
725 // white color is the default, if no colorize
727 switch (stat) {
728 case c_StatusTooFew:
729 color = c_ColorTooFew; // Magenta or Gray
730 break;
731 case c_StatusDefault:
732 color = c_ColorDefault; // default no colors
733 break;
734 case c_StatusGood:
735 color = c_ColorGood; // Good
736 break;
737 case c_StatusWarning:
738 color = c_ColorWarning; // Warning
739 break;
740 case c_StatusError:
741 color = c_ColorError; // Severe
742 break;
743 default:
744 color = c_ColorDefault; // default no colors
745 break;
746 }
747 return color;
748}
749
751{
752 if (!canvas) return;
753 auto color = getStatusColor(stat);
754
755 canvas->Pad()->SetFillColor(color);
756
757 canvas->Pad()->SetFrameFillColor(10); // White (kWhite is not used since it results in transparent!)
758 canvas->Pad()->SetFrameFillStyle(1001);// White
759 canvas->Pad()->Modified();
760 canvas->Pad()->Update();
761}
762
764{
765 B2INFO("Check PV Connections");
766
767 for (auto& it : m_epicsChID) {
768 printPVStatus(it);
769 }
770 B2INFO("Check PVs done");
771}
772
773void DQMHistAnalysisModule::printPVStatus(chid pv, bool onlyError)
774{
775 if (pv == nullptr) {
776 B2WARNING("PV chid was nullptr");
777 return;
778 }
779 auto state = ca_state(pv);
780 switch (state) {
781 case cs_never_conn: /* valid chid, server not found or unavailable */
782 B2WARNING("Channel never connected " << ca_name(pv));
783 break;
784 case cs_prev_conn: /* valid chid, previously connected to server */
785 B2WARNING("Channel was connected, but now is not " << ca_name(pv));
786 break;
787 case cs_closed: /* channel deleted by user */
788 B2WARNING("Channel deleted already " << ca_name(pv));
789 break;
790 case cs_conn: /* valid chid, connected to server */
791 if (!onlyError) B2INFO("Channel connected and OK " << ca_name(pv));
792 break;
793 default:
794 B2WARNING("Undefined status for channel " << ca_name(pv));
795 break;
796 }
797}
798
799void DQMHistAnalysisModule::CheckEpicsError(int state, const std::string& message, const std::string& name)
800{
801 if (state != ECA_NORMAL) {
802 B2WARNING(message << ": " << name);
803 printPVStatus(m_epicsNameToChID[name], false);
804 }
805}
806
807void DQMHistAnalysisModule::CheckEpicsError(int state, const std::string& message, chid id = nullptr)
808{
809 if (state != ECA_NORMAL) {
810 std::string name;
811 if (id) name = ca_name(id);
812 B2WARNING(message << ": " << name);
813 printPVStatus(id, false);
814 }
815}
816
static MonObjList s_monObjList
The list of MonitoringObjects.
TCanvas * findCanvas(TString cname)
Find canvas by name.
void printPVStatus(chid pv, bool onlyError=true)
check the status of a PVs and report if disconnected or not found
static TH1 * findRefHist(const std::string &histname, ERefScaling scaling=ERefScaling::c_RefScaleNone, const TH1 *hist=nullptr)
Get referencehistogram from list (no other search).
bool hasDeltaPar(const std::string &dirname, const std::string &histname)
Check if Delta histogram parameters exist for histogram.
std::map< std::string, HistObject > HistList
The type of list of histograms.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
std::map< std::string, MonitoringObject > MonObjList
The type of list of MonitoringObjects.
static TH1 * scaleReference(ERefScaling scaling, const TH1 *hist, TH1 *ref)
Using the original and reference, create scaled version.
void addDeltaPar(const std::string &dirname, const std::string &histname, HistDelta::EDeltaType t, int p, unsigned int a=1)
Add Delta histogram parameters.
static TH1 * findHistInFile(TFile *file, const std::string &histname)
Find histogram in specific TFile (e.g.
EStatusColor getStatusColor(EStatus status)
Return color for canvas state.
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).
static MonitoringObject * findMonitoringObject(const std::string &objName)
Find MonitoringObject.
void clearlist(void)
Clear all static global lists.
double getSigma68(TH1 *h) const
Helper function to compute half of the central interval covering 68% of a distribution.
EStatusColor
Status colors of histogram/canvas (corresponding to status)
@ c_ColorWarning
Analysis result: Warning, there may be minor issues.
@ c_ColorError
Analysis result: Severe issue found.
@ c_ColorTooFew
Not enough entries/event to judge.
@ c_ColorGood
Analysis result: Good.
@ c_ColorDefault
default for non-coloring
double getEpicsPV(std::string keyname)
Read value from a EPICS PV.
static int s_eventProcessed
Number of Events processed to fill histograms.
std::map< std::string, bool > CanvasUpdatedList
The type of list of canvas updated status.
std::string getEpicsStringPV(std::string keyname, bool &status)
Read value from a EPICS PV.
int registerExternalEpicsPV(std::string pvname, std::string keyname="")
Register a PV with its name and a key name.
static HistList s_histList
The list of Histograms.
static RefList s_refList
The list of references.
static std::string s_runType
The Run type.
static void clearHistList(void)
Clears the list of histograms.
TH1 * getDelta(const std::string &fullname, int n=0, bool onlyIfUpdated=true)
Get Delta histogram.
std::vector< std::string > StringSplit(const std::string &s, const char delim)
Helper function for string token split.
static void clearRefList(void)
Clears the list of ref histograms.
std::map< std::string, HistDelta > DeltaList
The type of list of delta settings and histograms.
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
static DeltaList s_deltaList
The list of Delta Histograms and settings.
DQMHistAnalysisModule()
Constructor / Destructor.
void checkPVStatus(void)
Check the status of all PVs and report if disconnected or not found.
std::map< std::string, RefHistObject > RefList
The type of list of references.
static bool m_epicsReadOnly
Flag if to use EPICS in ReadOnly mode (for reading limits) do not set by yourself,...
EStatus
Status flag of histogram/canvas.
@ c_StatusDefault
default for non-coloring
@ c_StatusTooFew
Not enough entries/event to judge.
@ c_StatusError
Analysis result: Severe issue found.
@ c_StatusWarning
Analysis result: Warning, there may be minor issues.
@ c_StatusGood
Analysis result: Good.
static bool addHist(const std::string &dirname, const std::string &histname, TH1 *h)
Add histogram.
bool getUseEpics(void)
Getter for EPICS usage.
void ExtractRunType(std::vector< TH1 * > &hs)
Extract Run Type from histogram title, called from input module.
void CheckEpicsError(int state, const std::string &message, const std::string &name)
check the return status and check PV in case of error
TH1 * findHistInCanvas(const std::string &hname, TCanvas **canvas=nullptr)
Find histogram in corresponding canvas.
static std::string m_PVPrefix
The Prefix for EPICS PVs.
void cleanupEpicsPVs(void)
Unsubscribe from EPICS PVs on terminate.
void clearCanvases(void)
Clear content of all Canvases.
ERefScaling
Reference plot scaling type.
@ c_RefScaleEntries
to number of entries (integral)
@ c_RefScaleMax
to maximum (bin entry)
EStatus makeStatus(bool enough, bool warn_flag, bool error_flag)
Helper function to judge the status for coloring and EPICS.
static bool m_useEpics
Flag if to use EPICS do not set by yourself, use EpicsEnable module to set.
static void initHistListBeforeEvent(void)
Reset the list of histograms.
void ExtractEvent(std::vector< TH1 * > &hs)
Extract event processed from daq histogram, called from input module.
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 resetDeltaList(void)
Reset Delta.
chid getEpicsPVChID(std::string keyname)
Get EPICS PV Channel Id.
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
int registerEpicsPVwithPrefix(std::string prefix, std::string pvname, std::string keyname="")
Register a PV with its name and a key name.
void setEpicsStringPV(std::string keyname, std::string value)
Write string to a EPICS PV.
int updateEpicsPVs(float timeout)
Update all EPICS PV (flush to network)
static CanvasUpdatedList s_canvasUpdatedList
The list of canvas updated status.
EDeltaType
enum definition for delta algo Disabled: nothing Entries: use nr histogram entries Underflow: use ent...
Definition: HistDelta.h:32
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
MonitoringObject is a basic object to hold data for the run-dependency monitoring Run summary TCanvas...
#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.