CMS 3D CMS Logo

CastorNumberingScheme.cc
Go to the documentation of this file.
1 // File: CastorNumberingScheme.cc
3 // Description: Numbering scheme for Castor
8 
9 #include "CLHEP/Units/GlobalSystemOfUnits.h"
10 #include "G4LogicalVolumeStore.hh"
11 #include <iostream>
12 
13 //#define EDM_ML_DEBUG
14 
16  : lvCASTFar(nullptr),
17  lvCASTNear(nullptr),
18  lvCAST(nullptr),
19  lvCAES(nullptr),
20  lvCEDS(nullptr),
21  lvCAHS(nullptr),
22  lvCHDS(nullptr),
23  lvCAER(nullptr),
24  lvCEDR(nullptr),
25  lvCAHR(nullptr),
26  lvCHDR(nullptr),
27  lvC3EF(nullptr),
28  lvC3HF(nullptr),
29  lvC4EF(nullptr),
30  lvC4HF(nullptr) {
31  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
32  std::vector<lvp>::const_iterator lvcite;
33  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
34  if (strcmp(((*lvcite)->GetName()).c_str(), "CASTFar") == 0)
35  lvCASTFar = (*lvcite);
36  if (strcmp(((*lvcite)->GetName()).c_str(), "CASTNear") == 0)
37  lvCASTNear = (*lvcite);
38  if (strcmp(((*lvcite)->GetName()).c_str(), "CAST") == 0)
39  lvCAST = (*lvcite);
40  if (strcmp(((*lvcite)->GetName()).c_str(), "CAES") == 0)
41  lvCAES = (*lvcite);
42  if (strcmp(((*lvcite)->GetName()).c_str(), "CEDS") == 0)
43  lvCEDS = (*lvcite);
44  if (strcmp(((*lvcite)->GetName()).c_str(), "CAHS") == 0)
45  lvCAHS = (*lvcite);
46  if (strcmp(((*lvcite)->GetName()).c_str(), "CHDS") == 0)
47  lvCHDS = (*lvcite);
48  if (strcmp(((*lvcite)->GetName()).c_str(), "CAER") == 0)
49  lvCAER = (*lvcite);
50  if (strcmp(((*lvcite)->GetName()).c_str(), "CEDR") == 0)
51  lvCEDR = (*lvcite);
52  if (strcmp(((*lvcite)->GetName()).c_str(), "CAHR") == 0)
53  lvCAHR = (*lvcite);
54  if (strcmp(((*lvcite)->GetName()).c_str(), "CHDR") == 0)
55  lvCHDR = (*lvcite);
56  if (strcmp(((*lvcite)->GetName()).c_str(), "C3EF") == 0)
57  lvC3EF = (*lvcite);
58  if (strcmp(((*lvcite)->GetName()).c_str(), "C3HF") == 0)
59  lvC3HF = (*lvcite);
60  if (strcmp(((*lvcite)->GetName()).c_str(), "C4EF") == 0)
61  lvC4EF = (*lvcite);
62  if (strcmp(((*lvcite)->GetName()).c_str(), "C4HF") == 0)
63  lvC4HF = (*lvcite);
64  }
65 #ifdef EDM_ML_DEBUG
66  edm::LogVerbatim("ForwardSim") << "Creating CastorNumberingScheme";
67  LogDebug("ForwardSim") << "CastorNumberingScheme:: LogicalVolume pointers\n"
68  << lvCASTFar << " for CASTFar; " << lvCASTNear << " for CASTNear; " << lvCAST << " for CAST; "
69  << lvCAES << " for CAES; " << lvCEDS << " for CEDS; " << lvCAHS << " for CAHS; " << lvCHDS
70  << " for CHDS; " << lvCAER << " for CAER; " << lvCEDR << " for CEDR; " << lvCAHR
71  << " for CAHR; " << lvCHDR << " for CHDR; " << lvC3EF << " for C3EF; " << lvC3HF
72  << " for C3HF; " << lvC4EF << " for C4EF; " << lvC4HF << " for C4HF.";
73 
74  LogDebug("ForwardSim") << "Call to init CastorNumberingScheme\n";
75  for (int mod = 0; mod < 15; mod++)
76  for (int sec = 0; sec < 17; sec++) {
77  HcalCastorDetId castorId = HcalCastorDetId(false, sec, mod);
78  LogDebug("ForwardSim") << "Mod: " << mod << " Sec: " << sec << " Id: " << castorId.rawId() << "\n";
79  }
80 
81 #endif
82 }
83 
85 
86 uint32_t CastorNumberingScheme::getUnitID(const G4Step* aStep) const {
87  uint32_t index = 0;
88  int level, copyno[20];
89  lvp lvs[20];
90  detectorLevel(aStep, level, copyno, lvs);
91 
92 #ifdef EDM_ML_DEBUG
93  LogDebug("ForwardSim") << "CastorNumberingScheme number of levels= " << level;
94 #endif
95 
96  if (level > 0) {
97  int zside = 0;
98  int sector = 0;
99  int module = 0;
100 
101  bool farSide = false;
102  int castorGeoVersion = 0; //0 = original // 1 = separated-halves geometry
103 
104  // HcalCastorDetId::Section section;
105  for (int ich = 0; ich < level; ich++) {
106  if (lvs[ich] == lvCAST) {
107  // Z index +Z = 1 ; -Z = 2
108  assert(1 <= copyno[ich] && copyno[ich] <= 3);
109  zside = copyno[ich] == 1 ? 1 : 2;
110  } // copyno 2 = Far : 3 = Near
111  else if (lvs[ich] == lvCASTFar || lvs[ich] == lvCASTNear) {
112  castorGeoVersion = 1; //detected separated-halves geometry
113  if (lvs[ich] == lvCASTFar)
114  farSide = true;
115  } else if (lvs[ich] == lvCAES || lvs[ich] == lvCEDS || lvs[ich] == lvCAHS || lvs[ich] == lvCHDS) {
116  // sector number for dead material 1 - 8
117  int copyn = copyno[ich];
118  if (castorGeoVersion == 1) {
119  //for separated-half geometry the copy numbers do not start at "3 o'clock" and go from 1-8.
120  //instead they start at "12 o'clock" for near side 1-4. and "6 o'clock" for far side 1-4 again
121  if (farSide) {
122  if (copyn < 3)
123  copyn += 6; //maps 1->7, 2->8
124  else
125  copyn -= 2; //maps 3->1 and 4->2
126  } else { //nearSide
127  copyn += 2; //maps 1->3, ...
128  }
129  } //endif separated-half geometry
130  if (copyn < 5)
131  sector = 5 - copyn;
132  else
133  sector = 13 - copyn;
134  } else if (lvs[ich] == lvCAER || lvs[ich] == lvCEDR) {
135  // zmodule number 1-2 for EM section (2 copies)
136  module = copyno[ich];
137  } else if (lvs[ich] == lvCAHR || lvs[ich] == lvCHDR) {
138  //zmodule number 3-14 for HAD section (12 copies)
139  module = copyno[ich] + 2;
140  } else if (lvs[ich] == lvC3EF || lvs[ich] == lvC3HF) {
141  // sector number for sensitive material 1 - 16
142  sector = sector * 2;
143  } else if (lvs[ich] == lvC4EF || lvs[ich] == lvC4HF) {
144  // sector number for sensitive material 1 - 16
145  sector = sector * 2 - 1;
146  }
147 
148 #ifdef EDM_ML_DEBUG
149  LogDebug("ForwardSim") << "CastorNumberingScheme "
150  << "ich = " << ich << "copyno = " << copyno[ich] << "name = " << lvs[ich]->GetName();
151 #endif
152  } //end for loop over levels
153 
154  // use for Castor number det = 9
155  //
156  // Z index +Z = 1 ; -Z = 2
157  // sector number 1 - 16
158  // zmodule number 1 - 18
159 
160  bool true_for_positive_eta = false;
161  if (zside == 1)
162  true_for_positive_eta = true;
163 
164  HcalCastorDetId castorId = HcalCastorDetId(true_for_positive_eta, sector, module);
165  index = castorId.rawId();
166 
167 #ifdef EDM_ML_DEBUG
168  uint32_t intindex = 0;
169  intindex = packIndex(zside, sector, module);
170  LogDebug("ForwardSim") << "CastorNumberingScheme: "
171  << " zside " << zside << " module " << module << " sector " << sector << " UnitID 0x"
172  << std::hex << intindex << std::dec << " index: " << index;
173 #endif
174  }
175  return index;
176 }
177 
179  /*
180  uint32_t idx=(section&31)<<28; //bits 28-31 (21-27 are free for now)
181  idx+=((z-1)&1)<<20; //bits 20 (1...2)
182  idx+=(sector&15)<<6; //bits 6-9 (1...16)
183  idx+=(module&63); //bits 0-5 (1...18)
184  return idx;
185  */
186 
187  uint32_t idx = ((z - 1) & 1) << 8; //bit 8
188  idx += (sector & 15) << 4; //bits 4-7 (1...16)
189  idx += (module & 15); //bits 0-3 (1...14)
190  return idx;
191 }
192 
193 void CastorNumberingScheme::unpackIndex(const uint32_t& idx, int& z, int& sector, int& module) {
194  /*
195  section = (idx>>28)&31;
196  z = (idx>>20)&1;
197  z += 1;
198  sector = (idx>>6)&15;
199  module= (idx&63);
200  */
201  z = (idx >> 8) & 1;
202  z += 1;
203  sector = (idx >> 4) & 15;
204  module = (idx & 15);
205 }
206 
207 void CastorNumberingScheme::detectorLevel(const G4Step* aStep, int& level, int* copyno, lvp* lvs) const {
208  //Get name and copy numbers
209  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
210  level = 0;
211  if (touch)
212  level = ((touch->GetHistoryDepth()) + 1);
213  if (level > 0) {
214  for (int ii = 0; ii < level; ii++) {
215  int i = level - ii - 1;
216  lvs[ii] = touch->GetVolume(i)->GetLogicalVolume();
217  copyno[ii] = touch->GetReplicaNumber(i);
218  }
219  }
220 }
Log< level::Info, true > LogVerbatim
virtual uint32_t getUnitID(const G4Step *aStep) const
void detectorLevel(const G4Step *, int &, int *, lvp *) const
int zside(DetId const &)
assert(be >=bs)
static void unpackIndex(const uint32_t &idx, int &z, int &sector, int &zmodule)
float float float z
ii
Definition: cuy.py:589
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
static uint32_t packIndex(int z, int sector, int zmodule)
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4
#define LogDebug(id)