CMS 3D CMS Logo

FP420NumberingScheme.cc
Go to the documentation of this file.
1 // File: FP420NumberingScheme.cc
3 // Date: 02.2006
4 // Description: Numbering scheme for FP420
5 // Modifications: 08.2008 mside and fside added
9 //
10 #include "CLHEP/Units/GlobalSystemOfUnits.h"
11 #include "globals.hh"
12 #include "G4Step.hh"
13 #include <iostream>
14 
15 //#define EDM_ML_DEBUG
16 
18  // sn0=3, pn0=6, rn0=7;
19 }
20 
22 #ifdef EDM_ML_DEBUG
23  edm::LogVerbatim("FP420") << " Deleting FP420NumberingScheme";
24 #endif
25 }
26 
27 int FP420NumberingScheme::detectorLevel(const G4Step* aStep) const {
28  //Find number of levels
29  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
30  int level = 0;
31  if (touch)
32  level = ((touch->GetHistoryDepth()) + 1);
33  return level;
34 }
35 
36 void FP420NumberingScheme::detectorLevel(const G4Step* aStep, int& level, int* copyno, G4String* name) const {
37  //Get name and copy numbers
38  if (level > 0) {
39  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
40  for (int ii = 0; ii < level; ii++) {
41  int i = level - ii - 1;
42  name[ii] = touch->GetVolume(i)->GetName();
43  copyno[ii] = touch->GetReplicaNumber(i);
44  }
45  }
46 }
47 
48 unsigned int FP420NumberingScheme::getUnitID(const G4Step* aStep) const {
49  unsigned intindex = 0;
50  int level = detectorLevel(aStep);
51 
52 #ifdef EDM_ML_DEBUG
53  edm::LogVerbatim("FP420") << "FP420NumberingScheme number of levels= " << level;
54 #endif
55  // unsigned int intIndex = 0;
56  if (level > 0) {
57  int* copyno = new int[level];
58  G4String* name = new G4String[level];
59  detectorLevel(aStep, level, copyno, name);
60 
61  // int det = static_cast<int>(FP420);;
62 
63  int det = 1;
64  int stationgen = 0;
65  int zside = 0;
66  int station = 0;
67  int plane = 0;
68  for (int ich = 0; ich < level; ich++) {
69  /*
70  // old set up configuration with equidistant stations
71  if(name[ich] == "FP420Ex") {
72  stationgen = copyno[ich];
73  } else if(name[ich] == "SISTATION") {
74  station = stationgen;
75  } else if(name[ich] == "SIPLANE") {
76  plane = copyno[ich];
77  } else if(name[ich] == "SIDETL") {
78  zside = 1;
79  } else if(name[ich] == "SIDETR") {
80  zside = 2;
81  }
82  */
83  // new and old set up configurations are possible:
84  if (name[ich] == "FP420E") {
85  det = copyno[ich];
86  } else if (name[ich] == "HPS240E") {
87  det = copyno[ich] + 2;
88  } else if (name[ich] == "FP420Ex1" || name[ich] == "HPS240Ex1") {
89  stationgen = 1;
90  // } else if(name[ich] == "FP420Ex2") {
91  // stationgen = 2;
92  } else if (name[ich] == "FP420Ex3" || name[ich] == "HPS240Ex3") {
93  stationgen = 2; // was =3
94  } else if (name[ich] == "SISTATION" || name[ich] == "HPS240SISTATION") {
95  station = stationgen;
96  } else if (name[ich] == "SIPLANE" || name[ich] == "HPS240SIPLANE") {
97  plane = copyno[ich];
98  // SIDETL (or R) can be ether X or Y type in next schemes of readout
99  // !!! (=...) zside
100  //
101  // 1 2 <---copyno
102  // Front(=2) Empty(=4) Back(=6) <--SIDETR OR SENSOR2
103  // 1 2 <---copyno
104  // Front(=1) Back(=3) Empty(=5) <--SIDETL OR SENSOR1
105  //
106  } else if (name[ich] == "SENSOR2" || name[ich] == "HPS240SENSOR2") {
107  // } else if(name[ich] == "SIDETR") {
108  // zside = 4 * copyno[ich] - 2 ;//= 2 6 (copyno=1,2)
109  zside = 3 * copyno[ich]; //= 3 6 (copyno=1,2)
110  } else if (name[ich] == "SENSOR1" || name[ich] == "HPS240SENSOR1") {
111  // } else if(name[ich] == "SIDETL") {
112  // zside = 2 * copyno[ich] - 1 ;//= 1 3 (copyno=1,2)
113  zside = copyno[ich]; //= 1 2 (copyno=1,2)
114  }
115  //
116 #ifdef EDM_ML_DEBUG
117  edm::LogVerbatim("FP420") << "FP420NumberingScheme "
118  << "ich=" << ich << "copyno" << copyno[ich] << "name=" << name[ich];
119 #endif
120  }
121  // det = 1 for +FP420 , = 2 for -FP420 / (det-1) = 0,1
122  // det = 3 for +HPS240 , = 4 for -HPS240 / (det-1) = 2,3
123  // 0 is as default for every below:
124  // Z index
125  // station number 1 - 8 (in reality just 2 ones)
126  // superplane(superlayer) number 1 - 16 (in reality just 5 ones)
127 
128  // intindex = myPacker.packEcalIndex (det, zside, station, plane);// see examples
129  // intindex = myPacker.packCastorIndex (det, zside, station, plane);// see examples
130  intindex = packFP420Index(det, zside, station, plane);
131  //
132 #ifdef EDM_ML_DEBUG
133  edm::LogVerbatim("FP420") << "FP420NumberingScheme det=" << det << " zside=" << zside << " station=" << station
134  << " plane=" << plane;
135  for (int ich = 0; ich < level; ich++) {
136  edm::LogVerbatim("FP420") << " name = " << name[ich] << " copy = " << copyno[ich];
137  edm::LogVerbatim("FP420") << " packed index = intindex" << intindex;
138  }
139 #endif
140 
141  delete[] copyno;
142  delete[] name;
143  }
144 
145  return intindex;
146 }
147 
148 unsigned FP420NumberingScheme::packFP420Index(int det, int zside, int station, int plane) {
149  //
150  unsigned int idx = ((det - 1) & 3) << 19; //bit 19-20 (det-1):0- 3 = 4-->2**2 = 4 -> 4-1 ->((det-1)&3) 2 bit: 0-1
151  // unsigned int idx = ((det-1)&1)<<20;//bit 20-20 (det-1):0- 1 = 2-->2**1 = 2 -> 2-1 ->((det-1)&1) 1 bit: 0
152  idx += (zside & 7) << 7; //bits 7- 9 zside:0- 7 = 8-->2**3 = 8 -> 8-1 -> (zside&7) 3 bits:0-2
153  // idx += (zside&3)<<7; //bits 7- 8 zside:0- 3 = 4-->2**2 = 4 -> 4-1 -> (zside&3) 2 bits:0-1
154  idx += (station & 7) << 4; //bits 4- 6 station:0- 7 = 8-->2**3 = 8 -> 8-1 ->(station&7) 3 bits:0-2
155  idx += (plane & 15); //bits 0- 3 plane:0-15 =16-->2**4 =16 ->16-1 -> (plane&15) 4 bits:0-3
156 
157  //
158 
159 #ifdef EDM_ML_DEBUG
160  edm::LogVerbatim("FP420") << "FP420 packing: det " << det << " zside " << zside << " station " << station
161  << " plane " << plane << " idx " << idx;
162 #endif
163  // int newdet, newzside, newstation,newplane;
164  // unpackFP420Index(idx, newdet, newzside, newstation,newplane);
165 
166  //
167 
168  return idx;
169 }
170 
171 void FP420NumberingScheme::unpackFP420Index(const unsigned int& idx, int& det, int& zside, int& station, int& plane) {
172  det = (idx >> 19) & 3;
173  //det = (idx>>20)&1;
174  det += 1;
175  zside = (idx >> 7) & 7;
176  // zside = (idx>>7)&3;
177  station = (idx >> 4) & 7;
178  plane = idx & 15;
179  //
180 
181 #ifdef EDM_ML_DEBUG
182  edm::LogVerbatim("FP420") << " FP420unpacking: idx=" << idx << " zside " << zside << " station " << station
183  << " plane " << plane;
184 #endif
185  //
186 }
Log< level::Info, true > LogVerbatim
static void unpackFP420Index(const unsigned int &idx, int &det, int &zside, int &station, int &superplane)
virtual int detectorLevel(const G4Step *) const
int zside(DetId const &)
static unsigned int packFP420Index(int det, int zside, int station, int superplane)
ii
Definition: cuy.py:589
virtual unsigned int getUnitID(const G4Step *aStep) const