CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 //---------------------------------------------------------------------------
18 //---------------------------------------------------------------------------
19 
20 #include <vector>
22 
23 #ifndef CMS_NOCXX11
24 #include <cstdint>
25 #else
26 #include "boost/cstdint.hpp"
27 #endif
28 
29 #include <cassert>
30 
31 class PixelDigi;
32 
34 public:
35 
36  class Pixel {
37  public:
38  constexpr Pixel() : x(0), y(0), adc(0){} // for root
39  constexpr Pixel(int pix_x, int pix_y, int pix_adc) :
40  x(pix_x), y(pix_y), adc(pix_adc) {}
41  uint16_t x;
42  uint16_t y;
43  uint16_t adc;
44  };
45 
46  //--- Integer shift in x and y directions.
47  class Shift {
48  public:
49  constexpr Shift( int dx, int dy) : dx_(dx), dy_(dy) {}
50  constexpr Shift() : dx_(0), dy_(0) {}
51  constexpr int dx() const { return dx_;}
52  constexpr int dy() const { return dy_;}
53  private:
54  int dx_;
55  int dy_;
56  };
57 
58  //--- Position of a SiPixel
59  class PixelPos {
60  public:
61  constexpr PixelPos() : row_(0), col_(0) {}
62  constexpr PixelPos(int row, int col) : row_(row) , col_(col) {}
63  constexpr int row() const { return row_;}
64  constexpr int col() const { return col_;}
66  return PixelPos( row() + shift.dx(), col() + shift.dy());
67  }
68  private:
69  int row_;
70  int col_;
71  };
72 
73  typedef std::vector<PixelDigi>::const_iterator PixelDigiIter;
74  typedef std::pair<PixelDigiIter,PixelDigiIter> PixelDigiRange;
75 
76 
77 #ifndef CMS_NOCXX11
78  static constexpr unsigned int POSBITS=10;
79  static constexpr unsigned int SPANBITS=6;
80  static constexpr unsigned int MAXSPAN=63;
81  static constexpr unsigned int MAXPOS=1023;
82 #else
83  static const unsigned int POSBITS=10;
84  static const unsigned int SPANBITS=6;
85  static const unsigned int MAXSPAN=63;
86  static const unsigned int MAXPOS=1023;
87 #endif
88 
93  SiPixelCluster() : thePixelRow(MAXPOS), thePixelCol(MAXPOS), err_x(-99999.9), err_y(-99999.9) {} // needed by many....
94 
95  SiPixelCluster(unsigned int isize, uint16_t const * adcs,
96  uint16_t const * xpos, uint16_t const * ypos,
97  uint16_t const xmin, uint16_t const ymin) :
98  thePixelOffset(2*isize), thePixelADC(adcs,adcs+isize), err_x(-99999.9), err_y(-99999.9) {
99  uint16_t maxCol = 0;
100  uint16_t maxRow = 0;
101  for (unsigned int i=0; i!=isize; ++i) {
102  uint16_t xoffset = xpos[i]-xmin;
103  uint16_t yoffset = ypos[i]-ymin;
104  thePixelOffset[i*2] = std::min(uint16_t(MAXSPAN),xoffset);
105  thePixelOffset[i*2+1] = std::min(uint16_t(MAXSPAN),yoffset);
106  if (xoffset > maxRow) maxRow = xoffset;
107  if (yoffset > maxCol) maxCol = yoffset;
108  }
109  packRow(xmin,maxRow);
110  packCol(ymin,maxCol);
111  }
112 
113 
114  // obsolete (only for regression tests)
115  SiPixelCluster( const PixelPos& pix, int adc);
116  void add( const PixelPos& pix, int adc);
117 
118  // Analog linear average position (barycenter)
119  float x() const {
120  float qm = 0.0;
121  int isize = thePixelADC.size();
122  for (int i=0; i<isize; ++i)
123  qm += float(thePixelADC[i]) * (thePixelOffset[i*2] + minPixelRow() + 0.5f);
124  return qm/charge();
125  }
126 
127  float y() const {
128  float qm = 0.0;
129  int isize = thePixelADC.size();
130  for (int i=0; i<isize; ++i)
131  qm += float(thePixelADC[i]) * (thePixelOffset[i*2+1] + minPixelCol() + 0.5f);
132  return qm/charge();
133  }
134 
135  // Return number of pixels.
136  int size() const { return thePixelADC.size();}
137 
138  // Return cluster dimension in the x direction.
139  int sizeX() const {verifyVersion(); return rowSpan() +1;}
140 
141  // Return cluster dimension in the y direction.
142  int sizeY() const {verifyVersion(); return colSpan() +1;}
143 
144 
145  inline float charge() const {
146  float qm = 0.0;
147  int isize = thePixelADC.size();
148  for (int i=0; i<isize; ++i)
149  qm += float(thePixelADC[i]);
150  return qm;
151  } // Return total cluster charge.
152 
153  inline int minPixelRow() const { return thePixelRow&MAXPOS;} // The min x index.
154  inline int maxPixelRow() const { verifyVersion(); return minPixelRow() + rowSpan();} // The max x index.
155  inline int minPixelCol() const { return thePixelCol&MAXPOS;} // The min y index.
156  inline int maxPixelCol() const { verifyVersion(); return minPixelCol() + colSpan();} // The max y index.
157 
158 
159  const std::vector<uint8_t> & pixelOffset() const { return thePixelOffset;}
160  const std::vector<uint16_t> & pixelADC() const { return thePixelADC;}
161 
162  // obsolete, use single pixel access below
163  const std::vector<Pixel> pixels() const {
164  std::vector<Pixel> oldPixVector;
165  int isize = thePixelADC.size();
166  oldPixVector.reserve(isize);
167  for(int i=0; i<isize; ++i) {
168  oldPixVector.push_back(pixel(i));
169  }
170  return oldPixVector;
171  }
172 
173  // infinite faster than above...
174  Pixel pixel(int i) const {
175  return Pixel(minPixelRow() + thePixelOffset[i*2],
176  minPixelCol() + thePixelOffset[i*2+1],
177  thePixelADC[i]
178  );
179  }
180 
181 private:
182 
183  static int span_(uint16_t packed) { return packed >> POSBITS;}
184  static int overflow_(uint16_t packed) { return span_(packed)==uint16_t(MAXSPAN);}
185  static uint16_t pack_(uint16_t zmin, unsigned short zspan) {
186  zspan = std::min(zspan, uint16_t(MAXSPAN));
187  return (zspan<<POSBITS) | zmin;
188  }
189 public:
190 
191  int colSpan() const {return span_(thePixelCol); }
192 
193  int rowSpan() const { return span_(thePixelRow); }
194 
195 
196  bool overflowCol() const { return overflow_(thePixelCol); }
197 
198  bool overflowRow() const { return overflow_(thePixelRow); }
199 
200  bool overflow() const { return overflowCol() || overflowRow(); }
201 
202  void packCol(uint16_t ymin, uint16_t yspan) {
203  thePixelCol = pack_(ymin,yspan);
204  }
205  void packRow(uint16_t xmin, uint16_t xspan) {
206  thePixelRow = pack_(xmin,xspan);
207  }
208 
209 
210 
212  void verifyVersion() const {
214  const_cast<SiPixelCluster*>(this)->computeMax();
215  }
216 
218  void computeMax() {
219  int maxRow = 0;
220  int maxCol = 0;
221  int isize = thePixelADC.size();
222  for (int i=0; i!=isize; ++i) {
223  int xsize = thePixelOffset[i*2];
224  if (xsize > maxRow) maxRow = xsize;
225  int ysize = thePixelOffset[i*2+1] ;
226  if (ysize > maxCol) maxCol = ysize;
227  }
228  // assume minimum is correct
229  int minCol= minPixelCol();
230  packCol(minCol,maxCol);
231  int minRow= minPixelRow();
232  packRow(minRow,maxRow);
233  }
234 
235  // ggiurgiu@fnal.gov, 01/05/12
236  // Getters and setters for the newly added data members (err_x and err_y). See below.
237  void setSplitClusterErrorX( float errx ) { err_x = errx; }
238  void setSplitClusterErrorY( float erry ) { err_y = erry; }
239  float getSplitClusterErrorX() const { return err_x; }
240  float getSplitClusterErrorY() const { return err_y; }
241 
242 
243 private:
244 
245  std::vector<uint8_t> thePixelOffset;
246  std::vector<uint16_t> thePixelADC;
247 
248 
249  uint16_t thePixelRow; // Minimum and span pixel index in the x direction (low edge).
250  uint16_t thePixelCol; // Minimum and span pixel index in the y direction (left edge).
251  // Need 10 bits for Postion information. the other 6 used for span
252 
253  // ggiurgiu@fnal.gov, 01/05/12
254  // Add cluster errors to be used by rechits from split clusters.
255  // A rechit from a split cluster has larger errors than rechits from normal clusters.
256  // However, when presented with a cluster, the CPE does not know if the cluster comes
257  // from a splitting procedure or not. That's why we have to instruct the CPE to use
258  // appropriate errors for split clusters.
259  // To avoid increase of data size on disk,these new data members are set as transient in:
260  // DataFormats/SiPixelCluster/src/classes_def.xml
261  float err_x;
262  float err_y;
263 
264 };
265 
266 // Comparison operators (no clue...)
267 inline bool operator<( const SiPixelCluster& one, const SiPixelCluster& other) {
268  if ( one.minPixelRow() < other.minPixelRow() ) {
269  return true;
270  } else if ( one.minPixelRow() > other.minPixelRow() ) {
271  return false;
272  } else if ( one.minPixelCol() < other.minPixelCol() ) {
273  return true;
274  } else {
275  return false;
276  }
277 }
278 
283 
288 
291 #endif
int adc(sample_type sample)
get the ADC sample (12 bits)
int i
Definition: DBlmapReader.cc:9
float charge() const
int minPixelCol() const
void setSplitClusterErrorY(float erry)
static uint16_t pack_(uint16_t zmin, unsigned short zspan)
bool overflowCol() const
PixelPos(int row, int col)
void verifyVersion() const
mostly to be compatible for &lt;610
std::vector< PixelDigi >::const_iterator PixelDigiIter
int rowSpan() const
static const unsigned int SPANBITS
void computeMax()
moslty to be compatible for &lt;610
const Int_t ysize
edm::RefProd< SiPixelClusterCollection > SiPixelClusterRefProd
#define constexpr
int maxPixelRow() const
uint16_t thePixelCol
#define unlikely(x)
void packRow(uint16_t xmin, uint16_t xspan)
void packCol(uint16_t ymin, uint16_t yspan)
std::vector< uint16_t > thePixelADC
bool operator<(const FedChannelConnection &, const FedChannelConnection &)
int minPixelRow() const
static const unsigned int POSBITS
PixelPos operator+(const Shift &shift) const
bool overflow() const
int colSpan() const
SiPixelCluster(unsigned int isize, uint16_t const *adcs, uint16_t const *xpos, uint16_t const *ypos, uint16_t const xmin, uint16_t const ymin)
static int span_(uint16_t packed)
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
static const unsigned int MAXPOS
void setSplitClusterErrorX(float errx)
float getSplitClusterErrorY() const
edm::DetSetVector< SiPixelCluster > SiPixelClusterCollection
Shift(int dx, int dy)
edm::Ref< SiPixelClusterCollectionNew, SiPixelCluster > SiPixelClusterRefNew
Pixel(int pix_x, int pix_y, int pix_adc)
int maxPixelCol() const
edmNew::DetSetVector< SiPixelCluster > SiPixelClusterCollectionNew
int sizeY() const
Pixel cluster – collection of neighboring pixels above threshold.
Pixel pixel(int i) const
static int overflow_(uint16_t packed)
float getSplitClusterErrorX() const
static unsigned int const shift
float y() const
const Int_t xsize
int sizeX() const
std::pair< PixelDigiIter, PixelDigiIter > PixelDigiRange
float x() const
uint16_t thePixelRow
const std::vector< Pixel > pixels() const
static const unsigned int MAXSPAN
int size() const