Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-17 22:58:45

0001 #include "Geometry/ForwardGeometry/interface/ZdcTopology.h"
0002 #include "DataFormats/HcalDetId/interface/HcalZDCDetId.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include <cmath>
0005 #include <iostream>
0006 #include <algorithm>
0007 
0008 ZdcTopology::ZdcTopology(const HcalDDDRecConstants* hcons)
0009     : hcons_(hcons),
0010       excludeEM_(false),
0011       excludeHAD_(false),
0012       excludeLUM_(false),
0013       excludeRPD_(false),
0014       excludeZP_(false),
0015       excludeZN_(false),
0016       firstEMModule_(1),
0017       lastEMModule_(HcalZDCDetId::kDepEM),
0018       firstHADModule_(1),
0019       lastHADModule_(HcalZDCDetId::kDepHAD),
0020       firstLUMModule_(1),
0021       lastLUMModule_(HcalZDCDetId::kDepLUM),
0022       firstRPDModule_(1),
0023       lastRPDModule_(HcalZDCDetId::kDepRPD) {
0024   mode_ = (HcalTopologyMode::Mode)(hcons_->getTopoMode());
0025   excludeRPD_ = ((mode_ != HcalTopologyMode::Mode::Run3) && (mode_ != HcalTopologyMode::Mode::Run4));
0026   edm::LogVerbatim("ForwardGeom") << "ZdcTopology : Mode " << mode_ << ":" << HcalTopologyMode::Mode::Run3
0027                                   << " ExcludeRPD " << excludeRPD_;
0028 }
0029 
0030 bool ZdcTopology::valid(const HcalZDCDetId& id) const {
0031   // check the raw rules
0032   bool ok = validRaw(id);
0033   ok = ok && !isExcluded(id);
0034   return ok;
0035 }
0036 
0037 bool ZdcTopology::isExcluded(const HcalZDCDetId& id) const {
0038   bool exed = false;
0039 
0040   // check for section exclutions
0041   switch (id.section()) {
0042     case (HcalZDCDetId::EM):
0043       exed = excludeEM_;
0044       break;
0045     case (HcalZDCDetId::HAD):
0046       exed = excludeHAD_;
0047       break;
0048     case (HcalZDCDetId::LUM):
0049       exed = excludeLUM_;
0050       break;
0051     case (HcalZDCDetId::RPD):
0052       exed = excludeRPD_;
0053       break;
0054     default:
0055       exed = true;
0056   }
0057 
0058   // check the entire list
0059   if (!exed && !exclusionList_.empty()) {
0060     std::vector<HcalZDCDetId>::const_iterator i = std::lower_bound(exclusionList_.begin(), exclusionList_.end(), id);
0061     if (i != exclusionList_.end() && *i == id)
0062       exed = true;
0063   }
0064   return exed;
0065 }
0066 
0067 void ZdcTopology::exclude(const HcalZDCDetId& id) {
0068   std::vector<HcalZDCDetId>::iterator i = std::lower_bound(exclusionList_.begin(), exclusionList_.end(), id);
0069   if (i == exclusionList_.end() || *i != id) {
0070     exclusionList_.insert(i, id);
0071   }
0072 }
0073 
0074 void ZdcTopology::exclude(int zside) {
0075   switch (zside) {
0076     case (1):
0077       excludeZP_ = true;
0078       break;
0079     case (-1):
0080       excludeZN_ = true;
0081       break;
0082     default:
0083       break;
0084   }
0085 }
0086 
0087 void ZdcTopology::exclude(int zside, HcalZDCDetId::Section section) {
0088   switch (zside) {
0089     case (1):
0090       excludeZP_ = true;
0091       break;
0092     case (-1):
0093       excludeZN_ = true;
0094       break;
0095     default:
0096       break;
0097   }
0098   switch (section) {
0099     case (HcalZDCDetId::EM):
0100       excludeEM_ = true;
0101       break;
0102     case (HcalZDCDetId::HAD):
0103       excludeHAD_ = true;
0104       break;
0105     case (HcalZDCDetId::LUM):
0106       excludeLUM_ = true;
0107       break;
0108     case (HcalZDCDetId::RPD):
0109       excludeRPD_ = true;
0110       break;
0111     default:
0112       break;
0113   }
0114 }
0115 
0116 int ZdcTopology::exclude(int zside, HcalZDCDetId::Section section, int ich1, int ich2) {
0117   bool exed = false;
0118   switch (zside) {
0119     case (1):
0120       exed = excludeZP_;
0121       break;
0122     case (-1):
0123       exed = excludeZN_;
0124       break;
0125     default:
0126       exed = false;
0127   }
0128   if (exed)
0129     return 0;
0130 
0131   switch (section) {
0132     case (HcalZDCDetId::EM):
0133       exed = excludeEM_;
0134       break;
0135     case (HcalZDCDetId::HAD):
0136       exed = excludeHAD_;
0137       break;
0138     case (HcalZDCDetId::LUM):
0139       exed = excludeLUM_;
0140       break;
0141     case (HcalZDCDetId::RPD):
0142       exed = excludeRPD_;
0143       break;
0144     default:
0145       exed = false;
0146   }
0147   if (exed)
0148     return 0;
0149 
0150   bool isPositive = false;
0151   if (zside == 1)
0152     isPositive = true;
0153 
0154   int n = 0;
0155   for (int ich = ich1; ich < ich2; ich++) {
0156     HcalZDCDetId id(section, isPositive, ich);
0157     if (validRaw(id))
0158       exclude(id);
0159     n++;
0160   }
0161   return n;
0162 }
0163 
0164 bool ZdcTopology::validRaw(const HcalZDCDetId& id) const {
0165   bool ok = true;
0166   if (abs(id.zside()) != 1)
0167     ok = false;
0168   else if (id.channel() <= 0)
0169     ok = false;
0170   else if (!(id.section() == HcalZDCDetId::EM || id.section() == HcalZDCDetId::HAD ||
0171              id.section() == HcalZDCDetId::LUM || id.section() == HcalZDCDetId::RPD))
0172     ok = false;
0173   else if (id.section() == HcalZDCDetId::EM && id.channel() > HcalZDCDetId::kDepEM)
0174     ok = false;
0175   else if (id.section() == HcalZDCDetId::HAD && id.channel() > HcalZDCDetId::kDepHAD)
0176     ok = false;
0177   else if (id.section() == HcalZDCDetId::LUM && id.channel() > HcalZDCDetId::kDepLUM)
0178     ok = false;
0179   else if (id.section() == HcalZDCDetId::RPD && id.channel() > HcalZDCDetId::kDepRPD)
0180     ok = false;
0181   else if (id.section() == HcalZDCDetId::Unknown)
0182     ok = false;
0183   return ok;
0184 }
0185 
0186 std::vector<DetId> ZdcTopology::transverse(const DetId& id) const {
0187   std::vector<DetId> vNeighborsDetId;
0188   HcalZDCDetId zdcId = HcalZDCDetId(id);
0189   HcalZDCDetId zdcDetId;
0190   if (validRaw(zdcId) && zdcId.section() == HcalZDCDetId::EM) {
0191     bool isPositive = false;
0192     if (zdcId.zside() == 1)
0193       isPositive = true;
0194     if (zdcId.channel() == 1) {
0195       zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() + 1);
0196       vNeighborsDetId.emplace_back(zdcDetId.rawId());
0197       return vNeighborsDetId;
0198     }
0199     if (zdcId.channel() == HcalZDCDetId::kDepEM) {
0200       zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() - 1);
0201       vNeighborsDetId.emplace_back(zdcDetId.rawId());
0202       return vNeighborsDetId;
0203     }
0204     zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() - 1);
0205     vNeighborsDetId.emplace_back(zdcDetId.rawId());
0206     zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() + 1);
0207     vNeighborsDetId.emplace_back(zdcDetId.rawId());
0208   }
0209   return vNeighborsDetId;
0210 }
0211 
0212 std::vector<DetId> ZdcTopology::longitudinal(const DetId& id) const {
0213   std::vector<DetId> vNeighborsDetId;
0214   HcalZDCDetId zdcId = HcalZDCDetId(id);
0215   HcalZDCDetId zdcDetId;
0216   if (validRaw(zdcId) && zdcId.section() == HcalZDCDetId::HAD) {
0217     bool isPositive = false;
0218     if (zdcId.zside() == 1)
0219       isPositive = true;
0220     if (zdcId.channel() == 1) {
0221       zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() + 1);
0222       vNeighborsDetId.emplace_back(zdcDetId.rawId());
0223       return vNeighborsDetId;
0224     }
0225     if (zdcId.channel() == HcalZDCDetId::kDepHAD) {
0226       zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() - 1);
0227       vNeighborsDetId.emplace_back(zdcDetId.rawId());
0228       return vNeighborsDetId;
0229     }
0230     zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() - 1);
0231     vNeighborsDetId.emplace_back(zdcDetId.rawId());
0232     zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() + 1);
0233     vNeighborsDetId.emplace_back(zdcDetId.rawId());
0234   }
0235   if (validRaw(zdcId) && zdcId.section() == HcalZDCDetId::LUM) {
0236     bool isPositive = false;
0237     if (zdcId.zside() == 1)
0238       isPositive = true;
0239     if (zdcId.channel() == 1) {
0240       zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() + 1);
0241       vNeighborsDetId.emplace_back(zdcDetId.rawId());
0242       return vNeighborsDetId;
0243     }
0244     if (zdcId.channel() == HcalZDCDetId::kDepLUM) {
0245       zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() - 1);
0246       vNeighborsDetId.emplace_back(zdcDetId.rawId());
0247       return vNeighborsDetId;
0248     }
0249   }
0250   if (validRaw(zdcId) && zdcId.section() == HcalZDCDetId::RPD) {
0251     bool isPositive = false;
0252     if (zdcId.zside() == 1)
0253       isPositive = true;
0254     if (zdcId.channel() == 1) {
0255       zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() + 1);
0256       vNeighborsDetId.emplace_back(zdcDetId.rawId());
0257       return vNeighborsDetId;
0258     }
0259     if (zdcId.channel() == HcalZDCDetId::kDepRPD) {
0260       zdcDetId = HcalZDCDetId(zdcId.section(), isPositive, zdcId.channel() - 1);
0261       vNeighborsDetId.emplace_back(zdcDetId.rawId());
0262       return vNeighborsDetId;
0263     }
0264   }
0265   return vNeighborsDetId;
0266 }
0267 
0268 std::vector<DetId> ZdcTopology::east(const DetId& /*id*/) const {
0269   edm::LogVerbatim("ForwardGeom") << "ZdcTopology::east() not yet implemented";
0270   std::vector<DetId> vNeighborsDetId;
0271   return vNeighborsDetId;
0272 }
0273 
0274 std::vector<DetId> ZdcTopology::west(const DetId& /*id*/) const {
0275   edm::LogVerbatim("ForwardGeom") << "ZdcTopology::west() not yet implemented";
0276   std::vector<DetId> vNeighborsDetId;
0277   return vNeighborsDetId;
0278 }
0279 
0280 std::vector<DetId> ZdcTopology::north(const DetId& /*id*/) const {
0281   edm::LogVerbatim("ForwardGeom") << "ZdcTopology::north() not yet implemented";
0282   std::vector<DetId> vNeighborsDetId;
0283   return vNeighborsDetId;
0284 }
0285 std::vector<DetId> ZdcTopology::south(const DetId& /*id*/) const {
0286   edm::LogVerbatim("ForwardGeom") << "ZdcTopology::south() not yet implemented";
0287   std::vector<DetId> vNeighborsDetId;
0288   return vNeighborsDetId;
0289 }
0290 std::vector<DetId> ZdcTopology::up(const DetId& /*id*/) const {
0291   edm::LogVerbatim("ForwardGeom") << "ZdcTopology::up() not yet implemented";
0292   std::vector<DetId> vNeighborsDetId;
0293   return vNeighborsDetId;
0294 }
0295 std::vector<DetId> ZdcTopology::down(const DetId& /*id*/) const {
0296   edm::LogVerbatim("ForwardGeom") << "ZdcTopology::down() not yet implemented";
0297   std::vector<DetId> vNeighborsDetId;
0298   return vNeighborsDetId;
0299 }
0300 
0301 int ZdcTopology::ncells(HcalZDCDetId::Section section) const {
0302   int ncells = 0;
0303   switch (section) {
0304     case (HcalZDCDetId::EM):
0305       ncells = HcalZDCDetId::kDepEM;
0306       break;
0307     case (HcalZDCDetId::HAD):
0308       ncells = HcalZDCDetId::kDepHAD;
0309       break;
0310     case (HcalZDCDetId::LUM):
0311       ncells = HcalZDCDetId::kDepLUM;
0312       break;
0313     case (HcalZDCDetId::RPD):
0314       ncells = HcalZDCDetId::kDepRPD;
0315       break;
0316     case (HcalZDCDetId::Unknown):
0317       ncells = 0;
0318       break;
0319   }
0320   return ncells;
0321 }
0322 
0323 int ZdcTopology::firstCell(HcalZDCDetId::Section section) const {
0324   int firstCell = 0;
0325   switch (section) {
0326     case (HcalZDCDetId::EM):
0327       firstCell = firstEMModule_;
0328       break;
0329     case (HcalZDCDetId::HAD):
0330       firstCell = firstHADModule_;
0331       break;
0332     case (HcalZDCDetId::LUM):
0333       firstCell = firstLUMModule_;
0334       break;
0335     case (HcalZDCDetId::RPD):
0336       firstCell = firstRPDModule_;
0337       break;
0338     case (HcalZDCDetId::Unknown):
0339       firstCell = 0;
0340       break;
0341   }
0342   return firstCell;
0343 }
0344 
0345 int ZdcTopology::lastCell(HcalZDCDetId::Section section) const {
0346   int lastCell = 0;
0347   switch (section) {
0348     case (HcalZDCDetId::EM):
0349       lastCell = lastEMModule_;
0350       break;
0351     case (HcalZDCDetId::HAD):
0352       lastCell = lastHADModule_;
0353       break;
0354     case (HcalZDCDetId::LUM):
0355       lastCell = lastLUMModule_;
0356       break;
0357     case (HcalZDCDetId::RPD):
0358       lastCell = lastRPDModule_;
0359       break;
0360     case (HcalZDCDetId::Unknown):
0361       lastCell = 0;
0362       break;
0363   }
0364   return lastCell;
0365 }
0366 
0367 uint32_t ZdcTopology::kSizeForDenseIndexing() const {
0368   return (mode_ >= HcalTopologyMode::Mode::Run3 ? HcalZDCDetId::kSizeForDenseIndexingRun3
0369                                                 : HcalZDCDetId::kSizeForDenseIndexingRun1);
0370 }
0371 
0372 DetId ZdcTopology::denseId2detId(uint32_t di) const {
0373   if (validDenseIndex(di)) {
0374     bool lz(false);
0375     uint32_t dp(0);
0376     HcalZDCDetId::Section se(HcalZDCDetId::Unknown);
0377     if (di >= 2 * HcalZDCDetId::kDepRun1) {
0378       lz = (di >= (HcalZDCDetId::kDepRun1 + HcalZDCDetId::kDepTot));
0379       se = HcalZDCDetId::RPD;
0380       dp = 1 + ((di - 2 * HcalZDCDetId::kDepRun1) % HcalZDCDetId::kDepRPD);
0381     } else {
0382       lz = (di >= HcalZDCDetId::kDepRun1);
0383       uint32_t in = (di % HcalZDCDetId::kDepRun1);
0384       se = (in < HcalZDCDetId::kDepEM
0385                 ? HcalZDCDetId::EM
0386                 : (in < HcalZDCDetId::kDepEM + HcalZDCDetId::kDepHAD ? HcalZDCDetId::HAD : HcalZDCDetId::LUM));
0387       dp = (se == HcalZDCDetId::EM ? in + 1
0388                                    : (se == HcalZDCDetId::HAD ? in - HcalZDCDetId::kDepEM + 1
0389                                                               : in - HcalZDCDetId::kDepEM - HcalZDCDetId::kDepHAD + 1));
0390     }
0391     return static_cast<DetId>(HcalZDCDetId(se, lz, dp));
0392   }
0393   return DetId();
0394 }
0395 
0396 uint32_t ZdcTopology::detId2DenseIndex(const DetId& id) const {
0397   HcalZDCDetId detId(id);
0398   const int32_t se(detId.section());
0399   uint32_t di = (detId.channel() - 1 +
0400                  (se == HcalZDCDetId::RPD
0401                       ? 2 * HcalZDCDetId::kDepRun1 + (detId.zside() < 0 ? 0 : HcalZDCDetId::kDepRPD)
0402                       : ((detId.zside() < 0 ? 0 : HcalZDCDetId::kDepRun1) +
0403                          (se == HcalZDCDetId::HAD
0404                               ? HcalZDCDetId::kDepEM
0405                               : (se == HcalZDCDetId::LUM ? HcalZDCDetId::kDepEM + HcalZDCDetId::kDepHAD : 0)))));
0406   return di;
0407 }