CMS 3D CMS Logo

PixelThresholdClusterizer.cc
Go to the documentation of this file.
1 //----------------------------------------------------------------------------
20 //----------------------------------------------------------------------------
21 
22 // Our own includes
24 #include "SiPixelArrayBuffer.h"
26 // Geometry
29 //#include "Geometry/CommonTopologies/RectangularPixelTopology.h"
30 
31 // STL
32 #include <stack>
33 #include <vector>
34 #include <iostream>
35 #include <atomic>
36 #include <algorithm>
37 #include <limits>
38 
39 //----------------------------------------------------------------------------
43 //----------------------------------------------------------------------------
45  : // Get thresholds in electrons
46  thePixelThreshold(conf.getParameter<int>("ChannelThreshold")),
47  theSeedThreshold(conf.getParameter<int>("SeedThreshold")),
48  theClusterThreshold(conf.getParameter<int>("ClusterThreshold")),
49  theClusterThreshold_L1(conf.getParameter<int>("ClusterThreshold_L1")),
50  theConversionFactor(conf.getParameter<int>("VCaltoElectronGain")),
51  theConversionFactor_L1(conf.getParameter<int>("VCaltoElectronGain_L1")),
52  theOffset(conf.getParameter<int>("VCaltoElectronOffset")),
53  theOffset_L1(conf.getParameter<int>("VCaltoElectronOffset_L1")),
54  theElectronPerADCGain(conf.getParameter<double>("ElectronPerADCGain")),
55  doPhase2Calibration(conf.getParameter<bool>("Phase2Calibration")),
56  dropDuplicates(conf.getParameter<bool>("DropDuplicates")),
57  thePhase2ReadoutMode(conf.getParameter<int>("Phase2ReadoutMode")),
58  thePhase2DigiBaseline(conf.getParameter<double>("Phase2DigiBaseline")),
59  thePhase2KinkADC(conf.getParameter<int>("Phase2KinkADC")),
60  theNumOfRows(0),
61  theNumOfCols(0),
62  theDetid(0),
63  // Get the constants for the miss-calibration studies
64  doMissCalibrate(conf.getParameter<bool>("MissCalibrate")),
65  doSplitClusters(conf.getParameter<bool>("SplitClusters")) {
67  theFakePixels.clear();
68  thePixelOccurrence.clear();
69 }
72 
73 // Configuration descriptions
75  desc.add<int>("ChannelThreshold", 1000);
76  desc.add<bool>("MissCalibrate", true);
77  desc.add<bool>("SplitClusters", false);
78  desc.add<int>("VCaltoElectronGain", 65);
79  desc.add<int>("VCaltoElectronGain_L1", 65);
80  desc.add<int>("VCaltoElectronOffset", -414);
81  desc.add<int>("VCaltoElectronOffset_L1", -414);
82  desc.add<int>("SeedThreshold", 1000);
83  desc.add<int>("ClusterThreshold_L1", 4000);
84  desc.add<int>("ClusterThreshold", 4000);
85  desc.add<double>("ElectronPerADCGain", 135.);
86  desc.add<bool>("DropDuplicates", true);
87  desc.add<bool>("Phase2Calibration", false);
88  desc.add<int>("Phase2ReadoutMode", -1);
89  desc.add<double>("Phase2DigiBaseline", 1200.);
90  desc.add<int>("Phase2KinkADC", 8);
91 }
92 
93 //----------------------------------------------------------------------------
96 //----------------------------------------------------------------------------
98  // Cache the topology.
99  const PixelTopology& topol = pixDet->specificTopology();
100 
101  // Get the new sizes.
102  int nrows = topol.nrows(); // rows in x
103  int ncols = topol.ncolumns(); // cols in y
104 
105  theNumOfRows = nrows; // Set new sizes
107 
108  if (nrows > theBuffer.rows() || ncols > theBuffer.columns()) { // change only when a larger is needed
109  if (nrows != theNumOfRows || ncols != theNumOfCols)
110  edm::LogWarning("setup()") << "pixel buffer redefined to" << nrows << " * " << ncols;
111  //theNumOfRows = nrows; // Set new sizes
112  //theNumOfCols = ncols;
113  // Resize the buffer
114  theBuffer.setSize(nrows, ncols); // Modify
115  }
116 
117  theFakePixels.resize(nrows * ncols, false);
118 
119  thePixelOccurrence.resize(nrows * ncols, 0);
120 
121  return true;
122 }
123 
124 #include "PixelThresholdClusterizer.icc"
125 
126 //----------------------------------------------------------------------------
134 //----------------------------------------------------------------------------
136  for (DigiIterator di = begin; di != end; ++di) {
137  theBuffer.set_adc(di->row(), di->column(), 0); // reset pixel adc to 0
138  }
139 }
140 
142  for (ClusterIterator ci = begin; ci != end; ++ci) {
143  for (int i = 0; i < ci->size(); ++i) {
144  const SiPixelCluster::Pixel pixel = ci->pixel(i);
145 
146  theBuffer.set_adc(pixel.x, pixel.y, 0); // reset pixel adc to 0
147  }
148  }
149 }
150 
151 //----------------------------------------------------------------------------
153 //----------------------------------------------------------------------------
155 #ifdef PIXELREGRESSION
156  static std::atomic<int> s_ic = 0;
157  in ic = ++s_ic;
158  if (ic == 1) {
159  // std::cout << (doMissCalibrate ? "VI from db" : "VI linear") << std::endl;
160  }
161 #endif
162 
163  //If called with empty/invalid DetSet, warn the user
164  if (end <= begin) {
165  edm::LogWarning("PixelThresholdClusterizer") << " copy_to_buffer called with empty or invalid range" << std::endl;
166  return;
167  }
168 
169  int electron[end - begin]; // pixel charge in electrons
170  memset(electron, 0, (end - begin) * sizeof(int));
171 
172  if (doPhase2Calibration) {
173  int i = 0;
174  for (DigiIterator di = begin; di != end; ++di) {
175  electron[i] = calibrate(di->adc(), di->column(), di->row());
176  i++;
177  }
178  assert(i == (end - begin));
179  }
180 
181  else {
182  if (doMissCalibrate) {
183  if (theLayer == 1) {
184  (*theSiPixelGainCalibrationService_)
186  } else {
187  (*theSiPixelGainCalibrationService_).calibrate(theDetid, begin, end, theConversionFactor, theOffset, electron);
188  }
189  } else {
190  int i = 0;
191  const float gain = theElectronPerADCGain; // default: 1 ADC = 135 electrons
192  for (DigiIterator di = begin; di != end; ++di) {
193  auto adc = di->adc();
194  const float pedestal = 0.; //
195  electron[i] = int(adc * gain + pedestal);
196  ++i;
197  }
198  assert(i == (end - begin));
199  }
200  }
201 
202  int i = 0;
203 #ifdef PIXELREGRESSION
204  static std::atomic<int> eqD = 0;
205 #endif
206  for (DigiIterator di = begin; di != end; ++di) {
207  int row = di->row();
208  int col = di->column();
209  // VV: do not calibrate a fake pixel, it already has a unit of 10e-:
210  int adc = (di->flag() != 0) ? di->adc() * 10 : electron[i]; // this is in electrons
211  i++;
212 
213 #ifdef PIXELREGRESSION
214  int adcOld = calibrate(di->adc(), col, row);
215  //assert(adc==adcOld);
216  if (adc != adcOld)
217  std::cout << "VI " << eqD << ' ' << ic << ' ' << end - begin << ' ' << i << ' ' << di->adc() << ' ' << adc << ' '
218  << adcOld << std::endl;
219  else
220  ++eqD;
221 #endif
222 
223  if (adc < 100)
224  adc = 100; // put all negative pixel charges into the 100 elec bin
225  /* This is semi-random good number. The exact number (in place of 100) is irrelevant from the point
226  of view of the final cluster charge since these are typically >= 20000.
227  */
228 
229  thePixelOccurrence[theBuffer.index(row, col)]++; // increment the occurrence counter
230  uint8_t occurrence =
231  (!dropDuplicates) ? 1 : thePixelOccurrence[theBuffer.index(row, col)]; // get the occurrence counter
232 
233  switch (occurrence) {
234  // the 1st occurrence (standard treatment)
235  case 1:
236  if (adc >= thePixelThreshold) {
237  theBuffer.set_adc(row, col, adc);
238  // VV: add pixel to the fake list. Only when running on digi collection
239  if (di->flag() != 0)
240  theFakePixels[row * theNumOfCols + col] = true;
241  if (adc >= theSeedThreshold)
242  theSeeds.push_back(SiPixelCluster::PixelPos(row, col));
243  }
244  break;
245 
246  // the 2nd occurrence (duplicate pixel: reset the buffer to 0 and remove from the list of seed pixels)
247  case 2:
248  theBuffer.set_adc(row, col, 0);
250  break;
251 
252  // in case a pixel appears more than twice, nothing needs to be done because it was already removed at the 2nd occurrence
253  }
254  }
255  assert(i == (end - begin));
256 }
257 
259  // loop over clusters
260  for (ClusterIterator ci = begin; ci != end; ++ci) {
261  // loop over pixels
262  for (int i = 0; i < ci->size(); ++i) {
263  const SiPixelCluster::Pixel pixel = ci->pixel(i);
264 
265  int row = pixel.x;
266  int col = pixel.y;
267  int adc = pixel.adc;
268  if (adc >= thePixelThreshold) {
269  theBuffer.add_adc(row, col, adc);
270  if (adc >= theSeedThreshold)
271  theSeeds.push_back(SiPixelCluster::PixelPos(row, col));
272  }
273  }
274  }
275 }
276 
277 //----------------------------------------------------------------------------
278 // Calibrate adc counts to electrons
279 //-----------------------------------------------------------------
281  int electrons = 0;
282 
283  if (doPhase2Calibration) {
284  const float gain = theElectronPerADCGain;
285  int p2rm = (thePhase2ReadoutMode < -1 ? -1 : thePhase2ReadoutMode);
286 
287  if (p2rm == -1) {
288  electrons = int(adc * gain);
289  } else {
290  if (adc < thePhase2KinkADC) {
291  electrons = int((adc + 0.5) * gain);
292  } else {
293  const int dualslopeparam = (thePhase2ReadoutMode < 10 ? thePhase2ReadoutMode : 10);
294  const int dualslope = int(dualslopeparam <= 1 ? 1. : pow(2, dualslopeparam - 1));
296  adc *= dualslope;
298  electrons = int((adc + 0.5 * dualslope) * gain);
299  }
301  }
302 
303  return electrons;
304  }
305 
306  if (doMissCalibrate) {
307  // do not perform calibration if pixel is dead!
308 
311  // Linear approximation of the TANH response
312  // Pixel(0,0,0)
313  //const float gain = 2.95; // 1 ADC = 2.95 VCALs (1/0.339)
314  //const float pedestal = -83.; // -28/0.339
315  // Roc-0 average
316  //const float gain = 1./0.357; // 1 ADC = 2.80 VCALs
317  //const float pedestal = -28.2 * gain; // -79.
318 
319  float DBgain = theSiPixelGainCalibrationService_->getGain(theDetid, col, row);
321  float DBpedestal = pedestal * DBgain;
322 
323  // Roc-6 average
324  //const float gain = 1./0.313; // 1 ADC = 3.19 VCALs
325  //const float pedestal = -6.2 * gain; // -19.8
326  //
327  float vcal = adc * DBgain - DBpedestal;
328 
329  // atanh calibration
330  // Roc-6 average
331  //const float p0 = 0.00492;
332  //const float p1 = 1.998;
333  //const float p2 = 90.6;
334  //const float p3 = 134.1;
335  // Roc-6 average
336  //const float p0 = 0.00382;
337  //const float p1 = 0.886;
338  //const float p2 = 112.7;
339  //const float p3 = 113.0;
340  //float vcal = ( atanh( (adc-p3)/p2) + p1)/p0;
341 
342  if (theLayer == 1) {
344  } else {
346  }
347  }
348  } else { // No misscalibration in the digitizer
349  // Simple (default) linear gain
350  const float gain = theElectronPerADCGain; // default: 1 ADC = 135 electrons
351  const float pedestal = 0.; //
352  electrons = int(adc * gain + pedestal);
353  }
354 
355  return electrons;
356 }
357 
358 //----------------------------------------------------------------------------
360 //----------------------------------------------------------------------------
363  //First we acquire the seeds for the clusters
364  uint16_t seed_adc;
365  std::stack<SiPixelCluster::PixelPos, std::vector<SiPixelCluster::PixelPos> > dead_pixel_stack;
366 
367  //The individual modules have been loaded into a buffer.
368  //After each pixel has been considered by the clusterizer, we set the adc count to 1
369  //to mark that we have already considered it.
370  //The only difference between dead/noisy pixels and standard ones is that for dead/noisy pixels,
371  //We consider the charge of the pixel to always be zero.
372 
373  /* this is not possible as dead and noisy pixel cannot make it into a seed...
374  if ( doMissCalibrate &&
375  (theSiPixelGainCalibrationService_->isDead(theDetid,pix.col(),pix.row()) ||
376  theSiPixelGainCalibrationService_->isNoisy(theDetid,pix.col(),pix.row())) )
377  {
378  std::cout << "IMPOSSIBLE" << std::endl;
379  seed_adc = 0;
380  theBuffer.set_adc(pix, 1);
381  }
382  else {
383  */
384  // Note: each ADC value is limited here to 65535 (std::numeric_limits<uint16_t>::max),
385  // as it is later stored as uint16_t in SiPixelCluster and PixelClusterizerBase/AccretionCluster
386  // (reminder: ADC values here may be expressed in number of electrons)
387  seed_adc = std::min(theBuffer(pix.row(), pix.col()), int(std::numeric_limits<uint16_t>::max()));
388  theBuffer.set_adc(pix, 1);
389  // }
390 
391  AccretionCluster acluster, cldata;
392  acluster.add(pix, seed_adc);
393  cldata.add(pix, seed_adc);
394 
395  //Here we search all pixels adjacent to all pixels in the cluster.
396  bool dead_flag = false;
397  while (!acluster.empty()) {
398  //This is the standard algorithm to find and add a pixel
399  auto curInd = acluster.top();
400  acluster.pop();
401  for (auto c = std::max(0, int(acluster.y[curInd]) - 1);
402  c < std::min(int(acluster.y[curInd]) + 2, theBuffer.columns());
403  ++c) {
404  for (auto r = std::max(0, int(acluster.x[curInd]) - 1);
405  r < std::min(int(acluster.x[curInd]) + 2, theBuffer.rows());
406  ++r) {
407  if (theBuffer(r, c) >= thePixelThreshold) {
408  SiPixelCluster::PixelPos newpix(r, c);
409  auto const newpix_adc = std::min(theBuffer(r, c), int(std::numeric_limits<uint16_t>::max()));
410  if (!acluster.add(newpix, newpix_adc))
411  goto endClus;
412  // VV: no fake pixels in cluster, leads to non-contiguous clusters
413  if (!theFakePixels[r * theNumOfCols + c]) {
414  cldata.add(newpix, newpix_adc);
415  }
416  theBuffer.set_adc(newpix, 1);
417  }
418 
419  /* //Commenting out the addition of dead pixels to the cluster until further testing -- dfehling 06/09
420  //Check on the bounds of the module; this is to keep the isDead and isNoisy modules from returning errors
421  else if(r>= 0 && c >= 0 && (r <= (theNumOfRows-1.)) && (c <= (theNumOfCols-1.))){
422  //Check for dead/noisy pixels check that the buffer is not -1 (already considered). Check whether we want to split clusters separated by dead pixels or not.
423  if((theSiPixelGainCalibrationService_->isDead(theDetid,c,r) || theSiPixelGainCalibrationService_->isNoisy(theDetid,c,r)) && theBuffer(r,c) != 1){
424 
425  //If a pixel is dead or noisy, check to see if we want to split the clusters or not.
426  //Push it into a dead pixel stack in case we want to split the clusters. Otherwise add it to the cluster.
427  //If we are splitting the clusters, we will iterate over the dead pixel stack later.
428 
429  SiPixelCluster::PixelPos newpix(r,c);
430  if(!doSplitClusters){
431 
432  cluster.add(newpix, std::min(theBuffer(r, c), int(std::numeric_limits<uint16_t>::max())));}
433  else if(doSplitClusters){
434  dead_pixel_stack.push(newpix);
435  dead_flag = true;}
436 
437  theBuffer.set_adc(newpix, 1);
438  }
439 
440  }
441  */
442  }
443  }
444 
445  } // while accretion
446 endClus:
447  SiPixelCluster cluster(cldata.isize, cldata.adc, cldata.x, cldata.y, cldata.xmin, cldata.ymin);
448  //Here we split the cluster, if the flag to do so is set and we have found a dead or noisy pixel.
449 
450  if (dead_flag && doSplitClusters) {
451  // Set separate cluster threshold for L1 (needed for phase1)
453  if (theLayer == 1)
455 
456  //Set the first cluster equal to the existing cluster.
457  SiPixelCluster first_cluster = cluster;
458  bool have_second_cluster = false;
459  while (!dead_pixel_stack.empty()) {
460  //consider each found dead pixel
461  SiPixelCluster::PixelPos deadpix = dead_pixel_stack.top();
462  dead_pixel_stack.pop();
463  theBuffer.set_adc(deadpix, 1);
464 
465  //Clusterize the split cluster using the dead pixel as a seed
466  SiPixelCluster second_cluster = make_cluster(deadpix, output);
467 
468  //If both clusters would normally have been found by the clusterizer, put them into output
469  if (second_cluster.charge() >= clusterThreshold && first_cluster.charge() >= clusterThreshold) {
470  output.push_back(second_cluster);
471  have_second_cluster = true;
472  }
473 
474  //We also want to keep the merged cluster in data and let the RecHit algorithm decide which set to keep
475  //This loop adds the second cluster to the first.
476  const std::vector<SiPixelCluster::Pixel>& branch_pixels = second_cluster.pixels();
477  for (unsigned int i = 0; i < branch_pixels.size(); i++) {
478  auto const temp_x = branch_pixels[i].x;
479  auto const temp_y = branch_pixels[i].y;
480  auto const temp_adc = branch_pixels[i].adc;
481  SiPixelCluster::PixelPos newpix(temp_x, temp_y);
482  cluster.add(newpix, temp_adc);
483  }
484  }
485 
486  //Remember to also add the first cluster if we added the second one.
487  if (first_cluster.charge() >= clusterThreshold && have_second_cluster) {
488  output.push_back(first_cluster);
489  std::push_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
490  return cl1.minPixelRow() < cl2.minPixelRow();
491  });
492  }
493  }
494 
495  return cluster;
496 }
static void fillPSetDescription(edm::ParameterSetDescription &desc)
SiPixelArrayBuffer theBuffer
Data storage.
PixelThresholdClusterizer(edm::ParameterSet const &conf)
SiPixelCluster make_cluster(const SiPixelCluster::PixelPos &pix, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
The actual clustering algorithm: group the neighboring pixels around the seed.
virtual int ncolumns() const =0
virtual int nrows() const =0
virtual bool isDead(const uint32_t &detID, const int &col, const int &row)=0
assert(be >=bs)
void set_adc(int row, int col, int adc)
int minPixelRow() const
int charge() const
const std::vector< Pixel > pixels() const
void add_adc(int row, int col, int adc)
int index(int row, int col) const
Definition of indexing within the buffer.
edm::DetSet< PixelDigi >::const_iterator DigiIterator
void add(const PixelPos &pix, int adc)
virtual float getPedestal(const uint32_t &detID, const int &col, const int &row)=0
int theNumOfRows
Geometry-related information.
std::vector< uint8_t > thePixelOccurrence
edmNew::DetSet< SiPixelCluster >::const_iterator ClusterIterator
void clear_buffer(DigiIterator begin, DigiIterator end)
Clear the internal buffer array.
void setSize(int rows, int cols)
SiPixelGainCalibrationServiceBase * theSiPixelGainCalibrationService_
def remove(d, key, TELL=False)
Definition: MatrixUtil.py:233
Pixel cluster – collection of neighboring pixels above threshold.
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
col
Definition: cuy.py:1009
constexpr int row() const
virtual float getGain(const uint32_t &detID, const int &col, const int &row)=0
bool setup(const PixelGeomDetUnit *pixDet)
Private helper methods:
constexpr int col() const
Definition: output.py:1
bool add(SiPixelCluster::PixelPos const &p, uint16_t const iadc)
virtual bool isNoisy(const uint32_t &detID, const int &col, const int &row)=0
int calibrate(int adc, int col, int row)
Log< level::Warning, false > LogWarning
std::vector< SiPixelCluster::PixelPos > theSeeds
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
void copy_to_buffer(DigiIterator begin, DigiIterator end)
Copy adc counts from PixelDigis into the buffer, identify seeds.
uint16_t *__restrict__ uint16_t const *__restrict__ adc