CMS 3D CMS Logo

FP420NumberingScheme.h
Go to the documentation of this file.
1 // File: FP420NumberingScheme.h
3 // Date: 02.2006
4 // Description: Numbering scheme for FP420
5 // Modifications:
7 #ifndef FP420NumberingScheme_h
8 #define FP420NumberingScheme_h
9 
10 #include <map>
11 
12 class G4Step;
13 class G4String;
14 
16 public:
18  virtual ~FP420NumberingScheme();
19 
20  virtual unsigned int getUnitID(const G4Step* aStep) const;
21 
22  // Utilities to get detector levels during a step
23  virtual int detectorLevel(const G4Step*) const;
24  virtual void detectorLevel(const G4Step*, int&, int*, G4String*) const;
25 
26  //protected:
27 
28  static unsigned int packFP420Index(int det, int zside, int station, int superplane);
29 
30  static void unpackFP420Index(const unsigned int& idx, int& det, int& zside, int& station, int& superplane);
31 
32  // private:
33  //
34 
35  static unsigned packMYIndex(int rn0, int pn0, int sn0, int det, int zside, int sector, int zmodule) {
36  int zScale = (rn0 - 1); // rn0=3 - current --> for update rn0=7
37  int sScale = zScale * (pn0 - 1); //pn0=6
38  int dScale = sScale * (sn0 - 1); //sn0=3
39 
40  unsigned int intindex = dScale * (det - 1) + sScale * (sector - 1) + zScale * (zmodule - 1) + zside;
41  //
42  return intindex;
43  }
44 
45  static void unpackMYIndex(const int& idx, int rn0, int pn0, int sn0, int& det, int& zside, int& sector, int& zmodule) {
46  int zScale = (rn0 - 1), sScale = (rn0 - 1) * (pn0 - 1), dScale = (rn0 - 1) * (pn0 - 1) * (sn0 - 1);
47 
48  det = (idx - 1) / dScale + 1;
49  sector = (idx - 1 - dScale * (det - 1)) / sScale + 1;
50  zmodule = (idx - 1 - dScale * (det - 1) - sScale * (sector - 1)) / zScale + 1;
51  zside = idx - dScale * (det - 1) - sScale * (sector - 1) - zScale * (zmodule - 1);
52  }
53 
54  static int unpackLayerIndex(int rn0, int zside) {
55  // 1,2
56  int layerIndex = 1, b;
57  float a = (zside + 1) / 2.;
58  b = (int)a;
59  if (a - b != 0.)
60  layerIndex = 2;
61  //
62  if (zside > (rn0 - 1) || zside < 1)
63  layerIndex = 0;
64  return layerIndex;
65  }
66 
67  static int unpackCopyIndex(int rn0, int zside) {
68  // 1,2,3
69  int copyIndex = 0;
70  if (zside > (rn0 - 1) || zside < 1) {
71  } else {
72  int layerIndex = 1, b;
73  float a = (zside + 1) / 2.;
74  b = (int)a;
75  if (a - b != 0.)
76  layerIndex = 2;
77  if (layerIndex == 2)
78  copyIndex = zside / 2;
79  if (layerIndex == 1)
80  copyIndex = (zside + 1) / 2;
81  }
82 
83  return copyIndex;
84  }
85 
86  static int unpackOrientation(int rn0, int zside) {
87  // Front: Orientation= 1; Back: Orientation= 2
88  int Orientation = 2;
89  if (zside > (rn0 - 1) || zside < 1)
90  Orientation = 0;
91  if (zside == 1 || zside == 2)
92  Orientation = 1;
93  //
94  return Orientation;
95  }
96 
97  static int realzside(int rn0, int zsideinorder) {
98  // zsideinorder:1 2 3 4 5 6
99  //sensorsold 1 0 0 2 0 0
100  //sensorsnew 1 0 5 2 4 0 ???
101  //sensorsnew 1 3 0 2 0 6
102  //zside 1 3 5 2 4 6 over layers 1 and 2
103  int zside, zsidereal;
104  if (zsideinorder < 0) {
105  zside = 0;
106  } else if (zsideinorder < 4) {
107  zside = 2 * zsideinorder - 1;
108  } else if (zsideinorder < 7) {
109  zside = 2 * zsideinorder - 6;
110  } else {
111  zside = 0;
112  }
113  zsidereal = zside;
114  //
115  //old:
116  if (rn0 == 3) {
117  if (zside > 2)
118  zsidereal = 0;
119  }
120  //new:
121  if (rn0 == 7) {
122  // if(zside==3 || zside==6) zsidereal=0; // ????
123  if (zside == 4 || zside == 5)
124  zsidereal = 0;
125  }
126  //
127  return zsidereal;
128  }
129 };
130 
131 #endif
static int unpackCopyIndex(int rn0, int zside)
static int realzside(int rn0, int zsideinorder)
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 int unpackLayerIndex(int rn0, int zside)
static void unpackMYIndex(const int &idx, int rn0, int pn0, int sn0, int &det, int &zside, int &sector, int &zmodule)
static unsigned packMYIndex(int rn0, int pn0, int sn0, int det, int zside, int sector, int zmodule)
static unsigned int packFP420Index(int det, int zside, int station, int superplane)
double b
Definition: hdecay.h:118
virtual unsigned int getUnitID(const G4Step *aStep) const
double a
Definition: hdecay.h:119
static int unpackOrientation(int rn0, int zside)