CMS 3D CMS Logo

RPixDetClusterizer.cc
Go to the documentation of this file.
1 #include <iostream>
3 
6 
7 namespace {
8  constexpr int maxCol = CTPPSPixelCluster::MAXCOL;
9  constexpr int maxRow = CTPPSPixelCluster::MAXROW;
10  constexpr double highRangeCal = 1800.;
11  constexpr double lowRangeCal = 260.;
12  constexpr int rocMask = 0xE000;
13  constexpr int rocOffset = 13;
14 } // namespace
15 
16 RPixDetClusterizer::RPixDetClusterizer(edm::ParameterSet const &conf) : params_(conf), SeedVector_(0) {
17  verbosity_ = conf.getUntrackedParameter<int>("RPixVerbosity");
18  SeedADCThreshold_ = conf.getParameter<int>("SeedADCThreshold");
19  ADCThreshold_ = conf.getParameter<int>("ADCThreshold");
20  ElectronADCGain_ = conf.getParameter<double>("ElectronADCGain");
21  VcaltoElectronGain_ = conf.getParameter<int>("VCaltoElectronGain");
22  VcaltoElectronOffset_ = conf.getParameter<int>("VCaltoElectronOffset");
23  doSingleCalibration_ = conf.getParameter<bool>("doSingleCalibration");
24 }
25 
27 
28 void RPixDetClusterizer::buildClusters(unsigned int detId,
29  const std::vector<CTPPSPixelDigi> &digi,
30  std::vector<CTPPSPixelCluster> &clusters,
31  const CTPPSPixelGainCalibrations *pcalibrations,
32  const CTPPSPixelAnalysisMask *maskera) {
33  std::map<uint32_t, CTPPSPixelROCAnalysisMask> const &mask = maskera->analysisMask;
34  std::set<std::pair<unsigned char, unsigned char> > maskedPixels;
35 
36  // read and store masked pixels after converting ROC numbering into module numbering
38  for (auto const &det : mask) {
39  uint32_t planeId = det.first & ~rocMask;
40 
41  if (planeId == detId) {
42  unsigned int rocNum = (det.first & rocMask) >> rocOffset;
43  if (rocNum > 5)
44  throw cms::Exception("InvalidRocNumber") << "roc number from mask: " << rocNum;
45  for (auto const &paio : det.second.maskedPixels) {
46  std::pair<unsigned char, unsigned char> pa = paio;
47  int modCol, modRow;
48  pI.transformToModule(pa.second, pa.first, rocNum, modCol, modRow);
49  std::pair<int, int> modPix(modRow, modCol);
50  maskedPixels.insert(modPix);
51  }
52  }
53  }
54 
55  if (verbosity_)
56  edm::LogInfo("RPixDetClusterizer") << detId << " received digi.size()=" << digi.size();
57 
58  // creating a set of CTPPSPixelDigi's and fill it
59 
60  rpix_digi_set_.clear();
61  rpix_digi_set_.insert(digi.begin(), digi.end());
62  SeedVector_.clear();
63 
64  // calibrate digis here
65  calib_rpix_digi_map_.clear();
66 
67  for (auto const &RPdit : rpix_digi_set_) {
68  uint8_t row = RPdit.row();
69  uint8_t column = RPdit.column();
70  if (row > maxRow || column > maxCol)
71  throw cms::Exception("CTPPSDigiOutofRange") << " row = " << row << " column = " << column;
72 
73  std::pair<unsigned char, unsigned char> pixel = std::make_pair(row, column);
74  unsigned short adc = RPdit.adc();
75  unsigned short electrons = 0;
76 
77  // check if pixel is noisy or dead (i.e. in the mask)
78  const bool is_in = maskedPixels.find(pixel) != maskedPixels.end();
79  if (!is_in) {
80  //calibrate digi and store the new ones
81  electrons = calibrate(detId, adc, row, column, pcalibrations);
84  RPixCalibDigi calibDigi(row, column, adc, electrons);
85  unsigned int index = column * maxRow + row;
86  calib_rpix_digi_map_.insert(std::pair<unsigned int, RPixCalibDigi>(index, calibDigi));
87  SeedVector_.push_back(calibDigi);
88  }
89  }
90  if (verbosity_)
91  edm::LogInfo("RPixDetClusterizer") << " RPix set size = " << calib_rpix_digi_map_.size();
92 
93  for (auto SeedIt : SeedVector_) {
94  make_cluster(SeedIt, clusters);
95  }
96 }
97 
98 void RPixDetClusterizer::make_cluster(RPixCalibDigi const &aSeed, std::vector<CTPPSPixelCluster> &clusters) {
100  unsigned int seedIndex = aSeed.column() * maxRow + aSeed.row();
101  if (calib_rpix_digi_map_.find(seedIndex) == calib_rpix_digi_map_.end()) {
102  return;
103  }
104  // creating a temporary cluster
105  RPixTempCluster atempCluster;
106 
107  // filling the cluster with the seed
108  atempCluster.addPixel(aSeed.row(), aSeed.column(), aSeed.electrons());
109  calib_rpix_digi_map_.erase(seedIndex);
110 
111  while (!atempCluster.empty()) {
112  //This is the standard algorithm to find and add a pixel
113  auto curInd = atempCluster.top();
114  atempCluster.pop();
115  for (auto c = std::max(0, int(atempCluster.col[curInd]) - 1);
116  c < std::min(int(atempCluster.col[curInd]) + 2, maxCol);
117  ++c) {
118  for (auto r = std::max(0, int(atempCluster.row[curInd]) - 1);
119  r < std::min(int(atempCluster.row[curInd]) + 2, maxRow);
120  ++r) {
121  unsigned int currIndex = c * maxRow + r;
122  if (calib_rpix_digi_map_.find(currIndex) != calib_rpix_digi_map_.end()) {
123  if (!atempCluster.addPixel(r, c, calib_rpix_digi_map_[currIndex].electrons())) {
124  clusters.emplace_back(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
125  return;
126  }
127  calib_rpix_digi_map_.erase(currIndex);
128  }
129  }
130  }
131 
132  } // while accretion
133 
134  clusters.emplace_back(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
135 }
136 
138  unsigned int detId, int adc, int row, int col, const CTPPSPixelGainCalibrations *pcalibrations) {
139  float gain = 0;
140  float pedestal = 0;
141  int electrons = 0;
142 
143  if (!doSingleCalibration_) {
144  const CTPPSPixelGainCalibration &DetCalibs = pcalibrations->getGainCalibration(detId);
145 
146  if (DetCalibs.getDetId() != 0) {
147  gain = DetCalibs.getGain(col, row) * highRangeCal / lowRangeCal;
148  pedestal = DetCalibs.getPed(col, row);
149  float vcal = (adc - pedestal) * gain;
151  } else {
153  pedestal = 0;
154  electrons = int(adc * gain - pedestal);
155  }
156 
157  if (verbosity_)
158  edm::LogInfo("RPixCalibration") << "calibrate: adc = " << adc << " electrons = " << electrons
159  << " gain = " << gain << " pedestal = " << pedestal;
160  } else {
162  pedestal = 0;
163  electrons = int(adc * gain - pedestal);
164  if (verbosity_)
165  edm::LogInfo("RPixCalibration") << "calibrate: adc = " << adc << " electrons = " << electrons
166  << " ElectronADCGain = " << gain << " pedestal = " << pedestal;
167  }
168  if (electrons < 0) {
169  edm::LogInfo("RPixCalibration") << "RPixDetClusterizer::calibrate: *** electrons < 0 *** --> " << electrons
170  << " --> electrons = 0";
171  electrons = 0;
172  }
173 
174  return electrons;
175 }
float getPed(const int &col, const int &row) const
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
unsigned short adc[MAXSIZE]
std::map< unsigned int, RPixCalibDigi > calib_rpix_digi_map_
std::map< uint32_t, CTPPSPixelROCAnalysisMask > analysisMask
unsigned short isize
RPixDetClusterizer(edm::ParameterSet const &conf)
constexpr uint32_t mask
Definition: gpuClustering.h:24
T getUntrackedParameter(std::string const &, T const &) const
float getGain(const int &col, const int &row) const
int electrons() const
std::vector< RPixCalibDigi > SeedVector_
uint8_t row[MAXSIZE]
void make_cluster(RPixCalibDigi const &aSeed, std::vector< CTPPSPixelCluster > &clusters)
std::set< CTPPSPixelDigi > rpix_digi_set_
Channel-mask mapping.
unsigned short SeedADCThreshold_
Log< level::Info, false > LogInfo
static constexpr uint8_t MAXROW
static constexpr uint8_t MAXCOL
unsigned short ADCThreshold_
const CTPPSPixelGainCalibration & getGainCalibration(const uint32_t &detid) const
bool addPixel(unsigned char myrow, unsigned char mycol, unsigned short const iadc)
uint8_t col[MAXSIZE]
int calibrate(unsigned int, int, int, int, const CTPPSPixelGainCalibrations *pcalibration)
int transformToModule(const int colROC, const int rowROC, const int rocId, int &col, int &row) const
void buildClusters(unsigned int detId, const std::vector< CTPPSPixelDigi > &digi, std::vector< CTPPSPixelCluster > &clusters, const CTPPSPixelGainCalibrations *pcalibration, const CTPPSPixelAnalysisMask *mask)
col
Definition: cuy.py:1009
int row() const
Access to digi information.
unsigned short top() const
int column() const
uint16_t *__restrict__ uint16_t const *__restrict__ adc