CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // $Id: HcalTB02HcalNumberingScheme.cc,v 1.2 2006/11/13 10:32:15 sunanda Exp $
12 //
13 
14 // system include files
15 
16 // user include files
19 
20 using namespace std;
21 
22 //
23 // constructors and destructor
24 //
25 
27  HcalTB02NumberingScheme(), phiScale(1000000), etaScale(10000) {
28  edm::LogInfo("HcalTBSim") << "Creating HcalTB02HcalNumberingScheme";
29 }
30 
32  edm::LogInfo("HcalTBSim") << "Deleting HcalTB02HcalNumberingScheme";
33 }
34 
35 //
36 // member functions
37 //
38 
39 int HcalTB02HcalNumberingScheme::getUnitID(const G4Step* aStep) const {
40 
41  int scintID = 0;
42 
43  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
44  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) ) return scintID=17;
53  if (hr > 3.830*m) 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(max(min(hz/hR,float(1.)),float(-1.))));
59  float hsintheta = sin(htheta);
60  float hphi = (hR*hsintheta == 0. ? 0. :acos( max(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  LogDebug("HcalTBSim") << "HcalTB02HcalNumberingScheme:: Layer "
68  << thePV->GetName() << " found at phi = " << phi
69  << " eta = " << eta << " lay = " << thePV->GetCopyNo()
70  << " " << ilayer;
71 
72  scintID = phiScale*phi + etaScale*eta + ilayer;
73  if (hy<0.) scintID = -scintID;
74 
75  return scintID;
76 }
77 
79 
80  sID = abs(sID);
81  int layerID = sID;
82  if ( (layerID != 17) && (layerID != 18) )
83  layerID = sID - int(float(sID)/float(etaScale))*etaScale;
84 
85  LogDebug("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID
86  << " layer = " << layerID;
87  return layerID;
88 }
89 
91 
92  float IDsign = 1.;
93  if (sID<0) IDsign = -1;
94  sID = abs(sID);
95  int phiID = int(float(sID)/float(phiScale));
96  LogDebug("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID
97  << " phi = " << phiID;
98  if (IDsign>0) {
99  phiID += 4;
100  } else {
101  phiID = abs(phiID-3);
102  }
103  return phiID;
104 }
105 
107 
108  sID = abs(sID);
109  int aux = sID - int(float(sID)/float(phiScale))*phiScale;
110  int etaID = int(float(aux)/float(etaScale));
111 
112  LogDebug("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID
113  << " eta = " << etaID;
114  return etaID;
115 
116 }
#define LogDebug(id)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
virtual int getUnitID(const G4Step *aStep) const
#define abs(x)
Definition: mlp_lapack.h:159
#define min(a, b)
Definition: mlp_lapack.h:161
T eta() const
const T & max(const T &a, const T &b)
T sqrt(T t)
Definition: SSEVec.h:48
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Definition: DDAxes.h:10