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