CMS 3D CMS Logo

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