File indexing completed on 2025-04-24 01:30:28
0001 #include <alpaka/alpaka.hpp>
0002
0003 #include "DataFormats/EcalDigi/interface/alpaka/EcalDigiPhase2DeviceCollection.h"
0004 #include "DataFormats/EcalDigi/interface/EcalDataFrame_Ph2.h"
0005 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
0006 #include "DataFormats/EcalRecHit/interface/alpaka/EcalUncalibratedRecHitDeviceCollection.h"
0007 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
0008
0009 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0010 #include "HeterogeneousCore/AlpakaInterface/interface/traits.h"
0011 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0012
0013 #include "EcalUncalibRecHitPhase2WeightsAlgoPortable.h"
0014 #include "EcalUncalibRecHitPhase2WeightsStruct.h"
0015
0016 namespace ALPAKA_ACCELERATOR_NAMESPACE::ecal::weights {
0017 using namespace cms::alpakatools;
0018
0019 class Phase2WeightsKernel {
0020 public:
0021 ALPAKA_FN_ACC void operator()(Acc1D const& acc,
0022 EcalUncalibRecHitPhase2Weights const* weightsObj,
0023 EcalDigiPhase2DeviceCollection::ConstView digisDev,
0024 EcalUncalibratedRecHitDeviceCollection::View uncalibratedRecHitsDev) const {
0025 constexpr int nsamples = ecalPh2::sampleSize;
0026 auto const nchannels = digisDev.size();
0027
0028 if (once_per_grid(acc)) {
0029 uncalibratedRecHitsDev.size() = digisDev.size();
0030 }
0031
0032 auto const* weightsdata = weightsObj->weights.data();
0033 auto const* timeWeightsdata = weightsObj->timeWeights.data();
0034
0035 for (auto tx : uniform_elements(acc, nchannels)) {
0036 bool g1 = false;
0037 const auto& digi = digisDev[tx].data();
0038 auto recHit = uncalibratedRecHitsDev[tx];
0039 recHit.amplitude() = 0;
0040 recHit.jitter() = 0;
0041 for (int s = 0; s < nsamples; ++s) {
0042 const auto sample = digi[s];
0043 const auto trace =
0044 (static_cast<float>(ecalLiteDTU::adc(sample))) * ecalPh2::gains[ecalLiteDTU::gainId(sample)];
0045 recHit.amplitude() += (trace * weightsdata[s]);
0046 recHit.jitter() += (trace * timeWeightsdata[s]);
0047 if (ecalLiteDTU::gainId(sample) == 1)
0048 g1 = true;
0049 recHit.outOfTimeAmplitudes()[s] = 0.;
0050 }
0051 recHit.amplitudeError() = 1.0f;
0052 recHit.id() = digisDev.id()[tx];
0053 recHit.flags() = 0;
0054 recHit.pedestal() = 0.;
0055 recHit.jitterError() = 0.;
0056 recHit.chi2() = 0.;
0057 recHit.aux() = 0;
0058 if (g1) {
0059 recHit.flags() = 0x1 << EcalUncalibratedRecHit::kHasSwitchToGain1;
0060 }
0061 }
0062 }
0063 };
0064
0065 void phase2Weights(EcalDigiPhase2DeviceCollection const& digis,
0066 EcalUncalibratedRecHitDeviceCollection& uncalibratedRecHits,
0067 EcalUncalibRecHitPhase2Weights const* weightsObj,
0068 Queue& queue) {
0069
0070 uint32_t items = 64;
0071
0072 uint32_t groups = divide_up_by(digis->metadata().size(), items);
0073
0074 auto workDiv = make_workdiv<Acc1D>(groups, items);
0075
0076 alpaka::exec<Acc1D>(
0077 queue, workDiv, Phase2WeightsKernel{}, weightsObj, digis.const_view(), uncalibratedRecHits.view());
0078 }
0079
0080 }