Belle II Software  release-05-01-25
PXDClusterPositionEstimatorPar.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Benjamin Schwenker *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 #pragma once
11 
12 #include <TObject.h>
13 #include <TH2F.h>
14 #include <map>
15 
16 #include <pxd/dbobjects/PXDClusterShapeClassifierPar.h>
17 #include <pxd/dbobjects/PXDClusterOffsetPar.h>
18 
19 namespace Belle2 {
29  class PXDClusterPositionEstimatorPar: public TObject {
30  public:
33 
36 
38  void addGrid(int clusterkind, const TH2F& grid)
39  {
40  m_gridmap[clusterkind] = grid;
41  m_shapeClassifiers[clusterkind] = std::vector<PXDClusterShapeClassifierPar>();
42 
43  for (auto uBin = 1; uBin <= m_gridmap[clusterkind].GetXaxis()->GetNbins(); uBin++) {
44  for (auto vBin = 1; vBin <= m_gridmap[clusterkind].GetYaxis()->GetNbins(); vBin++) {
45  auto size = m_shapeClassifiers[clusterkind].size();
46  m_gridmap[clusterkind].SetBinContent(uBin, vBin, size);
47  m_shapeClassifiers[clusterkind].push_back(PXDClusterShapeClassifierPar());
48  }
49  }
50  }
51 
53  const std::map<int, TH2F>& getGridMap() const
54  {
55  return m_gridmap;
56  }
57 
59  void setShapeClassifier(const PXDClusterShapeClassifierPar& classifier, int uBin, int vBin, int clusterkind)
60  {
61  int key = m_gridmap[clusterkind].GetBinContent(uBin, vBin);
62  m_shapeClassifiers[clusterkind][key] = classifier;
63  }
64 
66  const PXDClusterShapeClassifierPar& getShapeClassifier(int uBin, int vBin, int clusterkind) const
67  {
68  int key = m_gridmap.at(clusterkind).GetBinContent(uBin, vBin);
69  const PXDClusterShapeClassifierPar& classifier = m_shapeClassifiers.at(clusterkind)[key];
70  return classifier;
71  }
72 
74  const PXDClusterShapeClassifierPar& getShapeClassifier(double thetaU, double thetaV, int clusterkind) const
75  {
76  auto grid = m_gridmap.at(clusterkind);
77  int uBin = grid.GetXaxis()->FindBin(thetaU);
78  int vBin = grid.GetYaxis()->FindBin(thetaV);
79  int key = grid.GetBinContent(uBin, vBin);
80  return m_shapeClassifiers.at(clusterkind)[key];
81  }
82 
84  bool hasClassifier(double thetaU, double thetaV, int clusterkind) const
85  {
86  //Check clusterkind is valid
87  if (m_gridmap.find(clusterkind) == m_gridmap.end()) {
88  return false;
89  }
90 
91  // Check thetaU, thetaV are inside grid
92  auto grid = m_gridmap.at(clusterkind);
93  int uBin = grid.GetXaxis()->FindBin(thetaU);
94  int vBin = grid.GetYaxis()->FindBin(thetaV);
95  int uBins = grid.GetXaxis()->GetNbins();
96  int vBins = grid.GetYaxis()->GetNbins();
97  if ((uBin < 1) || (uBin > uBins) || (vBin < 1) || (vBin > vBins))
98  return false;
99 
100  return true;
101  }
102 
104  const PXDClusterOffsetPar* getOffset(int shape_index, float eta, double thetaU, double thetaV, int clusterkind) const
105  {
106  //Check if there is a classifier
107  if (not hasClassifier(thetaU, thetaV, clusterkind)) {
108  return nullptr;
109  }
110  const PXDClusterShapeClassifierPar& classifier = getShapeClassifier(thetaU, thetaV, clusterkind);
111  return classifier.getOffset(shape_index, eta);
112  }
113 
115  float getShapeLikelyhood(int shape_index, double thetaU, double thetaV, int clusterkind) const
116  {
117  // Return zero if there no classifier
118  if (not hasClassifier(thetaU, thetaV, clusterkind)) {
119  return 0;
120  }
121  auto likelyhoodMap = getShapeClassifier(thetaU, thetaV, clusterkind).getShapeLikelyhoodMap();
122 
123  // Return zero if no likelyhood was estimated
124  auto it = likelyhoodMap.find(shape_index);
125  if (it == likelyhoodMap.end())
126  return 0;
127 
128  return it->second;
129  }
130 
131  private:
132 
134  std::map<int, TH2F> m_gridmap;
136  std::map<int, std::vector<PXDClusterShapeClassifierPar> > m_shapeClassifiers;
137 
139  };
141 } // end of namespace Belle2
Belle2::PXDClusterPositionEstimatorPar::addGrid
void addGrid(int clusterkind, const TH2F &grid)
Add grid for clusterkind.
Definition: PXDClusterPositionEstimatorPar.h:46
Belle2::PXDClusterShapeClassifierPar
The class for PXD cluster shape classifier payload.
Definition: PXDClusterShapeClassifierPar.h:36
Belle2::PXDClusterOffsetPar
The class for PXD cluster position offset payload.
Definition: PXDClusterOffsetPar.h:32
Belle2::PXDClusterPositionEstimatorPar::getGridMap
const std::map< int, TH2F > & getGridMap() const
Return grid map.
Definition: PXDClusterPositionEstimatorPar.h:61
Belle2::PXDClusterPositionEstimatorPar::getShapeLikelyhood
float getShapeLikelyhood(int shape_index, double thetaU, double thetaV, int clusterkind) const
Returns shape likelyhood.
Definition: PXDClusterPositionEstimatorPar.h:123
Belle2::PXDClusterPositionEstimatorPar::hasClassifier
bool hasClassifier(double thetaU, double thetaV, int clusterkind) const
Returns True if there is a classifier available.
Definition: PXDClusterPositionEstimatorPar.h:92
Belle2::PXDClusterPositionEstimatorPar::setShapeClassifier
void setShapeClassifier(const PXDClusterShapeClassifierPar &classifier, int uBin, int vBin, int clusterkind)
Set shape classifier.
Definition: PXDClusterPositionEstimatorPar.h:67
Belle2::PXDClusterPositionEstimatorPar
The class for PXD cluster position lookup table payload.
Definition: PXDClusterPositionEstimatorPar.h:37
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::PXDClusterPositionEstimatorPar::getOffset
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...
Definition: PXDClusterPositionEstimatorPar.h:112
Belle2::PXDClusterShapeClassifierPar::getShapeLikelyhoodMap
const std::map< int, float > & getShapeLikelyhoodMap() const
Return shape likelyhood map
Definition: PXDClusterShapeClassifierPar.h:45
Belle2::PXDClusterPositionEstimatorPar::m_gridmap
std::map< int, TH2F > m_gridmap
Map of angular grids for different clusterkinds
Definition: PXDClusterPositionEstimatorPar.h:142
Belle2::PXDClusterPositionEstimatorPar::ClassDef
ClassDef(PXDClusterPositionEstimatorPar, 2)
ClassDef, must be the last term before the closing {}.
Belle2::PXDClusterPositionEstimatorPar::~ PXDClusterPositionEstimatorPar
~ PXDClusterPositionEstimatorPar()
Destructor.
Definition: PXDClusterPositionEstimatorPar.h:43
Belle2::PXDClusterPositionEstimatorPar::PXDClusterPositionEstimatorPar
PXDClusterPositionEstimatorPar()
Default constructor.
Definition: PXDClusterPositionEstimatorPar.h:40
Belle2::PXDClusterPositionEstimatorPar::getShapeClassifier
const PXDClusterShapeClassifierPar & getShapeClassifier(int uBin, int vBin, int clusterkind) const
Returns shape classifier.
Definition: PXDClusterPositionEstimatorPar.h:74
Belle2::PXDClusterPositionEstimatorPar::m_shapeClassifiers
std::map< int, std::vector< PXDClusterShapeClassifierPar > > m_shapeClassifiers
Map of cluster shape classifiers for different clusterkinds.
Definition: PXDClusterPositionEstimatorPar.h:144