Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "FWCore/Framework/interface/ESTransientHandle.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 
0004 #include "CondFormats/HcalObjects/interface/HcalRecoParamWithPulseShapeHost.h"
0005 #include "CondFormats/HcalObjects/interface/HcalRecoParamWithPulseShapeSoA.h"
0006 #include "CondFormats/DataRecord/interface/HcalRecoParamsRcd.h"
0007 
0008 #include "CalibCalorimetry/HcalAlgos/interface/HcalPulseShapes.h"
0009 #include "RecoLocalCalo/HcalRecAlgos/interface/PulseShapeFunctor.h"
0010 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalConstants.h"
0011 
0012 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESGetToken.h"
0013 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0014 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0015 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0016 #include "HeterogeneousCore/AlpakaInterface/interface/host.h"
0017 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0018 
0019 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0020   class HcalRecoParamWithPulseShapeESProducer : public ESProducer {
0021   public:
0022     HcalRecoParamWithPulseShapeESProducer(edm::ParameterSet const& iConfig) : ESProducer(iConfig) {
0023       auto cc = setWhatProduced(this);
0024       recoParamsToken_ = cc.consumes();
0025     }
0026 
0027     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0028       edm::ParameterSetDescription desc;
0029       descriptions.addWithDefaultLabel(desc);
0030     }
0031 
0032     std::unique_ptr<hcal::HcalRecoParamWithPulseShapeHost> produce(HcalRecoParamsRcd const& iRecord) {
0033       auto const& recoParams = iRecord.get(recoParamsToken_);
0034 
0035       auto const containers = recoParams.getAllContainers();
0036       size_t const totalChannels =
0037           recoParams.getAllContainers()[0].second.size() + recoParams.getAllContainers()[1].second.size();
0038 
0039       //Get unique ids
0040       HcalPulseShapes pulseShapes;
0041       std::unordered_map<unsigned int, uint32_t> idCache;  //<pulseShapeId,pulseShapeIdx>
0042 
0043       auto const& barrelValues = containers[0].second;
0044       for (uint64_t i = 0; i < barrelValues.size(); ++i) {
0045         auto const pulseShapeId = barrelValues[i].pulseShapeID();
0046         if (pulseShapeId == 0)
0047           continue;
0048         if (auto const iter = idCache.find(pulseShapeId); iter == idCache.end()) {
0049           idCache[pulseShapeId] = idCache.size();
0050         }
0051       }
0052       auto const& endcapValues = containers[1].second;
0053       for (uint64_t i = 0; i < endcapValues.size(); ++i) {
0054         auto const pulseShapeId = endcapValues[i].pulseShapeID();
0055         if (pulseShapeId == 0)
0056           continue;
0057         if (auto const iter = idCache.find(pulseShapeId); iter == idCache.end()) {
0058           idCache[pulseShapeId] = idCache.size();
0059         }
0060       }
0061 
0062       // Fill products
0063       auto product = std::make_unique<hcal::HcalRecoParamWithPulseShapeHost>(
0064           totalChannels, idCache.size(), cms::alpakatools::host());
0065       auto recoView = product->recoParamView();
0066       auto pulseShapeView = product->pulseShapeView();
0067       for (uint64_t i = 0; i < barrelValues.size(); ++i) {
0068         auto vi = recoView[i];
0069         vi.param1() = barrelValues[i].param1();
0070         vi.param2() = barrelValues[i].param2();
0071         vi.ids() = (barrelValues[i].pulseShapeID() == 0)
0072                        ? 0
0073                        : idCache.at(barrelValues[i].pulseShapeID());  //idx of the pulseShape of channel i
0074       }
0075       // fill in endcap
0076       auto const offset = barrelValues.size();
0077       for (uint64_t i = 0; i < endcapValues.size(); ++i) {
0078         auto vi = recoView[i + offset];
0079         vi.param1() = endcapValues[i].param1();
0080         vi.param2() = endcapValues[i].param2();
0081         vi.ids() = (endcapValues[i].pulseShapeID() == 0)
0082                        ? 0
0083                        : idCache.at(endcapValues[i].pulseShapeID());  //idx of the pulseShape of channel i
0084       }
0085 
0086       //fill pulseShape views
0087       for (auto& it : idCache) {
0088         auto const pulseShapeId = it.first;
0089         auto const arrId = it.second;
0090         auto const& pulseShape = pulseShapes.getShape(pulseShapeId);
0091         FitterFuncs::PulseShapeFunctor functor{pulseShape, false, false, false, 1, 0, 0, hcal::constants::maxSamples};
0092 
0093         for (int i = 0; i < hcal::constants::maxPSshapeBin; i++) {
0094           pulseShapeView[arrId].acc25nsVec()[i] = functor.acc25nsVec()[i];
0095           pulseShapeView[arrId].diff25nsItvlVec()[i] = functor.diff25nsItvlVec()[i];
0096         }
0097         for (int i = 0; i < hcal::constants::nsPerBX; i++) {
0098           pulseShapeView[arrId].accVarLenIdxMinusOneVec()[i] = functor.accVarLenIdxMinusOneVec()[i];
0099           pulseShapeView[arrId].diffVarItvlIdxMinusOneVec()[i] = functor.diffVarItvlIdxMinusOneVec()[i];
0100           pulseShapeView[arrId].accVarLenIdxZEROVec()[i] = functor.accVarLenIdxZEROVec()[i];
0101           pulseShapeView[arrId].diffVarItvlIdxZEROVec()[i] = functor.diffVarItvlIdxZEROVec()[i];
0102         }
0103       }
0104 
0105       return product;
0106     }
0107 
0108   private:
0109     edm::ESGetToken<HcalRecoParams, HcalRecoParamsRcd> recoParamsToken_;
0110   };
0111 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0112 
0113 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(HcalRecoParamWithPulseShapeESProducer);