CMS 3D CMS Logo

PixelThresholdClusterizer Class Reference

An explicit threshold-based clustering algorithm. More...

#include <RecoLocalTracker/SiPixelClusterizer/interface/PixelThresholdClusterizer.h>

Inheritance diagram for PixelThresholdClusterizer:

PixelClusterizerBase

List of all members.

Public Member Functions

void clusterizeDetUnit (const edm::DetSet< PixelDigi > &input, const PixelGeomDetUnit *pixDet, const std::vector< short > &badChannels, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
 Cluster pixels.
 PixelThresholdClusterizer (edm::ParameterSet const &conf)
 Constructor: Initilize the buffer to hold pixels from a detector module.
 ~PixelThresholdClusterizer ()

Private Member Functions

int calibrate (int adc, int col, int row)
void clear_buffer (DigiIterator begin, DigiIterator end)
 Clear the internal buffer array.
void copy_to_buffer (DigiIterator begin, DigiIterator end)
 Copy adc counts from PixelDigis into the buffer, identify seeds.
SiPixelCluster make_cluster (const SiPixelCluster::PixelPos &pix)
 The actual clustering algorithm: group the neighboring pixels around the seed.
bool setup (const PixelGeomDetUnit *pixDet)
 Private helper methods:.

Private Attributes

bool bufferAlreadySet
edm::ParameterSet conf_
uint32_t detid_
bool doMissCalibrate
SiPixelArrayBuffer theBuffer
 Data storage.
std::vector< SiPixelClustertheClusters
float theClusterThreshold
float theClusterThresholdInNoiseUnits
int theConversionFactor
int theNumOfCols
int theNumOfRows
 Geometry-related information.
int theOffset
int thePixelThreshold
float thePixelThresholdInNoiseUnits
 Clustering-related quantities:.
std::vector
< SiPixelCluster::PixelPos
theSeeds
int theSeedThreshold
float theSeedThresholdInNoiseUnits


Detailed Description

An explicit threshold-based clustering algorithm.

A specific threshold-based pixel clustering algorithm.

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.

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

Definition at line 59 of file PixelThresholdClusterizer.h.


Constructor & Destructor Documentation

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 43 of file PixelThresholdClusterizer.cc.

00043                                 :
00044     conf_(conf), bufferAlreadySet(false), theNumOfRows(0), theNumOfCols(0), detid_(0) {
00045 
00046    // Get thresholds in electrons
00047    thePixelThreshold   = 
00048      conf_.getParameter<int>("ChannelThreshold");
00049    theSeedThreshold    = 
00050      conf_.getParameter<int>("SeedThreshold");
00051    theClusterThreshold = 
00052      conf_.getParameter<double>("ClusterThreshold");
00053    theConversionFactor = 
00054      conf_.getParameter<int>("VCaltoElectronGain");
00055    theOffset = 
00056      conf_.getParameter<int>("VCaltoElectronOffset");
00057         
00058    
00059    // Get the constants for the miss-calibration studies
00060    doMissCalibrate=conf_.getUntrackedParameter<bool>("MissCalibrate",true); 
00061 
00062    theBuffer.setSize( theNumOfRows, theNumOfCols );
00063    //initTiming();
00064 }

PixelThresholdClusterizer::~PixelThresholdClusterizer (  ) 

Definition at line 66 of file PixelThresholdClusterizer.cc.

00066 {}


Member Function Documentation

int PixelThresholdClusterizer::calibrate ( int  adc,
int  col,
int  row 
) [private]

Definition at line 184 of file PixelThresholdClusterizer.cc.

References detid_, doMissCalibrate, SiPixelGainCalibrationServiceBase::getGain(), SiPixelGainCalibrationServiceBase::getPedestal(), int, SiPixelGainCalibrationServiceBase::isDead(), pedestal, theConversionFactor, theOffset, and PixelClusterizerBase::theSiPixelGainCalibrationService_.

Referenced by copy_to_buffer().

00184                                                                   {
00185   int electrons = 0;
00186 
00187   if(doMissCalibrate) {
00188     // do not perform calibration if pixel is dead!
00189 
00190     if(!theSiPixelGainCalibrationService_->isDead(detid_,col,row)){
00191 
00192       // Linear approximation of the TANH response
00193       // Pixel(0,0,0)
00194       //const float gain = 2.95; // 1 ADC = 2.95 VCALs (1/0.339)
00195       //const float pedestal = -83.; // -28/0.339
00196       // Roc-0 average
00197       //const float gain = 1./0.357; // 1 ADC = 2.80 VCALs 
00198       //const float pedestal = -28.2 * gain; // -79.
00199 
00200       float DBgain     = theSiPixelGainCalibrationService_->getGain(detid_, col, row);
00201       float DBpedestal = theSiPixelGainCalibrationService_->getPedestal(detid_, col, row) * DBgain;
00202 
00203       
00204       // Roc-6 average
00205       //const float gain = 1./0.313; // 1 ADC = 3.19 VCALs 
00206       //const float pedestal = -6.2 * gain; // -19.8
00207       // 
00208       float vcal = adc * DBgain - DBpedestal;
00209       
00210       // atanh calibration 
00211       // Roc-6 average
00212       //const float p0 = 0.00492;
00213       //const float p1 = 1.998;
00214       //const float p2 = 90.6;
00215       //const float p3 = 134.1; 
00216       // Roc-6 average
00217       //const float p0 = 0.00382;
00218       //const float p1 = 0.886;
00219       //const float p2 = 112.7;
00220       //const float p3 = 113.0; 
00221       //float vcal = ( atanh( (adc-p3)/p2) + p1)/p0;
00222   
00223       electrons = int( vcal * theConversionFactor + theOffset); 
00224     }
00225   } 
00226   else { // No misscalibration in the digitizer
00227     // Simple (default) linear gain 
00228     const float gain = 135.; // 1 ADC = 135 electrons
00229     const float pedestal = 0.; //
00230     electrons = int(adc * gain + pedestal);
00231   }
00232   
00233   return electrons;
00234 }

void PixelThresholdClusterizer::clear_buffer ( DigiIterator  begin,
DigiIterator  end 
) [private]

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 156 of file PixelThresholdClusterizer.cc.

References SiPixelArrayBuffer::set_adc(), and theBuffer.

Referenced by clusterizeDetUnit().

00156                                                                                    {
00157   // TimeMe tm1( *theClearTimer, false);
00158   for(DigiIterator di = begin; di != end; ++di ) {
00159     theBuffer.set_adc( di->row(), di->column(), 0 );   // reset pixel adc to 0
00160   }
00161 }

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

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

Implements PixelClusterizerBase.

Definition at line 101 of file PixelThresholdClusterizer.cc.

References edm::DetSet< T >::begin(), begin, SiPixelCluster::charge(), clear_buffer(), copy_to_buffer(), edm::DetSet< T >::detId(), detid_, edm::DetSet< T >::end(), end, i, make_cluster(), edmNew::DetSetVector< T >::push_back(), setup(), theBuffer, theClusterThreshold, and theSeeds.

00104                                                                                                          {
00105 
00106   //TimeMe tm1( *theClustersTimer, false);
00107 
00108    DigiIterator begin = input.begin();
00109    DigiIterator end   = input.end();
00110 
00111    // Do not bother for empty detectors
00112    //if (begin == end) cout << " PixelThresholdClusterizer::clusterizeDetUnit - No digis to clusterize";
00113 
00114    //  Set up the clusterization on this DetId.
00115    if (!setup(pixDet)) return;
00116    detid_ = input.detId();
00117    //  Copy PixelDigis to the buffer array; select the seed pixels
00118    //  on the way, and store them in theSeeds.
00119    copy_to_buffer(begin, end);
00120   //  At this point we know the number of seeds on this DetUnit, and thus
00121   //  also the maximal number of possible clusters, so resize theClusters
00122   //  in order to make vector<>::push_back() efficient.
00123   // output.reserve ( theSeeds.size() ); //GPetruc: It is better *not* to reserve, with the new DetSetVector!
00124     
00125   //  Loop over all seeds.  TO DO: wouldn't using iterators be faster?
00126   for (unsigned int i = 0; i < theSeeds.size(); i++) {
00127       
00128     if ( theBuffer(theSeeds[i]) != 0) {  // Is this seed still valid?
00129         
00130       //  Make a cluster around this seed
00131       SiPixelCluster cluster = make_cluster( theSeeds[i] );
00132         
00133       //  Check if the cluster is above threshold  
00134       // (TO DO: one is signed, other unsigned, gcc warns...)
00135       if ( cluster.charge() >= theClusterThreshold ) {
00136         output.push_back( cluster );
00137       }
00138     }
00139   }
00140   // Erase the seeds.
00141   theSeeds.clear();
00142 
00143   //  Need to clean unused pixels from the buffer array.
00144   clear_buffer(begin, end);
00145 
00146 }

void PixelThresholdClusterizer::copy_to_buffer ( DigiIterator  begin,
DigiIterator  end 
) [private]

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

Definition at line 166 of file PixelThresholdClusterizer.cc.

References ecalMGPA::adc(), calibrate(), parsecf::pyparsing::col(), row, SiPixelArrayBuffer::set_adc(), theBuffer, thePixelThreshold, theSeeds, and theSeedThreshold.

Referenced by clusterizeDetUnit().

00166                                                                                      {
00167   //TimeMe tm1( *theCopyTimer, false);
00168   for(DigiIterator di = begin; di != end; ++di) {
00169     int row = di->row();
00170     int col = di->column();
00171     int adc = calibrate(di->adc(),col,row); // convert ADC -> electrons
00172     if ( adc >= thePixelThreshold) {
00173       theBuffer.set_adc( row, col, adc);
00174       if ( adc >= theSeedThreshold) { 
00175         theSeeds.push_back( SiPixelCluster::PixelPos(row,col) );
00176       }
00177     }
00178   }
00179 }

SiPixelCluster PixelThresholdClusterizer::make_cluster ( const SiPixelCluster::PixelPos pix  )  [private]

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

Definition at line 239 of file PixelThresholdClusterizer.cc.

References SiPixelCluster::add(), c, SiPixelCluster::PixelPos::col(), r, SiPixelCluster::PixelPos::row(), SiPixelArrayBuffer::set_adc(), theBuffer, and thePixelThreshold.

Referenced by clusterizeDetUnit().

00240 {
00241   //TimeMe tm1( *theMakeClustTimer, false);
00242 
00243   // Make the cluster
00244   SiPixelCluster cluster( pix, theBuffer( pix.row(), pix.col()) );
00245   stack<SiPixelCluster::PixelPos, vector<SiPixelCluster::PixelPos> > pixel_stack;
00246 
00247   theBuffer.set_adc( pix, 0);
00248   pixel_stack.push( pix);
00249 
00250   while ( ! pixel_stack.empty()) {
00251     SiPixelCluster::PixelPos curpix = pixel_stack.top(); pixel_stack.pop();
00252     for ( int r = curpix.row()-1; r <= curpix.row()+1; ++r) {
00253       for ( int c = curpix.col()-1; c <= curpix.col()+1; ++c) {
00254           if ( theBuffer(r,c) >= thePixelThreshold) {
00255             SiPixelCluster::PixelPos newpix(r,c);
00256             cluster.add( newpix, theBuffer(r,c));
00257             theBuffer.set_adc( newpix, 0);
00258             pixel_stack.push( newpix);
00259           
00260         }
00261       }
00262     }
00263   }
00264   return cluster;
00265 }

bool PixelThresholdClusterizer::setup ( const PixelGeomDetUnit pixDet  )  [private]

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 72 of file PixelThresholdClusterizer.cc.

References bufferAlreadySet, PixelTopology::ncolumns(), PixelTopology::nrows(), SiPixelArrayBuffer::setSize(), PixelGeomDetUnit::specificTopology(), theBuffer, theNumOfCols, and theNumOfRows.

Referenced by clusterizeDetUnit().

00072                                                                      {
00073 
00074   // Cache the topology.
00075   const PixelTopology & topol = pixDet->specificTopology();
00076 
00077   // Get the new sizes.
00078   int nrows = topol.nrows();      // rows in x
00079   int ncols = topol.ncolumns();   // cols in y
00080 
00081   if( nrows > theNumOfRows || ncols > theNumOfCols ) { // change only when a larger is needed
00082     //if( nrows != theNumOfRows || ncols != theNumOfCols ) {
00083     //cout << " PixelThresholdClusterizer: pixel buffer redefined to " 
00084     // << nrows << " * " << ncols << endl;      
00085     theNumOfRows = nrows;  // Set new sizes
00086     theNumOfCols = ncols;
00087     // Resize the buffer
00088     theBuffer.setSize(nrows,ncols);  // Modify
00089     bufferAlreadySet = true;
00090   }
00091 
00092   return true;   
00093 }


Member Data Documentation

bool PixelThresholdClusterizer::bufferAlreadySet [private]

Definition at line 78 of file PixelThresholdClusterizer.h.

Referenced by setup().

edm::ParameterSet PixelThresholdClusterizer::conf_ [private]

Definition at line 74 of file PixelThresholdClusterizer.h.

uint32_t PixelThresholdClusterizer::detid_ [private]

Definition at line 96 of file PixelThresholdClusterizer.h.

Referenced by calibrate(), and clusterizeDetUnit().

bool PixelThresholdClusterizer::doMissCalibrate [private]

Definition at line 98 of file PixelThresholdClusterizer.h.

Referenced by calibrate().

SiPixelArrayBuffer PixelThresholdClusterizer::theBuffer [private]

Data storage.

Definition at line 77 of file PixelThresholdClusterizer.h.

Referenced by clear_buffer(), clusterizeDetUnit(), copy_to_buffer(), make_cluster(), and setup().

std::vector<SiPixelCluster> PixelThresholdClusterizer::theClusters [private]

Definition at line 80 of file PixelThresholdClusterizer.h.

float PixelThresholdClusterizer::theClusterThreshold [private]

Definition at line 89 of file PixelThresholdClusterizer.h.

Referenced by clusterizeDetUnit().

float PixelThresholdClusterizer::theClusterThresholdInNoiseUnits [private]

Definition at line 85 of file PixelThresholdClusterizer.h.

int PixelThresholdClusterizer::theConversionFactor [private]

Definition at line 90 of file PixelThresholdClusterizer.h.

Referenced by calibrate().

int PixelThresholdClusterizer::theNumOfCols [private]

Definition at line 95 of file PixelThresholdClusterizer.h.

Referenced by setup().

int PixelThresholdClusterizer::theNumOfRows [private]

Geometry-related information.

Definition at line 94 of file PixelThresholdClusterizer.h.

Referenced by setup().

int PixelThresholdClusterizer::theOffset [private]

Definition at line 91 of file PixelThresholdClusterizer.h.

Referenced by calibrate().

int PixelThresholdClusterizer::thePixelThreshold [private]

Definition at line 87 of file PixelThresholdClusterizer.h.

Referenced by copy_to_buffer(), and make_cluster().

float PixelThresholdClusterizer::thePixelThresholdInNoiseUnits [private]

Clustering-related quantities:.

Definition at line 83 of file PixelThresholdClusterizer.h.

std::vector<SiPixelCluster::PixelPos> PixelThresholdClusterizer::theSeeds [private]

Definition at line 79 of file PixelThresholdClusterizer.h.

Referenced by clusterizeDetUnit(), and copy_to_buffer().

int PixelThresholdClusterizer::theSeedThreshold [private]

Definition at line 88 of file PixelThresholdClusterizer.h.

Referenced by copy_to_buffer().

float PixelThresholdClusterizer::theSeedThresholdInNoiseUnits [private]

Definition at line 84 of file PixelThresholdClusterizer.h.


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:30:05 2009 for CMSSW by  doxygen 1.5.4