CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes
PixelThresholdClusterizer Class Reference

A specific threshold-based pixel clustering algorithm. More...

#include <PixelThresholdClusterizer.h>

Inheritance diagram for PixelThresholdClusterizer:
PixelClusterizerBase

Public Member Functions

void clusterizeDetUnit (const edm::DetSet< PixelDigi > &input, const PixelGeomDetUnit *pixDet, const TrackerTopology *tTopo, const std::vector< short > &badChannels, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output) override
 
void clusterizeDetUnit (const edmNew::DetSet< SiPixelCluster > &input, const PixelGeomDetUnit *pixDet, const TrackerTopology *tTopo, const std::vector< short > &badChannels, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output) override
 
 PixelThresholdClusterizer (edm::ParameterSet const &conf)
 
 ~PixelThresholdClusterizer () override
 
- Public Member Functions inherited from PixelClusterizerBase
void setSiPixelGainCalibrationService (SiPixelGainCalibrationServiceBase *in)
 
virtual ~PixelClusterizerBase ()
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &desc)
 

Protected Member Functions

int calibrate (int adc, int col, int row)
 
void clear_buffer (DigiIterator begin, DigiIterator end)
 Clear the internal buffer array. More...
 
void clear_buffer (ClusterIterator begin, ClusterIterator end)
 
template<typename T >
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 around each seed pixel. Input and output data stored in DetSet. More...
 
void copy_to_buffer (DigiIterator begin, DigiIterator end)
 Copy adc counts from PixelDigis into the buffer, identify seeds. More...
 
void copy_to_buffer (ClusterIterator begin, ClusterIterator end)
 
SiPixelCluster make_cluster (const SiPixelCluster::PixelPos &pix, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
 The actual clustering algorithm: group the neighboring pixels around the seed. More...
 
bool setup (const PixelGeomDetUnit *pixDet)
 Private helper methods: More...
 

Protected Attributes

const bool doMissCalibrate
 
const bool doPhase2Calibration
 
const bool doSplitClusters
 
const bool dropDuplicates
 
SiPixelArrayBuffer theBuffer
 Data storage. More...
 
std::vector< SiPixelClustertheClusters
 
const int theClusterThreshold
 
const int theClusterThreshold_L1
 
float theClusterThresholdInNoiseUnits
 
const int theConversionFactor
 
const int theConversionFactor_L1
 
uint32_t theDetid
 
const double theElectronPerADCGain
 
std::vector< bool > theFakePixels
 
int theLayer
 
int theNumOfCols
 
int theNumOfRows
 Geometry-related information. More...
 
const int theOffset
 
const int theOffset_L1
 
const double thePhase2DigiBaseline
 
const int thePhase2KinkADC
 
const int thePhase2ReadoutMode
 
std::vector< uint8_t > thePixelOccurrence
 
const int thePixelThreshold
 
float thePixelThresholdInNoiseUnits
 Clustering-related quantities: More...
 
std::vector< SiPixelCluster::PixelPostheSeeds
 
const int theSeedThreshold
 
float theSeedThresholdInNoiseUnits
 
- Protected Attributes inherited from PixelClusterizerBase
SiPixelGainCalibrationServiceBasetheSiPixelGainCalibrationService_
 

Additional Inherited Members

- Public Types inherited from PixelClusterizerBase
typedef edmNew::DetSet< SiPixelCluster >::const_iterator ClusterIterator
 
typedef edm::DetSet< PixelDigi >::const_iterator DigiIterator
 

Detailed Description

A specific threshold-based pixel clustering algorithm.

An explicit threshold-based clustering algorithm.

General logic of PixelThresholdClusterizer:

The clusterization is performed on a matrix with size equal to the size of the pixel detector, each cell containing the ADC count of the corresponding pixel. The matrix is reset after each clusterization.

The search starts from seed pixels, i.e. pixels with sufficiently large amplitudes, found at the time of filling of the matrix and stored in a SiPixelArrayBuffer.

Translate the pixel charge to electrons, we are suppose to do the calibrations ADC->electrons here. Modify the thresholds to be in electrons, convert adc to electrons. d.k. 20/3/06 Get rid of the noiseVector. d.k. 28/3/06

A threshold-based clustering algorithm which clusters SiPixelDigis into SiPixelClusters for each DetUnit. The algorithm is straightforward and purely topological: the clustering process starts with seed pixels and continues by adding adjacent pixels above the pixel threshold. Once the cluster is made, it has to be above the cluster threshold as well.

The clusterization is performed on a matrix with size equal to the size of the pixel detector, each cell containing the ADC count of the corresponding pixel. The matrix is reset after each clusterization.

The search starts from seed pixels, i.e. pixels with sufficiently large amplitudes, found at the time of filling of the matrix and stored in a

At this point the noise and dead channels are ignored, but soon they won't be.

SiPixelCluster contains a barrycenter, but it should be noted that that information is largely useless. One must use a PositionEstimator class to compute the RecHit position and its error for every given cluster.

Author
Largely copied from NewPixelClusterizer in ORCA written by Danek Kotlinski (PSI). Ported to CMSSW by Petar Maksimovic (JHU). DetSetVector data container implemented by V.Chiochia (Uni Zurich)

Sets the PixelArrayBuffer dimensions and pixel thresholds. Makes clusters and stores them in theCache if the option useCache has been set.

Definition at line 57 of file PixelThresholdClusterizer.h.

Constructor & Destructor Documentation

◆ PixelThresholdClusterizer()

PixelThresholdClusterizer::PixelThresholdClusterizer ( edm::ParameterSet const &  conf)

Constructor: Initilize the buffer to hold pixels from a detector module. This is a vector of 44k ints, stays valid all the time.

Definition at line 44 of file PixelThresholdClusterizer.cc.

References SiPixelArrayBuffer::setSize(), theBuffer, theFakePixels, theNumOfCols, theNumOfRows, and thePixelOccurrence.

45  : // Get thresholds in electrons
46  thePixelThreshold(conf.getParameter<int>("ChannelThreshold")),
47  theSeedThreshold(conf.getParameter<int>("SeedThreshold")),
48  theClusterThreshold(conf.getParameter<int>("ClusterThreshold")),
49  theClusterThreshold_L1(conf.getParameter<int>("ClusterThreshold_L1")),
50  theConversionFactor(conf.getParameter<int>("VCaltoElectronGain")),
51  theConversionFactor_L1(conf.getParameter<int>("VCaltoElectronGain_L1")),
52  theOffset(conf.getParameter<int>("VCaltoElectronOffset")),
53  theOffset_L1(conf.getParameter<int>("VCaltoElectronOffset_L1")),
54  theElectronPerADCGain(conf.getParameter<double>("ElectronPerADCGain")),
55  doPhase2Calibration(conf.getParameter<bool>("Phase2Calibration")),
56  dropDuplicates(conf.getParameter<bool>("DropDuplicates")),
57  thePhase2ReadoutMode(conf.getParameter<int>("Phase2ReadoutMode")),
58  thePhase2DigiBaseline(conf.getParameter<double>("Phase2DigiBaseline")),
59  thePhase2KinkADC(conf.getParameter<int>("Phase2KinkADC")),
60  theNumOfRows(0),
61  theNumOfCols(0),
62  theDetid(0),
63  // Get the constants for the miss-calibration studies
64  doMissCalibrate(conf.getParameter<bool>("MissCalibrate")),
65  doSplitClusters(conf.getParameter<bool>("SplitClusters")) {
67  theFakePixels.clear();
68  thePixelOccurrence.clear();
69 }
SiPixelArrayBuffer theBuffer
Data storage.
int theNumOfRows
Geometry-related information.
std::vector< uint8_t > thePixelOccurrence
void setSize(int rows, int cols)

◆ ~PixelThresholdClusterizer()

PixelThresholdClusterizer::~PixelThresholdClusterizer ( )
override

Definition at line 71 of file PixelThresholdClusterizer.cc.

71 {}

Member Function Documentation

◆ calibrate()

int PixelThresholdClusterizer::calibrate ( int  adc,
int  col,
int  row 
)
protected

Definition at line 352 of file PixelThresholdClusterizer.cc.

References gpuClustering::adc, cuy::col, doMissCalibrate, doPhase2Calibration, pwdgSkimBPark_cfi::electrons, PedestalClient_cfi::gain, SiPixelGainCalibrationServiceBase::getGain(), SiPixelGainCalibrationServiceBase::getPedestal(), createfilelist::int, SiPixelGainCalibrationServiceBase::isDead(), SiPixelGainCalibrationServiceBase::isNoisy(), EcalCondDBWriter_cfi::pedestal, funct::pow(), theConversionFactor, theConversionFactor_L1, theDetid, theElectronPerADCGain, theLayer, theOffset, theOffset_L1, thePhase2DigiBaseline, thePhase2KinkADC, thePhase2ReadoutMode, and PixelClusterizerBase::theSiPixelGainCalibrationService_.

Referenced by copy_to_buffer().

352  {
353  int electrons = 0;
354 
355  if (doPhase2Calibration) {
356  const float gain = theElectronPerADCGain;
357  int p2rm = (thePhase2ReadoutMode < -1 ? -1 : thePhase2ReadoutMode);
358 
359  if (p2rm == -1) {
360  electrons = int(adc * gain);
361  } else {
362  if (adc < thePhase2KinkADC) {
363  electrons = int((adc + 0.5) * gain);
364  } else {
365  const int dualslopeparam = (thePhase2ReadoutMode < 10 ? thePhase2ReadoutMode : 10);
366  const int dualslope = int(dualslopeparam <= 1 ? 1. : pow(2, dualslopeparam - 1));
368  adc *= dualslope;
370  electrons = int((adc + 0.5 * dualslope) * gain);
371  }
373  }
374 
375  return electrons;
376  }
377 
378  if (doMissCalibrate) {
379  // do not perform calibration if pixel is dead!
380 
383  // Linear approximation of the TANH response
384  // Pixel(0,0,0)
385  //const float gain = 2.95; // 1 ADC = 2.95 VCALs (1/0.339)
386  //const float pedestal = -83.; // -28/0.339
387  // Roc-0 average
388  //const float gain = 1./0.357; // 1 ADC = 2.80 VCALs
389  //const float pedestal = -28.2 * gain; // -79.
390 
391  float DBgain = theSiPixelGainCalibrationService_->getGain(theDetid, col, row);
393  float DBpedestal = pedestal * DBgain;
394 
395  // Roc-6 average
396  //const float gain = 1./0.313; // 1 ADC = 3.19 VCALs
397  //const float pedestal = -6.2 * gain; // -19.8
398  //
399  float vcal = adc * DBgain - DBpedestal;
400 
401  // atanh calibration
402  // Roc-6 average
403  //const float p0 = 0.00492;
404  //const float p1 = 1.998;
405  //const float p2 = 90.6;
406  //const float p3 = 134.1;
407  // Roc-6 average
408  //const float p0 = 0.00382;
409  //const float p1 = 0.886;
410  //const float p2 = 112.7;
411  //const float p3 = 113.0;
412  //float vcal = ( atanh( (adc-p3)/p2) + p1)/p0;
413 
414  if (theLayer == 1) {
416  } else {
418  }
419  }
420  } else { // No misscalibration in the digitizer
421  // Simple (default) linear gain
422  const float gain = theElectronPerADCGain; // default: 1 ADC = 135 electrons
423  const float pedestal = 0.; //
424  electrons = int(adc * gain + pedestal);
425  }
426 
427  return electrons;
428 }
virtual bool isDead(const uint32_t &detID, const int &col, const int &row)=0
virtual float getPedestal(const uint32_t &detID, const int &col, const int &row)=0
SiPixelGainCalibrationServiceBase * theSiPixelGainCalibrationService_
col
Definition: cuy.py:1009
virtual float getGain(const uint32_t &detID, const int &col, const int &row)=0
virtual bool isNoisy(const uint32_t &detID, const int &col, const int &row)=0
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ clear_buffer() [1/2]

void PixelThresholdClusterizer::clear_buffer ( DigiIterator  begin,
DigiIterator  end 
)
protected

Clear the internal buffer array.

Pixels which are not part of recognized clusters are NOT ERASED during the cluster finding. Erase them now.

TO DO: ask Danek... wouldn't it be faster to simply memcopy() zeros into the whole buffer array?

Definition at line 207 of file PixelThresholdClusterizer.cc.

References mps_fire::end, SiPixelArrayBuffer::set_adc(), and theBuffer.

Referenced by clusterizeDetUnitT().

207  {
208  for (DigiIterator di = begin; di != end; ++di) {
209  theBuffer.set_adc(di->row(), di->column(), 0); // reset pixel adc to 0
210  }
211 }
SiPixelArrayBuffer theBuffer
Data storage.
void set_adc(int row, int col, int adc)
edm::DetSet< PixelDigi >::const_iterator DigiIterator

◆ clear_buffer() [2/2]

void PixelThresholdClusterizer::clear_buffer ( ClusterIterator  begin,
ClusterIterator  end 
)
protected

Definition at line 213 of file PixelThresholdClusterizer.cc.

References mps_fire::end, mps_fire::i, muonClassificationByHits_cfi::pixel, SiPixelArrayBuffer::set_adc(), and theBuffer.

213  {
214  for (ClusterIterator ci = begin; ci != end; ++ci) {
215  for (int i = 0; i < ci->size(); ++i) {
216  const SiPixelCluster::Pixel pixel = ci->pixel(i);
217 
218  theBuffer.set_adc(pixel.x, pixel.y, 0); // reset pixel adc to 0
219  }
220  }
221 }
SiPixelArrayBuffer theBuffer
Data storage.
void set_adc(int row, int col, int adc)
edmNew::DetSet< SiPixelCluster >::const_iterator ClusterIterator

◆ clusterizeDetUnit() [1/2]

void PixelThresholdClusterizer::clusterizeDetUnit ( const edm::DetSet< PixelDigi > &  input,
const PixelGeomDetUnit pixDet,
const TrackerTopology tTopo,
const std::vector< short > &  badChannels,
edmNew::DetSetVector< SiPixelCluster >::FastFiller &  output 
)
inlineoverridevirtual

Implements PixelClusterizerBase.

Definition at line 63 of file PixelThresholdClusterizer.h.

References input, and convertSQLitetoXML_cfg::output.

67  {
68  clusterizeDetUnitT(input, pixDet, tTopo, badChannels, output);
69  }
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 std::string const input
Definition: EdmProvDump.cc:50

◆ clusterizeDetUnit() [2/2]

void PixelThresholdClusterizer::clusterizeDetUnit ( const edmNew::DetSet< SiPixelCluster > &  input,
const PixelGeomDetUnit pixDet,
const TrackerTopology tTopo,
const std::vector< short > &  badChannels,
edmNew::DetSetVector< SiPixelCluster >::FastFiller &  output 
)
inlineoverridevirtual

Implements PixelClusterizerBase.

Definition at line 70 of file PixelThresholdClusterizer.h.

References input, and convertSQLitetoXML_cfg::output.

74  {
75  clusterizeDetUnitT(input, pixDet, tTopo, badChannels, output);
76  }
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 std::string const input
Definition: EdmProvDump.cc:50

◆ clusterizeDetUnitT()

template<typename T >
void PixelThresholdClusterizer::clusterizeDetUnitT ( const T input,
const PixelGeomDetUnit pixDet,
const TrackerTopology tTopo,
const std::vector< short > &  badChannels,
edmNew::DetSetVector< SiPixelCluster >::FastFiller &  output 
)
protected

Cluster pixels. This method operates on a matrix of pixels and finds the largest contiguous cluster around each seed pixel. Input and output data stored in DetSet.

Definition at line 131 of file PixelThresholdClusterizer.cc.

References cms::cuda::assert(), SiPixelCluster::charge(), GetRecoTauVFromDQM_MC_cff::cl2, clear_buffer(), EMEnrichingFilter_cfi::clusterThreshold, copy_to_buffer(), mps_fire::end, mps_fire::i, input, make_cluster(), SiPixelCluster::minPixelRow(), eostools::move(), convertSQLitetoXML_cfg::output, TrackerTopology::pxbLayer(), setup(), DetId::subdetId(), theBuffer, theClusterThreshold, theClusterThreshold_L1, theDetid, theFakePixels, theLayer, thePixelOccurrence, theSeeds, and theSeedThreshold.

135  {
136  typename T::const_iterator begin = input.begin();
137  typename T::const_iterator end = input.end();
138 
139  // this should never happen and the raw2digi does not create empty detsets
140  if (begin == end) {
141  edm::LogError("PixelThresholdClusterizer") << "@SUB=PixelThresholdClusterizer::clusterizeDetUnitT()"
142  << " No digis to clusterize";
143  }
144 
145  // Set up the clusterization on this DetId.
146  if (!setup(pixDet))
147  return;
148 
149  theDetid = input.detId();
150 
151  // Set separate cluster threshold for L1 (needed for phase1)
153  theLayer = (DetId(theDetid).subdetId() == 1) ? tTopo->pxbLayer(theDetid) : 0;
154  if (theLayer == 1)
156 
157  // Copy PixelDigis to the buffer array; select the seed pixels
158  // on the way, and store them in theSeeds.
159  if (end > begin)
160  copy_to_buffer(begin, end);
161 
162  assert(output.empty());
163  // Loop over all seeds. TO DO: wouldn't using iterators be faster?
164  for (unsigned int i = 0; i < theSeeds.size(); i++) {
165  // Gavril : The charge of seeds that were already inlcuded in clusters is set to 1 electron
166  // so we don't want to call "make_cluster" for these cases
167  if (theBuffer(theSeeds[i]) >= theSeedThreshold) { // Is this seed still valid?
168  // Make a cluster around this seed
170 
171  // Check if the cluster is above threshold
172  // (TO DO: one is signed, other unsigned, gcc warns...)
173  if (cluster.charge() >= clusterThreshold) {
174  // sort by row (x)
175  output.push_back(std::move(cluster));
176  std::push_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
177  return cl1.minPixelRow() < cl2.minPixelRow();
178  });
179  }
180  }
181  }
182  // sort by row (x) maybe sorting the seed would suffice....
183  std::sort_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
184  return cl1.minPixelRow() < cl2.minPixelRow();
185  });
186 
187  // Erase the seeds.
188  theSeeds.clear();
189 
190  // Need to clean unused pixels from the buffer array.
191  clear_buffer(begin, end);
192 
193  theFakePixels.clear();
194 
195  thePixelOccurrence.clear();
196 }
unsigned int pxbLayer(const DetId &id) const
SiPixelArrayBuffer theBuffer
Data storage.
SiPixelCluster make_cluster(const SiPixelCluster::PixelPos &pix, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
The actual clustering algorithm: group the neighboring pixels around the seed.
Log< level::Error, false > LogError
assert(be >=bs)
static std::string const input
Definition: EdmProvDump.cc:50
int minPixelRow() const
int charge() const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
std::vector< uint8_t > thePixelOccurrence
Definition: DetId.h:17
void clear_buffer(DigiIterator begin, DigiIterator end)
Clear the internal buffer array.
Pixel cluster – collection of neighboring pixels above threshold.
bool setup(const PixelGeomDetUnit *pixDet)
Private helper methods:
std::vector< SiPixelCluster::PixelPos > theSeeds
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.

◆ copy_to_buffer() [1/2]

void PixelThresholdClusterizer::copy_to_buffer ( DigiIterator  begin,
DigiIterator  end 
)
protected

Copy adc counts from PixelDigis into the buffer, identify seeds.

Definition at line 226 of file PixelThresholdClusterizer.cc.

References gpuClustering::adc, cms::cuda::assert(), calibrate(), cuy::col, gather_cfg::cout, doMissCalibrate, doPhase2Calibration, dropDuplicates, HPSPFTauProducerPuppi_cfi::electron, mps_fire::end, PedestalClient_cfi::gain, mps_fire::i, recoMuon::in, SiPixelArrayBuffer::index(), createfilelist::int, EcalCondDBWriter_cfi::pedestal, MatrixUtil::remove(), SiPixelArrayBuffer::set_adc(), theBuffer, theConversionFactor, theConversionFactor_L1, theDetid, theElectronPerADCGain, theFakePixels, theLayer, theNumOfCols, theOffset, theOffset_L1, thePixelOccurrence, thePixelThreshold, theSeeds, and theSeedThreshold.

Referenced by clusterizeDetUnitT().

226  {
227 #ifdef PIXELREGRESSION
228  static std::atomic<int> s_ic = 0;
229  in ic = ++s_ic;
230  if (ic == 1) {
231  // std::cout << (doMissCalibrate ? "VI from db" : "VI linear") << std::endl;
232  }
233 #endif
234 
235  //If called with empty/invalid DetSet, warn the user
236  if (end <= begin) {
237  edm::LogWarning("PixelThresholdClusterizer") << " copy_to_buffer called with empty or invalid range" << std::endl;
238  return;
239  }
240 
241  int electron[end - begin]; // pixel charge in electrons
242  memset(electron, 0, (end - begin) * sizeof(int));
243 
244  if (doPhase2Calibration) {
245  int i = 0;
246  for (DigiIterator di = begin; di != end; ++di) {
247  electron[i] = calibrate(di->adc(), di->column(), di->row());
248  i++;
249  }
250  assert(i == (end - begin));
251  }
252 
253  else {
254  if (doMissCalibrate) {
255  if (theLayer == 1) {
256  (*theSiPixelGainCalibrationService_)
258  } else {
259  (*theSiPixelGainCalibrationService_).calibrate(theDetid, begin, end, theConversionFactor, theOffset, electron);
260  }
261  } else {
262  int i = 0;
263  const float gain = theElectronPerADCGain; // default: 1 ADC = 135 electrons
264  for (DigiIterator di = begin; di != end; ++di) {
265  auto adc = di->adc();
266  const float pedestal = 0.; //
267  electron[i] = int(adc * gain + pedestal);
268  ++i;
269  }
270  assert(i == (end - begin));
271  }
272  }
273 
274  int i = 0;
275 #ifdef PIXELREGRESSION
276  static std::atomic<int> eqD = 0;
277 #endif
278  for (DigiIterator di = begin; di != end; ++di) {
279  int row = di->row();
280  int col = di->column();
281  // VV: do not calibrate a fake pixel, it already has a unit of 10e-:
282  int adc = (di->flag() != 0) ? di->adc() * 10 : electron[i]; // this is in electrons
283  i++;
284 
285 #ifdef PIXELREGRESSION
286  int adcOld = calibrate(di->adc(), col, row);
287  //assert(adc==adcOld);
288  if (adc != adcOld)
289  std::cout << "VI " << eqD << ' ' << ic << ' ' << end - begin << ' ' << i << ' ' << di->adc() << ' ' << adc << ' '
290  << adcOld << std::endl;
291  else
292  ++eqD;
293 #endif
294 
295  if (adc < 100)
296  adc = 100; // put all negative pixel charges into the 100 elec bin
297  /* This is semi-random good number. The exact number (in place of 100) is irrelevant from the point
298  of view of the final cluster charge since these are typically >= 20000.
299  */
300 
301  thePixelOccurrence[theBuffer.index(row, col)]++; // increment the occurrence counter
302  uint8_t occurrence =
303  (!dropDuplicates) ? 1 : thePixelOccurrence[theBuffer.index(row, col)]; // get the occurrence counter
304 
305  switch (occurrence) {
306  // the 1st occurrence (standard treatment)
307  case 1:
308  if (adc >= thePixelThreshold) {
309  theBuffer.set_adc(row, col, adc);
310  // VV: add pixel to the fake list. Only when running on digi collection
311  if (di->flag() != 0)
312  theFakePixels[row * theNumOfCols + col] = true;
313  if (adc >= theSeedThreshold)
314  theSeeds.push_back(SiPixelCluster::PixelPos(row, col));
315  }
316  break;
317 
318  // the 2nd occurrence (duplicate pixel: reset the buffer to 0 and remove from the list of seed pixels)
319  case 2:
320  theBuffer.set_adc(row, col, 0);
322  break;
323 
324  // in case a pixel appears more than twice, nothing needs to be done because it was already removed at the 2nd occurrence
325  }
326  }
327  assert(i == (end - begin));
328 }
SiPixelArrayBuffer theBuffer
Data storage.
assert(be >=bs)
void set_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
std::vector< uint8_t > thePixelOccurrence
def remove(d, key, TELL=False)
Definition: MatrixUtil.py:223
col
Definition: cuy.py:1009
int calibrate(int adc, int col, int row)
Log< level::Warning, false > LogWarning
std::vector< SiPixelCluster::PixelPos > theSeeds
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ copy_to_buffer() [2/2]

void PixelThresholdClusterizer::copy_to_buffer ( ClusterIterator  begin,
ClusterIterator  end 
)
protected

Definition at line 330 of file PixelThresholdClusterizer.cc.

References gpuClustering::adc, SiPixelArrayBuffer::add_adc(), cuy::col, mps_fire::end, mps_fire::i, muonClassificationByHits_cfi::pixel, theBuffer, thePixelThreshold, theSeeds, and theSeedThreshold.

330  {
331  // loop over clusters
332  for (ClusterIterator ci = begin; ci != end; ++ci) {
333  // loop over pixels
334  for (int i = 0; i < ci->size(); ++i) {
335  const SiPixelCluster::Pixel pixel = ci->pixel(i);
336 
337  int row = pixel.x;
338  int col = pixel.y;
339  int adc = pixel.adc;
340  if (adc >= thePixelThreshold) {
341  theBuffer.add_adc(row, col, adc);
342  if (adc >= theSeedThreshold)
343  theSeeds.push_back(SiPixelCluster::PixelPos(row, col));
344  }
345  }
346  }
347 }
SiPixelArrayBuffer theBuffer
Data storage.
void add_adc(int row, int col, int adc)
edmNew::DetSet< SiPixelCluster >::const_iterator ClusterIterator
col
Definition: cuy.py:1009
std::vector< SiPixelCluster::PixelPos > theSeeds
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ fillPSetDescription()

void PixelThresholdClusterizer::fillPSetDescription ( edm::ParameterSetDescription desc)
static

Definition at line 74 of file PixelThresholdClusterizer.cc.

References submitPVResolutionJobs::desc.

Referenced by SiPixelClusterProducer::fillDescriptions().

74  {
75  desc.add<int>("ChannelThreshold", 1000);
76  desc.add<bool>("MissCalibrate", true);
77  desc.add<bool>("SplitClusters", false);
78  desc.add<int>("VCaltoElectronGain", 65);
79  desc.add<int>("VCaltoElectronGain_L1", 65);
80  desc.add<int>("VCaltoElectronOffset", -414);
81  desc.add<int>("VCaltoElectronOffset_L1", -414);
82  desc.add<int>("SeedThreshold", 1000);
83  desc.add<int>("ClusterThreshold_L1", 4000);
84  desc.add<int>("ClusterThreshold", 4000);
85  desc.add<double>("ElectronPerADCGain", 135.);
86  desc.add<bool>("DropDuplicates", true);
87  desc.add<bool>("Phase2Calibration", false);
88  desc.add<int>("Phase2ReadoutMode", -1);
89  desc.add<double>("Phase2DigiBaseline", 1200.);
90  desc.add<int>("Phase2KinkADC", 8);
91 }

◆ make_cluster()

SiPixelCluster PixelThresholdClusterizer::make_cluster ( const SiPixelCluster::PixelPos pix,
edmNew::DetSetVector< SiPixelCluster >::FastFiller &  output 
)
protected

The actual clustering algorithm: group the neighboring pixels around the seed.

Definition at line 433 of file PixelThresholdClusterizer.cc.

References PixelClusterizerBase::AccretionCluster::adc, PixelClusterizerBase::AccretionCluster::add(), SiPixelCluster::add(), c, SiPixelCluster::charge(), GetRecoTauVFromDQM_MC_cff::cl2, EMEnrichingFilter_cfi::clusterThreshold, SiPixelCluster::PixelPos::col(), SiPixelArrayBuffer::columns(), doSplitClusters, PixelClusterizerBase::AccretionCluster::empty(), mps_fire::i, PixelClusterizerBase::AccretionCluster::isize, SiStripPI::max, SiStripPI::min, SiPixelCluster::minPixelRow(), convertSQLitetoXML_cfg::output, SiPixelCluster::pixels(), PixelClusterizerBase::AccretionCluster::pop(), alignCSCRings::r, SiPixelCluster::PixelPos::row(), SiPixelArrayBuffer::rows(), SiPixelArrayBuffer::set_adc(), theBuffer, theClusterThreshold, theClusterThreshold_L1, theFakePixels, theLayer, theNumOfCols, thePixelThreshold, PixelClusterizerBase::AccretionCluster::top(), PixelClusterizerBase::AccretionCluster::x, PixelClusterizerBase::AccretionCluster::xmin, PixelClusterizerBase::AccretionCluster::y, and PixelClusterizerBase::AccretionCluster::ymin.

Referenced by clusterizeDetUnitT().

434  {
435  //First we acquire the seeds for the clusters
436  uint16_t seed_adc;
437  std::stack<SiPixelCluster::PixelPos, std::vector<SiPixelCluster::PixelPos> > dead_pixel_stack;
438 
439  //The individual modules have been loaded into a buffer.
440  //After each pixel has been considered by the clusterizer, we set the adc count to 1
441  //to mark that we have already considered it.
442  //The only difference between dead/noisy pixels and standard ones is that for dead/noisy pixels,
443  //We consider the charge of the pixel to always be zero.
444 
445  /* this is not possible as dead and noisy pixel cannot make it into a seed...
446  if ( doMissCalibrate &&
447  (theSiPixelGainCalibrationService_->isDead(theDetid,pix.col(),pix.row()) ||
448  theSiPixelGainCalibrationService_->isNoisy(theDetid,pix.col(),pix.row())) )
449  {
450  std::cout << "IMPOSSIBLE" << std::endl;
451  seed_adc = 0;
452  theBuffer.set_adc(pix, 1);
453  }
454  else {
455  */
456  // Note: each ADC value is limited here to 65535 (std::numeric_limits<uint16_t>::max),
457  // as it is later stored as uint16_t in SiPixelCluster and PixelClusterizerBase/AccretionCluster
458  // (reminder: ADC values here may be expressed in number of electrons)
459  seed_adc = std::min(theBuffer(pix.row(), pix.col()), int(std::numeric_limits<uint16_t>::max()));
460  theBuffer.set_adc(pix, 1);
461  // }
462 
463  AccretionCluster acluster, cldata;
464  acluster.add(pix, seed_adc);
465  cldata.add(pix, seed_adc);
466 
467  //Here we search all pixels adjacent to all pixels in the cluster.
468  bool dead_flag = false;
469  while (!acluster.empty()) {
470  //This is the standard algorithm to find and add a pixel
471  auto curInd = acluster.top();
472  acluster.pop();
473  for (auto c = std::max(0, int(acluster.y[curInd]) - 1);
474  c < std::min(int(acluster.y[curInd]) + 2, theBuffer.columns());
475  ++c) {
476  for (auto r = std::max(0, int(acluster.x[curInd]) - 1);
477  r < std::min(int(acluster.x[curInd]) + 2, theBuffer.rows());
478  ++r) {
479  if (theBuffer(r, c) >= thePixelThreshold) {
480  SiPixelCluster::PixelPos newpix(r, c);
481  auto const newpix_adc = std::min(theBuffer(r, c), int(std::numeric_limits<uint16_t>::max()));
482  if (!acluster.add(newpix, newpix_adc))
483  goto endClus;
484  // VV: no fake pixels in cluster, leads to non-contiguous clusters
485  if (!theFakePixels[r * theNumOfCols + c]) {
486  cldata.add(newpix, newpix_adc);
487  }
488  theBuffer.set_adc(newpix, 1);
489  }
490 
491  /* //Commenting out the addition of dead pixels to the cluster until further testing -- dfehling 06/09
492  //Check on the bounds of the module; this is to keep the isDead and isNoisy modules from returning errors
493  else if(r>= 0 && c >= 0 && (r <= (theNumOfRows-1.)) && (c <= (theNumOfCols-1.))){
494  //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.
495  if((theSiPixelGainCalibrationService_->isDead(theDetid,c,r) || theSiPixelGainCalibrationService_->isNoisy(theDetid,c,r)) && theBuffer(r,c) != 1){
496 
497  //If a pixel is dead or noisy, check to see if we want to split the clusters or not.
498  //Push it into a dead pixel stack in case we want to split the clusters. Otherwise add it to the cluster.
499  //If we are splitting the clusters, we will iterate over the dead pixel stack later.
500 
501  SiPixelCluster::PixelPos newpix(r,c);
502  if(!doSplitClusters){
503 
504  cluster.add(newpix, std::min(theBuffer(r, c), int(std::numeric_limits<uint16_t>::max())));}
505  else if(doSplitClusters){
506  dead_pixel_stack.push(newpix);
507  dead_flag = true;}
508 
509  theBuffer.set_adc(newpix, 1);
510  }
511 
512  }
513  */
514  }
515  }
516 
517  } // while accretion
518 endClus:
519  SiPixelCluster cluster(cldata.isize, cldata.adc, cldata.x, cldata.y, cldata.xmin, cldata.ymin);
520  //Here we split the cluster, if the flag to do so is set and we have found a dead or noisy pixel.
521 
522  if (dead_flag && doSplitClusters) {
523  // Set separate cluster threshold for L1 (needed for phase1)
525  if (theLayer == 1)
527 
528  //Set the first cluster equal to the existing cluster.
529  SiPixelCluster first_cluster = cluster;
530  bool have_second_cluster = false;
531  while (!dead_pixel_stack.empty()) {
532  //consider each found dead pixel
533  SiPixelCluster::PixelPos deadpix = dead_pixel_stack.top();
534  dead_pixel_stack.pop();
535  theBuffer.set_adc(deadpix, 1);
536 
537  //Clusterize the split cluster using the dead pixel as a seed
538  SiPixelCluster second_cluster = make_cluster(deadpix, output);
539 
540  //If both clusters would normally have been found by the clusterizer, put them into output
541  if (second_cluster.charge() >= clusterThreshold && first_cluster.charge() >= clusterThreshold) {
542  output.push_back(second_cluster);
543  have_second_cluster = true;
544  }
545 
546  //We also want to keep the merged cluster in data and let the RecHit algorithm decide which set to keep
547  //This loop adds the second cluster to the first.
548  const std::vector<SiPixelCluster::Pixel>& branch_pixels = second_cluster.pixels();
549  for (unsigned int i = 0; i < branch_pixels.size(); i++) {
550  auto const temp_x = branch_pixels[i].x;
551  auto const temp_y = branch_pixels[i].y;
552  auto const temp_adc = branch_pixels[i].adc;
553  SiPixelCluster::PixelPos newpix(temp_x, temp_y);
554  cluster.add(newpix, temp_adc);
555  }
556  }
557 
558  //Remember to also add the first cluster if we added the second one.
559  if (first_cluster.charge() >= clusterThreshold && have_second_cluster) {
560  output.push_back(first_cluster);
561  std::push_heap(output.begin(), output.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
562  return cl1.minPixelRow() < cl2.minPixelRow();
563  });
564  }
565  }
566 
567  return cluster;
568 }
SiPixelArrayBuffer theBuffer
Data storage.
SiPixelCluster make_cluster(const SiPixelCluster::PixelPos &pix, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
The actual clustering algorithm: group the neighboring pixels around the seed.
void set_adc(int row, int col, int adc)
int minPixelRow() const
int charge() const
const std::vector< Pixel > pixels() const
Pixel cluster – collection of neighboring pixels above threshold.
constexpr int row() const
constexpr int col() const

◆ setup()

bool PixelThresholdClusterizer::setup ( const PixelGeomDetUnit pixDet)
protected

Private helper methods:

Prepare the Clusterizer to work on a particular DetUnit. Re-init the size of the panel/plaquette (so update nrows and ncols),

Definition at line 97 of file PixelThresholdClusterizer.cc.

References SiPixelArrayBuffer::columns(), hgcalPlots::ncols, PixelTopology::ncolumns(), PixelTopology::nrows(), SiPixelArrayBuffer::rows(), SiPixelArrayBuffer::setSize(), PixelGeomDetUnit::specificTopology(), theBuffer, theFakePixels, theNumOfCols, theNumOfRows, and thePixelOccurrence.

Referenced by clusterizeDetUnitT().

97  {
98  // Cache the topology.
99  const PixelTopology& topol = pixDet->specificTopology();
100 
101  // Get the new sizes.
102  int nrows = topol.nrows(); // rows in x
103  int ncols = topol.ncolumns(); // cols in y
104 
105  theNumOfRows = nrows; // Set new sizes
107 
108  if (nrows > theBuffer.rows() || ncols > theBuffer.columns()) { // change only when a larger is needed
109  if (nrows != theNumOfRows || ncols != theNumOfCols)
110  edm::LogWarning("setup()") << "pixel buffer redefined to" << nrows << " * " << ncols;
111  //theNumOfRows = nrows; // Set new sizes
112  //theNumOfCols = ncols;
113  // Resize the buffer
114  theBuffer.setSize(nrows, ncols); // Modify
115  }
116 
117  theFakePixels.resize(nrows * ncols, false);
118 
119  thePixelOccurrence.resize(nrows * ncols, 0);
120 
121  return true;
122 }
SiPixelArrayBuffer theBuffer
Data storage.
virtual int ncolumns() const =0
virtual int nrows() const =0
int theNumOfRows
Geometry-related information.
std::vector< uint8_t > thePixelOccurrence
void setSize(int rows, int cols)
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
Log< level::Warning, false > LogWarning

Member Data Documentation

◆ doMissCalibrate

const bool PixelThresholdClusterizer::doMissCalibrate
protected

Definition at line 124 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), and copy_to_buffer().

◆ doPhase2Calibration

const bool PixelThresholdClusterizer::doPhase2Calibration
protected

Definition at line 113 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), and copy_to_buffer().

◆ doSplitClusters

const bool PixelThresholdClusterizer::doSplitClusters
protected

Definition at line 125 of file PixelThresholdClusterizer.h.

Referenced by make_cluster().

◆ dropDuplicates

const bool PixelThresholdClusterizer::dropDuplicates
protected

Definition at line 114 of file PixelThresholdClusterizer.h.

Referenced by copy_to_buffer().

◆ theBuffer

SiPixelArrayBuffer PixelThresholdClusterizer::theBuffer
protected

◆ theClusters

std::vector<SiPixelCluster> PixelThresholdClusterizer::theClusters
protected

Definition at line 91 of file PixelThresholdClusterizer.h.

◆ theClusterThreshold

const int PixelThresholdClusterizer::theClusterThreshold
protected

Definition at line 104 of file PixelThresholdClusterizer.h.

Referenced by clusterizeDetUnitT(), and make_cluster().

◆ theClusterThreshold_L1

const int PixelThresholdClusterizer::theClusterThreshold_L1
protected

Definition at line 105 of file PixelThresholdClusterizer.h.

Referenced by clusterizeDetUnitT(), and make_cluster().

◆ theClusterThresholdInNoiseUnits

float PixelThresholdClusterizer::theClusterThresholdInNoiseUnits
protected

Definition at line 100 of file PixelThresholdClusterizer.h.

◆ theConversionFactor

const int PixelThresholdClusterizer::theConversionFactor
protected

Definition at line 106 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), and copy_to_buffer().

◆ theConversionFactor_L1

const int PixelThresholdClusterizer::theConversionFactor_L1
protected

Definition at line 107 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), and copy_to_buffer().

◆ theDetid

uint32_t PixelThresholdClusterizer::theDetid
protected

Definition at line 122 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), clusterizeDetUnitT(), and copy_to_buffer().

◆ theElectronPerADCGain

const double PixelThresholdClusterizer::theElectronPerADCGain
protected

Definition at line 111 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), and copy_to_buffer().

◆ theFakePixels

std::vector<bool> PixelThresholdClusterizer::theFakePixels
protected

◆ theLayer

int PixelThresholdClusterizer::theLayer
protected

◆ theNumOfCols

int PixelThresholdClusterizer::theNumOfCols
protected

◆ theNumOfRows

int PixelThresholdClusterizer::theNumOfRows
protected

Geometry-related information.

Definition at line 120 of file PixelThresholdClusterizer.h.

Referenced by PixelThresholdClusterizer(), and setup().

◆ theOffset

const int PixelThresholdClusterizer::theOffset
protected

Definition at line 108 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), and copy_to_buffer().

◆ theOffset_L1

const int PixelThresholdClusterizer::theOffset_L1
protected

Definition at line 109 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), and copy_to_buffer().

◆ thePhase2DigiBaseline

const double PixelThresholdClusterizer::thePhase2DigiBaseline
protected

Definition at line 116 of file PixelThresholdClusterizer.h.

Referenced by calibrate().

◆ thePhase2KinkADC

const int PixelThresholdClusterizer::thePhase2KinkADC
protected

Definition at line 117 of file PixelThresholdClusterizer.h.

Referenced by calibrate().

◆ thePhase2ReadoutMode

const int PixelThresholdClusterizer::thePhase2ReadoutMode
protected

Definition at line 115 of file PixelThresholdClusterizer.h.

Referenced by calibrate().

◆ thePixelOccurrence

std::vector<uint8_t> PixelThresholdClusterizer::thePixelOccurrence
protected

◆ thePixelThreshold

const int PixelThresholdClusterizer::thePixelThreshold
protected

Definition at line 102 of file PixelThresholdClusterizer.h.

Referenced by copy_to_buffer(), and make_cluster().

◆ thePixelThresholdInNoiseUnits

float PixelThresholdClusterizer::thePixelThresholdInNoiseUnits
protected

Clustering-related quantities:

Definition at line 98 of file PixelThresholdClusterizer.h.

◆ theSeeds

std::vector<SiPixelCluster::PixelPos> PixelThresholdClusterizer::theSeeds
protected

Definition at line 90 of file PixelThresholdClusterizer.h.

Referenced by clusterizeDetUnitT(), and copy_to_buffer().

◆ theSeedThreshold

const int PixelThresholdClusterizer::theSeedThreshold
protected

Definition at line 103 of file PixelThresholdClusterizer.h.

Referenced by clusterizeDetUnitT(), and copy_to_buffer().

◆ theSeedThresholdInNoiseUnits

float PixelThresholdClusterizer::theSeedThresholdInNoiseUnits
protected

Definition at line 99 of file PixelThresholdClusterizer.h.