CMS 3D CMS Logo

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