Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-11-19 23:20:18

0001 #include "RecoLocalTracker/SiStripZeroSuppression/interface/SiStripRawProcessingAlgorithms.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 "DataFormats/SiStripDetId/interface/SiStripDetId.h"
0011 #include "RecoLocalTracker/SiStripZeroSuppression/interface/SiStripRawProcessingFactory.h"
0012 #include <memory>
0013 
0014 SiStripRawProcessingAlgorithms::SiStripRawProcessingAlgorithms(edm::ConsumesCollector iC,
0015                                                                std::unique_ptr<SiStripPedestalsSubtractor> ped,
0016                                                                std::unique_ptr<SiStripCommonModeNoiseSubtractor> cmn,
0017                                                                std::unique_ptr<SiStripFedZeroSuppression> zs,
0018                                                                std::unique_ptr<SiStripAPVRestorer> res,
0019                                                                bool doAPVRest,
0020                                                                bool useCMMap)
0021     : subtractorPed(std::move(ped)),
0022       subtractorCMN(std::move(cmn)),
0023       suppressor(std::move(zs)),
0024       restorer(std::move(res)),
0025       doAPVRestore(doAPVRest),
0026       useCMMeanMap(useCMMap),
0027       tkGeomToken_(iC.esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()) {}
0028 
0029 void SiStripRawProcessingAlgorithms::initialize(const edm::EventSetup& es) {
0030   subtractorPed->init(es);
0031   subtractorCMN->init(es);
0032   suppressor->init(es);
0033   if (restorer.get())
0034     restorer->init(es);
0035 
0036   trGeo = &es.getData(tkGeomToken_);
0037 }
0038 
0039 void SiStripRawProcessingAlgorithms::initialize(const edm::EventSetup& es, const edm::Event& e) {
0040   initialize(es);
0041   if (restorer.get() && doAPVRestore && useCMMeanMap)
0042     restorer->loadMeanCMMap(e);
0043 }
0044 
0045 /**
0046  * Zero-suppress "hybrid" raw data
0047  *
0048  * Subtracts common-mode noise, and inspects the digis then.
0049  * If flagged by the APV inspector, the zero-suppression is performed as usual.
0050  * Otherwise, the positive inputs are copied.
0051  *
0052  * @param hybridDigis input ADCs in ZS format (regular ZS or "hybrid", i.e. processed as x->(x+1024-ped)/2)
0053  * @param suppressedDigis zero-suppressed digis
0054  * @param firstAPV (optional) number of the first APV for which digis are should be handled (otherwise all present)
0055  *
0056  * @return number of restored APVs
0057  */
0058 uint16_t SiStripRawProcessingAlgorithms::suppressHybridData(const uint16_t maxNStrips,
0059                                                             const edm::DetSet<SiStripDigi>& hybridDigis,
0060                                                             edm::DetSet<SiStripDigi>& suppressedDigis,
0061                                                             uint16_t firstAPV) {
0062   uint16_t nAPVFlagged = 0;  // Count of flagged APVs
0063   auto currentDigi = hybridDigis.begin();
0064   const auto endDigi = hybridDigis.end();
0065   auto currentAPV = firstAPV;
0066 
0067   // Loop through the APVs in the range
0068   while (currentDigi != endDigi) {
0069     // Determine the range of digis belonging to the current APV
0070     const auto nextAPVBoundary = SiStripDigi((currentAPV + 1) * 128, 0);
0071 
0072     // Reject any APV larger than the max possible
0073     if (nextAPVBoundary.strip() > maxNStrips) {
0074       edm::LogError("SiStripRawProcessingAlgorithms")
0075           << "In DetId " << suppressedDigis.id << " encountered APV boundary with strip number "
0076           << nextAPVBoundary.strip() << ", which exceeds the maximum allowed value for this module (" << maxNStrips
0077           << "). Exiting loop.";
0078       break;
0079     }
0080 
0081     const auto nextAPVDigi = std::lower_bound(currentDigi, endDigi, nextAPVBoundary);
0082     const auto nDigisInAPV = std::distance(currentDigi, nextAPVDigi);
0083 
0084     // Handle based on the number of digis in the current APV
0085     if (nDigisInAPV > 64) {
0086       // Process hybrid data for noisy APV
0087       digivector_t workDigis(128, -1024);
0088 
0089       // Populate `workDigis` with values from `currentDigi`
0090       for (auto it = currentDigi; it != nextAPVDigi; ++it) {
0091         workDigis[it->strip() - 128 * currentAPV] = it->adc() * 2 - 1024;
0092       }
0093 
0094       // Perform pedestal subtraction
0095       digivector_t workDigisPedSubtracted(workDigis);
0096       subtractorCMN->subtract(hybridDigis.id, currentAPV, workDigis);
0097 
0098       // Inspect and restore digis
0099       const auto apvFlagged = restorer->inspectAndRestore(
0100           hybridDigis.id, currentAPV, workDigisPedSubtracted, workDigis, subtractorCMN->getAPVsCM());
0101       nAPVFlagged += apvFlagged;
0102 
0103       // Process based on the APV flag
0104       if (getAPVFlags()[currentAPV]) {
0105         // Suppress flagged APVs
0106         suppressor->suppress(workDigis, currentAPV, suppressedDigis);
0107       } else {
0108         // Handle bad APV: not flagged but exceeds threshold
0109         for (uint16_t i = 0; i < 128; ++i) {
0110           const auto digi = workDigisPedSubtracted[i];
0111           if (digi > 0) {
0112             suppressedDigis.push_back(SiStripDigi(currentAPV * 128 + i, suppressor->truncate(digi)));
0113           }
0114         }
0115       }
0116     } else {
0117       // Already zero-suppressed: copy and truncate
0118       std::transform(currentDigi, nextAPVDigi, std::back_inserter(suppressedDigis), [this](const SiStripDigi& inDigi) {
0119         return SiStripDigi(inDigi.strip(), suppressor->truncate(inDigi.adc()));
0120       });
0121     }
0122 
0123     // Move to the next APV
0124     currentDigi = nextAPVDigi;
0125     ++currentAPV;
0126   }
0127 
0128   return nAPVFlagged;
0129 }
0130 
0131 /**
0132  * Zero-suppress virgin raw data.
0133  *
0134  * Subtracts pedestals and common-mode noise, and (optionally, if doAPVRestore)
0135  * re-evaluates and subtracts the baseline.
0136  *
0137  * @param id module DetId
0138  * @param firstAPV index of the first APV to consider
0139  * @param procRawDigis input (virgin raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
0140  * @param output zero-suppressed digis
0141  * @return number of restored APVs
0142  */
0143 uint16_t SiStripRawProcessingAlgorithms::suppressVirginRawData(uint32_t id,
0144                                                                uint16_t firstAPV,
0145                                                                digivector_t& procRawDigis,
0146                                                                edm::DetSet<SiStripDigi>& output) {
0147   subtractorPed->subtract(id, firstAPV * 128, procRawDigis);
0148   return suppressProcessedRawData(id, firstAPV, procRawDigis, output);
0149 }
0150 
0151 /**
0152  * Zero-suppress virgin raw data.
0153  *
0154  * Subtracts pedestals and common-mode noise, and (optionally, if doAPVRestore)
0155  * re-evaluates and subtracts the baseline.
0156  *
0157  * @param rawDigis input (virgin) raw digis
0158  * @param output zero-suppressed digis
0159  * @return number of restored APVs
0160  */
0161 uint16_t SiStripRawProcessingAlgorithms::suppressVirginRawData(const edm::DetSet<SiStripRawDigi>& rawDigis,
0162                                                                edm::DetSet<SiStripDigi>& output) {
0163   digivector_t rawdigis;
0164   rawdigis.reserve(rawDigis.size());
0165   std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
0166     return digi.adc();
0167   });
0168   return suppressVirginRawData(rawDigis.id, 0, rawdigis, output);
0169 }
0170 
0171 /**
0172  * Zero-suppress processed (pedestals-subtracted) raw data.
0173  *
0174  * Subtracts common-mode noise and (optionally, if doAPVRestore)
0175  * re-evaluates and subtracts the baseline.
0176  *
0177  * @param id module DetId
0178  * @param firstAPV index of the first APV to consider
0179  * @param procRawDigis input (processed raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
0180  * @param output zero-suppressed digis
0181  * @return number of restored APVs
0182  */
0183 uint16_t SiStripRawProcessingAlgorithms::suppressProcessedRawData(uint32_t id,
0184                                                                   uint16_t firstAPV,
0185                                                                   digivector_t& procRawDigis,
0186                                                                   edm::DetSet<SiStripDigi>& output) {
0187   digivector_t procRawDigisPedSubtracted;
0188 
0189   int16_t nAPVFlagged = 0;
0190   if (doAPVRestore)
0191     procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
0192   subtractorCMN->subtract(id, firstAPV, procRawDigis);
0193   if (doAPVRestore)
0194     nAPVFlagged =
0195         restorer->inspectAndRestore(id, firstAPV, procRawDigisPedSubtracted, procRawDigis, subtractorCMN->getAPVsCM());
0196   suppressor->suppress(procRawDigis, firstAPV, output);
0197   return nAPVFlagged;
0198 }
0199 
0200 /**
0201  * Zero-suppress processed (pedestals-subtracted) raw data.
0202  *
0203  * Subtracts common-mode noise and (optionally, if doAPVRestore)
0204  * re-evaluates and subtracts the baseline.
0205  *
0206  * @param rawDigis input (processed) raw digis
0207  * @param output zero-suppressed digis
0208  * @return number of restored APVs
0209  */
0210 uint16_t SiStripRawProcessingAlgorithms::suppressProcessedRawData(const edm::DetSet<SiStripRawDigi>& rawDigis,
0211                                                                   edm::DetSet<SiStripDigi>& output) {
0212   digivector_t rawdigis;
0213   rawdigis.reserve(rawDigis.size());
0214   std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
0215     return digi.adc();
0216   });
0217   return suppressProcessedRawData(rawDigis.id, 0, rawdigis, output);
0218 }
0219 
0220 /**
0221  * Zero-suppress virgin raw data in "hybrid" mode
0222  *
0223  * Subtracts pedestals (in 11bit mode, x->(x+1024-ped)/2) and common-mode noise, and inspects the digis then.
0224  * If not flagged by the hybrid APV inspector, the zero-suppression is performed as usual
0225  * (evaluation and subtraction of the baseline, truncation).
0226  * Otherwise, the pedestal-subtracted digis (as above) are saved in one 128-strip cluster.
0227  * Note: the APV restorer is used, it must be configured with APVInspectMode='HybridEmulation' if this method is called.
0228  *
0229  * @param id module DetId
0230  * @param firstAPV index of the first APV considered
0231  * @param procRawDigis input (virgin raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
0232  * @param output zero-suppressed digis (or pedestal-subtracted digis, see above)
0233  * @return number of restored APVs
0234  */
0235 uint16_t SiStripRawProcessingAlgorithms::convertVirginRawToHybrid(uint32_t id,
0236                                                                   uint16_t firstAPV,
0237                                                                   digivector_t& procRawDigis,
0238                                                                   edm::DetSet<SiStripDigi>& output) {
0239   digivector_t procRawDigisPedSubtracted;
0240 
0241   for (auto& digi : procRawDigis) {
0242     digi += 1024;
0243   }  // adding one MSB
0244 
0245   subtractorPed->subtract(id, firstAPV * 128, procRawDigis);  // all strips are pedestals subtracted
0246 
0247   for (auto& digi : procRawDigis) {
0248     digi /= 2;
0249   }
0250 
0251   procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
0252 
0253   subtractorCMN->subtract(id, firstAPV, procRawDigis);
0254 
0255   const auto nAPVFlagged = restorer->inspect(id, firstAPV, procRawDigis, subtractorCMN->getAPVsCM());
0256 
0257   for (auto& digi : procRawDigis) {
0258     digi *= 2;
0259   }
0260 
0261   const std::vector<bool>& apvf = getAPVFlags();
0262   const std::size_t nAPVs = procRawDigis.size() / 128;
0263   for (uint16_t iAPV = firstAPV; iAPV < nAPVs + firstAPV; ++iAPV) {
0264     if (apvf[iAPV]) {
0265       //GB 23/6/08: truncation should be done at the very beginning
0266       for (uint16_t i = 0; i < 128; ++i) {
0267         const int16_t digi = procRawDigisPedSubtracted[128 * (iAPV - firstAPV) + i];
0268         output.push_back(SiStripDigi(128 * iAPV + i, (digi < 0 ? 0 : suppressor->truncate(digi))));
0269       }
0270     } else {
0271       const auto firstDigiIt = std::begin(procRawDigis) + 128 * (iAPV - firstAPV);
0272       std::vector<int16_t> singleAPVdigi(firstDigiIt, firstDigiIt + 128);
0273       suppressor->suppress(singleAPVdigi, iAPV, output);
0274     }
0275   }
0276 
0277   return nAPVFlagged;
0278 }
0279 
0280 /**
0281  * Zero-suppress virgin raw data in "hybrid" mode
0282  *
0283  * Subtracts pedestals (in 11bit mode, x->(x+1024-ped)/2) and common-mode noise, and inspects the digis then.
0284  * If flagged by the hybrid APV inspector, the zero-suppression is performed as usual
0285  * (evaluation and subtraction of the baseline, truncation).
0286  * Otherwise, the pedestal-subtracted digis are saved in one 128-strip cluster.
0287  *
0288  * @param rawDigis input (virgin) raw digis
0289  * @param output zero-suppressed digis (or pedestal-subtracted digis, see above)
0290  * @return number of restored APVs
0291  */
0292 uint16_t SiStripRawProcessingAlgorithms::convertVirginRawToHybrid(const edm::DetSet<SiStripRawDigi>& rawDigis,
0293                                                                   edm::DetSet<SiStripDigi>& suppressedDigis) {
0294   digivector_t rawdigis;
0295   rawdigis.reserve(rawDigis.size());
0296   std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
0297     return digi.adc();
0298   });
0299   return convertVirginRawToHybrid(rawDigis.id, 0, rawdigis, suppressedDigis);
0300 }