Belle II Software  release-08-01-10
DQMHistAnalysis.h
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 : DQMHistAnalysis.h
10 // Description : Baseclass for histogram analysis module for DQM
11 //-
12 
13 #pragma once
14 
15 #include <framework/core/Module.h>
16 #include <dqm/core/MonitoringObject.h>
17 #include <dqm/analysis/HistObject.h>
18 #include <dqm/analysis/HistDelta.h>
19 #include <TFile.h>
20 #include <TH1.h>
21 
22 #include <string>
23 #include <map>
24 
25 #ifdef _BELLE2_EPICS
26 // EPICS
27 #include "cadef.h"
28 #endif
29 
30 namespace Belle2 {
39  class DQMHistAnalysisModule : public Module {
40 
41  public:
45  enum EStatus {
50  c_StatusError = 4
51  };
52 
56  enum EStatusColor {
57  c_ColorTooFew = kGray,
58  c_ColorDefault = kWhite,
59  c_ColorGood = kGreen,
60  c_ColorWarning = kYellow,
61  c_ColorError = kRed
62  };
63 
67  typedef std::map<std::string, HistObject> HistList;
71  typedef std::map<std::string, MonitoringObject*> MonObjList;
72 
76  typedef std::map<std::string, HistDelta*> DeltaList;
77 
81  typedef std::map<std::string, bool> CanvasUpdatedList;
82 
83  private:
92 
97 
102 
110  inline static int s_eventProcessed = 0;
111 
115  inline static std::string s_runType = "";
116 
121  static bool m_useEpics;
122 
127  static bool m_epicsReadOnly;
128 
132  static std::string m_PVPrefix;
133 
134 
135 #ifdef _BELLE2_EPICS
137  std::vector <chid> m_epicsChID;
139  std::map <std::string, chid> m_epicsNameToChID;
140 #endif
141 
142  public:
147  static /*const*/ HistList& getHistList() { return s_histList;};
148 
153  static const MonObjList& getMonObjList() { return s_monObjList;};
154 
159  static const DeltaList& getDeltaList() { return s_deltaList;};
160 
166 
171  static const std::string& getRunType(void) { return s_runType;};
172 
177  static int getEventProcessed(void) { return s_eventProcessed;};
178 
183  void setRunType(std::string& t) {s_runType = t;};
184 
190 
196  TCanvas* findCanvas(TString cname);
197 
204  static TH1* findHist(const std::string& histname, bool onlyIfUpdated = false);
205 
213  static TH1* findHist(const std::string& dirname,
214  const std::string& histname, bool onlyIfUpdated = false);
215 
222  static TH1* findHistInFile(TFile* file, const std::string& histname);
223 
229  TH1* findHistInCanvas(const std::string& hname);
230 
236  static MonitoringObject* findMonitoringObject(const std::string& objName);
237 
244  double getSigma68(TH1* h) const;
245 
246  public:
254  static bool addHist(const std::string& dirname,
255  const std::string& histname, TH1* h);
256 
262  static MonitoringObject* getMonitoringObject(const std::string& histname);
263 
267  void clearCanvases(void);
268 
272  static void initHistListBeforeEvent(void);
273 
277  static void clearHistList(void);
278 
286  TH1* getDelta(const std::string& fullname, int n = 0, bool onlyIfUpdated = true);
287 
296  TH1* getDelta(const std::string& dirname, const std::string& histname, int n = 0, bool onlyIfUpdated = true);
297 
306  void addDeltaPar(const std::string& dirname, const std::string& histname, HistDelta::EDeltaType t, int p, unsigned int a = 1);
307 
314  bool hasDeltaPar(const std::string& dirname, const std::string& histname);
315 
321  void UpdateCanvas(std::string name, bool updated = true);
322 
328  void UpdateCanvas(TCanvas* canvas, bool updated = true);
329 
333  void ExtractRunType(std::vector <TH1*>& hs);
334 
338  void ExtractEvent(std::vector <TH1*>& hs);
339 
341 
353  int registerEpicsPV(std::string pvname, std::string keyname = "", bool update_pvs = true);
354 
360  void setEpicsPV(std::string keyname, double value);
361 
367  void setEpicsPV(std::string keyname, int value);
368 
374  void setEpicsStringPV(std::string keyname, std::string value);
375 
381  void setEpicsPV(int index, double value);
382 
388  void setEpicsPV(int index, int value);
389 
395  void setEpicsStringPV(int index, std::string value);
396 
402  double getEpicsPV(std::string keyname);
403 
409  double getEpicsPV(int index);
410 
417  std::string getEpicsStringPV(std::string keyname, bool& status);
418 
425  std::string getEpicsStringPV(int index, bool& status);
426 
432  int updateEpicsPVs(float timeout);
433 
439  chid getEpicsPVChID(std::string keyname);
440 
446  chid getEpicsPVChID(int index);
447 
457  bool requestLimitsFromEpicsPVs(chid id, double& lowerAlarm, double& lowerWarn, double& upperWarn, double& upperAlarm);
458 
468  bool requestLimitsFromEpicsPVs(std::string keyname, double& lowerAlarm, double& lowerWarn, double& upperWarn, double& upperAlarm);
469 
479  bool requestLimitsFromEpicsPVs(int index, double& lowerAlarm, double& lowerWarn, double& upperWarn, double& upperAlarm);
480 
485  void setUseEpics(bool flag) {m_useEpics = flag;};
486 
491  void setUseEpicsReadOnly(bool flag) {m_epicsReadOnly = flag;};
492 
497  bool getUseEpics(void) {return m_useEpics;};
498 
503  bool getUseEpicsReadOnly(void) {return m_epicsReadOnly;};
504 
508  void cleanupEpicsPVs(void);
509 
514  std::string& getPVPrefix(void) {return m_PVPrefix;};
515 
520  void setPVPrefix(std::string& prefix) { m_PVPrefix = prefix;};
521 
529  EStatus makeStatus(bool enough, bool warn_flag, bool error_flag);
530 
536  void colorizeCanvas(TCanvas* canvas, EStatus status);
537 
544 
548  void checkPVStatus(void);
549 
555  void printPVStatus(chid pv, bool onlyError = true);
556 
557  // Public functions
558  public:
559 
562  virtual ~DQMHistAnalysisModule() {};
563 
570  std::vector <std::string> StringSplit(const std::string& s, const char delim);
571 
572  // Data members
573  private:
574  };
576 } // end namespace Belle2
577 
The base class for the histogram analysis module.
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
bool hasDeltaPar(const std::string &dirname, const std::string &histname)
Check if Delta histogram parameters exist for histogram.
void setPVPrefix(std::string &prefix)
set global Prefix for EPICS PVs
std::map< std::string, HistObject > HistList
The type of list of histograms.
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.
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.
double getSigma68(TH1 *h) const
Helper function to compute half of the central interval covering 68% of a distribution.
static const DeltaList & getDeltaList()
Get the list of the delta histograms.
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.
static HistList s_histList
The list of Histograms.
bool getUseEpicsReadOnly(void)
Getter EPICS flag in read only mode.
static std::string s_runType
The Run type.
void setUseEpics(bool flag)
Setter for EPICS usage.
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 const CanvasUpdatedList & getCanvasUpdatedList()
Get the list of the canvas update status.
void setRunType(std::string &t)
Set the Run Type.
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
std::map< std::string, MonitoringObject * > MonObjList
The type of list of MonitoringObjects.
void setUseEpicsReadOnly(bool flag)
Setter EPICS flag in read only mode.
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.
void setEventProcessed(int e)
Set the number of processed events.
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.
static const std::string & getRunType(void)
Get the Run Type.
static const MonObjList & getMonObjList()
Get the list of MonitoringObjects.
std::string & getPVPrefix(void)
get global Prefix for EPICS PVs
static std::string m_PVPrefix
The Prefix for EPICS PVs.
TH1 * findHistInCanvas(const std::string &hname)
Find histogram in corresponding canvas.
void cleanupEpicsPVs(void)
Unsubscribe from EPICS PVs on terminate.
void clearCanvases(void)
Clear content of all Canvases.
static int getEventProcessed(void)
Get the number of processed events.
static HistList & getHistList()
Get the list of the histograms.
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.
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)
std::map< std::string, HistDelta * > DeltaList
The type of list of delta settings and histograms.
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.
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
MonitoringObject is a basic object to hold data for the run-dependency monitoring Run summary TCanvas...
Abstract base class for different kinds of events.