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