Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:19:51

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiStripMeanCMExtractor
0004 // Class:      SiStripMeanCMExtractor
0005 //
0006 /**\class SiStripMeanCMExtractor SiStripMeanCMExtractor.cc RecoLocalTracker/SiStripMeanCMExtractor/src/SiStripMeanCMExtractor.cc
0007 
0008  Description: [one line class summary]
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 //
0014 // Original Author:  Ivan Amos Cali,32 4-A08,+41227673039,
0015 //         Created:  Wed Oct 13 11:50:47 CEST 2010
0016 //
0017 //
0018 
0019 // system include files
0020 #include <sstream>
0021 #include <memory>
0022 #include <list>
0023 #include <algorithm>
0024 #include <cassert>
0025 
0026 // user include files
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/EventSetup.h"
0029 #include "DataFormats/Common/interface/Handle.h"
0030 #include "FWCore/Framework/interface/EventSetup.h"
0031 #include "FWCore/Framework/interface/Frameworkfwd.h"
0032 #include "FWCore/Framework/interface/one/EDProducer.h"
0033 
0034 #include "FWCore/Framework/interface/Event.h"
0035 #include "FWCore/Framework/interface/MakerMacros.h"
0036 #include "FWCore/Framework/interface/ESWatcher.h"
0037 
0038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0039 
0040 #include "DataFormats/Common/interface/DetSet.h"
0041 #include "DataFormats/Common/interface/DetSetVector.h"
0042 #include "DataFormats/DetId/interface/DetId.h"
0043 #include "DataFormats/SiStripDigi/interface/SiStripProcessedRawDigi.h"
0044 #include "DataFormats/SiStripDigi/interface/SiStripRawDigi.h"
0045 #include "CondFormats/SiStripObjects/interface/SiStripPedestals.h"
0046 #include "CondFormats/DataRecord/interface/SiStripPedestalsRcd.h"
0047 
0048 typedef std::map<uint32_t, std::vector<float>> CMMap;
0049 
0050 class SiStripMeanCMExtractor : public edm::one::EDProducer<> {
0051 public:
0052   explicit SiStripMeanCMExtractor(const edm::ParameterSet&);
0053   ~SiStripMeanCMExtractor() override;
0054   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0055 
0056 private:
0057   void beginJob() override;
0058   void produce(edm::Event&, const edm::EventSetup&) override;
0059 
0060   void init(const edm::EventSetup&);
0061 
0062   const edm::ESGetToken<SiStripPedestals, SiStripPedestalsRcd> pedestalsToken_;
0063   const SiStripPedestals* pedestalHandle_;
0064   edm::ESWatcher<SiStripPedestalsRcd> pedestalsWatcher_;
0065 
0066   void StoreMean(const edm::DetSetVector<SiStripProcessedRawDigi>&);
0067   void ConvertMeanMapToDetSetVector(std::vector<edm::DetSet<SiStripProcessedRawDigi>>&);
0068   void CMExtractorFromPedestals(const edm::DetSetVector<SiStripRawDigi>&,
0069                                 std::vector<edm::DetSet<SiStripProcessedRawDigi>>&);
0070   edm::InputTag _inputTag;
0071   const edm::EDGetTokenT<edm::DetSetVector<SiStripProcessedRawDigi>> _inputCMToken;
0072   const edm::EDGetTokenT<edm::DetSetVector<SiStripRawDigi>> _inputToken;
0073   std::string _Algorithm;
0074   uint16_t _nEventsToUse;
0075   uint16_t _actualEvent;
0076 
0077   CMMap _cmMap;  //it contains the sum of the CM calculated before. The normalization for the number of events it is done at the end when it is written in the DetSetVector.
0078 };
0079 
0080 SiStripMeanCMExtractor::SiStripMeanCMExtractor(const edm::ParameterSet& conf)
0081     : pedestalsToken_(esConsumes()),
0082       _inputTag(conf.getParameter<edm::InputTag>("CMCollection")),
0083       _inputCMToken(consumes<edm::DetSetVector<SiStripProcessedRawDigi>>(_inputTag)),
0084       _inputToken(consumes<edm::DetSetVector<SiStripRawDigi>>(_inputTag)),
0085       _Algorithm(conf.getParameter<std::string>("Algorithm")),
0086       _nEventsToUse(conf.getParameter<uint32_t>("NEvents")) {
0087   if (_nEventsToUse < 1)
0088     _nEventsToUse = 1;
0089   produces<edm::DetSetVector<SiStripProcessedRawDigi>>("MEANAPVCM");
0090 }
0091 
0092 SiStripMeanCMExtractor::~SiStripMeanCMExtractor() = default;
0093 
0094 void SiStripMeanCMExtractor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0095   edm::ParameterSetDescription desc;
0096   desc.add<edm::InputTag>("CMCollection", edm::InputTag("siStripZeroSuppression", "APVCM"));
0097   desc.add<std::string>("Algorithm", "Pedestals");
0098   desc.add<uint32_t>("NEvents", 100);
0099   descriptions.add("SiStripMeanCMExtractor", desc);
0100 }
0101 
0102 void SiStripMeanCMExtractor::init(const edm::EventSetup& es) {
0103   if (pedestalsWatcher_.check(es)) {
0104     pedestalHandle_ = &es.getData(pedestalsToken_);
0105   }
0106 }
0107 
0108 // ------------ method called to produce the data  ------------
0109 void SiStripMeanCMExtractor::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0110   using namespace edm;
0111 
0112   //if(_actualEvent > _nEventsToUse) return;
0113 
0114   std::vector<edm::DetSet<SiStripProcessedRawDigi>> meancm;
0115 
0116   if (_Algorithm == "StoredCM") {
0117     edm::Handle<edm::DetSetVector<SiStripProcessedRawDigi>> inputCM;
0118     iEvent.getByToken(_inputCMToken, inputCM);
0119 
0120     this->StoreMean(*inputCM);
0121     this->ConvertMeanMapToDetSetVector(meancm);
0122 
0123   } else if (_Algorithm == "Pedestals") {
0124     this->init(iSetup);
0125 
0126     edm::Handle<edm::DetSetVector<SiStripRawDigi>> input;
0127     iEvent.getByToken(_inputToken, input);
0128 
0129     this->CMExtractorFromPedestals(*input, meancm);
0130   }
0131 
0132   ++_actualEvent;
0133 
0134   iEvent.put(std::make_unique<edm::DetSetVector<SiStripProcessedRawDigi>>(meancm), "MEANAPVCM");
0135 }
0136 
0137 void SiStripMeanCMExtractor::CMExtractorFromPedestals(const edm::DetSetVector<SiStripRawDigi>& input,
0138                                                       std::vector<edm::DetSet<SiStripProcessedRawDigi>>& meancm) {
0139   meancm.clear();
0140   meancm.reserve(15000);
0141 
0142   for (edm::DetSetVector<SiStripRawDigi>::const_iterator rawDigis = input.begin(); rawDigis != input.end();
0143        rawDigis++) {
0144     SiStripPedestals::Range detPedestalRange = pedestalHandle_->getRange(rawDigis->id);
0145     edm::DetSet<SiStripProcessedRawDigi> MeanCMDetSet(rawDigis->id);
0146 
0147     for (uint16_t APV = 0; APV < rawDigis->size() / 128; ++APV) {
0148       uint16_t MinPed = 0;
0149       for (uint16_t strip = APV * 128; strip < (APV + 1) * 128; ++strip) {
0150         uint16_t ped = (uint16_t)pedestalHandle_->getPed(strip, detPedestalRange);
0151         if (ped < MinPed)
0152           MinPed = ped;
0153       }
0154       if (MinPed > 128)
0155         MinPed = 128;
0156       MeanCMDetSet.push_back(MinPed);
0157     }
0158 
0159     meancm.push_back(MeanCMDetSet);
0160   }
0161 }
0162 
0163 void SiStripMeanCMExtractor::StoreMean(const edm::DetSetVector<SiStripProcessedRawDigi>& Input) {
0164   uint32_t detId;
0165   CMMap::iterator itMap;
0166   edm::DetSetVector<SiStripProcessedRawDigi>::const_iterator itInput;
0167 
0168   for (itInput = Input.begin(); itInput != Input.end(); ++itInput) {
0169     detId = itInput->id;
0170     itMap = _cmMap.find(detId);
0171     edm::DetSet<SiStripProcessedRawDigi>::const_iterator itCM;
0172     std::vector<float> MeanCMNValue;
0173     MeanCMNValue.clear();
0174     if (itMap != _cmMap.end()) {  //the detId was already found
0175       std::vector<float>& MapContent = itMap->second;
0176       std::vector<float>::iterator itMapVector = MapContent.begin();
0177       for (itCM = itInput->begin(); itCM != itInput->end(); ++itCM, ++itMapVector) {
0178         MeanCMNValue.push_back(itCM->adc() + *itMapVector);
0179       }
0180       _cmMap.erase(itMap);
0181       _cmMap.insert(itMap, std::pair<uint32_t, std::vector<float>>(detId, MeanCMNValue));
0182     } else {  //no detId found
0183       for (itCM = itInput->begin(); itCM != itInput->end(); ++itCM)
0184         MeanCMNValue.push_back(itCM->adc());
0185       _cmMap.insert(std::pair<uint32_t, std::vector<float>>(detId, MeanCMNValue));
0186     }
0187   }
0188 }
0189 
0190 void SiStripMeanCMExtractor::ConvertMeanMapToDetSetVector(std::vector<edm::DetSet<SiStripProcessedRawDigi>>& meancm) {
0191   CMMap::iterator itMap;
0192   std::vector<float>::const_iterator itMapVector;
0193 
0194   meancm.clear();
0195   meancm.reserve(15000);
0196 
0197   for (itMap = _cmMap.begin(); itMap != _cmMap.end(); ++itMap) {
0198     edm::DetSet<SiStripProcessedRawDigi> MeanCMDetSet(itMap->first);
0199     for (itMapVector = (itMap->second).begin(); itMapVector != (itMap->second).end(); ++itMapVector)
0200       MeanCMDetSet.push_back(*itMapVector / (float)_actualEvent);
0201     meancm.push_back(MeanCMDetSet);
0202   }
0203 }
0204 void SiStripMeanCMExtractor::beginJob() {
0205   _actualEvent = 1;
0206 
0207   _cmMap.clear();
0208 }
0209 
0210 DEFINE_FWK_MODULE(SiStripMeanCMExtractor);