Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-05-19 07:20:05

0001 #include "FWCore/ParameterSet/interface/FileInPath.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "FWCore/Utilities/interface/ESGetToken.h"
0005 
0006 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0007 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0008 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0009 #include "HeterogeneousCore/AlpakaInterface/interface/host.h"
0010 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0011 
0012 #include "CondFormats/DataRecord/interface/HGCalElectronicsMappingRcd.h"
0013 #include "CondFormats/HGCalObjects/interface/HGCalMappingCellIndexer.h"
0014 #include "CondFormats/HGCalObjects/interface/HGCalMappingParameterHost.h"
0015 #include "CondFormats/HGCalObjects/interface/alpaka/HGCalMappingParameterDevice.h"
0016 #include "CondFormats/HGCalObjects/interface/HGCalMappingModuleIndexer.h"
0017 #include "DataFormats/HGCalDigi/interface/HGCalElectronicsId.h"
0018 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
0019 #include "DataFormats/ForwardDetId/interface/HGCScintillatorDetId.h"
0020 #include "Geometry/HGCalMapping/interface/HGCalMappingTools.h"
0021 
0022 #include <string>
0023 #include <iostream>
0024 #include <fstream>
0025 #include <sstream>
0026 #include <charconv>
0027 
0028 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0029 
0030   namespace hgcal {
0031 
0032     class HGCalMappingCellESProducer : public ESProducer {
0033     public:
0034       //
0035       HGCalMappingCellESProducer(const edm::ParameterSet& iConfig)
0036           : ESProducer(iConfig),
0037             filelist_(iConfig.getParameter<std::vector<std::string> >("filelist")),
0038             offsetfile_(iConfig.getParameter<edm::FileInPath>("offsetfile")) {
0039         auto cc = setWhatProduced(this);
0040         cellIndexTkn_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("cellindexer"));
0041         moduleIndexTkn_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("moduleindexer"));
0042       }
0043 
0044       //
0045       static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0046         edm::ParameterSetDescription desc;
0047         desc.add<std::vector<std::string> >("filelist", std::vector<std::string>({}))
0048             ->setComment("list of files with the readout cells of each module");
0049         desc.add<edm::ESInputTag>("cellindexer", edm::ESInputTag(""))->setComment("Dense cell index tool");
0050         desc.add<edm::ESInputTag>("moduleindexer", edm::ESInputTag(""))->setComment("Module index tool");
0051         desc.add<edm::FileInPath>(
0052                 "offsetfile",
0053                 edm::FileInPath("Geometry/HGCalMapping/data/CellMaps/calibration_to_surrounding_offsetMap.txt"))
0054             ->setComment("file containing the offsets between calibration and surrounding cells");
0055         descriptions.addWithDefaultLabel(desc);
0056       }
0057 
0058       std::map<int, int> makeOffsetMap(edm::FileInPath input_offsetfile,
0059                                        const HGCalMappingCellIndexer& cellIndexer,
0060                                        const HGCalMappingModuleIndexer& moduleIndexer) {
0061         std::map<int, int> offsetMap;
0062         const auto& offsetfile = input_offsetfile.fullPath();
0063         ::hgcal::mappingtools::HGCalEntityList omap;
0064         edm::FileInPath fip(offsetfile);
0065         omap.buildFrom(fip.fullPath());
0066         auto& mapEntries = omap.getEntries();
0067 
0068         for (auto row : mapEntries) {
0069           std::string typecode = omap.getAttr("Typecode", row);
0070           const auto& allTypecodes = moduleIndexer.getTypecodeMap();
0071           // Skip if typecode is not in the module indexer
0072           bool typecodeFound = false;
0073           for (const auto& key : allTypecodes) {
0074             if (key.first.find(typecode) != std::string::npos) {
0075               typecodeFound = true;
0076               break;
0077             }
0078           }
0079           if (!typecodeFound)
0080             continue;
0081           int roc = omap.getIntAttr("ROC", row);
0082           int halfroc = omap.getIntAttr("HalfROC", row);
0083           int readoutsequence = omap.getIntAttr("Seq", row);
0084           int offset = omap.getIntAttr("Offset", row);
0085           int idx = cellIndexer.denseIndex(typecode, roc, halfroc, readoutsequence);
0086           offsetMap[idx] = offset;
0087         }
0088         return offsetMap;
0089       }
0090 
0091       //
0092       std::optional<HGCalMappingCellParamHost> produce(const HGCalElectronicsMappingRcd& iRecord) {
0093         //get cell and module indexers
0094         const HGCalMappingCellIndexer& cellIndexer = iRecord.get(cellIndexTkn_);
0095         const HGCalMappingModuleIndexer& moduleIndexer = iRecord.get(moduleIndexTkn_);
0096         const uint32_t size = cellIndexer.maxDenseIndex();  // channel-level size
0097         HGCalMappingCellParamHost cellParams(size, cms::alpakatools::host());
0098         for (uint32_t i = 0; i < size; i++)
0099           cellParams.view()[i].valid() = false;
0100 
0101         auto offsetMap = makeOffsetMap(offsetfile_, cellIndexer, moduleIndexer);
0102 
0103         //loop over cell types and then over cells
0104         for (const auto& url : filelist_) {
0105           ::hgcal::mappingtools::HGCalEntityList pmap;
0106           edm::FileInPath fip(url);
0107           pmap.buildFrom(fip.fullPath());
0108           auto& entities = pmap.getEntries();
0109           for (auto row : entities) {
0110             //identify special cases (Si vs SiPM, calib vs normal)
0111             std::string typecode = pmap.getAttr("Typecode", row);
0112             auto typeidx = cellIndexer.getEnumFromTypecode(typecode);
0113             bool isSiPM = typecode.find("TM") != std::string::npos;
0114             int rocpin = pmap.getIntAttr("ROCpin", row);
0115             int celltype = pmap.getIntAttr("t", row);
0116             int i1(0), i2(0), sensorcell(0);
0117             bool isHD(false), iscalib(false);
0118             uint32_t detid(0);
0119             if (isSiPM) {
0120               i1 = pmap.getIntAttr("iring", row);
0121               i2 = pmap.getIntAttr("iphi", row);
0122             } else {
0123               i1 = pmap.getIntAttr("iu", row);
0124               i2 = pmap.getIntAttr("iv", row);
0125               isHD = {typecode.find("MH") != std::string::npos ? true : false};
0126               sensorcell = pmap.getIntAttr("SiCell", row);
0127               if (celltype == 0) {
0128                 iscalib = true;
0129                 std::string rocpinstr = pmap.getAttr("ROCpin", row);
0130                 rocpin = uint16_t(rocpinstr[rocpinstr.size() - 1]);
0131               }
0132 
0133               //detector id is initiated for a random sub-detector with Si wafers
0134               //we only care about the first 10bits where the cell u-v are stored
0135               DetId::Detector det(DetId::Detector::HGCalEE);
0136               detid = 0x3ff & HGCSiliconDetId(det, 0, 0, 0, 0, 0, i1, i2).rawId();
0137             }
0138 
0139             //fill cell info in the appopriate dense index
0140             int chip = pmap.getIntAttr("ROC", row);
0141             int half = pmap.getIntAttr("HalfROC", row);
0142             int seq = pmap.getIntAttr("Seq", row);
0143             int idx = cellIndexer.denseIndex(typecode, chip, half, seq);
0144             auto cell = cellParams.view()[idx];
0145             cell.valid() = true;
0146             cell.isHD() = isHD;
0147             cell.iscalib() = iscalib;
0148             cell.isSiPM() = isSiPM;
0149             cell.typeidx() = typeidx;
0150             cell.chip() = chip;
0151             cell.half() = half;
0152             cell.seq() = seq;
0153             cell.rocpin() = rocpin;
0154             cell.sensorcell() = sensorcell;
0155             cell.triglink() = pmap.getIntAttr("TrLink", row);
0156             cell.trigcell() = pmap.getIntAttr("TrCell", row);
0157             cell.i1() = i1;
0158             cell.i2() = i2;
0159             cell.t() = celltype;
0160             cell.trace() = pmap.getFloatAttr("trace", row);
0161             cell.eleid() = HGCalElectronicsId(false, 0, 0, 0, chip * 2 + half, seq).raw();
0162             cell.detid() = detid;
0163 
0164             int offset = (iscalib && offsetMap.find(idx) != offsetMap.end()) ? offsetMap[idx] : 0;
0165             cell.caliboffset() = offset;
0166 
0167           }  //end loop over entities
0168         }  //end loop over cell types
0169 
0170         return cellParams;
0171       }  // end of produce()
0172 
0173     private:
0174       edm::ESGetToken<HGCalMappingCellIndexer, HGCalElectronicsMappingRcd> cellIndexTkn_;
0175       edm::ESGetToken<HGCalMappingModuleIndexer, HGCalElectronicsMappingRcd> moduleIndexTkn_;
0176       const std::vector<std::string> filelist_;
0177       edm::FileInPath offsetfile_;
0178     };
0179 
0180   }  // namespace hgcal
0181 
0182 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0183 
0184 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(hgcal::HGCalMappingCellESProducer);