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 
25 namespace edm {
26  class ParameterSet;
28 } // namespace edm
29 
32 
33 protected:
35 
36 public:
37  virtual ~SiStripAPVRestorer() {}
38 
39  using digi_t = int16_t;
40  using digivector_t = std::vector<digi_t>;
41  using digimap_t = std::map<uint16_t, digi_t>;
42  using medians_t = std::vector<std::pair<short, float>>;
43  using baselinemap_t = std::map<uint16_t, digivector_t>;
44 
45  void init(const edm::EventSetup& es);
46 
47  uint16_t inspect(uint32_t detId, uint16_t firstAPV, const digivector_t& digis, const medians_t& vmedians);
48  void restore(uint16_t firstAPV, digivector_t& digis);
49 
50  uint16_t inspectAndRestore(uint32_t detId,
51  uint16_t firstAPV,
52  const digivector_t& rawDigisPedSubtracted,
53  digivector_t& processedRawDigi,
54  const medians_t& vmedians);
55 
56  void loadMeanCMMap(const edm::Event&);
57 
58  const baselinemap_t& getBaselineMap() const { return baselineMap_; }
59  const std::map<uint16_t, digimap_t>& getSmoothedPoints() const { return smoothedMaps_; }
60  const std::vector<bool>& getAPVFlags() const { return apvFlagsBool_; }
61 
62 private:
63  using CMMap = std::map<uint32_t, std::vector<float>>; //detId, Vector of MeanCM per detId
64  constexpr static uint16_t nTotStripsPerAPV = 128;
65 
66  uint16_t nullInspect(uint16_t firstAPV, const digivector_t& digis);
67  uint16_t abnormalBaselineInspect(uint16_t firstAPV, const digivector_t& digis);
68  uint16_t baselineFollowerInspect(uint16_t firstAPV, const digivector_t& digis);
69  uint16_t baselineAndSaturationInspect(uint16_t firstAPV, const digivector_t& digis);
70  uint16_t forceRestoreInspect(uint16_t firstAPV, const digivector_t& digis);
71  uint16_t hybridFormatInspect(uint16_t firstAPV, const digivector_t& digis);
72  uint16_t hybridEmulationInspect(uint16_t firstAPV, const digivector_t& digis);
73 
74  void flatRestore(uint16_t apvN, uint16_t firstAPV, digivector_t& digis);
75  bool checkBaseline(const std::vector<int16_t>& baseline) const;
76  void baselineFollowerRestore(uint16_t apvN, uint16_t firstAPV, float median, digivector_t& digis);
77  void derivativeFollowerRestore(uint16_t apvN, uint16_t firstAPV, digivector_t& digis);
78 
79  void baselineFollower(const digimap_t&, digivector_t& baseline, float median);
80  bool flatRegionsFinder(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
81 
82  void baselineCleaner(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
83  void cleaner_MonotonyChecker(digimap_t& smoothedpoints);
84  void cleaner_HighSlopeChecker(digimap_t& smoothedpoints);
85  void cleaner_LocalMinimumAdder(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
86 
89 
91 
92 private: // members
95 
105 
106  // event state
108  // state
109  uint32_t detId_;
110  std::vector<std::string> apvFlags_;
111  std::vector<bool> apvFlagsBool_;
112  std::vector<bool> apvFlagsBoolOverride_;
113  std::vector<float> median_;
114  std::vector<bool> badAPVs_;
115  std::map<uint16_t, digimap_t> smoothedMaps_;
117 
118  //--------------------------------------------------
119  // Configurable Parameters of Algorithm
120  //--------------------------------------------------
125  int32_t meanCM_;
126  uint32_t deltaCMThreshold_; // for BaselineFollower inspect
127  double fraction_; // fraction of strips deviating from nominal baseline
128  uint32_t deviation_; // ADC value of deviation from nominal baseline
129  double restoreThreshold_; // for Null inspect (fraction of adc=0 channels)
130  uint32_t nSaturatedStrip_; // for BaselineAndSaturation inspect
131  uint32_t nSigmaNoiseDerTh_; // threshold for rejecting hits strips
132  uint32_t consecThreshold_; // minimum length of flat region
133  uint32_t nSmooth_; // for smoothing and local minimum determination (odd number)
134  uint32_t distortionThreshold_; // (max-min) of flat regions to trigger baseline follower
137  int32_t slopeX_;
138  int32_t slopeY_;
139  uint32_t hitStripThreshold_; // height above median when strip is definitely a hit
140  uint32_t minStripsToFit_; // minimum strips to try spline algo (otherwise default to median)
148 };
149 #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)
static void fillDescriptions(edm::ParameterSetDescription &desc)
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:50
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
HLT enums.
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)