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