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 
30  class Pixel {
31  public:
32  constexpr Pixel() : x(0), y(0), adc(0){} // for root
33  constexpr Pixel(int pix_x, int pix_y, int pix_adc) :
34  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  class Shift {
42  public:
43  constexpr Shift( int dx, int dy) : dx_(dx), dy_(dy) {}
44  constexpr Shift() : dx_(0), dy_(0) {}
45  constexpr int dx() const { return dx_;}
46  constexpr int dy() const { return dy_;}
47  private:
48  int dx_;
49  int dy_;
50  };
51 
52  //--- Position of a SiPixel
53  class PixelPos {
54  public:
55  constexpr PixelPos() : row_(0), col_(0) {}
56  constexpr PixelPos(int row, int col) : row_(row) , col_(col) {}
57  constexpr int row() const { return row_;}
58  constexpr int col() const { return col_;}
60  return PixelPos( row() + shift.dx(), col() + shift.dy());
61  }
62  private:
63  int row_;
64  int col_;
65  };
66 
67  typedef std::vector<PixelDigi>::const_iterator PixelDigiIter;
68  typedef std::pair<PixelDigiIter,PixelDigiIter> PixelDigiRange;
69 
70 
71  static constexpr unsigned int MAXSPAN=255;
72  static constexpr unsigned int MAXPOS=2047;
73 
79 
80  SiPixelCluster(unsigned int isize, uint16_t const * adcs,
81  uint16_t const * xpos, uint16_t const * ypos,
82  uint16_t const xmin, 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) maxRow = xoffset;
92  if (yoffset > maxCol) maxCol = yoffset;
93  }
94  packRow(xmin,maxRow);
95  packCol(ymin,maxCol);
96  }
97 
98 
99  // obsolete (only for regression tests)
100  SiPixelCluster( const PixelPos& pix, int adc);
101  void add( const PixelPos& pix, int adc);
102 
103  // Analog linear average position (barycenter)
104  float x() const {
105  float qm = 0.0;
106  int isize = thePixelADC.size();
107  for (int i=0; i<isize; ++i)
108  qm += float(thePixelADC[i]) * (thePixelOffset[i*2] + minPixelRow() + 0.5f);
109  return qm/charge();
110  }
111 
112  float y() const {
113  float qm = 0.0;
114  int isize = thePixelADC.size();
115  for (int i=0; i<isize; ++i)
116  qm += float(thePixelADC[i]) * (thePixelOffset[i*2+1] + minPixelCol() + 0.5f);
117  return qm/charge();
118  }
119 
120  // Return number of pixels.
121  int size() const { return thePixelADC.size();}
122 
123  // Return cluster dimension in the x direction.
124  int sizeX() const { return rowSpan() +1;}
125 
126  // Return cluster dimension in the y direction.
127  int sizeY() const { return colSpan() +1;}
128 
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 
144  const std::vector<uint8_t> & pixelOffset() const { return thePixelOffset;}
145  const std::vector<uint16_t> & pixelADC() const { return thePixelADC;}
146 
147  // obsolete, use single pixel access below
148  const std::vector<Pixel> pixels() const {
149  std::vector<Pixel> oldPixVector;
150  int isize = thePixelADC.size();
151  oldPixVector.reserve(isize);
152  for(int i=0; i<isize; ++i) {
153  oldPixVector.push_back(pixel(i));
154  }
155  return oldPixVector;
156  }
157 
158  // infinite faster than above...
159  Pixel pixel(int i) const {
160  return Pixel(minPixelRow() + thePixelOffset[i*2],
161  minPixelCol() + thePixelOffset[i*2+1],
162  thePixelADC[i]
163  );
164  }
165 
166 private:
167 
168  static int overflow_(uint16_t span) { return span==uint16_t(MAXSPAN);}
169 
170 public:
171 
172  int colSpan() const {return thePixelColSpan; }
173 
174  int rowSpan() const { return thePixelRowSpan; }
175 
176 
177  bool overflowCol() const { return overflow_(thePixelColSpan); }
178 
179  bool overflowRow() const { return overflow_(thePixelRowSpan); }
180 
181  bool overflow() const { return overflowCol() || overflowRow(); }
182 
183  void packCol(uint16_t ymin, uint16_t yspan) {
185  thePixelColSpan = std::min(yspan, uint16_t(MAXSPAN));
186  }
187  void packRow(uint16_t xmin, uint16_t xspan) {
189  thePixelRowSpan = std::min(xspan, uint16_t(MAXSPAN));
190  }
191 
192  // ggiurgiu@fnal.gov, 01/05/12
193  // Getters and setters for the newly added data members (err_x and err_y). See below.
194  void setSplitClusterErrorX( float errx ) { err_x = errx; }
195  void setSplitClusterErrorY( float erry ) { err_y = erry; }
196  float getSplitClusterErrorX() const { return err_x; }
197  float getSplitClusterErrorY() const { return err_y; }
198 
199 
200 private:
201 
202  std::vector<uint8_t> thePixelOffset;
203  std::vector<uint16_t> thePixelADC;
204 
205 
206  uint16_t theMinPixelRow=MAXPOS; // Minimum pixel index in the x direction (low edge).
207  uint16_t theMinPixelCol=MAXPOS; // Minimum pixel index in the y direction (left edge).
208  uint8_t thePixelRowSpan=0; // Span pixel index in the x direction (low edge).
209  uint8_t thePixelColSpan=0; // Span pixel index in the y direction (left edge).
210 
211  float err_x=-99999.9f;
212  float err_y=-99999.9f;
213 
214 };
215 
216 
217 // Comparison operators (needed by DetSetVector)
218 inline bool operator<( const SiPixelCluster& one, const SiPixelCluster& other) {
219  if ( one.minPixelRow() < other.minPixelRow() ) {
220  return true;
221  } else if ( one.minPixelRow() > other.minPixelRow() ) {
222  return false;
223  } else if ( one.minPixelCol() < other.minPixelCol() ) {
224  return true;
225  } else {
226  return false;
227  }
228 }
229 
230 
235 
240 
243 #endif
int minPixelCol() const
void setSplitClusterErrorY(float erry)
constexpr int dx() const
bool overflowCol() const
std::vector< PixelDigi >::const_iterator PixelDigiIter
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
std::pair< PixelDigiIter, PixelDigiIter > PixelDigiRange
float x() const
constexpr int row() const
const std::vector< Pixel > pixels() const
#define constexpr
int size() const