CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Static Public Attributes | Static Private Member Functions | Private Attributes
Phase2ITPixelCluster Class Reference

Pixel cluster – collection of neighboring pixels above threshold. More...

#include <Phase2ITPixelCluster.h>

Classes

class  Pixel
 
class  PixelPos
 
class  Shift
 

Public Types

typedef std::vector< PixelDigi >::const_iterator PixelDigiIter
 
typedef std::pair< PixelDigiIter, PixelDigiIterPixelDigiRange
 

Public Member Functions

void add (const PixelPos &pix, uint32_t adc)
 
float charge () const
 
int colSpan () const
 
void computeMax ()
 moslty to be compatible for <610 More...
 
float getSplitClusterErrorX () const
 
float getSplitClusterErrorY () const
 
uint32_t maxPixelCol () const
 
uint32_t maxPixelRow () const
 
uint32_t minPixelCol () const
 
uint32_t minPixelRow () const
 
bool overflow () const
 
bool overflowCol () const
 
bool overflowRow () const
 
void packCol (uint32_t ymin, uint32_t yspan)
 
void packRow (uint32_t xmin, uint32_t xspan)
 
 Phase2ITPixelCluster ()
 
 Phase2ITPixelCluster (unsigned int isize, uint32_t const *adcs, uint32_t const *xpos, uint32_t const *ypos, uint32_t const xmin, uint32_t const ymin)
 
 Phase2ITPixelCluster (const PixelPos &pix, uint32_t adc)
 
Pixel pixel (int i) const
 
const std::vector< uint32_t > & pixelADC () const
 
const std::vector< uint16_t > & pixelOffset () const
 
const std::vector< Pixelpixels () const
 
int rowSpan () const
 
void setSplitClusterErrorX (float errx)
 
void setSplitClusterErrorY (float erry)
 
int size () const
 
int sizeX () const
 
int sizeY () const
 
void verifyVersion () const
 mostly to be compatible for <610 More...
 
float x () const
 
float y () const
 

Static Public Attributes

static unsigned int MAXPOS =2047
 
static unsigned int MAXSPAN =255
 
static unsigned int POSBITS =20
 
static unsigned int SPANBITS =12
 

Static Private Member Functions

static int overflow_ (uint32_t packed)
 
static uint32_t pack_ (uint32_t zmin, unsigned int zspan)
 
static int span_ (uint32_t packed)
 

Private Attributes

float err_x
 
float err_y
 
std::vector< uint32_t > thePixelADC
 
uint32_t thePixelCol
 
std::vector< uint16_t > thePixelOffset
 
uint32_t thePixelRow
 

Detailed Description

Pixel cluster – collection of neighboring pixels above threshold.

Pixel cluster – collection of pixels with ADC counts.

Class to contain and store all the topological information of pixel clusters: charge, global size, size and the barycenter in x and y local directions. It builds a vector of Pixel (which is an inner class) and a container of channels.

Class to contain and store all the topological information of pixel clusters: charge, global size, size and the barycenter in x and y local directions. It builds a vector of SiPixel (which is an inner class) and a container of channels.

March 2007: Edge pixel methods moved to RectangularPixelTopology (V.Chiochia) May 2008: Offset based packing (D.Fehling / A. Rizzi)

Author
Petar Maksimovic, JHU

Definition at line 23 of file Phase2ITPixelCluster.h.

Member Typedef Documentation

typedef std::vector<PixelDigi>::const_iterator Phase2ITPixelCluster::PixelDigiIter

Definition at line 63 of file Phase2ITPixelCluster.h.

Definition at line 64 of file Phase2ITPixelCluster.h.

Constructor & Destructor Documentation

Phase2ITPixelCluster::Phase2ITPixelCluster ( )
inline

Construct from a range of digis that form a cluster and from a DetID. The range is assumed to be non-empty.

Definition at line 75 of file Phase2ITPixelCluster.h.

Referenced by Phase2ITPixelCluster().

75 : thePixelRow(MAXPOS), thePixelCol(MAXPOS), err_x(-99999.9), err_y(-99999.9) {} // needed by many....
static unsigned int MAXPOS
Phase2ITPixelCluster::Phase2ITPixelCluster ( unsigned int  isize,
uint32_t const *  adcs,
uint32_t const *  xpos,
uint32_t const *  ypos,
uint32_t const  xmin,
uint32_t const  ymin 
)
inline

Definition at line 77 of file Phase2ITPixelCluster.h.

References Phase2ITPixelCluster::Pixel::adc, add(), i, min(), packCol(), packRow(), Phase2ITPixelCluster(), thePixelOffset, TrackerOfflineValidation_Dqm_cff::xmin, and Phase2TrackerMonitorDigi_cff::ymin.

79  :
80  thePixelOffset(2*isize), thePixelADC(adcs,adcs+isize), err_x(-99999.9), err_y(-99999.9) {
81  uint32_t maxCol = 0;
82  uint32_t maxRow = 0;
83  for (unsigned int i=0; i!=isize; ++i) {
84  uint32_t xoffset = xpos[i]-xmin;
85  uint32_t yoffset = ypos[i]-ymin;
86  thePixelOffset[i*2] = std::min(uint32_t(MAXSPAN),xoffset);
87  thePixelOffset[i*2+1] = std::min(uint32_t(MAXSPAN),yoffset);
88  if (xoffset > maxRow) maxRow = xoffset;
89  if (yoffset > maxCol) maxCol = yoffset;
90  }
91  packRow(xmin,maxRow);
92  packCol(ymin,maxCol);
93  }
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > thePixelADC
void packRow(uint32_t xmin, uint32_t xspan)
T min(T a, T b)
Definition: MathUtil.h:58
void packCol(uint32_t ymin, uint32_t yspan)
std::vector< uint16_t > thePixelOffset
static unsigned int MAXSPAN
Phase2ITPixelCluster::Phase2ITPixelCluster ( const PixelPos pix,
uint32_t  adc 
)

Definition at line 19 of file Phase2ITPixelCluster.cc.

References thePixelADC, and thePixelOffset.

19  :
20  thePixelRow(pix.row()),
21  thePixelCol(pix.col()),
22  // ggiurgiu@fnal.gov, 01/05/12
23  // Initialize the split cluster errors to un-physical values.
24  // The CPE will check these errors and if they are not un-physical,
25  // it will recognize the clusters as split and assign these (increased)
26  // errors to the corresponding rechit.
27  err_x(-99999.9),
28  err_y(-99999.9)
29 {
30  // First pixel in this cluster.
31  thePixelADC.push_back( adc );
32  thePixelOffset.push_back(0 );
33  thePixelOffset.push_back(0 );
34 }
int adc(sample_type sample)
get the ADC sample (12 bits)
std::vector< uint32_t > thePixelADC
std::vector< uint16_t > thePixelOffset

Member Function Documentation

void Phase2ITPixelCluster::add ( const PixelPos pix,
uint32_t  adc 
)

Definition at line 36 of file Phase2ITPixelCluster.cc.

References Phase2ITPixelCluster::PixelPos::col(), i, maxPixelCol(), maxPixelRow(), MAXSPAN, min(), minPixelCol(), minPixelRow(), overflowCol(), overflowRow(), packCol(), packRow(), Phase2ITPixelCluster::PixelPos::row(), thePixelADC, and thePixelOffset.

Referenced by Phase2ITPixelThresholdClusterizer::make_cluster(), Phase2ITPixelCluster(), and counter.Counter::register().

36  {
37 
38  uint32_t ominRow = minPixelRow();
39  uint32_t ominCol = minPixelCol();
40  bool recalculate = false;
41 
42  uint32_t minRow = ominRow;
43  uint32_t minCol = ominCol;
44 
45  if (pix.row() < minRow) {
46  minRow = pix.row();
47  recalculate = true;
48  }
49  if (pix.col() < minCol) {
50  minCol = pix.col();
51  recalculate = true;
52  }
53 
54  if (recalculate) {
55  int maxCol = 0;
56  int maxRow = 0;
57  int isize = thePixelADC.size();
58  for (int i=0; i<isize; ++i) {
59  int xoffset = thePixelOffset[i*2] + ominRow - minRow;
60  int yoffset = thePixelOffset[i*2+1] + ominCol -minCol;
61  thePixelOffset[i*2] = std::min(int(MAXSPAN),xoffset);
62  thePixelOffset[i*2+1] = std::min(int(MAXSPAN),yoffset);
63  if (xoffset > maxRow) maxRow = xoffset;
64  if (yoffset > maxCol) maxCol = yoffset;
65  }
66  packRow(minRow,maxRow);
67  packCol(minCol,maxCol);
68  }
69 
70  if ( (!overflowRow()) && pix.row() > maxPixelRow())
71  packRow(minRow,pix.row()-minRow);
72 
73  if ( (!overflowCol()) && pix.col() > maxPixelCol())
74  packCol(minCol,pix.col()-minCol);
75 
76  thePixelADC.push_back( adc );
77  thePixelOffset.push_back( std::min(MAXSPAN,pix.row() - minRow) );
78  thePixelOffset.push_back( std::min(MAXSPAN,pix.col() - minCol) );
79 }
int adc(sample_type sample)
get the ADC sample (12 bits)
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > thePixelADC
uint32_t minPixelRow() const
uint32_t maxPixelCol() const
void packRow(uint32_t xmin, uint32_t xspan)
uint32_t maxPixelRow() const
T min(T a, T b)
Definition: MathUtil.h:58
uint32_t minPixelCol() const
void packCol(uint32_t ymin, uint32_t yspan)
std::vector< uint16_t > thePixelOffset
static unsigned int MAXSPAN
float Phase2ITPixelCluster::charge ( ) const
inline

Definition at line 127 of file Phase2ITPixelCluster.h.

References i, and thePixelADC.

Referenced by Phase2ITPixelThresholdClusterizer::clusterizeDetUnit(), Phase2ITPixelThresholdClusterizer::make_cluster(), x(), and y().

127  {
128  float qm = 0.0;
129  int isize = thePixelADC.size();
130  for (int i=0; i<isize; ++i)
131  qm += float(thePixelADC[i]);
132  return qm;
133  } // Return total cluster charge.
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > thePixelADC
int Phase2ITPixelCluster::colSpan ( ) const
inline

Definition at line 173 of file Phase2ITPixelCluster.h.

References span_(), and thePixelCol.

Referenced by maxPixelCol(), and sizeY().

173 {return span_(thePixelCol); }
static int span_(uint32_t packed)
void Phase2ITPixelCluster::computeMax ( )
inline

moslty to be compatible for <610

Definition at line 200 of file Phase2ITPixelCluster.h.

References i, minPixelCol(), minPixelRow(), packCol(), packRow(), thePixelADC, thePixelOffset, xsize, and ysize.

Referenced by verifyVersion().

200  {
201  int maxRow = 0;
202  int maxCol = 0;
203  int isize = thePixelADC.size();
204  for (int i=0; i!=isize; ++i) {
205  int xsize = thePixelOffset[i*2];
206  if (xsize > maxRow) maxRow = xsize;
207  int ysize = thePixelOffset[i*2+1] ;
208  if (ysize > maxCol) maxCol = ysize;
209  }
210  // assume minimum is correct
211  uint32_t minCol= minPixelCol();
212  packCol(minCol,maxCol);
213  uint32_t minRow= minPixelRow();
214  packRow(minRow,maxRow);
215  }
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > thePixelADC
uint32_t minPixelRow() const
void packRow(uint32_t xmin, uint32_t xspan)
const Int_t ysize
uint32_t minPixelCol() const
void packCol(uint32_t ymin, uint32_t yspan)
std::vector< uint16_t > thePixelOffset
const Int_t xsize
float Phase2ITPixelCluster::getSplitClusterErrorX ( ) const
inline

Definition at line 220 of file Phase2ITPixelCluster.h.

References err_x.

220 { return err_x; }
float Phase2ITPixelCluster::getSplitClusterErrorY ( ) const
inline

Definition at line 221 of file Phase2ITPixelCluster.h.

References err_y.

221 { return err_y; }
uint32_t Phase2ITPixelCluster::maxPixelCol ( ) const
inline

Definition at line 138 of file Phase2ITPixelCluster.h.

References colSpan(), minPixelCol(), and verifyVersion().

Referenced by add().

138 { verifyVersion(); return minPixelCol() + colSpan();} // The max y index.
void verifyVersion() const
mostly to be compatible for <610
uint32_t minPixelCol() const
uint32_t Phase2ITPixelCluster::maxPixelRow ( ) const
inline

Definition at line 136 of file Phase2ITPixelCluster.h.

References minPixelRow(), rowSpan(), and verifyVersion().

Referenced by add().

136 { verifyVersion(); return minPixelRow() + rowSpan();} // The max x index.
uint32_t minPixelRow() const
void verifyVersion() const
mostly to be compatible for <610
uint32_t Phase2ITPixelCluster::minPixelCol ( ) const
inline

Definition at line 137 of file Phase2ITPixelCluster.h.

References MAXPOS, and thePixelCol.

Referenced by add(), computeMax(), maxPixelCol(), operator<(), pixel(), and y().

137 { return thePixelCol&MAXPOS;} // The min y index.
static unsigned int MAXPOS
uint32_t Phase2ITPixelCluster::minPixelRow ( ) const
inline

Definition at line 135 of file Phase2ITPixelCluster.h.

References MAXPOS, and thePixelRow.

Referenced by add(), computeMax(), maxPixelRow(), operator<(), pixel(), and x().

135 { return thePixelRow&MAXPOS;} // The min x index.
static unsigned int MAXPOS
bool Phase2ITPixelCluster::overflow ( ) const
inline

Definition at line 182 of file Phase2ITPixelCluster.h.

References overflowCol(), and overflowRow().

182 { return overflowCol() || overflowRow(); }
static int Phase2ITPixelCluster::overflow_ ( uint32_t  packed)
inlinestaticprivate

Definition at line 166 of file Phase2ITPixelCluster.h.

References span_().

Referenced by overflowCol(), and overflowRow().

166 { return span_(packed)==uint32_t(MAXSPAN);}
static int span_(uint32_t packed)
static unsigned int MAXSPAN
bool Phase2ITPixelCluster::overflowCol ( ) const
inline

Definition at line 178 of file Phase2ITPixelCluster.h.

References overflow_(), and thePixelCol.

Referenced by add(), and overflow().

178 { return overflow_(thePixelCol); }
static int overflow_(uint32_t packed)
bool Phase2ITPixelCluster::overflowRow ( ) const
inline

Definition at line 180 of file Phase2ITPixelCluster.h.

References overflow_(), and thePixelRow.

Referenced by add(), and overflow().

180 { return overflow_(thePixelRow); }
static int overflow_(uint32_t packed)
static uint32_t Phase2ITPixelCluster::pack_ ( uint32_t  zmin,
unsigned int  zspan 
)
inlinestaticprivate

Definition at line 167 of file Phase2ITPixelCluster.h.

References min().

Referenced by packCol(), and packRow().

167  {
168  zspan = std::min(zspan, uint32_t(MAXSPAN));
169  return (zspan<<POSBITS) | zmin;
170  }
static unsigned int POSBITS
T min(T a, T b)
Definition: MathUtil.h:58
static unsigned int MAXSPAN
void Phase2ITPixelCluster::packCol ( uint32_t  ymin,
uint32_t  yspan 
)
inline

Definition at line 184 of file Phase2ITPixelCluster.h.

References pack_(), and thePixelCol.

Referenced by add(), computeMax(), and Phase2ITPixelCluster().

184  {
185  thePixelCol = pack_(ymin,yspan);
186  }
static uint32_t pack_(uint32_t zmin, unsigned int zspan)
void Phase2ITPixelCluster::packRow ( uint32_t  xmin,
uint32_t  xspan 
)
inline

Definition at line 187 of file Phase2ITPixelCluster.h.

References pack_(), and thePixelRow.

Referenced by add(), computeMax(), and Phase2ITPixelCluster().

187  {
188  thePixelRow = pack_(xmin,xspan);
189  }
static uint32_t pack_(uint32_t zmin, unsigned int zspan)
Pixel Phase2ITPixelCluster::pixel ( int  i) const
inline

Definition at line 156 of file Phase2ITPixelCluster.h.

References minPixelCol(), minPixelRow(), Phase2ITPixelCluster::Pixel::Pixel(), thePixelADC, and thePixelOffset.

Referenced by pixels().

156  {
157  return Pixel(minPixelRow() + thePixelOffset[i*2],
158  minPixelCol() + thePixelOffset[i*2+1],
159  thePixelADC[i]
160  );
161  }
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > thePixelADC
uint32_t minPixelRow() const
uint32_t minPixelCol() const
std::vector< uint16_t > thePixelOffset
const std::vector<uint32_t>& Phase2ITPixelCluster::pixelADC ( ) const
inline

Definition at line 142 of file Phase2ITPixelCluster.h.

References thePixelADC.

142 { return thePixelADC;}
std::vector< uint32_t > thePixelADC
const std::vector<uint16_t>& Phase2ITPixelCluster::pixelOffset ( ) const
inline

Definition at line 141 of file Phase2ITPixelCluster.h.

References thePixelOffset.

141 { return thePixelOffset;}
std::vector< uint16_t > thePixelOffset
const std::vector<Pixel> Phase2ITPixelCluster::pixels ( ) const
inline

Definition at line 145 of file Phase2ITPixelCluster.h.

References i, pixel(), and thePixelADC.

Referenced by Phase2ITPixelThresholdClusterizer::make_cluster().

145  {
146  std::vector<Pixel> oldPixVector;
147  int isize = thePixelADC.size();
148  oldPixVector.reserve(isize);
149  for(int i=0; i<isize; ++i) {
150  oldPixVector.push_back(pixel(i));
151  }
152  return oldPixVector;
153  }
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > thePixelADC
Pixel pixel(int i) const
int Phase2ITPixelCluster::rowSpan ( ) const
inline

Definition at line 175 of file Phase2ITPixelCluster.h.

References span_(), and thePixelRow.

Referenced by maxPixelRow(), and sizeX().

175 { return span_(thePixelRow); }
static int span_(uint32_t packed)
void Phase2ITPixelCluster::setSplitClusterErrorX ( float  errx)
inline

Definition at line 218 of file Phase2ITPixelCluster.h.

References err_x.

218 { err_x = errx; }
void Phase2ITPixelCluster::setSplitClusterErrorY ( float  erry)
inline

Definition at line 219 of file Phase2ITPixelCluster.h.

References err_y.

219 { err_y = erry; }
int Phase2ITPixelCluster::size ( void  ) const
inline

Definition at line 118 of file Phase2ITPixelCluster.h.

References thePixelADC.

Referenced by ntuple._Collection::__iter__(), and ntuple._Collection::__len__().

118 { return thePixelADC.size();}
std::vector< uint32_t > thePixelADC
int Phase2ITPixelCluster::sizeX ( ) const
inline

Definition at line 121 of file Phase2ITPixelCluster.h.

References rowSpan(), and verifyVersion().

121 {verifyVersion(); return rowSpan() +1;}
void verifyVersion() const
mostly to be compatible for <610
int Phase2ITPixelCluster::sizeY ( ) const
inline

Definition at line 124 of file Phase2ITPixelCluster.h.

References colSpan(), and verifyVersion().

124 {verifyVersion(); return colSpan() +1;}
void verifyVersion() const
mostly to be compatible for <610
static int Phase2ITPixelCluster::span_ ( uint32_t  packed)
inlinestaticprivate

Definition at line 165 of file Phase2ITPixelCluster.h.

References POSBITS.

Referenced by colSpan(), overflow_(), and rowSpan().

165 { return packed >> POSBITS;}
static unsigned int POSBITS
void Phase2ITPixelCluster::verifyVersion ( ) const
inline

mostly to be compatible for <610

Definition at line 194 of file Phase2ITPixelCluster.h.

References computeMax(), thePixelCol, thePixelRow, and unlikely.

Referenced by maxPixelCol(), maxPixelRow(), sizeX(), and sizeY().

194  {
196  const_cast<Phase2ITPixelCluster*>(this)->computeMax();
197  }
void computeMax()
moslty to be compatible for <610
#define unlikely(x)
static unsigned int MAXPOS
Pixel cluster – collection of neighboring pixels above threshold.
float Phase2ITPixelCluster::x ( ) const
inline

Definition at line 101 of file Phase2ITPixelCluster.h.

References charge(), i, minPixelRow(), thePixelADC, and thePixelOffset.

Referenced by svgfig.Curve.Sample::__repr__(), svgfig.Ellipse::__repr__(), Vispa.Gui.WidgetContainer.WidgetContainer::autosize(), Vispa.Gui.VispaWidget.VispaWidget::boundingRect(), geometryXMLparser.Alignable::pos(), and Vispa.Gui.ConnectableWidget.ConnectableWidget::positionizeMenuWidget().

101  {
102  float qm = 0.0;
103  int isize = thePixelADC.size();
104  for (int i=0; i<isize; ++i)
105  qm += float(thePixelADC[i]) * (thePixelOffset[i*2] + minPixelRow() + 0.5f);
106  return qm/charge();
107  }
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > thePixelADC
uint32_t minPixelRow() const
std::vector< uint16_t > thePixelOffset
float Phase2ITPixelCluster::y ( ) const
inline

Definition at line 109 of file Phase2ITPixelCluster.h.

References charge(), i, minPixelCol(), thePixelADC, and thePixelOffset.

Referenced by svgfig.Ellipse::__repr__(), Vispa.Gui.WidgetContainer.WidgetContainer::autosize(), Vispa.Gui.VispaWidget.VispaWidget::boundingRect(), geometryXMLparser.Alignable::pos(), and Vispa.Gui.ConnectableWidget.ConnectableWidget::positionizeMenuWidget().

109  {
110  float qm = 0.0;
111  int isize = thePixelADC.size();
112  for (int i=0; i<isize; ++i)
113  qm += float(thePixelADC[i]) * (thePixelOffset[i*2+1] + minPixelCol() + 0.5f);
114  return qm/charge();
115  }
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > thePixelADC
uint32_t minPixelCol() const
std::vector< uint16_t > thePixelOffset

Member Data Documentation

float Phase2ITPixelCluster::err_x
private

Definition at line 238 of file Phase2ITPixelCluster.h.

Referenced by getSplitClusterErrorX(), and setSplitClusterErrorX().

float Phase2ITPixelCluster::err_y
private

Definition at line 239 of file Phase2ITPixelCluster.h.

Referenced by getSplitClusterErrorY(), and setSplitClusterErrorY().

unsigned int Phase2ITPixelCluster::MAXPOS =2047
static

Definition at line 69 of file Phase2ITPixelCluster.h.

Referenced by minPixelCol(), and minPixelRow().

unsigned int Phase2ITPixelCluster::MAXSPAN =255
static

Definition at line 68 of file Phase2ITPixelCluster.h.

Referenced by add().

unsigned int Phase2ITPixelCluster::POSBITS =20
static

Definition at line 66 of file Phase2ITPixelCluster.h.

Referenced by span_().

unsigned int Phase2ITPixelCluster::SPANBITS =12
static

Definition at line 67 of file Phase2ITPixelCluster.h.

std::vector<uint32_t> Phase2ITPixelCluster::thePixelADC
private
uint32_t Phase2ITPixelCluster::thePixelCol
private

Definition at line 231 of file Phase2ITPixelCluster.h.

Referenced by colSpan(), minPixelCol(), overflowCol(), packCol(), and verifyVersion().

std::vector<uint16_t> Phase2ITPixelCluster::thePixelOffset
private

Definition at line 226 of file Phase2ITPixelCluster.h.

Referenced by add(), computeMax(), Phase2ITPixelCluster(), pixel(), pixelOffset(), x(), and y().

uint32_t Phase2ITPixelCluster::thePixelRow
private

Definition at line 230 of file Phase2ITPixelCluster.h.

Referenced by minPixelRow(), overflowRow(), packRow(), rowSpan(), and verifyVersion().