CMS 3D CMS Logo

SiStripAPVRestorer.h
Go to the documentation of this file.
1 #ifndef RECOLOCALTRACKER_SISTRIPZEROSUPPRESSION_SISTRIPAPVRESTORER_H
2 #define RECOLOCALTRACKER_SISTRIPZEROSUPPRESSION_SISTRIPAPVRESTORER_H
3 
21 
22 #include <vector>
23 #include <cstdint>
24 
27 
28 protected:
30 
31 public:
32  virtual ~SiStripAPVRestorer(){};
33 
34  using digi_t = int16_t;
35  using digivector_t = std::vector<digi_t>;
36  using digimap_t = std::map<uint16_t, digi_t>;
37  using medians_t = std::vector<std::pair<short, float>>;
38  using baselinemap_t = std::map<uint16_t, digivector_t>;
39 
40  void init(const edm::EventSetup& es);
41 
42  uint16_t inspect(uint32_t detId, uint16_t firstAPV, const digivector_t& digis, const medians_t& vmedians);
43  void restore(uint16_t firstAPV, digivector_t& digis);
44 
45  uint16_t inspectAndRestore(uint32_t detId,
46  uint16_t firstAPV,
47  const digivector_t& rawDigisPedSubtracted,
48  digivector_t& processedRawDigi,
49  const medians_t& vmedians);
50 
51  void loadMeanCMMap(const edm::Event&);
52 
53  const baselinemap_t& getBaselineMap() const { return baselineMap_; }
54  const std::map<uint16_t, digimap_t>& getSmoothedPoints() const { return smoothedMaps_; }
55  const std::vector<bool>& getAPVFlags() const { return apvFlagsBool_; }
56 
57 private:
58  using CMMap = std::map<uint32_t, std::vector<float>>; //detId, Vector of MeanCM per detId
59  constexpr static uint16_t nTotStripsPerAPV = 128;
60 
61  uint16_t nullInspect(uint16_t firstAPV, const digivector_t& digis);
62  uint16_t abnormalBaselineInspect(uint16_t firstAPV, const digivector_t& digis);
63  uint16_t baselineFollowerInspect(uint16_t firstAPV, const digivector_t& digis);
64  uint16_t baselineAndSaturationInspect(uint16_t firstAPV, const digivector_t& digis);
65  uint16_t forceRestoreInspect(uint16_t firstAPV, const digivector_t& digis);
66  uint16_t hybridFormatInspect(uint16_t firstAPV, const digivector_t& digis);
67  uint16_t hybridEmulationInspect(uint16_t firstAPV, const digivector_t& digis);
68 
69  void flatRestore(uint16_t apvN, uint16_t firstAPV, digivector_t& digis);
70  bool checkBaseline(const std::vector<int16_t>& baseline) const;
71  void baselineFollowerRestore(uint16_t apvN, uint16_t firstAPV, float median, digivector_t& digis);
72  void derivativeFollowerRestore(uint16_t apvN, uint16_t firstAPV, digivector_t& digis);
73 
74  void baselineFollower(const digimap_t&, digivector_t& baseline, float median);
75  bool flatRegionsFinder(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
76 
77  void baselineCleaner(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
78  void cleaner_MonotonyChecker(digimap_t& smoothedpoints);
79  void cleaner_HighSlopeChecker(digimap_t& smoothedpoints);
80  void cleaner_LocalMinimumAdder(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
81 
84 
85 private: // members
88 
98 
99  // event state
101  // state
102  uint32_t detId_;
103  std::vector<std::string> apvFlags_;
104  std::vector<bool> apvFlagsBool_;
105  std::vector<bool> apvFlagsBoolOverride_;
106  std::vector<float> median_;
107  std::vector<bool> badAPVs_;
108  std::map<uint16_t, digimap_t> smoothedMaps_;
110 
111  //--------------------------------------------------
112  // Configurable Parameters of Algorithm
113  //--------------------------------------------------
118  int32_t meanCM_;
119  uint32_t deltaCMThreshold_; // for BaselineFollower inspect
120  double fraction_; // fraction of strips deviating from nominal baseline
121  uint32_t deviation_; // ADC value of deviation from nominal baseline
122  double restoreThreshold_; // for Null inspect (fraction of adc=0 channels)
123  uint32_t nSaturatedStrip_; // for BaselineAndSaturation inspect
124  uint32_t nSigmaNoiseDerTh_; // threshold for rejecting hits strips
125  uint32_t consecThreshold_; // minimum length of flat region
126  uint32_t nSmooth_; // for smoothing and local minimum determination (odd number)
127  uint32_t distortionThreshold_; // (max-min) of flat regions to trigger baseline follower
130  int32_t slopeX_;
131  int32_t slopeY_;
132  uint32_t hitStripThreshold_; // height above median when strip is definitely a hit
133  uint32_t minStripsToFit_; // minimum strips to try spline algo (otherwise default to median)
141 };
142 #endif
edm::ESWatcher< SiStripQualityRcd > qualityWatcher_
SiStripAPVRestorer(const edm::ParameterSet &conf, edm::ConsumesCollector)
const std::map< uint16_t, digimap_t > & getSmoothedPoints() const
void init(const edm::EventSetup &es)
void baselineFollowerRestore(uint16_t apvN, uint16_t firstAPV, float median, digivector_t &digis)
uint16_t baselineAndSaturationInspect(uint16_t firstAPV, const digivector_t &digis)
const SiStripPedestals * pedestalHandle
std::map< uint16_t, digimap_t > smoothedMaps_
uint16_t baselineFollowerInspect(uint16_t firstAPV, const digivector_t &digis)
double filteredBaselineDerivativeSumSquare_
void flatRestore(uint16_t apvN, uint16_t firstAPV, digivector_t &digis)
edm::ESGetToken< SiStripQuality, SiStripQualityRcd > qualityToken_
static std::string const input
Definition: EdmProvDump.cc:47
edm::ESWatcher< SiStripNoisesRcd > noiseWatcher_
void baselineFollower(const digimap_t &, digivector_t &baseline, float median)
void loadMeanCMMap(const edm::Event &)
edm::EDGetTokenT< edm::DetSetVector< SiStripRawDigi > > siStripRawDigiToken_
uint16_t hybridEmulationInspect(uint16_t firstAPV, const digivector_t &digis)
std::vector< std::string > apvFlags_
baselinemap_t baselineMap_
const SiStripQuality * qualityHandle
std::map< uint16_t, digivector_t > baselinemap_t
bool checkBaseline(const std::vector< int16_t > &baseline) const
void baselineCleaner(const digivector_t &adcs, digimap_t &smoothedpoints, uint16_t apvN)
void cleaner_LocalMinimumAdder(const digivector_t &adcs, digimap_t &smoothedpoints, uint16_t apvN)
void cleaner_HighSlopeChecker(digimap_t &smoothedpoints)
const SiStripNoises * noiseHandle
const baselinemap_t & getBaselineMap() const
std::map< uint32_t, std::vector< float > > CMMap
edm::EDGetTokenT< edm::DetSetVector< SiStripProcessedRawDigi > > siStripProcessedRawDigiToken_
bool flatRegionsFinder(const digivector_t &adcs, digimap_t &smoothedpoints, uint16_t apvN)
void cleaner_MonotonyChecker(digimap_t &smoothedpoints)
uint16_t inspectAndRestore(uint32_t detId, uint16_t firstAPV, const digivector_t &rawDigisPedSubtracted, digivector_t &processedRawDigi, const medians_t &vmedians)
void derivativeFollowerRestore(uint16_t apvN, uint16_t firstAPV, digivector_t &digis)
void createCMMapCMstored(const edm::DetSetVector< SiStripProcessedRawDigi > &input)
const std::vector< bool > & getAPVFlags() const
edm::ESGetToken< SiStripPedestals, SiStripPedestalsRcd > pedestalToken_
uint16_t hybridFormatInspect(uint16_t firstAPV, const digivector_t &digis)
edm::ESWatcher< SiStripPedestalsRcd > pedestalWatcher_
std::vector< bool > apvFlagsBoolOverride_
std::map< uint16_t, digi_t > digimap_t
edm::ESGetToken< SiStripNoises, SiStripNoisesRcd > noiseToken_
uint16_t nullInspect(uint16_t firstAPV, const digivector_t &digis)
std::vector< std::pair< short, float > > medians_t
std::vector< bool > badAPVs_
std::vector< bool > apvFlagsBool_
T median(std::vector< T > values)
Definition: median.h:16
uint16_t abnormalBaselineInspect(uint16_t firstAPV, const digivector_t &digis)
std::vector< digi_t > digivector_t
std::vector< float > median_
void restore(uint16_t firstAPV, digivector_t &digis)
uint16_t forceRestoreInspect(uint16_t firstAPV, const digivector_t &digis)
static constexpr uint16_t nTotStripsPerAPV
void createCMMapRealPed(const edm::DetSetVector< SiStripRawDigi > &input)
uint16_t inspect(uint32_t detId, uint16_t firstAPV, const digivector_t &digis, const medians_t &vmedians)