Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-03 04:18:11

0001 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0002 #include "FWCore/Utilities/interface/ESGetToken.h"
0003 
0004 #include "CondFormats/DataRecord/interface/HcalRecoParamsRcd.h"
0005 #include "CondFormats/DataRecord/interface/HcalPedestalsRcd.h"
0006 #include "CondFormats/DataRecord/interface/HcalGainsRcd.h"
0007 #include "CondFormats/DataRecord/interface/HcalLUTCorrsRcd.h"
0008 #include "CondFormats/DataRecord/interface/HcalRespCorrsRcd.h"
0009 #include "CondFormats/DataRecord/interface/HcalTimeCorrsRcd.h"
0010 #include "CondFormats/DataRecord/interface/HcalPedestalWidthsRcd.h"
0011 #include "CondFormats/DataRecord/interface/HcalGainWidthsRcd.h"
0012 #include "CondFormats/DataRecord/interface/HcalChannelQualityRcd.h"
0013 #include "CondFormats/DataRecord/interface/HcalQIETypesRcd.h"
0014 #include "CondFormats/DataRecord/interface/HcalQIEDataRcd.h"
0015 #include "CondFormats/DataRecord/interface/HcalSiPMParametersRcd.h"
0016 #include "CondFormats/HcalObjects/interface/HcalRecoParams.h"
0017 #include "CondFormats/HcalObjects/interface/HcalPedestals.h"
0018 #include "CondFormats/HcalObjects/interface/HcalGains.h"
0019 #include "CondFormats/HcalObjects/interface/HcalLUTCorrs.h"
0020 #include "CondFormats/HcalObjects/interface/HcalRespCorrs.h"
0021 #include "CondFormats/HcalObjects/interface/HcalTimeCorrs.h"
0022 #include "CondFormats/HcalObjects/interface/HcalPedestalWidths.h"
0023 #include "CondFormats/HcalObjects/interface/HcalGainWidths.h"
0024 #include "CondFormats/HcalObjects/interface/HcalChannelQuality.h"
0025 #include "CondFormats/HcalObjects/interface/HcalQIETypes.h"
0026 #include "CondFormats/HcalObjects/interface/HcalQIEData.h"
0027 #include "CondFormats/HcalObjects/interface/HcalSiPMParameters.h"
0028 
0029 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0030 #include "Geometry/HcalCommonData/interface/HcalDDDRecConstants.h"
0031 
0032 #include "CondFormats/HcalObjects/interface/alpaka/HcalMahiConditionsDevice.h"
0033 #include "CondFormats/HcalObjects/interface/HcalMahiConditionsSoA.h"
0034 #include "CondFormats/DataRecord/interface/HcalMahiConditionsRcd.h"
0035 
0036 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0037 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0038 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0039 #include "HeterogeneousCore/AlpakaInterface/interface/host.h"
0040 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0041 
0042 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0043   namespace {
0044     float convertPedWidths(
0045         float const value, float const width, int const i, HcalQIECoder const& coder, HcalQIEShape const& shape) {
0046       float const y = value;
0047       float const x = width;
0048       unsigned const x1 = static_cast<unsigned>(std::floor(y));
0049       unsigned const x2 = static_cast<unsigned>(std::floor(y + 1.));
0050       unsigned iun = static_cast<unsigned>(i);
0051       float const y1 = coder.charge(shape, x1, iun);
0052       float const y2 = coder.charge(shape, x2, iun);
0053       return (y2 - y1) * x;
0054     }
0055     float convertPed(float const x, int const i, HcalQIECoder const& coder, HcalQIEShape const& shape) {
0056       int const x1 = static_cast<int>(std::floor(x));
0057       int const x2 = static_cast<int>(std::floor(x + 1));
0058       float const y2 = coder.charge(shape, x2, i);
0059       float const y1 = coder.charge(shape, x1, i);
0060       return (y2 - y1) * (x - x1) + y1;
0061     }
0062   }  // namespace
0063 
0064   class HcalMahiConditionsESProducer : public ESProducer {
0065   public:
0066     HcalMahiConditionsESProducer(edm::ParameterSet const& iConfig) : ESProducer(iConfig) {
0067       auto cc = setWhatProduced(this);
0068       recoParamsToken_ = cc.consumes();
0069       pedestalsToken_ = cc.consumes();
0070       effectivePedestalsToken_ = cc.consumes(edm::ESInputTag{"", "withTopoEff"});
0071       gainsToken_ = cc.consumes();
0072       lutCorrsToken_ = cc.consumes();
0073       respCorrsToken_ = cc.consumes();
0074       timeCorrsToken_ = cc.consumes();
0075       pedestalWidthsToken_ = cc.consumes();
0076       effectivePedestalWidthsToken_ = cc.consumes(edm::ESInputTag{"", "withTopoEff"});
0077       gainWidthsToken_ = cc.consumes();
0078       channelQualityToken_ = cc.consumes();
0079       qieTypesToken_ = cc.consumes();
0080       qieDataToken_ = cc.consumes();
0081       sipmParametersToken_ = cc.consumes();
0082       topologyToken_ = cc.consumes();
0083       recConstantsToken_ = cc.consumes();
0084     }
0085 
0086     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0087       edm::ParameterSetDescription desc;
0088       descriptions.addWithDefaultLabel(desc);
0089     }
0090 
0091     std::unique_ptr<hcal::HcalMahiConditionsPortableHost> produce(HcalMahiConditionsRcd const& iRecord) {
0092       auto const& recoParams = iRecord.get(recoParamsToken_);
0093       auto const& pedestals = iRecord.get(pedestalsToken_);
0094       auto const& effectivePedestals = iRecord.get(effectivePedestalsToken_);
0095       auto const& gains = iRecord.get(gainsToken_);
0096       auto const& lutCorrs = iRecord.get(lutCorrsToken_);
0097       auto const& respCorrs = iRecord.get(respCorrsToken_);
0098       auto const& timeCorrs = iRecord.get(timeCorrsToken_);
0099       auto const& pedestalWidths = iRecord.get(pedestalWidthsToken_);
0100       auto const& effectivePedestalWidths = iRecord.get(effectivePedestalWidthsToken_);
0101       auto const& gainWidths = iRecord.get(gainWidthsToken_);
0102       auto const& channelQuality = iRecord.get(channelQualityToken_);
0103       auto const& qieTypes = iRecord.get(qieTypesToken_);
0104       auto const& qieData = iRecord.get(qieDataToken_);
0105       auto const& sipmParameters = iRecord.get(sipmParametersToken_);
0106       auto const& topology = iRecord.get(topologyToken_);
0107       auto const& recConstants = iRecord.get(recConstantsToken_);
0108 
0109       size_t const totalChannels =
0110           pedestals.getAllContainers()[0].second.size() + pedestals.getAllContainers()[1].second.size();
0111 
0112       auto product = std::make_unique<hcal::HcalMahiConditionsPortableHost>(totalChannels, cms::alpakatools::host());
0113 
0114       auto view = product->view();
0115 
0116       // convert pedestals
0117       auto const ped_unitIsADC = pedestals.isADC();
0118       auto const effPed_unitIsADC = effectivePedestals.isADC();
0119 
0120       // fill HB channels
0121       auto const recoParams_containers = recoParams.getAllContainers();
0122       auto const pedestals_containers = pedestals.getAllContainers();
0123       auto const effectivePedestals_containers = effectivePedestals.getAllContainers();
0124       auto const gains_containers = gains.getAllContainers();
0125       auto const lutCorrs_containers = lutCorrs.getAllContainers();
0126       auto const respCorrs_containers = respCorrs.getAllContainers();
0127       auto const timeCorrs_containers = timeCorrs.getAllContainers();
0128       auto const pedestalWidths_containers = pedestalWidths.getAllContainers();
0129       auto const effectivePedestalWidths_containers = effectivePedestalWidths.getAllContainers();
0130       auto const gainWidths_containers = gainWidths.getAllContainers();
0131       auto const channelQuality_containers = channelQuality.getAllContainers();
0132       auto const qieTypes_containers = qieTypes.getAllContainers();
0133       auto const qieData_containers = qieData.getAllContainers();
0134       auto const sipmParameters_containers = sipmParameters.getAllContainers();
0135 
0136       auto const& pedestals_barrel = pedestals_containers[0].second;
0137       auto const& effectivePedestals_barrel = effectivePedestals_containers[0].second;
0138       auto const& gains_barrel = gains_containers[0].second;
0139       auto const& lutCorrs_barrel = lutCorrs_containers[0].second;
0140       auto const& respCorrs_barrel = respCorrs_containers[0].second;
0141       auto const& timeCorrs_barrel = timeCorrs_containers[0].second;
0142       auto const& pedestalWidths_barrel = pedestalWidths_containers[0].second;
0143       auto const& effectivePedestalWidths_barrel = effectivePedestalWidths_containers[0].second;
0144       auto const& gainWidths_barrel = gainWidths_containers[0].second;
0145       auto const& channelQuality_barrel = channelQuality_containers[0].second;
0146       auto const& qieTypes_barrel = qieTypes_containers[0].second;
0147       auto const& qieData_barrel = qieData_containers[0].second;
0148       auto const& sipmParameters_barrel = sipmParameters_containers[0].second;
0149 
0150       for (uint64_t i = 0; i < pedestals_barrel.size(); ++i) {
0151         auto vi = view[i];
0152 
0153         // convert pedestals
0154         auto const& qieCoder = qieData_barrel[i];
0155         auto const qieType = qieTypes_barrel[i].getValue() > 1 ? 1 : 0;
0156         auto const& qieShape = qieData.getShape(qieType);
0157 
0158         // covert pedestal values if unit is ADC
0159         vi.pedestals_value()[0] = ped_unitIsADC ? convertPed(pedestals_barrel[i].getValue(0), 0, qieCoder, qieShape)
0160                                                 : pedestals_barrel[i].getValue(0);
0161         vi.pedestals_value()[1] = ped_unitIsADC ? convertPed(pedestals_barrel[i].getValue(1), 1, qieCoder, qieShape)
0162                                                 : pedestals_barrel[i].getValue(1);
0163         vi.pedestals_value()[2] = ped_unitIsADC ? convertPed(pedestals_barrel[i].getValue(2), 2, qieCoder, qieShape)
0164                                                 : pedestals_barrel[i].getValue(2);
0165         vi.pedestals_value()[3] = ped_unitIsADC ? convertPed(pedestals_barrel[i].getValue(3), 3, qieCoder, qieShape)
0166                                                 : pedestals_barrel[i].getValue(3);
0167 
0168         vi.pedestals_width()[0] =
0169             ped_unitIsADC
0170                 ? convertPedWidths(
0171                       pedestals_barrel[i].getValue(0), pedestalWidths_barrel[i].getWidth(0), 0, qieCoder, qieShape)
0172                 : pedestalWidths_barrel[i].getWidth(0);
0173         vi.pedestals_width()[1] =
0174             ped_unitIsADC
0175                 ? convertPedWidths(
0176                       pedestals_barrel[i].getValue(1), pedestalWidths_barrel[i].getWidth(1), 1, qieCoder, qieShape)
0177                 : pedestalWidths_barrel[i].getWidth(1);
0178         vi.pedestals_width()[2] =
0179             ped_unitIsADC
0180                 ? convertPedWidths(
0181                       pedestals_barrel[i].getValue(2), pedestalWidths_barrel[i].getWidth(2), 2, qieCoder, qieShape)
0182                 : pedestalWidths_barrel[i].getWidth(2);
0183         vi.pedestals_width()[3] =
0184             ped_unitIsADC
0185                 ? convertPedWidths(
0186                       pedestals_barrel[i].getValue(3), pedestalWidths_barrel[i].getWidth(3), 3, qieCoder, qieShape)
0187                 : pedestalWidths_barrel[i].getWidth(3);
0188 
0189         vi.effectivePedestals()[0] = effPed_unitIsADC
0190                                          ? convertPed(effectivePedestals_barrel[i].getValue(0), 0, qieCoder, qieShape)
0191                                          : effectivePedestals_barrel[i].getValue(0);
0192         vi.effectivePedestals()[1] = effPed_unitIsADC
0193                                          ? convertPed(effectivePedestals_barrel[i].getValue(1), 1, qieCoder, qieShape)
0194                                          : effectivePedestals_barrel[i].getValue(1);
0195         vi.effectivePedestals()[2] = effPed_unitIsADC
0196                                          ? convertPed(effectivePedestals_barrel[i].getValue(2), 2, qieCoder, qieShape)
0197                                          : effectivePedestals_barrel[i].getValue(2);
0198         vi.effectivePedestals()[3] = effPed_unitIsADC
0199                                          ? convertPed(effectivePedestals_barrel[i].getValue(3), 3, qieCoder, qieShape)
0200                                          : effectivePedestals_barrel[i].getValue(3);
0201 
0202         vi.effectivePedestalWidths()[0] = effPed_unitIsADC
0203                                               ? convertPedWidths(effectivePedestals_barrel[i].getValue(0),
0204                                                                  effectivePedestalWidths_barrel[i].getWidth(0),
0205                                                                  0,
0206                                                                  qieCoder,
0207                                                                  qieShape)
0208                                               : effectivePedestalWidths_barrel[i].getWidth(0);
0209         vi.effectivePedestalWidths()[1] = effPed_unitIsADC
0210                                               ? convertPedWidths(effectivePedestals_barrel[i].getValue(1),
0211                                                                  effectivePedestalWidths_barrel[i].getWidth(1),
0212                                                                  1,
0213                                                                  qieCoder,
0214                                                                  qieShape)
0215                                               : effectivePedestalWidths_barrel[i].getWidth(1);
0216         vi.effectivePedestalWidths()[2] = effPed_unitIsADC
0217                                               ? convertPedWidths(effectivePedestals_barrel[i].getValue(2),
0218                                                                  effectivePedestalWidths_barrel[i].getWidth(2),
0219                                                                  2,
0220                                                                  qieCoder,
0221                                                                  qieShape)
0222                                               : effectivePedestalWidths_barrel[i].getWidth(2);
0223         vi.effectivePedestalWidths()[3] = effPed_unitIsADC
0224                                               ? convertPedWidths(effectivePedestals_barrel[i].getValue(3),
0225                                                                  effectivePedestalWidths_barrel[i].getWidth(3),
0226                                                                  3,
0227                                                                  qieCoder,
0228                                                                  qieShape)
0229                                               : effectivePedestalWidths_barrel[i].getWidth(3);
0230 
0231         vi.gains_value()[0] = gains_barrel[i].getValue(0);
0232         vi.gains_value()[1] = gains_barrel[i].getValue(1);
0233         vi.gains_value()[2] = gains_barrel[i].getValue(2);
0234         vi.gains_value()[3] = gains_barrel[i].getValue(3);
0235 
0236         vi.lutCorrs_values() = lutCorrs_barrel[i].getValue();
0237         vi.respCorrs_values() = respCorrs_barrel[i].getValue();
0238         vi.timeCorrs_values() = timeCorrs_barrel[i].getValue();
0239 
0240         vi.pedestalWidths_sigma00() = *(pedestalWidths_barrel[i].getValues());
0241         vi.pedestalWidths_sigma01() = *(pedestalWidths_barrel[i].getValues() + 1);
0242         vi.pedestalWidths_sigma02() = *(pedestalWidths_barrel[i].getValues() + 2);
0243         vi.pedestalWidths_sigma03() = *(pedestalWidths_barrel[i].getValues() + 3);
0244         vi.pedestalWidths_sigma10() = *(pedestalWidths_barrel[i].getValues() + 4);
0245         vi.pedestalWidths_sigma11() = *(pedestalWidths_barrel[i].getValues() + 5);
0246         vi.pedestalWidths_sigma12() = *(pedestalWidths_barrel[i].getValues() + 6);
0247         vi.pedestalWidths_sigma13() = *(pedestalWidths_barrel[i].getValues() + 7);
0248         vi.pedestalWidths_sigma20() = *(pedestalWidths_barrel[i].getValues() + 8);
0249         vi.pedestalWidths_sigma21() = *(pedestalWidths_barrel[i].getValues() + 9);
0250         vi.pedestalWidths_sigma22() = *(pedestalWidths_barrel[i].getValues() + 10);
0251         vi.pedestalWidths_sigma23() = *(pedestalWidths_barrel[i].getValues() + 11);
0252         vi.pedestalWidths_sigma30() = *(pedestalWidths_barrel[i].getValues() + 12);
0253         vi.pedestalWidths_sigma31() = *(pedestalWidths_barrel[i].getValues() + 13);
0254         vi.pedestalWidths_sigma32() = *(pedestalWidths_barrel[i].getValues() + 14);
0255         vi.pedestalWidths_sigma33() = *(pedestalWidths_barrel[i].getValues() + 15);
0256 
0257         vi.gainWidths_value0() = gainWidths_barrel[i].getValue(0);
0258         vi.gainWidths_value1() = gainWidths_barrel[i].getValue(1);
0259         vi.gainWidths_value2() = gainWidths_barrel[i].getValue(2);
0260         vi.gainWidths_value3() = gainWidths_barrel[i].getValue(3);
0261 
0262         vi.channelQuality_status() = channelQuality_barrel[i].getValue();
0263         vi.qieTypes_values() = qieTypes_barrel[i].getValue();
0264 
0265         for (uint32_t k = 0; k < 4; k++)
0266           for (uint32_t l = 0; l < 4; l++) {
0267             auto const linear = k * 4 + l;
0268             vi.qieCoders_offsets()[linear] = qieData_barrel[i].offset(k, l);
0269             vi.qieCoders_slopes()[linear] = qieData_barrel[i].slope(k, l);
0270           }
0271 
0272         vi.sipmPar_type() = sipmParameters_barrel[i].getType();
0273         vi.sipmPar_auxi1() = sipmParameters_barrel[i].getauxi1();
0274         vi.sipmPar_fcByPE() = sipmParameters_barrel[i].getFCByPE();
0275         vi.sipmPar_darkCurrent() = sipmParameters_barrel[i].getDarkCurrent();
0276         vi.sipmPar_auxi2() = sipmParameters_barrel[i].getauxi2();
0277       }
0278 
0279       // fill HE channels
0280       auto const& pedestals_endcaps = pedestals_containers[1].second;
0281       auto const& effectivePedestals_endcaps = effectivePedestals_containers[1].second;
0282       auto const& gains_endcaps = gains_containers[1].second;
0283       auto const& lutCorrs_endcaps = lutCorrs_containers[1].second;
0284       auto const& respCorrs_endcaps = respCorrs_containers[1].second;
0285       auto const& timeCorrs_endcaps = timeCorrs_containers[1].second;
0286       auto const& pedestalWidths_endcaps = pedestalWidths_containers[1].second;
0287       auto const& effectivePedestalWidths_endcaps = effectivePedestalWidths_containers[1].second;
0288       auto const& gainWidths_endcaps = gainWidths_containers[1].second;
0289       auto const& channelQuality_endcaps = channelQuality_containers[1].second;
0290       auto const& qieTypes_endcaps = qieTypes_containers[1].second;
0291       auto const& qieData_endcaps = qieData_containers[1].second;
0292       auto const& sipmParameters_endcaps = sipmParameters_containers[1].second;
0293 
0294       auto const offset = pedestals_barrel.size();
0295 
0296       for (uint64_t i = 0; i < pedestals_endcaps.size(); ++i) {
0297         auto const& qieCoder = qieData_endcaps[i];
0298         auto const qieType = qieTypes_endcaps[i].getValue() > 1 ? 1 : 0;
0299         auto const& qieShape = qieData.getShape(qieType);
0300 
0301         auto vi = view[offset + i];
0302 
0303         vi.pedestals_value()[0] = ped_unitIsADC ? convertPed(pedestals_endcaps[i].getValue(0), 0, qieCoder, qieShape)
0304                                                 : pedestals_endcaps[i].getValue(0);
0305         vi.pedestals_value()[1] = ped_unitIsADC ? convertPed(pedestals_endcaps[i].getValue(1), 1, qieCoder, qieShape)
0306                                                 : pedestals_endcaps[i].getValue(1);
0307         vi.pedestals_value()[2] = ped_unitIsADC ? convertPed(pedestals_endcaps[i].getValue(2), 2, qieCoder, qieShape)
0308                                                 : pedestals_endcaps[i].getValue(2);
0309         vi.pedestals_value()[3] = ped_unitIsADC ? convertPed(pedestals_endcaps[i].getValue(3), 3, qieCoder, qieShape)
0310                                                 : pedestals_endcaps[i].getValue(3);
0311 
0312         vi.pedestals_width()[0] =
0313             ped_unitIsADC
0314                 ? convertPedWidths(
0315                       pedestals_endcaps[i].getValue(0), pedestalWidths_endcaps[i].getWidth(0), 0, qieCoder, qieShape)
0316                 : pedestalWidths_endcaps[i].getWidth(0);
0317         vi.pedestals_width()[1] =
0318             ped_unitIsADC
0319                 ? convertPedWidths(
0320                       pedestals_endcaps[i].getValue(1), pedestalWidths_endcaps[i].getWidth(1), 1, qieCoder, qieShape)
0321                 : pedestalWidths_endcaps[i].getWidth(1);
0322         vi.pedestals_width()[2] =
0323             ped_unitIsADC
0324                 ? convertPedWidths(
0325                       pedestals_endcaps[i].getValue(2), pedestalWidths_endcaps[i].getWidth(2), 2, qieCoder, qieShape)
0326                 : pedestalWidths_endcaps[i].getWidth(2);
0327         vi.pedestals_width()[3] =
0328             ped_unitIsADC
0329                 ? convertPedWidths(
0330                       pedestals_endcaps[i].getValue(3), pedestalWidths_endcaps[i].getWidth(3), 3, qieCoder, qieShape)
0331                 : pedestalWidths_endcaps[i].getWidth(3);
0332 
0333         vi.effectivePedestals()[0] = effPed_unitIsADC
0334                                          ? convertPed(effectivePedestals_endcaps[i].getValue(0), 0, qieCoder, qieShape)
0335                                          : effectivePedestals_endcaps[i].getValue(0);
0336         vi.effectivePedestals()[1] = effPed_unitIsADC
0337                                          ? convertPed(effectivePedestals_endcaps[i].getValue(1), 1, qieCoder, qieShape)
0338                                          : effectivePedestals_endcaps[i].getValue(1);
0339         vi.effectivePedestals()[2] = effPed_unitIsADC
0340                                          ? convertPed(effectivePedestals_endcaps[i].getValue(2), 2, qieCoder, qieShape)
0341                                          : effectivePedestals_endcaps[i].getValue(2);
0342         vi.effectivePedestals()[3] = effPed_unitIsADC
0343                                          ? convertPed(effectivePedestals_endcaps[i].getValue(3), 3, qieCoder, qieShape)
0344                                          : effectivePedestals_endcaps[i].getValue(3);
0345 
0346         vi.effectivePedestalWidths()[0] = effPed_unitIsADC
0347                                               ? convertPedWidths(effectivePedestals_endcaps[i].getValue(0),
0348                                                                  effectivePedestalWidths_endcaps[i].getWidth(0),
0349                                                                  0,
0350                                                                  qieCoder,
0351                                                                  qieShape)
0352                                               : effectivePedestalWidths_endcaps[i].getWidth(0);
0353         vi.effectivePedestalWidths()[1] = effPed_unitIsADC
0354                                               ? convertPedWidths(effectivePedestals_endcaps[i].getValue(1),
0355                                                                  effectivePedestalWidths_endcaps[i].getWidth(1),
0356                                                                  1,
0357                                                                  qieCoder,
0358                                                                  qieShape)
0359                                               : effectivePedestalWidths_endcaps[i].getWidth(1);
0360         vi.effectivePedestalWidths()[2] = effPed_unitIsADC
0361                                               ? convertPedWidths(effectivePedestals_endcaps[i].getValue(2),
0362                                                                  effectivePedestalWidths_endcaps[i].getWidth(2),
0363                                                                  2,
0364                                                                  qieCoder,
0365                                                                  qieShape)
0366                                               : effectivePedestalWidths_endcaps[i].getWidth(2);
0367         vi.effectivePedestalWidths()[3] = effPed_unitIsADC
0368                                               ? convertPedWidths(effectivePedestals_endcaps[i].getValue(3),
0369                                                                  effectivePedestalWidths_endcaps[i].getWidth(3),
0370                                                                  3,
0371                                                                  qieCoder,
0372                                                                  qieShape)
0373                                               : effectivePedestalWidths_endcaps[i].getWidth(3);
0374 
0375         vi.gains_value()[0] = gains_endcaps[i].getValue(0);
0376         vi.gains_value()[1] = gains_endcaps[i].getValue(1);
0377         vi.gains_value()[2] = gains_endcaps[i].getValue(2);
0378         vi.gains_value()[3] = gains_endcaps[i].getValue(3);
0379 
0380         vi.lutCorrs_values() = lutCorrs_endcaps[i].getValue();
0381         vi.respCorrs_values() = respCorrs_endcaps[i].getValue();
0382         vi.timeCorrs_values() = timeCorrs_endcaps[i].getValue();
0383 
0384         vi.pedestalWidths_sigma00() = *(pedestalWidths_endcaps[i].getValues());
0385         vi.pedestalWidths_sigma01() = *(pedestalWidths_endcaps[i].getValues() + 1);
0386         vi.pedestalWidths_sigma02() = *(pedestalWidths_endcaps[i].getValues() + 2);
0387         vi.pedestalWidths_sigma03() = *(pedestalWidths_endcaps[i].getValues() + 3);
0388         vi.pedestalWidths_sigma10() = *(pedestalWidths_endcaps[i].getValues() + 4);
0389         vi.pedestalWidths_sigma11() = *(pedestalWidths_endcaps[i].getValues() + 5);
0390         vi.pedestalWidths_sigma12() = *(pedestalWidths_endcaps[i].getValues() + 6);
0391         vi.pedestalWidths_sigma13() = *(pedestalWidths_endcaps[i].getValues() + 7);
0392         vi.pedestalWidths_sigma20() = *(pedestalWidths_endcaps[i].getValues() + 8);
0393         vi.pedestalWidths_sigma21() = *(pedestalWidths_endcaps[i].getValues() + 9);
0394         vi.pedestalWidths_sigma22() = *(pedestalWidths_endcaps[i].getValues() + 10);
0395         vi.pedestalWidths_sigma23() = *(pedestalWidths_endcaps[i].getValues() + 11);
0396         vi.pedestalWidths_sigma30() = *(pedestalWidths_endcaps[i].getValues() + 12);
0397         vi.pedestalWidths_sigma31() = *(pedestalWidths_endcaps[i].getValues() + 13);
0398         vi.pedestalWidths_sigma32() = *(pedestalWidths_endcaps[i].getValues() + 14);
0399         vi.pedestalWidths_sigma33() = *(pedestalWidths_endcaps[i].getValues() + 15);
0400 
0401         vi.gainWidths_value0() = gainWidths_endcaps[i].getValue(0);
0402         vi.gainWidths_value1() = gainWidths_endcaps[i].getValue(1);
0403         vi.gainWidths_value2() = gainWidths_endcaps[i].getValue(2);
0404         vi.gainWidths_value3() = gainWidths_endcaps[i].getValue(3);
0405 
0406         vi.channelQuality_status() = channelQuality_endcaps[i].getValue();
0407         vi.qieTypes_values() = qieTypes_endcaps[i].getValue();
0408 
0409         for (uint32_t k = 0; k < 4; k++)
0410           for (uint32_t l = 0; l < 4; l++) {
0411             auto const linear = k * 4u + l;
0412             vi.qieCoders_offsets()[linear] = qieData_endcaps[i].offset(k, l);
0413             vi.qieCoders_slopes()[linear] = qieData_endcaps[i].slope(k, l);
0414           }
0415 
0416         vi.sipmPar_type() = sipmParameters_endcaps[i].getType();
0417         vi.sipmPar_auxi1() = sipmParameters_endcaps[i].getauxi1();
0418         vi.sipmPar_fcByPE() = sipmParameters_endcaps[i].getFCByPE();
0419         vi.sipmPar_darkCurrent() = sipmParameters_endcaps[i].getDarkCurrent();
0420         vi.sipmPar_auxi2() = sipmParameters_endcaps[i].getauxi2();
0421       }
0422       //fill the scalars
0423       static const int IPHI_MAX = 72;  // private member of topology
0424 
0425       view.maxDepthHB() = topology.maxDepthHB();
0426       view.maxDepthHE() = topology.maxDepthHE();
0427       view.maxPhiHE() = recConstants.getNPhi(1) > IPHI_MAX ? recConstants.getNPhi(1) : IPHI_MAX;
0428       view.firstHBRing() = topology.firstHBRing();
0429       view.lastHBRing() = topology.lastHBRing();
0430       view.firstHERing() = topology.firstHERing();
0431       view.lastHERing() = topology.lastHERing();
0432       view.nEtaHB() = recConstants.getEtaRange(0).second - recConstants.getEtaRange(0).first + 1;
0433       view.nEtaHE() =
0434           topology.firstHERing() > topology.lastHERing() ? 0 : (topology.lastHERing() - topology.firstHERing() + 1);
0435       view.offsetForHashes() = offset;
0436 
0437       return product;
0438     }
0439 
0440   private:
0441     edm::ESGetToken<HcalRecoParams, HcalRecoParamsRcd> recoParamsToken_;
0442     edm::ESGetToken<HcalPedestals, HcalPedestalsRcd> pedestalsToken_;
0443     edm::ESGetToken<HcalPedestals, HcalPedestalsRcd> effectivePedestalsToken_;
0444     edm::ESGetToken<HcalGains, HcalGainsRcd> gainsToken_;
0445     edm::ESGetToken<HcalLUTCorrs, HcalLUTCorrsRcd> lutCorrsToken_;
0446     edm::ESGetToken<HcalRespCorrs, HcalRespCorrsRcd> respCorrsToken_;
0447     edm::ESGetToken<HcalTimeCorrs, HcalTimeCorrsRcd> timeCorrsToken_;
0448     edm::ESGetToken<HcalPedestalWidths, HcalPedestalWidthsRcd> pedestalWidthsToken_;
0449     edm::ESGetToken<HcalPedestalWidths, HcalPedestalWidthsRcd> effectivePedestalWidthsToken_;
0450     edm::ESGetToken<HcalGainWidths, HcalGainWidthsRcd> gainWidthsToken_;
0451     edm::ESGetToken<HcalChannelQuality, HcalChannelQualityRcd> channelQualityToken_;
0452     edm::ESGetToken<HcalQIETypes, HcalQIETypesRcd> qieTypesToken_;
0453     edm::ESGetToken<HcalQIEData, HcalQIEDataRcd> qieDataToken_;
0454     edm::ESGetToken<HcalSiPMParameters, HcalSiPMParametersRcd> sipmParametersToken_;
0455     edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> topologyToken_;
0456     edm::ESGetToken<HcalDDDRecConstants, HcalRecNumberingRecord> recConstantsToken_;
0457   };
0458 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0459 
0460 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(HcalMahiConditionsESProducer);