Belle II Software  release-08-00-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  typedef std::map<std::string, HistObject> HistList;
49  typedef std::map<std::string, MonitoringObject*> MonObjList;
50 
54  typedef std::map<std::string, HistDelta*> DeltaList;
55 
59  typedef std::map<std::string, bool> CanvasUpdatedList;
60 
61  private:
70 
75 
80 
88  inline static int s_eventProcessed = 0;
89 
93  inline static std::string s_runType = "";
94 
99  static bool m_useEpics;
100 
105  static bool m_epicsReadOnly;
106 
110  static std::string m_PVPrefix;
111 
112 
113 #ifdef _BELLE2_EPICS
115  std::vector <chid> m_epicsChID;
117  std::map <std::string, chid> m_epicsNameToChID;
118 #endif
119 
120  public:
125  static /*const*/ HistList& getHistList() { return s_histList;};
126 
131  static const MonObjList& getMonObjList() { return s_monObjList;};
132 
137  static const DeltaList& getDeltaList() { return s_deltaList;};
138 
144 
149  static const std::string& getRunType(void) { return s_runType;};
150 
155  static int getEventProcessed(void) { return s_eventProcessed;};
156 
161  void setRunType(std::string& t) {s_runType = t;};
162 
168 
174  TCanvas* findCanvas(TString cname);
175 
182  static TH1* findHist(const std::string& histname, bool onlyIfUpdated = false);
183 
191  static TH1* findHist(const std::string& dirname,
192  const std::string& histname, bool onlyIfUpdated = false);
193 
200  static TH1* findHistInFile(TFile* file, const std::string& histname);
201 
207  TH1* findHistInCanvas(const std::string& hname);
208 
214  static MonitoringObject* findMonitoringObject(const std::string& objName);
215 
222  double getSigma68(TH1* h) const;
223 
224  public:
232  static bool addHist(const std::string& dirname,
233  const std::string& histname, TH1* h);
234 
240  static MonitoringObject* getMonitoringObject(const std::string& histname);
241 
245  static void initHistListBeforeEvent(void);
246 
250  static void clearHistList(void);
251 
259  TH1* getDelta(const std::string& fullname, int n = 0, bool onlyIfUpdated = true);
260 
269  TH1* getDelta(const std::string& dirname, const std::string& histname, int n = 0, bool onlyIfUpdated = true);
270 
279  void addDeltaPar(const std::string& dirname, const std::string& histname, HistDelta::EDeltaType t, int p, unsigned int a = 1);
280 
287  bool hasDeltaPar(const std::string& dirname, const std::string& histname);
288 
294  void UpdateCanvas(std::string name, bool updated = true);
295 
301  void UpdateCanvas(TCanvas* canvas, bool updated = true);
302 
306  void ExtractRunType(std::vector <TH1*>& hs);
307 
311  void ExtractEvent(std::vector <TH1*>& hs);
312 
314 
326  int registerEpicsPV(std::string pvname, std::string keyname = "", bool update_pvs = true);
327 
333  void setEpicsPV(std::string keyname, double value);
334 
340  void setEpicsPV(std::string keyname, int value);
341 
347  void setEpicsStringPV(std::string keyname, std::string value);
348 
354  void setEpicsPV(int index, double value);
355 
361  void setEpicsPV(int index, int value);
362 
368  void setEpicsStringPV(int index, std::string value);
369 
375  double getEpicsPV(std::string keyname);
376 
382  double getEpicsPV(int index);
383 
390  std::string getEpicsStringPV(std::string keyname, bool& status);
391 
398  std::string getEpicsStringPV(int index, bool& status);
399 
404  void updateEpicsPVs(float timeout);
405 
411  chid getEpicsPVChID(std::string keyname);
412 
418  chid getEpicsPVChID(int index);
419 
429  bool requestLimitsFromEpicsPVs(chid id, double& lowerAlarm, double& lowerWarn, double& upperWarn, double& upperAlarm);
430 
440  bool requestLimitsFromEpicsPVs(std::string keyname, double& lowerAlarm, double& lowerWarn, double& upperWarn, double& upperAlarm);
441 
451  bool requestLimitsFromEpicsPVs(int index, double& lowerAlarm, double& lowerWarn, double& upperWarn, double& upperAlarm);
452 
457  void setUseEpics(bool flag) {m_useEpics = flag;};
458 
463  void setUseEpicsReadOnly(bool flag) {m_epicsReadOnly = flag;};
464 
469  bool getUseEpics(void) {return m_useEpics;};
470 
475  bool getUseEpicsReadOnly(void) {return m_epicsReadOnly;};
476 
480  void cleanupEpicsPVs(void);
481 
486  std::string& getPVPrefix(void) {return m_PVPrefix;};
487 
492  void setPVPrefix(std::string& prefix) { m_PVPrefix = prefix;};
493 
494  // Public functions
495  public:
496 
499  virtual ~DQMHistAnalysisModule() {};
500 
507  std::vector <std::string> StringSplit(const std::string& s, const char delim);
508 
509  // Data members
510  private:
511  };
513 } // end namespace Belle2
514 
The base class for the histogram analysis module.
static MonObjList s_monObjList
The list of MonitoringObjects.
TCanvas * findCanvas(TString cname)
Find canvas by name.
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.
void updateEpicsPVs(float timeout)
Update all EPICS PV (flush to network)
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.
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.
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 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,...
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.
static int getEventProcessed(void)
Get the number of processed events.
static HistList & getHistList()
Get the list of the histograms.
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.
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.