Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:37:44

0001 #ifndef RECOLOCALTRACKER_SISTRIPZEROSUPPRESSION_SISTRIPAPVRESTORER_H
0002 #define RECOLOCALTRACKER_SISTRIPZEROSUPPRESSION_SISTRIPAPVRESTORER_H
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "DataFormats/Common/interface/Handle.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/ESWatcher.h"
0009 #include "CondFormats/SiStripObjects/interface/SiStripNoises.h"
0010 #include "CondFormats/SiStripObjects/interface/SiStripPedestals.h"
0011 #include "CalibFormats/SiStripObjects/interface/SiStripQuality.h"
0012 #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
0013 #include "CondFormats/DataRecord/interface/SiStripPedestalsRcd.h"
0014 #include "CalibTracker/Records/interface/SiStripQualityRcd.h"
0015 #include "DataFormats/Common/interface/DetSet.h"
0016 #include "DataFormats/Common/interface/DetSetVector.h"
0017 #include "DataFormats/DetId/interface/DetId.h"
0018 #include "DataFormats/SiStripDigi/interface/SiStripRawDigi.h"
0019 #include "DataFormats/SiStripDigi/interface/SiStripProcessedRawDigi.h"
0020 #include "RecoLocalTracker/SiStripZeroSuppression/interface/SiStripCommonModeNoiseSubtractor.h"
0021 
0022 #include <vector>
0023 #include <cstdint>
0024 
0025 namespace edm {
0026   class ParameterSet;
0027   class ParameterSetDescription;
0028 }  // namespace edm
0029 
0030 class SiStripAPVRestorer {
0031   friend class SiStripRawProcessingFactory;
0032 
0033 protected:
0034   SiStripAPVRestorer(const edm::ParameterSet& conf, edm::ConsumesCollector);
0035 
0036 public:
0037   virtual ~SiStripAPVRestorer() {}
0038 
0039   using digi_t = int16_t;
0040   using digivector_t = std::vector<digi_t>;
0041   using digimap_t = std::map<uint16_t, digi_t>;
0042   using medians_t = std::vector<std::pair<short, float>>;
0043   using baselinemap_t = std::map<uint16_t, digivector_t>;
0044 
0045   void init(const edm::EventSetup& es);
0046 
0047   uint16_t inspect(uint32_t detId, uint16_t firstAPV, const digivector_t& digis, const medians_t& vmedians);
0048   void restore(uint16_t firstAPV, digivector_t& digis);
0049 
0050   uint16_t inspectAndRestore(uint32_t detId,
0051                              uint16_t firstAPV,
0052                              const digivector_t& rawDigisPedSubtracted,
0053                              digivector_t& processedRawDigi,
0054                              const medians_t& vmedians);
0055 
0056   void loadMeanCMMap(const edm::Event&);
0057 
0058   const baselinemap_t& getBaselineMap() const { return baselineMap_; }
0059   const std::map<uint16_t, digimap_t>& getSmoothedPoints() const { return smoothedMaps_; }
0060   const std::vector<bool>& getAPVFlags() const { return apvFlagsBool_; }
0061 
0062 private:
0063   using CMMap = std::map<uint32_t, std::vector<float>>;  //detId, Vector of MeanCM per detId
0064   constexpr static uint16_t nTotStripsPerAPV = 128;
0065 
0066   uint16_t nullInspect(uint16_t firstAPV, const digivector_t& digis);
0067   uint16_t abnormalBaselineInspect(uint16_t firstAPV, const digivector_t& digis);
0068   uint16_t baselineFollowerInspect(uint16_t firstAPV, const digivector_t& digis);
0069   uint16_t baselineAndSaturationInspect(uint16_t firstAPV, const digivector_t& digis);
0070   uint16_t forceRestoreInspect(uint16_t firstAPV, const digivector_t& digis);
0071   uint16_t hybridFormatInspect(uint16_t firstAPV, const digivector_t& digis);
0072   uint16_t hybridEmulationInspect(uint16_t firstAPV, const digivector_t& digis);
0073 
0074   void flatRestore(uint16_t apvN, uint16_t firstAPV, digivector_t& digis);
0075   bool checkBaseline(const std::vector<int16_t>& baseline) const;
0076   void baselineFollowerRestore(uint16_t apvN, uint16_t firstAPV, float median, digivector_t& digis);
0077   void derivativeFollowerRestore(uint16_t apvN, uint16_t firstAPV, digivector_t& digis);
0078 
0079   void baselineFollower(const digimap_t&, digivector_t& baseline, float median);
0080   bool flatRegionsFinder(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
0081 
0082   void baselineCleaner(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
0083   void cleaner_MonotonyChecker(digimap_t& smoothedpoints);
0084   void cleaner_HighSlopeChecker(digimap_t& smoothedpoints);
0085   void cleaner_LocalMinimumAdder(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN);
0086 
0087   void createCMMapRealPed(const edm::DetSetVector<SiStripRawDigi>& input);
0088   void createCMMapCMstored(const edm::DetSetVector<SiStripProcessedRawDigi>& input);
0089 
0090   static void fillDescriptions(edm::ParameterSetDescription& desc);
0091 
0092 private:  // members
0093   edm::EDGetTokenT<edm::DetSetVector<SiStripRawDigi>> siStripRawDigiToken_;
0094   edm::EDGetTokenT<edm::DetSetVector<SiStripProcessedRawDigi>> siStripProcessedRawDigiToken_;
0095 
0096   edm::ESGetToken<SiStripQuality, SiStripQualityRcd> qualityToken_;
0097   edm::ESGetToken<SiStripNoises, SiStripNoisesRcd> noiseToken_;
0098   edm::ESGetToken<SiStripPedestals, SiStripPedestalsRcd> pedestalToken_;
0099   const SiStripQuality* qualityHandle;
0100   const SiStripNoises* noiseHandle;
0101   const SiStripPedestals* pedestalHandle;
0102   edm::ESWatcher<SiStripQualityRcd> qualityWatcher_;
0103   edm::ESWatcher<SiStripNoisesRcd> noiseWatcher_;
0104   edm::ESWatcher<SiStripPedestalsRcd> pedestalWatcher_;
0105 
0106   // event state
0107   CMMap meanCMmap_;
0108   // state
0109   uint32_t detId_;
0110   std::vector<std::string> apvFlags_;
0111   std::vector<bool> apvFlagsBool_;
0112   std::vector<bool> apvFlagsBoolOverride_;
0113   std::vector<float> median_;
0114   std::vector<bool> badAPVs_;
0115   std::map<uint16_t, digimap_t> smoothedMaps_;
0116   baselinemap_t baselineMap_;
0117 
0118   //--------------------------------------------------
0119   // Configurable Parameters of Algorithm
0120   //--------------------------------------------------
0121   bool forceNoRestore_;
0122   std::string inspectAlgo_;
0123   std::string restoreAlgo_;
0124   bool useRealMeanCM_;
0125   int32_t meanCM_;
0126   uint32_t deltaCMThreshold_;     // for BaselineFollower inspect
0127   double fraction_;               // fraction of strips deviating from nominal baseline
0128   uint32_t deviation_;            // ADC value of deviation from nominal baseline
0129   double restoreThreshold_;       // for Null inspect (fraction of adc=0 channels)
0130   uint32_t nSaturatedStrip_;      // for BaselineAndSaturation inspect
0131   uint32_t nSigmaNoiseDerTh_;     // threshold for rejecting hits strips
0132   uint32_t consecThreshold_;      // minimum length of flat region
0133   uint32_t nSmooth_;              // for smoothing and local minimum determination (odd number)
0134   uint32_t distortionThreshold_;  // (max-min) of flat regions to trigger baseline follower
0135   bool applyBaselineCleaner_;
0136   uint32_t cleaningSequence_;
0137   int32_t slopeX_;
0138   int32_t slopeY_;
0139   uint32_t hitStripThreshold_;  // height above median when strip is definitely a hit
0140   uint32_t minStripsToFit_;     // minimum strips to try spline algo (otherwise default to median)
0141   bool applyBaselineRejection_;
0142   double filteredBaselineMax_;
0143   double filteredBaselineDerivativeSumSquare_;
0144   int gradient_threshold_;
0145   int last_gradient_;
0146   int size_window_;
0147   int width_cluster_;
0148 };
0149 #endif