Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:05:11

0001 #ifndef DataFormats_SiPixel_Cluster_SiPixelCluster_h
0002 #define DataFormats_SiPixel_Cluster_SiPixelCluster_h
0003 
0004 //---------------------------------------------------------------------------
0005 //!  \class SiPixelCluster
0006 //!  \brief Pixel cluster -- collection of neighboring pixels above threshold
0007 //!
0008 //!  Class to contain and store all the topological information of pixel clusters:
0009 //!  charge, global size, size and the barycenter in x and y
0010 //!  local directions. It builds a vector of SiPixel (which is
0011 //!  an inner class) and a container of channels.
0012 //!
0013 //!  March 2007: Edge methods moved to RectangularPixelTopology class (V.Chiochia)
0014 //!  Feb 2008: Modify the Pixel class from float to shorts
0015 //!  May   2008: Offset based packing (D.Fehling / A. Rizzi)
0016 //!  Sep 2012: added Max back, removed detId (V.I.)
0017 //!  sizeX and sizeY now clipped at 127
0018 //!  July 2017 make it compatible with PhaseII, remove errs....
0019 //---------------------------------------------------------------------------
0020 
0021 #include <vector>
0022 #include <cstdint>
0023 #include <cassert>
0024 #include <limits>
0025 
0026 class PixelDigi;
0027 
0028 class SiPixelCluster {
0029 public:
0030   // FIXME make it POD
0031   class Pixel {
0032   public:
0033     constexpr Pixel() : x(0), y(0), adc(0) {}  // for root
0034     constexpr Pixel(int pix_x, int pix_y, int pix_adc) : x(pix_x), y(pix_y), adc(pix_adc) {}
0035     uint16_t x;
0036     uint16_t y;
0037     uint16_t adc;
0038   };
0039 
0040   //--- Integer shift in x and y directions.
0041   // FIXME make    it POD
0042   class Shift {
0043   public:
0044     constexpr Shift(int dx, int dy) : dx_(dx), dy_(dy) {}
0045     constexpr Shift() : dx_(0), dy_(0) {}
0046     constexpr int dx() const { return dx_; }
0047     constexpr int dy() const { return dy_; }
0048 
0049   private:
0050     int dx_;
0051     int dy_;
0052   };
0053 
0054   //--- Position of a SiPixel
0055   // FIXME make    it POD
0056   class PixelPos {
0057   public:
0058     constexpr PixelPos() : row_(0), col_(0) {}
0059     constexpr PixelPos(int row, int col) : row_(row), col_(col) {}
0060     constexpr int row() const { return row_; }
0061     constexpr int col() const { return col_; }
0062     constexpr PixelPos operator+(const Shift& shift) const { return PixelPos(row() + shift.dx(), col() + shift.dy()); }
0063     constexpr bool operator==(const PixelPos& pos) const { return (row() == pos.row() && col() == pos.col()); }
0064 
0065   private:
0066     int row_;
0067     int col_;
0068   };
0069 
0070   typedef std::vector<PixelDigi>::const_iterator PixelDigiIter;
0071   typedef std::pair<PixelDigiIter, PixelDigiIter> PixelDigiRange;
0072 
0073   static constexpr unsigned int MAXSPAN = 255;
0074   static constexpr unsigned int MAXPOS = 2047;
0075 
0076   static constexpr uint16_t invalidClusterId = std::numeric_limits<uint16_t>::max();
0077 
0078   /** Construct from a range of digis that form a cluster and from 
0079    *  a DetID. The range is assumed to be non-empty.
0080    */
0081 
0082   SiPixelCluster() = default;
0083   ~SiPixelCluster() = default;
0084   SiPixelCluster(SiPixelCluster const&) = default;
0085   SiPixelCluster(SiPixelCluster&&) = default;
0086   SiPixelCluster& operator=(SiPixelCluster const&) = default;
0087   SiPixelCluster& operator=(SiPixelCluster&&) = default;
0088 
0089   SiPixelCluster(unsigned int isize,
0090                  uint16_t const* adcs,
0091                  uint16_t const* xpos,
0092                  uint16_t const* ypos,
0093                  uint16_t xmin,
0094                  uint16_t ymin,
0095                  uint16_t id = invalidClusterId)
0096       : thePixelOffset(2 * isize), thePixelADC(adcs, adcs + isize), theOriginalClusterId(id) {
0097     uint16_t maxCol = 0;
0098     uint16_t maxRow = 0;
0099     for (unsigned int i = 0; i < isize; ++i) {
0100       uint16_t xoffset = xpos[i] - xmin;
0101       uint16_t yoffset = ypos[i] - ymin;
0102       thePixelOffset[i * 2] = std::min(uint16_t(MAXSPAN), xoffset);
0103       thePixelOffset[i * 2 + 1] = std::min(uint16_t(MAXSPAN), yoffset);
0104       if (xoffset > maxRow)
0105         maxRow = xoffset;
0106       if (yoffset > maxCol)
0107         maxCol = yoffset;
0108     }
0109     packRow(xmin, maxRow);
0110     packCol(ymin, maxCol);
0111   }
0112 
0113   // obsolete (only for regression tests)
0114   SiPixelCluster(const PixelPos& pix, int adc);
0115   void add(const PixelPos& pix, int adc);
0116 
0117   // Analog linear average position (barycenter)
0118   float x() const {
0119     float qm = 0.0;
0120     int isize = thePixelADC.size();
0121     for (int i = 0; i < isize; ++i)
0122       qm += float(thePixelADC[i]) * (thePixelOffset[i * 2] + minPixelRow() + 0.5f);
0123     return qm / charge();
0124   }
0125 
0126   float y() const {
0127     float qm = 0.0;
0128     int isize = thePixelADC.size();
0129     for (int i = 0; i < isize; ++i)
0130       qm += float(thePixelADC[i]) * (thePixelOffset[i * 2 + 1] + minPixelCol() + 0.5f);
0131     return qm / charge();
0132   }
0133 
0134   // Return number of pixels.
0135   int size() const { return thePixelADC.size(); }
0136 
0137   // Return cluster dimension in the x direction.
0138   int sizeX() const { return rowSpan() + 1; }
0139 
0140   // Return cluster dimension in the y direction.
0141   int sizeY() const { return colSpan() + 1; }
0142 
0143   inline int charge() const {
0144     int qm = 0;
0145     int isize = thePixelADC.size();
0146     for (int i = 0; i < isize; ++i)
0147       qm += thePixelADC[i];
0148     return qm;
0149   }  // Return total cluster charge.
0150 
0151   inline int minPixelRow() const { return theMinPixelRow; }             // The min x index.
0152   inline int maxPixelRow() const { return minPixelRow() + rowSpan(); }  // The max x index.
0153   inline int minPixelCol() const { return theMinPixelCol; }             // The min y index.
0154   inline int maxPixelCol() const { return minPixelCol() + colSpan(); }  // The max y index.
0155 
0156   const std::vector<uint8_t>& pixelOffset() const { return thePixelOffset; }
0157   const std::vector<uint16_t>& pixelADC() const { return thePixelADC; }
0158 
0159   // obsolete, use single pixel access below
0160   const std::vector<Pixel> pixels() const {
0161     std::vector<Pixel> oldPixVector;
0162     int isize = thePixelADC.size();
0163     oldPixVector.reserve(isize);
0164     for (int i = 0; i < isize; ++i) {
0165       oldPixVector.push_back(pixel(i));
0166     }
0167     return oldPixVector;
0168   }
0169 
0170   // infinite faster than above...
0171   Pixel pixel(int i) const {
0172     return Pixel(minPixelRow() + thePixelOffset[i * 2], minPixelCol() + thePixelOffset[i * 2 + 1], thePixelADC[i]);
0173   }
0174 
0175 private:
0176   static int overflow_(uint16_t span) { return span == uint16_t(MAXSPAN); }
0177 
0178 public:
0179   int colSpan() const { return thePixelColSpan; }
0180 
0181   int rowSpan() const { return thePixelRowSpan; }
0182 
0183   bool overflowCol() const { return overflow_(thePixelColSpan); }
0184 
0185   bool overflowRow() const { return overflow_(thePixelRowSpan); }
0186 
0187   bool overflow() const { return overflowCol() || overflowRow(); }
0188 
0189   void packCol(uint16_t ymin, uint16_t yspan) {
0190     theMinPixelCol = ymin;
0191     thePixelColSpan = std::min(yspan, uint16_t(MAXSPAN));
0192   }
0193   void packRow(uint16_t xmin, uint16_t xspan) {
0194     theMinPixelRow = xmin;
0195     thePixelRowSpan = std::min(xspan, uint16_t(MAXSPAN));
0196   }
0197 
0198   // ggiurgiu@fnal.gov, 01/05/12
0199   // Getters and setters for the newly added data members (err_x and err_y). See below.
0200   void setSplitClusterErrorX(float errx) { err_x = errx; }
0201   void setSplitClusterErrorY(float erry) { err_y = erry; }
0202   float getSplitClusterErrorX() const { return err_x; }
0203   float getSplitClusterErrorY() const { return err_y; }
0204 
0205   // the original id (they get sorted)
0206   auto originalId() const { return theOriginalClusterId; }
0207   void setOriginalId(uint16_t id) { theOriginalClusterId = id; }
0208 
0209 private:
0210   std::vector<uint8_t> thePixelOffset;
0211   std::vector<uint16_t> thePixelADC;
0212 
0213   uint16_t theMinPixelRow = MAXPOS;  // Minimum pixel index in the x direction (low edge).
0214   uint16_t theMinPixelCol = MAXPOS;  // Minimum pixel index in the y direction (left edge).
0215   uint8_t thePixelRowSpan = 0;       // Span pixel index in the x direction (low edge).
0216   uint8_t thePixelColSpan = 0;       // Span pixel index in the y direction (left edge).
0217 
0218   uint16_t theOriginalClusterId = invalidClusterId;
0219 
0220   float err_x = -99999.9f;
0221   float err_y = -99999.9f;
0222 };
0223 
0224 // Comparison operators  (needed by DetSetVector)
0225 inline bool operator<(const SiPixelCluster& one, const SiPixelCluster& other) {
0226   if (one.minPixelRow() < other.minPixelRow()) {
0227     return true;
0228   } else if (one.minPixelRow() > other.minPixelRow()) {
0229     return false;
0230   } else if (one.minPixelCol() < other.minPixelCol()) {
0231     return true;
0232   } else {
0233     return false;
0234   }
0235 }
0236 
0237 #include "DataFormats/Common/interface/DetSetVector.h"
0238 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0239 #include "DataFormats/Common/interface/Ref.h"
0240 #include "DataFormats/Common/interface/DetSetRefVector.h"
0241 
0242 typedef edm::DetSetVector<SiPixelCluster> SiPixelClusterCollection;
0243 typedef edm::Ref<SiPixelClusterCollection, SiPixelCluster> SiPixelClusterRef;
0244 typedef edm::DetSetRefVector<SiPixelCluster> SiPixelClusterRefVector;
0245 typedef edm::RefProd<SiPixelClusterCollection> SiPixelClusterRefProd;
0246 
0247 typedef edmNew::DetSetVector<SiPixelCluster> SiPixelClusterCollectionNew;
0248 typedef edm::Ref<SiPixelClusterCollectionNew, SiPixelCluster> SiPixelClusterRefNew;
0249 #endif