Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:31

0001 #include "RecoLocalTracker/SiStripZeroSuppression/plugins/SiStripZeroSuppression.h"
0002 
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "DataFormats/Common/interface/DetSetVector.h"
0008 #include "DataFormats/SiStripDigi/interface/SiStripDigi.h"
0009 #include "DataFormats/SiStripDigi/interface/SiStripRawDigi.h"
0010 #include "RecoLocalTracker/SiStripZeroSuppression/interface/SiStripRawProcessingFactory.h"
0011 #include "FWCore/Utilities/interface/transform.h"
0012 #include <memory>
0013 
0014 SiStripZeroSuppression::SiStripZeroSuppression(edm::ParameterSet const& conf)
0015     : algorithms(
0016           SiStripRawProcessingFactory::create(conf.getParameter<edm::ParameterSet>("Algorithms"), consumesCollector())),
0017       produceRawDigis(conf.getParameter<bool>("produceRawDigis")),
0018       storeCM(conf.getParameter<bool>("storeCM")),
0019       fixCM(conf.getParameter<bool>("fixCM")),
0020       produceCalculatedBaseline(conf.getParameter<bool>("produceCalculatedBaseline")),
0021       produceBaselinePoints(conf.getParameter<bool>("produceBaselinePoints")),
0022       storeInZScollBadAPV(conf.getParameter<bool>("storeInZScollBadAPV")),
0023       produceHybridFormat(conf.getParameter<bool>("produceHybridFormat")) {
0024   for (const auto& inputTag : conf.getParameter<std::vector<edm::InputTag>>("RawDigiProducersList")) {
0025     const auto& tagName = inputTag.instance();
0026     produces<edm::DetSetVector<SiStripDigi>>(tagName);
0027     if (produceRawDigis)
0028       produces<edm::DetSetVector<SiStripRawDigi>>(tagName);
0029     if (storeCM)
0030       produces<edm::DetSetVector<SiStripProcessedRawDigi>>("APVCM" + tagName);
0031     if (produceCalculatedBaseline)
0032       produces<edm::DetSetVector<SiStripProcessedRawDigi>>("BADAPVBASELINE" + tagName);
0033     if (produceBaselinePoints)
0034       produces<edm::DetSetVector<SiStripDigi>>("BADAPVBASELINEPOINTS" + tagName);
0035 
0036     RawType inputType = RawType::Unknown;
0037     if (tagName == "ProcessedRaw") {
0038       inputType = RawType::ProcessedRaw;
0039       if (produceHybridFormat)
0040         throw cms::Exception("Processed Raw Cannot be converted in hybrid Format");
0041     } else if (tagName == "VirginRaw") {
0042       inputType = RawType::VirginRaw;
0043     } else if (tagName == "ScopeMode") {
0044       inputType = RawType::ScopeMode;
0045       if (produceHybridFormat)
0046         throw cms::Exception("Scope Mode cannot be converted in hybrid Format");
0047     }
0048     if (RawType::Unknown != inputType) {
0049       rawInputs.emplace_back(tagName, inputType, consumes<edm::DetSetVector<SiStripRawDigi>>(inputTag));
0050     } else if (tagName == "ZeroSuppressed") {
0051       hybridInputs.emplace_back(tagName, consumes<edm::DetSetVector<SiStripDigi>>(inputTag));
0052     } else {
0053       throw cms::Exception("Unknown input type")
0054           << tagName << " unknown.  "
0055           << "SiStripZeroZuppression can only process types \"VirginRaw\", \"ProcessedRaw\" and \"ZeroSuppressed\"";
0056     }
0057   }
0058 
0059   if (produceHybridFormat &&
0060       ("HybridEmulation" !=
0061        conf.getParameter<edm::ParameterSet>("Algorithms").getParameter<std::string>("APVInspectMode")))
0062     throw cms::Exception("Invalid option") << "When producing data in the hybrid format, the APV restorer must be "
0063                                               "configured with APVInspectMode='HybridEmulation'";
0064 
0065   if ((!hybridInputs.empty()) && produceRawDigis) {
0066     edm::LogInfo("SiStripZeroSuppression") << "Raw digis will not be saved for hybrid inputs";
0067   }
0068 
0069   if (!(rawInputs.empty() && hybridInputs.empty())) {
0070     output_base.reserve(16000);
0071     if (produceRawDigis && !rawInputs.empty())
0072       output_base_raw.reserve(16000);
0073     if (storeCM)
0074       output_apvcm.reserve(16000);
0075     if (produceCalculatedBaseline)
0076       output_baseline.reserve(16000);
0077     if (produceBaselinePoints)
0078       output_baseline_points.reserve(16000);
0079   }
0080 }
0081 
0082 void SiStripZeroSuppression::produce(edm::Event& e, const edm::EventSetup& es) {
0083   algorithms->initialize(es, e);
0084 
0085   for (const auto& input : rawInputs) {
0086     clearOutputs();
0087     edm::Handle<edm::DetSetVector<SiStripRawDigi>> inDigis;
0088     e.getByToken(std::get<rawtoken_t>(input), inDigis);
0089     if (!inDigis->empty())
0090       processRaw(*inDigis, std::get<RawType>(input));
0091     putOutputs(e, std::get<std::string>(input));
0092   }
0093   for (const auto& input : hybridInputs) {
0094     clearOutputs();
0095     edm::Handle<edm::DetSetVector<SiStripDigi>> inDigis;
0096     e.getByToken(std::get<zstoken_t>(input), inDigis);
0097     if (!inDigis->empty()) {
0098       processHybrid(*inDigis);
0099     }
0100     putOutputs(e, std::get<std::string>(input));
0101   }
0102 }
0103 
0104 inline void SiStripZeroSuppression::clearOutputs() {
0105   output_base.clear();
0106   output_base_raw.clear();
0107   output_baseline.clear();
0108   output_baseline_points.clear();
0109   output_apvcm.clear();
0110 }
0111 inline void SiStripZeroSuppression::putOutputs(edm::Event& evt, const std::string& tagName) {
0112   evt.put(std::make_unique<edm::DetSetVector<SiStripDigi>>(output_base), tagName);
0113   if (produceRawDigis && !rawInputs.empty())
0114     evt.put(std::make_unique<edm::DetSetVector<SiStripRawDigi>>(output_base_raw), tagName);
0115   if (produceCalculatedBaseline)
0116     evt.put(std::make_unique<edm::DetSetVector<SiStripProcessedRawDigi>>(output_baseline), "BADAPVBASELINE" + tagName);
0117   if (produceBaselinePoints)
0118     evt.put(std::make_unique<edm::DetSetVector<SiStripDigi>>(output_baseline_points), "BADAPVBASELINEPOINTS" + tagName);
0119   if (storeCM)
0120     evt.put(std::make_unique<edm::DetSetVector<SiStripProcessedRawDigi>>(output_apvcm), "APVCM" + tagName);
0121 }
0122 
0123 inline void SiStripZeroSuppression::processRaw(const edm::DetSetVector<SiStripRawDigi>& input, RawType inType) {
0124   for (const auto& rawDigis : input) {
0125     edm::DetSet<SiStripDigi> suppressedDigis(rawDigis.id);
0126 
0127     int16_t nAPVflagged = 0;
0128     if (RawType::ProcessedRaw == inType) {
0129       nAPVflagged = algorithms->suppressProcessedRawData(rawDigis, suppressedDigis);
0130     } else if (RawType::ScopeMode == inType) {
0131       nAPVflagged = algorithms->suppressVirginRawData(rawDigis, suppressedDigis);
0132     } else if (RawType::VirginRaw == inType) {
0133       if (produceHybridFormat) {
0134         nAPVflagged = algorithms->convertVirginRawToHybrid(rawDigis, suppressedDigis);
0135       } else {
0136         nAPVflagged = algorithms->suppressVirginRawData(rawDigis, suppressedDigis);
0137       }
0138     }
0139 
0140     storeExtraOutput(rawDigis.id, nAPVflagged);
0141     if (!suppressedDigis.empty() && (storeInZScollBadAPV || nAPVflagged == 0))
0142       output_base.push_back(std::move(suppressedDigis));
0143 
0144     if (produceRawDigis && nAPVflagged > 0) {
0145       output_base_raw.push_back(formatRawDigis(rawDigis));
0146     }
0147   }
0148 }
0149 
0150 inline void SiStripZeroSuppression::processHybrid(const edm::DetSetVector<SiStripDigi>& input) {
0151   for (const auto& inDigis : input) {
0152     edm::DetSet<SiStripDigi> suppressedDigis(inDigis.id);
0153 
0154     uint16_t nAPVflagged = 0;
0155     nAPVflagged = algorithms->suppressHybridData(inDigis, suppressedDigis);
0156 
0157     storeExtraOutput(inDigis.id, nAPVflagged);
0158     if (!suppressedDigis.empty())
0159       output_base.push_back(std::move(suppressedDigis));
0160   }
0161 }
0162 
0163 inline edm::DetSet<SiStripRawDigi> SiStripZeroSuppression::formatRawDigis(const edm::DetSet<SiStripRawDigi>& rawDigis) {
0164   edm::DetSet<SiStripRawDigi> outRawDigis(rawDigis.id);
0165   outRawDigis.reserve(rawDigis.size());
0166   const std::vector<bool>& apvf = algorithms->getAPVFlags();
0167   uint32_t strip = 0;
0168   for (const auto rawDigi : rawDigis) {
0169     int16_t apvN = strip / 128;
0170     if (apvf[apvN])
0171       outRawDigis.push_back(rawDigi);
0172     else
0173       outRawDigis.push_back(SiStripRawDigi(0));
0174     ++strip;
0175   }
0176   return outRawDigis;
0177 }
0178 
0179 inline edm::DetSet<SiStripRawDigi> SiStripZeroSuppression::formatRawDigis(uint32_t detId,
0180                                                                           const std::vector<int16_t>& rawDigis) {
0181   edm::DetSet<SiStripRawDigi> outRawDigis(detId);
0182   outRawDigis.reserve(rawDigis.size());
0183   const std::vector<bool>& apvf = algorithms->getAPVFlags();
0184   uint32_t strip = 0;
0185   for (const auto rawDigi : rawDigis) {
0186     int16_t apvN = strip / 128;
0187     if (apvf[apvN])
0188       outRawDigis.push_back(SiStripRawDigi(rawDigi));
0189     else
0190       outRawDigis.push_back(SiStripRawDigi(0));
0191     ++strip;
0192   }
0193   return outRawDigis;
0194 }
0195 
0196 inline void SiStripZeroSuppression::storeExtraOutput(uint32_t id, int16_t nAPVflagged) {
0197   const auto& vmedians = algorithms->getAPVsCM();
0198   if (storeCM)
0199     storeCMN(id, vmedians);
0200   if (nAPVflagged > 0) {
0201     if (produceCalculatedBaseline)
0202       storeBaseline(id, vmedians);
0203     if (produceBaselinePoints)
0204       storeBaselinePoints(id);
0205   }
0206 }
0207 
0208 inline void SiStripZeroSuppression::storeBaseline(uint32_t id, const medians_t& vmedians) {
0209   const auto& baselinemap = algorithms->getBaselineMap();
0210 
0211   edm::DetSet<SiStripProcessedRawDigi> baselineDetSet(id);
0212   baselineDetSet.reserve(vmedians.size() * 128);
0213   for (const auto& vmed : vmedians) {
0214     const uint16_t apvN = vmed.first;
0215     const float median = vmed.second;
0216     auto itBaselineMap = baselinemap.find(apvN);
0217     if (baselinemap.end() == itBaselineMap) {
0218       for (size_t strip = 0; strip < 128; ++strip)
0219         baselineDetSet.push_back(SiStripProcessedRawDigi(median));
0220     } else {
0221       for (size_t strip = 0; strip < 128; ++strip)
0222         baselineDetSet.push_back(SiStripProcessedRawDigi((itBaselineMap->second)[strip]));
0223     }
0224   }
0225 
0226   if (!baselineDetSet.empty())
0227     output_baseline.push_back(baselineDetSet);
0228 }
0229 
0230 inline void SiStripZeroSuppression::storeBaselinePoints(uint32_t id) {
0231   edm::DetSet<SiStripDigi> baspointDetSet(id);
0232   for (const auto& itBaselinePointVect : algorithms->getSmoothedPoints()) {
0233     const uint16_t apvN = itBaselinePointVect.first;
0234     for (const auto& itBaselinePointMap : itBaselinePointVect.second) {
0235       const uint16_t bpstrip = itBaselinePointMap.first + apvN * 128;
0236       const int16_t bp = itBaselinePointMap.second;
0237       baspointDetSet.push_back(SiStripDigi(bpstrip, bp + 128));
0238     }
0239   }
0240 
0241   if (!baspointDetSet.empty())
0242     output_baseline_points.push_back(std::move(baspointDetSet));
0243 }
0244 
0245 inline void SiStripZeroSuppression::storeCMN(uint32_t id, const medians_t& vmedians) {
0246   std::vector<bool> apvf(6, false);
0247   if (fixCM) {
0248     const auto& apvFlagged = algorithms->getAPVFlags();
0249     std::copy(std::begin(apvFlagged), std::end(apvFlagged), std::begin(apvf));
0250   }
0251 
0252   edm::DetSet<SiStripProcessedRawDigi> apvDetSet(id);
0253   short apvNb = 0;
0254   for (const auto& vmed : vmedians) {
0255     if (vmed.first > apvNb) {
0256       for (int i{0}; i < vmed.first - apvNb; ++i) {
0257         apvDetSet.push_back(SiStripProcessedRawDigi(-999.));
0258         apvNb++;
0259       }
0260     }
0261 
0262     if (fixCM && apvf[vmed.first]) {
0263       apvDetSet.push_back(SiStripProcessedRawDigi(-999.));
0264     } else {
0265       apvDetSet.push_back(SiStripProcessedRawDigi(vmed.second));
0266     }
0267     apvNb++;
0268   }
0269 
0270   if (!apvDetSet.empty())
0271     output_apvcm.push_back(std::move(apvDetSet));
0272 }