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 "G4SystemOfUnits.hh"
20 //#define EDM_ML_DEBUG
21 //
22 // constructors and destructor
23 //
24 
26  : HcalTB02NumberingScheme(), phiScale(1000000), etaScale(10000) {
27  edm::LogVerbatim("HcalTBSim") << "Creating HcalTB02HcalNumberingScheme";
28 }
29 
31 #ifdef EDM_ML_DEBUG
32  edm::LogVerbatim("HcalTBSim") << "Deleting HcalTB02HcalNumberingScheme";
33 #endif
34 }
35 
36 //
37 // member functions
38 //
39 
40 int HcalTB02HcalNumberingScheme::getUnitID(const G4Step* aStep) const {
41  int scintID = 0;
42 
43  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
44  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
45  float hx = hitPoint.x();
46  float hy = hitPoint.y();
47  float hz = hitPoint.z();
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 scintID = 17;
54  if (hr > 3.830 * m)
55  return scintID = 18;
56 
57  // Compute the scintID in the HB.
58 
59  float hR = hitPoint.mag(); //sqrt( pow(hx,2)+pow(hy,2)+pow(hz,2) );
60  float htheta = (hR == 0. ? 0. : acos(std::max(std::min(hz / hR, float(1.)), float(-1.))));
61  float hsintheta = sin(htheta);
62  float hphi = (hR * hsintheta == 0. ? 0. : acos(std::max(std::min(hx / (hR * hsintheta), float(1.)), float(-1.))));
63  float heta = (std::fabs(hsintheta) == 1. ? 0. : -std::log(std::fabs(tan(htheta / 2.))));
64  int eta = int(heta / 0.087);
65  int phi = int(hphi / (5. * degree));
66 
67  G4VPhysicalVolume* thePV = preStepPoint->GetPhysicalVolume();
68  int ilayer = ((thePV->GetCopyNo()) / 10) % 100;
69 #ifdef EDM_ML_DEBUG
70  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: Layer " << thePV->GetName()
71  << " found at phi = " << phi << " eta = " << eta << " lay = " << thePV->GetCopyNo()
72  << " " << ilayer;
73 #endif
74  scintID = phiScale * phi + etaScale * eta + ilayer;
75  if (hy < 0.)
76  scintID = -scintID;
77 
78  return scintID;
79 }
80 
82  sID = std::abs(sID);
83  int layerID = sID;
84  if ((layerID != 17) && (layerID != 18))
85  layerID = sID - int(float(sID) / float(etaScale)) * etaScale;
86 
87 #ifdef EDM_ML_DEBUG
88  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " layer = " << layerID;
89 #endif
90  return layerID;
91 }
92 
94  float IDsign = 1.;
95  if (sID < 0)
96  IDsign = -1;
97  sID = std::abs(sID);
98  int phiID = int(float(sID) / float(phiScale));
99 #ifdef EDM_ML_DEBUG
100  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " phi = " << phiID;
101 #endif
102  if (IDsign > 0) {
103  phiID += 4;
104  } else {
105  phiID = std::abs(phiID - 3);
106  }
107  return phiID;
108 }
109 
111  sID = std::abs(sID);
112  int aux = sID - int(float(sID) / float(phiScale)) * phiScale;
113  int etaID = int(float(aux) / float(etaScale));
114 
115 #ifdef EDM_ML_DEBUG
116  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " eta = " << etaID;
117 #endif
118  return etaID;
119 }
Log< level::Info, true > LogVerbatim
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:19
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