Belle II Software  release-08-01-10
PXDClusterPositionEstimatorPar.h
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 #pragma once
9 
10 #include <TObject.h>
11 #include <TH2F.h>
12 #include <map>
13 
14 #include <pxd/dbobjects/PXDClusterShapeClassifierPar.h>
15 #include <pxd/dbobjects/PXDClusterOffsetPar.h>
16 
17 namespace Belle2 {
27  class PXDClusterPositionEstimatorPar: public TObject {
28  public:
31 
34 
36  void addGrid(int clusterkind, const TH2F& grid)
37  {
38  m_gridmap[clusterkind] = grid;
39  m_shapeClassifiers[clusterkind] = std::vector<PXDClusterShapeClassifierPar>();
40 
41  for (auto uBin = 1; uBin <= m_gridmap[clusterkind].GetXaxis()->GetNbins(); uBin++) {
42  for (auto vBin = 1; vBin <= m_gridmap[clusterkind].GetYaxis()->GetNbins(); vBin++) {
43  auto size = m_shapeClassifiers[clusterkind].size();
44  m_gridmap[clusterkind].SetBinContent(uBin, vBin, size);
45  m_shapeClassifiers[clusterkind].push_back(PXDClusterShapeClassifierPar());
46  }
47  }
48  }
49 
51  const std::map<int, TH2F>& getGridMap() const
52  {
53  return m_gridmap;
54  }
55 
57  void setShapeClassifier(const PXDClusterShapeClassifierPar& classifier, int uBin, int vBin, int clusterkind)
58  {
59  int key = m_gridmap[clusterkind].GetBinContent(uBin, vBin);
60  m_shapeClassifiers[clusterkind][key] = classifier;
61  }
62 
64  const PXDClusterShapeClassifierPar& getShapeClassifier(int uBin, int vBin, int clusterkind) const
65  {
66  int key = m_gridmap.at(clusterkind).GetBinContent(uBin, vBin);
67  const PXDClusterShapeClassifierPar& classifier = m_shapeClassifiers.at(clusterkind)[key];
68  return classifier;
69  }
70 
72  const PXDClusterShapeClassifierPar& getShapeClassifier(double thetaU, double thetaV, int clusterkind) const
73  {
74  auto grid = m_gridmap.at(clusterkind);
75  int uBin = grid.GetXaxis()->FindBin(thetaU);
76  int vBin = grid.GetYaxis()->FindBin(thetaV);
77  int key = grid.GetBinContent(uBin, vBin);
78  return m_shapeClassifiers.at(clusterkind)[key];
79  }
80 
82  bool hasClassifier(double thetaU, double thetaV, int clusterkind) const
83  {
84  //Check clusterkind is valid
85  if (m_gridmap.find(clusterkind) == m_gridmap.end()) {
86  return false;
87  }
88 
89  // Check thetaU, thetaV are inside grid
90  auto grid = m_gridmap.at(clusterkind);
91  int uBin = grid.GetXaxis()->FindBin(thetaU);
92  int vBin = grid.GetYaxis()->FindBin(thetaV);
93  int uBins = grid.GetXaxis()->GetNbins();
94  int vBins = grid.GetYaxis()->GetNbins();
95  if ((uBin < 1) || (uBin > uBins) || (vBin < 1) || (vBin > vBins))
96  return false;
97 
98  return true;
99  }
100 
102  const PXDClusterOffsetPar* getOffset(int shape_index, float eta, double thetaU, double thetaV, int clusterkind) const
103  {
104  //Check if there is a classifier
105  if (not hasClassifier(thetaU, thetaV, clusterkind)) {
106  return nullptr;
107  }
108  const PXDClusterShapeClassifierPar& classifier = getShapeClassifier(thetaU, thetaV, clusterkind);
109  return classifier.getOffset(shape_index, eta);
110  }
111 
113  float getShapeLikelyhood(int shape_index, double thetaU, double thetaV, int clusterkind) const
114  {
115  // Return zero if there no classifier
116  if (not hasClassifier(thetaU, thetaV, clusterkind)) {
117  return 0;
118  }
119  auto likelyhoodMap = getShapeClassifier(thetaU, thetaV, clusterkind).getShapeLikelyhoodMap();
120 
121  // Return zero if no likelyhood was estimated
122  auto it = likelyhoodMap.find(shape_index);
123  if (it == likelyhoodMap.end())
124  return 0;
125 
126  return it->second;
127  }
128 
129  private:
130 
132  std::map<int, TH2F> m_gridmap;
134  std::map<int, std::vector<PXDClusterShapeClassifierPar> > m_shapeClassifiers;
135 
137  };
139 } // end of namespace Belle2
The class for PXD cluster position offset payload.
The class for PXD cluster position lookup table payload.
void setShapeClassifier(const PXDClusterShapeClassifierPar &classifier, int uBin, int vBin, int clusterkind)
Set shape classifier.
ClassDef(PXDClusterPositionEstimatorPar, 2)
ClassDef, must be the last term before the closing {}.
float getShapeLikelyhood(int shape_index, double thetaU, double thetaV, int clusterkind) const
Returns shape likelyhood.
std::map< int, TH2F > m_gridmap
Map of angular grids for different clusterkinds
const std::map< int, TH2F > & getGridMap() const
Return grid map.
const PXDClusterShapeClassifierPar & getShapeClassifier(int uBin, int vBin, int clusterkind) const
Returns shape classifier.
const PXDClusterOffsetPar * getOffset(int shape_index, float eta, double thetaU, double thetaV, int clusterkind) const
Returns correction (offset) for cluster shape relative to center of pixel (startU/startV) if availabl...
bool hasClassifier(double thetaU, double thetaV, int clusterkind) const
Returns True if there is a classifier available.
void addGrid(int clusterkind, const TH2F &grid)
Add grid for clusterkind.
std::map< int, std::vector< PXDClusterShapeClassifierPar > > m_shapeClassifiers
Map of cluster shape classifiers for different clusterkinds.
const PXDClusterShapeClassifierPar & getShapeClassifier(double thetaU, double thetaV, int clusterkind) const
Returns shape classifier for incidence angles and clusterkind.
The class for PXD cluster shape classifier payload.
const std::map< int, float > & getShapeLikelyhoodMap() const
Return shape likelyhood map
Abstract base class for different kinds of events.