CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/DataFormats/SiPixelCluster/interface/SiPixelCluster.h

Go to the documentation of this file.
00001 #ifndef DataFormats_SiPixel_Cluster_SiPixelCluster_h
00002 #define DataFormats_SiPixel_Cluster_SiPixelCluster_h
00003 
00004 //---------------------------------------------------------------------------
00018 //---------------------------------------------------------------------------
00019 
00020 #include <vector>
00021 #include "FWCore/Utilities/interface/GCC11Compatibility.h"
00022 
00023 #ifndef CMS_NOCXX11
00024 #include <cstdint>
00025 #else
00026 #include "boost/cstdint.hpp"
00027 #endif
00028 
00029 #include <cassert>
00030 
00031 class PixelDigi;
00032 
00033 class SiPixelCluster {
00034 public:
00035   
00036   class Pixel {
00037   public:
00038     constexpr Pixel() : x(0), y(0), adc(0){} // for root
00039     constexpr Pixel(int pix_x, int pix_y, int pix_adc) :
00040       x(pix_x), y(pix_y), adc(pix_adc) {}
00041     uint16_t  x;
00042     uint16_t y;
00043     uint16_t adc;
00044   };
00045   
00046   //--- Integer shift in x and y directions.
00047   class Shift {
00048   public:
00049     constexpr Shift( int dx, int dy) : dx_(dx), dy_(dy) {}
00050     constexpr Shift() : dx_(0), dy_(0) {}
00051     constexpr int dx() const { return dx_;}
00052     constexpr int dy() const { return dy_;}
00053   private:
00054     int dx_;
00055     int dy_;
00056   };
00057   
00058   //--- Position of a SiPixel
00059   class PixelPos {
00060   public:
00061     constexpr PixelPos() : row_(0), col_(0) {}
00062     constexpr PixelPos(int row, int col) : row_(row) , col_(col) {}
00063     constexpr int row() const { return row_;}
00064     constexpr int col() const { return col_;}
00065     constexpr PixelPos operator+( const Shift& shift) {
00066       return PixelPos( row() + shift.dx(), col() + shift.dy());
00067     }
00068   private:
00069     int row_;
00070     int col_;
00071   };
00072   
00073   typedef std::vector<PixelDigi>::const_iterator   PixelDigiIter;
00074   typedef std::pair<PixelDigiIter,PixelDigiIter>   PixelDigiRange;
00075   
00076   
00077 #ifndef CMS_NOCXX11
00078   static constexpr unsigned int POSBITS=10;
00079   static constexpr unsigned int SPANBITS=6;
00080   static constexpr unsigned int MAXSPAN=63;
00081   static constexpr unsigned int MAXPOS=1023;
00082 #else
00083   static const unsigned int POSBITS=10;
00084   static const unsigned int SPANBITS=6;
00085   static const unsigned int MAXSPAN=63;
00086   static const unsigned int MAXPOS=1023;
00087 #endif  
00088   
00093   SiPixelCluster() : thePixelRow(MAXPOS), thePixelCol(MAXPOS), err_x(-99999.9), err_y(-99999.9) {}  // needed by many....
00094   
00095   SiPixelCluster(unsigned int isize, uint16_t const * adcs,
00096                  uint16_t const * xpos,  uint16_t const * ypos, 
00097                  uint16_t const  xmin,  uint16_t const  ymin) :   
00098     thePixelOffset(2*isize), thePixelADC(adcs,adcs+isize), err_x(-99999.9), err_y(-99999.9) {
00099     uint16_t maxCol = 0;
00100     uint16_t maxRow = 0;
00101     for (unsigned int i=0; i!=isize; ++i) {
00102       uint16_t xoffset = xpos[i]-xmin;
00103       uint16_t yoffset = ypos[i]-ymin;
00104       thePixelOffset[i*2] = std::min(uint16_t(MAXSPAN),xoffset);
00105       thePixelOffset[i*2+1] = std::min(uint16_t(MAXSPAN),yoffset);
00106       if (xoffset > maxRow) maxRow = xoffset; 
00107       if (yoffset > maxCol) maxCol = yoffset; 
00108     }
00109     packRow(xmin,maxRow);
00110     packCol(ymin,maxCol);
00111   }
00112   
00113   
00114   // obsolete (only for regression tests)
00115   SiPixelCluster( const PixelPos& pix, int adc);
00116   void add( const PixelPos& pix, int adc);
00117   
00118   // Analog linear average position (barycenter) 
00119   float x() const {
00120     float qm = 0.0;
00121     int isize = thePixelADC.size();
00122     for (int i=0; i<isize; ++i)
00123       qm += float(thePixelADC[i]) * (thePixelOffset[i*2] + minPixelRow() + 0.5f);
00124     return qm/charge();
00125   }
00126   
00127   float y() const {
00128     float qm = 0.0;
00129     int isize = thePixelADC.size();
00130     for (int i=0; i<isize; ++i)
00131       qm += float(thePixelADC[i]) * (thePixelOffset[i*2+1]  + minPixelCol() + 0.5f);
00132     return qm/charge();
00133   }
00134   
00135   // Return number of pixels.
00136   int size() const { return thePixelADC.size();}
00137   
00138   // Return cluster dimension in the x direction.
00139   int sizeX() const {verifyVersion(); return rowSpan() +1;}
00140   
00141   // Return cluster dimension in the y direction.
00142   int sizeY() const {verifyVersion(); return colSpan() +1;}
00143   
00144   
00145   inline float charge() const {
00146     float qm = 0.0;
00147     int isize = thePixelADC.size();
00148     for (int i=0; i<isize; ++i) 
00149       qm += float(thePixelADC[i]);
00150     return qm;
00151   } // Return total cluster charge.
00152   
00153   inline int minPixelRow() const { return thePixelRow&MAXPOS;} // The min x index.
00154   inline int maxPixelRow() const { verifyVersion(); return minPixelRow() + rowSpan();} // The max x index.
00155   inline int minPixelCol() const { return thePixelCol&MAXPOS;} // The min y index.
00156   inline int maxPixelCol() const { verifyVersion(); return minPixelCol() + colSpan();} // The max y index.
00157   
00158   
00159   const std::vector<uint8_t> & pixelOffset() const { return thePixelOffset;}
00160   const std::vector<uint16_t> & pixelADC() const { return thePixelADC;}
00161   
00162   // obsolete, use single pixel access below
00163   const std::vector<Pixel> pixels() const {
00164     std::vector<Pixel> oldPixVector;
00165     int isize = thePixelADC.size();
00166     oldPixVector.reserve(isize); 
00167     for(int i=0; i<isize; ++i) {
00168       oldPixVector.push_back(pixel(i));
00169     }
00170     return oldPixVector;
00171   }
00172   
00173   // infinite faster than above...
00174   Pixel pixel(int i) const {
00175     return Pixel(minPixelRow() + thePixelOffset[i*2],
00176                  minPixelCol() + thePixelOffset[i*2+1],
00177                  thePixelADC[i]
00178                  );
00179   }
00180   
00181 private:
00182   
00183   static int span_(uint16_t packed) { return packed >> POSBITS;}
00184   static int overflow_(uint16_t packed) { return span_(packed)==uint16_t(MAXSPAN);}
00185   static uint16_t pack_(uint16_t zmin, unsigned  short zspan) {
00186     zspan = std::min(zspan, uint16_t(MAXSPAN));
00187     return (zspan<<POSBITS) | zmin;
00188   }
00189 public:
00190   
00191   int colSpan() const {return span_(thePixelCol); }
00192   
00193   int rowSpan() const { return span_(thePixelRow); }
00194   
00195   
00196   bool overflowCol() const { return overflow_(thePixelCol); }
00197   
00198   bool overflowRow() const { return overflow_(thePixelRow); }
00199   
00200   bool overflow() const { return  overflowCol() || overflowRow(); }
00201   
00202   void packCol(uint16_t ymin, uint16_t yspan) {
00203     thePixelCol = pack_(ymin,yspan);
00204   }
00205   void packRow(uint16_t xmin, uint16_t xspan) {
00206     thePixelRow = pack_(xmin,xspan);
00207   }
00208   
00209   
00210   
00212   void verifyVersion() const {
00213     if unlikely( thePixelRow<MAXPOS && thePixelCol<MAXPOS)
00214                  const_cast<SiPixelCluster*>(this)->computeMax();
00215   }
00216   
00218   void computeMax()  {
00219     int maxRow = 0;
00220     int maxCol = 0;
00221     int isize  = thePixelADC.size();
00222     for (int i=0; i!=isize; ++i) {
00223       int xsize  = thePixelOffset[i*2];
00224       if (xsize > maxRow) maxRow = xsize;
00225       int ysize = thePixelOffset[i*2+1] ;
00226       if (ysize > maxCol) maxCol = ysize;
00227     }
00228     // assume minimum is correct
00229     int minCol= minPixelCol();
00230     packCol(minCol,maxCol);
00231     int minRow= minPixelRow();
00232     packRow(minRow,maxRow);
00233   }
00234   
00235   // ggiurgiu@fnal.gov, 01/05/12 
00236   // Getters and setters for the newly added data members (err_x and err_y). See below. 
00237   void setSplitClusterErrorX( float errx ) { err_x = errx; }
00238   void setSplitClusterErrorY( float erry ) { err_y = erry; }
00239   float getSplitClusterErrorX() const { return err_x; }
00240   float getSplitClusterErrorY() const { return err_y; }
00241   
00242   
00243 private:
00244   
00245   std::vector<uint8_t>  thePixelOffset;
00246   std::vector<uint16_t> thePixelADC;
00247   
00248   
00249   uint16_t thePixelRow; // Minimum and span pixel index in the x direction (low edge).
00250   uint16_t thePixelCol; // Minimum and span pixel index in the y direction (left edge).
00251   // Need 10 bits for Postion information. the other 6 used for span
00252   
00253   // ggiurgiu@fnal.gov, 01/05/12
00254   // Add cluster errors to be used by rechits from split clusters. 
00255   // A rechit from a split cluster has larger errors than rechits from normal clusters. 
00256   // However, when presented with a cluster, the CPE does not know if the cluster comes 
00257   // from a splitting procedure or not. That's why we have to instruct the CPE to use 
00258   // appropriate errors for split clusters.
00259   // To avoid increase of data size on disk,these new data members are set as transient in: 
00260   // DataFormats/SiPixelCluster/src/classes_def.xml
00261   float err_x;
00262   float err_y;
00263   
00264 };
00265 
00266 // Comparison operators  (no clue...)
00267 inline bool operator<( const SiPixelCluster& one, const SiPixelCluster& other) {
00268   if ( one.minPixelRow() < other.minPixelRow() ) {
00269     return true;
00270   } else if ( one.minPixelRow() > other.minPixelRow() ) {
00271     return false;
00272   } else if ( one.minPixelCol() < other.minPixelCol() ) {
00273     return true;
00274   } else {
00275     return false;
00276   }
00277 }
00278 
00279 #include "DataFormats/Common/interface/DetSetVector.h"
00280 #include "DataFormats/Common/interface/DetSetVectorNew.h"
00281 #include "DataFormats/Common/interface/Ref.h"
00282 #include "DataFormats/Common/interface/DetSetRefVector.h"
00283 
00284 typedef edm::DetSetVector<SiPixelCluster> SiPixelClusterCollection;
00285 typedef edm::Ref<SiPixelClusterCollection, SiPixelCluster> SiPixelClusterRef;
00286 typedef edm::DetSetRefVector<SiPixelCluster> SiPixelClusterRefVector;
00287 typedef edm::RefProd<SiPixelClusterCollection> SiPixelClusterRefProd;
00288 
00289 typedef edmNew::DetSetVector<SiPixelCluster> SiPixelClusterCollectionNew;
00290 typedef edm::Ref<SiPixelClusterCollectionNew, SiPixelCluster> SiPixelClusterRefNew;
00291 #endif