Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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 edm::DetSet<SiStripDigi>& hybridDigis,
0059                                                             edm::DetSet<SiStripDigi>& suppressedDigis,
0060                                                             uint16_t firstAPV) {
0061   uint16_t nAPVFlagged = 0;
0062   auto beginAPV = hybridDigis.begin();
0063   const auto indigis_end = hybridDigis.end();
0064   auto iAPV = firstAPV;
0065   while (beginAPV != indigis_end) {
0066     const auto endAPV = std::lower_bound(beginAPV, indigis_end, SiStripDigi((iAPV + 1) * 128, 0));
0067     const auto nDigisInAPV = std::distance(beginAPV, endAPV);
0068     if (nDigisInAPV > 64) {
0069       digivector_t workDigis(128, -1024);
0070       for (auto it = beginAPV; it != endAPV; ++it) {
0071         workDigis[it->strip() - 128 * iAPV] = it->adc() * 2 - 1024;
0072       }
0073       digivector_t workDigisPedSubtracted(workDigis);
0074       subtractorCMN->subtract(hybridDigis.id, iAPV, workDigis);
0075       const auto apvFlagged = restorer->inspectAndRestore(
0076           hybridDigis.id, iAPV, workDigisPedSubtracted, workDigis, subtractorCMN->getAPVsCM());
0077       nAPVFlagged += apvFlagged;
0078       if (getAPVFlags()[iAPV]) {
0079         suppressor->suppress(workDigis, iAPV, suppressedDigis);
0080       } else {  // bad APV: more than 64 but not flagged
0081         for (uint16_t i = 0; i != 128; ++i) {
0082           const auto digi = workDigisPedSubtracted[i];
0083           if (digi > 0) {
0084             suppressedDigis.push_back(SiStripDigi(iAPV * 128 + i, suppressor->truncate(digi)));
0085           }
0086         }
0087       }
0088     } else {  // already zero-suppressed, copy and truncate
0089       std::transform(beginAPV, endAPV, std::back_inserter(suppressedDigis), [this](const SiStripDigi inDigi) {
0090         return SiStripDigi(inDigi.strip(), suppressor->truncate(inDigi.adc()));
0091       });
0092     }
0093     beginAPV = endAPV;
0094     ++iAPV;
0095   }
0096   return nAPVFlagged;
0097 }
0098 
0099 /**
0100  * Zero-suppress virgin raw data.
0101  *
0102  * Subtracts pedestals and common-mode noise, and (optionally, if doAPVRestore)
0103  * re-evaluates and subtracts the baseline.
0104  *
0105  * @param id module DetId
0106  * @param firstAPV index of the first APV to consider
0107  * @param procRawDigis input (virgin raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
0108  * @param output zero-suppressed digis
0109  * @return number of restored APVs
0110  */
0111 uint16_t SiStripRawProcessingAlgorithms::suppressVirginRawData(uint32_t id,
0112                                                                uint16_t firstAPV,
0113                                                                digivector_t& procRawDigis,
0114                                                                edm::DetSet<SiStripDigi>& output) {
0115   subtractorPed->subtract(id, firstAPV * 128, procRawDigis);
0116   return suppressProcessedRawData(id, firstAPV, procRawDigis, output);
0117 }
0118 
0119 /**
0120  * Zero-suppress virgin raw data.
0121  *
0122  * Subtracts pedestals and common-mode noise, and (optionally, if doAPVRestore)
0123  * re-evaluates and subtracts the baseline.
0124  *
0125  * @param rawDigis input (virgin) raw digis
0126  * @param output zero-suppressed digis
0127  * @return number of restored APVs
0128  */
0129 uint16_t SiStripRawProcessingAlgorithms::suppressVirginRawData(const edm::DetSet<SiStripRawDigi>& rawDigis,
0130                                                                edm::DetSet<SiStripDigi>& output) {
0131   digivector_t rawdigis;
0132   rawdigis.reserve(rawDigis.size());
0133   std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
0134     return digi.adc();
0135   });
0136   return suppressVirginRawData(rawDigis.id, 0, rawdigis, output);
0137 }
0138 
0139 /**
0140  * Zero-suppress processed (pedestals-subtracted) raw data.
0141  *
0142  * Subtracts common-mode noise and (optionally, if doAPVRestore)
0143  * re-evaluates and subtracts the baseline.
0144  *
0145  * @param id module DetId
0146  * @param firstAPV index of the first APV to consider
0147  * @param procRawDigis input (processed raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
0148  * @param output zero-suppressed digis
0149  * @return number of restored APVs
0150  */
0151 uint16_t SiStripRawProcessingAlgorithms::suppressProcessedRawData(uint32_t id,
0152                                                                   uint16_t firstAPV,
0153                                                                   digivector_t& procRawDigis,
0154                                                                   edm::DetSet<SiStripDigi>& output) {
0155   digivector_t procRawDigisPedSubtracted;
0156 
0157   int16_t nAPVFlagged = 0;
0158   if (doAPVRestore)
0159     procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
0160   subtractorCMN->subtract(id, firstAPV, procRawDigis);
0161   if (doAPVRestore)
0162     nAPVFlagged =
0163         restorer->inspectAndRestore(id, firstAPV, procRawDigisPedSubtracted, procRawDigis, subtractorCMN->getAPVsCM());
0164   suppressor->suppress(procRawDigis, firstAPV, output);
0165   return nAPVFlagged;
0166 }
0167 
0168 /**
0169  * Zero-suppress processed (pedestals-subtracted) raw data.
0170  *
0171  * Subtracts common-mode noise and (optionally, if doAPVRestore)
0172  * re-evaluates and subtracts the baseline.
0173  *
0174  * @param rawDigis input (processed) raw digis
0175  * @param output zero-suppressed digis
0176  * @return number of restored APVs
0177  */
0178 uint16_t SiStripRawProcessingAlgorithms::suppressProcessedRawData(const edm::DetSet<SiStripRawDigi>& rawDigis,
0179                                                                   edm::DetSet<SiStripDigi>& output) {
0180   digivector_t rawdigis;
0181   rawdigis.reserve(rawDigis.size());
0182   std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
0183     return digi.adc();
0184   });
0185   return suppressProcessedRawData(rawDigis.id, 0, rawdigis, output);
0186 }
0187 
0188 /**
0189  * Zero-suppress virgin raw data in "hybrid" mode
0190  *
0191  * Subtracts pedestals (in 11bit mode, x->(x+1024-ped)/2) and common-mode noise, and inspects the digis then.
0192  * If not flagged by the hybrid APV inspector, the zero-suppression is performed as usual
0193  * (evaluation and subtraction of the baseline, truncation).
0194  * Otherwise, the pedestal-subtracted digis (as above) are saved in one 128-strip cluster.
0195  * Note: the APV restorer is used, it must be configured with APVInspectMode='HybridEmulation' if this method is called.
0196  *
0197  * @param id module DetId
0198  * @param firstAPV index of the first APV considered
0199  * @param procRawDigis input (virgin raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
0200  * @param output zero-suppressed digis (or pedestal-subtracted digis, see above)
0201  * @return number of restored APVs
0202  */
0203 uint16_t SiStripRawProcessingAlgorithms::convertVirginRawToHybrid(uint32_t id,
0204                                                                   uint16_t firstAPV,
0205                                                                   digivector_t& procRawDigis,
0206                                                                   edm::DetSet<SiStripDigi>& output) {
0207   digivector_t procRawDigisPedSubtracted;
0208 
0209   for (auto& digi : procRawDigis) {
0210     digi += 1024;
0211   }  // adding one MSB
0212 
0213   subtractorPed->subtract(id, firstAPV * 128, procRawDigis);  // all strips are pedestals subtracted
0214 
0215   for (auto& digi : procRawDigis) {
0216     digi /= 2;
0217   }
0218 
0219   procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
0220 
0221   subtractorCMN->subtract(id, firstAPV, procRawDigis);
0222 
0223   const auto nAPVFlagged = restorer->inspect(id, firstAPV, procRawDigis, subtractorCMN->getAPVsCM());
0224 
0225   for (auto& digi : procRawDigis) {
0226     digi *= 2;
0227   }
0228 
0229   const std::vector<bool>& apvf = getAPVFlags();
0230   const std::size_t nAPVs = procRawDigis.size() / 128;
0231   for (uint16_t iAPV = firstAPV; iAPV < nAPVs + firstAPV; ++iAPV) {
0232     if (apvf[iAPV]) {
0233       //GB 23/6/08: truncation should be done at the very beginning
0234       for (uint16_t i = 0; i < 128; ++i) {
0235         const int16_t digi = procRawDigisPedSubtracted[128 * (iAPV - firstAPV) + i];
0236         output.push_back(SiStripDigi(128 * iAPV + i, (digi < 0 ? 0 : suppressor->truncate(digi))));
0237       }
0238     } else {
0239       const auto firstDigiIt = std::begin(procRawDigis) + 128 * (iAPV - firstAPV);
0240       std::vector<int16_t> singleAPVdigi(firstDigiIt, firstDigiIt + 128);
0241       suppressor->suppress(singleAPVdigi, iAPV, output);
0242     }
0243   }
0244 
0245   return nAPVFlagged;
0246 }
0247 
0248 /**
0249  * Zero-suppress virgin raw data in "hybrid" mode
0250  *
0251  * Subtracts pedestals (in 11bit mode, x->(x+1024-ped)/2) and common-mode noise, and inspects the digis then.
0252  * If flagged by the hybrid APV inspector, the zero-suppression is performed as usual
0253  * (evaluation and subtraction of the baseline, truncation).
0254  * Otherwise, the pedestal-subtracted digis are saved in one 128-strip cluster.
0255  *
0256  * @param rawDigis input (virgin) raw digis
0257  * @param output zero-suppressed digis (or pedestal-subtracted digis, see above)
0258  * @return number of restored APVs
0259  */
0260 uint16_t SiStripRawProcessingAlgorithms::convertVirginRawToHybrid(const edm::DetSet<SiStripRawDigi>& rawDigis,
0261                                                                   edm::DetSet<SiStripDigi>& suppressedDigis) {
0262   digivector_t rawdigis;
0263   rawdigis.reserve(rawDigis.size());
0264   std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
0265     return digi.adc();
0266   });
0267   return convertVirginRawToHybrid(rawDigis.id, 0, rawdigis, suppressedDigis);
0268 }