Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:33

0001 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibrationForHLTHost.h"
0002 #include "CalibTracker/Records/interface/SiPixelGainCalibrationForHLTSoARcd.h"
0003 #include "CondFormats/DataRecord/interface/SiPixelGainCalibrationForHLTRcd.h"
0004 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibrationForHLT.h"
0005 #include "DataFormats/Portable/interface/alpaka/PortableCollection.h"
0006 #include "FWCore/Framework/interface/ESProducer.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/ESHandle.h"
0009 #include "FWCore/Framework/interface/ModuleFactory.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0012 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0013 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0014 
0015 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0016 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0017 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0018 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020 
0021 #include <memory>
0022 
0023 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0024 
0025   class SiPixelGainCalibrationForHLTSoAESProducer : public ESProducer {
0026   public:
0027     explicit SiPixelGainCalibrationForHLTSoAESProducer(const edm::ParameterSet& iConfig);
0028     std::unique_ptr<SiPixelGainCalibrationForHLTHost> produce(const SiPixelGainCalibrationForHLTSoARcd& iRecord);
0029 
0030     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0031 
0032   private:
0033     edm::ESGetToken<SiPixelGainCalibrationForHLT, SiPixelGainCalibrationForHLTRcd> gainsToken_;
0034     edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geometryToken_;
0035   };
0036 
0037   SiPixelGainCalibrationForHLTSoAESProducer::SiPixelGainCalibrationForHLTSoAESProducer(const edm::ParameterSet& iConfig)
0038       : ESProducer(iConfig) {
0039     auto cc = setWhatProduced(this);
0040     gainsToken_ = cc.consumes();
0041     geometryToken_ = cc.consumes();
0042   }
0043 
0044   void SiPixelGainCalibrationForHLTSoAESProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0045     edm::ParameterSetDescription desc;
0046     descriptions.addWithDefaultLabel(desc);
0047   }
0048 
0049   std::unique_ptr<SiPixelGainCalibrationForHLTHost> SiPixelGainCalibrationForHLTSoAESProducer::produce(
0050       const SiPixelGainCalibrationForHLTSoARcd& iRecord) {
0051     auto const& gains = iRecord.get(gainsToken_);
0052     auto const& geom = iRecord.get(geometryToken_);
0053 
0054     auto product = std::make_unique<SiPixelGainCalibrationForHLTHost>(gains.data().size(), cms::alpakatools::host());
0055 
0056     // bizzarre logic (looking for fist strip-det) don't ask
0057     auto const& dus = geom.detUnits();
0058     unsigned int n_detectors = dus.size();
0059     for (unsigned int i = 1; i < 7; ++i) {
0060       const auto offset = geom.offsetDU(GeomDetEnumerators::tkDetEnum[i]);
0061       if (offset != dus.size() && dus[offset]->type().isTrackerStrip()) {
0062         if (n_detectors > offset)
0063           n_detectors = offset;
0064       }
0065     }
0066 
0067     LogDebug("SiPixelGainCalibrationForHLTSoA")
0068         << "caching calibs for " << n_detectors << " pixel detectors of size " << gains.data().size() << '\n'
0069         << "sizes " << sizeof(char) << ' ' << sizeof(uint8_t) << ' ' << sizeof(siPixelGainsSoA::DecodingStructure);
0070 
0071     for (size_t i = 0; i < gains.data().size(); i = i + 2) {
0072       product->view().v_pedestals()[i / 2].gain = gains.data()[i];
0073       product->view().v_pedestals()[i / 2].ped = gains.data()[i + 1];
0074     }
0075 
0076     //std::copy here
0077     // do not read back from the (possibly write-combined) memory buffer
0078     auto minPed = gains.getPedLow();
0079     auto maxPed = gains.getPedHigh();
0080     auto minGain = gains.getGainLow();
0081     auto maxGain = gains.getGainHigh();
0082     auto nBinsToUseForEncoding = 253;
0083 
0084     // we will simplify later (not everything is needed....)
0085     product->view().minPed() = minPed;
0086     product->view().maxPed() = maxPed;
0087     product->view().minGain() = minGain;
0088     product->view().maxGain() = maxGain;
0089 
0090     product->view().numberOfRowsAveragedOver() = 80;
0091     product->view().nBinsToUseForEncoding() = nBinsToUseForEncoding;
0092     product->view().deadFlag() = 255;
0093     product->view().noisyFlag() = 254;
0094 
0095     product->view().pedPrecision() = static_cast<float>(maxPed - minPed) / nBinsToUseForEncoding;
0096     product->view().gainPrecision() = static_cast<float>(maxGain - minGain) / nBinsToUseForEncoding;
0097 
0098     LogDebug("SiPixelGainCalibrationForHLTSoA")
0099         << "precisions g " << product->view().pedPrecision() << ' ' << product->view().gainPrecision();
0100 
0101     // fill the index map
0102     auto const& ind = gains.getIndexes();
0103     LogDebug("SiPixelGainCalibrationForHLTSoA") << ind.size() << " " << n_detectors;
0104 
0105     for (auto i = 0U; i < n_detectors; ++i) {
0106       auto p = std::lower_bound(
0107           ind.begin(), ind.end(), dus[i]->geographicalId().rawId(), SiPixelGainCalibrationForHLT::StrictWeakOrdering());
0108       assert(p != ind.end() && p->detid == dus[i]->geographicalId());
0109       assert(p->iend <= gains.data().size());
0110       assert(p->iend >= p->ibegin);
0111       assert(0 == p->ibegin % 2);
0112       assert(0 == p->iend % 2);
0113       assert(p->ibegin != p->iend);
0114       assert(p->ncols > 0);
0115 
0116       product->view().modStarts()[i] = p->ibegin;
0117       product->view().modEnds()[i] = p->iend;
0118       product->view().modCols()[i] = p->ncols;
0119 
0120       if (ind[i].detid != dus[i]->geographicalId())
0121         LogDebug("SiPixelGainCalibrationForHLTSoA") << ind[i].detid << "!=" << dus[i]->geographicalId();
0122     }
0123 
0124     return product;
0125   }
0126 
0127 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0128 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(SiPixelGainCalibrationForHLTSoAESProducer);