CMS 3D CMS Logo

FastTimeTopology.cc
Go to the documentation of this file.
3 
4 //#define EDM_ML_DEBUG
5 
7  : hdcons_(hdcons), subdet_(sub), type_(type) {
10  kHGhalf_ = nEtaZ_ * nPhi_;
11  kHGeomHalf_ = 1;
12  kSizeForDenseIndexing = (unsigned int)(2 * kHGhalf_);
13 #ifdef EDM_ML_DEBUG
14  std::cout << "FastTimeTopology initialized for subDetetcor " << subdet_ << " Type " << type_ << " with " << nEtaZ_
15  << " cells along Z|Eta and " << nPhi_ << " cells along phi: total channels " << kSizeForDenseIndexing << ":"
16  << (2 * kHGeomHalf_) << std::endl;
17 #endif
18 }
19 
20 uint32_t FastTimeTopology::detId2denseId(const DetId& id) const {
22  uint32_t idx = (uint32_t)(((id_.zside > 0) ? kHGhalf_ : 0) + ((id_.iEtaZ - 1) * nPhi_ + id_.iPhi - 1));
23  return idx;
24 }
25 
27  if (validHashIndex(hi)) {
29  id_.iType = type_;
30  id_.zside = ((int)(hi) < kHGhalf_ ? -1 : 1);
31  int di = ((int)(hi) % kHGhalf_);
32  int iPhi = (di % nPhi_);
33  id_.iPhi = iPhi + 1;
34  id_.iEtaZ = (((di - iPhi) / nPhi_) + 1);
35  return encode(id_);
36  } else {
37  return DetId(0);
38  }
39 }
40 
41 uint32_t FastTimeTopology::detId2denseGeomId(const DetId& id) const {
43  uint32_t idx = (uint32_t)((id_.zside > 0) ? kHGeomHalf_ : 0);
44  return idx;
45 }
46 
47 bool FastTimeTopology::valid(const DetId& id) const {
49  bool flag = hdcons_.isValidXY(id_.iType, id_.iEtaZ, id_.iPhi);
50  return flag;
51 }
52 
53 DetId FastTimeTopology::offsetBy(const DetId startId, int nrStepsX, int nrStepsY) const {
54  if (startId.det() == DetId::Forward && startId.subdetId() == (int)(subdet_)) {
55  DetId id = changeXY(startId, nrStepsX, nrStepsY);
56  if (valid(id))
57  return id;
58  }
59  return DetId(0);
60 }
61 
63  if (startId.det() == DetId::Forward && startId.subdetId() == (int)(subdet_)) {
65  id_.zside = -id_.zside;
66  DetId id = encode(id_);
67  if (valid(id))
68  return id;
69  }
70  return DetId(0);
71 }
72 
75  if (hi < totalGeomModules()) {
76  id_.zside = ((int)(hi) < kHGeomHalf_ ? -1 : 1);
77  }
78  return id_;
79 }
80 
83  FastTimeDetId id(startId);
84  id_.iPhi = id.iphi();
85  id_.iEtaZ = id.ieta();
86  id_.iType = id.type();
87  id_.zside = id.zside();
88  id_.subdet = id.subdetId();
89  return id_;
90 }
91 
93  DetId id = FastTimeDetId(id_.iType, id_.iEtaZ, id_.iPhi, id_.zside);
94  return id;
95 }
96 
97 DetId FastTimeTopology::changeXY(const DetId& id, int nrStepsX, int nrStepsY) const {
99  int iEtaZ = id_.iEtaZ + nrStepsX;
100  int iPhi = id_.iPhi + nrStepsY;
101  if (iPhi < 1)
102  iPhi += nPhi_;
103  else if (iPhi > nPhi_)
104  iPhi -= nPhi_;
105  if (id_.iType == 1 && iEtaZ < 0) {
106  iEtaZ = -iEtaZ;
107  id_.zside = -id_.zside;
108  }
109  id_.iPhi = iPhi;
110  id_.iEtaZ = iEtaZ;
111  DetId nextPoint = encode(id_);
112  if (valid(nextPoint))
113  return nextPoint;
114  else
115  return DetId(0);
116 }
117 
119 
FastTimeTopology.h
FastTimeTopology::validHashIndex
bool validHashIndex(uint32_t ix) const
Definition: FastTimeTopology.h:76
FastTimeTopology::DecodedDetId::iPhi
int iPhi
Definition: FastTimeTopology.h:89
FastTimeTopology::offsetBy
DetId offsetBy(const DetId startId, int nrStepsX, int nrStepsY) const
Definition: FastTimeTopology.cc:53
FastTimeTopology::nPhi_
int nPhi_
Definition: FastTimeTopology.h:105
FastTimeTopology::DecodedDetId::subdet
int subdet
Definition: FastTimeTopology.h:89
ForwardSubdetector
ForwardSubdetector
Definition: ForwardSubdetector.h:4
DetId::det
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
gather_cfg.cout
cout
Definition: gather_cfg.py:144
typelookup.h
FastTimeDDDConstants
Definition: FastTimeDDDConstants.h:19
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
DetId
Definition: DetId.h:17
FastTimeTopology::decode
DecodedDetId decode(const DetId &id) const
Definition: FastTimeTopology.cc:81
FastTimeTopology::kHGhalf_
int kHGhalf_
Definition: FastTimeTopology.h:105
FastTimeTopology::kHGeomHalf_
int kHGeomHalf_
Definition: FastTimeTopology.h:105
FastTimeTopology::FastTimeTopology
FastTimeTopology(const FastTimeDDDConstants &hdcons, ForwardSubdetector subdet, int type)
create a new Topology
Definition: FastTimeTopology.cc:6
FastTimeTopology::DecodedDetId::zside
int zside
Definition: FastTimeTopology.h:89
FastTimeTopology::valid
bool valid(const DetId &id) const override
Is this a valid cell id.
Definition: FastTimeTopology.cc:47
FastTimeTopology::geomDenseId2decId
DecodedDetId geomDenseId2decId(const uint32_t &hi) const
Definition: FastTimeTopology.cc:73
FastTimeTopology::detId2denseGeomId
virtual uint32_t detId2denseGeomId(const DetId &id) const
Definition: FastTimeTopology.cc:41
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
FastTimeTopology::totalGeomModules
unsigned int totalGeomModules() const
Definition: FastTimeTopology.h:79
FastTimeTopology::subdet_
ForwardSubdetector subdet_
Definition: FastTimeTopology.h:104
CaloSubdetectorGeometry.h
FastTimeDDDConstants::numberPhi
int numberPhi(int type) const
Definition: FastTimeDDDConstants.cc:184
type
type
Definition: SiPixelVCal_PayloadInspector.cc:39
FastTimeTopology::switchZSide
DetId switchZSide(const DetId startId) const
Definition: FastTimeTopology.cc:62
createfilelist.int
int
Definition: createfilelist.py:10
FastTimeDDDConstants::numberEtaZ
int numberEtaZ(int type) const
Definition: FastTimeDDDConstants.cc:174
TYPELOOKUP_DATA_REG
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:102
FastTimeTopology::DecodedDetId
Definition: FastTimeTopology.h:87
FastTimeTopology::kSizeForDenseIndexing
unsigned int kSizeForDenseIndexing
Definition: FastTimeTopology.h:106
FastTimeTopology::nEtaZ_
int nEtaZ_
Definition: FastTimeTopology.h:105
FastTimeTopology::DecodedDetId::iEtaZ
int iEtaZ
Definition: FastTimeTopology.h:89
hi
Definition: EPCuts.h:4
FastTimeTopology::encode
DetId encode(const DecodedDetId &id_) const
Definition: FastTimeTopology.cc:92
FastTimeTopology
Definition: FastTimeTopology.h:11
FastTimeTopology::changeXY
DetId changeXY(const DetId &id, int nrStepsX, int nrStepsY) const
move the nagivator along x, y
Definition: FastTimeTopology.cc:97
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:29
FastTimeDetId
Definition: FastTimeDetId.h:8
FastTimeDDDConstants::isValidXY
bool isValidXY(int type, int izeta, int iphi) const
Definition: FastTimeDDDConstants.cc:164
FastTimeTopology::detId2denseId
uint32_t detId2denseId(const DetId &id) const override
Dense indexing.
Definition: FastTimeTopology.cc:20
FastTimeTopology::hdcons_
const FastTimeDDDConstants & hdcons_
Definition: FastTimeTopology.h:103
DetId::Forward
Definition: DetId.h:30
FastTimeTopology::DecodedDetId::iType
int iType
Definition: FastTimeTopology.h:89
FastTimeTopology::type_
int type_
Definition: FastTimeTopology.h:105
FastTimeTopology::denseId2detId
DetId denseId2detId(uint32_t denseId) const override
Definition: FastTimeTopology.cc:26
RemoveAddSevLevel.flag
flag
Definition: RemoveAddSevLevel.py:117