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 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
103 
104  if (nrows > theBuffer.rows() || ncols > theBuffer.columns()) { // change only when a larger is needed
105  //if( nrows != theNumOfRows || ncols != theNumOfCols ) {
106  //cout << " PixelThresholdClusterizer: pixel buffer redefined to "
107  // << nrows << " * " << ncols << endl;
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  // Do not bother for empty detectors
135  //if (begin == end) cout << " PixelThresholdClusterizer::clusterizeDetUnit - No digis to clusterize";
136 
137  // Set up the clusterization on this DetId.
138  if (!setup(pixDet))
139  return;
140 
141  theDetid = input.detId();
142 
143  // Set separate cluster threshold for L1 (needed for phase1)
145  theLayer = (DetId(theDetid).subdetId() == 1) ? tTopo->pxbLayer(theDetid) : 0;
146  if (theLayer == 1)
148 
149  // Copy PixelDigis to the buffer array; select the seed pixels
150  // on the way, and store them in theSeeds.
151  copy_to_buffer(begin, end);
152 
153  assert(output.empty());
154  // Loop over all seeds. TO DO: wouldn't using iterators be faster?
155  // edm::LogError("PixelThresholdClusterizer") << "Starting clusterizing" << endl;
156  for (unsigned int i = 0; i < theSeeds.size(); i++) {
157  // Gavril : The charge of seeds that were already inlcuded in clusters is set to 1 electron
158  // so we don't want to call "make_cluster" for these cases
159  if (theBuffer(theSeeds[i]) >= theSeedThreshold) { // Is this seed still valid?
160  // Make a cluster around this seed
162 
163  // Check if the cluster is above threshold
164  // (TO DO: one is signed, other unsigned, gcc warns...)
165  if (cluster.charge() >= clusterThreshold) {
166  // std::cout << "putting in this cluster " << i << " " << cluster.charge() << " " << cluster.pixelADC().size() << endl;
167  // sort by row (x)
168  output.push_back(std::move(cluster));
169  std::push_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
170  return cl1.minPixelRow() < cl2.minPixelRow();
171  });
172  }
173  }
174  }
175  // sort by row (x) maybe sorting the seed would suffice....
176  std::sort_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
177  return cl1.minPixelRow() < cl2.minPixelRow();
178  });
179 
180  // Erase the seeds.
181  theSeeds.clear();
182 
183  // Need to clean unused pixels from the buffer array.
184  clear_buffer(begin, end);
185 
186  theFakePixels.clear();
187 }
188 
189 //----------------------------------------------------------------------------
197 //----------------------------------------------------------------------------
199  for (DigiIterator di = begin; di != end; ++di) {
200  theBuffer.set_adc(di->row(), di->column(), 0); // reset pixel adc to 0
201  }
202 }
203 
205  for (ClusterIterator ci = begin; ci != end; ++ci) {
206  for (int i = 0; i < ci->size(); ++i) {
207  const SiPixelCluster::Pixel pixel = ci->pixel(i);
208 
209  theBuffer.set_adc(pixel.x, pixel.y, 0); // reset pixel adc to 0
210  }
211  }
212 }
213 
214 //----------------------------------------------------------------------------
216 //----------------------------------------------------------------------------
218 #ifdef PIXELREGRESSION
219  static std::atomic<int> s_ic = 0;
220  in ic = ++s_ic;
221  if (ic == 1) {
222  // std::cout << (doMissCalibrate ? "VI from db" : "VI linear") << std::endl;
223  }
224 #endif
225  int electron[end - begin]; // pixel charge in electrons
226  memset(electron, 0, (end - begin) * sizeof(int));
227 
228  if (doPhase2Calibration) {
229  int i = 0;
230  for (DigiIterator di = begin; di != end; ++di) {
231  electron[i] = calibrate(di->adc(), di->column(), di->row());
232  i++;
233  }
234  assert(i == (end - begin));
235  }
236 
237  else {
238  if (doMissCalibrate) {
239  if (theLayer == 1) {
240  (*theSiPixelGainCalibrationService_)
242  } else {
243  (*theSiPixelGainCalibrationService_).calibrate(theDetid, begin, end, theConversionFactor, theOffset, electron);
244  }
245  } else {
246  int i = 0;
247  const float gain = theElectronPerADCGain; // default: 1 ADC = 135 electrons
248  for (DigiIterator di = begin; di != end; ++di) {
249  auto adc = di->adc();
250  const float pedestal = 0.; //
251  electron[i] = int(adc * gain + pedestal);
252  ++i;
253  }
254  assert(i == (end - begin));
255  }
256  }
257 
258  int i = 0;
259 #ifdef PIXELREGRESSION
260  static std::atomic<int> eqD = 0;
261 #endif
262  for (DigiIterator di = begin; di != end; ++di) {
263  int row = di->row();
264  int col = di->column();
265  // VV: do not calibrate a fake pixel, it already has a unit of 10e-:
266  int adc = (di->flag() != 0) ? di->adc() * 10 : electron[i]; // this is in electrons
267  i++;
268 
269 #ifdef PIXELREGRESSION
270  int adcOld = calibrate(di->adc(), col, row);
271  //assert(adc==adcOld);
272  if (adc != adcOld)
273  std::cout << "VI " << eqD << ' ' << ic << ' ' << end - begin << ' ' << i << ' ' << di->adc() << ' ' << adc << ' '
274  << adcOld << std::endl;
275  else
276  ++eqD;
277 #endif
278 
279  if (adc < 100)
280  adc = 100; // put all negative pixel charges into the 100 elec bin
281  /* This is semi-random good number. The exact number (in place of 100) is irrelevant from the point
282  of view of the final cluster charge since these are typically >= 20000.
283  */
284 
285  if (adc >= thePixelThreshold) {
286  theBuffer.set_adc(row, col, adc);
287  // VV: add pixel to the fake list. Only when running on digi collection
288  if (di->flag() != 0)
289  theFakePixels[row * theNumOfCols + col] = true;
290  if (adc >= theSeedThreshold)
291  theSeeds.push_back(SiPixelCluster::PixelPos(row, col));
292  }
293  }
294  assert(i == (end - begin));
295 }
296 
298  // loop over clusters
299  for (ClusterIterator ci = begin; ci != end; ++ci) {
300  // loop over pixels
301  for (int i = 0; i < ci->size(); ++i) {
302  const SiPixelCluster::Pixel pixel = ci->pixel(i);
303 
304  int row = pixel.x;
305  int col = pixel.y;
306  int adc = pixel.adc;
307  if (adc >= thePixelThreshold) {
308  theBuffer.add_adc(row, col, adc);
309  if (adc >= theSeedThreshold)
310  theSeeds.push_back(SiPixelCluster::PixelPos(row, col));
311  }
312  }
313  }
314 }
315 
316 //----------------------------------------------------------------------------
317 // Calibrate adc counts to electrons
318 //-----------------------------------------------------------------
320  int electrons = 0;
321 
322  if (doPhase2Calibration) {
323  const float gain = theElectronPerADCGain;
324  int p2rm = (thePhase2ReadoutMode < -1 ? -1 : thePhase2ReadoutMode);
325 
326  if (p2rm == -1) {
327  electrons = int(adc * gain);
328  } else {
329  if (adc < thePhase2KinkADC) {
330  electrons = int((adc + 0.5) * gain);
331  } else {
332  const int dualslopeparam = (thePhase2ReadoutMode < 10 ? thePhase2ReadoutMode : 10);
333  const int dualslope = int(dualslopeparam <= 1 ? 1. : pow(2, dualslopeparam - 1));
335  adc *= dualslope;
337  electrons = int((adc + 0.5 * dualslope) * gain);
338  }
340  }
341 
342  return electrons;
343  }
344 
345  if (doMissCalibrate) {
346  // do not perform calibration if pixel is dead!
347 
350  // Linear approximation of the TANH response
351  // Pixel(0,0,0)
352  //const float gain = 2.95; // 1 ADC = 2.95 VCALs (1/0.339)
353  //const float pedestal = -83.; // -28/0.339
354  // Roc-0 average
355  //const float gain = 1./0.357; // 1 ADC = 2.80 VCALs
356  //const float pedestal = -28.2 * gain; // -79.
357 
358  float DBgain = theSiPixelGainCalibrationService_->getGain(theDetid, col, row);
360  float DBpedestal = pedestal * DBgain;
361 
362  // Roc-6 average
363  //const float gain = 1./0.313; // 1 ADC = 3.19 VCALs
364  //const float pedestal = -6.2 * gain; // -19.8
365  //
366  float vcal = adc * DBgain - DBpedestal;
367 
368  // atanh calibration
369  // Roc-6 average
370  //const float p0 = 0.00492;
371  //const float p1 = 1.998;
372  //const float p2 = 90.6;
373  //const float p3 = 134.1;
374  // Roc-6 average
375  //const float p0 = 0.00382;
376  //const float p1 = 0.886;
377  //const float p2 = 112.7;
378  //const float p3 = 113.0;
379  //float vcal = ( atanh( (adc-p3)/p2) + p1)/p0;
380 
381  if (theLayer == 1) {
383  } else {
385  }
386  }
387  } else { // No misscalibration in the digitizer
388  // Simple (default) linear gain
389  const float gain = theElectronPerADCGain; // default: 1 ADC = 135 electrons
390  const float pedestal = 0.; //
391  electrons = int(adc * gain + pedestal);
392  }
393 
394  return electrons;
395 }
396 
397 //----------------------------------------------------------------------------
399 //----------------------------------------------------------------------------
402  //First we acquire the seeds for the clusters
403  int seed_adc;
404  stack<SiPixelCluster::PixelPos, vector<SiPixelCluster::PixelPos> > dead_pixel_stack;
405 
406  //The individual modules have been loaded into a buffer.
407  //After each pixel has been considered by the clusterizer, we set the adc count to 1
408  //to mark that we have already considered it.
409  //The only difference between dead/noisy pixels and standard ones is that for dead/noisy pixels,
410  //We consider the charge of the pixel to always be zero.
411 
412  /* this is not possible as dead and noisy pixel cannot make it into a seed...
413  if ( doMissCalibrate &&
414  (theSiPixelGainCalibrationService_->isDead(theDetid,pix.col(),pix.row()) ||
415  theSiPixelGainCalibrationService_->isNoisy(theDetid,pix.col(),pix.row())) )
416  {
417  std::cout << "IMPOSSIBLE" << std::endl;
418  seed_adc = 0;
419  theBuffer.set_adc(pix, 1);
420  }
421  else {
422  */
423  seed_adc = theBuffer(pix.row(), pix.col());
424  theBuffer.set_adc(pix, 1);
425  // }
426 
427  AccretionCluster acluster, cldata;
428  acluster.add(pix, seed_adc);
429  cldata.add(pix, seed_adc);
430 
431  //Here we search all pixels adjacent to all pixels in the cluster.
432  bool dead_flag = false;
433  while (!acluster.empty()) {
434  //This is the standard algorithm to find and add a pixel
435  auto curInd = acluster.top();
436  acluster.pop();
437  for (auto c = std::max(0, int(acluster.y[curInd]) - 1);
438  c < std::min(int(acluster.y[curInd]) + 2, theBuffer.columns());
439  ++c) {
440  for (auto r = std::max(0, int(acluster.x[curInd]) - 1);
441  r < std::min(int(acluster.x[curInd]) + 2, theBuffer.rows());
442  ++r) {
443  if (theBuffer(r, c) >= thePixelThreshold) {
444  SiPixelCluster::PixelPos newpix(r, c);
445  if (!acluster.add(newpix, theBuffer(r, c)))
446  goto endClus;
447  // VV: no fake pixels in cluster, leads to non-contiguous clusters
448  if (!theFakePixels[r * theNumOfCols + c]) {
449  cldata.add(newpix, theBuffer(r, c));
450  }
451  theBuffer.set_adc(newpix, 1);
452  }
453 
454  /* //Commenting out the addition of dead pixels to the cluster until further testing -- dfehling 06/09
455  //Check on the bounds of the module; this is to keep the isDead and isNoisy modules from returning errors
456  else if(r>= 0 && c >= 0 && (r <= (theNumOfRows-1.)) && (c <= (theNumOfCols-1.))){
457  //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.
458  if((theSiPixelGainCalibrationService_->isDead(theDetid,c,r) || theSiPixelGainCalibrationService_->isNoisy(theDetid,c,r)) && theBuffer(r,c) != 1){
459 
460  //If a pixel is dead or noisy, check to see if we want to split the clusters or not.
461  //Push it into a dead pixel stack in case we want to split the clusters. Otherwise add it to the cluster.
462  //If we are splitting the clusters, we will iterate over the dead pixel stack later.
463 
464  SiPixelCluster::PixelPos newpix(r,c);
465  if(!doSplitClusters){
466 
467  cluster.add(newpix, theBuffer(r,c));}
468  else if(doSplitClusters){
469  dead_pixel_stack.push(newpix);
470  dead_flag = true;}
471 
472  theBuffer.set_adc(newpix, 1);
473  }
474 
475  }
476  */
477  }
478  }
479 
480  } // while accretion
481 endClus:
482  SiPixelCluster cluster(cldata.isize, cldata.adc, cldata.x, cldata.y, cldata.xmin, cldata.ymin);
483  //Here we split the cluster, if the flag to do so is set and we have found a dead or noisy pixel.
484 
485  if (dead_flag && doSplitClusters) {
486  // Set separate cluster threshold for L1 (needed for phase1)
488  if (theLayer == 1)
490 
491  //Set the first cluster equal to the existing cluster.
492  SiPixelCluster first_cluster = cluster;
493  bool have_second_cluster = false;
494  while (!dead_pixel_stack.empty()) {
495  //consider each found dead pixel
496  SiPixelCluster::PixelPos deadpix = dead_pixel_stack.top();
497  dead_pixel_stack.pop();
498  theBuffer.set_adc(deadpix, 1);
499 
500  //Clusterize the split cluster using the dead pixel as a seed
501  SiPixelCluster second_cluster = make_cluster(deadpix, output);
502 
503  //If both clusters would normally have been found by the clusterizer, put them into output
504  if (second_cluster.charge() >= clusterThreshold && first_cluster.charge() >= clusterThreshold) {
505  output.push_back(second_cluster);
506  have_second_cluster = true;
507  }
508 
509  //We also want to keep the merged cluster in data and let the RecHit algorithm decide which set to keep
510  //This loop adds the second cluster to the first.
511  const std::vector<SiPixelCluster::Pixel>& branch_pixels = second_cluster.pixels();
512  for (unsigned int i = 0; i < branch_pixels.size(); i++) {
513  int temp_x = branch_pixels[i].x;
514  int temp_y = branch_pixels[i].y;
515  int temp_adc = branch_pixels[i].adc;
516  SiPixelCluster::PixelPos newpix(temp_x, temp_y);
517  cluster.add(newpix, temp_adc);
518  }
519  }
520 
521  //Remember to also add the first cluster if we added the second one.
522  if (first_cluster.charge() >= clusterThreshold && have_second_cluster) {
523  output.push_back(first_cluster);
524  std::push_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
525  return cl1.minPixelRow() < cl2.minPixelRow();
526  });
527  }
528  }
529 
530  return cluster;
531 }
PixelClusterizerBase::DigiIterator
edm::DetSet< PixelDigi >::const_iterator DigiIterator
Definition: PixelClusterizerBase.h:20
hgcalPlots.ncols
ncols
Definition: hgcalPlots.py:104
electrons_cff.bool
bool
Definition: electrons_cff.py:366
mps_fire.i
i
Definition: mps_fire.py:428
input
static const std::string input
Definition: EdmProvDump.cc:48
PixelClusterizerBase::ClusterIterator
edmNew::DetSet< SiPixelCluster >::const_iterator ClusterIterator
Definition: PixelClusterizerBase.h:21
PixelClusterizerBase::AccretionCluster::isize
unsigned int isize
Definition: PixelClusterizerBase.h:30
PixelTopology.h
SiPixelArrayBuffer::rows
int rows() const
Definition: SiPixelArrayBuffer.h:33
PixelClusterizerBase::AccretionCluster::adc
uint16_t adc[MAXSIZE]
Definition: PixelClusterizerBase.h:25
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
PixelThresholdClusterizer::setup
bool setup(const PixelGeomDetUnit *pixDet)
Private helper methods:
Definition: PixelThresholdClusterizer.cc:93
SiPixelArrayBuffer::set_adc
void set_adc(int row, int col, int adc)
Definition: SiPixelArrayBuffer.h:67
min
T min(T a, T b)
Definition: MathUtil.h:58
PixelThresholdClusterizer::theNumOfRows
int theNumOfRows
Geometry-related information.
Definition: PixelThresholdClusterizer.h:117
SiPixelGainCalibrationServiceBase::isDead
virtual bool isDead(const uint32_t &detID, const int &col, const int &row)=0
gpuClustering::adc
uint16_t *__restrict__ uint16_t const *__restrict__ adc
Definition: gpuClusterChargeCut.h:20
TrackerTopology
Definition: TrackerTopology.h:16
PixelThresholdClusterizer::theOffset
const int theOffset
Definition: PixelThresholdClusterizer.h:106
SiPixelGainCalibrationServiceBase::getGain
virtual float getGain(const uint32_t &detID, const int &col, const int &row)=0
PixelThresholdClusterizer::theClusterThreshold_L1
const int theClusterThreshold_L1
Definition: PixelThresholdClusterizer.h:103
cuy.col
col
Definition: cuy.py:1009
gather_cfg.cout
cout
Definition: gather_cfg.py:144
PixelClusterizerBase::theSiPixelGainCalibrationService_
SiPixelGainCalibrationServiceBase * theSiPixelGainCalibrationService_
Definition: PixelClusterizerBase.h:83
PixelThresholdClusterizer::theClusterThreshold
const int theClusterThreshold
Definition: PixelThresholdClusterizer.h:102
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
muonClassificationByHits_cfi.pixel
pixel
Definition: muonClassificationByHits_cfi.py:9
PixelThresholdClusterizer::theSeeds
std::vector< SiPixelCluster::PixelPos > theSeeds
Definition: PixelThresholdClusterizer.h:90
cms::cuda::assert
assert(be >=bs)
PixelThresholdClusterizer::theNumOfCols
int theNumOfCols
Definition: PixelThresholdClusterizer.h:118
PixelThresholdClusterizer::make_cluster
SiPixelCluster make_cluster(const SiPixelCluster::PixelPos &pix, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
The actual clustering algorithm: group the neighboring pixels around the seed.
Definition: PixelThresholdClusterizer.cc:400
PixelThresholdClusterizer::thePhase2KinkADC
const int thePhase2KinkADC
Definition: PixelThresholdClusterizer.h:114
GetRecoTauVFromDQM_MC_cff.cl2
cl2
Definition: GetRecoTauVFromDQM_MC_cff.py:44
SiPixelCluster
Pixel cluster – collection of neighboring pixels above threshold.
Definition: SiPixelCluster.h:28
SiPixelCluster::Pixel
Definition: SiPixelCluster.h:31
PixelClusterizerBase::AccretionCluster::ymin
uint16_t ymin
Definition: PixelClusterizerBase.h:29
PixelThresholdClusterizer::copy_to_buffer
void copy_to_buffer(DigiIterator begin, DigiIterator end)
Copy adc counts from PixelDigis into the buffer, identify seeds.
Definition: PixelThresholdClusterizer.cc:217
TrackerTopology::pxbLayer
unsigned int pxbLayer(const DetId &id) const
Definition: TrackerTopology.h:144
DetId
Definition: DetId.h:17
PixelThresholdClusterizer::~PixelThresholdClusterizer
~PixelThresholdClusterizer() override
Definition: PixelThresholdClusterizer.cc:68
SiPixelArrayBuffer::columns
int columns() const
Definition: SiPixelArrayBuffer.h:34
SiPixelCluster::add
void add(const PixelPos &pix, int adc)
Definition: SiPixelCluster.cc:25
PixelClusterizerBase::AccretionCluster::xmin
uint16_t xmin
Definition: PixelClusterizerBase.h:28
SiPixelArrayBuffer.h
PixelClusterizerBase::AccretionCluster::x
uint16_t x[MAXSIZE]
Definition: PixelClusterizerBase.h:26
SiPixelCluster::minPixelRow
int minPixelRow() const
Definition: SiPixelCluster.h:150
PixelGeomDetUnit
Definition: PixelGeomDetUnit.h:15
mps_fire.end
end
Definition: mps_fire.py:242
PixelTopology::ncolumns
virtual int ncolumns() const =0
PixelClusterizerBase::AccretionCluster::pop
void pop()
Definition: PixelClusterizerBase.h:36
PixelThresholdClusterizer::PixelThresholdClusterizer
PixelThresholdClusterizer(edm::ParameterSet const &conf)
Definition: PixelThresholdClusterizer.cc:43
EMEnrichingFilter_cfi.clusterThreshold
clusterThreshold
Definition: EMEnrichingFilter_cfi.py:12
PixelThresholdClusterizer::thePhase2DigiBaseline
const double thePhase2DigiBaseline
Definition: PixelThresholdClusterizer.h:113
PixelTopology
Definition: PixelTopology.h:10
PixelThresholdClusterizer::thePhase2ReadoutMode
const int thePhase2ReadoutMode
Definition: PixelThresholdClusterizer.h:112
PixelThresholdClusterizer::clusterizeDetUnitT
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...
Definition: PixelThresholdClusterizer.cc:126
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
SiPixelGainCalibrationOffline.h
PixelThresholdClusterizer::theLayer
int theLayer
Definition: PixelThresholdClusterizer.h:120
PixelThresholdClusterizer::thePixelThreshold
const int thePixelThreshold
Definition: PixelThresholdClusterizer.h:100
edm::ParameterSet
Definition: ParameterSet.h:47
PixelClusterizerBase::AccretionCluster
Definition: PixelClusterizerBase.h:23
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
PixelGeomDetUnit::specificTopology
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
Definition: PixelGeomDetUnit.cc:17
PixelThresholdClusterizer::clear_buffer
void clear_buffer(DigiIterator begin, DigiIterator end)
Clear the internal buffer array.
Definition: PixelThresholdClusterizer.cc:198
recoMuon::in
Definition: RecoMuonEnumerators.h:6
PixelClusterizerBase::AccretionCluster::empty
bool empty()
Definition: PixelClusterizerBase.h:37
PixelThresholdClusterizer::theFakePixels
std::vector< bool > theFakePixels
Definition: PixelThresholdClusterizer.h:93
createfilelist.int
int
Definition: createfilelist.py:10
SiPixelArrayBuffer::add_adc
void add_adc(int row, int col, int adc)
Definition: SiPixelArrayBuffer.h:71
PixelThresholdClusterizer.h
PixelThresholdClusterizer::theOffset_L1
const int theOffset_L1
Definition: PixelThresholdClusterizer.h:107
EcalCondDBWriter_cfi.pedestal
pedestal
Definition: EcalCondDBWriter_cfi.py:49
PixelThresholdClusterizer::theConversionFactor
const int theConversionFactor
Definition: PixelThresholdClusterizer.h:104
PixelThresholdClusterizer::theConversionFactor_L1
const int theConversionFactor_L1
Definition: PixelThresholdClusterizer.h:105
HPSPFTauProducerPuppi_cfi.electron
electron
Definition: HPSPFTauProducerPuppi_cfi.py:13
alignCSCRings.r
r
Definition: alignCSCRings.py:93
PixelThresholdClusterizer::doSplitClusters
const bool doSplitClusters
Definition: PixelThresholdClusterizer.h:122
PedestalClient_cfi.gain
gain
Definition: PedestalClient_cfi.py:37
PixelClusterizerBase::AccretionCluster::add
bool add(SiPixelCluster::PixelPos const &p, uint16_t const iadc)
Definition: PixelClusterizerBase.h:47
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
T
long double T
Definition: Basic3DVectorLD.h:48
SiPixelCluster::PixelPos::row
constexpr int row() const
Definition: SiPixelCluster.h:60
PixelGeomDetUnit.h
PixelClusterizerBase::AccretionCluster::top
uint16_t top() const
Definition: PixelClusterizerBase.h:35
SiPixelCluster::charge
int charge() const
Definition: SiPixelCluster.h:142
PixelThresholdClusterizer::doPhase2Calibration
const bool doPhase2Calibration
Definition: PixelThresholdClusterizer.h:111
edmNew::DetSetVector::FastFiller
Definition: DetSetVectorNew.h:202
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
SiPixelCluster::pixels
const std::vector< Pixel > pixels() const
Definition: SiPixelCluster.h:159
PixelClusterizerBase::AccretionCluster::y
uint16_t y[MAXSIZE]
Definition: PixelClusterizerBase.h:27
SiPixelGainCalibrationServiceBase::getPedestal
virtual float getPedestal(const uint32_t &detID, const int &col, const int &row)=0
SiPixelArrayBuffer::setSize
void setSize(int rows, int cols)
Definition: SiPixelArrayBuffer.h:54
PixelThresholdClusterizer::doMissCalibrate
const bool doMissCalibrate
Definition: PixelThresholdClusterizer.h:121
c
auto & c
Definition: CAHitNtupletGeneratorKernelsImpl.h:56
PixelThresholdClusterizer::theSeedThreshold
const int theSeedThreshold
Definition: PixelThresholdClusterizer.h:101
SiPixelCluster::PixelPos
Definition: SiPixelCluster.h:56
PixelTopology::nrows
virtual int nrows() const =0
PixelThresholdClusterizer::theBuffer
SiPixelArrayBuffer theBuffer
Data storage.
Definition: PixelThresholdClusterizer.h:89
PixelThresholdClusterizer::calibrate
int calibrate(int adc, int col, int row)
Definition: PixelThresholdClusterizer.cc:319
PixelThresholdClusterizer::fillPSetDescription
static void fillPSetDescription(edm::ParameterSetDescription &desc)
Definition: PixelThresholdClusterizer.cc:71
PixelThresholdClusterizer::theElectronPerADCGain
const double theElectronPerADCGain
Definition: PixelThresholdClusterizer.h:109
SiPixelCluster::PixelPos::col
constexpr int col() const
Definition: SiPixelCluster.h:61
PixelThresholdClusterizer::theDetid
uint32_t theDetid
Definition: PixelThresholdClusterizer.h:119
SiPixelGainCalibrationServiceBase::isNoisy
virtual bool isNoisy(const uint32_t &detID, const int &col, const int &row)=0