Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:23:51

0001 #ifndef HcalSimAlgos_HcalTriggerPrimitiveAlgo_h
0002 #define HcalSimAlgos_HcalTriggerPrimitiveAlgo_h
0003 
0004 #include "DataFormats/HcalDetId/interface/HcalTrigTowerDetId.h"
0005 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0006 
0007 #include "CalibCalorimetry/HcalTPGAlgos/interface/HcaluLUTTPGCoder.h"
0008 #include "CalibFormats/CaloTPG/interface/HcalTPGCompressor.h"
0009 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
0010 #include "CondFormats/HcalObjects/interface/HcalElectronicsMap.h"
0011 
0012 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0013 
0014 #include "Geometry/HcalTowerAlgo/interface/HcalTrigTowerGeometry.h"
0015 
0016 #include "SimCalorimetry/HcalTrigPrimAlgos/interface/HcalFeatureHFEMBit.h"
0017 #include "SimCalorimetry/HcalTrigPrimAlgos/interface/HcalFinegrainBit.h"
0018 
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 
0021 #include <map>
0022 #include <vector>
0023 
0024 class CaloGeometry;
0025 class IntegerCaloSamples;
0026 
0027 class HcalTriggerPrimitiveAlgo {
0028 public:
0029   HcalTriggerPrimitiveAlgo(bool pf,
0030                            const std::vector<double>& w,
0031                            int latency,
0032                            uint32_t FG_threshold,
0033                            const std::vector<uint32_t>& FG_HF_thresholds,
0034                            uint32_t ZS_threshold,
0035                            int numberOfSamples,
0036                            int numberOfPresamples,
0037                            int numberOfFilterPresamplesHBQIE11,
0038                            int numberOfFilterPresamplesHEQIE11,
0039                            int numberOfSamplesHF,
0040                            int numberOfPresamplesHF,
0041                            bool useTDCInMinBiasBits,
0042                            uint32_t minSignalThreshold = 0,
0043                            uint32_t PMT_NoiseThreshold = 0);
0044   ~HcalTriggerPrimitiveAlgo();
0045 
0046   template <typename... Digis>
0047   void run(const HcalTPGCoder* incoder,
0048            const HcalTPGCompressor* outcoder,
0049            const HcalDbService* conditions,
0050            HcalTrigPrimDigiCollection& result,
0051            const HcalTrigTowerGeometry* trigTowerGeometry,
0052            float rctlsb,
0053            const HcalFeatureBit* LongvrsShortCut,
0054            const Digis&... digis);
0055 
0056   template <typename T, typename... Args>
0057   void addDigis(const T& collection, const Args&... digis) {
0058     addDigis(collection);
0059     addDigis(digis...);
0060   };
0061 
0062   template <typename T>
0063   void addDigis(const T& collection) {
0064     for (const auto& digi : collection) {
0065       addSignal(digi);
0066     }
0067   };
0068 
0069   template <typename D>
0070   void addDigis(const HcalDataFrameContainer<D>& collection) {
0071     for (auto i = collection.begin(); i != collection.end(); ++i) {
0072       D digi(*i);
0073       addSignal(digi);
0074     }
0075   };
0076 
0077   void runZS(HcalTrigPrimDigiCollection& tp);
0078   void runFEFormatError(const FEDRawDataCollection* rawraw,
0079                         const HcalElectronicsMap* emap,
0080                         HcalTrigPrimDigiCollection& result);
0081   void setPeakFinderAlgorithm(int algo);
0082   void setWeightsQIE11(const edm::ParameterSet& weightsQIE11);
0083   void setWeightQIE11(int aieta, int weight);
0084   void setNCTScaleShift(int);
0085   void setRCTScaleShift(int);
0086 
0087   void setNumFilterPresamplesHBQIE11(int presamples) { numberOfFilterPresamplesHBQIE11_ = presamples; }
0088 
0089   void setNumFilterPresamplesHEQIE11(int presamples) { numberOfFilterPresamplesHEQIE11_ = presamples; }
0090 
0091   void setUpgradeFlags(bool hb, bool he, bool hf);
0092   void setFixSaturationFlag(bool fix_saturation);
0093   void overrideParameters(const edm::ParameterSet& ps);
0094 
0095 private:
0096   /// adds the signal to the map
0097   void addSignal(const HBHEDataFrame& frame);
0098   void addSignal(const HFDataFrame& frame);
0099   void addSignal(const QIE10DataFrame& frame);
0100   void addSignal(const QIE11DataFrame& frame);
0101   void addSignal(const IntegerCaloSamples& samples);
0102   void addFG(const HcalTrigTowerDetId& id, std::vector<bool>& msb);
0103   void addUpgradeFG(const HcalTrigTowerDetId& id, int depth, const std::vector<std::bitset<2>>& bits);
0104   void addUpgradeTDCFG(const HcalTrigTowerDetId& id, const QIE11DataFrame& frame);
0105 
0106   bool passTDC(const QIE10DataFrame& digi, int ts) const;
0107   bool validUpgradeFG(const HcalTrigTowerDetId& id, int depth) const;
0108   bool validChannel(const QIE10DataFrame& digi, int ts) const;
0109   bool needLegacyFG(const HcalTrigTowerDetId& id) const;
0110   bool needUpgradeID(const HcalTrigTowerDetId& id, int depth) const;
0111 
0112   /// adds the actual digis
0113   void analyze(IntegerCaloSamples& samples, HcalTriggerPrimitiveDigi& result);
0114   // 2017 and later: QIE11
0115   void analyzeQIE11(IntegerCaloSamples& samples,
0116                     std::vector<bool> sample_saturation,
0117                     HcalTriggerPrimitiveDigi& result,
0118                     const HcalFinegrainBit& fg_algo);
0119   // Version 0: RCT
0120   void analyzeHF(IntegerCaloSamples& samples, HcalTriggerPrimitiveDigi& result, const int hf_lumi_shift);
0121   // Version 1: 1x1
0122   void analyzeHF2016(const IntegerCaloSamples& SAMPLES,
0123                      HcalTriggerPrimitiveDigi& result,
0124                      const int HF_LUMI_SHIFT,
0125                      const HcalFeatureBit* HCALFEM);
0126   // With dual anode readout
0127   void analyzeHFQIE10(const IntegerCaloSamples& SAMPLES,
0128                       HcalTriggerPrimitiveDigi& result,
0129                       const int HF_LUMI_SHIFT,
0130                       const HcalFeatureBit* HCALFEM);
0131 
0132   // Member initialized by constructor
0133   const HcaluLUTTPGCoder* incoder_;
0134   const HcalTPGCompressor* outcoder_;
0135   const HcalDbService* conditions_;
0136   double theThreshold;
0137   bool peakfind_;
0138   std::vector<double> weights_;
0139   std::array<std::array<int, 2>, 29> weightsQIE11_;
0140   int latency_;
0141   uint32_t FG_threshold_;
0142   std::vector<uint32_t> FG_HF_thresholds_;
0143   uint32_t ZS_threshold_;
0144   int ZS_threshold_I_;
0145   int numberOfSamples_;
0146   int numberOfPresamples_;
0147   int numberOfFilterPresamplesHBQIE11_;
0148   int numberOfFilterPresamplesHEQIE11_;
0149   int numberOfSamplesHF_;
0150   int numberOfPresamplesHF_;
0151   bool useTDCInMinBiasBits_;
0152   uint32_t minSignalThreshold_;
0153   uint32_t PMT_NoiseThreshold_;
0154   int NCTScaleShift;
0155   int RCTScaleShift;
0156 
0157   // Algo1: isPeak = TS[i-1] < TS[i] && TS[i] >= TS[i+1]
0158   // Algo2: isPeak = TSS[i-1] < TSS[i] && TSS[i] >= TSS[i+1],
0159   // TSS[i] = TS[i] + TS[i+1]
0160   // Default: Algo2
0161   int peak_finder_algorithm_;
0162 
0163   // Member not initialzed
0164   //std::vector<HcalTrigTowerDetId> towerIds(const HcalDetId & id) const;
0165 
0166   const HcalTrigTowerGeometry* theTrigTowerGeometry;
0167 
0168   typedef std::map<HcalTrigTowerDetId, IntegerCaloSamples> SumMap;
0169   SumMap theSumMap;
0170 
0171   typedef std::map<HcalTrigTowerDetId, std::vector<bool>> SatMap;
0172   SatMap theSatMap;
0173 
0174   struct HFDetails {
0175     IntegerCaloSamples long_fiber;
0176     IntegerCaloSamples short_fiber;
0177     HFDataFrame ShortDigi;
0178     HFDataFrame LongDigi;
0179   };
0180   typedef std::map<HcalTrigTowerDetId, std::map<uint32_t, HFDetails>> HFDetailMap;
0181   HFDetailMap theHFDetailMap;
0182 
0183   struct HFUpgradeDetails {
0184     IntegerCaloSamples samples;
0185     QIE10DataFrame digi;
0186     std::vector<bool> validity;
0187     std::vector<std::bitset<2>> fgbits;
0188     std::vector<bool> passTDC;
0189   };
0190   typedef std::map<HcalTrigTowerDetId, std::map<uint32_t, std::array<HFUpgradeDetails, 4>>> HFUpgradeDetailMap;
0191   HFUpgradeDetailMap theHFUpgradeDetailMap;
0192 
0193   typedef std::vector<IntegerCaloSamples> SumFGContainer;
0194   typedef std::map<HcalTrigTowerDetId, SumFGContainer> TowerMapFGSum;
0195   TowerMapFGSum theTowerMapFGSum;
0196 
0197   // ==============================
0198   // =  HF Veto
0199   // ==============================
0200   // Sum = Long + Short;" // intermediate calculation.
0201   //  if ((Short < MinSignalThresholdET OR Long  < MinSignalThresholdET)
0202   //     AND Sum > PMTNoiseThresholdET) VetoedSum = 0;
0203   //  else VetoedSum = Sum;
0204   // ==============================
0205   // Map from FG id to veto booleans
0206   HcalFeatureBit* LongvrsShortCut;
0207   typedef std::map<uint32_t, std::vector<bool>> TowerMapVeto;
0208   TowerMapVeto HF_Veto;
0209 
0210   typedef std::map<HcalTrigTowerDetId, std::vector<bool>> FGbitMap;
0211   FGbitMap fgMap_;
0212 
0213   typedef std::vector<HcalFinegrainBit::Tower> FGUpgradeContainer;
0214   typedef std::map<HcalTrigTowerDetId, FGUpgradeContainer> FGUpgradeMap;
0215   FGUpgradeMap fgUpgradeMap_;
0216 
0217   typedef std::vector<HcalFinegrainBit::TowerTDC> FGUpgradeTDCContainer;
0218   typedef std::map<HcalTrigTowerDetId, FGUpgradeTDCContainer> FGUpgradeTDCMap;
0219   FGUpgradeTDCMap fgUpgradeTDCMap_;
0220 
0221   bool upgrade_hb_ = false;
0222   bool upgrade_he_ = false;
0223   bool upgrade_hf_ = false;
0224 
0225   bool fix_saturation_ = false;
0226 
0227   edm::ParameterSet override_parameters_;
0228 
0229   bool override_adc_hf_ = false;
0230   uint32_t override_adc_hf_value_;
0231   bool override_tdc_hf_ = false;
0232   unsigned long long override_tdc_hf_value_;
0233 
0234   // HE constants
0235   static const int HBHE_OVERLAP_TOWER = 16;
0236   static const int FIRST_DEPTH7_TOWER = 26;
0237   static const int LAST_FINEGRAIN_DEPTH = 6;
0238   static const int LAST_FINEGRAIN_TOWER = 28;
0239 
0240   // Fine-grain in HF ignores tower 29, and starts with 30
0241   static const int FIRST_FINEGRAIN_TOWER = 30;
0242 
0243   static const int QIE8_LINEARIZATION_ET = HcaluLUTTPGCoder::QIE8_LUT_BITMASK;
0244   static const int QIE10_LINEARIZATION_ET = HcaluLUTTPGCoder::QIE10_LUT_BITMASK;
0245   static const int QIE11_LINEARIZATION_ET = HcaluLUTTPGCoder::QIE11_LUT_BITMASK;
0246   // Consider CaloTPGTranscoderULUT.h for values
0247   static const int QIE10_MAX_LINEARIZATION_ET = 0x7FF;
0248   static const int QIE11_MAX_LINEARIZATION_ET = 0x7FF;
0249 };
0250 
0251 template <typename... Digis>
0252 void HcalTriggerPrimitiveAlgo::run(const HcalTPGCoder* incoder,
0253                                    const HcalTPGCompressor* outcoder,
0254                                    const HcalDbService* conditions,
0255                                    HcalTrigPrimDigiCollection& result,
0256                                    const HcalTrigTowerGeometry* trigTowerGeometry,
0257                                    float rctlsb,
0258                                    const HcalFeatureBit* LongvrsShortCut,
0259                                    const Digis&... digis) {
0260   theTrigTowerGeometry = trigTowerGeometry;
0261 
0262   incoder_ = dynamic_cast<const HcaluLUTTPGCoder*>(incoder);
0263   outcoder_ = outcoder;
0264   conditions_ = conditions;
0265 
0266   theSumMap.clear();
0267   theSatMap.clear();
0268   theTowerMapFGSum.clear();
0269   HF_Veto.clear();
0270   fgMap_.clear();
0271   fgUpgradeMap_.clear();
0272   fgUpgradeTDCMap_.clear();
0273   theHFDetailMap.clear();
0274   theHFUpgradeDetailMap.clear();
0275 
0276   // Add all digi collections
0277   addDigis(digis...);
0278 
0279   // Prepare the fine-grain calculation algorithm for HB/HE
0280   int version = 0;
0281   if (upgrade_he_ or upgrade_hb_)
0282     version = conditions_->getHcalTPParameters()->getFGVersionHBHE();
0283   if (override_parameters_.exists("FGVersionHBHE"))
0284     version = override_parameters_.getParameter<uint32_t>("FGVersionHBHE");
0285   HcalFinegrainBit fg_algo(version);
0286 
0287   // VME produces additional bits on the front used by lumi but not the
0288   // trigger, this shift corrects those out by right shifting over them.
0289   for (auto& item : theSumMap) {
0290     result.push_back(HcalTriggerPrimitiveDigi(item.first));
0291     HcalTrigTowerDetId detId(item.second.id());
0292     if (detId.ietaAbs() >= theTrigTowerGeometry->firstHFTower(detId.version())) {
0293       if (detId.version() == 0) {
0294         analyzeHF(item.second, result.back(), RCTScaleShift);
0295       } else if (detId.version() == 1) {
0296         if (upgrade_hf_)
0297           analyzeHFQIE10(item.second, result.back(), NCTScaleShift, LongvrsShortCut);
0298         else
0299           analyzeHF2016(item.second, result.back(), NCTScaleShift, LongvrsShortCut);
0300       } else {
0301         // Things are going to go poorly
0302       }
0303     } else {
0304       // Determine which energy reconstruction path to take based on the
0305       // fine-grain availability:
0306       //  * QIE8 TP add entries into fgMap_
0307       //  * QIE11 TP add entries into fgUpgradeMap_
0308       //    (not for tower 16 unless HB is upgraded, too)
0309       if (fgMap_.find(item.first) != fgMap_.end()) {
0310         analyze(item.second, result.back());
0311       } else if (fgUpgradeMap_.find(item.first) != fgUpgradeMap_.end()) {
0312         SatMap::iterator item_sat = theSatMap.find(detId);
0313         if (item_sat == theSatMap.end())
0314           analyzeQIE11(item.second, std::vector<bool>(), result.back(), fg_algo);
0315         else
0316           analyzeQIE11(item.second, item_sat->second, result.back(), fg_algo);
0317       }
0318     }
0319   }
0320 
0321   // Free up some memory
0322   theSumMap.clear();
0323   theSatMap.clear();
0324   theTowerMapFGSum.clear();
0325   HF_Veto.clear();
0326   fgMap_.clear();
0327   fgUpgradeMap_.clear();
0328   fgUpgradeTDCMap_.clear();
0329   theHFDetailMap.clear();
0330   theHFUpgradeDetailMap.clear();
0331 
0332   return;
0333 }
0334 
0335 #endif