CMS 3D CMS Logo

CSCHitFromStripOnly.h
Go to the documentation of this file.
1 #ifndef CSCRecHitD_CSCHitFromStripOnly_h
2 #define CSCRecHitD_CSCHitFromStripOnly_h
3 
16 
21 
23 
26 
27 #include <vector>
28 #include <array>
29 
30 class CSCLayer;
31 class CSCStripDigi;
32 class CSCPedestalChoice;
33 
35 public:
36  typedef std::array<CSCStripData, 100> PulseHeightMap;
37 
38  explicit CSCHitFromStripOnly(const edm::ParameterSet& ps);
39 
41 
42  std::vector<CSCStripHit> runStrip(const CSCDetId& id,
43  const CSCLayer* layer,
44  const CSCStripDigiCollection::Range& rstripd);
45 
47 
48  bool ganged() { return ganged_; }
49  void setGanged(bool ig) { ganged_ = ig; }
50 
51 private:
54 
56  void findMaxima(const CSCDetId& id);
57  // What we call a peak
58  bool isPeakOK(int iStrip, float heightCluster);
59 
61  float makeCluster(int centerStrip);
62 
64  CSCStripHitData makeStripData(int centerStrip, int offset);
65 
67  bool isNearDeadStrip(const CSCDetId& id, int centralStrip, int nstrips);
68 
70  bool isDeadStrip(const CSCDetId& id, int centralStrip, int nstrips);
71 
73  float findHitOnStripPosition(const std::vector<CSCStripHitData>& data, const int& centerStrip);
74 
75  // MEMBER DATA
76 
77  // Hold pointers to current layer, conditions data
79  const CSCLayer* layer_;
81  // Number of strips in layer
82  unsigned nstrips_;
83  // gain correction weights and crosstalks read in from conditions database.
84  float gainWeight[80];
85 
86  // The specific pedestal calculator
88 
89  // The cuts for forming the strip hits are described in the config file
90  bool useCalib;
91  static const int theClusterSize = 3;
94 
95  // working buffer for sca pulseheights
97 
98  std::vector<int> theMaxima;
99  std::vector<int> theConsecutiveStrips; //... with charge for a given maximum
100  std::vector<int> theClosestMaximum; // this is number of strips to the closest other maximum
101 
102  // Variables entering the CSCStripHit construction:
103  int tmax_cluster; // Peaking time for strip hit, in time bin units
105  std::vector<float> strips_adc;
106  std::vector<float> strips_adcRaw;
107  std::vector<int> theStrips;
108 
109  bool ganged_; // only True if ME1/1A AND it is ganged
110 };
111 
112 #endif
const CSCLayer * layer_
CSCPedestalChoice * calcped_
bool isNearDeadStrip(const CSCDetId &id, int centralStrip, int nstrips)
Is either neighbour &#39;bad&#39;?
void findMaxima(const CSCDetId &id)
Find local maxima.
std::vector< int > theConsecutiveStrips
bool isDeadStrip(const CSCDetId &id, int centralStrip, int nstrips)
Is the strip &#39;bad&#39;?
std::vector< CSCStripHit > runStrip(const CSCDetId &id, const CSCLayer *layer, const CSCStripDigiCollection::Range &rstripd)
std::vector< int > theStrips
std::vector< int > theClosestMaximum
float findHitOnStripPosition(const std::vector< CSCStripHitData > &data, const int &centerStrip)
Find position of hit in strip cluster in terms of strip #.
CSCHitFromStripOnly(const edm::ParameterSet &ps)
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
static const int theClusterSize
PulseHeightMap thePulseHeightMap
const CSCRecoConditions * recoConditions_
float makeCluster(int centerStrip)
Make clusters using local maxima.
std::pair< const_iterator, const_iterator > Range
std::vector< float > strips_adcRaw
std::vector< int > theMaxima
bool isPeakOK(int iStrip, float heightCluster)
void setConditions(const CSCRecoConditions *reco)
fixed size matrix
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
std::vector< float > strips_adc
std::array< CSCStripData, 100 > PulseHeightMap
void fillPulseHeights(const CSCStripDigiCollection::Range &rstripd)
Store SCA pulseheight information from strips in digis of one layer.
CSCStripHitData makeStripData(int centerStrip, int offset)