Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-04 00:29:53

0001 #include "FWCore/Framework/interface/ESTransientHandle.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 
0004 #include "CondFormats/DataRecord/interface/EcalADCToGeVConstantRcd.h"
0005 #include "CondFormats/DataRecord/interface/EcalChannelStatusRcd.h"
0006 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
0007 #include "CondFormats/DataRecord/interface/EcalLaserAPDPNRatiosRcd.h"
0008 #include "CondFormats/DataRecord/interface/EcalLaserAPDPNRatiosRefRcd.h"
0009 #include "CondFormats/DataRecord/interface/EcalLaserAlphasRcd.h"
0010 #include "CondFormats/DataRecord/interface/EcalLinearCorrectionsRcd.h"
0011 #include "CondFormats/DataRecord/interface/EcalTimeCalibConstantsRcd.h"
0012 #include "CondFormats/DataRecord/interface/EcalTimeOffsetConstantRcd.h"
0013 #include "CondFormats/EcalObjects/interface/EcalADCToGeVConstant.h"
0014 #include "CondFormats/EcalObjects/interface/EcalChannelStatus.h"
0015 #include "CondFormats/EcalObjects/interface/EcalIntercalibConstants.h"
0016 #include "CondFormats/EcalObjects/interface/EcalLaserAPDPNRatios.h"
0017 #include "CondFormats/EcalObjects/interface/EcalLaserAPDPNRatiosRef.h"
0018 #include "CondFormats/EcalObjects/interface/EcalLaserAlphas.h"
0019 #include "CondFormats/EcalObjects/interface/EcalLinearCorrections.h"
0020 #include "CondFormats/EcalObjects/interface/EcalTimeCalibConstants.h"
0021 #include "CondFormats/EcalObjects/interface/EcalTimeOffsetConstant.h"
0022 
0023 #include "CondFormats/EcalObjects/interface/alpaka/EcalRecHitConditionsDevice.h"
0024 #include "CondFormats/EcalObjects/interface/EcalRecHitConditionsSoA.h"
0025 #include "CondFormats/DataRecord/interface/EcalRecHitConditionsRcd.h"
0026 
0027 #include "DataFormats/EcalDigi/interface/EcalConstants.h"
0028 
0029 #include "DataFormats/EcalDetId/interface/EcalElectronicsId.h"
0030 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESGetToken.h"
0031 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0032 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0033 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0034 #include "HeterogeneousCore/AlpakaInterface/interface/host.h"
0035 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0036 
0037 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0038 
0039   class EcalRecHitConditionsESProducer : public ESProducer {
0040   public:
0041     EcalRecHitConditionsESProducer(edm::ParameterSet const& iConfig)
0042         : ESProducer(iConfig), isPhase2_{iConfig.getParameter<bool>("isPhase2")} {
0043       auto cc = setWhatProduced(this);
0044       adcToGeVConstantToken_ = cc.consumes();
0045       intercalibConstantsToken_ = cc.consumes();
0046       channelStatusToken_ = cc.consumes();
0047       laserAPDPNRatiosToken_ = cc.consumes();
0048       laserAPDPNRatiosRefToken_ = cc.consumes();
0049       laserAlphasToken_ = cc.consumes();
0050       linearCorrectionsToken_ = cc.consumes();
0051       timeCalibConstantsToken_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("timeCalibTag"));
0052       timeOffsetConstantToken_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("timeOffsetTag"));
0053     }
0054 
0055     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0056       edm::ParameterSetDescription desc;
0057       desc.add<edm::ESInputTag>("timeCalibTag", edm::ESInputTag());
0058       desc.add<edm::ESInputTag>("timeOffsetTag", edm::ESInputTag());
0059       desc.add<bool>("isPhase2", false);
0060       descriptions.addWithDefaultLabel(desc);
0061     }
0062 
0063     std::unique_ptr<EcalRecHitConditionsHost> produce(EcalRecHitConditionsRcd const& iRecord) {
0064       auto const& adcToGeVConstantData = iRecord.get(adcToGeVConstantToken_);
0065       auto const& intercalibConstantsData = iRecord.get(intercalibConstantsToken_);
0066       auto const& channelStatusData = iRecord.get(channelStatusToken_);
0067       auto const& laserAPDPNRatiosData = iRecord.get(laserAPDPNRatiosToken_);
0068       auto const& laserAPDPNRatiosRefData = iRecord.get(laserAPDPNRatiosRefToken_);
0069       auto const& laserAlphasData = iRecord.get(laserAlphasToken_);
0070       auto const& linearCorrectionsData = iRecord.get(linearCorrectionsToken_);
0071       auto const& timeCalibConstantsData = iRecord.get(timeCalibConstantsToken_);
0072       auto const& timeOffsetConstantData = iRecord.get(timeOffsetConstantToken_);
0073 
0074       const auto barrelSize = channelStatusData.barrelItems().size();
0075       auto numberOfXtals = barrelSize;
0076       if (!isPhase2_) {
0077         numberOfXtals += channelStatusData.endcapItems().size();
0078       }
0079 
0080       auto product = std::make_unique<EcalRecHitConditionsHost>(numberOfXtals, cms::alpakatools::host());
0081       auto view = product->view();
0082 
0083       // Filling crystal level conditions
0084       auto const& intercalibConstantsEB = intercalibConstantsData.barrelItems();
0085       auto const& channelStatusEB = channelStatusData.barrelItems();
0086       auto const& laserAPDPNRatiosLaserEB = laserAPDPNRatiosData.getLaserMap().barrelItems();
0087       auto const& laserAPDPNRatiosTime = laserAPDPNRatiosData.getTimeMap();
0088       auto const& laserAPDPNRatiosRefEB = laserAPDPNRatiosRefData.barrelItems();
0089       auto const& laserAlphasEB = laserAlphasData.barrelItems();
0090       auto const& linearCorrectionsEB = linearCorrectionsData.getValueMap().barrelItems();
0091       auto const& linearCorrectionsTime = linearCorrectionsData.getTimeMap();
0092       auto const& timeCalibConstantsEB = timeCalibConstantsData.barrelItems();
0093 
0094       // barrel conditions
0095       for (unsigned int i = 0; i < barrelSize; ++i) {
0096         auto vi = view[i];
0097 
0098         vi.intercalibConstants() = intercalibConstantsEB[i];
0099         vi.channelStatus() = channelStatusEB[i].getEncodedStatusCode();
0100 
0101         vi.laserAPDPNRatios_p1() = laserAPDPNRatiosLaserEB[i].p1;
0102         vi.laserAPDPNRatios_p2() = laserAPDPNRatiosLaserEB[i].p2;
0103         vi.laserAPDPNRatios_p3() = laserAPDPNRatiosLaserEB[i].p3;
0104 
0105         vi.laserAPDPNref() = laserAPDPNRatiosRefEB[i];
0106         vi.laserAlpha() = laserAlphasEB[i];
0107 
0108         vi.linearCorrections_p1() = linearCorrectionsEB[i].p1;
0109         vi.linearCorrections_p2() = linearCorrectionsEB[i].p2;
0110         vi.linearCorrections_p3() = linearCorrectionsEB[i].p3;
0111 
0112         vi.timeCalibConstants() = timeCalibConstantsEB[i];
0113       }  // end Barrel loop
0114 
0115       // time maps
0116       for (unsigned int i = 0; i < laserAPDPNRatiosData.getTimeMap().size(); ++i) {
0117         auto vi = view[i];
0118         vi.laserAPDPNRatios_t1() = laserAPDPNRatiosTime[i].t1.value();
0119         vi.laserAPDPNRatios_t2() = laserAPDPNRatiosTime[i].t2.value();
0120         vi.laserAPDPNRatios_t3() = laserAPDPNRatiosTime[i].t3.value();
0121       }
0122 
0123       for (unsigned int i = 0; i < linearCorrectionsData.getTimeMap().size(); ++i) {
0124         auto vi = view[i];
0125         vi.linearCorrections_t1() = linearCorrectionsTime[i].t1.value();
0126         vi.linearCorrections_t2() = linearCorrectionsTime[i].t2.value();
0127         vi.linearCorrections_t3() = linearCorrectionsTime[i].t3.value();
0128       }
0129 
0130       // scalar data
0131       // ADC to GeV constants
0132       view.adcToGeVConstantEB() = adcToGeVConstantData.getEBValue();
0133 
0134       // time offset constants
0135       view.timeOffsetConstantEB() = timeOffsetConstantData.getEBValue();
0136 
0137       // endcap conditions
0138       if (!isPhase2_) {
0139         auto const& intercalibConstantsEE = intercalibConstantsData.endcapItems();
0140         auto const& channelStatusEE = channelStatusData.endcapItems();
0141         auto const& laserAPDPNRatiosLaserEE = laserAPDPNRatiosData.getLaserMap().endcapItems();
0142         auto const& laserAPDPNRatiosRefEE = laserAPDPNRatiosRefData.endcapItems();
0143         auto const& laserAlphasEE = laserAlphasData.endcapItems();
0144         auto const& linearCorrectionsEE = linearCorrectionsData.getValueMap().endcapItems();
0145         auto const& timeCalibConstantsEE = timeCalibConstantsData.endcapItems();
0146 
0147         const auto endcapSize = channelStatusData.endcapItems().size();
0148         for (unsigned int i = 0; i < endcapSize; ++i) {
0149           auto vi = view[barrelSize + i];
0150 
0151           vi.intercalibConstants() = intercalibConstantsEE[i];
0152           vi.channelStatus() = channelStatusEE[i].getEncodedStatusCode();
0153 
0154           vi.laserAPDPNRatios_p1() = laserAPDPNRatiosLaserEE[i].p1;
0155           vi.laserAPDPNRatios_p2() = laserAPDPNRatiosLaserEE[i].p2;
0156           vi.laserAPDPNRatios_p3() = laserAPDPNRatiosLaserEE[i].p3;
0157 
0158           vi.laserAPDPNref() = laserAPDPNRatiosRefEE[i];
0159           vi.laserAlpha() = laserAlphasEE[i];
0160 
0161           vi.linearCorrections_p1() = linearCorrectionsEE[i].p1;
0162           vi.linearCorrections_p2() = linearCorrectionsEE[i].p2;
0163           vi.linearCorrections_p3() = linearCorrectionsEE[i].p3;
0164 
0165           vi.timeCalibConstants() = timeCalibConstantsEE[i];
0166         }  // end Endcap loop
0167 
0168         // scalar data
0169         // ADC to GeV constants
0170         view.adcToGeVConstantEE() = adcToGeVConstantData.getEEValue();
0171 
0172         // time offset constants
0173         view.timeOffsetConstantEE() = timeOffsetConstantData.getEEValue();
0174       }
0175 
0176       // number of barrel items as offset for hashed ID access to EE items of columns
0177       view.offsetEE() = barrelSize;
0178 
0179       return product;
0180     }
0181 
0182   private:
0183     edm::ESGetToken<EcalADCToGeVConstant, EcalADCToGeVConstantRcd> adcToGeVConstantToken_;
0184     edm::ESGetToken<EcalIntercalibConstants, EcalIntercalibConstantsRcd> intercalibConstantsToken_;
0185     edm::ESGetToken<EcalChannelStatus, EcalChannelStatusRcd> channelStatusToken_;
0186     edm::ESGetToken<EcalLaserAPDPNRatios, EcalLaserAPDPNRatiosRcd> laserAPDPNRatiosToken_;
0187     edm::ESGetToken<EcalLaserAPDPNRatiosRef, EcalLaserAPDPNRatiosRefRcd> laserAPDPNRatiosRefToken_;
0188     edm::ESGetToken<EcalLaserAlphas, EcalLaserAlphasRcd> laserAlphasToken_;
0189     edm::ESGetToken<EcalLinearCorrections, EcalLinearCorrectionsRcd> linearCorrectionsToken_;
0190     edm::ESGetToken<EcalTimeCalibConstants, EcalTimeCalibConstantsRcd> timeCalibConstantsToken_;
0191     edm::ESGetToken<EcalTimeOffsetConstant, EcalTimeOffsetConstantRcd> timeOffsetConstantToken_;
0192 
0193     bool const isPhase2_;
0194   };
0195 
0196 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0197 
0198 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(EcalRecHitConditionsESProducer);