Belle II Software  release-08-01-10
SVDChannelMappingModule.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 <svd/modules/svdCalibration/SVDChannelMappingModule.h>
10 #include <vxd/geometry/GeoCache.h>
11 
12 #include <TH2F.h>
13 
14 using namespace std;
15 using namespace Belle2;
16 
17 //-----------------------------------------------------------------
18 // Register the Module
19 //-----------------------------------------------------------------
20 REG_MODULE(SVDChannelMapping);
21 
22 //-----------------------------------------------------------------
23 // Implementation
24 //-----------------------------------------------------------------
25 
26 SVDChannelMappingModule::SVDChannelMappingModule()
27  : Module()
28  , hInterDictionary(172, [](const Belle2::VxdID & vxdid) {return (size_t)vxdid.getID(); })
29 , n_events(0)
30 {
31  //Set module properties
32  setDescription("SVD Channel Mapping Verification Module");
34 
35  addParam("outputFileName", m_rootFileName, "Name of output root file.", std::string("SVDChannelMapping.root"));
36 
37  addParam("ShaperDigitsName", m_SVDShaperDigitsName,
38  "name of the list of SVDShaperDigits", std::string(""));
39 
40  addParam("ClustersName", m_SVDClustersName,
41  "name of the list of SVDClusters", std::string(""));
42 
43  addParam("InterceptsName", m_InterceptsName,
44  "name of the list of interceptions", std::string(""));
45 
46 }
47 
48 
49 
51 {
52 
53  m_histoList_digits = new TList;
54  m_histoList_clusters = new TList;
55 
57 
58  m_shapers.isRequired(m_SVDShaperDigitsName);
60  m_Intercepts.isRequired(m_InterceptsName);
61 
62  n_events = 0;
63 
64  m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
65 
66 }
67 
69 {
70 
71  n_events++;
72 
73 
74 
75  for (auto& it : m_Intercepts)
77 }
78 
80 {
81 
82  // VXD::GeoCache& aGeometry = VXD::GeoCache::getInstance();
83 
84  string name; //name of the histogram
85  string title; //title of the histogram
86  TH2F* tmp2D; //temporary 2D histo used to set axis title
87  TH1F* tmp1D; //temporary 1D histo used to set axis title
88 
89 
90  std::set<Belle2::VxdID> svdLayers = m_aGeometry.getLayers(VXD::SensorInfoBase::SVD);
91  std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
92 
93  while (itSvdLayers != svdLayers.end()) {
94 
95  std::set<Belle2::VxdID> svdLadders = m_aGeometry.getLadders(*itSvdLayers);
96  std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
97 
98  while (itSvdLadders != svdLadders.end()) {
99 
100  std::set<Belle2::VxdID> svdSensors = m_aGeometry.getSensors(*itSvdLadders);
101  std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
102 
103  while (itSvdSensors != svdSensors.end()) {
104 
105  string sensorid = std::to_string(itSvdSensors->getLayerNumber()) + "_" + std::to_string(itSvdSensors->getLadderNumber()) + "_" +
106  std::to_string(itSvdSensors->getSensorNumber());
107 
108 
109  // ------ HISTOGRAMS WITH A FILL PER INTERCEPT -------
110 
111  // coor U and V
112  name = "hCoorU_" + sensorid;
113  title = "U coordinate of the extrapolation in U for sensor " + sensorid;
114  tmp1D = new TH1F(name.c_str(), title.c_str(), 100, -2.9, 2.9);
115  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
116  (
117  (Belle2::VxdID)*itSvdSensors,
119  tmp1D,
120  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorU()); }
121  )
122  )
123  );
124  m_histoList_digits->Add(tmp1D);
125 
126  name = "hCoorV_" + sensorid;
127  title = "V coordinate of the extrapolation in V for sensor " + sensorid;
128  tmp1D = new TH1F(name.c_str(), title.c_str(), 100, -6.15, 6.15);
129  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
130  (
131  (Belle2::VxdID)*itSvdSensors,
133  tmp1D,
134  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorV()); }
135  )
136  )
137  );
138  m_histoList_digits->Add(tmp1D);
139 
140  // Intercept U vs V coordinate
141  name = "hCoorU_vs_CoorV_" + sensorid;
142  title = "U vs V intercept (cm) " + sensorid;
143  tmp2D = new TH2F(name.c_str(), title.c_str(), 100, -2.9, 2.9, 100, -6.15, 6.15);
144  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
145  tmp2D->GetYaxis()->SetTitle("intercept V coor (cm)");
146  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
147  (
148  (Belle2::VxdID)*itSvdSensors,
150  tmp2D,
151  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorU(), inter->getCoorV()); }
152  )
153  )
154  );
155  m_histoList_digits->Add(tmp2D);
156 
157  // sigma U and V
158  name = "hStatErrU_" + sensorid;
159  title = "stat error of the extrapolation in U for sensor " + sensorid;
160  tmp1D = new TH1F(name.c_str(), title.c_str(), 100, 0, 0.35);
161  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
162  (
163  (Belle2::VxdID)*itSvdSensors,
165  tmp1D,
166  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getSigmaU()); }
167  )
168  )
169  );
170  m_histoList_digits->Add(tmp1D);
171 
172  name = "hStatErrV_" + sensorid;
173  title = "stat error of the extrapolation in V for sensor " + sensorid;
174  tmp1D = new TH1F(name.c_str(), title.c_str(), 100, 0, 0.35);
175  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
176  (
177  (Belle2::VxdID)*itSvdSensors,
179  tmp1D,
180  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getSigmaV()); }
181  )
182  )
183  );
184  m_histoList_digits->Add(tmp1D);
185 
186  //---------------
187  //digits
188  //1D residuals
189  name = "hDigitResidU_" + sensorid;
190  title = "U residuals = intercept - digit, for sensor " + sensorid;
191  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -2.9, 2.9);
192  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
193  (
194  (Belle2::VxdID)*itSvdSensors,
196  tmp1D,
197  [this](TH1 * hPtr, const SVDIntercept * inter) {
198  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
199 
200  for (auto& it : SVDShaperDigits)
201  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
202  if (it.isUStrip()) {
203  const VXD::SensorInfoBase& aSensorInfo = m_aGeometry.getSensorInfo(it.getSensorID());
204  hPtr->Fill(inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID()));
205  }
206  }
207  }
208  )
209  )
210  );
211  m_histoList_digits->Add(tmp1D);
212 
213 
214  name = "hDigitResidV_" + sensorid;
215  title = "V residuals = intercept - digit, for sensor " + sensorid;
216  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -6.15, 6.15);
217  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
218  (
219  (Belle2::VxdID)*itSvdSensors,
221  tmp1D,
222  [this](TH1 * hPtr, const SVDIntercept * inter) {
223  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
224 
225  for (auto& it : SVDShaperDigits)
226  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
227  if (!it.isUStrip()) {
228  const VXD::SensorInfoBase& aSensorInfo = m_aGeometry.getSensorInfo(it.getSensorID());
229  hPtr->Fill(inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID()));
230  }
231  }
232  }
233  )
234  )
235  );
236  m_histoList_digits->Add(tmp1D);
237 
238  //residual U,V vs coordinate U,V
239  name = "hDigitResidU_vs_DigitU_" + sensorid;
240  title = "U residual (cm) vs digit U (cm) " + sensorid;
241  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -2.9, 2.9, 1000, -2.9, 2.9);
242  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
243  tmp2D->GetXaxis()->SetTitle("U digit (cm)");
244  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
245  (
246  (Belle2::VxdID)*itSvdSensors,
248  tmp2D,
249  [this](TH1 * hPtr, const SVDIntercept * inter) {
250  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
251 
252  for (auto& it : SVDShaperDigits)
253  if (((int)it.getSensorID() == (int)inter->getSensorID()) && it.isUStrip()) {
254  const VXD::SensorInfoBase& aSensorInfo = m_aGeometry.getSensorInfo(it.getSensorID());
255  double resid = inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID());
256  hPtr->Fill(aSensorInfo.getUCellPosition(it.getCellID()), resid);
257  }
258  }
259  )
260  )
261  );
262  m_histoList_digits->Add(tmp2D);
263 
264  name = "hDigitResidV_vs_DigitV_" + sensorid;
265  title = "V residual (cm) vs digit V (cm) " + sensorid;
266  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -6.15, 6.15, 1000, -6.15, 6.15);
267  tmp2D->GetYaxis()->SetTitle("V resid (cm)");
268  tmp2D->GetXaxis()->SetTitle("V digit (cm)");
269  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
270  (
271  (Belle2::VxdID)*itSvdSensors,
273  tmp2D,
274  [this](TH1 * hPtr, const SVDIntercept * inter) {
275  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
276 
277  for (auto& it : SVDShaperDigits)
278  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
279  const VXD::SensorInfoBase& aSensorInfo = m_aGeometry.getSensorInfo(it.getSensorID());
280  double resid = inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID());
281  hPtr->Fill(aSensorInfo.getVCellPosition(it.getCellID()), resid);
282  }
283  }
284  )
285  )
286  );
287  m_histoList_digits->Add(tmp2D);
288 
289 
290  // scatter plot: U,V intercept in cm VS U,V cell position
291  name = "hCoorU_vs_DigitU_" + sensorid;
292  title = "U intercept (cm) vs U digit (cm) " + sensorid;
293  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -2.9, 2.9, 1000, -2.9, 2.9);
294  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
295  tmp2D->GetYaxis()->SetTitle("digit U coor (cm)");
296  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
297  (
298  (Belle2::VxdID)*itSvdSensors,
300  tmp2D,
301  [this](TH1 * hPtr, const SVDIntercept * inter) {
302  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
303 
304  for (auto& it : SVDShaperDigits)
305  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUStrip())) {
306  const VXD::SensorInfoBase& aSensorInfo = m_aGeometry.getSensorInfo(it.getSensorID());
307  hPtr->Fill(inter->getCoorU(), aSensorInfo.getUCellPosition(it.getCellID()));
308  }
309  }
310  )
311  )
312  );
313  m_histoList_digits->Add(tmp2D);
314 
315  name = "hCoorV_vs_DigitV_" + sensorid;
316  title = "V intercept (cm) vs V digit (ID) " + sensorid;
317  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -6.15, 6.15, 1000, -6.15, 6.15);
318  tmp2D->GetXaxis()->SetTitle("intercept V coor (cm)");
319  tmp2D->GetYaxis()->SetTitle("digi V coor (cm)");
320  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
321  (
322  (Belle2::VxdID)*itSvdSensors,
324  tmp2D,
325  [this](TH1 * hPtr, const SVDIntercept * inter) {
326  StoreArray<SVDShaperDigit> SVDShaperDigits(this->m_SVDShaperDigitsName);
327 
328  for (auto& it : SVDShaperDigits) {
329  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
330  const VXD::SensorInfoBase& aSensorInfo = m_aGeometry.getSensorInfo(it.getSensorID());
331  hPtr->Fill(inter->getCoorV(), aSensorInfo.getVCellPosition(it.getCellID()));
332  }
333  }
334  }
335  )
336  )
337  );
338  m_histoList_digits->Add(tmp2D);
339 
340  //-----------------
341  //clusters
342  //1D residuals
343  name = "hClusterResidU_" + sensorid;
344  title = "U residuals = intercept - cluster, for sensor " + sensorid;
345  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -2.9, 2.9);
346  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
347  (
348  (Belle2::VxdID)*itSvdSensors,
350  tmp1D,
351  [this](TH1 * hPtr, const SVDIntercept * inter) {
352  StoreArray<SVDCluster> SVDClusters(this->m_SVDClustersName);
353 
354  for (auto& it : SVDClusters)
355  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
356  if (it.isUCluster()) {
357  hPtr->Fill(inter->getCoorU() - it.getPosition());
358  }
359  }
360  }
361  )
362  )
363  );
364  m_histoList_clusters->Add(tmp1D);
365 
366 
367  name = "hClusterResidV_" + sensorid;
368  title = "V residuals = intercept - cluster, for sensor " + sensorid;
369  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -6.15, 6.15);
370  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
371  (
372  (Belle2::VxdID)*itSvdSensors,
374  tmp1D,
375  [this](TH1 * hPtr, const SVDIntercept * inter) {
376  StoreArray<SVDCluster> SVDClusters(this->m_SVDClustersName);
377 
378  for (auto& it : SVDClusters)
379  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
380  if (!it.isUCluster()) {
381  hPtr->Fill(inter->getCoorV() - it.getPosition());
382  }
383  }
384  }
385  )
386  )
387  );
388  m_histoList_clusters->Add(tmp1D);
389 
390  //residual U,V vs coordinate U,V
391  name = "hClusterResidU_vs_ClusterU_" + sensorid;
392  title = "U residual (cm) vs cluster U (cm) " + sensorid;
393  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -2.9, 2.9, 1000, -2.9, 2.9);
394  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
395  tmp2D->GetXaxis()->SetTitle("U cluster (cm)");
396  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
397  (
398  (Belle2::VxdID)*itSvdSensors,
400  tmp2D,
401  [this](TH1 * hPtr, const SVDIntercept * inter) {
402  StoreArray<SVDCluster> SVDClusters(this->m_SVDClustersName);
403 
404  for (auto& it : SVDClusters)
405  if (((int)it.getSensorID() == (int)inter->getSensorID()) && it.isUCluster()) {
406  double resid = inter->getCoorU() - it.getPosition();
407  hPtr->Fill(it.getPosition(), resid);
408  }
409  }
410  )
411  )
412  );
413  m_histoList_clusters->Add(tmp2D);
414 
415  name = "hClusterResidV_vs_ClusterV_" + sensorid;
416  title = "V residual (cm) vs cluster V (cm) " + sensorid;
417  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -6.15, 6.15, 1000, -6.15, 6.15);
418  tmp2D->GetYaxis()->SetTitle("V resid (cm)");
419  tmp2D->GetXaxis()->SetTitle("V cluster (cm)");
420  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
421  (
422  (Belle2::VxdID)*itSvdSensors,
424  tmp2D,
425  [this](TH1 * hPtr, const SVDIntercept * inter) {
426  StoreArray<SVDCluster> SVDClusters(this->m_SVDClustersName);
427 
428  for (auto& it : SVDClusters)
429  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUCluster())) {
430  double resid = inter->getCoorV() - it.getPosition();
431  hPtr->Fill(it.getPosition(), resid);
432  }
433  }
434  )
435  )
436  );
437  m_histoList_clusters->Add(tmp2D);
438 
439 
440  // scatter plot: U,V intercept in cm VS U,V cell position
441  name = "hCoorU_vs_ClusterU_" + sensorid;
442  title = "U intercept (cm) vs U cluster (cm) " + sensorid;
443  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -2.9, 2.9, 1000, -2.9, 2.9);
444  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
445  tmp2D->GetYaxis()->SetTitle("cluster U coor (cm)");
446  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
447  (
448  (Belle2::VxdID)*itSvdSensors,
450  tmp2D,
451  [this](TH1 * hPtr, const SVDIntercept * inter) {
452  StoreArray<SVDCluster> SVDClusters(this->m_SVDClustersName);
453 
454  for (auto& it : SVDClusters)
455  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUCluster())) {
456  hPtr->Fill(inter->getCoorU(), it.getPosition());
457  }
458  }
459  )
460  )
461  );
462  m_histoList_clusters->Add(tmp2D);
463 
464  name = "hCoorV_vs_ClusterV_" + sensorid;
465  title = "V intercept (cm) vs V cluster (ID) " + sensorid;
466  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -6.15, 6.15, 1000, -6.15, 6.15);
467  tmp2D->GetXaxis()->SetTitle("intercept V coor (cm)");
468  tmp2D->GetYaxis()->SetTitle("cluster V coor (cm)");
469  hInterDictionary.insert(pair< Belle2::VxdID, InterHistoAndFill >
470  (
471  (Belle2::VxdID)*itSvdSensors,
473  tmp2D,
474  [this](TH1 * hPtr, const SVDIntercept * inter) {
475  StoreArray<SVDCluster> SVDClusters(this->m_SVDClustersName);
476 
477  for (auto& it : SVDClusters) {
478  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUCluster())) {
479  hPtr->Fill(inter->getCoorV(), it.getPosition());
480  }
481  }
482  }
483  )
484  )
485  );
486  m_histoList_clusters->Add(tmp2D);
487 
488  ++itSvdSensors;
489  }
490  ++itSvdLadders;
491  }
492  ++itSvdLayers;
493  }
494 
495 }
496 
498 {
499 
500  auto its = hInterDictionary.equal_range(inter->getSensorID());
501 
502  for (auto it = its.first; it != its.second; ++it) {
503  InterHistoAndFill aInterHistoAndFill = it->second;
504  aInterHistoAndFill.second(aInterHistoAndFill.first, inter);
505  }
506 
507 }
508 
510 {
511 
512  if (m_rootFilePtr != nullptr) {
513  m_rootFilePtr->cd();
514 
515 
516  TDirectory* oldDir = gDirectory;
517  TObject* obj;
518 
519  TDirectory* dir_digits = oldDir->mkdir("digits");
520  dir_digits->cd();
521  TIter nextH_digits(m_histoList_digits);
522  while ((obj = nextH_digits()))
523  obj->Write();
524 
525  TDirectory* dir_clusters = oldDir->mkdir("clusters");
526  dir_clusters->cd();
527  TIter nextH_clusters(m_histoList_clusters);
528  while ((obj = nextH_clusters()))
529  obj->Write();
530 
531  }
532 }
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
std::string m_rootFileName
root file name
std::unordered_multimap< Belle2::VxdID, InterHistoAndFill, std::function< size_t(const Belle2::VxdID &)> > hInterDictionary
map of histograms to be filled once per intercept
std::pair< TH1 *, std::function< void(TH1 *, const SVDIntercept *) > > InterHistoAndFill
typedef: histograms to be filled once per intercept + filling function
std::string m_InterceptsName
Name of the SVDIntercept StoreArray.
VXD::GeoCache & m_aGeometry
the geometry
void initialize() override
initialize
void fillSensorInterHistos(const SVDIntercept *inter)
fill histograms per sensor, filled once per intercept
void createHistosDictionaries()
create the dictionary
std::string m_SVDClustersName
reco digit list name
StoreArray< SVDShaperDigit > m_shapers
SVDSHaperDigits StoreArray.
TList * m_histoList_clusters
list of histograms for the clusters
StoreArray< SVDCluster > m_clusters
SVDClusters StoreArray.
TList * m_histoList_digits
list of histograms for the strips
StoreArray< SVDIntercept > m_Intercepts
SVDintercept StoreArray.
TFile * m_rootFilePtr
pointer at root file used for storing histograms
std::string m_SVDShaperDigitsName
shaper digit list name
SVDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an SVD ...
Definition: SVDIntercept.h:22
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
double getSigmaV() const
return the statistical error on the V coordinate of the intercept
Definition: VXDIntercept.h:62
double getCoorV() const
return the V coordinate of the intercept
Definition: VXDIntercept.h:60
double getSigmaU() const
return the statistical error on the U coordinate of the intercept
Definition: VXDIntercept.h:61
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
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
Definition: GeoCache.cc:176
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
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:204
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:193
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.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#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.