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 
21 //
22 // constructors and destructor
23 //
24 
26  : HcalTB02NumberingScheme(), phiScale(1000000), etaScale(10000) {
27  edm::LogVerbatim("HcalTBSim") << "Creating HcalTB02HcalNumberingScheme";
28 }
29 
31  edm::LogVerbatim("HcalTBSim") << "Deleting HcalTB02HcalNumberingScheme";
32 }
33 
34 //
35 // member functions
36 //
37 
38 int HcalTB02HcalNumberingScheme::getUnitID(const G4Step* aStep) const {
39  int scintID = 0;
40 
41  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
42  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
43  float hx = hitPoint.x();
44  float hy = hitPoint.y();
45  float hz = hitPoint.z();
46  float hr = std::sqrt(pow(hx, 2) + pow(hy, 2));
47 
48  // Check if hit happened in first HO layer or second.
49 
50  if ((hr > 3. * m) && (hr < 3.830 * m))
51  return scintID = 17;
52  if (hr > 3.830 * m)
53  return scintID = 18;
54 
55  // Compute the scintID in the HB.
56 
57  float hR = hitPoint.mag(); //sqrt( pow(hx,2)+pow(hy,2)+pow(hz,2) );
58  float htheta = (hR == 0. ? 0. : acos(std::max(std::min(hz / hR, float(1.)), float(-1.))));
59  float hsintheta = sin(htheta);
60  float hphi = (hR * hsintheta == 0. ? 0. : acos(std::max(std::min(hx / (hR * hsintheta), float(1.)), float(-1.))));
61  float heta = (fabs(hsintheta) == 1. ? 0. : -log(fabs(tan(htheta / 2.))));
62  int eta = int(heta / 0.087);
63  int phi = int(hphi / (5. * degree));
64 
65  G4VPhysicalVolume* thePV = preStepPoint->GetPhysicalVolume();
66  int ilayer = ((thePV->GetCopyNo()) / 10) % 100;
67  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: Layer " << thePV->GetName() << " found at phi = "
68  << phi << " eta = " << eta << " lay = " << thePV->GetCopyNo() << " " << ilayer;
69 
70  scintID = phiScale * phi + etaScale * eta + ilayer;
71  if (hy < 0.)
72  scintID = -scintID;
73 
74  return scintID;
75 }
76 
78  sID = abs(sID);
79  int layerID = sID;
80  if ((layerID != 17) && (layerID != 18))
81  layerID = sID - int(float(sID) / float(etaScale)) * etaScale;
82 
83  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " layer = " << layerID;
84  return layerID;
85 }
86 
88  float IDsign = 1.;
89  if (sID < 0)
90  IDsign = -1;
91  sID = abs(sID);
92  int phiID = int(float(sID) / float(phiScale));
93  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " phi = " << phiID;
94  if (IDsign > 0) {
95  phiID += 4;
96  } else {
97  phiID = abs(phiID - 3);
98  }
99  return phiID;
100 }
101 
103  sID = abs(sID);
104  int aux = sID - int(float(sID) / float(phiScale)) * phiScale;
105  int etaID = int(float(aux) / float(etaScale));
106 
107  edm::LogVerbatim("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID << " eta = " << etaID;
108  return etaID;
109 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:19
susybsm::HSCParticleRef hr
Definition: classes.h:26
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
int getUnitID(const G4Step *aStep) const override
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30