Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:33

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