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  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
45  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
46  float hx = hitPoint.x();
47  float hy = hitPoint.y();
48  float hr = std::sqrt(pow(hx, 2) + pow(hy, 2));
49 
50  // Check if hit happened in first HO layer or second.
51 
52  if ((hr > 3. * m) && (hr < 3.830 * m))
53  return 17;
54  if (hr > 3.830 * m)
55  return 18;
56 
57  // Compute the scintID in the HB.
58  int scintID = 0;
59  float hz = hitPoint.z();
60  float hR = hitPoint.mag(); //sqrt( pow(hx,2)+pow(hy,2)+pow(hz,2) );
61  float htheta = (hR == 0. ? 0. : acos(std::max(std::min(hz / hR, float(1.)), float(-1.))));
62  float hsintheta = sin(htheta);
63  float hphi = (hR * hsintheta == 0. ? 0. : acos(std::max(std::min(hx / (hR * hsintheta), float(1.)), float(-1.))));
64  float heta = (std::fabs(hsintheta) == 1. ? 0. : -std::log(std::fabs(tan(htheta / 2.))));
65  int eta = int(heta / 0.087);
66  int phi = int(hphi / (5. * degree));
67 
68  G4VPhysicalVolume* thePV = preStepPoint->GetPhysicalVolume();
69  int ilayer = ((thePV->GetCopyNo()) / 10) % 100;
70 #ifdef EDM_ML_DEBUG
71  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: Layer " << thePV->GetName()
72  << " found at phi = " << phi << " eta = " << eta << " lay = " << thePV->GetCopyNo()
73  << " " << ilayer;
74 #endif
75  scintID = phiScale * phi + etaScale * eta + ilayer;
76  if (hy < 0.)
77  scintID = -scintID;
78 
79  return scintID;
80 }
81 
83  sID = std::abs(sID);
84  int layerID = sID;
85  if ((layerID != 17) && (layerID != 18))
86  layerID = sID - int(float(sID) / float(etaScale)) * etaScale;
87 
88 #ifdef EDM_ML_DEBUG
89  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " layer = " << layerID;
90 #endif
91  return layerID;
92 }
93 
95  float IDsign = 1.;
96  if (sID < 0)
97  IDsign = -1;
98  sID = std::abs(sID);
99  int phiID = int(float(sID) / float(phiScale));
100 #ifdef EDM_ML_DEBUG
101  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " phi = " << phiID;
102 #endif
103  if (IDsign > 0) {
104  phiID += 4;
105  } else {
106  phiID = std::abs(phiID - 3);
107  }
108  return phiID;
109 }
110 
112  sID = std::abs(sID);
113  int aux = sID - int(float(sID) / float(phiScale)) * phiScale;
114  int etaID = int(float(aux) / float(etaScale));
115 
116 #ifdef EDM_ML_DEBUG
117  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " eta = " << etaID;
118 #endif
119  return etaID;
120 }
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