Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-18 03:42:04

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/DataRecord/interface/HGCalDenseIndexInfoRcd.h"
0012 #include "CondFormats/HGCalObjects/interface/HGCalMappingModuleIndexer.h"
0013 #include "CondFormats/HGCalObjects/interface/HGCalMappingParameterHost.h"
0014 #include "CondFormats/HGCalObjects/interface/alpaka/HGCalMappingParameterDevice.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 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0020 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0021 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
0022 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0023 
0024 #include <string>
0025 #include <iostream>
0026 #include <fstream>
0027 #include <sstream>
0028 
0029 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0030 
0031   namespace hgcal {
0032 
0033     class HGCalDenseIndexInfoESProducer : public ESProducer {
0034     public:
0035       //
0036       HGCalDenseIndexInfoESProducer(const edm::ParameterSet& iConfig) : ESProducer(iConfig) {
0037         auto cc = setWhatProduced(this);
0038         moduleIndexTkn_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("moduleindexer"));
0039         cellIndexTkn_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("cellindexer"));
0040         moduleInfoTkn_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("moduleinfo"));
0041         cellInfoTkn_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("cellinfo"));
0042         caloGeomToken_ = cc.consumes();
0043       }
0044 
0045       //
0046       static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0047         edm::ParameterSetDescription desc;
0048         desc.add<edm::ESInputTag>("moduleindexer", edm::ESInputTag(""))->setComment("Module index tool");
0049         desc.add<edm::ESInputTag>("cellindexer", edm::ESInputTag(""))->setComment("Dense cell index tool");
0050         desc.add<edm::ESInputTag>("moduleinfo", edm::ESInputTag(""))->setComment("Module info table");
0051         desc.add<edm::ESInputTag>("cellinfo", edm::ESInputTag(""))->setComment("Cell info table");
0052         descriptions.addWithDefaultLabel(desc);
0053       }
0054 
0055       //
0056       std::optional<HGCalDenseIndexInfoHost> produce(const HGCalDenseIndexInfoRcd& iRecord) {
0057         //geometry
0058         auto const& geo = iRecord.get(caloGeomToken_);
0059 
0060         //get cell and module indexer
0061         auto const& modIndexer = iRecord.get(moduleIndexTkn_);
0062         auto const& cellIndexer = iRecord.get(cellIndexTkn_);
0063 
0064         //get cell and module info
0065         auto const& moduleInfo = iRecord.get(moduleInfoTkn_);
0066         auto const& cellInfo = iRecord.get(cellInfoTkn_);
0067 
0068         //declare the dense index info collection to be produced
0069         //the size is determined by the module indexer
0070         const uint32_t nIndices = modIndexer.getMaxDataSize();
0071         HGCalDenseIndexInfoHost denseIdxInfo(nIndices, cms::alpakatools::host());
0072         for (auto fedRS : modIndexer.getFEDReadoutSequences()) {
0073           uint32_t fedId = fedRS.id;
0074           for (size_t imod = 0; imod < fedRS.readoutTypes_.size(); imod++) {
0075             //the number of words expected, the first channel dense index
0076             int modTypeIdx = fedRS.readoutTypes_[imod];
0077             uint32_t nch = modIndexer.getGlobalTypesNWords()[modTypeIdx];
0078             int off = fedRS.chDataOffsets_[imod];
0079 
0080             //get additional necessary module info
0081             int modIdx = modIndexer.getIndexForModule(fedId, imod);
0082             auto module_row = moduleInfo.view()[modIdx];
0083             bool isSiPM = module_row.isSiPM();
0084             uint16_t typeidx = module_row.typeidx();
0085             uint32_t eleid = module_row.eleid();
0086             uint32_t detid = module_row.detid();
0087 
0088             //get the appropriate geometry
0089             DetId::Detector det = DetId(detid).det();
0090             int subdet = ForwardSubdetector::ForwardEmpty;
0091             const HGCalGeometry* hgcal_geom =
0092                 static_cast<const HGCalGeometry*>(geo.getSubdetectorGeometry(det, subdet));
0093 
0094             //get the offset to start reading the cell info from sequential
0095             uint32_t cellInfoOffset = cellIndexer.offsets_[typeidx];
0096 
0097             //now fill the information sequentially on the cells of this module
0098             for (uint32_t ich = 0; ich < nch; ich++) {
0099               //finalize assigning the dense index
0100               uint32_t denseIdx = off + ich;
0101 
0102               auto row = denseIdxInfo.view()[denseIdx];
0103 
0104               //fill the fields
0105               row.fedId() = fedId;
0106               row.fedReadoutSeq() = imod;
0107               row.chNumber() = ich;
0108               row.modInfoIdx() = modIdx;
0109               uint32_t cellIdx = cellInfoOffset + ich;
0110               row.cellInfoIdx() = cellIdx;
0111 
0112               auto cell_row = cellInfo.view()[cellIdx];
0113               row.eleid() = eleid + cell_row.eleid();
0114 
0115               //assign det id only for full and calibration cells
0116               row.detid() = 0;
0117               if (cell_row.t() == 1 || cell_row.t() == 0) {
0118                 if (isSiPM) {
0119                   row.detid() = ::hgcal::mappingtools::getSiPMDetId(module_row.zside(),
0120                                                                     module_row.plane(),
0121                                                                     module_row.i2(),
0122                                                                     module_row.celltype(),
0123                                                                     cell_row.i1(),
0124                                                                     cell_row.i2());
0125                 } else {
0126                   row.detid() = module_row.detid() + cell_row.detid();
0127                 }
0128 
0129                 //assign position from geometry
0130                 row.x() = 0;
0131                 row.y() = 0;
0132                 row.z() = 0;
0133                 if (hgcal_geom != nullptr) {
0134                   GlobalPoint position = hgcal_geom->getPosition(row.detid());
0135                   row.x() = position.x();
0136                   row.y() = position.y();
0137                   row.z() = position.z();
0138                 }
0139               }
0140             }  // end cell loop
0141 
0142           }  // end module loop
0143 
0144         }  // end fed readout sequence loop
0145 
0146         return denseIdxInfo;
0147       }  // end of produce()
0148 
0149     private:
0150       edm::ESGetToken<HGCalMappingModuleIndexer, HGCalElectronicsMappingRcd> moduleIndexTkn_;
0151       edm::ESGetToken<HGCalMappingCellIndexer, HGCalElectronicsMappingRcd> cellIndexTkn_;
0152       edm::ESGetToken<hgcal::HGCalMappingModuleParamHost, HGCalElectronicsMappingRcd> moduleInfoTkn_;
0153       edm::ESGetToken<hgcal::HGCalMappingCellParamHost, HGCalElectronicsMappingRcd> cellInfoTkn_;
0154       edm::ESGetToken<CaloGeometry, CaloGeometryRecord> caloGeomToken_;
0155     };
0156 
0157   }  // namespace hgcal
0158 
0159 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0160 
0161 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(hgcal::HGCalDenseIndexInfoESProducer);