Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-02 00:53:55

0001 #include "Geometry/HGCalCommonData/interface/HGCalGeomRotation.h"
0002 #include "DataFormats/ForwardDetId/interface/ForwardSubdetector.h"
0003 #include "DataFormats/ForwardDetId/interface/HFNoseTriggerDetId.h"
0004 #include "DataFormats/ForwardDetId/interface/HGCScintillatorDetId.h"
0005 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetIdToROC.h"
0006 #include "DataFormats/ForwardDetId/interface/HGCalDetId.h"
0007 #include "DataFormats/ForwardDetId/interface/HGCalTriggerDetId.h"
0008 #include "FWCore/ParameterSet/interface/FileInPath.h"
0009 #include "L1Trigger/L1THGCal/interface/HGCalTriggerGeometryBase.h"
0010 #include "DataFormats/ForwardDetId/interface/HGCalTriggerModuleDetId.h"
0011 #include "DataFormats/ForwardDetId/interface/HGCalTriggerBackendDetId.h"
0012 #include "DataFormats/ForwardDetId/interface/HFNoseDetIdToModule.h"
0013 
0014 #include <fstream>
0015 #include <vector>
0016 
0017 #include <nlohmann/json.hpp>
0018 using json = nlohmann::json;
0019 
0020 class HGCalTriggerGeometryV9Imp3 : public HGCalTriggerGeometryBase {
0021 public:
0022   HGCalTriggerGeometryV9Imp3(const edm::ParameterSet& conf);
0023 
0024   void initialize(const HGCalGeometry*, const HGCalGeometry*, const HGCalGeometry*) final;
0025   void initialize(const HGCalGeometry*, const HGCalGeometry*, const HGCalGeometry*, const HGCalGeometry*) final;
0026   void reset() final;
0027 
0028   unsigned getTriggerCellFromCell(const unsigned) const final;
0029   unsigned getModuleFromCell(const unsigned) const final;
0030   unsigned getModuleFromTriggerCell(const unsigned) const final;
0031 
0032   geom_set getCellsFromTriggerCell(const unsigned) const final;
0033   geom_set getCellsFromModule(const unsigned) const final;
0034   geom_set getTriggerCellsFromModule(const unsigned) const final;
0035 
0036   geom_ordered_set getOrderedCellsFromModule(const unsigned) const final;
0037   geom_ordered_set getOrderedTriggerCellsFromModule(const unsigned) const final;
0038 
0039   geom_set getNeighborsFromTriggerCell(const unsigned) const final;
0040 
0041   geom_set getStage1FpgasFromStage2Fpga(const unsigned) const final;
0042   geom_set getStage2FpgasFromStage1Fpga(const unsigned) const final;
0043 
0044   geom_set getStage1LinksFromStage2Fpga(const unsigned) const final;
0045   unsigned getStage1FpgaFromStage1Link(const unsigned) const final;
0046   unsigned getStage2FpgaFromStage1Link(const unsigned) const final;
0047   geom_set getStage1LinksFromStage1Fpga(const unsigned) const final;
0048   std::vector<unsigned> getLpgbtsFromStage1Fpga(const unsigned) const final;
0049   unsigned getStage1FpgaFromLpgbt(const unsigned) const final;
0050   geom_set getModulesFromLpgbt(const unsigned) const final;
0051   geom_set getLpgbtsFromModule(const unsigned) const final;
0052   unsigned getStage1FpgaFromModule(const unsigned module_id) const final;
0053 
0054   unsigned getLinksInModule(const unsigned module_id) const final;
0055   unsigned getModuleSize(const unsigned module_id) const final;
0056 
0057   GlobalPoint getTriggerCellPosition(const unsigned) const final;
0058   GlobalPoint getModulePosition(const unsigned) const final;
0059 
0060   bool validCell(const unsigned) const final;
0061   bool validTriggerCell(const unsigned) const final;
0062   bool disconnectedModule(const unsigned) const final;
0063   unsigned lastTriggerLayer() const final { return last_trigger_layer_; }
0064   unsigned triggerLayer(const unsigned) const final;
0065 
0066 private:
0067   // HSc trigger cell grouping
0068   unsigned hSc_triggercell_size_ = 2;
0069   static constexpr unsigned hSc_num_panels_per_sector_ = 12;
0070   static constexpr unsigned hSc_tcs_per_module_phi_ = 4;
0071   static constexpr unsigned hSc_front_layers_split_ = 12;
0072   static constexpr unsigned hSc_back_layers_split_ = 8;
0073   static constexpr unsigned hSc_layer_for_split_ = 40;
0074   static constexpr int hSc_tc_layer0_min_ = 24;
0075   static constexpr int ntc_per_wafer_ = 48;
0076   static constexpr int nSectors_ = 3;
0077 
0078   edm::FileInPath jsonMappingFile_;
0079 
0080   // rotation class
0081   HGCalGeomRotation geom_rotation_120_ = {HGCalGeomRotation::SectorType::Sector120Degrees};
0082 
0083   // module related maps
0084   std::unordered_map<unsigned, unsigned> links_per_module_;
0085 
0086   std::unordered_multimap<unsigned, unsigned> stage2_to_stage1links_;
0087   std::unordered_map<unsigned, bool> stage1links_samesector_;
0088   std::unordered_map<unsigned, unsigned> stage1link_to_stage2_;
0089   std::unordered_map<unsigned, unsigned> stage1link_to_stage1_;
0090   std::unordered_multimap<unsigned, unsigned> stage1_to_stage1links_;
0091   std::unordered_map<unsigned, std::vector<unsigned>> stage1_to_lpgbts_;
0092   std::unordered_map<unsigned, unsigned> lpgbt_to_stage1_;
0093   std::unordered_multimap<unsigned, unsigned> lpgbt_to_modules_;
0094   std::unordered_multimap<unsigned, unsigned> module_to_lpgbts_;
0095   std::unordered_map<unsigned, unsigned> module_to_stage1_;
0096 
0097   // Disconnected modules and layers
0098   std::unordered_set<unsigned> disconnected_layers_;
0099   std::vector<unsigned> trigger_layers_;
0100   std::vector<unsigned> trigger_nose_layers_;
0101   unsigned last_trigger_layer_ = 0;
0102 
0103   // layer offsets
0104   unsigned heOffset_ = 0;
0105   unsigned noseLayers_ = 0;
0106   unsigned totalLayers_ = 0;
0107 
0108   void fillMaps();
0109   bool validCellId(unsigned det, unsigned cell_id) const;
0110   bool validTriggerCellFromCells(const unsigned) const;
0111 
0112   int detIdWaferType(unsigned det, unsigned layer, short waferU, short waferV) const;
0113   void layerWithoutOffsetAndSubdetId(unsigned& layer, int& subdetId, bool isSilicon) const;
0114   unsigned packLayerSubdetWaferId(unsigned layer, int subdet, int waferU, int waferV) const;
0115   void unpackLayerSubdetWaferId(unsigned wafer, unsigned& layer, int& subdet, int& waferU, int& waferV) const;
0116   HGCalGeomRotation::WaferCentring getWaferCentring(unsigned layer, int subdet) const;
0117   void etaphiMappingFromSector0(int& ieta, int& iphi, unsigned sector) const;
0118   unsigned tcEtaphiMappingToSector0(int& tc_ieta, int& tc_iphi) const;
0119   void getScintillatoriEtaiPhi(int& ieta, int& iphi, int tc_eta, int tc_phi, unsigned layer) const;
0120   unsigned layerWithOffset(unsigned) const;
0121   unsigned getNextSector(const unsigned sector) const;
0122   unsigned getPreviousSector(const unsigned sector) const;
0123 };
0124 
0125 HGCalTriggerGeometryV9Imp3::HGCalTriggerGeometryV9Imp3(const edm::ParameterSet& conf)
0126     : HGCalTriggerGeometryBase(conf),
0127       hSc_triggercell_size_(conf.getParameter<unsigned>("ScintillatorTriggerCellSize")),
0128       jsonMappingFile_(conf.getParameter<edm::FileInPath>("JsonMappingFile")) {
0129   std::vector<unsigned> tmp_vector = conf.getParameter<std::vector<unsigned>>("DisconnectedLayers");
0130   std::move(tmp_vector.begin(), tmp_vector.end(), std::inserter(disconnected_layers_, disconnected_layers_.end()));
0131 }
0132 
0133 void HGCalTriggerGeometryV9Imp3::reset() {
0134   stage2_to_stage1links_.clear();
0135   stage1links_samesector_.clear();
0136   stage1link_to_stage2_.clear();
0137   stage1link_to_stage1_.clear();
0138   stage1_to_stage1links_.clear();
0139   stage1_to_lpgbts_.clear();
0140   lpgbt_to_stage1_.clear();
0141   lpgbt_to_modules_.clear();
0142   module_to_lpgbts_.clear();
0143   module_to_stage1_.clear();
0144 }
0145 
0146 void HGCalTriggerGeometryV9Imp3::initialize(const HGCalGeometry* hgc_ee_geometry,
0147                                             const HGCalGeometry* hgc_hsi_geometry,
0148                                             const HGCalGeometry* hgc_hsc_geometry) {
0149   setEEGeometry(hgc_ee_geometry);
0150   setHSiGeometry(hgc_hsi_geometry);
0151   setHScGeometry(hgc_hsc_geometry);
0152   heOffset_ = eeTopology().dddConstants().layers(true);
0153   totalLayers_ = heOffset_ + hsiTopology().dddConstants().layers(true);
0154   trigger_layers_.resize(totalLayers_ + 1);
0155   trigger_layers_[0] = 0;  // layer number 0 doesn't exist
0156   unsigned trigger_layer = 1;
0157   for (unsigned layer = 1; layer < trigger_layers_.size(); layer++) {
0158     if (disconnected_layers_.find(layer) == disconnected_layers_.end()) {
0159       // Increase trigger layer number if the layer is not disconnected
0160       trigger_layers_[layer] = trigger_layer;
0161       trigger_layer++;
0162     } else {
0163       trigger_layers_[layer] = 0;
0164     }
0165   }
0166   last_trigger_layer_ = trigger_layer - 1;
0167   fillMaps();
0168 }
0169 
0170 void HGCalTriggerGeometryV9Imp3::initialize(const HGCalGeometry* hgc_ee_geometry,
0171                                             const HGCalGeometry* hgc_hsi_geometry,
0172                                             const HGCalGeometry* hgc_hsc_geometry,
0173                                             const HGCalGeometry* hgc_nose_geometry) {
0174   setEEGeometry(hgc_ee_geometry);
0175   setHSiGeometry(hgc_hsi_geometry);
0176   setHScGeometry(hgc_hsc_geometry);
0177   setNoseGeometry(hgc_nose_geometry);
0178 
0179   heOffset_ = eeTopology().dddConstants().layers(true);
0180   totalLayers_ = heOffset_ + hsiTopology().dddConstants().layers(true);
0181 
0182   trigger_layers_.resize(totalLayers_ + 1);
0183   trigger_layers_[0] = 0;  // layer number 0 doesn't exist
0184   unsigned trigger_layer = 1;
0185   for (unsigned layer = 1; layer < trigger_layers_.size(); layer++) {
0186     if (disconnected_layers_.find(layer) == disconnected_layers_.end()) {
0187       // Increase trigger layer number if the layer is not disconnected
0188       trigger_layers_[layer] = trigger_layer;
0189       trigger_layer++;
0190     } else {
0191       trigger_layers_[layer] = 0;
0192     }
0193   }
0194   last_trigger_layer_ = trigger_layer - 1;
0195   fillMaps();
0196 
0197   noseLayers_ = noseTopology().dddConstants().layers(true);
0198 
0199   trigger_nose_layers_.resize(noseLayers_ + 1);
0200   trigger_nose_layers_[0] = 0;  // layer number 0 doesn't exist
0201   unsigned trigger_nose_layer = 1;
0202   for (unsigned layer = 1; layer < trigger_nose_layers_.size(); layer++) {
0203     trigger_nose_layers_[layer] = trigger_nose_layer;
0204     trigger_nose_layer++;
0205   }
0206 }
0207 
0208 unsigned HGCalTriggerGeometryV9Imp3::getTriggerCellFromCell(const unsigned cell_id) const {
0209   unsigned det = DetId(cell_id).det();
0210   unsigned trigger_cell_id = 0;
0211   // Scintillator
0212   if (det == DetId::HGCalHSc) {
0213     // Very rough mapping from cells to TC
0214     HGCScintillatorDetId cell_sc_id(cell_id);
0215     int ieta = ((cell_sc_id.ietaAbs() - 1) / hSc_triggercell_size_ + 1) * cell_sc_id.zside();
0216     int iphi = (cell_sc_id.iphi() - 1) / hSc_triggercell_size_ + 1;
0217     trigger_cell_id = HGCScintillatorDetId(cell_sc_id.type(), cell_sc_id.layer(), ieta, iphi);
0218   }
0219   // HFNose
0220   else if (det == DetId::Forward && DetId(cell_id).subdetId() == ForwardSubdetector::HFNose) {
0221     HFNoseDetId cell_nose_id(cell_id);
0222     trigger_cell_id = HFNoseTriggerDetId(HGCalTriggerSubdetector::HFNoseTrigger,
0223                                          cell_nose_id.zside(),
0224                                          cell_nose_id.type(),
0225                                          cell_nose_id.layer(),
0226                                          cell_nose_id.waferU(),
0227                                          cell_nose_id.waferV(),
0228                                          cell_nose_id.triggerCellU(),
0229                                          cell_nose_id.triggerCellV());
0230   }
0231   // Silicon
0232   else if (det == DetId::HGCalEE || det == DetId::HGCalHSi) {
0233     HGCSiliconDetId cell_si_id(cell_id);
0234     trigger_cell_id = HGCalTriggerDetId(
0235         det == DetId::HGCalEE ? HGCalTriggerSubdetector::HGCalEETrigger : HGCalTriggerSubdetector::HGCalHSiTrigger,
0236         cell_si_id.zside(),
0237         cell_si_id.type(),
0238         cell_si_id.layer(),
0239         cell_si_id.waferU(),
0240         cell_si_id.waferV(),
0241         cell_si_id.triggerCellU(),
0242         cell_si_id.triggerCellV());
0243   }
0244   return trigger_cell_id;
0245 }
0246 
0247 unsigned HGCalTriggerGeometryV9Imp3::getModuleFromCell(const unsigned cell_id) const {
0248   return getModuleFromTriggerCell(getTriggerCellFromCell(cell_id));
0249 }
0250 
0251 unsigned HGCalTriggerGeometryV9Imp3::getModuleFromTriggerCell(const unsigned trigger_cell_id) const {
0252   unsigned det = DetId(trigger_cell_id).det();
0253   int zside = 0;
0254   unsigned tc_type = 1;
0255   unsigned layer = 0;
0256   unsigned module_id = 0;
0257 
0258   // Scintillator
0259   if (det == DetId::HGCalHSc) {
0260     HGCScintillatorDetId trigger_cell_sc_id(trigger_cell_id);
0261     tc_type = trigger_cell_sc_id.type();
0262     layer = trigger_cell_sc_id.layer();
0263     zside = trigger_cell_sc_id.zside();
0264     int tc_eta = trigger_cell_sc_id.ietaAbs();
0265     int tc_phi = trigger_cell_sc_id.iphi();
0266     unsigned sector = tcEtaphiMappingToSector0(tc_eta, tc_phi);
0267     int ieta = 0;
0268     int iphi = 0;
0269     getScintillatoriEtaiPhi(ieta, iphi, tc_eta, tc_phi, layer);
0270     module_id =
0271         HGCalTriggerModuleDetId(HGCalTriggerSubdetector::HGCalHScTrigger, zside, tc_type, layer, sector, ieta, iphi);
0272   }
0273   // HFNose
0274   else if (det == DetId::HGCalTrigger and
0275            HGCalTriggerDetId(trigger_cell_id).subdet() == HGCalTriggerSubdetector::HFNoseTrigger) {
0276     HFNoseTriggerDetId trigger_cell_trig_id(trigger_cell_id);
0277     tc_type = trigger_cell_trig_id.type();
0278     layer = trigger_cell_trig_id.layer();
0279     zside = trigger_cell_trig_id.zside();
0280     int waferu = trigger_cell_trig_id.waferU();
0281     int waferv = trigger_cell_trig_id.waferV();
0282     unsigned sector = geom_rotation_120_.uvMappingToSector0(
0283         getWaferCentring(layer, HGCalTriggerSubdetector::HFNoseTrigger), waferu, waferv);
0284     module_id =
0285         HGCalTriggerModuleDetId(HGCalTriggerSubdetector::HFNoseTrigger, zside, tc_type, layer, sector, waferu, waferv);
0286   }
0287   // Silicon
0288   else {
0289     HGCalTriggerDetId trigger_cell_trig_id(trigger_cell_id);
0290     unsigned subdet = trigger_cell_trig_id.subdet();
0291     tc_type = trigger_cell_trig_id.type();
0292     layer = trigger_cell_trig_id.layer();
0293     zside = trigger_cell_trig_id.zside();
0294     int waferu = trigger_cell_trig_id.waferU();
0295     int waferv = trigger_cell_trig_id.waferV();
0296     unsigned sector = geom_rotation_120_.uvMappingToSector0(getWaferCentring(layer, subdet), waferu, waferv);
0297     module_id = HGCalTriggerModuleDetId(HGCalTriggerSubdetector(subdet), zside, tc_type, layer, sector, waferu, waferv);
0298   }
0299   return module_id;
0300 }
0301 
0302 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getCellsFromTriggerCell(
0303     const unsigned trigger_cell_id) const {
0304   DetId trigger_cell_det_id(trigger_cell_id);
0305   unsigned det = trigger_cell_det_id.det();
0306   geom_set cell_det_ids;
0307   // Scintillator
0308   if (det == DetId::HGCalHSc) {
0309     HGCScintillatorDetId trigger_cell_sc_id(trigger_cell_id);
0310     int ieta0 = (trigger_cell_sc_id.ietaAbs() - 1) * hSc_triggercell_size_ + 1;
0311     int iphi0 = (trigger_cell_sc_id.iphi() - 1) * hSc_triggercell_size_ + 1;
0312     for (int ietaAbs = ieta0; ietaAbs < ieta0 + (int)hSc_triggercell_size_; ietaAbs++) {
0313       int ieta = ietaAbs * trigger_cell_sc_id.zside();
0314       for (int iphi = iphi0; iphi < iphi0 + (int)hSc_triggercell_size_; iphi++) {
0315         unsigned cell_id = HGCScintillatorDetId(trigger_cell_sc_id.type(), trigger_cell_sc_id.layer(), ieta, iphi);
0316         if (validCellId(DetId::HGCalHSc, cell_id))
0317           cell_det_ids.emplace(cell_id);
0318       }
0319     }
0320   }
0321   // HFNose
0322   else if (det == DetId::HGCalTrigger and
0323            HGCalTriggerDetId(trigger_cell_id).subdet() == HGCalTriggerSubdetector::HFNoseTrigger) {
0324     HFNoseTriggerDetId trigger_cell_nose_id(trigger_cell_id);
0325     int layer = trigger_cell_nose_id.layer();
0326     int zside = trigger_cell_nose_id.zside();
0327     int type = trigger_cell_nose_id.type();
0328     int waferu = trigger_cell_nose_id.waferU();
0329     int waferv = trigger_cell_nose_id.waferV();
0330     std::vector<int> cellus = trigger_cell_nose_id.cellU();
0331     std::vector<int> cellvs = trigger_cell_nose_id.cellV();
0332     for (unsigned ic = 0; ic < cellus.size(); ic++) {
0333       HFNoseDetId cell_det_id(zside, type, layer, waferu, waferv, cellus[ic], cellvs[ic]);
0334       cell_det_ids.emplace(cell_det_id);
0335     }
0336   }
0337   // Silicon
0338   else {
0339     HGCalTriggerDetId trigger_cell_trig_id(trigger_cell_id);
0340     unsigned subdet = trigger_cell_trig_id.subdet();
0341     if (subdet == HGCalTriggerSubdetector::HGCalEETrigger || subdet == HGCalTriggerSubdetector::HGCalHSiTrigger) {
0342       DetId::Detector cell_det = (subdet == HGCalTriggerSubdetector::HGCalEETrigger ? DetId::HGCalEE : DetId::HGCalHSi);
0343       int layer = trigger_cell_trig_id.layer();
0344       int zside = trigger_cell_trig_id.zside();
0345       int type = trigger_cell_trig_id.type();
0346       int waferu = trigger_cell_trig_id.waferU();
0347       int waferv = trigger_cell_trig_id.waferV();
0348       std::vector<int> cellus = trigger_cell_trig_id.cellU();
0349       std::vector<int> cellvs = trigger_cell_trig_id.cellV();
0350       for (unsigned ic = 0; ic < cellus.size(); ic++) {
0351         HGCSiliconDetId cell_det_id(cell_det, zside, type, layer, waferu, waferv, cellus[ic], cellvs[ic]);
0352         cell_det_ids.emplace(cell_det_id);
0353       }
0354     }
0355   }
0356   return cell_det_ids;
0357 }
0358 
0359 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getCellsFromModule(const unsigned module_id) const {
0360   geom_set cell_det_ids;
0361   geom_set trigger_cells = getTriggerCellsFromModule(module_id);
0362 
0363   for (auto trigger_cell_id : trigger_cells) {
0364     geom_set cells = getCellsFromTriggerCell(trigger_cell_id);
0365     cell_det_ids.insert(cells.begin(), cells.end());
0366   }
0367   return cell_det_ids;
0368 }
0369 
0370 HGCalTriggerGeometryBase::geom_ordered_set HGCalTriggerGeometryV9Imp3::getOrderedCellsFromModule(
0371     const unsigned module_id) const {
0372   geom_ordered_set cell_det_ids;
0373   geom_ordered_set trigger_cells = getOrderedTriggerCellsFromModule(module_id);
0374   for (auto trigger_cell_id : trigger_cells) {
0375     geom_set cells = getCellsFromTriggerCell(trigger_cell_id);
0376     cell_det_ids.insert(cells.begin(), cells.end());
0377   }
0378   return cell_det_ids;
0379 }
0380 
0381 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getTriggerCellsFromModule(
0382     const unsigned module_id) const {
0383   HGCalTriggerModuleDetId hgc_module_id(module_id);
0384   unsigned subdet = hgc_module_id.triggerSubdetId();
0385 
0386   geom_set trigger_cell_det_ids;
0387   // Scintillator
0388   if (subdet == HGCalTriggerSubdetector::HGCalHScTrigger) {
0389     int ieta0 = hgc_module_id.eta();
0390     int iphi0 = hgc_module_id.phi();
0391 
0392     unsigned layer = hgc_module_id.layer();
0393     etaphiMappingFromSector0(ieta0, iphi0, hgc_module_id.sector());
0394     int split = (layer > hSc_layer_for_split_) ? hSc_back_layers_split_ : hSc_front_layers_split_;
0395     if (ieta0 == 1) {
0396       ieta0 = ieta0 + split;
0397     } else {
0398       ieta0 = ieta0 + 1;
0399     }
0400     iphi0 = (iphi0 * hSc_tcs_per_module_phi_) + hSc_tc_layer0_min_ + 1;
0401     int total_tcs = hSc_num_panels_per_sector_ * hSc_tcs_per_module_phi_ * nSectors_;
0402     if (iphi0 > total_tcs) {
0403       iphi0 = iphi0 - total_tcs;
0404     }
0405 
0406     int hSc_tcs_per_module_eta = (layer > hSc_layer_for_split_) ? hSc_back_layers_split_ : hSc_front_layers_split_;
0407 
0408     for (int ietaAbs = ieta0; ietaAbs < ieta0 + (int)hSc_tcs_per_module_eta; ietaAbs++) {
0409       int ieta = ietaAbs * hgc_module_id.zside();
0410       for (int iphi = iphi0; iphi < iphi0 + (int)hSc_tcs_per_module_phi_; iphi++) {
0411         unsigned trigger_cell_id = HGCScintillatorDetId(hgc_module_id.type(), hgc_module_id.layer(), ieta, iphi);
0412         if (validTriggerCellFromCells(trigger_cell_id))
0413           trigger_cell_det_ids.emplace(trigger_cell_id);
0414       }
0415     }
0416   }
0417   // HFNose
0418   else if (subdet == HGCalTriggerSubdetector::HFNoseTrigger) {
0419     HFNoseDetId module_nose_id(module_id);
0420     HFNoseDetIdToModule hfn;
0421     std::vector<HFNoseTriggerDetId> ids = hfn.getTriggerDetIds(module_nose_id);
0422     for (auto const& idx : ids) {
0423       if (validTriggerCellFromCells(idx.rawId()))
0424         trigger_cell_det_ids.emplace(idx);
0425     }
0426   }
0427   // Silicon
0428   else {
0429     HGCSiliconDetIdToROC tc2roc;
0430     int moduleU = hgc_module_id.moduleU();
0431     int moduleV = hgc_module_id.moduleV();
0432     unsigned layer = hgc_module_id.layer();
0433 
0434     //Rotate to sector
0435     geom_rotation_120_.uvMappingFromSector0(getWaferCentring(layer, subdet), moduleU, moduleV, hgc_module_id.sector());
0436 
0437     DetId::Detector det = (subdet == HGCalTriggerSubdetector::HGCalEETrigger ? DetId::HGCalEE : DetId::HGCalHSi);
0438 
0439     unsigned wafer_type = detIdWaferType(det, layer, moduleU, moduleV);
0440     int nroc = (((wafer_type == HGCSiliconDetId::HGCalHD120) || (wafer_type == HGCSiliconDetId::HGCalHD200)) ? 6 : 3);
0441     // Loop on ROCs in wafer
0442     for (int roc = 1; roc <= nroc; roc++) {
0443       // loop on TCs in ROC
0444       auto tc_uvs = tc2roc.getTriggerId(roc, wafer_type);
0445       for (const auto& tc_uv : tc_uvs) {
0446         HGCalTriggerDetId trigger_cell_id(
0447             subdet, hgc_module_id.zside(), wafer_type, layer, moduleU, moduleV, tc_uv.first, tc_uv.second);
0448         if (validTriggerCellFromCells(trigger_cell_id.rawId()))
0449           trigger_cell_det_ids.emplace(trigger_cell_id);
0450       }
0451     }
0452   }
0453 
0454   return trigger_cell_det_ids;
0455 }
0456 
0457 HGCalTriggerGeometryBase::geom_ordered_set HGCalTriggerGeometryV9Imp3::getOrderedTriggerCellsFromModule(
0458     const unsigned module_id) const {
0459   HGCalTriggerModuleDetId hgc_module_id(module_id);
0460   unsigned subdet = hgc_module_id.triggerSubdetId();
0461 
0462   geom_ordered_set trigger_cell_det_ids;
0463 
0464   // Scintillator
0465   if (subdet == HGCalTriggerSubdetector::HGCalHScTrigger) {
0466     int ieta0 = hgc_module_id.eta();
0467     int iphi0 = hgc_module_id.phi();
0468 
0469     unsigned layer = hgc_module_id.layer();
0470     etaphiMappingFromSector0(ieta0, iphi0, hgc_module_id.sector());
0471     int split = (layer > hSc_layer_for_split_) ? hSc_back_layers_split_ : hSc_front_layers_split_;
0472     if (ieta0 == 1) {
0473       ieta0 = ieta0 + split;
0474     } else {
0475       ieta0 = ieta0 + 1;
0476     }
0477     iphi0 = (iphi0 * hSc_tcs_per_module_phi_) + hSc_tc_layer0_min_ + 1;
0478     int total_tcs = hSc_num_panels_per_sector_ * hSc_tcs_per_module_phi_ * nSectors_;
0479     if (iphi0 > total_tcs) {
0480       iphi0 = iphi0 - total_tcs;
0481     }
0482 
0483     int hSc_tcs_per_module_eta = (layer > hSc_layer_for_split_) ? hSc_back_layers_split_ : hSc_front_layers_split_;
0484 
0485     for (int ietaAbs = ieta0; ietaAbs < ieta0 + (int)hSc_tcs_per_module_eta; ietaAbs++) {
0486       int ieta = ietaAbs * hgc_module_id.zside();
0487       for (int iphi = iphi0; iphi < iphi0 + (int)hSc_tcs_per_module_phi_; iphi++) {
0488         unsigned trigger_cell_id = HGCScintillatorDetId(hgc_module_id.type(), hgc_module_id.layer(), ieta, iphi);
0489         if (validTriggerCellFromCells(trigger_cell_id))
0490           trigger_cell_det_ids.emplace(trigger_cell_id);
0491       }
0492     }
0493   }
0494 
0495   // HFNose
0496   else if (subdet == HGCalTriggerSubdetector::HFNoseTrigger) {
0497     HFNoseDetId module_nose_id(module_id);
0498     HFNoseDetIdToModule hfn;
0499     std::vector<HFNoseTriggerDetId> ids = hfn.getTriggerDetIds(module_nose_id);
0500     for (auto const& idx : ids) {
0501       if (validTriggerCellFromCells(idx.rawId()))
0502         trigger_cell_det_ids.emplace(idx);
0503     }
0504   }
0505   // Silicon
0506   else {
0507     HGCSiliconDetIdToROC tc2roc;
0508     int moduleU = hgc_module_id.moduleU();
0509     int moduleV = hgc_module_id.moduleV();
0510     unsigned layer = hgc_module_id.layer();
0511 
0512     //Rotate to sector
0513     geom_rotation_120_.uvMappingFromSector0(getWaferCentring(layer, subdet), moduleU, moduleV, hgc_module_id.sector());
0514 
0515     DetId::Detector det = (subdet == HGCalTriggerSubdetector::HGCalEETrigger ? DetId::HGCalEE : DetId::HGCalHSi);
0516 
0517     unsigned wafer_type = detIdWaferType(det, layer, moduleU, moduleV);
0518     int nroc = (((wafer_type == HGCSiliconDetId::HGCalHD120) || (wafer_type == HGCSiliconDetId::HGCalHD200)) ? 6 : 3);
0519     // Loop on ROCs in wafer
0520     for (int roc = 1; roc <= nroc; roc++) {
0521       // loop on TCs in ROC
0522       auto tc_uvs = tc2roc.getTriggerId(roc, wafer_type);
0523       for (const auto& tc_uv : tc_uvs) {
0524         HGCalTriggerDetId trigger_cell_id(
0525             subdet, hgc_module_id.zside(), wafer_type, layer, moduleU, moduleV, tc_uv.first, tc_uv.second);
0526         if (validTriggerCellFromCells(trigger_cell_id.rawId()))
0527           trigger_cell_det_ids.emplace(trigger_cell_id);
0528       }
0529     }
0530   }
0531 
0532   return trigger_cell_det_ids;
0533 }
0534 
0535 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getNeighborsFromTriggerCell(
0536     const unsigned trigger_cell_id) const {
0537   throw cms::Exception("FeatureNotImplemented") << "Neighbor search is not implemented in HGCalTriggerGeometryV9Imp3";
0538 }
0539 
0540 unsigned HGCalTriggerGeometryV9Imp3::getLinksInModule(const unsigned module_id) const {
0541   HGCalTriggerModuleDetId module_det_id(module_id);
0542   unsigned subdet = module_det_id.triggerSubdetId();
0543 
0544   unsigned links = 0;
0545   // HF Nose
0546   if (subdet == HGCalTriggerSubdetector::HFNoseTrigger) {
0547     links = 1;
0548   }
0549   // TO ADD HFNOSE : getLinksInModule
0550   // Silicon and Scintillator
0551   else {
0552     int packed_module =
0553         packLayerSubdetWaferId(module_det_id.layer(), subdet, module_det_id.moduleU(), module_det_id.moduleV());
0554     links = links_per_module_.at(packed_module);
0555   }
0556   return links;
0557 }
0558 
0559 unsigned HGCalTriggerGeometryV9Imp3::getModuleSize(const unsigned module_id) const {
0560   unsigned nWafers = 1;
0561   return nWafers;
0562 }
0563 
0564 unsigned HGCalTriggerGeometryV9Imp3::getNextSector(const unsigned sector) const {
0565   unsigned next_sector = 0;
0566   if (sector < 2) {
0567     next_sector = sector + 1;
0568   }
0569   return next_sector;
0570 }
0571 
0572 unsigned HGCalTriggerGeometryV9Imp3::getPreviousSector(const unsigned sector) const {
0573   unsigned previous_sector = 2;
0574   if (sector > 0) {
0575     previous_sector = sector - 1;
0576   }
0577   return previous_sector;
0578 }
0579 
0580 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getStage1FpgasFromStage2Fpga(
0581     const unsigned stage2_id) const {
0582   geom_set stage1_ids;
0583 
0584   geom_set stage1_links = getStage1LinksFromStage2Fpga(stage2_id);
0585   for (const auto& stage1_link : stage1_links) {
0586     stage1_ids.emplace(getStage1FpgaFromStage1Link(stage1_link));
0587   }
0588 
0589   return stage1_ids;
0590 }
0591 
0592 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getStage2FpgasFromStage1Fpga(
0593     const unsigned stage1_id) const {
0594   geom_set stage2_ids;
0595 
0596   geom_set stage1_links = getStage1LinksFromStage1Fpga(stage1_id);
0597   for (const auto& stage1_link : stage1_links) {
0598     stage2_ids.emplace(getStage2FpgaFromStage1Link(stage1_link));
0599   }
0600 
0601   return stage2_ids;
0602 }
0603 
0604 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getStage1LinksFromStage2Fpga(
0605     const unsigned stage2_id) const {
0606   geom_set stage1link_ids;
0607   HGCalTriggerBackendDetId id(stage2_id);
0608   auto stage2_itrs = stage2_to_stage1links_.equal_range(id.label());
0609   for (auto stage2_itr = stage2_itrs.first; stage2_itr != stage2_itrs.second; stage2_itr++) {
0610     unsigned label = stage2_itr->second;
0611     if (stage1links_samesector_.at(label) == true) {  //link and stage2 FPGA are the same sector
0612       stage1link_ids.emplace(
0613           HGCalTriggerBackendDetId(id.zside(), HGCalTriggerBackendDetId::BackendType::Stage1Link, id.sector(), label));
0614     } else {  //link is from the next sector (anti-clockwise)
0615       stage1link_ids.emplace(HGCalTriggerBackendDetId(
0616           id.zside(), HGCalTriggerBackendDetId::BackendType::Stage1Link, getNextSector(id.sector()), label));
0617     }
0618   }
0619 
0620   return stage1link_ids;
0621 }
0622 
0623 unsigned HGCalTriggerGeometryV9Imp3::getStage1FpgaFromStage1Link(const unsigned link_id) const {
0624   HGCalTriggerBackendDetId id(link_id);
0625   unsigned stage1_label = stage1link_to_stage1_.at(id.label());
0626 
0627   return HGCalTriggerBackendDetId(
0628       id.zside(), HGCalTriggerBackendDetId::BackendType::Stage1FPGA, id.sector(), stage1_label);
0629 }
0630 
0631 unsigned HGCalTriggerGeometryV9Imp3::getStage2FpgaFromStage1Link(const unsigned link_id) const {
0632   HGCalTriggerBackendDetId id(link_id);
0633   bool same_sector = stage1link_to_stage2_.at(id.label());
0634   unsigned sector = id.sector();
0635 
0636   if (!same_sector) {
0637     sector = getPreviousSector(sector);
0638   }
0639 
0640   return HGCalTriggerBackendDetId(id.zside(), HGCalTriggerBackendDetId::BackendType::Stage2FPGA, sector, 0);
0641 }
0642 
0643 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getStage1LinksFromStage1Fpga(
0644     const unsigned stage1_id) const {
0645   geom_set stage1link_ids;
0646   HGCalTriggerBackendDetId id(stage1_id);
0647 
0648   auto stage1_itrs = stage1_to_stage1links_.equal_range(id.label());
0649   for (auto stage1_itr = stage1_itrs.first; stage1_itr != stage1_itrs.second; stage1_itr++) {
0650     stage1link_ids.emplace(HGCalTriggerBackendDetId(
0651         id.zside(), HGCalTriggerBackendDetId::BackendType::Stage1Link, id.sector(), stage1_itr->second));
0652   }
0653 
0654   return stage1link_ids;
0655 }
0656 
0657 std::vector<unsigned> HGCalTriggerGeometryV9Imp3::getLpgbtsFromStage1Fpga(const unsigned stage1_id) const {
0658   std::vector<unsigned> lpgbt_ids;
0659   HGCalTriggerBackendDetId id(stage1_id);
0660 
0661   const auto stage1_lpgbts = stage1_to_lpgbts_.at(id.label());
0662   lpgbt_ids.reserve(stage1_lpgbts.size());
0663   for (const auto& stage1_lpgbt : stage1_lpgbts) {
0664     lpgbt_ids.emplace_back(
0665         HGCalTriggerBackendDetId(id.zside(), HGCalTriggerBackendDetId::BackendType::LpGBT, id.sector(), stage1_lpgbt));
0666   }
0667 
0668   return lpgbt_ids;
0669 }
0670 
0671 unsigned HGCalTriggerGeometryV9Imp3::getStage1FpgaFromLpgbt(const unsigned lpgbt_id) const {
0672   HGCalTriggerBackendDetId id(lpgbt_id);
0673   unsigned stage1_label = lpgbt_to_stage1_.at(id.label());
0674 
0675   return HGCalTriggerBackendDetId(
0676       id.zside(), HGCalTriggerBackendDetId::BackendType::Stage1FPGA, id.sector(), stage1_label);
0677 }
0678 
0679 HGCalTriggerGeometryBase::geom_set HGCalTriggerGeometryV9Imp3::getModulesFromLpgbt(const unsigned lpgbt_id) const {
0680   geom_set modules;
0681   HGCalTriggerBackendDetId id(lpgbt_id);
0682 
0683   auto lpgbt_itrs = lpgbt_to_modules_.equal_range(id.label());
0684   for (auto lpgbt_itr = lpgbt_itrs.first; lpgbt_itr != lpgbt_itrs.second; lpgbt_itr++) {
0685     unsigned layer = 0;
0686     int moduleU = 0;
0687     int moduleV = 0;
0688     int subdet = 0;
0689     unpackLayerSubdetWaferId(lpgbt_itr->second, layer, subdet, moduleU, moduleV);
0690     unsigned det = 0;
0691     switch (subdet) {
0692       case HGCalTriggerSubdetector::HGCalEETrigger:
0693         det = DetId::HGCalEE;
0694         break;
0695       case HGCalTriggerSubdetector::HGCalHSiTrigger:
0696         det = DetId::HGCalHSi;
0697         break;
0698       case HGCalTriggerSubdetector::HGCalHScTrigger:
0699         det = DetId::HGCalHSc;
0700         break;
0701       default:
0702         det = DetId::HGCalEE;
0703         break;
0704     }
0705 
0706     int type = detIdWaferType(det, layer, moduleU, moduleV);
0707     modules.emplace(HGCalTriggerModuleDetId(
0708         HGCalTriggerSubdetector(subdet), id.zside(), type, layer, id.sector(), moduleU, moduleV));
0709   }
0710 
0711   return modules;
0712 }
0713 
0714 HGCalTriggerGeometryV9Imp3::geom_set HGCalTriggerGeometryV9Imp3::getLpgbtsFromModule(const unsigned module_id) const {
0715   geom_set lpgbt_ids;
0716   HGCalTriggerModuleDetId id(module_id);
0717 
0718   auto module_itrs = module_to_lpgbts_.equal_range(
0719       packLayerSubdetWaferId(id.layer(), id.triggerSubdetId(), id.moduleU(), id.moduleV()));
0720   for (auto module_itr = module_itrs.first; module_itr != module_itrs.second; module_itr++) {
0721     lpgbt_ids.emplace(HGCalTriggerBackendDetId(
0722         id.zside(), HGCalTriggerBackendDetId::BackendType::LpGBT, id.sector(), module_itr->second));
0723   }
0724 
0725   return lpgbt_ids;
0726 }
0727 
0728 unsigned HGCalTriggerGeometryV9Imp3::getStage1FpgaFromModule(const unsigned module_id) const {
0729   HGCalTriggerModuleDetId id(module_id);
0730 
0731   unsigned stage1_label =
0732       module_to_stage1_.at(packLayerSubdetWaferId(id.layer(), id.triggerSubdetId(), id.moduleU(), id.moduleV()));
0733 
0734   return HGCalTriggerBackendDetId(
0735       id.zside(), HGCalTriggerBackendDetId::BackendType::Stage1FPGA, id.sector(), stage1_label);
0736 }
0737 
0738 GlobalPoint HGCalTriggerGeometryV9Imp3::getTriggerCellPosition(const unsigned trigger_cell_det_id) const {
0739   unsigned det = DetId(trigger_cell_det_id).det();
0740 
0741   // Position: barycenter of the trigger cell.
0742   Basic3DVector<float> triggerCellVector(0., 0., 0.);
0743   const auto cell_ids = getCellsFromTriggerCell(trigger_cell_det_id);
0744   // Scintillator
0745   if (det == DetId::HGCalHSc) {
0746     for (const auto& cell : cell_ids) {
0747       triggerCellVector += hscGeometry()->getPosition(cell).basicVector();
0748     }
0749   }
0750   // HFNose
0751   else if (det == DetId::HGCalTrigger and
0752            HGCalTriggerDetId(trigger_cell_det_id).subdet() == HGCalTriggerSubdetector::HFNoseTrigger) {
0753     for (const auto& cell : cell_ids) {
0754       HFNoseDetId cellDetId(cell);
0755       triggerCellVector += noseGeometry()->getPosition(cellDetId).basicVector();
0756     }
0757   }
0758   // Silicon
0759   else {
0760     for (const auto& cell : cell_ids) {
0761       HGCSiliconDetId cellDetId(cell);
0762       triggerCellVector += (cellDetId.det() == DetId::HGCalEE ? eeGeometry()->getPosition(cellDetId)
0763                                                               : hsiGeometry()->getPosition(cellDetId))
0764                                .basicVector();
0765     }
0766   }
0767   return GlobalPoint(triggerCellVector / cell_ids.size());
0768 }
0769 
0770 GlobalPoint HGCalTriggerGeometryV9Imp3::getModulePosition(const unsigned module_det_id) const {
0771   unsigned subdet = HGCalTriggerModuleDetId(module_det_id).triggerSubdetId();
0772   // Position: barycenter of the module.
0773   Basic3DVector<float> moduleVector(0., 0., 0.);
0774   const auto cell_ids = getCellsFromModule(module_det_id);
0775   // Scintillator
0776   if (subdet == HGCalTriggerSubdetector::HGCalHScTrigger) {
0777     for (const auto& cell : cell_ids) {
0778       moduleVector += hscGeometry()->getPosition(cell).basicVector();
0779     }
0780   }
0781   // HFNose
0782   else if (subdet == HGCalTriggerSubdetector::HFNoseTrigger) {
0783     for (const auto& cell : cell_ids) {
0784       HFNoseDetId cellDetId(cell);
0785       moduleVector += noseGeometry()->getPosition(cellDetId).basicVector();
0786     }
0787   }  // Silicon
0788   else {
0789     for (const auto& cell : cell_ids) {
0790       HGCSiliconDetId cellDetId(cell);
0791       moduleVector += (cellDetId.det() == DetId::HGCalEE ? eeGeometry()->getPosition(cellDetId)
0792                                                          : hsiGeometry()->getPosition(cellDetId))
0793                           .basicVector();
0794     }
0795   }
0796 
0797   return GlobalPoint(moduleVector / cell_ids.size());
0798 }
0799 
0800 void HGCalTriggerGeometryV9Imp3::fillMaps() {
0801   // read json mapping file
0802   json mapping_config;
0803   std::ifstream json_input_file(jsonMappingFile_.fullPath());
0804   if (!json_input_file.is_open()) {
0805     throw cms::Exception("MissingDataFile") << "Cannot open HGCalTriggerGeometry L1TMapping file\n";
0806   }
0807   json_input_file >> mapping_config;
0808 
0809   try {
0810     //Stage 2 to Stage 1 links mapping
0811     for (unsigned stage2_id = 0; stage2_id < mapping_config.at("Stage2").size(); stage2_id++) {
0812       for (unsigned link_id = 0; link_id < mapping_config.at("Stage2").at(stage2_id).at("Stage1Links").size();
0813            link_id++) {
0814         stage2_to_stage1links_.emplace(stage2_id, link_id);
0815         stage1links_samesector_.emplace(
0816             link_id, mapping_config.at("Stage2").at(stage2_id).at("Stage1Links").at(link_id).at("SameSector"));
0817       }
0818     }
0819   } catch (const json::exception& e) {
0820     edm::LogError("HGCalTriggerGeometryV9Imp3")
0821         << "The mapping input json file does not have the expected structure for the Stage2 block";
0822   }
0823 
0824   try {
0825     for (unsigned link_id = 0; link_id < mapping_config.at("Stage1Links").size(); link_id++) {
0826       //Stage 1 links to Stage 1 FPGAs mapping
0827       stage1link_to_stage1_.emplace(link_id, mapping_config.at("Stage1Links").at(link_id).at("Stage1"));
0828 
0829       //Stage 1 links to Stage 2 mapping
0830       stage1link_to_stage2_.emplace(link_id, mapping_config.at("Stage1Links").at(link_id).at("Stage2SameSector"));
0831     }
0832   } catch (const json::exception& e) {
0833     edm::LogError("HGCalTriggerGeometryV9Imp3")
0834         << "The mapping input json file does not have the expected structure for the Stage1Links block";
0835   }
0836 
0837   try {
0838     for (unsigned stage1_id = 0; stage1_id < mapping_config.at("Stage1").size(); stage1_id++) {
0839       //Stage 1 to Stage 1 links mapping
0840       for (auto& link_id : mapping_config.at("Stage1").at(stage1_id).at("Stage1Links")) {
0841         stage1_to_stage1links_.emplace(stage1_id, link_id);
0842       }
0843 
0844       //Stage 1 to lpgbt mapping
0845       std::vector<unsigned> lpgbt_id_vec;
0846       for (auto& lpgbt_id : mapping_config.at("Stage1").at(stage1_id).at("lpgbts")) {
0847         lpgbt_id_vec.push_back(lpgbt_id);
0848       }
0849       stage1_to_lpgbts_.emplace(stage1_id, lpgbt_id_vec);
0850     }
0851 
0852   } catch (const json::exception& e) {
0853     edm::LogError("HGCalTriggerGeometryV9Imp3")
0854         << "The mapping input json file does not have the expected structure for the Stage1 block";
0855   }
0856 
0857   try {
0858     for (unsigned lpgbt_id = 0; lpgbt_id < mapping_config.at("lpgbt").size(); lpgbt_id++) {
0859       //lpgbt to Stage 1 mapping
0860       unsigned stage1_id = mapping_config.at("lpgbt").at(lpgbt_id).at("Stage1");
0861       lpgbt_to_stage1_.emplace(lpgbt_id, stage1_id);
0862 
0863       //lpgbt to module mapping
0864       for (auto& modules : mapping_config.at("lpgbt").at(lpgbt_id).at("Modules")) {
0865         unsigned layer = modules.at("layer");
0866         int subdetId = 0;
0867         bool isSilicon = modules.at("isSilicon");
0868         layerWithoutOffsetAndSubdetId(layer, subdetId, isSilicon);
0869         unsigned packed_value = packLayerSubdetWaferId(layer, subdetId, modules.at("u"), modules.at("v"));
0870         lpgbt_to_modules_.emplace(lpgbt_id, packed_value);
0871 
0872         //fill subsiduary module to stage 1 mapping
0873         auto result = module_to_stage1_.emplace(packed_value, stage1_id);
0874         if (result.second == false &&
0875             stage1_id != result.first->second) {  //check that the stage1_id is the same as in the existing map
0876           edm::LogError("HGCalTriggerGeometryV9Imp3") << "One module is connected to two separate Stage1 FPGAs";
0877         }
0878       }
0879     }
0880 
0881   } catch (const json::exception& e) {
0882     edm::LogError("HGCalTriggerGeometryV9Imp3")
0883         << "The mapping input json file does not have the expected structure for the lpGBT block";
0884   }
0885 
0886   try {
0887     //module to lpgbt mapping
0888     for (unsigned module = 0; module < mapping_config.at("Module").size(); module++) {
0889       unsigned num_elinks = 0;  //Sum number of e-links in each module over lpGBTs
0890       unsigned layer = mapping_config.at("Module").at(module).at("layer");
0891       unsigned moduleU = mapping_config.at("Module").at(module).at("u");
0892       unsigned moduleV = mapping_config.at("Module").at(module).at("v");
0893       bool isSilicon = mapping_config.at("Module").at(module).at("isSilicon");
0894       int subdetId = 0;
0895       layerWithoutOffsetAndSubdetId(layer, subdetId, isSilicon);
0896 
0897       for (auto& lpgbt : mapping_config.at("Module").at(module).at("lpgbts")) {
0898         module_to_lpgbts_.emplace(packLayerSubdetWaferId(layer, subdetId, moduleU, moduleV), lpgbt.at("id"));
0899         num_elinks += unsigned(lpgbt.at("nElinks"));
0900       }
0901       int packed_module = packLayerSubdetWaferId(layer, subdetId, moduleU, moduleV);
0902       links_per_module_.emplace(packed_module, num_elinks);
0903     }
0904   } catch (const json::exception& e) {
0905     edm::LogError("HGCalTriggerGeometryV9Imp3")
0906         << "The mapping input json file does not have the expected structure for the Module block";
0907   }
0908 
0909   json_input_file.close();
0910 }
0911 
0912 unsigned HGCalTriggerGeometryV9Imp3::packLayerSubdetWaferId(unsigned layer, int subdet, int waferU, int waferV) const {
0913   unsigned packed_value = 0;
0914 
0915   packed_value |=
0916       ((waferU & HGCalTriggerModuleDetId::kHGCalModuleUMask) << HGCalTriggerModuleDetId::kHGCalModuleUOffset);
0917   packed_value |=
0918       ((waferV & HGCalTriggerModuleDetId::kHGCalModuleVMask) << HGCalTriggerModuleDetId::kHGCalModuleVOffset);
0919   packed_value |= ((subdet & HGCalTriggerModuleDetId::kHGCalTriggerSubdetMask)
0920                    << HGCalTriggerModuleDetId::kHGCalTriggerSubdetOffset);
0921   packed_value |= ((layer & HGCalTriggerModuleDetId::kHGCalLayerMask) << HGCalTriggerModuleDetId::kHGCalLayerOffset);
0922   return packed_value;
0923 }
0924 
0925 void HGCalTriggerGeometryV9Imp3::unpackLayerSubdetWaferId(
0926     unsigned wafer, unsigned& layer, int& subdet, int& waferU, int& waferV) const {
0927   waferU = (wafer >> HGCalTriggerModuleDetId::kHGCalModuleUOffset) & HGCalTriggerModuleDetId::kHGCalModuleUMask;
0928   waferV = (wafer >> HGCalTriggerModuleDetId::kHGCalModuleVOffset) & HGCalTriggerModuleDetId::kHGCalModuleVMask;
0929   subdet =
0930       (wafer >> HGCalTriggerModuleDetId::kHGCalTriggerSubdetOffset) & HGCalTriggerModuleDetId::kHGCalTriggerSubdetMask;
0931   layer = (wafer >> HGCalTriggerModuleDetId::kHGCalLayerOffset) & HGCalTriggerModuleDetId::kHGCalLayerMask;
0932 }
0933 
0934 void HGCalTriggerGeometryV9Imp3::etaphiMappingFromSector0(int& ieta, int& iphi, unsigned sector) const {
0935   if (sector == 0) {
0936     return;
0937   }
0938   if (sector == 2) {
0939     iphi = iphi + hSc_num_panels_per_sector_;
0940   } else if (sector == 1) {
0941     iphi = iphi + (2 * hSc_num_panels_per_sector_);
0942   }
0943 }
0944 
0945 HGCalGeomRotation::WaferCentring HGCalTriggerGeometryV9Imp3::getWaferCentring(unsigned layer, int subdet) const {
0946   if (subdet == HGCalTriggerSubdetector::HGCalEETrigger) {  // CE-E
0947     return HGCalGeomRotation::WaferCentring::WaferCentred;
0948   } else if (subdet == HGCalTriggerSubdetector::HGCalHSiTrigger) {
0949     if ((layer % 2) == 1) {  // CE-H Odd
0950       return HGCalGeomRotation::WaferCentring::CornerCentredY;
0951     } else {  // CE-H Even
0952       return HGCalGeomRotation::WaferCentring::CornerCentredMercedes;
0953     }
0954   } else if (subdet == HGCalTriggerSubdetector::HFNoseTrigger) {  //HFNose
0955     return HGCalGeomRotation::WaferCentring::WaferCentred;
0956   } else {
0957     edm::LogError("HGCalTriggerGeometryV9Imp3")
0958         << "HGCalTriggerGeometryV9Imp3: trigger sub-detector expected to be silicon";
0959     return HGCalGeomRotation::WaferCentring::WaferCentred;
0960   }
0961 }
0962 
0963 unsigned HGCalTriggerGeometryV9Imp3::tcEtaphiMappingToSector0(int& tc_ieta, int& tc_iphi) const {
0964   unsigned sector = 0;
0965 
0966   if (tc_iphi > hSc_tc_layer0_min_ && tc_iphi <= hSc_tc_layer0_min_ + ntc_per_wafer_) {
0967     sector = 0;
0968   } else if (tc_iphi > hSc_tc_layer0_min_ + ntc_per_wafer_ && tc_iphi <= hSc_tc_layer0_min_ + 2 * ntc_per_wafer_) {
0969     sector = 2;
0970   } else {
0971     sector = 1;
0972   }
0973 
0974   if (sector == 0) {
0975     tc_iphi = tc_iphi - hSc_tc_layer0_min_;
0976   } else if (sector == 2) {
0977     tc_iphi = tc_iphi - (hSc_tc_layer0_min_ + ntc_per_wafer_);
0978   } else if (sector == 1) {
0979     if (tc_iphi <= hSc_tc_layer0_min_) {
0980       tc_iphi = tc_iphi + nSectors_ * ntc_per_wafer_;
0981     }
0982     tc_iphi = tc_iphi - (nSectors_ * ntc_per_wafer_ - hSc_tc_layer0_min_);
0983   }
0984 
0985   return sector;
0986 }
0987 
0988 void HGCalTriggerGeometryV9Imp3::getScintillatoriEtaiPhi(
0989     int& ieta, int& iphi, int tc_eta, int tc_phi, unsigned layer) const {
0990   iphi = (tc_phi - 1) / hSc_tcs_per_module_phi_;  //Phi index 1-12
0991 
0992   int split = hSc_front_layers_split_;
0993   if (layer > hSc_layer_for_split_) {
0994     split = hSc_back_layers_split_;
0995   }
0996   if (tc_eta <= split) {
0997     ieta = 0;
0998   } else {
0999     ieta = 1;
1000   }
1001 }
1002 
1003 bool HGCalTriggerGeometryV9Imp3::validTriggerCell(const unsigned trigger_cell_id) const {
1004   return validTriggerCellFromCells(trigger_cell_id);
1005 }
1006 
1007 bool HGCalTriggerGeometryV9Imp3::disconnectedModule(const unsigned module_id) const {
1008   bool disconnected = false;
1009   HGCalTriggerModuleDetId id(module_id);
1010   if (module_to_stage1_.find(packLayerSubdetWaferId(id.layer(), id.triggerSubdetId(), id.moduleU(), id.moduleV())) ==
1011       module_to_stage1_.end()) {
1012     disconnected = true;
1013   }
1014   if (disconnected_layers_.find(layerWithOffset(module_id)) != disconnected_layers_.end()) {
1015     disconnected = true;
1016   }
1017   return disconnected;
1018 }
1019 
1020 unsigned HGCalTriggerGeometryV9Imp3::triggerLayer(const unsigned id) const {
1021   unsigned layer = layerWithOffset(id);
1022 
1023   if (DetId(id).det() == DetId::HGCalTrigger and
1024       HGCalTriggerDetId(id).subdet() == HGCalTriggerSubdetector::HFNoseTrigger) {
1025     if (layer >= trigger_nose_layers_.size())
1026       return 0;
1027     return trigger_nose_layers_[layer];
1028   }
1029   if (layer >= trigger_layers_.size())
1030     return 0;
1031   return trigger_layers_[layer];
1032 }
1033 
1034 bool HGCalTriggerGeometryV9Imp3::validCell(unsigned cell_id) const {
1035   bool is_valid = false;
1036   unsigned det = DetId(cell_id).det();
1037   switch (det) {
1038     case DetId::HGCalEE:
1039       is_valid = eeTopology().valid(cell_id);
1040       break;
1041     case DetId::HGCalHSi:
1042       is_valid = hsiTopology().valid(cell_id);
1043       break;
1044     case DetId::HGCalHSc:
1045       is_valid = hscTopology().valid(cell_id);
1046       break;
1047     case DetId::Forward:
1048       is_valid = noseTopology().valid(cell_id);
1049       break;
1050     default:
1051       is_valid = false;
1052       break;
1053   }
1054   return is_valid;
1055 }
1056 
1057 bool HGCalTriggerGeometryV9Imp3::validTriggerCellFromCells(const unsigned trigger_cell_id) const {
1058   // Check the validity of a trigger cell with the
1059   // validity of the cells. One valid cell in the
1060   // trigger cell is enough to make the trigger cell
1061   // valid.
1062   const geom_set cells = getCellsFromTriggerCell(trigger_cell_id);
1063   bool is_valid = false;
1064   for (const auto cell_id : cells) {
1065     unsigned det = DetId(cell_id).det();
1066     is_valid |= validCellId(det, cell_id);
1067     if (is_valid)
1068       break;
1069   }
1070   return is_valid;
1071 }
1072 
1073 bool HGCalTriggerGeometryV9Imp3::validCellId(unsigned subdet, unsigned cell_id) const {
1074   bool is_valid = false;
1075   switch (subdet) {
1076     case DetId::HGCalEE:
1077       is_valid = eeTopology().valid(cell_id);
1078       break;
1079     case DetId::HGCalHSi:
1080       is_valid = hsiTopology().valid(cell_id);
1081       break;
1082     case DetId::HGCalHSc:
1083       is_valid = hscTopology().valid(cell_id);
1084       break;
1085     case DetId::Forward:
1086       is_valid = noseTopology().valid(cell_id);
1087       break;
1088     default:
1089       is_valid = false;
1090       break;
1091   }
1092   return is_valid;
1093 }
1094 
1095 int HGCalTriggerGeometryV9Imp3::detIdWaferType(unsigned det, unsigned layer, short waferU, short waferV) const {
1096   int wafer_type = 0;
1097   switch (det) {
1098     case DetId::HGCalEE:
1099       wafer_type = eeTopology().dddConstants().getTypeHex(layer, waferU, waferV);
1100       break;
1101     case DetId::HGCalHSi:
1102       wafer_type = hsiTopology().dddConstants().getTypeHex(layer, waferU, waferV);
1103       break;
1104     case DetId::HGCalHSc:
1105       wafer_type = hscTopology().dddConstants().getTypeTrap(layer);
1106       break;
1107     default:
1108       break;
1109   };
1110   return wafer_type;
1111 }
1112 
1113 void HGCalTriggerGeometryV9Imp3::layerWithoutOffsetAndSubdetId(unsigned& layer, int& subdetId, bool isSilicon) const {
1114   if (!isSilicon) {
1115     layer = layer - heOffset_;
1116     subdetId = HGCalTriggerSubdetector::HGCalHScTrigger;
1117   } else {
1118     if (layer > heOffset_) {
1119       subdetId = HGCalTriggerSubdetector::HGCalHSiTrigger;
1120       layer = layer - heOffset_;
1121     } else {
1122       subdetId = HGCalTriggerSubdetector::HGCalEETrigger;
1123     }
1124   }
1125 }
1126 
1127 unsigned HGCalTriggerGeometryV9Imp3::layerWithOffset(unsigned id) const {
1128   unsigned det = DetId(id).det();
1129   unsigned layer = 0;
1130 
1131   if (det == DetId::HGCalTrigger) {
1132     unsigned subdet = HGCalTriggerDetId(id).subdet();
1133     if (subdet == HGCalTriggerSubdetector::HGCalEETrigger) {
1134       layer = HGCalTriggerDetId(id).layer();
1135     } else if (subdet == HGCalTriggerSubdetector::HGCalHSiTrigger) {
1136       layer = heOffset_ + HGCalTriggerDetId(id).layer();
1137     } else if (subdet == HGCalTriggerSubdetector::HFNoseTrigger) {
1138       layer = HFNoseTriggerDetId(id).layer();
1139     }
1140   } else if (det == DetId::HGCalHSc) {
1141     layer = heOffset_ + HGCScintillatorDetId(id).layer();
1142   } else if (det == DetId::Forward) {
1143     unsigned subdet = HGCalTriggerModuleDetId(id).triggerSubdetId();
1144     if (subdet == HGCalTriggerSubdetector::HGCalEETrigger) {
1145       layer = HGCalTriggerModuleDetId(id).layer();
1146     } else if (subdet == HGCalTriggerSubdetector::HGCalHSiTrigger ||
1147                subdet == HGCalTriggerSubdetector::HGCalHScTrigger) {
1148       layer = heOffset_ + HGCalDetId(id).layer();
1149     } else if (subdet == HGCalTriggerSubdetector::HFNoseTrigger) {
1150       layer = HGCalTriggerModuleDetId(id).layer();
1151     }
1152   }
1153   return layer;
1154 }
1155 
1156 DEFINE_EDM_PLUGIN(HGCalTriggerGeometryFactory, HGCalTriggerGeometryV9Imp3, "HGCalTriggerGeometryV9Imp3");