Belle II Software  release-08-00-10
DQMHistAnalysisECLOutOfTimeDigits.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 //THIS MODULE
10 #include <dqm/analysis/modules/DQMHistAnalysisECLOutOfTimeDigits.h>
11 
12 //ROOT
13 #include <TF1.h>
14 #include <TH1F.h>
15 
16 using namespace Belle2;
17 
18 REG_MODULE(DQMHistAnalysisECLOutOfTimeDigits);
19 
22 {
23  B2DEBUG(20, "DQMHistAnalysisECLOutOfTimeDigits: Constructor done.");
24  setDescription("Module to collect and process 'out of time' ECL digits");
25  addParam("pvPrefix", m_pvPrefix, "Prefix to use for PVs registered by this module",
26  std::string("ECL:out_of_time_digits:"));
27  addParam("onlyIfUpdated", m_onlyIfUpdated, "If true (default), update EPICS PVs only if there histograms were updated.",
28  true);
29 }
30 
31 
33 {
34  // Register EPICS PVs
35  for (auto& event_type : {"rand", "dphy", "physics"}) {
36  for (auto& ecl_part : {"All", "FWDEndcap", "Barrel", "BWDEndcap"}) {
37  std::string pv_name = event_type + std::string(":") + ecl_part;
38  registerEpicsPV(m_pvPrefix + pv_name, pv_name);
39  }
40  }
41  updateEpicsPVs(5.0);
42 
44 
45  B2DEBUG(20, "DQMHistAnalysisECLOutOfTimeDigits: initialized.");
46 }
47 
49 {
50  //== Get DQM info
51  for (auto& event_type : {"rand", "dphy", "physics"}) {
52  for (auto& ecl_part : {"All", "FWDEndcap", "Barrel", "BWDEndcap"}) {
53  std::string pv_name = event_type + std::string(":") + ecl_part;
54  std::string var_name = pv_name;
55  std::replace(var_name.begin(), var_name.end(), ':', '_');
56 
57  m_out_of_time_digits[pv_name] = 0;
58 
59  std::string hist_name = "ECL/out_of_time_" + var_name;
60  auto hist = (TH1F*)findHist(hist_name, m_onlyIfUpdated);
61 
62  if (!hist) continue;
63 
64  m_out_of_time_digits[pv_name] = hist->GetMean();
65 
66  //== Set EPICS PVs
67 
68  double selected_value = m_out_of_time_digits[pv_name];
69  setEpicsPV(pv_name, selected_value);
70  }
71  }
72  updateEpicsPVs(5.0);
73 }
74 
76 {
77  B2DEBUG(20, "DQMHistAnalysisECLOutOfTimeDigits: endRun called");
78 
79  auto main_hist = (TH1F*)findHist("ECL/out_of_time_physics_All");
80 
81  if (main_hist == nullptr) {
82  m_monObj->setVariable("comment_out_of_time_digits", "No ECL out-of-time ECLCalDigits histograms available");
83  B2INFO("Histogram named ECL/out_of_time_physics_All is not found.");
84  return;
85  }
86 
87  TF1 gaus("fit_func", "gaus");
88 
89  for (auto& event_type : {"rand", "dphy", "physics"}) {
90  for (auto& ecl_part : {"All", "FWDEndcap", "Barrel", "BWDEndcap"}) {
91  std::string pv_name = event_type + std::string(":") + ecl_part;
92  std::string hist_name = "ECL/out_of_time_" + pv_name;
93  std::string var_name = "out_of_time_digits_" + pv_name;
94 
95  std::replace(hist_name.begin(), hist_name.end(), ':', '_');
96  std::replace(var_name.begin(), var_name.end(), ':', '_');
97 
98  // If enough statistics, obtain more detailed information for MiraBelle
99  auto hist = (TH1F*)findHist(hist_name);
100  if (hist && hist->GetEntries() > 1000) {
101  // Fit the histogram to get the peak of a distribution
102  hist->Fit(&gaus);
103  m_monObj->setVariable(var_name, gaus.GetParameter(1));
104  m_monObj->setVariable(var_name + "_stddev", gaus.GetParameter(2));
105  } else {
106  // Use simple mean from the histogram
107  m_monObj->setVariable(var_name, m_out_of_time_digits[pv_name]);
108  m_monObj->setVariable(var_name + "_stddev", 0);
109  }
110  }
111  }
112 }
113 
114 
void initialize() override final
Initialize the module.
std::map< std::string, double > m_out_of_time_digits
Out-of-time ECLCalDigits for several cases.
DQMHistAnalysisECLOutOfTimeDigitsModule()
< derived from DQMHistAnalysisModule class.
std::string m_pvPrefix
Prefix to use for PVs registered by this module.
bool m_onlyIfUpdated
If true (default), update EPICS PVs only if there were changes in the histograms.
The base class for the histogram analysis module.
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.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
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)
REG_MODULE(arichBtest)
Register the Module.
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
Abstract base class for different kinds of events.