Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-04 04:04:30

0001 #include "FWCore/ParameterSet/interface/FileInPath.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 #include "FWCore/Utilities/interface/ESGetToken.h"
0004 
0005 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0006 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0007 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0008 #include "HeterogeneousCore/AlpakaInterface/interface/host.h"
0009 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0010 #include "CondFormats/DataRecord/interface/HGCalElectronicsMappingRcd.h"
0011 #include "CondFormats/HGCalObjects/interface/HGCalMappingModuleIndexer.h"
0012 #include "CondFormats/HGCalObjects/interface/HGCalMappingParameterHostCollection.h"
0013 #include "CondFormats/HGCalObjects/interface/alpaka/HGCalMappingParameterDeviceCollection.h"
0014 #include "DataFormats/HGCalDigi/interface/HGCalElectronicsId.h"
0015 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
0016 #include "DataFormats/ForwardDetId/interface/HGCScintillatorDetId.h"
0017 #include "Geometry/HGCalMapping/interface/HGCalMappingTools.h"
0018 
0019 #include <string>
0020 #include <iostream>
0021 #include <fstream>
0022 #include <sstream>
0023 
0024 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0025 
0026   namespace hgcal {
0027 
0028     class HGCalMappingModuleESProducer : public ESProducer {
0029     public:
0030       //
0031       HGCalMappingModuleESProducer(const edm::ParameterSet& iConfig)
0032           : ESProducer(iConfig), filename_(iConfig.getParameter<edm::FileInPath>("filename")) {
0033         auto cc = setWhatProduced(this);
0034         moduleIndexTkn_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("moduleindexer"));
0035       }
0036 
0037       //
0038       static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0039         edm::ParameterSetDescription desc;
0040         desc.add<edm::FileInPath>("filename")->setComment("module locator file");
0041         desc.add<edm::ESInputTag>("moduleindexer", edm::ESInputTag(""))->setComment("Dense module index tool");
0042         descriptions.addWithDefaultLabel(desc);
0043       }
0044 
0045       //
0046       std::optional<HGCalMappingModuleParamHostCollection> produce(const HGCalElectronicsMappingRcd& iRecord) {
0047         //get cell and module indexer
0048         auto modIndexer = iRecord.get(moduleIndexTkn_);
0049 
0050         // load dense indexing
0051         const uint32_t size = modIndexer.maxModulesIdx_;
0052         HGCalMappingModuleParamHostCollection moduleParams(size, cms::alpakatools::host());
0053         for (size_t i = 0; i < size; i++)
0054           moduleParams.view()[i].valid() = false;
0055 
0056         ::hgcal::mappingtools::HGCalEntityList pmap;
0057         pmap.buildFrom(filename_.fullPath());
0058         auto& entities = pmap.getEntries();
0059         for (auto row : entities) {
0060           int fedid = pmap.getIntAttr("fedid", row);
0061           int captureblockidx = pmap.getIntAttr("captureblockidx", row);
0062           int econdidx = pmap.getIntAttr("econdidx", row);
0063           int idx = modIndexer.getIndexForModule(fedid, captureblockidx, econdidx);
0064           int typeidx = modIndexer.getTypeForModule(fedid, captureblockidx, econdidx);
0065           std::string typecode = pmap.getAttr("typecode", row);
0066           auto celltypes = modIndexer.convertTypeCode(typecode);
0067           bool isSiPM = celltypes.first;
0068           int celltype = celltypes.second;
0069           int zside = pmap.getIntAttr("zside", row);
0070           int plane = pmap.getIntAttr("plane", row);
0071           int i1 = pmap.getIntAttr("u", row);
0072           int i2 = pmap.getIntAttr("v", row);
0073           uint32_t eleid = HGCalElectronicsId((zside > 0), fedid, captureblockidx, econdidx, 0, 0).raw();
0074           uint32_t detid(0);
0075           if (!isSiPM) {
0076             int zp(zside > 0 ? 1 : -1);
0077             DetId::Detector det = plane <= 26 ? DetId::Detector::HGCalEE : DetId::Detector::HGCalHSi;
0078             detid = HGCSiliconDetId(det, zp, celltype, plane, i1, i2, 0, 0).rawId();
0079           }
0080 
0081           auto module = moduleParams.view()[idx];
0082           module.valid() = true;
0083           module.zside() = (zside > 0);
0084           module.isSiPM() = isSiPM;
0085           module.celltype() = celltype;
0086           module.plane() = plane;
0087           module.i1() = i1;
0088           module.i2() = i2;
0089           module.typeidx() = typeidx;
0090           module.fedid() = fedid;
0091           module.slinkidx() = pmap.getIntAttr("slinkidx", row);
0092           module.captureblock() = pmap.getIntAttr("captureblock", row);
0093           ;
0094           module.econdidx() = econdidx;
0095           module.captureblockidx() = captureblockidx;
0096           module.eleid() = eleid;
0097           module.detid() = detid;
0098         }
0099 
0100         return moduleParams;
0101 
0102       }  // end of produce()
0103 
0104     private:
0105       edm::ESGetToken<HGCalMappingModuleIndexer, HGCalElectronicsMappingRcd> moduleIndexTkn_;
0106       const edm::FileInPath filename_;
0107     };
0108 
0109   }  // namespace hgcal
0110 
0111 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0112 
0113 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(hgcal::HGCalMappingModuleESProducer);