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  thePhase2ReadoutMode(conf.getParameter<int>("Phase2ReadoutMode")),
57  thePhase2DigiBaseline(conf.getParameter<double>("Phase2DigiBaseline")),
58  thePhase2KinkADC(conf.getParameter<int>("Phase2KinkADC")),
59  theNumOfRows(0),
60  theNumOfCols(0),
61  theDetid(0),
62  // Get the constants for the miss-calibration studies
63  doMissCalibrate(conf.getParameter<bool>("MissCalibrate")),
64  doSplitClusters(conf.getParameter<bool>("SplitClusters")) {
66  theFakePixels.clear();
67  thePixelOccurrence.clear();
68 }
71 
72 // Configuration descriptions
74  desc.add<int>("ChannelThreshold", 1000);
75  desc.add<bool>("MissCalibrate", true);
76  desc.add<bool>("SplitClusters", false);
77  desc.add<int>("VCaltoElectronGain", 65);
78  desc.add<int>("VCaltoElectronGain_L1", 65);
79  desc.add<int>("VCaltoElectronOffset", -414);
80  desc.add<int>("VCaltoElectronOffset_L1", -414);
81  desc.add<int>("SeedThreshold", 1000);
82  desc.add<int>("ClusterThreshold_L1", 4000);
83  desc.add<int>("ClusterThreshold", 4000);
84  desc.add<double>("ElectronPerADCGain", 135.);
85  desc.add<bool>("Phase2Calibration", false);
86  desc.add<int>("Phase2ReadoutMode", -1);
87  desc.add<double>("Phase2DigiBaseline", 1200.);
88  desc.add<int>("Phase2KinkADC", 8);
89 }
90 
91 //----------------------------------------------------------------------------
94 //----------------------------------------------------------------------------
96  // Cache the topology.
97  const PixelTopology& topol = pixDet->specificTopology();
98 
99  // Get the new sizes.
100  int nrows = topol.nrows(); // rows in x
101  int ncols = topol.ncolumns(); // cols in y
102 
103  theNumOfRows = nrows; // Set new sizes
105 
106  if (nrows > theBuffer.rows() || ncols > theBuffer.columns()) { // change only when a larger is needed
107  if (nrows != theNumOfRows || ncols != theNumOfCols)
108  edm::LogWarning("setup()") << "pixel buffer redefined to" << nrows << " * " << ncols;
109  //theNumOfRows = nrows; // Set new sizes
110  //theNumOfCols = ncols;
111  // Resize the buffer
112  theBuffer.setSize(nrows, ncols); // Modify
113  }
114 
115  theFakePixels.resize(nrows * ncols, false);
116 
117  thePixelOccurrence.resize(nrows * ncols, 0);
118 
119  return true;
120 }
121 //----------------------------------------------------------------------------
127 //----------------------------------------------------------------------------
128 template <typename T>
130  const PixelGeomDetUnit* pixDet,
131  const TrackerTopology* tTopo,
132  const std::vector<short>& badChannels,
134  typename T::const_iterator begin = input.begin();
135  typename T::const_iterator end = input.end();
136 
137  // this should never happen and the raw2digi does not create empty detsets
138  if (begin == end) {
139  edm::LogError("PixelThresholdClusterizer") << "@SUB=PixelThresholdClusterizer::clusterizeDetUnitT()"
140  << " No digis to clusterize";
141  }
142 
143  // Set up the clusterization on this DetId.
144  if (!setup(pixDet))
145  return;
146 
147  theDetid = input.detId();
148 
149  // Set separate cluster threshold for L1 (needed for phase1)
151  theLayer = (DetId(theDetid).subdetId() == 1) ? tTopo->pxbLayer(theDetid) : 0;
152  if (theLayer == 1)
154 
155  // Copy PixelDigis to the buffer array; select the seed pixels
156  // on the way, and store them in theSeeds.
157  if (end > begin)
158  copy_to_buffer(begin, end);
159 
160  assert(output.empty());
161  // Loop over all seeds. TO DO: wouldn't using iterators be faster?
162  for (unsigned int i = 0; i < theSeeds.size(); i++) {
163  // Gavril : The charge of seeds that were already inlcuded in clusters is set to 1 electron
164  // so we don't want to call "make_cluster" for these cases
165  if (theBuffer(theSeeds[i]) >= theSeedThreshold) { // Is this seed still valid?
166  // Make a cluster around this seed
168 
169  // Check if the cluster is above threshold
170  // (TO DO: one is signed, other unsigned, gcc warns...)
171  if (cluster.charge() >= clusterThreshold) {
172  // sort by row (x)
173  output.push_back(std::move(cluster));
174  std::push_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
175  return cl1.minPixelRow() < cl2.minPixelRow();
176  });
177  }
178  }
179  }
180  // sort by row (x) maybe sorting the seed would suffice....
181  std::sort_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
182  return cl1.minPixelRow() < cl2.minPixelRow();
183  });
184 
185  // Erase the seeds.
186  theSeeds.clear();
187 
188  // Need to clean unused pixels from the buffer array.
189  clear_buffer(begin, end);
190 
191  theFakePixels.clear();
192 
193  thePixelOccurrence.clear();
194 }
195 
196 //----------------------------------------------------------------------------
204 //----------------------------------------------------------------------------
206  for (DigiIterator di = begin; di != end; ++di) {
207  theBuffer.set_adc(di->row(), di->column(), 0); // reset pixel adc to 0
208  }
209 }
210 
212  for (ClusterIterator ci = begin; ci != end; ++ci) {
213  for (int i = 0; i < ci->size(); ++i) {
214  const SiPixelCluster::Pixel pixel = ci->pixel(i);
215 
216  theBuffer.set_adc(pixel.x, pixel.y, 0); // reset pixel adc to 0
217  }
218  }
219 }
220 
221 //----------------------------------------------------------------------------
223 //----------------------------------------------------------------------------
225 #ifdef PIXELREGRESSION
226  static std::atomic<int> s_ic = 0;
227  in ic = ++s_ic;
228  if (ic == 1) {
229  // std::cout << (doMissCalibrate ? "VI from db" : "VI linear") << std::endl;
230  }
231 #endif
232 
233  //If called with empty/invalid DetSet, warn the user
234  if (end <= begin) {
235  edm::LogWarning("PixelThresholdClusterizer") << " copy_to_buffer called with empty or invalid range" << std::endl;
236  return;
237  }
238 
239  int electron[end - begin]; // pixel charge in electrons
240  memset(electron, 0, (end - begin) * sizeof(int));
241 
242  if (doPhase2Calibration) {
243  int i = 0;
244  for (DigiIterator di = begin; di != end; ++di) {
245  electron[i] = calibrate(di->adc(), di->column(), di->row());
246  i++;
247  }
248  assert(i == (end - begin));
249  }
250 
251  else {
252  if (doMissCalibrate) {
253  if (theLayer == 1) {
254  (*theSiPixelGainCalibrationService_)
256  } else {
257  (*theSiPixelGainCalibrationService_).calibrate(theDetid, begin, end, theConversionFactor, theOffset, electron);
258  }
259  } else {
260  int i = 0;
261  const float gain = theElectronPerADCGain; // default: 1 ADC = 135 electrons
262  for (DigiIterator di = begin; di != end; ++di) {
263  auto adc = di->adc();
264  const float pedestal = 0.; //
265  electron[i] = int(adc * gain + pedestal);
266  ++i;
267  }
268  assert(i == (end - begin));
269  }
270  }
271 
272  int i = 0;
273 #ifdef PIXELREGRESSION
274  static std::atomic<int> eqD = 0;
275 #endif
276  for (DigiIterator di = begin; di != end; ++di) {
277  int row = di->row();
278  int col = di->column();
279  // VV: do not calibrate a fake pixel, it already has a unit of 10e-:
280  int adc = (di->flag() != 0) ? di->adc() * 10 : electron[i]; // this is in electrons
281  i++;
282 
283 #ifdef PIXELREGRESSION
284  int adcOld = calibrate(di->adc(), col, row);
285  //assert(adc==adcOld);
286  if (adc != adcOld)
287  std::cout << "VI " << eqD << ' ' << ic << ' ' << end - begin << ' ' << i << ' ' << di->adc() << ' ' << adc << ' '
288  << adcOld << std::endl;
289  else
290  ++eqD;
291 #endif
292 
293  if (adc < 100)
294  adc = 100; // put all negative pixel charges into the 100 elec bin
295  /* This is semi-random good number. The exact number (in place of 100) is irrelevant from the point
296  of view of the final cluster charge since these are typically >= 20000.
297  */
298 
299  thePixelOccurrence[theBuffer.index(row, col)]++; // increment the occurrence counter
300  uint8_t occurrence = thePixelOccurrence[theBuffer.index(row, col)]; // get the occurrence counter
301 
302  switch (occurrence) {
303  // the 1st occurrence (standard treatment)
304  case 1:
305  if (adc >= thePixelThreshold) {
306  theBuffer.set_adc(row, col, adc);
307  // VV: add pixel to the fake list. Only when running on digi collection
308  if (di->flag() != 0)
309  theFakePixels[row * theNumOfCols + col] = true;
310  if (adc >= theSeedThreshold)
311  theSeeds.push_back(SiPixelCluster::PixelPos(row, col));
312  }
313  break;
314 
315  // the 2nd occurrence (duplicate pixel: reset the buffer to 0 and remove from the list of seed pixels)
316  case 2:
317  theBuffer.set_adc(row, col, 0);
319  break;
320 
321  // in case a pixel appears more than twice, nothing needs to be done because it was already removed at the 2nd occurrence
322  }
323  }
324  assert(i == (end - begin));
325 }
326 
328  // loop over clusters
329  for (ClusterIterator ci = begin; ci != end; ++ci) {
330  // loop over pixels
331  for (int i = 0; i < ci->size(); ++i) {
332  const SiPixelCluster::Pixel pixel = ci->pixel(i);
333 
334  int row = pixel.x;
335  int col = pixel.y;
336  int adc = pixel.adc;
337  if (adc >= thePixelThreshold) {
338  theBuffer.add_adc(row, col, adc);
339  if (adc >= theSeedThreshold)
340  theSeeds.push_back(SiPixelCluster::PixelPos(row, col));
341  }
342  }
343  }
344 }
345 
346 //----------------------------------------------------------------------------
347 // Calibrate adc counts to electrons
348 //-----------------------------------------------------------------
350  int electrons = 0;
351 
352  if (doPhase2Calibration) {
353  const float gain = theElectronPerADCGain;
354  int p2rm = (thePhase2ReadoutMode < -1 ? -1 : thePhase2ReadoutMode);
355 
356  if (p2rm == -1) {
357  electrons = int(adc * gain);
358  } else {
359  if (adc < thePhase2KinkADC) {
360  electrons = int((adc + 0.5) * gain);
361  } else {
362  const int dualslopeparam = (thePhase2ReadoutMode < 10 ? thePhase2ReadoutMode : 10);
363  const int dualslope = int(dualslopeparam <= 1 ? 1. : pow(2, dualslopeparam - 1));
365  adc *= dualslope;
367  electrons = int((adc + 0.5 * dualslope) * gain);
368  }
370  }
371 
372  return electrons;
373  }
374 
375  if (doMissCalibrate) {
376  // do not perform calibration if pixel is dead!
377 
380  // Linear approximation of the TANH response
381  // Pixel(0,0,0)
382  //const float gain = 2.95; // 1 ADC = 2.95 VCALs (1/0.339)
383  //const float pedestal = -83.; // -28/0.339
384  // Roc-0 average
385  //const float gain = 1./0.357; // 1 ADC = 2.80 VCALs
386  //const float pedestal = -28.2 * gain; // -79.
387 
388  float DBgain = theSiPixelGainCalibrationService_->getGain(theDetid, col, row);
390  float DBpedestal = pedestal * DBgain;
391 
392  // Roc-6 average
393  //const float gain = 1./0.313; // 1 ADC = 3.19 VCALs
394  //const float pedestal = -6.2 * gain; // -19.8
395  //
396  float vcal = adc * DBgain - DBpedestal;
397 
398  // atanh calibration
399  // Roc-6 average
400  //const float p0 = 0.00492;
401  //const float p1 = 1.998;
402  //const float p2 = 90.6;
403  //const float p3 = 134.1;
404  // Roc-6 average
405  //const float p0 = 0.00382;
406  //const float p1 = 0.886;
407  //const float p2 = 112.7;
408  //const float p3 = 113.0;
409  //float vcal = ( atanh( (adc-p3)/p2) + p1)/p0;
410 
411  if (theLayer == 1) {
413  } else {
415  }
416  }
417  } else { // No misscalibration in the digitizer
418  // Simple (default) linear gain
419  const float gain = theElectronPerADCGain; // default: 1 ADC = 135 electrons
420  const float pedestal = 0.; //
421  electrons = int(adc * gain + pedestal);
422  }
423 
424  return electrons;
425 }
426 
427 //----------------------------------------------------------------------------
429 //----------------------------------------------------------------------------
432  //First we acquire the seeds for the clusters
433  uint16_t seed_adc;
434  std::stack<SiPixelCluster::PixelPos, std::vector<SiPixelCluster::PixelPos> > dead_pixel_stack;
435 
436  //The individual modules have been loaded into a buffer.
437  //After each pixel has been considered by the clusterizer, we set the adc count to 1
438  //to mark that we have already considered it.
439  //The only difference between dead/noisy pixels and standard ones is that for dead/noisy pixels,
440  //We consider the charge of the pixel to always be zero.
441 
442  /* this is not possible as dead and noisy pixel cannot make it into a seed...
443  if ( doMissCalibrate &&
444  (theSiPixelGainCalibrationService_->isDead(theDetid,pix.col(),pix.row()) ||
445  theSiPixelGainCalibrationService_->isNoisy(theDetid,pix.col(),pix.row())) )
446  {
447  std::cout << "IMPOSSIBLE" << std::endl;
448  seed_adc = 0;
449  theBuffer.set_adc(pix, 1);
450  }
451  else {
452  */
453  // Note: each ADC value is limited here to 65535 (std::numeric_limits<uint16_t>::max),
454  // as it is later stored as uint16_t in SiPixelCluster and PixelClusterizerBase/AccretionCluster
455  // (reminder: ADC values here may be expressed in number of electrons)
456  seed_adc = std::min(theBuffer(pix.row(), pix.col()), int(std::numeric_limits<uint16_t>::max()));
457  theBuffer.set_adc(pix, 1);
458  // }
459 
460  AccretionCluster acluster, cldata;
461  acluster.add(pix, seed_adc);
462  cldata.add(pix, seed_adc);
463 
464  //Here we search all pixels adjacent to all pixels in the cluster.
465  bool dead_flag = false;
466  while (!acluster.empty()) {
467  //This is the standard algorithm to find and add a pixel
468  auto curInd = acluster.top();
469  acluster.pop();
470  for (auto c = std::max(0, int(acluster.y[curInd]) - 1);
471  c < std::min(int(acluster.y[curInd]) + 2, theBuffer.columns());
472  ++c) {
473  for (auto r = std::max(0, int(acluster.x[curInd]) - 1);
474  r < std::min(int(acluster.x[curInd]) + 2, theBuffer.rows());
475  ++r) {
476  if (theBuffer(r, c) >= thePixelThreshold) {
477  SiPixelCluster::PixelPos newpix(r, c);
478  auto const newpix_adc = std::min(theBuffer(r, c), int(std::numeric_limits<uint16_t>::max()));
479  if (!acluster.add(newpix, newpix_adc))
480  goto endClus;
481  // VV: no fake pixels in cluster, leads to non-contiguous clusters
482  if (!theFakePixels[r * theNumOfCols + c]) {
483  cldata.add(newpix, newpix_adc);
484  }
485  theBuffer.set_adc(newpix, 1);
486  }
487 
488  /* //Commenting out the addition of dead pixels to the cluster until further testing -- dfehling 06/09
489  //Check on the bounds of the module; this is to keep the isDead and isNoisy modules from returning errors
490  else if(r>= 0 && c >= 0 && (r <= (theNumOfRows-1.)) && (c <= (theNumOfCols-1.))){
491  //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.
492  if((theSiPixelGainCalibrationService_->isDead(theDetid,c,r) || theSiPixelGainCalibrationService_->isNoisy(theDetid,c,r)) && theBuffer(r,c) != 1){
493 
494  //If a pixel is dead or noisy, check to see if we want to split the clusters or not.
495  //Push it into a dead pixel stack in case we want to split the clusters. Otherwise add it to the cluster.
496  //If we are splitting the clusters, we will iterate over the dead pixel stack later.
497 
498  SiPixelCluster::PixelPos newpix(r,c);
499  if(!doSplitClusters){
500 
501  cluster.add(newpix, std::min(theBuffer(r, c), int(std::numeric_limits<uint16_t>::max())));}
502  else if(doSplitClusters){
503  dead_pixel_stack.push(newpix);
504  dead_flag = true;}
505 
506  theBuffer.set_adc(newpix, 1);
507  }
508 
509  }
510  */
511  }
512  }
513 
514  } // while accretion
515 endClus:
516  SiPixelCluster cluster(cldata.isize, cldata.adc, cldata.x, cldata.y, cldata.xmin, cldata.ymin);
517  //Here we split the cluster, if the flag to do so is set and we have found a dead or noisy pixel.
518 
519  if (dead_flag && doSplitClusters) {
520  // Set separate cluster threshold for L1 (needed for phase1)
522  if (theLayer == 1)
524 
525  //Set the first cluster equal to the existing cluster.
526  SiPixelCluster first_cluster = cluster;
527  bool have_second_cluster = false;
528  while (!dead_pixel_stack.empty()) {
529  //consider each found dead pixel
530  SiPixelCluster::PixelPos deadpix = dead_pixel_stack.top();
531  dead_pixel_stack.pop();
532  theBuffer.set_adc(deadpix, 1);
533 
534  //Clusterize the split cluster using the dead pixel as a seed
535  SiPixelCluster second_cluster = make_cluster(deadpix, output);
536 
537  //If both clusters would normally have been found by the clusterizer, put them into output
538  if (second_cluster.charge() >= clusterThreshold && first_cluster.charge() >= clusterThreshold) {
539  output.push_back(second_cluster);
540  have_second_cluster = true;
541  }
542 
543  //We also want to keep the merged cluster in data and let the RecHit algorithm decide which set to keep
544  //This loop adds the second cluster to the first.
545  const std::vector<SiPixelCluster::Pixel>& branch_pixels = second_cluster.pixels();
546  for (unsigned int i = 0; i < branch_pixels.size(); i++) {
547  auto const temp_x = branch_pixels[i].x;
548  auto const temp_y = branch_pixels[i].y;
549  auto const temp_adc = branch_pixels[i].adc;
550  SiPixelCluster::PixelPos newpix(temp_x, temp_y);
551  cluster.add(newpix, temp_adc);
552  }
553  }
554 
555  //Remember to also add the first cluster if we added the second one.
556  if (first_cluster.charge() >= clusterThreshold && have_second_cluster) {
557  output.push_back(first_cluster);
558  std::push_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
559  return cl1.minPixelRow() < cl2.minPixelRow();
560  });
561  }
562  }
563 
564  return cluster;
565 }
void clusterizeDetUnitT(const T &input, const PixelGeomDetUnit *pixDet, const TrackerTopology *tTopo, const std::vector< short > &badChannels, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
Cluster pixels. This method operates on a matrix of pixels and finds the largest contiguous cluster a...
static void fillPSetDescription(edm::ParameterSetDescription &desc)
unsigned int pxbLayer(const DetId &id) const
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
Log< level::Error, false > LogError
assert(be >=bs)
void set_adc(int row, int col, int adc)
static std::string const input
Definition: EdmProvDump.cc:47
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
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
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
Definition: DetId.h:17
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:223
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
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
long double T
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511
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