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>
31 if ((*m_shapeIndexFromDB).isValid() && (*m_positionEstimatorFromDB).isValid()) {
42 m_shapeIndexPar = **m_shapeIndexFromDB;
47 m_positionEstimatorPar = **m_positionEstimatorFromDB;
61 double thetaU = TMath::ATan2(tu, 1.0) * 180.0 / M_PI;
62 double thetaV = TMath::ATan2(tv, 1.0) * 180.0 / M_PI;
63 int sector_index = getSectorIndex(thetaU, thetaV);
65 int clusterkind = cluster.getKind();
66 int shape_index = cluster.getSectorShapeIndices().at(sector_index);
67 float eta = cluster.getSectorEtaValues().at(sector_index);
68 return m_positionEstimatorPar.getOffset(shape_index, eta, thetaU, thetaV, clusterkind);
74 double thetaU = TMath::ATan2(tu, 1.0) * 180.0 / M_PI;
75 double thetaV = TMath::ATan2(tv, 1.0) * 180.0 / M_PI;
76 int clusterkind = cluster.getKind();
77 int sector_index = getSectorIndex(thetaU, thetaV);
78 int shape_index = cluster.getSectorShapeIndices().at(sector_index);
79 return m_positionEstimatorPar.getShapeLikelyhood(shape_index, thetaU, thetaV, clusterkind);
84 double thetaU,
double thetaV)
const
86 const Belle2::PXD::Pixel& headPixel = getHeadPixel(pixels, vStart, vSize, thetaU, thetaV);
87 const Belle2::PXD::Pixel& tailPixel = getTailPixel(pixels, vStart, vSize, thetaU, thetaV);
99 int vStart,
int vSize,
double thetaU,
104 return getLastPixelWithVOffset(pixels, vStart, vSize - 1);
106 return getFirstPixelWithVOffset(pixels, vStart, vSize - 1);
110 return getLastPixelWithVOffset(pixels, vStart, 0);
112 return getFirstPixelWithVOffset(pixels, vStart, 0);
118 int vStart,
int vSize,
double thetaU,
123 return getFirstPixelWithVOffset(pixels, vStart, 0);
125 return getLastPixelWithVOffset(pixels, vStart, 0);
129 return getFirstPixelWithVOffset(pixels, vStart, vSize - 1);
131 return getLastPixelWithVOffset(pixels, vStart, vSize - 1);
138 int vStart,
int vOffset)
const
140 for (
auto pxit = pixels.cbegin(); pxit != pixels.cend(); ++pxit) {
141 int v = pxit->
getV() - vStart;
143 if (pxit == pixels.cbegin()) {
151 B2FATAL(
"Found cluster with empty pixel set. ");
153 auto pxit = --pixels.cend();
158 const std::set<Belle2::PXD::Pixel>& pixels,
159 int vStart,
int vOffset)
const
162 int v = px.
getV() - vStart;
168 B2FATAL(
"Found cluster with empty pixel set. ");
170 return *pixels.cbegin();
174 int vStart,
int vSize,
double thetaU,
177 const Belle2::PXD::Pixel& headPixel = getHeadPixel(pixels, vStart, vSize, thetaU, thetaV);
178 const Belle2::PXD::Pixel& tailPixel = getTailPixel(pixels, vStart, vSize, thetaU, thetaV);
179 std::string name =
"S";
181 name +=
"D" + std::to_string(tailPixel.
getV() - vStart) +
'.' + std::to_string(tailPixel.
getU() - uStart);
184 name +=
"D" + std::to_string(headPixel.
getV() - vStart) +
'.' + std::to_string(headPixel.
getU() - uStart);
190 int vSize,
double thetaU,
194 auto shape_name = getShortName(pixels, uStart, vStart, vSize, thetaU, thetaV);
196 return m_shapeIndexPar.getShapeIndex(shape_name);
202 int vStart,
int vSize,
double thetaU,
205 const Belle2::PXD::Pixel& headPixel = getHeadPixel(pixels, vStart, vSize, thetaU, thetaV);
206 const Belle2::PXD::Pixel& tailPixel = getTailPixel(pixels, vStart, vSize, thetaU, thetaV);
207 int vmax = vSize - 1;
209 std::string name =
"S";
210 name +=
"D" + std::to_string(vmax - tailPixel.
getV() + vStart) +
'.' + std::to_string(tailPixel.
getU() - uStart);
213 name +=
"D" + std::to_string(vmax - headPixel.
getV() + vStart) +
'.' + std::to_string(headPixel.
getU() - uStart);
221 return std::accumulate(pixels.begin(), pixels.end(), std::string(
"F"),
222 [uStart, vStart](
auto name,
auto px) {
223 return name +
"D" + std::to_string(px.getV() - vStart) +
"." + std::to_string(px.getU() - uStart);
229 std::set<int> pixelkinds;
239 pixelkinds.insert(pixelkind);
242 if (digit.getVCellID() == 0 or digit.getVCellID() >= 767)
245 if (digit.getUCellID() == 0 or digit.getUCellID() >= 249)
250 int clusterkind = *pixelkinds.begin();
254 if (pixelkinds.size() > 1 || uEdge || vEdge)
263 std::set<int> pixelkinds;
272 pixelkinds.insert(pixelkind);
275 if (pix.getV() == 0 or pix.getV() >= 767)
278 if (pix.getU() == 0 or pix.getU() >= 249)
283 int clusterkind = *pixelkinds.begin();
287 if (pixelkinds.size() > 1 || uEdge || vEdge)
Class for accessing objects in the database.
The class for PXD cluster position offset payload.
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
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.
std::unique_ptr< DBObjPtr< PXDClusterPositionEstimatorPar > > m_positionEstimatorFromDB
PXDClusterPositionEstimatorPar retrieved from DB.
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.
std::unique_ptr< DBObjPtr< PXDClusterShapeIndexPar > > m_shapeIndexFromDB
PXDClusterShapeIndex retrieved from DB.
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.
void setShapeIndexFromDB()
Set ShapeIndex from DB.
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.
unsigned short getU() const
Return the CellID in u.
unsigned int getIndex() const
Return the Index of the digit.
float getCharge() const
Return the Charge of the Pixel.
unsigned short getV() const
Return the CellID in v.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
int getPixelKindNew(const VxdID &sensorID, int vID) const
Return pixel kind ID.
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
static GeoCache & getInstance()
Return a reference to the singleton instance.
Class to uniquely identify a any structure of the PXD and SVD.