CMS 3D CMS Logo

SiPixelCluster.h
Go to the documentation of this file.
1 #ifndef DataFormats_SiPixel_Cluster_SiPixelCluster_h
2 #define DataFormats_SiPixel_Cluster_SiPixelCluster_h
3 
4 //---------------------------------------------------------------------------
19 //---------------------------------------------------------------------------
20 
21 #include <vector>
22 #include <cstdint>
23 #include <cassert>
24 
25 class PixelDigi;
26 
28 public:
29  class Pixel {
30  public:
31  constexpr Pixel() : x(0), y(0), adc(0) {} // for root
32  constexpr Pixel(int pix_x, int pix_y, int pix_adc) : x(pix_x), y(pix_y), adc(pix_adc) {}
33  uint16_t x;
34  uint16_t y;
35  uint16_t adc;
36  };
37 
38  //--- Integer shift in x and y directions.
39  class Shift {
40  public:
41  constexpr Shift(int dx, int dy) : dx_(dx), dy_(dy) {}
42  constexpr Shift() : dx_(0), dy_(0) {}
43  constexpr int dx() const { return dx_; }
44  constexpr int dy() const { return dy_; }
45 
46  private:
47  int dx_;
48  int dy_;
49  };
50 
51  //--- Position of a SiPixel
52  class PixelPos {
53  public:
54  constexpr PixelPos() : row_(0), col_(0) {}
55  constexpr PixelPos(int row, int col) : row_(row), col_(col) {}
56  constexpr int row() const { return row_; }
57  constexpr int col() const { return col_; }
58  constexpr PixelPos operator+(const Shift& shift) const { return PixelPos(row() + shift.dx(), col() + shift.dy()); }
59 
60  private:
61  int row_;
62  int col_;
63  };
64 
65  typedef std::vector<PixelDigi>::const_iterator PixelDigiIter;
66  typedef std::pair<PixelDigiIter, PixelDigiIter> PixelDigiRange;
67 
68  static constexpr unsigned int MAXSPAN = 255;
69  static constexpr unsigned int MAXPOS = 2047;
70 
76 
77  SiPixelCluster(unsigned int isize,
78  uint16_t const* adcs,
79  uint16_t const* xpos,
80  uint16_t const* ypos,
81  uint16_t const xmin,
82  uint16_t const ymin)
83  : thePixelOffset(2 * isize), thePixelADC(adcs, adcs + isize) {
84  uint16_t maxCol = 0;
85  uint16_t maxRow = 0;
86  for (unsigned int i = 0; i != isize; ++i) {
87  uint16_t xoffset = xpos[i] - xmin;
88  uint16_t yoffset = ypos[i] - ymin;
89  thePixelOffset[i * 2] = std::min(uint16_t(MAXSPAN), xoffset);
90  thePixelOffset[i * 2 + 1] = std::min(uint16_t(MAXSPAN), yoffset);
91  if (xoffset > maxRow)
92  maxRow = xoffset;
93  if (yoffset > maxCol)
94  maxCol = yoffset;
95  }
96  packRow(xmin, maxRow);
97  packCol(ymin, maxCol);
98  }
99 
100  // obsolete (only for regression tests)
101  SiPixelCluster(const PixelPos& pix, int adc);
102  void add(const PixelPos& pix, int adc);
103 
104  // Analog linear average position (barycenter)
105  float x() const {
106  float qm = 0.0;
107  int isize = thePixelADC.size();
108  for (int i = 0; i < isize; ++i)
109  qm += float(thePixelADC[i]) * (thePixelOffset[i * 2] + minPixelRow() + 0.5f);
110  return qm / charge();
111  }
112 
113  float y() const {
114  float qm = 0.0;
115  int isize = thePixelADC.size();
116  for (int i = 0; i < isize; ++i)
117  qm += float(thePixelADC[i]) * (thePixelOffset[i * 2 + 1] + minPixelCol() + 0.5f);
118  return qm / charge();
119  }
120 
121  // Return number of pixels.
122  int size() const { return thePixelADC.size(); }
123 
124  // Return cluster dimension in the x direction.
125  int sizeX() const { return rowSpan() + 1; }
126 
127  // Return cluster dimension in the y direction.
128  int sizeY() const { return colSpan() + 1; }
129 
130  inline int charge() const {
131  int qm = 0;
132  int isize = thePixelADC.size();
133  for (int i = 0; i < isize; ++i)
134  qm += thePixelADC[i];
135  return qm;
136  } // Return total cluster charge.
137 
138  inline int minPixelRow() const { return theMinPixelRow; } // The min x index.
139  inline int maxPixelRow() const { return minPixelRow() + rowSpan(); } // The max x index.
140  inline int minPixelCol() const { return theMinPixelCol; } // The min y index.
141  inline int maxPixelCol() const { return minPixelCol() + colSpan(); } // The max y index.
142 
143  const std::vector<uint8_t>& pixelOffset() const { return thePixelOffset; }
144  const std::vector<uint16_t>& pixelADC() const { return thePixelADC; }
145 
146  // obsolete, use single pixel access below
147  const std::vector<Pixel> pixels() const {
148  std::vector<Pixel> oldPixVector;
149  int isize = thePixelADC.size();
150  oldPixVector.reserve(isize);
151  for (int i = 0; i < isize; ++i) {
152  oldPixVector.push_back(pixel(i));
153  }
154  return oldPixVector;
155  }
156 
157  // infinite faster than above...
158  Pixel pixel(int i) const {
159  return Pixel(minPixelRow() + thePixelOffset[i * 2], minPixelCol() + thePixelOffset[i * 2 + 1], thePixelADC[i]);
160  }
161 
162 private:
163  static int overflow_(uint16_t span) { return span == uint16_t(MAXSPAN); }
164 
165 public:
166  int colSpan() const { return thePixelColSpan; }
167 
168  int rowSpan() const { return thePixelRowSpan; }
169 
170  bool overflowCol() const { return overflow_(thePixelColSpan); }
171 
172  bool overflowRow() const { return overflow_(thePixelRowSpan); }
173 
174  bool overflow() const { return overflowCol() || overflowRow(); }
175 
176  void packCol(uint16_t ymin, uint16_t yspan) {
178  thePixelColSpan = std::min(yspan, uint16_t(MAXSPAN));
179  }
180  void packRow(uint16_t xmin, uint16_t xspan) {
182  thePixelRowSpan = std::min(xspan, uint16_t(MAXSPAN));
183  }
184 
185  // ggiurgiu@fnal.gov, 01/05/12
186  // Getters and setters for the newly added data members (err_x and err_y). See below.
187  void setSplitClusterErrorX(float errx) { err_x = errx; }
188  void setSplitClusterErrorY(float erry) { err_y = erry; }
189  float getSplitClusterErrorX() const { return err_x; }
190  float getSplitClusterErrorY() const { return err_y; }
191 
192 private:
193  std::vector<uint8_t> thePixelOffset;
194  std::vector<uint16_t> thePixelADC;
195 
196  uint16_t theMinPixelRow = MAXPOS; // Minimum pixel index in the x direction (low edge).
197  uint16_t theMinPixelCol = MAXPOS; // Minimum pixel index in the y direction (left edge).
198  uint8_t thePixelRowSpan = 0; // Span pixel index in the x direction (low edge).
199  uint8_t thePixelColSpan = 0; // Span pixel index in the y direction (left edge).
200 
201  float err_x = -99999.9f;
202  float err_y = -99999.9f;
203 };
204 
205 // Comparison operators (needed by DetSetVector)
206 inline bool operator<(const SiPixelCluster& one, const SiPixelCluster& other) {
207  if (one.minPixelRow() < other.minPixelRow()) {
208  return true;
209  } else if (one.minPixelRow() > other.minPixelRow()) {
210  return false;
211  } else if (one.minPixelCol() < other.minPixelCol()) {
212  return true;
213  } else {
214  return false;
215  }
216 }
217 
222 
227 
230 #endif
int minPixelCol() const
void setSplitClusterErrorY(float erry)
constexpr int dx() const
bool overflowCol() const
std::vector< PixelDigi >::const_iterator PixelDigiIter
std::pair< PixelDigiIter, PixelDigiIter > PixelDigiRange
static int overflow_(uint16_t span)
int rowSpan() const
int charge() const
constexpr Shift(int dx, int dy)
edm::RefProd< SiPixelClusterCollection > SiPixelClusterRefProd
int maxPixelRow() const
void packRow(uint16_t xmin, uint16_t xspan)
static constexpr unsigned int MAXPOS
void packCol(uint16_t ymin, uint16_t yspan)
constexpr PixelPos operator+(const Shift &shift) const
std::vector< uint16_t > thePixelADC
int minPixelRow() const
bool overflow() const
int colSpan() const
constexpr Pixel(int pix_x, int pix_y, int pix_adc)
uint16_t theMinPixelRow
SiPixelCluster(unsigned int isize, uint16_t const *adcs, uint16_t const *xpos, uint16_t const *ypos, uint16_t const xmin, uint16_t const ymin)
edm::Ref< SiPixelClusterCollection, SiPixelCluster > SiPixelClusterRef
void add(const PixelPos &pix, int adc)
T min(T a, T b)
Definition: MathUtil.h:58
bool overflowRow() const
edm::DetSetRefVector< SiPixelCluster > SiPixelClusterRefVector
const std::vector< uint8_t > & pixelOffset() const
bool operator<(const SiPixelCluster &one, const SiPixelCluster &other)
const std::vector< uint16_t > & pixelADC() const
std::vector< uint8_t > thePixelOffset
constexpr PixelPos(int row, int col)
constexpr int col() const
void setSplitClusterErrorX(float errx)
float getSplitClusterErrorY() const
edm::DetSetVector< SiPixelCluster > SiPixelClusterCollection
edm::Ref< SiPixelClusterCollectionNew, SiPixelCluster > SiPixelClusterRefNew
int maxPixelCol() const
constexpr int dy() const
edmNew::DetSetVector< SiPixelCluster > SiPixelClusterCollectionNew
int sizeY() const
uint16_t theMinPixelCol
Pixel cluster – collection of neighboring pixels above threshold.
Pixel pixel(int i) const
float getSplitClusterErrorX() const
col
Definition: cuy.py:1010
static unsigned int const shift
float y() const
static constexpr unsigned int MAXSPAN
uint8_t thePixelColSpan
uint8_t thePixelRowSpan
int sizeX() const
float x() const
constexpr int row() const
const std::vector< Pixel > pixels() const
#define constexpr
int size() const