CMS 3D CMS Logo

PixelThresholdClusterizer.h
Go to the documentation of this file.
1 #ifndef RecoLocalTracker_SiPixelClusterizer_PixelThresholdClusterizer_H
2 #define RecoLocalTracker_SiPixelClusterizer_PixelThresholdClusterizer_H
3 
4 //-----------------------------------------------------------------------
39 //-----------------------------------------------------------------------
40 
41 // Base class, defines SiPixelDigi and SiPixelCluster. The latter includes
42 // Pixel, PixelPos and Shift as inner classes.
43 //
45 #include "PixelClusterizerBase.h"
46 
47 // The private pixel buffer
48 #include "SiPixelArrayBuffer.h"
49 
50 // Parameter Set:
52 
55 
56 #include <vector>
57 
58 
60  public:
61 
64 
65  // Full I/O in DetSet
67  const PixelGeomDetUnit * pixDet,
68  const TrackerTopology* tTopo,
69  const std::vector<short>& badChannels,
70  edmNew::DetSetVector<SiPixelCluster>::FastFiller& output) { clusterizeDetUnitT(input, pixDet, tTopo, badChannels, output); }
72  const PixelGeomDetUnit * pixDet,
73  const TrackerTopology* tTopo,
74  const std::vector<short>& badChannels,
75  edmNew::DetSetVector<SiPixelCluster>::FastFiller& output) { clusterizeDetUnitT(input, pixDet, tTopo, badChannels, output); }
76 
77  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
78 
79  private:
80 
81  template<typename T>
82  void clusterizeDetUnitT( const T & input,
83  const PixelGeomDetUnit * pixDet,
84  const TrackerTopology* tTopo,
85  const std::vector<short>& badChannels,
87 
89  SiPixelArrayBuffer theBuffer; // internal nrow * ncol matrix
90  bool bufferAlreadySet; // status of the buffer array
91  std::vector<SiPixelCluster::PixelPos> theSeeds; // cached seed pixels
92  std::vector<SiPixelCluster> theClusters; // resulting clusters
93 
95  float thePixelThresholdInNoiseUnits; // Pixel threshold in units of noise
96  float theSeedThresholdInNoiseUnits; // Pixel cluster seed in units of noise
97  float theClusterThresholdInNoiseUnits; // Cluster threshold in units of noise
98 
99  const int thePixelThreshold; // Pixel threshold in electrons
100  const int theSeedThreshold; // Seed threshold in electrons
101  const int theClusterThreshold; // Cluster threshold in electrons
102  const int theClusterThreshold_L1; // Cluster threshold in electrons for Layer 1
103  const int theConversionFactor; // adc to electron conversion factor
104  const int theConversionFactor_L1; // adc to electron conversion factor for Layer 1
105  const int theOffset; // adc to electron conversion offset
106  const int theOffset_L1; // adc to electron conversion offset for Layer 1
107 
108  const int theStackADC_; // The maximum ADC count for the stack layers
109  const int theFirstStack_; // The index of the first stack layer
110  const double theElectronPerADCGain_; // ADC to electrons conversion
111 
115  uint32_t detid_;
116  int layer_;
117  bool dead_flag;
118  const bool doMissCalibrate; // Use calibration or not
119  const bool doSplitClusters;
121  bool setup(const PixelGeomDetUnit * pixDet);
122  void copy_to_buffer( DigiIterator begin, DigiIterator end );
123  void copy_to_buffer( ClusterIterator begin, ClusterIterator end );
124  void clear_buffer( DigiIterator begin, DigiIterator end );
125  void clear_buffer( ClusterIterator begin, ClusterIterator end );
127  // Calibrate the ADC charge to electrons
128  int calibrate(int adc, int col, int row);
129 
130 };
131 
132 #endif
int adc(sample_type sample)
get the ADC sample (12 bits)
#define dso_hidden
SiPixelArrayBuffer theBuffer
Data storage.
void clusterizeDetUnit(const edmNew::DetSet< SiPixelCluster > &input, const PixelGeomDetUnit *pixDet, const TrackerTopology *tTopo, const std::vector< short > &badChannels, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
static std::string const input
Definition: EdmProvDump.cc:44
edm::DetSet< PixelDigi >::const_iterator DigiIterator
#define end
Definition: vmac.h:37
int theNumOfRows
Geometry-related information.
edmNew::DetSet< SiPixelCluster >::const_iterator ClusterIterator
std::vector< SiPixelCluster > theClusters
Pixel cluster – collection of neighboring pixels above threshold.
float thePixelThresholdInNoiseUnits
Clustering-related quantities:
#define begin
Definition: vmac.h:30
col
Definition: cuy.py:1008
void clusterizeDetUnit(const edm::DetSet< PixelDigi > &input, const PixelGeomDetUnit *pixDet, const TrackerTopology *tTopo, const std::vector< short > &badChannels, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)
std::vector< SiPixelCluster::PixelPos > theSeeds
long double T
A specific threshold-based pixel clustering algorithm.
Class to store ADC counts during clustering.