Belle II Software  release-06-00-14
DQMHistAnalysisEpicsExample.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 #include <dqm/analysis/modules/DQMHistAnalysisEpicsExample.h>
10 #include <TROOT.h>
11 #include <TClass.h>
12 
13 using namespace std;
14 using namespace Belle2;
15 
16 //-----------------------------------------------------------------
17 // Register the Module
18 //-----------------------------------------------------------------
19 REG_MODULE(DQMHistAnalysisEpicsExample)
20 
21 //-----------------------------------------------------------------
22 // Implementation
23 //-----------------------------------------------------------------
24 
27 {
28  // This module CAN NOT be run in parallel!
29 
30  //Parameter definition
31  addParam("HistoName", m_histoname, "Name of Histogram (incl dir)", std::string(""));
32  addParam("Function", m_function, "Fit function definition", std::string("gaus"));
33  addParam("Parameters", m_parameters, "Fit function parameters for EPICS", 3);
34  addParam("PVName", m_pvPrefix, "PV Prefix", std::string("DQM:TEST:hist:"));
35  B2DEBUG(20, "DQMHistAnalysisEpicsExample: Constructor done.");
36 }
37 
38 
39 DQMHistAnalysisEpicsExampleModule::~DQMHistAnalysisEpicsExampleModule()
40 {
41 #ifdef _BELLE2_EPICS
42  if (ca_current_context()) ca_context_destroy();
43 #endif
44 }
45 
46 void DQMHistAnalysisEpicsExampleModule::initialize()
47 {
48  B2DEBUG(20, "DQMHistAnalysisEpicsExample: initialized.");
49 
50  TString a;
51  a = m_histoname;
52  a.ReplaceAll("/", "_");
53  m_c1 = new TCanvas("c_" + a);
54  m_f1 = new TF1("f_" + a, TString(m_function), -30, 300);
55  m_f1->SetParameter(0, 1000);
56  m_f1->SetParameter(1, 0);
57  m_f1->SetParameter(2, 10);
58  m_f1->SetLineColor(4);
59  m_f1->SetNpx(512);
60  m_f1->SetNumberFitPoints(512);
61 
62  m_line = new TLine(0, 10, 0, 0);
63  m_line->SetVertical(true);
64  m_line->SetLineColor(8);
65  m_line->SetLineWidth(3);
66 
67  m_line_lo = new TLine(0, 10, 0, 0);
68  m_line_lo->SetVertical(true);
69  m_line_lo->SetLineColor(2);
70  m_line_lo->SetLineWidth(3);
71 
72  m_line_hi = new TLine(0, 10, 0, 0);
73  m_line_hi->SetVertical(true);
74  m_line_hi->SetLineColor(2);
75  m_line_hi->SetLineWidth(3);
76 
77  m_line_lo->SetX1(5);// get from epics
78  m_line_lo->SetX2(5);
79 
80  m_line_hi->SetX1(50);// get from epics
81  m_line_hi->SetX2(50);
82 
83  // need the function to get parameter names
84  if (m_parameters > 0) {
85  if (m_parameters > 10) m_parameters = 10; // hard limit
86 #ifdef _BELLE2_EPICS
87  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
88 #endif
89  for (auto i = 0; i < m_parameters; i++) {
90  std::string aa;
91  aa = m_f1->GetParName(i);
92  if (aa == "") aa = string("par") + string(TString::Itoa(i, 10).Data());
93  aa = m_pvPrefix + aa;
94 #ifdef _BELLE2_EPICS
95  SEVCHK(ca_create_channel(aa.c_str(), NULL, NULL, 10, &mychid[i]), "ca_create_channel failure");
96  // Read LO and HI limits from EPICS, seems this needs additional channels?
97  // SEVCHK(ca_get(DBR_DOUBLE,mychid[i],(void*)&data),"ca_get failure"); // data is only valid after ca_pend_io!!
98 #endif
99  }
100 #ifdef _BELLE2_EPICS
101  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
102 #endif
103  } else {
104  m_parameters = 0;
105  }
106 }
107 
108 
109 void DQMHistAnalysisEpicsExampleModule::beginRun()
110 {
111  //m_serv->SetTimer(100, kFALSE);
112  B2DEBUG(20, "DQMHistAnalysisEpicsExample: beginRun called.");
113  m_c1->Clear();
114 
115  TH1* hh1;
116  hh1 = findHist(m_histoname.c_str());
117 
118  if (hh1 == NULL) {
119  B2DEBUG(20, "Histo " << m_histoname << " not in memfile");
120  TDirectory* d = gROOT;
121  TString myl = m_histoname;
122  TString tok;
123  Ssiz_t from = 0;
124  while (myl.Tokenize(tok, from, "/")) {
125  TString dummy;
126  Ssiz_t f;
127  f = from;
128  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
129  auto e = d->GetDirectory(tok);
130  if (e) {
131  B2DEBUG(20, "Cd Dir " << tok);
132  d = e;
133  }
134  d->cd();
135  } else {
136  break;
137  }
138  }
139  TObject* obj = d->FindObject(tok);
140  if (obj != NULL) {
141  if (obj->IsA()->InheritsFrom("TH1")) {
142  B2DEBUG(20, "Histo " << m_histoname << " found in mem");
143  hh1 = (TH1*)obj;
144  }
145  } else {
146  B2DEBUG(20, "Histo " << m_histoname << " NOT found in mem");
147  }
148  }
149 
150  if (hh1 != NULL) {
151  m_c1->cd();
152  hh1->Draw();
153  m_line->Draw();
154  m_line_lo->Draw();
155  m_line_hi->Draw();
156  } else {
157  B2DEBUG(20, "Histo " << m_histoname << " not found");
158  }
159 }
160 
161 void DQMHistAnalysisEpicsExampleModule::event()
162 {
163  TH1* hh1;
164  bool flag = false;
165 
166  hh1 = findHist(m_histoname.c_str());
167  if (hh1 == NULL) {
168  B2DEBUG(20, "Histo " << m_histoname << " not in memfile");
169  TDirectory* d = gROOT;
170  TString myl = m_histoname;
171  TString tok;
172  Ssiz_t from = 0;
173  while (myl.Tokenize(tok, from, "/")) {
174  TString dummy;
175  Ssiz_t f;
176  f = from;
177  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
178  auto e = d->GetDirectory(tok);
179  if (e) {
180  B2DEBUG(20, "Cd Dir " << tok);
181  d = e;
182  }
183  d->cd();
184  } else {
185  break;
186  }
187  }
188  TObject* obj = d->FindObject(tok);
189  if (obj != NULL) {
190  if (obj->IsA()->InheritsFrom("TH1")) {
191  B2DEBUG(20, "Histo " << m_histoname << " found in mem");
192  hh1 = (TH1*)obj;
193  flag = true;
194  }
195  } else {
196  B2DEBUG(20, "Histo " << m_histoname << " NOT found in mem");
197  }
198  }
199  if (hh1 != NULL) {
200  m_c1->cd();// necessary!
201  hh1->Fit(m_f1, "");
202  double y1 = hh1->GetMaximum();
203  double y2 = hh1->GetMinimum();
204  m_line->SetY1(y1 + (y1 - y2) * 0.05);
205  m_line_lo->SetY1(y1 + (y1 - y2) * 0.05);
206  m_line_hi->SetY1(y1 + (y1 - y2) * 0.05);
207 // m_line->SetY2(y2-(y1-y2)*0.05);
208 // m_line_lo->SetY2(y2-(y1-y2)*0.05);
209 // m_line_hi->SetY2(y2-(y1-y2)*0.05);
210  double x = m_f1->GetParameter(1);
211  m_line->SetX1(x);
212  m_line->SetX2(x);
213  if (!flag) {
214  // dont add another line...
215  m_line->Draw();
216  m_line_lo->Draw();
217  m_line_hi->Draw();
218  }
219  m_c1->Modified();
220  m_c1->Update();
221  } else {
222  B2DEBUG(20, "Histo " << m_histoname << " not found");
223  }
224 
225 #ifdef _BELLE2_EPICS
226  if (m_parameters > 0) {
227  for (auto i = 0; i < m_parameters; i++) {
228  double data;
229  data = m_f1->GetParameter(i);
230  if (mychid[i]) SEVCHK(ca_put(DBR_DOUBLE, mychid[i], (void*)&data), "ca_set failure");
231  }
232  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
233  }
234 #endif
235 }
236 
237 void DQMHistAnalysisEpicsExampleModule::endRun()
238 {
239  B2DEBUG(20, "DQMHistAnalysisEpicsExample : endRun called");
240 }
241 
242 
243 void DQMHistAnalysisEpicsExampleModule::terminate()
244 {
245 #ifdef _BELLE2_EPICS
246  if (m_parameters > 0) {
247  for (auto i = 0; i < m_parameters; i++) {
248  if (mychid[i]) SEVCHK(ca_clear_channel(mychid[i]), "ca_clear_channel failure");
249  }
250  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
251  }
252 #endif
253  B2DEBUG(20, "DQMHistAnalysisEpicsExample: terminate called");
254 }
255 
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.