CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
64  private:
65  int row_;
66  int col_;
67  };
68 
69  typedef std::vector<PixelDigi>::const_iterator PixelDigiIter;
70  typedef std::pair<PixelDigiIter, PixelDigiIter> PixelDigiRange;
71 
72  static constexpr unsigned int MAXSPAN = 255;
73  static constexpr unsigned int MAXPOS = 2047;
74 
76 
81  SiPixelCluster() = default;
82  ~SiPixelCluster() = default;
83  SiPixelCluster(SiPixelCluster const&) = default;
84  SiPixelCluster(SiPixelCluster&&) = default;
85  SiPixelCluster& operator=(SiPixelCluster const&) = default;
87 
88  SiPixelCluster(unsigned int isize,
89  uint16_t const* adcs,
90  uint16_t const* xpos,
91  uint16_t const* ypos,
92  uint16_t xmin,
93  uint16_t ymin,
94  uint16_t id = invalidClusterId)
95  : thePixelOffset(2 * isize), thePixelADC(adcs, adcs + isize), theOriginalClusterId(id) {
96  uint16_t maxCol = 0;
97  uint16_t maxRow = 0;
98  for (unsigned int i = 0; i < isize; ++i) {
99  uint16_t xoffset = xpos[i] - xmin;
100  uint16_t yoffset = ypos[i] - ymin;
101  thePixelOffset[i * 2] = std::min(uint16_t(MAXSPAN), xoffset);
102  thePixelOffset[i * 2 + 1] = std::min(uint16_t(MAXSPAN), yoffset);
103  if (xoffset > maxRow)
104  maxRow = xoffset;
105  if (yoffset > maxCol)
106  maxCol = yoffset;
107  }
108  packRow(xmin, maxRow);
109  packCol(ymin, maxCol);
110  }
111 
112  // obsolete (only for regression tests)
113  SiPixelCluster(const PixelPos& pix, int adc);
114  void add(const PixelPos& pix, int adc);
115 
116  // Analog linear average position (barycenter)
117  float x() const {
118  float qm = 0.0;
119  int isize = thePixelADC.size();
120  for (int i = 0; i < isize; ++i)
121  qm += float(thePixelADC[i]) * (thePixelOffset[i * 2] + minPixelRow() + 0.5f);
122  return qm / charge();
123  }
124 
125  float y() const {
126  float qm = 0.0;
127  int isize = thePixelADC.size();
128  for (int i = 0; i < isize; ++i)
129  qm += float(thePixelADC[i]) * (thePixelOffset[i * 2 + 1] + minPixelCol() + 0.5f);
130  return qm / charge();
131  }
132 
133  // Return number of pixels.
134  int size() const { return thePixelADC.size(); }
135 
136  // Return cluster dimension in the x direction.
137  int sizeX() const { return rowSpan() + 1; }
138 
139  // Return cluster dimension in the y direction.
140  int sizeY() const { return colSpan() + 1; }
141 
142  inline int charge() const {
143  int qm = 0;
144  int isize = thePixelADC.size();
145  for (int i = 0; i < isize; ++i)
146  qm += thePixelADC[i];
147  return qm;
148  } // Return total cluster charge.
149 
150  inline int minPixelRow() const { return theMinPixelRow; } // The min x index.
151  inline int maxPixelRow() const { return minPixelRow() + rowSpan(); } // The max x index.
152  inline int minPixelCol() const { return theMinPixelCol; } // The min y index.
153  inline int maxPixelCol() const { return minPixelCol() + colSpan(); } // The max y index.
154 
155  const std::vector<uint8_t>& pixelOffset() const { return thePixelOffset; }
156  const std::vector<uint16_t>& pixelADC() const { return thePixelADC; }
157 
158  // obsolete, use single pixel access below
159  const std::vector<Pixel> pixels() const {
160  std::vector<Pixel> oldPixVector;
161  int isize = thePixelADC.size();
162  oldPixVector.reserve(isize);
163  for (int i = 0; i < isize; ++i) {
164  oldPixVector.push_back(pixel(i));
165  }
166  return oldPixVector;
167  }
168 
169  // infinite faster than above...
170  Pixel pixel(int i) const {
171  return Pixel(minPixelRow() + thePixelOffset[i * 2], minPixelCol() + thePixelOffset[i * 2 + 1], thePixelADC[i]);
172  }
173 
174 private:
175  static int overflow_(uint16_t span) { return span == uint16_t(MAXSPAN); }
176 
177 public:
178  int colSpan() const { return thePixelColSpan; }
179 
180  int rowSpan() const { return thePixelRowSpan; }
181 
182  bool overflowCol() const { return overflow_(thePixelColSpan); }
183 
184  bool overflowRow() const { return overflow_(thePixelRowSpan); }
185 
186  bool overflow() const { return overflowCol() || overflowRow(); }
187 
188  void packCol(uint16_t ymin, uint16_t yspan) {
190  thePixelColSpan = std::min(yspan, uint16_t(MAXSPAN));
191  }
192  void packRow(uint16_t xmin, uint16_t xspan) {
194  thePixelRowSpan = std::min(xspan, uint16_t(MAXSPAN));
195  }
196 
197  // ggiurgiu@fnal.gov, 01/05/12
198  // Getters and setters for the newly added data members (err_x and err_y). See below.
199  void setSplitClusterErrorX(float errx) { err_x = errx; }
200  void setSplitClusterErrorY(float erry) { err_y = erry; }
201  float getSplitClusterErrorX() const { return err_x; }
202  float getSplitClusterErrorY() const { return err_y; }
203 
204  // the original id (they get sorted)
205  auto originalId() const { return theOriginalClusterId; }
206  void setOriginalId(uint16_t id) { theOriginalClusterId = id; }
207 
208 private:
209  std::vector<uint8_t> thePixelOffset;
210  std::vector<uint16_t> thePixelADC;
211 
212  uint16_t theMinPixelRow = MAXPOS; // Minimum pixel index in the x direction (low edge).
213  uint16_t theMinPixelCol = MAXPOS; // Minimum pixel index in the y direction (left edge).
214  uint8_t thePixelRowSpan = 0; // Span pixel index in the x direction (low edge).
215  uint8_t thePixelColSpan = 0; // Span pixel index in the y direction (left edge).
216 
218 
219  float err_x = -99999.9f;
220  float err_y = -99999.9f;
221 };
222 
223 // Comparison operators (needed by DetSetVector)
224 inline bool operator<(const SiPixelCluster& one, const SiPixelCluster& other) {
225  if (one.minPixelRow() < other.minPixelRow()) {
226  return true;
227  } else if (one.minPixelRow() > other.minPixelRow()) {
228  return false;
229  } else if (one.minPixelCol() < other.minPixelCol()) {
230  return true;
231  } else {
232  return false;
233  }
234 }
235 
240 
245 
248 #endif
int minPixelCol() const
void setSplitClusterErrorY(float erry)
constexpr int dx() const
uint16_t *__restrict__ id
bool overflowCol() const
SiPixelCluster()=default
std::vector< PixelDigi >::const_iterator PixelDigiIter
std::pair< PixelDigiIter, PixelDigiIter > PixelDigiRange
static int overflow_(uint16_t span)
int rowSpan() const
uint16_t theOriginalClusterId
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
static constexpr uint16_t invalidClusterId
void setOriginalId(uint16_t id)
bool overflow() const
int colSpan() const
constexpr Pixel(int pix_x, int pix_y, int pix_adc)
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)
uint16_t theMinPixelRow
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
const std::vector< uint16_t > & pixelADC() const
std::vector< uint8_t > thePixelOffset
constexpr PixelPos(int row, int col)
auto originalId() const
constexpr int col() const
~SiPixelCluster()=default
void setSplitClusterErrorX(float errx)
float getSplitClusterErrorY() const
SiPixelCluster & operator=(SiPixelCluster const &)=default
edm::DetSetVector< SiPixelCluster > SiPixelClusterCollection
edm::Ref< SiPixelClusterCollectionNew, SiPixelCluster > SiPixelClusterRefNew
int maxPixelCol() const
constexpr int dy() const
edmNew::DetSetVector< SiPixelCluster > SiPixelClusterCollectionNew
int sizeY() const
bool operator<(DTCELinkId const &lhs, DTCELinkId const &rhs)
Definition: DTCELinkId.h:70
uint16_t theMinPixelCol
Pixel cluster – collection of neighboring pixels above threshold.
Pixel pixel(int i) const
float getSplitClusterErrorX() const
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
uint16_t *__restrict__ uint16_t const *__restrict__ adc
int size() const