Belle II Software  release-06-01-15
PXDClusterPositionEstimator.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 <framework/logging/Logger.h>
10 #include <pxd/reconstruction/PXDClusterPositionEstimator.h>
11 #include <vxd/dataobjects/VxdID.h>
12 #include <pxd/geometry/SensorInfo.h>
13 #include <vxd/geometry/GeoCache.h>
14 
15 #include <cmath>
16 #include <numeric>
17 
18 using namespace std;
19 
20 //namespace Belle2 {
21 // namespace PXD {
22 
24 {
25  m_shapeIndexFromDB = unique_ptr<Belle2::DBObjPtr<Belle2::PXDClusterShapeIndexPar>>(new
27  m_positionEstimatorFromDB = unique_ptr<Belle2::DBObjPtr<Belle2::PXDClusterPositionEstimatorPar>>
29 
30  if ((*m_shapeIndexFromDB).isValid() && (*m_positionEstimatorFromDB).isValid()) {
31  setShapeIndexFromDB();
32  (*m_shapeIndexFromDB).addCallback(this, &Belle2::PXD::PXDClusterPositionEstimator::setShapeIndexFromDB);
33 
34  setPositionEstimatorFromDB();
35  (*m_positionEstimatorFromDB).addCallback(this, &Belle2::PXD::PXDClusterPositionEstimator::setPositionEstimatorFromDB);
36  }
37 }
38 
40 {
41  m_shapeIndexPar = **m_shapeIndexFromDB;
42 }
43 
45 {
46  m_positionEstimatorPar = **m_positionEstimatorFromDB;
47 }
48 
50 {
51  static std::unique_ptr<Belle2::PXD::PXDClusterPositionEstimator> instance(new Belle2::PXD::PXDClusterPositionEstimator());
52  return *instance;
53 }
54 
55 
57  double tu,
58  double tv) const
59 {
60  double thetaU = TMath::ATan2(tu, 1.0) * 180.0 / M_PI;
61  double thetaV = TMath::ATan2(tv, 1.0) * 180.0 / M_PI;
62  int sector_index = getSectorIndex(thetaU, thetaV);
63 
64  int clusterkind = cluster.getKind();
65  int shape_index = cluster.getSectorShapeIndices().at(sector_index);
66  float eta = cluster.getSectorEtaValues().at(sector_index);
67  return m_positionEstimatorPar.getOffset(shape_index, eta, thetaU, thetaV, clusterkind);
68 }
69 
70 
72 {
73  double thetaU = TMath::ATan2(tu, 1.0) * 180.0 / M_PI;
74  double thetaV = TMath::ATan2(tv, 1.0) * 180.0 / M_PI;
75  int clusterkind = cluster.getKind();
76  int sector_index = getSectorIndex(thetaU, thetaV);
77  int shape_index = cluster.getSectorShapeIndices().at(sector_index);
78  return m_positionEstimatorPar.getShapeLikelyhood(shape_index, thetaU, thetaV, clusterkind);
79 }
80 
81 
82 float Belle2::PXD::PXDClusterPositionEstimator::computeEta(const std::set<Belle2::PXD::Pixel>& pixels, int vStart, int vSize,
83  double thetaU, double thetaV) const
84 {
85  const Belle2::PXD::Pixel& headPixel = getHeadPixel(pixels, vStart, vSize, thetaU, thetaV);
86  const Belle2::PXD::Pixel& tailPixel = getTailPixel(pixels, vStart, vSize, thetaU, thetaV);
87  float eta = 0;
88  if (headPixel.getIndex() != tailPixel.getIndex()) {
89  eta = tailPixel.getCharge();
90  eta /= (tailPixel.getCharge() + headPixel.getCharge());
91  } else {
92  eta = tailPixel.getCharge();
93  }
94  return eta;
95 }
96 
97 const Belle2::PXD::Pixel& Belle2::PXD::PXDClusterPositionEstimator::getHeadPixel(const std::set<Belle2::PXD::Pixel>& pixels,
98  int vStart, int vSize, double thetaU,
99  double thetaV) const
100 {
101  if (thetaV >= 0) {
102  if (thetaU >= 0) {
103  return getLastPixelWithVOffset(pixels, vStart, vSize - 1);
104  } else {
105  return getFirstPixelWithVOffset(pixels, vStart, vSize - 1);
106  }
107  } else {
108  if (thetaU >= 0) {
109  return getLastPixelWithVOffset(pixels, vStart, 0);
110  } else {
111  return getFirstPixelWithVOffset(pixels, vStart, 0);
112  }
113  }
114 }
115 
116 const Belle2::PXD::Pixel& Belle2::PXD::PXDClusterPositionEstimator::getTailPixel(const std::set<Belle2::PXD::Pixel>& pixels,
117  int vStart, int vSize, double thetaU,
118  double thetaV) const
119 {
120  if (thetaV >= 0) {
121  if (thetaU >= 0) {
122  return getFirstPixelWithVOffset(pixels, vStart, 0);
123  } else {
124  return getLastPixelWithVOffset(pixels, vStart, 0);
125  }
126  } else {
127  if (thetaU >= 0) {
128  return getFirstPixelWithVOffset(pixels, vStart, vSize - 1);
129  } else {
130  return getLastPixelWithVOffset(pixels, vStart, vSize - 1);
131  }
132  }
133 }
134 
136  pixels,
137  int vStart, int vOffset) const
138 {
139  for (auto pxit = pixels.cbegin(); pxit != pixels.cend(); ++pxit) {
140  int v = pxit->getV() - vStart;
141  if (vOffset < v) {
142  if (pxit == pixels.cbegin()) {
143  } else {
144  pxit--;
145  return *pxit;
146  }
147  }
148  }
149  if (pixels.empty())
150  B2FATAL("Found cluster with empty pixel set. ");
151 
152  auto pxit = --pixels.cend();
153  return *pxit;
154 }
155 
157  const std::set<Belle2::PXD::Pixel>& pixels,
158  int vStart, int vOffset) const
159 {
160  for (const Belle2::PXD::Pixel& px : pixels) {
161  int v = px.getV() - vStart;
162  if (vOffset == v) {
163  return px;
164  }
165  }
166  if (pixels.empty())
167  B2FATAL("Found cluster with empty pixel set. ");
168 
169  return *pixels.cbegin();
170 }
171 
172 const std::string Belle2::PXD::PXDClusterPositionEstimator::getShortName(const std::set<Belle2::PXD::Pixel>& pixels, int uStart,
173  int vStart, int vSize, double thetaU,
174  double thetaV) const
175 {
176  const Belle2::PXD::Pixel& headPixel = getHeadPixel(pixels, vStart, vSize, thetaU, thetaV);
177  const Belle2::PXD::Pixel& tailPixel = getTailPixel(pixels, vStart, vSize, thetaU, thetaV);
178  std::string name = "S";
179 
180  name += "D" + std::to_string(tailPixel.getV() - vStart) + '.' + std::to_string(tailPixel.getU() - uStart);
181 
182  if (headPixel.getIndex() != tailPixel.getIndex()) {
183  name += "D" + std::to_string(headPixel.getV() - vStart) + '.' + std::to_string(headPixel.getU() - uStart);
184  }
185  return name;
186 }
187 
188 int Belle2::PXD::PXDClusterPositionEstimator::computeShapeIndex(const std::set<Belle2::PXD::Pixel>& pixels, int uStart, int vStart,
189  int vSize, double thetaU,
190  double thetaV) const
191 {
192  // Compute shape name
193  auto shape_name = getShortName(pixels, uStart, vStart, vSize, thetaU, thetaV);
194  // Return shape index
195  return m_shapeIndexPar.getShapeIndex(shape_name);
196 }
197 
198 
199 const std::string Belle2::PXD::PXDClusterPositionEstimator::getMirroredShortName(const std::set<Belle2::PXD::Pixel>& pixels,
200  int uStart,
201  int vStart, int vSize, double thetaU,
202  double thetaV) const
203 {
204  const Belle2::PXD::Pixel& headPixel = getHeadPixel(pixels, vStart, vSize, thetaU, thetaV);
205  const Belle2::PXD::Pixel& tailPixel = getTailPixel(pixels, vStart, vSize, thetaU, thetaV);
206  int vmax = vSize - 1;
207 
208  std::string name = "S";
209  name += "D" + std::to_string(vmax - tailPixel.getV() + vStart) + '.' + std::to_string(tailPixel.getU() - uStart);
210 
211  if (headPixel.getIndex() != tailPixel.getIndex()) {
212  name += "D" + std::to_string(vmax - headPixel.getV() + vStart) + '.' + std::to_string(headPixel.getU() - uStart);
213  }
214  return name;
215 }
216 
217 const std::string Belle2::PXD::PXDClusterPositionEstimator::getFullName(const std::set<Belle2::PXD::Pixel>& pixels, int uStart,
218  int vStart) const
219 {
220  return std::accumulate(pixels.begin(), pixels.end(), std::string("F"),
221  [uStart, vStart](auto name, auto px) {
222  return name + "D" + std::to_string(px.getV() - vStart) + "." + std::to_string(px.getU() - uStart);
223  });
224 }
225 
227 {
228  std::set<int> pixelkinds;
229  bool uEdge = false;
230  bool vEdge = false;
231 
232  Belle2::VxdID sensorID = cluster.getSensorID();
233  const Belle2::PXD::SensorInfo& Info = dynamic_cast<const Belle2::PXD::SensorInfo&>(Belle2::VXD::GeoCache::get(sensorID));
234 
235  for (const Belle2::PXDDigit& digit : cluster.getRelationsTo<Belle2::PXDDigit>("PXDDigits")) {
236  int pixelkind = Info.getPixelKindNew(sensorID, digit.getVCellID());
237  pixelkinds.insert(pixelkind);
238 
239  // Cluster at v sensor edge
240  if (digit.getVCellID() == 0 or digit.getVCellID() >= 767)
241  vEdge = true;
242  // Cluster at u sensor edge
243  if (digit.getUCellID() == 0 or digit.getUCellID() >= 249)
244  uEdge = true;
245  }
246 
247  // In most cases, clusterkind is just pixelkind of first digit
248  int clusterkind = *pixelkinds.begin();
249 
250  // Clusters with different pixelkinds or edge digits are special
251  // TODO: At the moment, clusterkind >3 will not be corrected
252  if (pixelkinds.size() > 1 || uEdge || vEdge)
253  clusterkind = 4;
254 
255  return clusterkind;
256 }
257 
258 int Belle2::PXD::PXDClusterPositionEstimator::getClusterkind(const std::vector<Belle2::PXD::Pixel>& pixels,
259  const Belle2::VxdID& sensorID) const
260 {
261  std::set<int> pixelkinds;
262  bool uEdge = false;
263  bool vEdge = false;
264 
265  const Belle2::PXD::SensorInfo& Info = dynamic_cast<const Belle2::PXD::SensorInfo&>(Belle2::VXD::GeoCache::get(sensorID));
266 
267  for (const Belle2::PXD::Pixel& pix : pixels) {
268  int pixelkind = Info.getPixelKindNew(sensorID, pix.getV());
269  pixelkinds.insert(pixelkind);
270 
271  // Cluster at v sensor edge
272  if (pix.getV() == 0 or pix.getV() >= 767)
273  vEdge = true;
274  // Cluster at u sensor edge
275  if (pix.getU() == 0 or pix.getU() >= 249)
276  uEdge = true;
277  }
278 
279  // In most cases, clusterkind is just pixelkind of first digit
280  int clusterkind = *pixelkinds.begin();
281 
282  // Clusters with different pixelkinds or edge digits are special
283  // TODO: At the moment, clusterkind >3 will not be corrected
284  if (pixelkinds.size() > 1 || uEdge || vEdge)
285  clusterkind = 4;
286 
287  return clusterkind;
288 }
289 
290 int Belle2::PXD::PXDClusterPositionEstimator::getSectorIndex(double thetaU, double thetaV) const
291 {
292  int sectorIndex = 0;
293  if (thetaU >= 0) {
294  if (thetaV >= 0) {
295  sectorIndex = 0;
296  } else {
297  sectorIndex = 3;
298  }
299  } else {
300  if (thetaV >= 0) {
301  sectorIndex = 1;
302  } else {
303  sectorIndex = 2;
304  }
305  }
306  return sectorIndex;
307 }
308 
309 
310 
311 // }
312 //}
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
The class for PXD cluster position offset payload.
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
The PXD digit class.
Definition: PXDDigit.h:27
Singleton class that estimates cluster positions taking into account the estimated track incidence an...
void initialize()
Initialize PXDClusterPositionEstimator from DB.
const std::string getFullName(const std::set< Pixel > &pixels, int uStart, int vStart) const
Return a name for the pixel set.
static PXDClusterPositionEstimator & getInstance()
Main (and only) way to access the PXDClusterPositionEstimator.
int getSectorIndex(double thetaU, double thetaV) const
Get sector index from angles.
const Pixel & getTailPixel(const std::set< Pixel > &pixels, int vStart, int vSize, double thetaU, double thetaV) const
Return reference to the tail pixel in pixel set.
const std::string getMirroredShortName(const std::set< Pixel > &pixels, int uStart, int vStart, int vSize, double thetaU, double thetaV) const
Return the mirrored name for the pixel set.
const Pixel & getLastPixelWithVOffset(const std::set< Pixel > &pixels, int vStart, int vOffset) const
Return reference to the last pixel in pixel set with given vOffset from vStart.
const PXDClusterOffsetPar * getClusterOffset(const PXDCluster &cluster, double tu, double tv) const
Return pointer to cluster offsets, can be nullptr.
float computeEta(const std::set< Pixel > &pixels, int vStart, int vSize, double thetaU, double thetaV) const
Return the normed charge ratio between head and tail pixels (size>=2) or the charge of the seed (size...
float getShapeLikelyhood(const PXDCluster &cluster, double tu, double tv) const
Return cluster shape likelyhood.
int computeShapeIndex(const std::set< Pixel > &pixels, int uStart, int vStart, int vSize, double thetaU, double thetaV) const
Return the shape index of the pixels.
void setPositionEstimatorFromDB()
Set PositionEstimator from DB.
const std::string getShortName(const std::set< Pixel > &pixels, int uStart, int vStart, int vSize, double thetaU, double thetaV) const
Return the name for the pixel set.
int getClusterkind(const PXDCluster &cluster) const
Return kind of cluster needed to find cluster position correction.
const Pixel & getFirstPixelWithVOffset(const std::set< Pixel > &pixels, int vStart, int vOffset) const
Return reference to the first pixel in pixel set with given vOffset from vStart.
const Pixel & getHeadPixel(const std::set< Pixel > &pixels, int vStart, int vSize, double thetaU, double thetaV) const
Return reference to the head pixel in pixel set.
Class to represent one pixel, used in clustering for fast access.
Definition: Pixel.h:36
unsigned short getU() const
Return the CellID in u.
Definition: Pixel.h:66
unsigned int getIndex() const
Return the Index of the digit.
Definition: Pixel.h:72
float getCharge() const
Return the Charge of the Pixel.
Definition: Pixel.h:70
unsigned short getV() const
Return the CellID in v.
Definition: Pixel.h:68
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
int getPixelKindNew(const VxdID &sensorID, int vID) const
Return pixel kind ID.
Definition: SensorInfo.cc:77
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33