File indexing completed on 2024-04-06 12:27:34
0001 #include <memory>
0002
0003 #include "CondFormats/EcalObjects/interface/EcalPFRecHitThresholds.h"
0004 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0007 #include "RecoParticleFlow/PFRecHitProducer/interface/PFRecHitParamsHostCollection.h"
0008 #include "RecoParticleFlow/PFRecHitProducer/interface/PFRecHitParamsRecord.h"
0009 #include "CalorimeterDefinitions.h"
0010
0011 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0012 using namespace particleFlowRecHitProducer;
0013
0014 class PFRecHitECALParamsESProducer : public ESProducer {
0015 public:
0016 PFRecHitECALParamsESProducer(edm::ParameterSet const& iConfig)
0017 : ESProducer(iConfig), cleaningThreshold_(iConfig.getParameter<double>("cleaningThreshold")) {
0018 auto cc = setWhatProduced(this);
0019 thresholdsToken_ = cc.consumes();
0020 }
0021
0022 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0023 edm::ParameterSetDescription desc;
0024 desc.add<double>("cleaningThreshold", 2);
0025 descriptions.addWithDefaultLabel(desc);
0026 }
0027
0028 std::unique_ptr<reco::PFRecHitECALParamsHostCollection> produce(const EcalPFRecHitThresholdsRcd& iRecord) {
0029 const auto& thresholds = iRecord.get(thresholdsToken_);
0030 auto product = std::make_unique<reco::PFRecHitECALParamsHostCollection>(ECAL::kSize, cms::alpakatools::host());
0031 for (uint32_t denseId = 0; denseId < ECAL::Barrel::kSize; denseId++)
0032 product->view().energyThresholds()[denseId] = thresholds.barrel(denseId);
0033 for (uint32_t denseId = 0; denseId < ECAL::Endcap::kSize; denseId++)
0034 product->view().energyThresholds()[denseId + ECAL::Barrel::kSize] = thresholds.endcap(denseId);
0035 product->view().cleaningThreshold() = cleaningThreshold_;
0036 return product;
0037 }
0038
0039 private:
0040 const double cleaningThreshold_;
0041 edm::ESGetToken<EcalPFRecHitThresholds, EcalPFRecHitThresholdsRcd> thresholdsToken_;
0042 };
0043
0044 }
0045
0046 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0047 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(PFRecHitECALParamsESProducer);