Belle II Software  release-06-02-00
DQMHistAnalysisPXDER.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/DQMHistAnalysisPXDER.h"
10 
11 #include <pxd/geometry/SensorInfo.h>
12 #include <vxd/geometry/SensorInfoBase.h>
13 #include <vxd/geometry/GeoCache.h>
14 
15 #include <boost/format.hpp>
16 
17 #include <TClass.h>
18 #include <TKey.h>
19 #include <TDirectory.h>
20 #include <TFile.h>
21 #include <TROOT.h>
22 
23 #include <memory>
24 
25 using namespace std;
26 using boost::format;
27 using namespace Belle2;
28 
29 //-----------------------------------------------------------------
30 // Register the Module
31 //-----------------------------------------------------------------
32 REG_MODULE(DQMHistAnalysisPXDER)
33 
34 
35 //-----------------------------------------------------------------
36 // Implementation
37 //-----------------------------------------------------------------
38 
40 {
41  //Set module properties
42  setDescription("PXD DQM analysis module for Express Reco ");
43  addParam("RefHistoFile", m_refFileName, "Reference histrogram file name", std::string("refHisto.root"));
44 
45  // NO parallel processing
46 }
47 
48 
49 DQMHistAnalysisPXDERModule::~DQMHistAnalysisPXDERModule()
50 {
51 }
52 
53 void DQMHistAnalysisPXDERModule::initialize()
54 {
55  m_refFile = NULL;
56  if (m_refFileName != "") {
57  m_refFile = new TFile(m_refFileName.data());
58  }
59 
60  gROOT->cd(); // this seems to be important, or strange things happen
61  // basic constants presets:
62  VXD::GeoCache& geo = VXD::GeoCache::getInstance();
63  c_nVXDLayers = geo.getLayers().size();
64  c_firstVXDLayer = 1; // counting start from 1...
65  c_lastVXDLayer = c_nVXDLayers;
66  c_nPXDLayers = geo.getLayers(VXD::SensorInfoBase::SensorType::PXD).size();
67  c_firstPXDLayer = c_firstVXDLayer;
68  c_lastPXDLayer = c_nPXDLayers;
69  c_nSVDLayers = geo.getLayers(VXD::SensorInfoBase::SensorType::SVD).size();
70  c_firstSVDLayer = c_nPXDLayers + c_firstPXDLayer;
71  c_lastSVDLayer = c_firstSVDLayer + c_nSVDLayers;
72 
73  c_nPXDSensors = 0;
74  for (VxdID layer : geo.getLayers()) {
75  for (VxdID ladder : geo.getLadders(layer)) {
76  if (layer.getLayerNumber() <= c_lastPXDLayer) { // PXD
77  c_nPXDSensors += geo.getLadders(layer).size() * geo.getSensors(ladder).size();
78  }
79  break;
80  }
81  }
82 
83 // m_hitMapCounts = "DQMER_PXD_PixelHitmapCounts";
84 // m_hitMapClCounts = "DQMER_PXD_ClusterHitmapCounts";
85 
86  for (int i = 0; i < c_nPXDSensors; i++) {
87  int iLayer = 0;
88  int iLadder = 0;
89  int iSensor = 0;
90  getIDsFromIndex(i, iLayer, iLadder, iSensor);
91  VxdID sensorID(iLayer, iLadder, iSensor);
92  PXD::SensorInfo SensorInfo = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
93  string sensorDescr = str(format("%1%_%2%_%3%") % iLayer % iLadder % iSensor);
94  //----------------------------------------------------------------
95  // Number of fired pixels per frame
96  //----------------------------------------------------------------
97  m_fired.emplace_back(str(format("DQMER_PXD_%1%_Fired") % sensorDescr));
98  m_ref_fired.emplace_back("ref/" + m_fired.back());
99  //----------------------------------------------------------------
100  // Number of clusters per frame
101  //----------------------------------------------------------------
102  m_clusters.emplace_back(str(format("DQMER_PXD_%1%_Clusters") % sensorDescr));
103  m_ref_clusters.emplace_back("ref/" + m_clusters.back());
104  //----------------------------------------------------------------
105  // Start row distribution
106  //----------------------------------------------------------------
107  m_startRow.emplace_back(str(format("DQMER_PXD_%1%_StartRow") % sensorDescr));
108  m_ref_startRow.emplace_back("ref/" + m_startRow.back());
109  //----------------------------------------------------------------
110  // Cluster seed charge by distance from the start row
111  //----------------------------------------------------------------
112  m_chargStartRow.emplace_back(str(format("DQMER_PXD_%1%_AverageSeedByStartRow") % sensorDescr));
113  m_ref_chargStartRow.emplace_back("ref/" + m_chargStartRow.back());
114 
115 
116  m_startRowCount.emplace_back(str(format("DQMER_PXD_%1%_SeedCountsByStartRow") % sensorDescr));
117  m_ref_startRowCount.emplace_back("ref/" + m_startRowCount.back());
118  //----------------------------------------------------------------
119  // Cluster Charge
120  //----------------------------------------------------------------
121  m_clusterCharge.emplace_back(str(format("DQMER_PXD_%1%_ClusterCharge") % sensorDescr));
122  m_ref_clusterCharge.emplace_back("ref/" + m_clusterCharge.back());
123  //----------------------------------------------------------------
124  // Pixel Signal
125  //----------------------------------------------------------------
126  m_pixelSignal.emplace_back(str(format("DQMER_PXD_%1%_PixelSignal") % sensorDescr));
127  m_ref_pixelSignal.emplace_back("ref/" + m_pixelSignal.back());
128  //----------------------------------------------------------------
129  // Cluster Size in U
130  //----------------------------------------------------------------
131  m_clusterSizeU.emplace_back(str(format("DQMER_PXD_%1%_ClusterSizeU") % sensorDescr));
132  m_ref_clusterSizeU.emplace_back("ref/" + m_clusterSizeU.back());
133  //----------------------------------------------------------------
134  // Cluster Size in V
135  //----------------------------------------------------------------
136  m_clusterSizeV.emplace_back(str(format("DQMER_PXD_%1%_ClusterSizeV") % sensorDescr));
137  m_ref_clusterSizeV.emplace_back("ref/" + m_clusterSizeV.back());
138  //----------------------------------------------------------------
139  // Cluster Size in U+V
140  //----------------------------------------------------------------
141  m_clusterSizeUV.emplace_back(str(format("DQMER_PXD_%1%_ClusterSizeUV") % sensorDescr));
142  m_ref_clusterSizeUV.emplace_back("ref/" + m_clusterSizeUV.back());
143  }
144 // m_fHitMapCountsFlag = NULL;
145 // m_fHitMapClCountsFlag = NULL;
146 // m_hitMapCounts = NULL;
147 // m_hitMapClCounts = NULL;
148 
149 // Create flag histograms:
150 // DirPXDFlags->cd();
151  m_fFiredFlag = new TH1I("DQMER_PXD_FiredFlag", "DQM ER PXD Fired Flag",
152  c_nPXDSensors, 0, c_nPXDSensors);
153  m_fFiredFlag->GetXaxis()->SetTitle("Sensor ID");
154  m_fFiredFlag->GetYaxis()->SetTitle("flag");
155  m_fClustersFlag = new TH1I("DQMER_PXD_ClustersFlag", "DQM ER PXD Clusters Flag",
156  c_nPXDSensors, 0, c_nPXDSensors);
157  m_fClustersFlag->GetXaxis()->SetTitle("Sensor ID");
158  m_fClustersFlag->GetYaxis()->SetTitle("flag");
159  m_fStartRowFlag = new TH1I("DQMER_PXD_StartRowFlag", "DQM ER PXD Start Row Flag",
160  c_nPXDSensors, 0, c_nPXDSensors);
161  m_fStartRowFlag->GetXaxis()->SetTitle("Sensor ID");
162  m_fStartRowFlag->GetYaxis()->SetTitle("flag");
163  m_fChargStartRowFlag = new TH1I("DQMER_PXD_ChargStartRowFlag", "DQM ER PXD Charg Start Row Flag",
164  c_nPXDSensors, 0, c_nPXDSensors);
165  m_fChargStartRowFlag->GetXaxis()->SetTitle("Sensor ID");
166  m_fChargStartRowFlag->GetYaxis()->SetTitle("flag");
167  m_fStartRowCountFlag = new TH1I("DQMER_PXD_StartRowCountFlag", "DQM ER PXD Row Count Flag",
168  c_nPXDSensors, 0, c_nPXDSensors);
169  m_fStartRowCountFlag->GetXaxis()->SetTitle("Sensor ID");
170  m_fStartRowCountFlag->GetYaxis()->SetTitle("flag");
171 // m_fHitMapCountsFlag = new TH1I("DQMER_PXD_PixelHitmapCountsFlag", "DQM ER PXD Pixel Hitmaps Counts Flag",
172 // c_nPXDSensors, 0, c_nPXDSensors);
173 // m_fHitMapCountsFlag->GetXaxis()->SetTitle("Sensor ID");
174 // m_fHitMapCountsFlag->GetYaxis()->SetTitle("flag");
175 // m_fHitMapClCountsFlag = new TH1I("DQMER_PXD_ClusterHitmapCountsFlag", "DQM ER PXD Cluster Hitmaps Counts Flag",
176 // c_nPXDSensors, 0, c_nPXDSensors);
177 // m_fHitMapClCountsFlag->GetXaxis()->SetTitle("Sensor ID");
178 // m_fHitMapClCountsFlag->GetYaxis()->SetTitle("flag");
179  m_fClusterChargeFlag = new TH1I("DQMER_PXD_ClusterChargeFlag", "DQM ER PXD Cluster Charge Flag",
180  c_nPXDSensors, 0, c_nPXDSensors);
181  m_fClusterChargeFlag->GetXaxis()->SetTitle("Sensor ID");
182  m_fClusterChargeFlag->GetYaxis()->SetTitle("flag");
183  m_fPixelSignalFlag = new TH1I("DQMER_PXD_PixelSignalFlag", "DQM ER PXD Pixel Signal Flag",
184  c_nPXDSensors, 0, c_nPXDSensors);
185  m_fPixelSignalFlag->GetXaxis()->SetTitle("Sensor ID");
186  m_fPixelSignalFlag->GetYaxis()->SetTitle("flag");
187  m_fClusterSizeUFlag = new TH1I("DQMER_PXD_ClasterSizeUFlag", "DQM ER PXD Cluster Size U Flag",
188  c_nPXDSensors, 0, c_nPXDSensors);
189  m_fClusterSizeUFlag->GetXaxis()->SetTitle("Sensor ID");
190  m_fClusterSizeUFlag->GetYaxis()->SetTitle("flag");
191  m_fClusterSizeVFlag = new TH1I("DQMER_PXD_ClasterSizeVFlag", "DQM ER PXD Cluster Size V Flag",
192  c_nPXDSensors, 0, c_nPXDSensors);
193  m_fClusterSizeVFlag->GetXaxis()->SetTitle("Sensor ID");
194  m_fClusterSizeVFlag->GetYaxis()->SetTitle("flag");
195  m_fClusterSizeUVFlag = new TH1I("DQMER_PXD_ClasterSizeUVFlag", "DQM ER PXD Cluster Size UV Flag",
196  c_nPXDSensors, 0, c_nPXDSensors);
197  m_fClusterSizeUVFlag->GetXaxis()->SetTitle("Sensor ID");
198  m_fClusterSizeUVFlag->GetYaxis()->SetTitle("flag");
199 
200  for (int i = 0; i < c_nPXDSensors; i++) {
201  int iLayer = 0;
202  int iLadder = 0;
203  int iSensor = 0;
204  getIDsFromIndex(i, iLayer, iLadder, iSensor);
205  TString AxisTicks = Form("%i_%i_%i", iLayer, iLadder, iSensor);
206 // m_hitMapCounts->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
207 // m_hitMapClCounts->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
208  m_fFiredFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
209  m_fClustersFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
210  m_fStartRowFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
211  m_fChargStartRowFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
212  m_fStartRowCountFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
213 // m_fHitMapCountsFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
214 // m_fHitMapClCountsFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
215  m_fClusterChargeFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
216  m_fPixelSignalFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
217  m_fClusterSizeUFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
218  m_fClusterSizeVFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
219  m_fClusterSizeUVFlag->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
220  }
221 // m_oldDir->cd();
222 }
223 
224 void DQMHistAnalysisPXDERModule::beginRun()
225 {
226  // Just to make sure, reset all the histograms.
227 // if (m_fHitMapCountsFlag != NULL) m_fHitMapCountsFlag->Reset();
228 // if (m_fHitMapClCountsFlag != NULL) m_fHitMapClCountsFlag->Reset();
229  if (m_fFiredFlag != NULL) m_fFiredFlag->Reset();
230  if (m_fClustersFlag != NULL) m_fClustersFlag->Reset();
231  if (m_fStartRowFlag != NULL) m_fStartRowFlag->Reset();
232  if (m_fChargStartRowFlag != NULL) m_fChargStartRowFlag->Reset();
233  if (m_fStartRowCountFlag != NULL) m_fStartRowCountFlag->Reset();
234  if (m_fClusterChargeFlag != NULL) m_fClusterChargeFlag->Reset();
235  if (m_fPixelSignalFlag != NULL) m_fPixelSignalFlag->Reset();
236  if (m_fClusterSizeUFlag != NULL) m_fClusterSizeUFlag->Reset();
237  if (m_fClusterSizeVFlag != NULL) m_fClusterSizeVFlag->Reset();
238  if (m_fClusterSizeUVFlag != NULL) m_fClusterSizeUVFlag->Reset();
239 
240 // if (m_hitMapCounts != NULL) m_hitMapCounts->Reset();
241 // if (m_hitMapClCounts != NULL) m_hitMapClCounts->Reset();
242 
243 }
244 
245 
246 void DQMHistAnalysisPXDERModule::event()
247 {
248 
249  // Dont sum up!
250 // if (m_fHitMapCountsFlag != NULL) m_fHitMapCountsFlag->Reset();
251 // if (m_fHitMapClCountsFlag != NULL) m_fHitMapClCountsFlag->Reset();
252  if (m_fFiredFlag != NULL) m_fFiredFlag->Reset();
253  if (m_fClustersFlag != NULL) m_fClustersFlag->Reset();
254  if (m_fStartRowFlag != NULL) m_fStartRowFlag->Reset();
255  if (m_fChargStartRowFlag != NULL) m_fChargStartRowFlag->Reset();
256  if (m_fStartRowCountFlag != NULL) m_fStartRowCountFlag->Reset();
257  if (m_fClusterChargeFlag != NULL) m_fClusterChargeFlag->Reset();
258  if (m_fPixelSignalFlag != NULL) m_fPixelSignalFlag->Reset();
259  if (m_fClusterSizeUFlag != NULL) m_fClusterSizeUFlag->Reset();
260  if (m_fClusterSizeVFlag != NULL) m_fClusterSizeVFlag->Reset();
261  if (m_fClusterSizeUVFlag != NULL) m_fClusterSizeUVFlag->Reset();
262 
263  // Compare histograms with reference histograms and create flags:
264  for (int i = 0; i < c_nPXDSensors; i++) {
265  double pars[2];
266  pars[0] = 0.01;// Probabilty value error?
267  pars[1] = 0.05;// Probabilty value warning?
268 
269  double m_NoOfEvents = 1., m_NoOfEventsRef = 1.; // workaround
270 
271 // SetFlag(9, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
272 // m_hitMapCounts, r_hitMapCounts, m_fHitMapCountsFlag);
273 // SetFlag(9, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
274 // m_hitMapClCounts, r_hitMapClCounts, m_fHitMapClCountsFlag);
275  SetFlag(2, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
276  m_fired.at(i), m_ref_fired.at(i), m_fFiredFlag);
277  SetFlag(2, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
278  m_clusters.at(i), m_ref_clusters.at(i), m_fClustersFlag);
279  SetFlag(100, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
280  m_startRow.at(i), m_ref_startRow.at(i), m_fStartRowFlag);
281  SetFlag(100, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
282  m_chargStartRow.at(i), m_ref_chargStartRow.at(i), m_fChargStartRowFlag);
283  SetFlag(100, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
284  m_startRowCount.at(i), m_ref_startRowCount.at(i), m_fStartRowCountFlag);
285  SetFlag(5, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
286  m_clusterCharge.at(i), m_ref_clusterCharge.at(i), m_fClusterChargeFlag);
287  SetFlag(5, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
288  m_pixelSignal.at(i), m_ref_pixelSignal.at(i), m_fPixelSignalFlag);
289  SetFlag(2, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
290  m_clusterSizeU.at(i), m_ref_clusterSizeU.at(i), m_fClusterSizeUFlag);
291  SetFlag(2, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
292  m_clusterSizeV.at(i), m_ref_clusterSizeV.at(i), m_fClusterSizeVFlag);
293  SetFlag(2, i, pars, (double)m_NoOfEvents / m_NoOfEventsRef,
294  m_clusterSizeUV.at(i), m_ref_clusterSizeUV.at(i), m_fClusterSizeUVFlag);
295 
296  }
297 }
298 
299 void DQMHistAnalysisPXDERModule::endRun()
300 {
301 }
302 
303 void DQMHistAnalysisPXDERModule::terminate()
304 {
305 }
306 
307 void DQMHistAnalysisPXDERModule::getIDsFromIndex(const int Index, int& Layer, int& Ladder, int& Sensor) const
308 {
309  VXD::GeoCache& geo = VXD::GeoCache::getInstance();
310  int tempcounter = 0;
311  for (VxdID layer : geo.getLayers()) {
312  if (layer.getLayerNumber() > c_lastPXDLayer) continue; // need PXD
313  for (VxdID ladder : geo.getLadders(layer)) {
314  for (VxdID sensor : geo.getSensors(ladder)) {
315  if (tempcounter == Index) {
316  Layer = layer.getLayerNumber();
317  Ladder = ladder.getLadderNumber();
318  Sensor = sensor.getSensorNumber();
319  return;
320  }
321  tempcounter++;
322  }
323  }
324  }
325 }
326 
327 int DQMHistAnalysisPXDERModule::SetFlag(int Type, int bin, double* pars, double ratio, const std::string& name_hist,
328  const std::string& name_refhist, TH1I* flaghist)
329 {
330  int iret = 0;
331  float WarningLevel = 6.0;
332  float ErrorLevel = 10.0;
333 
334  TH1* hist, *refhist;
335 
336  hist = GetHisto(name_hist);
337  if (!hist) return -1;
338  refhist = GetHisto(name_refhist);
339  if (!refhist) return -1;
340 
341  // What happens if they are TH1I, TH1D and not TH1F
342 
343  auto temp = std::unique_ptr<TH1F>(new TH1F("temp", "temp", hist->GetNbinsX(), hist->GetXaxis()->GetXmin(),
344  hist->GetXaxis()->GetXmax()));
345  double NEvents = 0;
346  double flagInt = 0;
347  double flagrInt = 0;
348  for (int j = 0; j < hist->GetNbinsX(); j++) {
349  double val = hist->GetBinContent(j + 1);
350  NEvents += val;
351  val = val / ratio;
352  temp->SetBinContent(j + 1, val);
353  flagInt += temp->GetBinContent(j + 1);
354  flagrInt += refhist->GetBinContent(j + 1);
355  }
356  if (NEvents < 100) { // not enough information for comparition
357  iret = -1;
358  flaghist->SetBinContent(bin + 1, -1);
359  return iret;
360  }
361  double flag = temp->GetMean();
362  double flagErr = temp->GetMeanError();
363  double flagRMS = temp->GetRMS();
364  double flagRMSErr = temp->GetRMSError();
365  double flagr = refhist->GetMean();
366  double flagrErr = refhist->GetMeanError();
367  double flagrRMS = refhist->GetRMS();
368  double flagrRMSErr = refhist->GetRMSError();
369  TString strDebugInfo = Form("Conditions for Flag--->\n source %f %f+-%f %f+-%f\n referen %f %f+-%f %f+-%f\n",
370  flagInt, flag, flagErr, flagRMS, flagRMSErr,
371  flagrInt, flagr, flagrErr, flagrRMS, flagrRMSErr
372  );
373  B2DEBUG(130, strDebugInfo.Data());
374  if (Type == 1) { // counts+mean+RMS use
375  if ((fabs(flag - flagr) > ErrorLevel * (flagErr + flagrErr)) ||
376  (fabs(flagRMS - flagrRMS) > ErrorLevel * (flagRMSErr + flagrRMSErr)) ||
377  (fabs(flagInt - flagrInt) > ErrorLevel * (sqrt(flagInt) + sqrt(flagrInt)))
378  ) {
379  flaghist->SetBinContent(bin + 1, 2);
380  } else if ((fabs(flag - flagr) > WarningLevel * (flagErr + flagrErr)) ||
381  (fabs(flagRMS - flagrRMS) > WarningLevel * (flagRMSErr + flagrRMSErr)) ||
382  (fabs(flagInt - flagrInt) > WarningLevel * (sqrt(flagInt) + sqrt(flagrInt)))
383  ) {
384  flaghist->SetBinContent(bin + 1, 1);
385  } else {
386  flaghist->SetBinContent(bin + 1, 0);
387  }
388  iret = 1;
389  } else if (Type == 2) { // counts use
390  if (fabs(flagInt - flagrInt) > ErrorLevel * (sqrt(flagInt) + sqrt(flagrInt))) {
391  flaghist->SetBinContent(bin + 1, 2);
392  } else if (fabs(flagInt - flagrInt) > WarningLevel * (sqrt(flagInt) + sqrt(flagrInt))) {
393  flaghist->SetBinContent(bin + 1, 1);
394  } else {
395  flaghist->SetBinContent(bin + 1, 0);
396  }
397  iret = 1;
398  } else if (Type == 3) { // mean use
399  if (fabs(flag - flagr) > ErrorLevel * (flagErr + flagrErr)) {
400  flaghist->SetBinContent(bin + 1, 2);
401  } else if (fabs(flag - flagr) > WarningLevel * (flagErr + flagrErr)) {
402  flaghist->SetBinContent(bin + 1, 1);
403  } else {
404  flaghist->SetBinContent(bin + 1, 0);
405  }
406  iret = 1;
407  } else if (Type == 4) { // RMS use
408  if (fabs(flagRMS - flagrRMS) > ErrorLevel * (flagRMSErr + flagrRMSErr)) {
409  flaghist->SetBinContent(bin + 1, 2);
410  } else if (fabs(flagRMS - flagrRMS) > WarningLevel * (flagRMSErr + flagrRMSErr)) {
411  flaghist->SetBinContent(bin + 1, 1);
412  } else {
413  flaghist->SetBinContent(bin + 1, 0);
414  }
415  iret = 1;
416  } else if (Type == 5) { // counts+mean use
417  if ((fabs(flag - flagr) > ErrorLevel * (flagErr + flagrErr)) ||
418  (fabs(flagInt - flagrInt) > ErrorLevel * (sqrt(flagInt) + sqrt(flagrInt)))
419  ) {
420  flaghist->SetBinContent(bin + 1, 2);
421  } else if ((fabs(flag - flagr) > WarningLevel * (flagErr + flagrErr)) ||
422  (fabs(flagInt - flagrInt) > WarningLevel * (sqrt(flagInt) + sqrt(flagrInt)))
423  ) {
424  flaghist->SetBinContent(bin + 1, 1);
425  } else {
426  flaghist->SetBinContent(bin + 1, 0);
427  }
428  iret = 1;
429  } else if (Type == 9) { // bin content use
430  flagInt = temp->GetBinContent(bin + 1);
431  flagrInt = refhist->GetBinContent(bin + 1);
432  if (fabs(flagInt - flagrInt) > ErrorLevel * (sqrt(flagInt) + sqrt(flagrInt))) {
433  flaghist->SetBinContent(bin + 1, 2);
434  } else if (fabs(flagInt - flagrInt) > WarningLevel * (sqrt(flagInt) + sqrt(flagrInt))) {
435  flaghist->SetBinContent(bin + 1, 1);
436  } else {
437  flaghist->SetBinContent(bin + 1, 0);
438  }
439  iret = 1;
440  } else if (Type == 10) {
441  float flag2 = refhist->Chi2Test(temp.get());
442  flaghist->SetBinContent(bin + 1, 0);
443  if (flag2 > pars[1])
444  flaghist->SetBinContent(bin + 1, 2);
445  if (flag2 > pars[0])
446  flaghist->SetBinContent(bin + 1, 1);
447  iret = 1;
448  } else if (Type == 100) {
449  flaghist->SetBinContent(bin + 1, 0);
450  iret = 1;
451  } else {
452  flaghist->SetBinContent(bin + 1, -3);
453  iret = -1;
454  }
455  strDebugInfo = Form("SetFlag---> %f, type %i\n", flaghist->GetBinContent(bin + 1), Type);
456  B2DEBUG(130, strDebugInfo.Data());
457  return iret;
458 }
459 
460 // int DQMHistAnalysisPXDERModule::SetFlag(int Type, int bin, double* pars, double ratio, std::string name_hist, std::string name_refhist, TH1I* flaghist)
461 // {
462 //
463 // TH1F* histF = new TH1F("histF", "histF", hist->GetNbinsX(), hist->GetXaxis()->GetXmin(), hist->GetXaxis()->GetXmax());
464 // TH1F* refhistF = new TH1F("refhistF", "refhistF", refhist->GetNbinsX(), refhist->GetXaxis()->GetXmin(),
465 // refhist->GetXaxis()->GetXmax());
466 // for (int j = 0; j < hist->GetNbinsX(); j++) {
467 // histF->SetBinContent(j + 1, hist->GetBinContent(j + 1));
468 // refhistF->SetBinContent(j + 1, refhist->GetBinContent(j + 1));
469 // }
470 // int ret = SetFlag(Type, bin, pars, ratio, histF, refhistF, flaghist);
471 // delete histF;
472 // delete refhistF;
473 // return ret;
474 // }
475 
476 TH1* DQMHistAnalysisPXDERModule::GetHisto(TString histoname)
477 {
478  TH1* hh1;
479  hh1 = findHist(histoname.Data());
480  if (hh1 == NULL) {
481  B2INFO("Histo " << histoname << " not in memfile");
482 
483  // first search reference root file ... if ther is one
484  if (m_refFile && m_refFile->IsOpen()) {
485  TDirectory* d = m_refFile;
486  TString myl = histoname;
487  TString tok;
488  Ssiz_t from = 0;
489  B2INFO(myl);
490  while (myl.Tokenize(tok, from, "/")) {
491  TString dummy;
492  Ssiz_t f;
493  f = from;
494  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
495  auto e = d->GetDirectory(tok);
496  if (e) {
497  B2INFO("Cd Dir " << tok);
498  d = e;
499  } else {
500  B2INFO("cd failed " << tok);
501  }
502  } else {
503  break;
504  }
505  }
506  TObject* obj = d->FindObject(tok);
507  if (obj != NULL) {
508  if (obj->IsA()->InheritsFrom("TH1")) {
509  B2INFO("Histo " << histoname << " found in ref file");
510  hh1 = (TH1*)obj;
511  } else {
512  B2INFO("Histo " << histoname << " found in ref file but wrong type");
513  }
514  } else {
515  // seems find will only find objects, not keys, thus get the object on first access
516  TIter next(d->GetListOfKeys());
517  TKey* key;
518  while ((key = (TKey*)next())) {
519  TObject* obj2 = key->ReadObj() ;
520  if (obj2->InheritsFrom("TH1")) {
521  if (obj2->GetName() == tok) {
522  hh1 = (TH1*)obj2;
523  B2INFO("Histo " << histoname << " found as key -> readobj");
524  break;
525  }
526  }
527  }
528  if (hh1 == NULL) B2INFO("Histo " << histoname << " NOT found in ref file " << tok);
529  }
530  }
531 
532  if (hh1 == NULL) {
533  B2INFO("Histo " << histoname << " not in memfile or ref file");
534 
535  TDirectory* d = gROOT;
536  TString myl = histoname;
537  TString tok;
538  Ssiz_t from = 0;
539  while (myl.Tokenize(tok, from, "/")) {
540  TString dummy;
541  Ssiz_t f;
542  f = from;
543  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
544  auto e = d->GetDirectory(tok);
545  if (e) {
546  B2INFO("Cd Dir " << tok);
547  d = e;
548  } else B2INFO("cd failed " << tok);
549  d->cd();
550  } else {
551  break;
552  }
553  }
554  TObject* obj = d->FindObject(tok);
555  if (obj != NULL) {
556  if (obj->IsA()->InheritsFrom("TH1")) {
557  B2INFO("Histo " << histoname << " found in mem");
558  hh1 = (TH1*)obj;
559  }
560  } else {
561  B2INFO("Histo " << histoname << " NOT found in mem");
562  }
563  }
564  }
565 
566  if (hh1 == NULL) {
567  B2INFO("Histo " << histoname << " not found");
568  }
569 
570  return hh1;
571 }
The base class for the histogram analysis module.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
Definition: GeoCache.cc:175
const std::set< Belle2::VxdID > & getSensors(Belle2::VxdID ladder) const
Return a set of all sensor IDs belonging to a given ladder.
Definition: GeoCache.cc:203
const std::set< Belle2::VxdID > & getLadders(Belle2::VxdID layer) const
Return a set of all ladder IDs belonging to a given layer.
Definition: GeoCache.cc:192
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
std::map< Digit, DigitValue > Sensor
Map of all hits in one Sensor.
Abstract base class for different kinds of events.