Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-03-30 23:47:24

0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "Geometry/CaloTopology/interface/FastTimeTopology.h"
0003 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
0004 
0005 //#define EDM_ML_DEBUG
0006 
0007 FastTimeTopology::FastTimeTopology(const FastTimeDDDConstants& hdcons, ForwardSubdetector sub, int type)
0008     : hdcons_(hdcons), subdet_(sub), type_(type) {
0009   nEtaZ_ = hdcons_.numberEtaZ(type_);
0010   nPhi_ = hdcons_.numberPhi(type_);
0011   kHGhalf_ = nEtaZ_ * nPhi_;
0012   kHGeomHalf_ = 1;
0013   kSizeForDenseIndexing = (unsigned int)(2 * kHGhalf_);
0014 #ifdef EDM_ML_DEBUG
0015   edm::LogVerbatim("FastTime") << "FastTimeTopology initialized for subDetetcor " << subdet_ << " Type " << type_
0016                                << "  with " << nEtaZ_ << " cells along Z|Eta and " << nPhi_
0017                                << " cells along phi: total channels " << kSizeForDenseIndexing << ":"
0018                                << (2 * kHGeomHalf_);
0019 #endif
0020 }
0021 
0022 uint32_t FastTimeTopology::detId2denseId(const DetId& id) const {
0023   FastTimeTopology::DecodedDetId id_ = decode(id);
0024   uint32_t idx = (uint32_t)(((id_.zside > 0) ? kHGhalf_ : 0) + ((id_.iEtaZ - 1) * nPhi_ + id_.iPhi - 1));
0025   return idx;
0026 }
0027 
0028 DetId FastTimeTopology::denseId2detId(uint32_t hi) const {
0029   if (validHashIndex(hi)) {
0030     FastTimeTopology::DecodedDetId id_;
0031     id_.iType = type_;
0032     id_.zside = ((int)(hi) < kHGhalf_ ? -1 : 1);
0033     int di = ((int)(hi) % kHGhalf_);
0034     int iPhi = (di % nPhi_);
0035     id_.iPhi = iPhi + 1;
0036     id_.iEtaZ = (((di - iPhi) / nPhi_) + 1);
0037     return encode(id_);
0038   } else {
0039     return DetId(0);
0040   }
0041 }
0042 
0043 uint32_t FastTimeTopology::detId2denseGeomId(const DetId& id) const {
0044   FastTimeTopology::DecodedDetId id_ = decode(id);
0045   uint32_t idx = (uint32_t)((id_.zside > 0) ? kHGeomHalf_ : 0);
0046   return idx;
0047 }
0048 
0049 bool FastTimeTopology::valid(const DetId& id) const {
0050   FastTimeTopology::DecodedDetId id_ = decode(id);
0051   bool flag = hdcons_.isValidXY(id_.iType, id_.iEtaZ, id_.iPhi);
0052   return flag;
0053 }
0054 
0055 DetId FastTimeTopology::offsetBy(const DetId startId, int nrStepsX, int nrStepsY) const {
0056   if (startId.det() == DetId::Forward && startId.subdetId() == (int)(subdet_)) {
0057     DetId id = changeXY(startId, nrStepsX, nrStepsY);
0058     if (valid(id))
0059       return id;
0060   }
0061   return DetId(0);
0062 }
0063 
0064 DetId FastTimeTopology::switchZSide(const DetId startId) const {
0065   if (startId.det() == DetId::Forward && startId.subdetId() == (int)(subdet_)) {
0066     FastTimeTopology::DecodedDetId id_ = decode(startId);
0067     id_.zside = -id_.zside;
0068     DetId id = encode(id_);
0069     if (valid(id))
0070       return id;
0071   }
0072   return DetId(0);
0073 }
0074 
0075 FastTimeTopology::DecodedDetId FastTimeTopology::geomDenseId2decId(const uint32_t& hi) const {
0076   FastTimeTopology::DecodedDetId id_;
0077   if (hi < totalGeomModules()) {
0078     id_.zside = ((int)(hi) < kHGeomHalf_ ? -1 : 1);
0079   }
0080   return id_;
0081 }
0082 
0083 FastTimeTopology::DecodedDetId FastTimeTopology::decode(const DetId& startId) const {
0084   FastTimeTopology::DecodedDetId id_;
0085   FastTimeDetId id(startId);
0086   id_.iPhi = id.iphi();
0087   id_.iEtaZ = id.ieta();
0088   id_.iType = id.type();
0089   id_.zside = id.zside();
0090   id_.subdet = id.subdetId();
0091   return id_;
0092 }
0093 
0094 DetId FastTimeTopology::encode(const FastTimeTopology::DecodedDetId& id_) const {
0095   DetId id = FastTimeDetId(id_.iType, id_.iEtaZ, id_.iPhi, id_.zside);
0096   return id;
0097 }
0098 
0099 DetId FastTimeTopology::changeXY(const DetId& id, int nrStepsX, int nrStepsY) const {
0100   FastTimeTopology::DecodedDetId id_ = decode(id);
0101   int iEtaZ = id_.iEtaZ + nrStepsX;
0102   int iPhi = id_.iPhi + nrStepsY;
0103   if (iPhi < 1)
0104     iPhi += nPhi_;
0105   else if (iPhi > nPhi_)
0106     iPhi -= nPhi_;
0107   if (id_.iType == 1 && iEtaZ < 0) {
0108     iEtaZ = -iEtaZ;
0109     id_.zside = -id_.zside;
0110   }
0111   id_.iPhi = iPhi;
0112   id_.iEtaZ = iEtaZ;
0113   DetId nextPoint = encode(id_);
0114   if (valid(nextPoint))
0115     return nextPoint;
0116   else
0117     return DetId(0);
0118 }
0119 
0120 #include "FWCore/Utilities/interface/typelookup.h"
0121 
0122 TYPELOOKUP_DATA_REG(FastTimeTopology);