Belle II Software  release-06-01-15
SVDROIDQMModule.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 <tracking/modules/svdROIFinder/SVDROIDQMModule.h>
10 #include <vxd/geometry/GeoCache.h>
11 
12 #include <TDirectory.h>
13 #include <TH2F.h>
14 
15 using namespace std;
16 using namespace Belle2;
17 
18 //-----------------------------------------------------------------
19 // Register the Module
20 //-----------------------------------------------------------------
21 REG_MODULE(SVDROIDQM)
22 
23 //-----------------------------------------------------------------
24 // Implementation
25 //-----------------------------------------------------------------
26 
28  : HistoModule()
29  , m_InterDir(nullptr)
30  , m_ROIDir(nullptr)
31  , hInterDictionary(172, [](const Belle2::VxdID & vxdid) {return (size_t)vxdid.getID(); })
32 , hROIDictionary(172, [](const Belle2::VxdID& vxdid) {return (size_t)vxdid.getID(); })
33 , hROIDictionaryEvt(172, [](const Belle2::VxdID& vxdid) {return (size_t)vxdid.getID(); })
34 , m_numModules(0)
35 , hnROIs(nullptr)
36 , hnInter(nullptr)
37 , harea(nullptr)
38 , hredFactor(nullptr)
39 , hCellU(nullptr)
40 , hCellV(nullptr)
41 , n_events(0)
42 {
43  //Set module properties
44  setDescription("Monitor of the ROIs creation on HLT");
45  setPropertyFlags(c_ParallelProcessingCertified);
46 
47  addParam("SVDShaperDigitsName", m_SVDShaperDigitsName,
48  "name of the list of SVDShaperDigits", std::string(""));
49  addParam("SVDRecoDigitsName", m_SVDRecoDigitsName,
50  "name of the list of SVDRecoDigits", std::string(""));
51 
52  addParam("InterceptsName", m_InterceptsName,
53  "name of the list of interceptions", std::string(""));
54 
55  addParam("ROIsName", m_ROIsName,
56  "name of the list of ROIs", std::string(""));
57 
58 }
59 
60 void SVDROIDQMModule::defineHisto()
61 {
62 
63  // Create a separate histogram directory and cd into it.
64  TDirectory* oldDir = gDirectory;
65  TDirectory* roiDir = oldDir->mkdir("SVDROIs");
66  m_InterDir = roiDir->mkdir("intercept");
67  m_ROIDir = roiDir->mkdir("roi");
68 
69  hCellU = new TH1F("hCellU", "CellID U", 769, -0.5, 768.5);
70  hCellU->GetXaxis()->SetTitle("U cell ID");
71  hCellV = new TH1F("hCellV", "CellID V", 769, -0.5, 768.5);
72  hCellV->GetXaxis()->SetTitle("V cell ID");
73 
74  m_InterDir->cd();
75  hnInter = new TH1F("hnInter", "number of intercepts", 100, 0, 100);
76 
77  m_ROIDir->cd();
78  hnROIs = new TH1F("hnROIs", "number of ROIs", 100, 0, 100);
79  harea = new TH1F("harea", "ROIs area", 100, 0, 100000);
80  hredFactor = new TH1F("hredFactor", "ROI reduction factor", 1000, 0, 1);
81 
82 
83  createHistosDictionaries();
84 
85  oldDir->cd();
86 
87 }
88 
89 void SVDROIDQMModule::initialize()
90 {
91  REG_HISTOGRAM
92 
93  m_SVDShaperDigits.isOptional(m_SVDShaperDigitsName);
94  m_SVDRecoDigits.isOptional(m_SVDRecoDigitsName);
95  m_ROIs.isRequired(m_ROIsName);
96  m_Intercepts.isRequired(m_InterceptsName);
97 
98  n_events = 0;
99 
100 
101 }
102 
103 void SVDROIDQMModule::event()
104 {
105 
106  n_events++;
107 
108  for (auto& it : m_SVDShaperDigits)
109  if (it.isUStrip())
110  hCellU->Fill(it.getCellID());
111  else
112  hCellV->Fill(it.getCellID());
113 
114  hnInter->Fill(m_Intercepts.getEntries());
115 
116  for (auto& it : m_Intercepts)
117  fillSensorInterHistos(&it);
118 
119 
120  for (auto it = hROIDictionaryEvt.begin(); it != hROIDictionaryEvt.end(); ++it)
121  (it->second).value = 0;
122 
123  int ROIarea = 0;
124  double redFactor = 0;
125 
126  for (auto& it : m_ROIs) {
127 
128  fillSensorROIHistos(&it);
129 
130  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
131  const int nPixelsU = aSensorInfo.getUCells();
132  const int nPixelsV = aSensorInfo.getVCells();
133 
134  int minU = it.getMinUid();
135  int minV = it.getMinVid();
136  int maxU = it.getMaxUid();
137  int maxV = it.getMaxVid();
138 
139  int tmpROIarea = (maxU - minU) * (maxV - minV);
140  ROIarea += tmpROIarea;
141  redFactor += (double)tmpROIarea / (nPixelsU * nPixelsV * m_numModules);
142 
143  }
144 
145  hnROIs->Fill(m_ROIs.getEntries());
146 
147  harea->Fill((double)ROIarea);
148 
149  hredFactor->Fill((double)redFactor);
150 
151 
152  for (auto it = hROIDictionaryEvt.begin(); it != hROIDictionaryEvt.end(); ++it) {
153  ROIHistoAccumulateAndFill aROIHistoAccumulateAndFill = it->second;
154  aROIHistoAccumulateAndFill.fill(aROIHistoAccumulateAndFill.hPtr, aROIHistoAccumulateAndFill.value);
155  }
156 
157 }
158 
159 void SVDROIDQMModule::createHistosDictionaries()
160 {
161 
162  // VXD::GeoCache& aGeometry = VXD::GeoCache::getInstance();
163 
164  string name; //name of the histogram
165  string title; //title of the histogram
166  TH2F* tmp2D; //temporary 2D histo used to set axis title
167  TH1F* tmp1D; //temporary 1D histo used to set axis title
168 
169  m_numModules = 0;
170 
171  std::set<Belle2::VxdID> svdLayers = m_geoCache.getLayers(VXD::SensorInfoBase::SVD);
172  std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
173 
174  while (itSvdLayers != svdLayers.end()) {
175 
176  std::set<Belle2::VxdID> svdLadders = m_geoCache.getLadders(*itSvdLayers);
177  std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
178 
179  while (itSvdLadders != svdLadders.end()) {
180 
181  std::set<Belle2::VxdID> svdSensors = m_geoCache.getSensors(*itSvdLadders);
182  std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
183 
184  while (itSvdSensors != svdSensors.end()) {
185 
186  m_numModules++; //counting the total number of modules
187 
188  const VXD::SensorInfoBase& wSensorInfo = m_geoCache.getSensorInfo(*itSvdSensors);
189 
190  const int nPixelsU = wSensorInfo.getUCells();
191  const int nPixelsV = wSensorInfo.getVCells();
192  string sensorid = std::to_string(itSvdSensors->getLayerNumber()) + "_" + std::to_string(itSvdSensors->getLadderNumber()) + "_" +
193  std::to_string(itSvdSensors->getSensorNumber());
194 
195 
196  // ------ HISTOGRAMS WITH AN ACCUMULATE PER ROI AND A FILL PER EVENT -------
197  m_ROIDir->cd();
198 
199  name = "hNROIs_" + sensorid;
200  title = "number of ROIs for sensor " + sensorid;
201  double value = 0;
203  new TH1F(name.c_str(), title.c_str(), 25, 0, 25),
204  [](const ROIid*, double & val) {val++;},
205  [](TH1 * hPtr, double & val) { hPtr->Fill(val); },
206  value
207  };
208  hROIDictionaryEvt.insert(pair< Belle2::VxdID, ROIHistoAccumulateAndFill& > ((Belle2::VxdID)*itSvdSensors, *aHAAF));
209 
210 
211 
212 
213  // ------ HISTOGRAMS WITH A FILL PER INTERCEPT -------
214  m_InterDir->cd();
215 
216  // coor U and V
217  name = "hCoorU_" + sensorid;
218  title = "U coordinate of the extrapolation in U for sensor " + sensorid;
219  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
220  (
221  (Belle2::VxdID)*itSvdSensors,
223  new TH1F(name.c_str(), title.c_str(), 100, -5, 5),
224  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorU()); }
225  )
226  )
227  );
228 
229  name = "hCoorV_" + sensorid;
230  title = "V coordinate of the extrapolation in V for sensor " + sensorid;
231  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
232  (
233  (Belle2::VxdID)*itSvdSensors,
235  new TH1F(name.c_str(), title.c_str(), 100, -5, 5),
236  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorV()); }
237  )
238  )
239  );
240 
241  // Intercept U vs V coordinate
242  name = "hCoorU_vs_CoorV_" + sensorid;
243  title = "U vs V intercept (cm) " + sensorid;
244  tmp2D = new TH2F(name.c_str(), title.c_str(), 100, -5, 5, 100, -5, 5);
245  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
246  tmp2D->GetYaxis()->SetTitle("intercept V coor (cm)");
247  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
248  (
249  (Belle2::VxdID)*itSvdSensors,
251  tmp2D,
252  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorU(), inter->getCoorV()); }
253  )
254  )
255  );
256 
257 
258  // sigma U and V
259  name = "hStatErrU_" + sensorid;
260  title = "stat error of the extrapolation in U for sensor " + sensorid;
261  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
262  (
263  (Belle2::VxdID)*itSvdSensors,
265  new TH1F(name.c_str(), title.c_str(), 100, 0, 0.35),
266  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getSigmaU()); }
267  )
268  )
269  );
270  name = "hStatErrV_" + sensorid;
271  title = "stat error of the extrapolation in V for sensor " + sensorid;
272  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
273  (
274  (Belle2::VxdID)*itSvdSensors,
276  new TH1F(name.c_str(), title.c_str(), 100, 0, 0.35),
277  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getSigmaV()); }
278  )
279  )
280  );
281 
282  //1D residuals
283  name = "hResidU_" + sensorid;
284  title = "U residuals = intercept - digit, for sensor " + sensorid;
285  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
286  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
287  (
288  (Belle2::VxdID)*itSvdSensors,
290  tmp1D,
291  [this](TH1 * hPtr, const SVDIntercept * inter) {
292  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
293 
294  for (auto& it : SVDShaperDigits)
295  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
296  if (it.isUStrip()) {
297  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
298  hPtr->Fill(inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID()));
299  }
300  }
301  }
302  )
303  )
304  );
305 
306  name = "hResidV_" + sensorid;
307  title = "V residuals = intercept - digit, for sensor " + sensorid;
308  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
309  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
310  (
311  (Belle2::VxdID)*itSvdSensors,
313  tmp1D,
314  [this](TH1 * hPtr, const SVDIntercept * inter) {
315  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
316 
317  for (auto& it : SVDShaperDigits)
318  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
319  if (!it.isUStrip()) {
320  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
321  hPtr->Fill(inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID()));
322  }
323  }
324  }
325  )
326  )
327  );
328 
329  //residual U,V vs coordinate U,V
330  name = "hResidU_vs_CoorU_" + sensorid;
331  title = "U residual (cm) vs coor U (cm) " + sensorid;
332  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
333  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
334  tmp2D->GetXaxis()->SetTitle("U coor (cm)");
335  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
336  (
337  (Belle2::VxdID)*itSvdSensors,
339  tmp2D,
340  [this](TH1 * hPtr, const SVDIntercept * inter) {
341  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
342 
343  for (auto& it : SVDShaperDigits)
344  if (((int)it.getSensorID() == (int)inter->getSensorID()) && it.isUStrip()) {
345  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
346  double resid = inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID());
347  hPtr->Fill(inter->getCoorU(), resid);
348  }
349  }
350  )
351  )
352  );
353 
354  name = "hResidV_vs_CoorV_" + sensorid;
355  title = "V residual (cm) vs coor V (cm) " + sensorid;
356  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
357  tmp2D->GetYaxis()->SetTitle("V resid (cm)");
358  tmp2D->GetXaxis()->SetTitle("V coor (cm)");
359  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
360  (
361  (Belle2::VxdID)*itSvdSensors,
363  tmp2D,
364  [this](TH1 * hPtr, const SVDIntercept * inter) {
365  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
366 
367  for (auto& it : SVDShaperDigits)
368  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
369  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
370  double resid = inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID());
371  hPtr->Fill(inter->getCoorV(), resid);
372  }
373  }
374  )
375  )
376  );
377 
378 
379  //residual vs charge
380  name = "hResidU_vs_charge_" + sensorid;
381  title = "U residual (cm) vs charge " + sensorid;
382  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
383  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
384  tmp2D->GetXaxis()->SetTitle("charge");
385  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
386  (
387  (Belle2::VxdID)*itSvdSensors,
389  tmp2D,
390  [this](TH1 * hPtr, const SVDIntercept * inter) {
391  StoreArray<SVDRecoDigit> SVDRecoDigits(this->m_SVDRecoDigitsName);
392 
393  for (auto& it : SVDRecoDigits)
394  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUStrip())) {
395  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
396  double resid = inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID());
397  hPtr->Fill(it.getCharge(), resid);
398  }
399  }
400  )
401  )
402  );
403 
404  name = "hResidV_vs_charge_" + sensorid;
405  title = "V residual (cm) vs charge " + sensorid;
406  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
407  tmp2D->GetYaxis()->SetTitle("V resid (cm)");
408  tmp2D->GetXaxis()->SetTitle("charge");
409  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
410  (
411  (Belle2::VxdID)*itSvdSensors,
413  tmp2D,
414  [this](TH1 * hPtr, const SVDIntercept * inter) {
415  StoreArray<SVDRecoDigit> SVDRecoDigits(this->m_SVDRecoDigitsName);
416 
417  for (auto& it : SVDRecoDigits)
418  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
419  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
420  double resid = inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID());
421  hPtr->Fill(it.getCharge(), resid);
422  }
423  }
424  )
425  )
426  );
427 
428 
429  // scatter plot: U,V intercept in cm VS U,V cell position
430  name = "hCoorU_vs_UDigit_" + sensorid;
431  title = "U intercept (cm) vs U Digit (ID) " + sensorid;
432  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
433  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
434  tmp2D->GetYaxis()->SetTitle("digit U coor (cm)");
435  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
436  (
437  (Belle2::VxdID)*itSvdSensors,
439  tmp2D,
440  [this](TH1 * hPtr, const SVDIntercept * inter) {
441  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
442 
443  for (auto& it : SVDShaperDigits)
444  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUStrip())) {
445  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
446  hPtr->Fill(inter->getCoorU(), aSensorInfo.getUCellPosition(it.getCellID()));
447  // hPtr->Fill( inter->getCoorU(), it.getVCellID()*75e-4 );
448  }
449  }
450  )
451  )
452  );
453 
454  name = "hCoorV_vs_VDigit_" + sensorid;
455  title = "V intercept (cm) vs V Digit (ID) " + sensorid;
456  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
457  tmp2D->GetXaxis()->SetTitle("intercept V coor (cm)");
458  tmp2D->GetYaxis()->SetTitle("digi V coor (cm)");
459  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
460  (
461  (Belle2::VxdID)*itSvdSensors,
463  tmp2D,
464  [this](TH1 * hPtr, const SVDIntercept * inter) {
465  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
466 
467  for (auto& it : SVDShaperDigits) {
468  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
469  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
470  hPtr->Fill(inter->getCoorV(), aSensorInfo.getVCellPosition(it.getCellID()));
471  // hPtr->Fill( inter->getCoorV(), it.getUCellID()*50e-4 );
472  }
473  }
474  }
475  )
476  )
477  );
478 
479 
480 
481 
482 
483 
484  // ------ HISTOGRAMS WITH A FILL PER ROI -------
485  m_ROIDir->cd();
486 
487  // MIN in U and V
488  name = "hminU_" + sensorid;
489  title = "ROI min in U for sensor " + sensorid;
490  hROIDictionary.insert(pair< Belle2::VxdID, ROIHistoAndFill >
491  (
492  (Belle2::VxdID)*itSvdSensors,
494  new TH1F(name.c_str(), title.c_str(), nPixelsU, 0, nPixelsU),
495  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMinUid()); }
496  )
497  )
498  );
499  name = "hminV_" + sensorid;
500  title = "ROI min in V for sensor " + sensorid;
501  hROIDictionary.insert(pair< Belle2::VxdID, ROIHistoAndFill >
502  (
503  (Belle2::VxdID)*itSvdSensors,
505  new TH1F(name.c_str(), title.c_str(), nPixelsV, 0, nPixelsV),
506  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMinVid()); }
507  )
508  )
509  );
510  //--------------------------
511  // MAX in U and V
512  name = "hmaxU_" + sensorid;
513  title = "ROI max in U for sensor " + sensorid;
514  hROIDictionary.insert(pair< Belle2::VxdID, ROIHistoAndFill >
515  (
516  (Belle2::VxdID)*itSvdSensors,
518  new TH1F(name.c_str(), title.c_str(), nPixelsU, 0, nPixelsU),
519  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMaxUid()); }
520  )
521  )
522  );
523  name = "hmaxV_" + sensorid;
524  title = "ROI max in V for sensor " + sensorid;
525  hROIDictionary.insert(pair< Belle2::VxdID, ROIHistoAndFill >
526  (
527  (Belle2::VxdID)*itSvdSensors,
529  new TH1F(name.c_str(), title.c_str(), nPixelsV, 0, nPixelsV),
530  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMaxVid()); }
531  )
532  )
533  );
534  //--------------------------
535 
536  // WIDTH in U and V
537  name = "hwidthU_" + sensorid;
538  title = "ROI width in U for sensor " + sensorid;
539  hROIDictionary.insert(pair< Belle2::VxdID, ROIHistoAndFill >
540  (
541  (Belle2::VxdID)*itSvdSensors,
543  new TH1F(name.c_str(), title.c_str(), nPixelsU, 0, nPixelsU),
544  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMaxUid() - roi->getMinUid()); }
545  )
546  )
547  );
548  name = "hwidthV_" + sensorid;
549  title = "ROI width in V for sensor " + sensorid;
550  hROIDictionary.insert(pair< Belle2::VxdID, ROIHistoAndFill >
551  (
552  (Belle2::VxdID)*itSvdSensors,
554  new TH1F(name.c_str(), title.c_str(), nPixelsV, 0, nPixelsV),
555  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMaxVid() - roi->getMinVid()); }
556  )
557  )
558  );
559 
560  // ROI center
561  name = "hROIcenter_" + sensorid;
562  title = "ROI center " + sensorid;
563  tmp2D = new TH2F(name.c_str(), title.c_str(), nPixelsU, 0, nPixelsU, nPixelsV, 0, nPixelsV);
564  tmp2D->GetXaxis()->SetTitle(" U (ID)");
565  tmp2D->GetYaxis()->SetTitle(" V (ID)");
566  hROIDictionary.insert(pair< Belle2::VxdID, ROIHistoAndFill >
567  (
568  (Belle2::VxdID)*itSvdSensors,
570  tmp2D,
571  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill((roi->getMaxUid() + roi->getMinUid()) / 2, (roi->getMaxVid() + roi->getMinVid()) / 2); }
572  )
573  )
574  );
575 
576  //--------------------------
577 
578  ++itSvdSensors;
579  }
580  ++itSvdLadders;
581  }
582  ++itSvdLayers;
583  }
584 
585 }
586 
587 void SVDROIDQMModule::fillSensorInterHistos(const SVDIntercept* inter)
588 {
589 
590  auto its = hInterDictionary.equal_range(inter->getSensorID());
591 
592  for (auto it = its.first; it != its.second; ++it) {
593  InterHistoAndFill aInterHistoAndFill = it->second;
594  aInterHistoAndFill.second(aInterHistoAndFill.first, inter);
595  }
596 
597 }
598 
599 void SVDROIDQMModule::fillSensorROIHistos(const ROIid* roi)
600 {
601 
602  auto its = hROIDictionary.equal_range(roi->getSensorID());
603 
604  for (auto it = its.first; it != its.second; ++it) {
605  ROIHistoAndFill aROIHistoAndFill = it->second;
606  aROIHistoAndFill.second(aROIHistoAndFill.first, roi);
607  }
608 
609  auto itsEvt = hROIDictionaryEvt.equal_range(roi->getSensorID());
610  for (auto it = itsEvt.first; it != itsEvt.second; ++it)
611  (it->second).accumulate(roi, (it->second).value);
612 }
613 
614 void SVDROIDQMModule::endRun()
615 {
616 
617  hCellU->Scale((double)1 / n_events);
618  hCellV->Scale((double)1 / n_events);
619 
620  for (auto it = hROIDictionaryEvt.begin(); it != hROIDictionaryEvt.end(); ++it)
621  delete &(it->second);
622 }
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
ROIid stores the U and V ids and the sensor id of the Region Of Interest.
Definition: ROIid.h:25
SVDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an SVD ...
Definition: SVDIntercept.h:22
Creates basic DQM for ROI creation on ExpressReco
std::pair< TH1 *, std::function< void(TH1 *, const SVDIntercept *) > > InterHistoAndFill
typedef: histograms to be filled once per intercept + filling function
std::pair< TH1 *, std::function< void(TH1 *, const ROIid *) > > ROIHistoAndFill
typedef: histograms to be filled once per roi + filling function
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
double getCoorV() const
return the V coordinate of the intercept
Definition: VXDIntercept.h:60
VxdID::baseType getSensorID() const
return the sensor ID
Definition: VXDIntercept.h:66
double getCoorU() const
return the U coordinate of the intercept
Definition: VXDIntercept.h:59
Base class to provide Sensor Information for PXD and SVD.
double getVCellPosition(int vID) const
Return the position of a specific strip/pixel in v direction.
double getUCellPosition(int uID, int vID=-1) const
Return the position of a specific strip/pixel in u direction.
int getVCells() const
Return number of pixel/strips in v direction.
int getUCells() const
Return number of pixel/strips in u direction.
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
Abstract base class for different kinds of events.
struct: histograms to be filled once per event + filling fucntion + accumulate function
std::function< void(TH1 *, double &) > fill
fill function