Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-05 23:51:43

0001 // CMSSW imports
0002 #include "FWCore/Framework/interface/ESWatcher.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "FWCore/Utilities/interface/InputTag.h"
0009 #include "FWCore/Utilities/interface/StreamID.h"
0010 #include <iomanip>  // for std::setw
0011 #include <future>
0012 #include <chrono>
0013 
0014 // Alpaka imports
0015 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/stream/EDProducer.h"
0016 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/EDPutToken.h"
0017 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESGetToken.h"
0018 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/Event.h"
0019 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0020 #include "HeterogeneousCore/AlpakaInterface/interface/CopyToDevice.h"
0021 
0022 // includes for data formats
0023 #include "DataFormats/HGCalDigi/interface/HGCalDigiHost.h"
0024 #include "DataFormats/HGCalDigi/interface/alpaka/HGCalDigiDevice.h"
0025 #include "DataFormats/HGCalRecHit/interface/HGCalRecHitHost.h"
0026 #include "DataFormats/HGCalRecHit/interface/alpaka/HGCalRecHitDevice.h"
0027 
0028 // includes for size, calibration, and configuration parameters
0029 #include "CondFormats/DataRecord/interface/HGCalElectronicsMappingRcd.h"
0030 #include "CondFormats/DataRecord/interface/HGCalModuleConfigurationRcd.h"
0031 #include "CondFormats/HGCalObjects/interface/HGCalMappingModuleIndexer.h"
0032 #include "CondFormats/HGCalObjects/interface/HGCalCalibrationParameterHost.h"
0033 #include "CondFormats/HGCalObjects/interface/alpaka/HGCalCalibrationParameterDevice.h"
0034 #include "RecoLocalCalo/HGCalRecAlgos/interface/alpaka/HGCalRecHitCalibrationAlgorithms.h"
0035 
0036 // flag to assist the computational performance test
0037 // #define HGCAL_PERF_TEST
0038 
0039 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0040 
0041   using namespace cms::alpakatools;
0042 
0043   class HGCalRecHitsProducer : public stream::EDProducer<> {
0044   public:
0045     explicit HGCalRecHitsProducer(const edm::ParameterSet&);
0046     static void fillDescriptions(edm::ConfigurationDescriptions&);
0047 
0048   private:
0049     void produce(device::Event&, device::EventSetup const&) override;
0050     edm::ESWatcher<HGCalElectronicsMappingRcd> calibWatcher_;
0051     edm::ESWatcher<HGCalModuleConfigurationRcd> configWatcher_;
0052     const edm::EDGetTokenT<hgcaldigi::HGCalDigiHost> digisToken_;
0053     edm::ESGetToken<hgcalrechit::HGCalCalibParamHost, HGCalModuleConfigurationRcd> calibToken_;
0054     device::ESGetToken<hgcalrechit::HGCalConfigParamDevice, HGCalModuleConfigurationRcd> configToken_;
0055     const device::EDPutToken<hgcalrechit::HGCalRecHitDevice> recHitsToken_;
0056     const HGCalRecHitCalibrationAlgorithms calibrator_;
0057     const int n_hits_scale;
0058   };
0059 
0060   HGCalRecHitsProducer::HGCalRecHitsProducer(const edm::ParameterSet& iConfig)
0061       : EDProducer(iConfig),
0062         digisToken_{consumes<hgcaldigi::HGCalDigiHost>(iConfig.getParameter<edm::InputTag>("digis"))},
0063         calibToken_{esConsumes(iConfig.getParameter<edm::ESInputTag>("calibSource"))},
0064         configToken_{esConsumes(iConfig.getParameter<edm::ESInputTag>("configSource"))},
0065         recHitsToken_{produces()},
0066         calibrator_{iConfig.getParameter<int>("n_blocks"), iConfig.getParameter<int>("n_threads")},
0067         n_hits_scale{iConfig.getParameter<int>("n_hits_scale")} {
0068 #ifndef HGCAL_PERF_TEST
0069     if (n_hits_scale > 1) {
0070       throw cms::Exception("RuntimeError") << "Build with `HGCAL_PERF_TEST` flag to activate `n_hits_scale`.";
0071     }
0072 #endif
0073   }
0074 
0075   void HGCalRecHitsProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0076     edm::ParameterSetDescription desc;
0077     desc.add<edm::InputTag>("digis", edm::InputTag("hgcalDigis", "DIGI", "TEST"));
0078     desc.add("calibSource", edm::ESInputTag{})->setComment("Label for calibration parameters");
0079     desc.add("configSource", edm::ESInputTag{})->setComment("Label for ROC configuration parameters");
0080     desc.add<int>("n_blocks", -1);
0081     desc.add<int>("n_threads", -1);
0082     desc.add<int>("n_hits_scale", -1);
0083     descriptions.addWithDefaultLabel(desc);
0084   }
0085 
0086   void HGCalRecHitsProducer::produce(device::Event& iEvent, device::EventSetup const& iSetup) {
0087     auto& queue = iEvent.queue();
0088 
0089     // Read digis
0090     auto const& hostCalibParamProvider = iSetup.getData(calibToken_);
0091     auto const& deviceConfigParamProvider = iSetup.getData(configToken_);
0092     auto const& hostDigisIn = iEvent.get(digisToken_);
0093 
0094     //printout new conditions if available
0095     LogDebug("HGCalCalibrationParameter").log([&](auto& log) {
0096       if (calibWatcher_.check(iSetup)) {
0097         for (int i = 0; i < deviceConfigParamProvider.view().metadata().size(); i++) {
0098           log << "idx = " << i << ", "
0099               << "gain = " << deviceConfigParamProvider.view()[i].gain() << ", ";
0100         }
0101         for (int i = 0; i < hostCalibParamProvider.view().metadata().size(); i++) {
0102           log << "idx = " << i << ", "
0103               << "ADC_ped = " << hostCalibParamProvider.view()[i].ADC_ped() << ", "
0104               << "CM_slope = " << hostCalibParamProvider.view()[i].CM_slope() << ", "
0105               << "CM_ped = " << hostCalibParamProvider.view()[i].CM_ped() << ", "
0106               << "BXm1_slope = " << hostCalibParamProvider.view()[i].BXm1_slope() << ", ";
0107         }
0108       }
0109     });
0110 
0111 #ifdef HGCAL_PERF_TEST
0112     uint32_t oldSize = hostDigisIn.view().metadata().size();
0113     uint32_t newSize = oldSize * (n_hits_scale > 0 ? (unsigned)n_hits_scale : 1);
0114     auto hostDigis = HGCalDigiHost(newSize, queue);
0115     auto hostCalibParam = HGCalCalibParamHost(newSize, queue);
0116     // TODO: replace with memcp ?
0117     for (uint32_t i = 0; i < newSize; i++) {
0118       hostDigis.view()[i].tctp() = hostDigisIn.view()[i % oldSize].tctp();
0119       hostDigis.view()[i].adcm1() = hostDigisIn.view()[i % oldSize].adcm1();
0120       hostDigis.view()[i].adc() = hostDigisIn.view()[i % oldSize].adc();
0121       hostDigis.view()[i].tot() = hostDigisIn.view()[i % oldSize].tot();
0122       hostDigis.view()[i].toa() = hostDigisIn.view()[i % oldSize].toa();
0123       hostDigis.view()[i].cm() = hostDigisIn.view()[i % oldSize].cm();
0124       hostDigis.view()[i].flags() = hostDigisIn.view()[i % oldSize].flags();
0125 
0126       hostCalibParam.view()[i].ADC_ped() = hostCalibParamProvider.view()[i % oldSize].ADC_ped();
0127       hostCalibParam.view()[i].Noise() = hostCalibParamProvider.view()[i % oldSize].Noise();
0128       hostCalibParam.view()[i].CM_slope() = hostCalibParamProvider.view()[i % oldSize].CM_slope();
0129       hostCalibParam.view()[i].CM_ped() = hostCalibParamProvider.view()[i % oldSize].CM_ped();
0130       hostCalibParam.view()[i].BXm1_slope() = hostCalibParamProvider.view()[i % oldSize].BXm1_slope();
0131       hostCalibParam.view()[i].TOTtoADC() = hostCalibParamProvider.view()[i % oldSize].TOTtoADC();
0132       hostCalibParam.view()[i].TOT_ped() = hostCalibParamProvider.view()[i % oldSize].TOT_ped();
0133       hostCalibParam.view()[i].TOT_lin() = hostCalibParamProvider.view()[i % oldSize].TOT_lin();
0134       hostCalibParam.view()[i].TOT_P0() = hostCalibParamProvider.view()[i % oldSize].TOT_P0();
0135       hostCalibParam.view()[i].TOT_P1() = hostCalibParamProvider.view()[i % oldSize].TOT_P1();
0136       hostCalibParam.view()[i].TOT_P2() = hostCalibParamProvider.view()[i % oldSize].TOT_P2();
0137       hostCalibParam.view()[i].TOAtops() = hostCalibParamProvider.view()[i % oldSize].TOAtops();
0138       hostCalibParam.view()[i].MIPS_scale() = hostCalibParamProvider.view()[i % oldSize].MIPS_scale();
0139       hostCalibParam.view()[i].valid() = hostCalibParamProvider.view()[i % oldSize].valid();
0140     }
0141 #else
0142     const auto& hostDigis = hostDigisIn;
0143     const auto& hostCalibParam = hostCalibParamProvider;
0144 #endif
0145 
0146     LogDebug("HGCalRecHitsProducer") << "Loaded host digis: " << hostDigis.view().metadata().size();  //<< std::endl;
0147 
0148     LogDebug("HGCalRecHitsProducer") << "\n\nINFO -- calling calibrate method";  //<< std::endl;
0149 
0150 #ifdef EDM_ML_DEBUG
0151     alpaka::wait(queue);
0152     auto start = std::chrono::steady_clock::now();
0153 #endif
0154 
0155     LogDebug("HGCalRecHitsProducer") << "\n\nINFO -- Copying the calib to the device\n\n" << std::endl;
0156     HGCalCalibParamDevice deviceCalibParam(hostCalibParam.view().metadata().size(), queue);
0157     alpaka::memcpy(queue, deviceCalibParam.buffer(), hostCalibParam.const_buffer());
0158 
0159 #ifdef HGCAL_PERF_TEST
0160     auto tmpRecHits = calibrator_.calibrate(queue, hostDigis, deviceCalibParam, deviceConfigParamProvider);
0161     HGCalRecHitDevice recHits(oldSize, queue);
0162     alpaka::memcpy(queue, recHits.buffer(), tmpRecHits.const_buffer(), oldSize);
0163 #else
0164     auto recHits = calibrator_.calibrate(queue, hostDigis, deviceCalibParam, deviceConfigParamProvider);
0165 #endif
0166 
0167 #ifdef EDM_ML_DEBUG
0168     alpaka::wait(queue);
0169     auto stop = std::chrono::steady_clock::now();
0170     std::chrono::duration<float> elapsed = stop - start;
0171     LogDebug("HGCalRecHitsProducer") << "Time spent calibrating " << hostDigis.view().metadata().size()
0172                                      << " digis: " << elapsed.count();  //<< std::endl;
0173 #endif
0174 
0175     LogDebug("HGCalRecHitsProducer") << "\n\nINFO -- storing rec hits in the event";  //<< std::endl;
0176     iEvent.emplace(recHitsToken_, std::move(recHits));
0177   }
0178 
0179 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0180 
0181 // define this as a plug-in
0182 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0183 DEFINE_FWK_ALPAKA_MODULE(HGCalRecHitsProducer);