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 
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  bool fullyMaskedRoc[6][6] = {{false}};
37  CTPPSPixelDetId currentId(detId);
38 
39  // read and store masked pixels after converting ROC numbering into module numbering
41  for (auto const &det : mask) {
42  uint32_t planeId = det.first & ~rocMask;
43 
44  if (planeId == detId) {
45  unsigned int rocNum = (det.first & rocMask) >> rocOffset;
46  if (rocNum > 5)
47  throw cms::Exception("InvalidRocNumber") << "roc number from mask: " << rocNum;
48  if (det.second.fullMask)
49  fullyMaskedRoc[currentId.plane()][rocNum] = true;
50  for (auto const &paio : det.second.maskedPixels) {
51  std::pair<unsigned char, unsigned char> pa = paio;
52  int modCol, modRow;
53  pI.transformToModule(pa.second, pa.first, rocNum, modCol, modRow);
54  std::pair<int, int> modPix(modRow, modCol);
55  maskedPixels.insert(modPix);
56  }
57  }
58  }
59 
60  if (verbosity_)
61  edm::LogInfo("RPixDetClusterizer") << detId << " received digi.size()=" << digi.size();
62 
63  // creating a set of CTPPSPixelDigi's and fill it
64 
65  rpix_digi_set_.clear();
66  rpix_digi_set_.insert(digi.begin(), digi.end());
67  SeedVector_.clear();
68 
69  // calibrate digis here
70  calib_rpix_digi_map_.clear();
71 
72  for (auto const &RPdit : rpix_digi_set_) {
73  uint8_t row = RPdit.row();
74  uint8_t column = RPdit.column();
75  if (row > maxRow || column > maxCol)
76  throw cms::Exception("CTPPSDigiOutofRange") << " row = " << row << " column = " << column;
77 
78  std::pair<unsigned char, unsigned char> pixel = std::make_pair(row, column);
79  unsigned short adc = RPdit.adc();
80  unsigned short electrons = 0;
81 
82  // check if pixel is noisy or dead (i.e. in the mask)
83  const bool is_in = maskedPixels.find(pixel) != maskedPixels.end();
84 
85  // check if pixel inside a fully masked roc
86  int piano = currentId.plane();
87  int rocId = pI.getROCId(column, row);
88 
89  if (!is_in && !fullyMaskedRoc[piano][rocId]) {
90  //calibrate digi and store the new ones
91  electrons = calibrate(detId, adc, row, column, pcalibrations);
94  RPixCalibDigi calibDigi(row, column, adc, electrons);
95  unsigned int index = column * maxRow + row;
96  calib_rpix_digi_map_.insert(std::pair<unsigned int, RPixCalibDigi>(index, calibDigi));
97  SeedVector_.push_back(calibDigi);
98  }
99  }
100  if (verbosity_)
101  edm::LogInfo("RPixDetClusterizer") << " RPix set size = " << calib_rpix_digi_map_.size();
102 
103  for (auto SeedIt : SeedVector_) {
104  make_cluster(SeedIt, clusters);
105  }
106 }
107 
108 void RPixDetClusterizer::make_cluster(RPixCalibDigi const &aSeed, std::vector<CTPPSPixelCluster> &clusters) {
110  unsigned int seedIndex = aSeed.column() * maxRow + aSeed.row();
111  if (calib_rpix_digi_map_.find(seedIndex) == calib_rpix_digi_map_.end()) {
112  return;
113  }
114  // creating a temporary cluster
115  RPixTempCluster atempCluster;
116 
117  // filling the cluster with the seed
118  atempCluster.addPixel(aSeed.row(), aSeed.column(), aSeed.electrons());
119  calib_rpix_digi_map_.erase(seedIndex);
120 
121  while (!atempCluster.empty()) {
122  //This is the standard algorithm to find and add a pixel
123  auto curInd = atempCluster.top();
124  atempCluster.pop();
125  for (auto c = std::max(0, int(atempCluster.col[curInd]) - 1);
126  c < std::min(int(atempCluster.col[curInd]) + 2, maxCol);
127  ++c) {
128  for (auto r = std::max(0, int(atempCluster.row[curInd]) - 1);
129  r < std::min(int(atempCluster.row[curInd]) + 2, maxRow);
130  ++r) {
131  unsigned int currIndex = c * maxRow + r;
132  if (calib_rpix_digi_map_.find(currIndex) != calib_rpix_digi_map_.end()) {
133  if (!atempCluster.addPixel(r, c, calib_rpix_digi_map_[currIndex].electrons())) {
134  clusters.emplace_back(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
135  return;
136  }
137  calib_rpix_digi_map_.erase(currIndex);
138  }
139  }
140  }
141 
142  } // while accretion
143 
144  clusters.emplace_back(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
145 }
146 
148  unsigned int detId, int adc, int row, int col, const CTPPSPixelGainCalibrations *pcalibrations) {
149  float gain = 0;
150  float pedestal = 0;
151  int electrons = 0;
152 
153  if (!doSingleCalibration_) {
154  const CTPPSPixelGainCalibration &DetCalibs = pcalibrations->getGainCalibration(detId);
155 
156  if (DetCalibs.getDetId() != 0) {
157  gain = DetCalibs.getGain(col, row) * highRangeCal / lowRangeCal;
158  pedestal = DetCalibs.getPed(col, row);
159  float vcal = (adc - pedestal) * gain;
161  } else {
163  pedestal = 0;
164  electrons = int(adc * gain - pedestal);
165  }
166 
167  if (verbosity_)
168  edm::LogInfo("RPixCalibration") << "calibrate: adc = " << adc << " electrons = " << electrons
169  << " gain = " << gain << " pedestal = " << pedestal;
170  } else {
172  pedestal = 0;
173  electrons = int(adc * gain - pedestal);
174  if (verbosity_)
175  edm::LogInfo("RPixCalibration") << "calibrate: adc = " << adc << " electrons = " << electrons
176  << " ElectronADCGain = " << gain << " pedestal = " << pedestal;
177  }
178  if (electrons < 0) {
179  edm::LogInfo("RPixCalibration") << "RPixDetClusterizer::calibrate: *** electrons < 0 *** --> " << electrons
180  << " --> electrons = 0";
181  electrons = 0;
182  }
183 
184  return electrons;
185 }
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:26
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_
uint32_t plane() const
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)
int getROCId(const int col, const int row) const
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