CMS 3D CMS Logo

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