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 
0011 #include "CondFormats/DataRecord/interface/HGCalElectronicsMappingRcd.h"
0012 #include "CondFormats/HGCalObjects/interface/HGCalMappingCellIndexer.h"
0013 #include "CondFormats/HGCalObjects/interface/HGCalMappingParameterHostCollection.h"
0014 #include "CondFormats/HGCalObjects/interface/alpaka/HGCalMappingParameterDeviceCollection.h"
0015 #include "DataFormats/HGCalDigi/interface/HGCalElectronicsId.h"
0016 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
0017 #include "DataFormats/ForwardDetId/interface/HGCScintillatorDetId.h"
0018 #include "Geometry/HGCalMapping/interface/HGCalMappingTools.h"
0019 
0020 #include <string>
0021 #include <iostream>
0022 #include <fstream>
0023 #include <sstream>
0024 
0025 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0026 
0027   namespace hgcal {
0028 
0029     class HGCalMappingCellESProducer : public ESProducer {
0030     public:
0031       //
0032       HGCalMappingCellESProducer(const edm::ParameterSet& iConfig)
0033           : ESProducer(iConfig), filelist_(iConfig.getParameter<std::vector<std::string> >("filelist")) {
0034         auto cc = setWhatProduced(this);
0035         cellIndexTkn_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("cellindexer"));
0036       }
0037 
0038       //
0039       static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0040         edm::ParameterSetDescription desc;
0041         desc.add<std::vector<std::string> >("filelist", std::vector<std::string>({}))
0042             ->setComment("list of files with the readout cells of each module");
0043         desc.add<edm::ESInputTag>("cellindexer", edm::ESInputTag(""))->setComment("Dense cell index tool");
0044         descriptions.addWithDefaultLabel(desc);
0045       }
0046 
0047       //
0048       std::optional<HGCalMappingCellParamHostCollection> produce(const HGCalElectronicsMappingRcd& iRecord) {
0049         //get cell indexer
0050         const HGCalMappingCellIndexer& cellIndexer = iRecord.get(cellIndexTkn_);
0051         const uint32_t size = cellIndexer.maxDenseIndex();  // channel-level size
0052         HGCalMappingCellParamHostCollection cellParams(size, cms::alpakatools::host());
0053         for (uint32_t i = 0; i < size; i++)
0054           cellParams.view()[i].valid() = false;
0055 
0056         //loop over cell types and then over cells
0057         for (auto url : filelist_) {
0058           ::hgcal::mappingtools::HGCalEntityList pmap;
0059           edm::FileInPath fip(url);
0060           pmap.buildFrom(fip.fullPath());
0061           auto& entities = pmap.getEntries();
0062           for (auto row : entities) {
0063             //identify special cases (Si vs SiPM, calib vs normal)
0064             std::string typecode = pmap.getAttr("Typecode", row);
0065             auto typeidx = cellIndexer.getEnumFromTypecode(typecode);
0066             bool isSiPM = typecode.find("TM") != std::string::npos;
0067             int rocpin = pmap.getIntAttr("ROCpin", row);
0068             int celltype = pmap.getIntAttr("t", row);
0069             int i1(0), i2(0), sensorcell(0);
0070             bool isHD(false), iscalib(false);
0071             uint32_t detid(0);
0072             if (isSiPM) {
0073               i1 = pmap.getIntAttr("iring", row);
0074               i2 = pmap.getIntAttr("iphi", row);
0075             } else {
0076               i1 = pmap.getIntAttr("iu", row);
0077               i2 = pmap.getIntAttr("iv", row);
0078               isHD = {typecode.find("MH") != std::string::npos ? true : false};
0079               sensorcell = pmap.getIntAttr("SiCell", row);
0080               if (celltype == 0) {
0081                 iscalib = true;
0082                 std::string rocpinstr = pmap.getAttr("ROCpin", row);
0083                 rocpin = uint16_t(rocpinstr[rocpinstr.size() - 1]);
0084               }
0085 
0086               //detector id is initiated for a random sub-detector with Si wafers
0087               //we only care about the first 10bits where the cell u-v are stored
0088               DetId::Detector det(DetId::Detector::HGCalEE);
0089               detid = 0x3ff & HGCSiliconDetId(det, 0, 0, 0, 0, 0, i1, i2).rawId();
0090             }
0091 
0092             //fill cell info in the appopriate dense index
0093             int chip = pmap.getIntAttr("ROC", row);
0094             int half = pmap.getIntAttr("HalfROC", row);
0095             int seq = pmap.getIntAttr("Seq", row);
0096             int idx = cellIndexer.denseIndex(typecode, chip, half, seq);
0097             auto cell = cellParams.view()[idx];
0098             cell.valid() = true;
0099             cell.isHD() = isHD;
0100             cell.iscalib() = iscalib;
0101             cell.isSiPM() = isSiPM;
0102             cell.typeidx() = typeidx;
0103             cell.chip() = chip;
0104             cell.half() = half;
0105             cell.seq() = seq;
0106             cell.rocpin() = rocpin;
0107             cell.sensorcell() = sensorcell;
0108             cell.triglink() = pmap.getIntAttr("TrLink", row);
0109             cell.trigcell() = pmap.getIntAttr("TrCell", row);
0110             cell.i1() = i1;
0111             cell.i2() = i2;
0112             cell.t() = celltype;
0113             cell.trace() = pmap.getFloatAttr("trace", row);
0114             cell.eleid() = HGCalElectronicsId(false, 0, 0, 0, chip * 2 + half, seq).raw();
0115             cell.detid() = detid;
0116           }  //end loop over entities
0117         }    //end loop over cell types
0118 
0119         return cellParams;
0120       }  // end of produce()
0121 
0122     private:
0123       edm::ESGetToken<HGCalMappingCellIndexer, HGCalElectronicsMappingRcd> cellIndexTkn_;
0124       const std::vector<std::string> filelist_;
0125     };
0126 
0127   }  // namespace hgcal
0128 
0129 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0130 
0131 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(hgcal::HGCalMappingCellESProducer);