Belle II Software  release-06-00-14
DQMHistAnalysisIP.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 : DQMHistAnalysisIP.cc
10 // Description : Mean for IP position with delta histogramming
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisIP.h>
15 #include <TROOT.h>
16 #include <TClass.h>
17 
18 using namespace std;
19 using namespace Belle2;
20 
21 //-----------------------------------------------------------------
22 // Register the Module
23 //-----------------------------------------------------------------
24 REG_MODULE(DQMHistAnalysisIP)
25 
26 //-----------------------------------------------------------------
27 // Implementation
28 //-----------------------------------------------------------------
29 
32 {
33  // This module CAN NOT be run in parallel!
34 
35  //Parameter definition
36  addParam("HistoName", m_histoname, "Name of Histogram (incl dir)", std::string(""));
37  addParam("PVName", m_pvPrefix, "PV Prefix", std::string("DQM:TEST:hist:"));
38  addParam("MonitorPrefix", m_monPrefix, "Monitor Prefix");// force to be set!
39  addParam("useEpics", m_useEpics, "useEpics", true);
40  addParam("minEntries", m_minEntries, "minimum number of new Entries for a fit", 1000);
41  B2DEBUG(20, "DQMHistAnalysisIP: Constructor done.");
42 }
43 
44 
45 DQMHistAnalysisIPModule::~DQMHistAnalysisIPModule()
46 {
47 #ifdef _BELLE2_EPICS
48  if (m_useEpics) {
49  if (ca_current_context()) ca_context_destroy();
50  }
51 #endif
52 }
53 
54 void DQMHistAnalysisIPModule::initialize()
55 {
56  B2DEBUG(20, "DQMHistAnalysisIP: initialized.");
57 
58  m_monObj = getMonitoringObject("ip");
59 
60  TString a;
61  a = m_histoname;
62  m_c1 = new TCanvas(a + "_fit"); // prefer to change canvas name to monitorPrefix, but then chnages on the web gui are needed :-(
63 
64  m_line = new TLine(0, 10, 0, 0);
65  m_line->SetVertical(true);
66  m_line->SetLineColor(8);
67  m_line->SetLineWidth(3);
68 
69  m_monObj->addCanvas(m_c1);
70 
71  // need the function to get parameter names
72 #ifdef _BELLE2_EPICS
73  if (m_useEpics) {
74  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
75  std::string aa;
76  aa = m_pvPrefix + "Mean";
77  SEVCHK(ca_create_channel(aa.c_str(), NULL, NULL, 10, &mychid[0]), "ca_create_channel failure");
78  aa = m_pvPrefix + "RMS";
79  SEVCHK(ca_create_channel(aa.c_str(), NULL, NULL, 10, &mychid[1]), "ca_create_channel failure");
80  // Read LO and HI limits from EPICS, seems this needs additional channels?
81  // SEVCHK(ca_get(DBR_DOUBLE,mychid[i],(void*)&data),"ca_get failure"); // data is only valid after ca_pend_io!!
82  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
83  }
84 #endif
85 }
86 
87 
88 void DQMHistAnalysisIPModule::beginRun()
89 {
90  //m_serv->SetTimer(100, kFALSE);
91  B2DEBUG(20, "DQMHistAnalysisIP: beginRun called.");
92  m_c1->Clear();
93 
94  TH1* hh1;
95  hh1 = findHist(m_histoname.c_str());
96 
97  if (hh1 == NULL) {
98  B2DEBUG(20, "Histo " << m_histoname << " not in memfile");
99  TDirectory* d = gROOT;
100  TString myl = m_histoname;
101  TString tok;
102  Ssiz_t from = 0;
103  while (myl.Tokenize(tok, from, "/")) {
104  TString dummy;
105  Ssiz_t f;
106  f = from;
107  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
108  auto e = d->GetDirectory(tok);
109  if (e) {
110  B2DEBUG(20, "Cd Dir " << tok);
111  d = e;
112  }
113  d->cd();
114  } else {
115  break;
116  }
117  }
118  TObject* obj = d->FindObject(tok);
119  if (obj != NULL) {
120  if (obj->IsA()->InheritsFrom("TH1")) {
121  B2DEBUG(20, "Histo " << m_histoname << " found in mem");
122  hh1 = (TH1*)obj;
123  }
124  } else {
125  B2DEBUG(20, "Histo " << m_histoname << " NOT found in mem");
126  }
127  }
128 
129  if (hh1 != NULL) {
130  m_c1->cd();
131  hh1->Draw();
132  m_line->Draw();
133  } else {
134  B2DEBUG(20, "Histo " << m_histoname << " not found");
135  }
136 
137  if (m_h_last) m_h_last->Reset();
138 }
139 
140 void DQMHistAnalysisIPModule::event()
141 {
142  TH1* hh1;
143  bool flag = false;
144 
145  hh1 = findHist(m_histoname.c_str());
146  if (hh1 == NULL) {
147  B2DEBUG(20, "Histo " << m_histoname << " not in memfile");
148  TDirectory* d = gROOT;
149  TString myl = m_histoname;
150  TString tok;
151  Ssiz_t from = 0;
152  while (myl.Tokenize(tok, from, "/")) {
153  TString dummy;
154  Ssiz_t f;
155  f = from;
156  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
157  auto e = d->GetDirectory(tok);
158  if (e) {
159  B2DEBUG(20, "Cd Dir " << tok);
160  d = e;
161  }
162  d->cd();
163  } else {
164  break;
165  }
166  }
167  TObject* obj = d->FindObject(tok);
168  if (obj != NULL) {
169  if (obj->IsA()->InheritsFrom("TH1")) {
170  B2DEBUG(20, "Histo " << m_histoname << " found in mem");
171  hh1 = (TH1*)obj;
172  flag = true;
173  }
174  } else {
175  B2DEBUG(20, "Histo " << m_histoname << " NOT found in mem");
176  }
177  }
178  if (hh1 != NULL) {
179 
180  if (!m_h_last) {
181  B2DEBUG(20, "Initial Clone");
182  gROOT->cd();
183  m_h_last = (TH1*)hh1->Clone();
184  m_h_last->SetName(TString(hh1->GetName()) + "_last");
185  m_h_last->Reset();
186  }
187 
188  double last = m_h_last->GetEntries();
189  double current = hh1->GetEntries();
190  B2DEBUG(20, "Entries: " << last << "," << current);
191  if (current - last >= m_minEntries) {
192  B2DEBUG(20, "Update Delta & Fit");
193  gROOT->cd();
194  TH1* delta = (TH1*)hh1->Clone();
195  delta->Add(m_h_last, -1.);
196 
197  m_c1->Clear();
198  m_c1->cd();// necessary!
199  delta->Draw("hist");
200  m_h_last->Reset();
201  m_h_last->Add(hh1);
202 
203  double x = delta->GetMean();// must be double bc of EPICS below
204  double w = delta->GetRMS();// must be double bc of EPICS below
205  double y1 = delta->GetMaximum();
206  double y2 = delta->GetMinimum();
207  B2DEBUG(20, "Fit " << x << "," << w << "," << y1 << "," << y2);
208  m_line->SetY1(y1 + (y1 - y2) * 0.05);
209  m_line->SetX1(x);
210  m_line->SetX2(x);
211  delta->GetXaxis()->SetRangeUser(x - 3 * w, x + 3 * w);
212  if (!flag) {
213  // dont add another line...
214  m_line->Draw();
215  }
216  m_c1->Modified();
217  m_c1->Update();
218 
219  m_monObj->setVariable(m_monPrefix + "_mean", x);
220  m_monObj->setVariable(m_monPrefix + "_width", w);
221 
222 #ifdef _BELLE2_EPICS
223  if (m_useEpics) {
224  B2INFO("Update EPICS");
225  if (mychid[0]) SEVCHK(ca_put(DBR_DOUBLE, mychid[0], (void*)&x), "ca_set failure");
226  if (mychid[1]) SEVCHK(ca_put(DBR_DOUBLE, mychid[1], (void*)&w), "ca_set failure");
227  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
228  }
229 #endif
230  }
231  } else {
232  B2DEBUG(20, "Histo " << m_histoname << " not found");
233  }
234 
235 }
236 
237 void DQMHistAnalysisIPModule::terminate()
238 {
239 #ifdef _BELLE2_EPICS
240  if (m_useEpics) {
241  for (auto i = 0; i < m_parameters; i++) {
242  if (mychid[i]) SEVCHK(ca_clear_channel(mychid[i]), "ca_clear_channel failure");
243  }
244  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
245  }
246 #endif
247  B2DEBUG(20, "DQMHistAnalysisIP: terminate called");
248 }
249 
Class definition for the output module of Sequential ROOT I/O.
The base class for the histogram analysis module.
#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.