CMS 3D CMS Logo

EcalMultifitParametersHostESProducer.cc
Go to the documentation of this file.
1 #include <cassert>
2 
4 
8 
10 
17 
20  public:
22  ~EcalMultifitParametersHostESProducer() override = default;
23 
25  std::unique_ptr<EcalMultifitParametersHost> produce(EcalMultifitParametersRcd const&);
26 
27  private:
28  std::vector<float> ebTimeFitParameters_;
29  std::vector<float> eeTimeFitParameters_;
30  std::vector<float> ebAmplitudeFitParameters_;
31  std::vector<float> eeAmplitudeFitParameters_;
32  };
33 
35  : ESProducer(iConfig) {
36  setWhatProduced(this);
37 
38  auto const ebTimeFitParamsFromPSet = iConfig.getParameter<std::vector<double>>("EBtimeFitParameters");
39  auto const eeTimeFitParamsFromPSet = iConfig.getParameter<std::vector<double>>("EEtimeFitParameters");
40  // Assert that there are as many parameters as the EcalMultiFitParametersSoA expects
41  assert(ebTimeFitParamsFromPSet.size() == kNTimeFitParams);
42  assert(eeTimeFitParamsFromPSet.size() == kNTimeFitParams);
43  ebTimeFitParameters_.assign(ebTimeFitParamsFromPSet.begin(), ebTimeFitParamsFromPSet.end());
44  eeTimeFitParameters_.assign(eeTimeFitParamsFromPSet.begin(), eeTimeFitParamsFromPSet.end());
45 
46  auto const ebAmplFitParamsFromPSet = iConfig.getParameter<std::vector<double>>("EBamplitudeFitParameters");
47  auto const eeAmplFitParamsFromPSet = iConfig.getParameter<std::vector<double>>("EEamplitudeFitParameters");
48  // Assert that there are as many parameters as the EcalMultiFitParametersSoA expects
49  assert(ebAmplFitParamsFromPSet.size() == kNAmplitudeFitParams);
50  assert(eeAmplFitParamsFromPSet.size() == kNAmplitudeFitParams);
51  ebAmplitudeFitParameters_.assign(ebAmplFitParamsFromPSet.begin(), ebAmplFitParamsFromPSet.end());
52  eeAmplitudeFitParameters_.assign(eeAmplFitParamsFromPSet.begin(), eeAmplFitParamsFromPSet.end());
53  }
54 
57  desc.add<std::vector<double>>("EBtimeFitParameters",
58  {-2.015452e+00,
59  3.130702e+00,
60  -1.234730e+01,
61  4.188921e+01,
62  -8.283944e+01,
63  9.101147e+01,
64  -5.035761e+01,
65  1.105621e+01});
66  desc.add<std::vector<double>>("EEtimeFitParameters",
67  {-2.390548e+00,
68  3.553628e+00,
69  -1.762341e+01,
70  6.767538e+01,
71  -1.332130e+02,
72  1.407432e+02,
73  -7.541106e+01,
74  1.620277e+01});
75  desc.add<std::vector<double>>("EBamplitudeFitParameters", {1.138, 1.652});
76  desc.add<std::vector<double>>("EEamplitudeFitParameters", {1.890, 1.400});
77  descriptions.addWithDefaultLabel(desc);
78  }
79 
80  std::unique_ptr<EcalMultifitParametersHost> EcalMultifitParametersHostESProducer::produce(
81  EcalMultifitParametersRcd const& iRecord) {
82  size_t const sizeone = 1;
83  auto product = std::make_unique<EcalMultifitParametersHost>(sizeone, cms::alpakatools::host());
84  auto view = product->view();
85 
86  std::memcpy(view.timeFitParamsEB().data(), ebTimeFitParameters_.data(), sizeof(float) * kNTimeFitParams);
87  std::memcpy(view.timeFitParamsEE().data(), eeTimeFitParameters_.data(), sizeof(float) * kNTimeFitParams);
88 
89  std::memcpy(
90  view.amplitudeFitParamsEB().data(), ebAmplitudeFitParameters_.data(), sizeof(float) * kNAmplitudeFitParams);
91  std::memcpy(
92  view.amplitudeFitParamsEE().data(), eeAmplitudeFitParameters_.data(), sizeof(float) * kNAmplitudeFitParams);
93 
94  return product;
95  }
96 
97 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
98 
99 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(EcalMultifitParametersHostESProducer);
constexpr size_t kNTimeFitParams
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:166
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
assert(be >=bs)
constexpr size_t kNAmplitudeFitParams
alpaka::DevCpu const & host()
Definition: host.h:14
std::unique_ptr< EcalMultifitParametersHost > produce(EcalMultifitParametersRcd const &)
#define DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(type)
Definition: ModuleFactory.h:17