File indexing completed on 2024-04-06 12:25:41
0001 #include <cassert>
0002
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004
0005 #include "CondFormats/EcalObjects/interface/alpaka/EcalMultifitParametersDevice.h"
0006 #include "CondFormats/EcalObjects/interface/EcalMultifitParametersSoA.h"
0007 #include "CondFormats/DataRecord/interface/EcalMultifitParametersRcd.h"
0008
0009 #include "DataFormats/EcalDigi/interface/EcalConstants.h"
0010
0011 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESGetToken.h"
0012 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0013 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0014 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0015 #include "HeterogeneousCore/AlpakaInterface/interface/host.h"
0016 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0017
0018 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0019 class EcalMultifitParametersHostESProducer : public ESProducer {
0020 public:
0021 EcalMultifitParametersHostESProducer(edm::ParameterSet const&);
0022 ~EcalMultifitParametersHostESProducer() override = default;
0023
0024 static void fillDescriptions(edm::ConfigurationDescriptions&);
0025 std::unique_ptr<EcalMultifitParametersHost> produce(EcalMultifitParametersRcd const&);
0026
0027 private:
0028 std::vector<float> ebTimeFitParameters_;
0029 std::vector<float> eeTimeFitParameters_;
0030 std::vector<float> ebAmplitudeFitParameters_;
0031 std::vector<float> eeAmplitudeFitParameters_;
0032 };
0033
0034 EcalMultifitParametersHostESProducer::EcalMultifitParametersHostESProducer(edm::ParameterSet const& iConfig)
0035 : ESProducer(iConfig) {
0036 setWhatProduced(this);
0037
0038 auto const ebTimeFitParamsFromPSet = iConfig.getParameter<std::vector<double>>("EBtimeFitParameters");
0039 auto const eeTimeFitParamsFromPSet = iConfig.getParameter<std::vector<double>>("EEtimeFitParameters");
0040
0041 assert(ebTimeFitParamsFromPSet.size() == kNTimeFitParams);
0042 assert(eeTimeFitParamsFromPSet.size() == kNTimeFitParams);
0043 ebTimeFitParameters_.assign(ebTimeFitParamsFromPSet.begin(), ebTimeFitParamsFromPSet.end());
0044 eeTimeFitParameters_.assign(eeTimeFitParamsFromPSet.begin(), eeTimeFitParamsFromPSet.end());
0045
0046 auto const ebAmplFitParamsFromPSet = iConfig.getParameter<std::vector<double>>("EBamplitudeFitParameters");
0047 auto const eeAmplFitParamsFromPSet = iConfig.getParameter<std::vector<double>>("EEamplitudeFitParameters");
0048
0049 assert(ebAmplFitParamsFromPSet.size() == kNAmplitudeFitParams);
0050 assert(eeAmplFitParamsFromPSet.size() == kNAmplitudeFitParams);
0051 ebAmplitudeFitParameters_.assign(ebAmplFitParamsFromPSet.begin(), ebAmplFitParamsFromPSet.end());
0052 eeAmplitudeFitParameters_.assign(eeAmplFitParamsFromPSet.begin(), eeAmplFitParamsFromPSet.end());
0053 }
0054
0055 void EcalMultifitParametersHostESProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0056 edm::ParameterSetDescription desc;
0057 desc.add<std::vector<double>>("EBtimeFitParameters",
0058 {-2.015452e+00,
0059 3.130702e+00,
0060 -1.234730e+01,
0061 4.188921e+01,
0062 -8.283944e+01,
0063 9.101147e+01,
0064 -5.035761e+01,
0065 1.105621e+01});
0066 desc.add<std::vector<double>>("EEtimeFitParameters",
0067 {-2.390548e+00,
0068 3.553628e+00,
0069 -1.762341e+01,
0070 6.767538e+01,
0071 -1.332130e+02,
0072 1.407432e+02,
0073 -7.541106e+01,
0074 1.620277e+01});
0075 desc.add<std::vector<double>>("EBamplitudeFitParameters", {1.138, 1.652});
0076 desc.add<std::vector<double>>("EEamplitudeFitParameters", {1.890, 1.400});
0077 descriptions.addWithDefaultLabel(desc);
0078 }
0079
0080 std::unique_ptr<EcalMultifitParametersHost> EcalMultifitParametersHostESProducer::produce(
0081 EcalMultifitParametersRcd const& iRecord) {
0082 size_t const sizeone = 1;
0083 auto product = std::make_unique<EcalMultifitParametersHost>(sizeone, cms::alpakatools::host());
0084 auto view = product->view();
0085
0086 std::memcpy(view.timeFitParamsEB().data(), ebTimeFitParameters_.data(), sizeof(float) * kNTimeFitParams);
0087 std::memcpy(view.timeFitParamsEE().data(), eeTimeFitParameters_.data(), sizeof(float) * kNTimeFitParams);
0088
0089 std::memcpy(
0090 view.amplitudeFitParamsEB().data(), ebAmplitudeFitParameters_.data(), sizeof(float) * kNAmplitudeFitParams);
0091 std::memcpy(
0092 view.amplitudeFitParamsEE().data(), eeAmplitudeFitParameters_.data(), sizeof(float) * kNAmplitudeFitParams);
0093
0094 return product;
0095 }
0096
0097 }
0098
0099 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(EcalMultifitParametersHostESProducer);