Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:18:57

0001 #include "RecoLocalCalo/HGCalRecAlgos/interface/RecHitTools.h"
0002 
0003 #include "DataFormats/ForwardDetId/interface/HGCalDetId.h"
0004 #include "DataFormats/ForwardDetId/interface/HGCScintillatorDetId.h"
0005 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
0006 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0007 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0008 #include "Geometry/HcalTowerAlgo/interface/HcalGeometry.h"
0009 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0010 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0011 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0012 
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 
0018 using namespace hgcal;
0019 
0020 namespace {
0021   template <typename DDD>
0022   inline void check_ddd(const DDD* ddd) {
0023     if (nullptr == ddd) {
0024       throw cms::Exception("hgcal::RecHitTools") << "DDDConstants not accessibl to hgcal::RecHitTools!";
0025     }
0026   }
0027 
0028   template <typename GEOM>
0029   inline void check_geom(const GEOM* geom) {
0030     if (nullptr == geom) {
0031       throw cms::Exception("hgcal::RecHitTools") << "Geometry not provided yet to hgcal::RecHitTools!";
0032     }
0033   }
0034 
0035   inline const HGCalDDDConstants* get_ddd(const CaloSubdetectorGeometry* geom, const HGCalDetId& detid) {
0036     const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
0037     const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
0038     check_ddd(ddd);
0039     return ddd;
0040   }
0041 
0042   inline const HGCalDDDConstants* get_ddd(const CaloSubdetectorGeometry* geom, const HGCSiliconDetId& detid) {
0043     const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
0044     const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
0045     check_ddd(ddd);
0046     return ddd;
0047   }
0048 
0049   inline const HGCalDDDConstants* get_ddd(const CaloSubdetectorGeometry* geom, const HFNoseDetId& detid) {
0050     const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
0051     const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
0052     check_ddd(ddd);
0053     return ddd;
0054   }
0055 
0056   inline const HGCalDDDConstants* get_ddd(const CaloGeometry* geom, int type, int maxLayerEE, int layer) {
0057     DetId::Detector det = ((type == 0) ? DetId::Forward : ((layer > maxLayerEE) ? DetId::HGCalHSi : DetId::HGCalEE));
0058     int subdet = ((type != 0) ? ForwardSubdetector::ForwardEmpty
0059                               : ((layer > maxLayerEE) ? ForwardSubdetector::HGCEE : ForwardSubdetector::HGCHEF));
0060     const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(det, subdet));
0061     const HGCalDDDConstants* ddd = &(hg->topology().dddConstants());
0062     check_ddd(ddd);
0063     return ddd;
0064   }
0065 
0066 }  // namespace
0067 
0068 void RecHitTools::setGeometry(const CaloGeometry& geom) {
0069   geom_ = &geom;
0070   unsigned int wmaxEE(0), wmaxFH(0);
0071   auto geomEE = static_cast<const HGCalGeometry*>(
0072       geom_->getSubdetectorGeometry(DetId::HGCalEE, ForwardSubdetector::ForwardEmpty));
0073   //check if it's the new geometry
0074   if (geomEE) {
0075     geometryType_ = 1;
0076     eeOffset_ = (geomEE->topology().dddConstants()).getLayerOffset();
0077     wmaxEE = (geomEE->topology().dddConstants()).waferCount(0);
0078     auto geomFH = static_cast<const HGCalGeometry*>(
0079         geom_->getSubdetectorGeometry(DetId::HGCalHSi, ForwardSubdetector::ForwardEmpty));
0080     fhOffset_ = (geomFH->topology().dddConstants()).getLayerOffset();
0081     wmaxFH = (geomFH->topology().dddConstants()).waferCount(0);
0082     fhLastLayer_ = fhOffset_ + (geomFH->topology().dddConstants()).lastLayer(true);
0083     auto geomBH = static_cast<const HGCalGeometry*>(
0084         geom_->getSubdetectorGeometry(DetId::HGCalHSc, ForwardSubdetector::ForwardEmpty));
0085     bhOffset_ = (geomBH->topology().dddConstants()).getLayerOffset();
0086     bhFirstLayer_ = bhOffset_ + (geomBH->topology().dddConstants()).firstLayer();
0087     bhLastLayer_ = bhOffset_ + (geomBH->topology().dddConstants()).lastLayer(true);
0088     bhMaxIphi_ = geomBH->topology().dddConstants().maxCells(true);
0089   } else {
0090     geometryType_ = 0;
0091     geomEE =
0092         static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HGCEE));
0093     eeOffset_ = (geomEE->topology().dddConstants()).getLayerOffset();
0094     wmaxEE = 1 + (geomEE->topology().dddConstants()).waferMax();
0095     auto geomFH =
0096         static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HGCHEF));
0097     fhOffset_ = (geomFH->topology().dddConstants()).getLayerOffset();
0098     fhLastLayer_ = fhOffset_ + (geomFH->topology().dddConstants()).layers(true);
0099     bhOffset_ = fhLastLayer_;
0100     bhFirstLayer_ = bhOffset_ + 1;
0101     wmaxFH = 1 + (geomFH->topology().dddConstants()).waferMax();
0102     auto geomBH =
0103         static_cast<const HcalGeometry*>(geom_->getSubdetectorGeometry(DetId::Hcal, HcalSubdetector::HcalEndcap));
0104     bhLastLayer_ = bhOffset_ + (geomBH->topology().dddConstants())->getMaxDepth(1);
0105   }
0106   maxNumberOfWafersPerLayer_ = std::max(wmaxEE, wmaxFH);
0107   // For nose geometry
0108   auto geomNose =
0109       static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HFNose));
0110   if (geomNose) {
0111     maxNumberOfWafersNose_ = (geomNose->topology().dddConstants()).waferCount(0);
0112     noseLastLayer_ = (geomNose->topology().dddConstants()).layers(true);
0113   } else {
0114     maxNumberOfWafersNose_ = 0;
0115     noseLastLayer_ = 0;
0116   }
0117 }
0118 
0119 const CaloSubdetectorGeometry* RecHitTools::getSubdetectorGeometry(const DetId& id) const {
0120   DetId::Detector det = id.det();
0121   int subdet = (det == DetId::HGCalEE || det == DetId::HGCalHSi || det == DetId::HGCalHSc)
0122                    ? ForwardSubdetector::ForwardEmpty
0123                    : id.subdetId();
0124   auto geom = geom_->getSubdetectorGeometry(det, subdet);
0125   check_geom(geom);
0126   return geom;
0127 }
0128 
0129 GlobalPoint RecHitTools::getPosition(const DetId& id) const {
0130   auto geom = getSubdetectorGeometry(id);
0131   GlobalPoint position;
0132   if (id.det() == DetId::Hcal) {
0133     position = geom->getGeometry(id)->getPosition();
0134   } else {
0135     auto hg = static_cast<const HGCalGeometry*>(geom);
0136     position = hg->getPosition(id);
0137   }
0138   return position;
0139 }
0140 
0141 GlobalPoint RecHitTools::getPositionLayer(int layer, bool nose) const {
0142   unsigned int lay = std::abs(layer);
0143   double z(0);
0144   if (nose) {
0145     auto geomNose =
0146         static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HFNose));
0147     if (geomNose) {
0148       const HGCalDDDConstants* ddd = &(geomNose->topology().dddConstants());
0149       if (ddd)
0150         z = (layer > 0) ? ddd->waferZ(lay, true) : -ddd->waferZ(lay, true);
0151     }
0152   } else {
0153     const HGCalDDDConstants* ddd = get_ddd(geom_, geometryType_, fhOffset_, lay);
0154     if (geometryType_ == 1) {
0155       if (lay > fhOffset_)
0156         lay -= fhOffset_;
0157     }
0158     z = (layer > 0) ? ddd->waferZ(lay, true) : -ddd->waferZ(lay, true);
0159   }
0160   return GlobalPoint(0, 0, z);
0161 }
0162 
0163 int RecHitTools::zside(const DetId& id) const {
0164   int zside = 0;
0165   if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
0166     zside = HGCSiliconDetId(id).zside();
0167   } else if (id.det() == DetId::HGCalHSc) {
0168     zside = HGCScintillatorDetId(id).zside();
0169   } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
0170     zside = HFNoseDetId(id).zside();
0171   } else if (id.det() == DetId::Forward) {
0172     zside = HGCalDetId(id).zside();
0173   } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
0174     zside = HcalDetId(id).zside();
0175   }
0176   return zside;
0177 }
0178 
0179 std::float_t RecHitTools::getSiThickness(const DetId& id) const {
0180   auto geom = getSubdetectorGeometry(id);
0181   std::float_t thick(0.37);
0182   if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
0183     const HGCSiliconDetId hid(id);
0184     auto ddd = get_ddd(geom, hid);
0185     thick = ddd->cellThickness(hid.layer(), hid.waferU(), hid.waferV());
0186   } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
0187     const HFNoseDetId hid(id);
0188     auto ddd = get_ddd(geom, hid);
0189     thick = ddd->cellThickness(hid.layer(), hid.waferU(), hid.waferV());
0190   } else if (id.det() == DetId::Forward) {
0191     const HGCalDetId hid(id);
0192     auto ddd = get_ddd(geom, hid);
0193     thick = ddd->cellThickness(hid.layer(), hid.wafer(), 0);
0194   } else {
0195     LogDebug("getSiThickness::InvalidSiliconDetid")
0196         << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
0197     // It does not make any sense to return 0.37 as thickness for a detector
0198     // that is not Silicon(does it?). Mark it as 0. to be able to distinguish
0199     // the case.
0200     thick = 0.f;
0201   }
0202   return thick;
0203 }
0204 
0205 int RecHitTools::getSiThickIndex(const DetId& id) const {
0206   int thickIndex = -1;
0207   if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
0208     thickIndex = HGCSiliconDetId(id).type();
0209   } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
0210     thickIndex = HFNoseDetId(id).type();
0211   } else if (id.det() == DetId::Forward) {
0212     float thickness = getSiThickness(id);
0213     if (thickness > 99. && thickness < 121.)
0214       thickIndex = 0;
0215     else if (thickness > 199. && thickness < 201.)
0216       thickIndex = 1;
0217     else if (thickness > 299. && thickness < 301.)
0218       thickIndex = 2;
0219     else
0220       assert(thickIndex > 0 && "ERROR - silicon thickness has a nonsensical value");
0221   }
0222   return thickIndex;
0223 }
0224 
0225 std::pair<float, float> RecHitTools::getScintDEtaDPhi(const DetId& id) const {
0226   if (!isScintillator(id)) {
0227     LogDebug("getScintDEtaDPhi::InvalidScintDetid")
0228         << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal scintillator!";
0229     return {0.f, 0.f};
0230   }
0231   auto cellGeom = getSubdetectorGeometry(id)->getGeometry(id);
0232   return {cellGeom->etaSpan(), cellGeom->phiSpan()};
0233 }
0234 
0235 std::float_t RecHitTools::getRadiusToSide(const DetId& id) const {
0236   auto geom = getSubdetectorGeometry(id);
0237   std::float_t size(std::numeric_limits<std::float_t>::max());
0238   if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
0239     const HGCSiliconDetId hid(id);
0240     auto ddd = get_ddd(geom, hid);
0241     size = ddd->cellSizeHex(hid.type());
0242   } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
0243     const HFNoseDetId hid(id);
0244     auto ddd = get_ddd(geom, hid);
0245     size = ddd->cellSizeHex(hid.type());
0246   } else if (id.det() == DetId::Forward) {
0247     const HGCalDetId hid(id);
0248     auto ddd = get_ddd(geom, hid);
0249     size = ddd->cellSizeHex(hid.waferType());
0250   } else {
0251     edm::LogError("getRadiusToSide::InvalidSiliconDetid")
0252         << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
0253   }
0254   return size;
0255 }
0256 
0257 unsigned int RecHitTools::getLayer(const ForwardSubdetector type) const {
0258   int layer(0);
0259   switch (type) {
0260     case (ForwardSubdetector::HGCEE): {
0261       auto geomEE =
0262           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HGCEE));
0263       layer = (geomEE->topology().dddConstants()).layers(true);
0264       break;
0265     }
0266     case (ForwardSubdetector::HGCHEF): {
0267       auto geomFH =
0268           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HGCHEF));
0269       layer = (geomFH->topology().dddConstants()).layers(true);
0270       break;
0271     }
0272     case (ForwardSubdetector::HGCHEB): {
0273       auto geomBH =
0274           static_cast<const HcalGeometry*>(geom_->getSubdetectorGeometry(DetId::Hcal, HcalSubdetector::HcalEndcap));
0275       layer = (geomBH->topology().dddConstants())->getMaxDepth(1);
0276       break;
0277     }
0278     case (ForwardSubdetector::ForwardEmpty): {
0279       auto geomEE =
0280           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HGCEE));
0281       layer = (geomEE->topology().dddConstants()).layers(true);
0282       auto geomFH =
0283           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HGCHEF));
0284       layer += (geomFH->topology().dddConstants()).layers(true);
0285       auto geomBH =
0286           static_cast<const HcalGeometry*>(geom_->getSubdetectorGeometry(DetId::Hcal, HcalSubdetector::HcalEndcap));
0287       if (geomBH)
0288         layer += (geomBH->topology().dddConstants())->getMaxDepth(1);
0289       auto geomBH2 =
0290           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HGCHEB));
0291       if (geomBH2)
0292         layer += (geomBH2->topology().dddConstants()).layers(true);
0293       break;
0294     }
0295     case (ForwardSubdetector::HFNose): {
0296       auto geomHFN =
0297           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HFNose));
0298       layer = (geomHFN->topology().dddConstants()).layers(true);
0299       break;
0300     }
0301     default:
0302       layer = 0;
0303   }
0304   return (unsigned int)(layer);
0305 }
0306 
0307 unsigned int RecHitTools::getLayer(const DetId::Detector type, bool nose) const {
0308   int layer;
0309   switch (type) {
0310     case (DetId::HGCalEE): {
0311       auto geomEE =
0312           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(type, ForwardSubdetector::ForwardEmpty));
0313       layer = (geomEE->topology().dddConstants()).layers(true);
0314       break;
0315     }
0316     case (DetId::HGCalHSi): {
0317       auto geomFH =
0318           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(type, ForwardSubdetector::ForwardEmpty));
0319       layer = (geomFH->topology().dddConstants()).layers(true);
0320       break;
0321     }
0322     case (DetId::HGCalHSc): {
0323       auto geomBH =
0324           static_cast<const HGCalGeometry*>(geom_->getSubdetectorGeometry(type, ForwardSubdetector::ForwardEmpty));
0325       layer = (geomBH->topology().dddConstants()).layers(true);
0326       break;
0327     }
0328     case (DetId::Forward): {
0329       if (nose) {
0330         auto geomHFN = static_cast<const HGCalGeometry*>(
0331             geom_->getSubdetectorGeometry(DetId::Forward, ForwardSubdetector::HFNose));
0332         layer = (geomHFN->topology().dddConstants()).layers(true);
0333       } else {
0334         auto geomEE = static_cast<const HGCalGeometry*>(
0335             geom_->getSubdetectorGeometry(DetId::HGCalEE, ForwardSubdetector::ForwardEmpty));
0336         layer = (geomEE->topology().dddConstants()).layers(true);
0337         auto geomFH = static_cast<const HGCalGeometry*>(
0338             geom_->getSubdetectorGeometry(DetId::HGCalHSi, ForwardSubdetector::ForwardEmpty));
0339         layer += (geomFH->topology().dddConstants()).layers(true);
0340       }
0341       break;
0342     }
0343     default:
0344       layer = 0;
0345   }
0346   return (unsigned int)(layer);
0347 }
0348 
0349 unsigned int RecHitTools::getLayer(const DetId& id) const {
0350   unsigned int layer = std::numeric_limits<unsigned int>::max();
0351   if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
0352     layer = HGCSiliconDetId(id).layer();
0353   } else if (id.det() == DetId::HGCalHSc) {
0354     layer = HGCScintillatorDetId(id).layer();
0355   } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
0356     layer = HFNoseDetId(id).layer();
0357   } else if (id.det() == DetId::Forward) {
0358     layer = HGCalDetId(id).layer();
0359   } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
0360     layer = HcalDetId(id).depth();
0361   }
0362   return layer;
0363 }
0364 
0365 unsigned int RecHitTools::getLayerWithOffset(const DetId& id) const {
0366   unsigned int layer = getLayer(id);
0367   if (id.det() == DetId::Forward && id.subdetId() == HGCHEF) {
0368     layer += fhOffset_;
0369   } else if (id.det() == DetId::HGCalHSi || id.det() == DetId::HGCalHSc) {
0370     // DetId::HGCalHSc hits include the offset w.r.t. EE already
0371     layer += fhOffset_;
0372   } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
0373     layer += bhOffset_;
0374   }
0375   // no need to add offset for HFnose
0376   return layer;
0377 }
0378 
0379 std::pair<int, int> RecHitTools::getWafer(const DetId& id) const {
0380   int waferU = std::numeric_limits<int>::max();
0381   int waferV = 0;
0382   if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
0383     waferU = HGCSiliconDetId(id).waferU();
0384     waferV = HGCSiliconDetId(id).waferV();
0385   } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
0386     waferU = HFNoseDetId(id).waferU();
0387     waferV = HFNoseDetId(id).waferV();
0388   } else if (id.det() == DetId::Forward) {
0389     waferU = HGCalDetId(id).wafer();
0390   } else {
0391     edm::LogError("getWafer::InvalidSiliconDetid")
0392         << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
0393   }
0394   return std::pair<int, int>(waferU, waferV);
0395 }
0396 
0397 std::pair<int, int> RecHitTools::getCell(const DetId& id) const {
0398   int cellU = std::numeric_limits<int>::max();
0399   int cellV = 0;
0400   if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
0401     cellU = HGCSiliconDetId(id).cellU();
0402     cellV = HGCSiliconDetId(id).cellV();
0403   } else if (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)) {
0404     cellU = HFNoseDetId(id).cellU();
0405     cellV = HFNoseDetId(id).cellV();
0406   } else if (id.det() == DetId::Forward) {
0407     cellU = HGCalDetId(id).cell();
0408   } else {
0409     edm::LogError("getCell::InvalidSiliconDetid")
0410         << "det id: " << std::hex << id.rawId() << std::dec << ":" << id.det() << " is not HGCal silicon!";
0411   }
0412   return std::pair<int, int>(cellU, cellV);
0413 }
0414 
0415 bool RecHitTools::isHalfCell(const DetId& id) const {
0416   bool ishalf = false;
0417   if (id.det() == DetId::Forward) {
0418     HGCalDetId hid(id);
0419     auto geom = getSubdetectorGeometry(hid);
0420     auto ddd = get_ddd(geom, hid);
0421     const int waferType = ddd->waferTypeT(hid.waferType());
0422     return ddd->isHalfCell(waferType, hid.cell());
0423   }
0424   //new geometry is always false
0425   return ishalf;
0426 }
0427 
0428 bool RecHitTools::isSilicon(const DetId& id) const {
0429   return (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi ||
0430           (id.det() == DetId::Forward && id.subdetId() == static_cast<int>(HFNose)));
0431 }
0432 
0433 bool RecHitTools::isScintillator(const DetId& id) const { return id.det() == DetId::HGCalHSc; }
0434 
0435 bool RecHitTools::isOnlySilicon(const unsigned int layer) const {
0436   // HFnose TODO
0437   bool isonlysilicon = (layer % bhLastLayer_) < bhOffset_;
0438   return isonlysilicon;
0439 }
0440 
0441 float RecHitTools::getEta(const GlobalPoint& position, const float& vertex_z) const {
0442   GlobalPoint corrected_position = GlobalPoint(position.x(), position.y(), position.z() - vertex_z);
0443   return corrected_position.eta();
0444 }
0445 
0446 float RecHitTools::getEta(const DetId& id, const float& vertex_z) const {
0447   GlobalPoint position = getPosition(id);
0448   float eta = getEta(position, vertex_z);
0449   return eta;
0450 }
0451 
0452 float RecHitTools::getPhi(const GlobalPoint& position) const {
0453   float phi = atan2(position.y(), position.x());
0454   return phi;
0455 }
0456 
0457 float RecHitTools::getPhi(const DetId& id) const {
0458   GlobalPoint position = getPosition(id);
0459   float phi = atan2(position.y(), position.x());
0460   return phi;
0461 }
0462 
0463 float RecHitTools::getPt(const GlobalPoint& position, const float& hitEnergy, const float& vertex_z) const {
0464   float eta = getEta(position, vertex_z);
0465   float pt = hitEnergy / cosh(eta);
0466   return pt;
0467 }
0468 
0469 float RecHitTools::getPt(const DetId& id, const float& hitEnergy, const float& vertex_z) const {
0470   GlobalPoint position = getPosition(id);
0471   float eta = getEta(position, vertex_z);
0472   float pt = hitEnergy / cosh(eta);
0473   return pt;
0474 }
0475 
0476 std::pair<uint32_t, uint32_t> RecHitTools::firstAndLastLayer(DetId::Detector det, int subdet) const {
0477   if ((det == DetId::HGCalEE) || ((det == DetId::Forward) && (subdet == HGCEE))) {
0478     return std::make_pair(eeOffset_ + 1, fhOffset_);
0479   } else if ((det == DetId::HGCalHSi) || ((det == DetId::Forward) && (subdet == HGCHEF))) {
0480     return std::make_pair(fhOffset_ + 1, fhLastLayer_);
0481   } else if ((det == DetId::Forward) && (subdet == HFNose)) {
0482     return std::make_pair(1, noseLastLayer_);
0483   } else {
0484     return std::make_pair(bhFirstLayer_, bhLastLayer_);
0485   }
0486 }
0487 
0488 bool RecHitTools::maskCell(const DetId& id, int corners) const {
0489   if (id.det() == DetId::Hcal) {
0490     return false;
0491   } else {
0492     auto hg = static_cast<const HGCalGeometry*>(getSubdetectorGeometry(id));
0493     return hg->topology().dddConstants().maskCell(id, corners);
0494   }
0495 }