CMS 3D CMS Logo

HcalTB02HcalNumberingScheme.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalTestBeam
4 // Class : HcalTB02HcalNumberingScheme
5 //
6 // Implementation:
7 // Numbering scheme for hadron calorimeter in 2002 test beam
8 //
9 // Original Author:
10 // Created: Sun 21 10:14:34 CEST 2006
11 //
12 
13 // system include files
14 
15 // user include files
18 
19 #include <CLHEP/Units/SystemOfUnits.h>
20 using CLHEP::degree;
21 using CLHEP::m;
22 
23 //#define EDM_ML_DEBUG
24 //
25 // constructors and destructor
26 //
27 
29  : HcalTB02NumberingScheme(), phiScale(1000000), etaScale(10000) {
30  edm::LogVerbatim("HcalTBSim") << "Creating HcalTB02HcalNumberingScheme";
31 }
32 
34 #ifdef EDM_ML_DEBUG
35  edm::LogVerbatim("HcalTBSim") << "Deleting HcalTB02HcalNumberingScheme";
36 #endif
37 }
38 
39 //
40 // member functions
41 //
42 
43 int HcalTB02HcalNumberingScheme::getUnitID(const G4Step* aStep) const {
44  int scintID = 0;
45 
46  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
47  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
48  float hx = hitPoint.x();
49  float hy = hitPoint.y();
50  float hz = hitPoint.z();
51  float hr = std::sqrt(pow(hx, 2) + pow(hy, 2));
52 
53  // Check if hit happened in first HO layer or second.
54 
55  if ((hr > 3. * m) && (hr < 3.830 * m))
56  return scintID = 17;
57  if (hr > 3.830 * m)
58  return scintID = 18;
59 
60  // Compute the scintID in the HB.
61 
62  float hR = hitPoint.mag(); //sqrt( pow(hx,2)+pow(hy,2)+pow(hz,2) );
63  float htheta = (hR == 0. ? 0. : acos(std::max(std::min(hz / hR, float(1.)), float(-1.))));
64  float hsintheta = sin(htheta);
65  float hphi = (hR * hsintheta == 0. ? 0. : acos(std::max(std::min(hx / (hR * hsintheta), float(1.)), float(-1.))));
66  float heta = (std::fabs(hsintheta) == 1. ? 0. : -std::log(std::fabs(tan(htheta / 2.))));
67  int eta = int(heta / 0.087);
68  int phi = int(hphi / (5. * degree));
69 
70  G4VPhysicalVolume* thePV = preStepPoint->GetPhysicalVolume();
71  int ilayer = ((thePV->GetCopyNo()) / 10) % 100;
72 #ifdef EDM_ML_DEBUG
73  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: Layer " << thePV->GetName()
74  << " found at phi = " << phi << " eta = " << eta << " lay = " << thePV->GetCopyNo()
75  << " " << ilayer;
76 #endif
77  scintID = phiScale * phi + etaScale * eta + ilayer;
78  if (hy < 0.)
79  scintID = -scintID;
80 
81  return scintID;
82 }
83 
85  sID = std::abs(sID);
86  int layerID = sID;
87  if ((layerID != 17) && (layerID != 18))
88  layerID = sID - int(float(sID) / float(etaScale)) * etaScale;
89 
90 #ifdef EDM_ML_DEBUG
91  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " layer = " << layerID;
92 #endif
93  return layerID;
94 }
95 
97  float IDsign = 1.;
98  if (sID < 0)
99  IDsign = -1;
100  sID = std::abs(sID);
101  int phiID = int(float(sID) / float(phiScale));
102 #ifdef EDM_ML_DEBUG
103  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " phi = " << phiID;
104 #endif
105  if (IDsign > 0) {
106  phiID += 4;
107  } else {
108  phiID = std::abs(phiID - 3);
109  }
110  return phiID;
111 }
112 
114  sID = std::abs(sID);
115  int aux = sID - int(float(sID) / float(phiScale)) * phiScale;
116  int etaID = int(float(aux) / float(etaScale));
117 
118 #ifdef EDM_ML_DEBUG
119  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " eta = " << etaID;
120 #endif
121  return etaID;
122 }
Log< level::Info, true > LogVerbatim
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:23
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int getUnitID(const G4Step *aStep) const override
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29