CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/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 //---------------------------------------------------------------------------
00016 //---------------------------------------------------------------------------
00017 
00018 #include <vector>
00019 #include "boost/cstdint.hpp"
00020 
00021 class PixelDigi;
00022 
00023 class SiPixelCluster {
00024  public:
00025   
00026   class Pixel {
00027   public:
00028     Pixel() {} // for root
00029     Pixel(int pix_x, int pix_y, int pix_adc) :
00030       x(pix_x), y(pix_y), adc(pix_adc) {}
00031     unsigned char  x;
00032     unsigned short y;
00033     unsigned short adc;
00034   };
00035   
00036   //--- Integer shift in x and y directions.
00037   class Shift {
00038   public:
00039     Shift( int dx, int dy) : dx_(dx), dy_(dy) {}
00040     Shift() : dx_(0), dy_(0) {}
00041     int dx() const { return dx_;}
00042     int dy() const { return dy_;}
00043   private:
00044     int dx_;
00045     int dy_;
00046   };
00047   
00048   //--- Position of a SiPixel
00049   class PixelPos {
00050   public:
00051     PixelPos() : row_(0), col_(0) {}
00052     PixelPos(int row, int col) : row_(row) , col_(col) {}
00053     int row() const { return row_;}
00054     int col() const { return col_;}
00055     PixelPos operator+( const Shift& shift) {
00056       return PixelPos( row() + shift.dx(), col() + shift.dy());
00057     }
00058   private:
00059     int row_;
00060     int col_;
00061   };
00062   
00063   typedef std::vector<PixelDigi>::const_iterator   PixelDigiIter;
00064   typedef std::pair<PixelDigiIter,PixelDigiIter>   PixelDigiRange;
00065   
00070   // &&& Decide the fate of the two strip-like constructors below:
00071   SiPixelCluster() : detId_(0) {}  // needed by vector::push_back()!
00072   // SiPixelCluster( unsigned int detid, const PixelDigiRange& range)
00073     
00074   SiPixelCluster( const PixelPos& pix, int adc);
00075   
00076   void add( const PixelPos& pix, int adc);
00077   
00078   // Analog linear average position (barycenter) 
00079   float x() const {
00080                 float qm = 0.0;
00081                 int isize = thePixelADC.size();
00082                 for (int i=0; i<isize; ++i)
00083                         qm += float(thePixelADC[i]) * (thePixelOffset[i*2] + theMinPixelRow + 0.5);
00084                 return qm/charge();
00085                         }
00086   float y() const {
00087                 float qm = 0.0;
00088                 int isize = thePixelADC.size();
00089                 for (int i=0; i<isize; ++i)
00090                         qm += float(thePixelADC[i]) * (thePixelOffset[i*2+1]  + theMinPixelCol + 0.5);
00091                 return qm/charge();
00092         }
00093 
00094   // Return number of pixels.
00095   int size() const { return thePixelADC.size();}
00096 
00097   // Return cluster dimension in the x direction.
00098   int sizeX() const {return maxPixelRow() - theMinPixelRow +1;}
00099 
00100   // Return cluster dimension in the y direction.
00101   int sizeY() const {return maxPixelCol() - theMinPixelCol +1;}
00102 
00103   // Detect clusters at the edge of the detector.
00104   // NOTE: Moved to RectangularPixelTopology class
00105   // bool edgeHitX() const;
00106   // bool edgeHitY() const;
00107 
00108         inline float charge() const {
00109                 float qm = 0.0;
00110                 int isize = thePixelADC.size();
00111                 for (int i=0; i<isize; ++i) 
00112                         qm += float(thePixelADC[i]);
00113                 return qm;
00114         } // Return total cluster charge.
00115 
00116         inline int minPixelRow() const { return theMinPixelRow;} // The min x index.
00117   inline int minPixelCol() const { return theMinPixelCol & 511;} // The min y index.
00118         
00119   inline int maxPixelRow() const {
00120                 int maxRow = 0;
00121                 int isize  = thePixelADC.size();
00122                 for (int i=0; i<isize; ++i) {
00123                         int xsize  = thePixelOffset[i*2];
00124                         if (xsize > maxRow) maxRow = xsize;
00125                 }
00126         return maxRow + theMinPixelRow; // The max x index.
00127         }
00128         
00129         inline int maxPixelCol() const {
00130                 int maxCol = 0;
00131                 int isize = thePixelADC.size();
00132                 for (int i=0; i<isize; ++i) {
00133                         int ysize = thePixelOffset[i*2+1] ;
00134                         if (ysize > maxCol) maxCol = ysize;
00135                 }
00136                 return maxCol + theMinPixelCol; // The max y index.
00137         }
00138   
00139   const std::vector<uint8_t> & pixelOffset() const { return thePixelOffset;}
00140         const std::vector<uint16_t> & pixelADC() const { return thePixelADC;}
00141 
00142         const std::vector<Pixel> pixels() const {
00143                 std::vector<Pixel> oldPixVector;
00144                 int isize = thePixelADC.size();
00145                 oldPixVector.reserve(isize); 
00146                 for(int i=0; i<isize; ++i) {
00147                         int x = theMinPixelRow + (thePixelOffset[i*2]  );
00148                         int y = theMinPixelCol + (thePixelOffset[i*2+1] );
00149                         oldPixVector.push_back(Pixel(x,y,thePixelADC[i]));
00150                 }
00151                 return oldPixVector;
00152         }
00153   //--- Cloned fom Strips:
00154   
00158   unsigned int geographicalId() const {return detId_;}
00159   
00160   // &&& Decide if we still need these two:
00161   // typedef vector<Digi::ChannelType>    ChannelContainer;
00162   // ChannelContainer  channels() const;
00163   
00164  private:
00165   unsigned int         detId_;
00166   
00167   std::vector<uint8_t>  thePixelOffset;
00168         std::vector<uint16_t> thePixelADC;
00169 
00170         /*  float theSumX;  // Sum of charge weighted pixel positions.
00171   float theSumY;
00172   float theCharge;  // Total charge
00173         uint8_t  theMaxPixelRow; // Maximum pixel index in the x direction (top edge).
00174         uint16_t theMaxPixelCol; // Maximum pixel index in the y direction (right edge).
00175         */
00176         uint8_t  theMinPixelRow; // Minimum pixel index in the x direction (low edge).
00177         uint16_t theMinPixelCol; // Minimum pixel index in the y direction (left edge).
00178                                  // Need 9 bits for Col information. Use 1 bit for whether larger
00179                                  // cluster than 9x33. Other 6 bits for quality information.
00180   
00181 };
00182 
00183 // Comparison operators
00184 inline bool operator<( const SiPixelCluster& one, const SiPixelCluster& other) {
00185   if ( one.geographicalId() < other.geographicalId() ) {
00186     return true;
00187   } else if ( one.geographicalId() > other.geographicalId() ) {
00188     return false;
00189   } else if ( one.minPixelRow() < other.minPixelRow() ) {
00190     return true;
00191   } else if ( one.minPixelRow() > other.minPixelRow() ) {
00192     return false;
00193   } else if ( one.minPixelCol() < other.minPixelCol() ) {
00194     return true;
00195   } else {
00196     return false;
00197   }
00198 }
00199 
00200 #include "DataFormats/Common/interface/DetSetVector.h"
00201 #include "DataFormats/Common/interface/DetSetVectorNew.h"
00202 #include "DataFormats/Common/interface/Ref.h"
00203 #include "DataFormats/Common/interface/DetSetRefVector.h"
00204 
00205 typedef edm::DetSetVector<SiPixelCluster> SiPixelClusterCollection;
00206 typedef edm::Ref<SiPixelClusterCollection, SiPixelCluster> SiPixelClusterRef;
00207 typedef edm::DetSetRefVector<SiPixelCluster> SiPixelClusterRefVector;
00208 typedef edm::RefProd<SiPixelClusterCollection> SiPixelClusterRefProd;
00209 
00210 typedef edmNew::DetSetVector<SiPixelCluster> SiPixelClusterCollectionNew;
00211 typedef edm::Ref<SiPixelClusterCollectionNew, SiPixelCluster> SiPixelClusterRefNew;
00212 #endif