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