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