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