CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 debug
14 
15 CastorNumberingScheme::CastorNumberingScheme(): lvCAST(0),lvCAES(0),lvCEDS(0),
16  lvCAHS(0),lvCHDS(0),lvCAER(0),
17  lvCEDR(0),lvCAHR(0),lvCHDR(0),
18  lvC3EF(0),lvC3HF(0),lvC4EF(0),
19  lvC4HF(0) {
20  edm::LogInfo("ForwardSim") << "Creating CastorNumberingScheme";
21  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
22  std::vector<lvp>::const_iterator lvcite;
23  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
24  if (strcmp(((*lvcite)->GetName()).c_str(),"CAST") == 0) lvCAST = (*lvcite);
25  if (strcmp(((*lvcite)->GetName()).c_str(),"CAES") == 0) lvCAES = (*lvcite);
26  if (strcmp(((*lvcite)->GetName()).c_str(),"CEDS") == 0) lvCEDS = (*lvcite);
27  if (strcmp(((*lvcite)->GetName()).c_str(),"CAHS") == 0) lvCAHS = (*lvcite);
28  if (strcmp(((*lvcite)->GetName()).c_str(),"CHDS") == 0) lvCHDS = (*lvcite);
29  if (strcmp(((*lvcite)->GetName()).c_str(),"CAER") == 0) lvCAER = (*lvcite);
30  if (strcmp(((*lvcite)->GetName()).c_str(),"CEDR") == 0) lvCEDR = (*lvcite);
31  if (strcmp(((*lvcite)->GetName()).c_str(),"CAHR") == 0) lvCAHR = (*lvcite);
32  if (strcmp(((*lvcite)->GetName()).c_str(),"CHDR") == 0) lvCHDR = (*lvcite);
33  if (strcmp(((*lvcite)->GetName()).c_str(),"C3EF") == 0) lvC3EF = (*lvcite);
34  if (strcmp(((*lvcite)->GetName()).c_str(),"C3HF") == 0) lvC3HF = (*lvcite);
35  if (strcmp(((*lvcite)->GetName()).c_str(),"C4EF") == 0) lvC4EF = (*lvcite);
36  if (strcmp(((*lvcite)->GetName()).c_str(),"C4HF") == 0) lvC4HF = (*lvcite);
37  }
38 #ifdef debug
39  LogDebug("ForwardSim") << "CastorNumberingScheme:: LogicalVolume pointers\n"
40  << lvCAST << " for CAST; " << lvCAES << " for CAES; "
41  << lvCEDS << " for CEDS; " << lvCAHS << " for CAHS; "
42  << lvCHDS << " for CHDS; " << lvCAER << " for CAER; "
43  << lvCEDR << " for CEDR; " << lvCAHR << " for CAHR; "
44  << lvCHDR << " for CHDR; " << lvC3EF << " for C3EF; "
45  << lvC3HF << " for C3HF; " << lvC4EF << " for C4EF; "
46  << lvC4HF << " for C4HF.";
47 #endif
48 }
49 
51  edm::LogInfo("ForwardSim") << "Deleting CastorNumberingScheme";
52 }
53 
54 uint32_t CastorNumberingScheme::getUnitID(const G4Step* aStep) const {
55 
56  uint32_t intindex = 0;
57 
58  uint32_t index=0;
59  int level, copyno[20];
60  lvp lvs[20];
61  detectorLevel(aStep, level, copyno, lvs);
62 
63 #ifdef debug
64  LogDebug("ForwardSim") << "CastorNumberingScheme number of levels= " <<level;
65 #endif
66 
67  if (level > 0) {
68 
69  int zside = 0;
70  int sector = 0;
71  int module = 0;
72 
73  // HcalCastorDetId::Section section;
74  for (int ich=0; ich < level; ich++) {
75  if(lvs[ich] == lvCAST) {
76  // Z index +Z = 1 ; -Z = 2
77  zside = copyno[ich];
78  } else if(lvs[ich] == lvCAES || lvs[ich] == lvCEDS) {
79  // sector number for dead material 1 - 8
80  //section = HcalCastorDetId::EM;
81  if (copyno[ich]<5) {sector = 5-copyno[ich] ;
82  }else{sector = 13-copyno[ich] ;}
83  } else if(lvs[ich] == lvCAHS || lvs[ich] == lvCHDS) {
84  // sector number for dead material 1 - 8
85  if (copyno[ich]<5) {sector = 5-copyno[ich] ;
86  }else{sector = 13-copyno[ich] ;}
87  //section = HcalCastorDetId::HAD;
88  } else if(lvs[ich] == lvCAER || lvs[ich] == lvCEDR) {
89  // zmodule number 1-2 for EM section (2 copies)
90  module = copyno[ich];
91  } else if(lvs[ich] == lvCAHR || lvs[ich] == lvCHDR) {
92  //zmodule number 3-14 for HAD section (12 copies)
93  module = copyno[ich] + 2;
94  } else if(lvs[ich] == lvC3EF || lvs[ich] == lvC3HF) {
95  // sector number for sensitive material 1 - 16
96  sector = sector*2 ;
97  } else if(lvs[ich] == lvC4EF || lvs[ich] == lvC4HF) {
98  // sector number for sensitive material 1 - 16
99  sector = sector*2 -1;
100  }
101 
102 
103 
104 #ifdef debug
105  LogDebug("ForwardSim") << "CastorNumberingScheme " << "ich = " << ich
106  << "copyno" << copyno[ich] << "name = "
107  << lvs[ich]->GetName();
108 #endif
109  }
110  // use for Castor number 9
111  //
112  // Z index +Z = 1 ; -Z = 2
113  // sector number 1 - 16
114  // zmodule number 1 - 18
115 
116 
117  // int det = 9;
118  // intindex = packIndex (det, zside, sector, zmodule);
119 
120  //intindex = packIndex (section, zside, sector, module);
121 
122  intindex = packIndex(zside, sector, module);
123 
124  bool true_for_positive_eta = false;
125  if(zside == 1)true_for_positive_eta = true;
126  if(zside == -1)true_for_positive_eta = false;
127 
128  HcalCastorDetId castorId = HcalCastorDetId(true_for_positive_eta, sector, module);
129  index = castorId.rawId();
130 
131 #ifdef debug
132  LogDebug("ForwardSim") << "CastorNumberingScheme :" <<" zside "
133  << zside << " sector " << sector << " module "
134  << module << " UnitID 0x" << std::hex << intindex
135  << std::dec;
136 #endif
137  }
138  return index;
139 
140 }
141 
142 //uint32_t CastorNumberingScheme::packIndex(int section, int z, int sector, int module ) {
143 
144 uint32_t CastorNumberingScheme::packIndex(int z, int sector, int module) {
145  /*
146  uint32_t idx=(section&31)<<28; //bits 28-31 (21-27 are free for now)
147  idx+=((z-1)&1)<<20; //bits 20 (1...2)
148  idx+=(sector&15)<<6; //bits 6-9 (1...16)
149  idx+=(module&63); //bits 0-5 (1...18)
150  return idx;
151  */
152 
153  uint32_t idx=((z-1)&1)<<8; //bit 8
154  idx+=(sector&15)<<4; //bits 4-7 (1...16)
155  idx+=(module&15); //bits 0-3 (1...14)
156  return idx;
157 
158 }
159 
160 //void CastorNumberingScheme::unpackIndex(const uint32_t& idx, int& section, int& z, int& sector, int& module) {
161 
162 void CastorNumberingScheme::unpackIndex(const uint32_t& idx, int& z, int& sector, int& module) {
163  /*
164  section = (idx>>28)&31;
165  z = (idx>>20)&1;
166  z += 1;
167  sector = (idx>>6)&15;
168  module= (idx&63);
169  */
170  z = (idx>>8)&1;
171  z += 1;
172  sector = (idx>>4)&15;
173  module= (idx&15);
174 }
175 
176 void CastorNumberingScheme::detectorLevel(const G4Step* aStep, int& level,
177  int* copyno, lvp* lvs) const {
178 
179  //Get name and copy numbers
180  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
181  level = 0;
182  if (touch) level = ((touch->GetHistoryDepth())+1);
183  if (level > 0) {
184  for (int ii = 0; ii < level; ii++) {
185  int i = level - ii - 1;
186  lvs[ii] = touch->GetVolume(i)->GetLogicalVolume();
187  copyno[ii] = touch->GetReplicaNumber(i);
188  }
189  }
190 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
static void unpackIndex(const uint32_t &idx, int &z, int &sector, int &zmodule)
double double double z
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
void detectorLevel(const G4Step *, int &, int *, lvp *) const
virtual uint32_t getUnitID(const G4Step *aStep) const
tuple level
Definition: testEve_cfg.py:34
static uint32_t packIndex(int z, int sector, int zmodule)
Definition: vlib.h:209