CMS 3D CMS Logo

EcalMultifitConditionsHostESProducer.cc
Go to the documentation of this file.
3 
4 #include <algorithm>
5 #include <cassert>
24 
28 
31 
39 
42  public:
44  auto cc = setWhatProduced(this);
45  pedestalsToken_ = cc.consumes();
46  gainRatiosToken_ = cc.consumes();
47  pulseShapesToken_ = cc.consumes();
48  pulseCovariancesToken_ = cc.consumes();
49  samplesCorrelationToken_ = cc.consumes();
50  timeBiasCorrectionsToken_ = cc.consumes();
51  timeCalibConstantsToken_ = cc.consumes();
52  sampleMaskToken_ = cc.consumes();
53  timeOffsetConstantToken_ = cc.consumes();
54  }
55 
56  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
58  descriptions.addWithDefaultLabel(desc);
59  }
60 
61  std::unique_ptr<EcalMultifitConditionsHost> produce(EcalMultifitConditionsRcd const& iRecord) {
62  auto const& pedestalsData = iRecord.get(pedestalsToken_);
63  auto const& gainRatiosData = iRecord.get(gainRatiosToken_);
64  auto const& pulseShapesData = iRecord.get(pulseShapesToken_);
65  auto const& pulseCovariancesData = iRecord.get(pulseCovariancesToken_);
66  auto const& samplesCorrelationData = iRecord.get(samplesCorrelationToken_);
67  auto const& timeBiasCorrectionsData = iRecord.get(timeBiasCorrectionsToken_);
68  auto const& timeCalibConstantsData = iRecord.get(timeCalibConstantsToken_);
69  auto const& sampleMaskData = iRecord.get(sampleMaskToken_);
70  auto const& timeOffsetConstantData = iRecord.get(timeOffsetConstantToken_);
71 
72  size_t numberOfXtals = pedestalsData.size();
73 
74  auto product = std::make_unique<EcalMultifitConditionsHost>(numberOfXtals, cms::alpakatools::host());
75  auto view = product->view();
76 
77  // Filling pedestals
78  const auto barrelSize = pedestalsData.barrelItems().size();
79  const auto endcapSize = pedestalsData.endcapItems().size();
80 
81  auto const& pedestalsEB = pedestalsData.barrelItems();
82  auto const& pedestalsEE = pedestalsData.endcapItems();
83  auto const& gainRatiosEB = gainRatiosData.barrelItems();
84  auto const& gainRatiosEE = gainRatiosData.endcapItems();
85  auto const& pulseShapesEB = pulseShapesData.barrelItems();
86  auto const& pulseShapesEE = pulseShapesData.endcapItems();
87  auto const& pulseCovariancesEB = pulseCovariancesData.barrelItems();
88  auto const& pulseCovariancesEE = pulseCovariancesData.endcapItems();
89  auto const& timeCalibConstantsEB = timeCalibConstantsData.barrelItems();
90  auto const& timeCalibConstantsEE = timeCalibConstantsData.endcapItems();
91 
92  for (unsigned int i = 0; i < barrelSize; ++i) {
93  auto vi = view[i];
94 
95  vi.pedestals_mean_x12() = pedestalsEB[i].mean_x12;
96  vi.pedestals_rms_x12() = pedestalsEB[i].rms_x12;
97  vi.pedestals_mean_x6() = pedestalsEB[i].mean_x6;
98  vi.pedestals_rms_x6() = pedestalsEB[i].rms_x6;
99  vi.pedestals_mean_x1() = pedestalsEB[i].mean_x1;
100  vi.pedestals_rms_x1() = pedestalsEB[i].rms_x1;
101 
102  vi.gain12Over6() = gainRatiosEB[i].gain12Over6();
103  vi.gain6Over1() = gainRatiosEB[i].gain6Over1();
104 
105  vi.timeCalibConstants() = timeCalibConstantsEB[i];
106 
107  std::memcpy(vi.pulseShapes().data(), pulseShapesEB[i].pdfval, sizeof(float) * EcalPulseShape::TEMPLATESAMPLES);
108  for (unsigned int j = 0; j < EcalPulseShape::TEMPLATESAMPLES; ++j) {
109  for (unsigned int k = 0; k < EcalPulseShape::TEMPLATESAMPLES; ++k) {
110  vi.pulseCovariance()(j, k) = pulseCovariancesEB[i].val(j, k);
111  }
112  }
113  } // end Barrel loop
114  for (unsigned int i = 0; i < endcapSize; ++i) {
115  auto vi = view[barrelSize + i];
116 
117  vi.pedestals_mean_x12() = pedestalsEE[i].mean_x12;
118  vi.pedestals_rms_x12() = pedestalsEE[i].rms_x12;
119  vi.pedestals_mean_x6() = pedestalsEE[i].mean_x6;
120  vi.pedestals_rms_x6() = pedestalsEE[i].rms_x6;
121  vi.pedestals_mean_x1() = pedestalsEE[i].mean_x1;
122  vi.pedestals_rms_x1() = pedestalsEE[i].rms_x1;
123 
124  vi.gain12Over6() = gainRatiosEE[i].gain12Over6();
125  vi.gain6Over1() = gainRatiosEE[i].gain6Over1();
126 
127  vi.timeCalibConstants() = timeCalibConstantsEE[i];
128 
129  std::memcpy(vi.pulseShapes().data(), pulseShapesEE[i].pdfval, sizeof(float) * EcalPulseShape::TEMPLATESAMPLES);
130 
131  for (unsigned int j = 0; j < EcalPulseShape::TEMPLATESAMPLES; ++j) {
132  for (unsigned int k = 0; k < EcalPulseShape::TEMPLATESAMPLES; ++k) {
133  vi.pulseCovariance()(j, k) = pulseCovariancesEE[i].val(j, k);
134  }
135  }
136  } // end Endcap loop
137 
138  // === Scalar data (not by xtal)
139  // TimeBiasCorrection
140  auto const timeCorrAmplBinsSizeEB = timeBiasCorrectionsData.EBTimeCorrAmplitudeBins.size();
141  auto const timeCorrShiftBinsSizeEB = timeBiasCorrectionsData.EBTimeCorrShiftBins.size();
142  // Assert that there are not more parameters than the EcalMultiFitConditionsSoA expects
143  assert(timeCorrAmplBinsSizeEB <= kMaxTimeBiasCorrectionBinsEB);
144  assert(timeCorrShiftBinsSizeEB <= kMaxTimeBiasCorrectionBinsEB);
145  std::memcpy(view.timeBiasCorrections_amplitude_EB().data(),
146  timeBiasCorrectionsData.EBTimeCorrAmplitudeBins.data(),
147  sizeof(float) * timeCorrAmplBinsSizeEB);
148  std::memcpy(view.timeBiasCorrections_shift_EB().data(),
149  timeBiasCorrectionsData.EBTimeCorrShiftBins.data(),
150  sizeof(float) * timeCorrShiftBinsSizeEB);
151 
152  auto const timeCorrAmplBinsSizeEE = timeBiasCorrectionsData.EETimeCorrAmplitudeBins.size();
153  auto const timeCorrShiftBinsSizeEE = timeBiasCorrectionsData.EETimeCorrShiftBins.size();
154  // Assert that there are not more parameters than the EcalMultiFitConditionsSoA expects
155  assert(timeCorrAmplBinsSizeEE <= kMaxTimeBiasCorrectionBinsEE);
156  assert(timeCorrShiftBinsSizeEE <= kMaxTimeBiasCorrectionBinsEE);
157  std::memcpy(view.timeBiasCorrections_amplitude_EE().data(),
158  timeBiasCorrectionsData.EETimeCorrAmplitudeBins.data(),
159  sizeof(float) * timeCorrAmplBinsSizeEE);
160  std::memcpy(view.timeBiasCorrections_shift_EE().data(),
161  timeBiasCorrectionsData.EETimeCorrShiftBins.data(),
162  sizeof(float) * timeCorrShiftBinsSizeEE);
163 
164  view.timeBiasCorrectionSizeEB() = std::min(timeCorrAmplBinsSizeEB, kMaxTimeBiasCorrectionBinsEB);
165  view.timeBiasCorrectionSizeEE() = std::min(timeCorrAmplBinsSizeEE, kMaxTimeBiasCorrectionBinsEE);
166 
167  // SampleCorrelation
168  std::memcpy(view.sampleCorrelation_EB_G12().data(),
169  samplesCorrelationData.EBG12SamplesCorrelation.data(),
170  sizeof(double) * ecalPh1::sampleSize);
171  std::memcpy(view.sampleCorrelation_EB_G6().data(),
172  samplesCorrelationData.EBG6SamplesCorrelation.data(),
173  sizeof(double) * ecalPh1::sampleSize);
174  std::memcpy(view.sampleCorrelation_EB_G1().data(),
175  samplesCorrelationData.EBG1SamplesCorrelation.data(),
176  sizeof(double) * ecalPh1::sampleSize);
177 
178  std::memcpy(view.sampleCorrelation_EE_G12().data(),
179  samplesCorrelationData.EEG12SamplesCorrelation.data(),
180  sizeof(double) * ecalPh1::sampleSize);
181  std::memcpy(view.sampleCorrelation_EE_G6().data(),
182  samplesCorrelationData.EEG6SamplesCorrelation.data(),
183  sizeof(double) * ecalPh1::sampleSize);
184  std::memcpy(view.sampleCorrelation_EE_G1().data(),
185  samplesCorrelationData.EEG1SamplesCorrelation.data(),
186  sizeof(double) * ecalPh1::sampleSize);
187 
188  // Sample masks
189  view.sampleMask_EB() = sampleMaskData.getEcalSampleMaskRecordEB();
190  view.sampleMask_EE() = sampleMaskData.getEcalSampleMaskRecordEE();
191 
192  // Time offsets
193  view.timeOffset_EB() = timeOffsetConstantData.getEBValue();
194  view.timeOffset_EE() = timeOffsetConstantData.getEEValue();
195 
196  // number of barrel items as offset for hashed ID access to EE items of columns
197  view.offsetEE() = barrelSize;
198 
199  return product;
200  }
201 
202  private:
212  };
213 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
214 
215 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(EcalMultifitConditionsHostESProducer);
edm::ESGetToken< EcalSamplesCorrelation, EcalSamplesCorrelationRcd > samplesCorrelationToken_
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:166
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
static const int TEMPLATESAMPLES
constexpr size_t kMaxTimeBiasCorrectionBinsEB
assert(be >=bs)
edm::ESGetToken< EcalTimeBiasCorrections, EcalTimeBiasCorrectionsRcd > timeBiasCorrectionsToken_
std::unique_ptr< EcalMultifitConditionsHost > produce(EcalMultifitConditionsRcd const &iRecord)
edm::ESGetToken< EcalTimeCalibConstants, EcalTimeCalibConstantsRcd > timeCalibConstantsToken_
constexpr size_t kMaxTimeBiasCorrectionBinsEE
edm::ESGetToken< EcalPulseCovariances, EcalPulseCovariancesRcd > pulseCovariancesToken_
alpaka::DevCpu const & host()
Definition: host.h:14
edm::ESGetToken< EcalTimeOffsetConstant, EcalTimeOffsetConstantRcd > timeOffsetConstantToken_
static constexpr unsigned int sampleSize
Definition: EcalConstants.h:53
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(type)
Definition: ModuleFactory.h:17
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const