Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:58

0001 #ifndef CSCRecHitD_CSCHitFromStripOnly_h
0002 #define CSCRecHitD_CSCHitFromStripOnly_h
0003 
0004 /** \class CSCHitFromStripOnly (old comment below)
0005  *
0006  * Search for strip with ADC output exceeding theThresholdForAPeak.  For each of these strips,
0007  * build a cluster of strip of size theClusterSize (typically 5 strips).  Finally, make
0008  * a Strip Hit out of these clusters by finding the center-of-mass position of the hit
0009  * The DetId, strip hit position, and peaking time are stored in a CSCStripHit collection.
0010  *
0011  * Be careful with the ME_1/a chambers: the 48 strips are ganged into 16 channels,
0012  * Each channel contains signals from three strips, each separated by 16 strips from the next.
0013  *
0014  */
0015 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
0016 
0017 #include "RecoLocalMuon/CSCRecHitD/src/CSCStripData.h"
0018 #include "RecoLocalMuon/CSCRecHitD/src/CSCStripHitData.h"
0019 #include "RecoLocalMuon/CSCRecHitD/src/CSCStripHit.h"
0020 #include "RecoLocalMuon/CSCRecHitD/src/CSCRecoConditions.h"
0021 
0022 #include "DataFormats/CSCDigi/interface/CSCStripDigiCollection.h"
0023 
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0026 
0027 #include <vector>
0028 #include <array>
0029 
0030 class CSCLayer;
0031 class CSCStripDigi;
0032 class CSCPedestalChoice;
0033 
0034 class CSCHitFromStripOnly {
0035 public:
0036   typedef std::array<CSCStripData, 100> PulseHeightMap;
0037 
0038   explicit CSCHitFromStripOnly(const edm::ParameterSet& ps);
0039 
0040   ~CSCHitFromStripOnly();
0041 
0042   std::vector<CSCStripHit> runStrip(const CSCDetId& id,
0043                                     const CSCLayer* layer,
0044                                     const CSCStripDigiCollection::Range& rstripd);
0045 
0046   void setConditions(const CSCRecoConditions* reco) { recoConditions_ = reco; }
0047 
0048   bool ganged() { return ganged_; }
0049   void setGanged(bool ig) { ganged_ = ig; }
0050 
0051 private:
0052   /// Store SCA pulseheight information from strips in digis of one layer
0053   void fillPulseHeights(const CSCStripDigiCollection::Range& rstripd);
0054 
0055   /// Find local maxima
0056   void findMaxima(const CSCDetId& id);
0057   // What we call a peak
0058   bool isPeakOK(int iStrip, float heightCluster);
0059 
0060   /// Make clusters using local maxima
0061   float makeCluster(int centerStrip);
0062 
0063   ///
0064   CSCStripHitData makeStripData(int centerStrip, int offset);
0065 
0066   /// Is either neighbour 'bad'?
0067   bool isNearDeadStrip(const CSCDetId& id, int centralStrip, int nstrips);
0068 
0069   /// Is the strip 'bad'?
0070   bool isDeadStrip(const CSCDetId& id, int centralStrip, int nstrips);
0071 
0072   /// Find position of hit in strip cluster in terms of strip #
0073   float findHitOnStripPosition(const std::vector<CSCStripHitData>& data, const int& centerStrip);
0074 
0075   // MEMBER DATA
0076 
0077   // Hold pointers to current layer, conditions data
0078   CSCDetId id_;
0079   const CSCLayer* layer_;
0080   const CSCRecoConditions* recoConditions_;
0081   // Number of strips in layer
0082   unsigned nstrips_;
0083   // gain correction weights and crosstalks read in from conditions database.
0084   float gainWeight[80];
0085 
0086   // The specific pedestal calculator
0087   CSCPedestalChoice* calcped_;
0088 
0089   // The cuts for forming the strip hits are described in the config file
0090   bool useCalib;
0091   static const int theClusterSize = 3;
0092   float theThresholdForAPeak;
0093   float theThresholdForCluster;
0094 
0095   // working buffer for sca pulseheights
0096   PulseHeightMap thePulseHeightMap;
0097 
0098   std::vector<int> theMaxima;
0099   std::vector<int> theConsecutiveStrips;  //... with charge for a given maximum
0100   std::vector<int> theClosestMaximum;     // this is number of strips to the closest other maximum
0101 
0102   // Variables entering the CSCStripHit construction:
0103   int tmax_cluster;  // Peaking time for strip hit, in time bin units
0104   int clusterSize;
0105   std::vector<float> strips_adc;
0106   std::vector<float> strips_adcRaw;
0107   std::vector<int> theStrips;
0108 
0109   bool ganged_;  // only True if ME1/1A AND it is ganged
0110 };
0111 
0112 #endif