CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/SimG4CMS/HcalTestBeam/src/HcalTB02HcalNumberingScheme.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     HcalTestBeam
00004 // Class  :     HcalTB02HcalNumberingScheme
00005 //
00006 // Implementation:
00007 //     Numbering scheme for hadron calorimeter in 2002 test beam
00008 //
00009 // Original Author:
00010 //         Created:  Sun 21 10:14:34 CEST 2006
00011 // $Id: HcalTB02HcalNumberingScheme.cc,v 1.2 2006/11/13 10:32:15 sunanda Exp $
00012 //
00013   
00014 // system include files
00015   
00016 // user include files
00017 #include "SimG4CMS/HcalTestBeam/interface/HcalTB02HcalNumberingScheme.h"
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 
00020 using namespace std;
00021 
00022 //
00023 // constructors and destructor
00024 //
00025 
00026 HcalTB02HcalNumberingScheme::HcalTB02HcalNumberingScheme() : 
00027   HcalTB02NumberingScheme(), phiScale(1000000), etaScale(10000) {
00028   edm::LogInfo("HcalTBSim") << "Creating HcalTB02HcalNumberingScheme";
00029 }
00030 
00031 HcalTB02HcalNumberingScheme::~HcalTB02HcalNumberingScheme() {
00032   edm::LogInfo("HcalTBSim") << "Deleting HcalTB02HcalNumberingScheme";
00033 }
00034 
00035 //
00036 // member functions
00037 //
00038  
00039 int HcalTB02HcalNumberingScheme::getUnitID(const G4Step* aStep) const {
00040 
00041   int scintID = 0;
00042 
00043   G4StepPoint* preStepPoint = aStep->GetPreStepPoint(); 
00044   G4ThreeVector    hitPoint = preStepPoint->GetPosition();
00045   float hx = hitPoint.x();
00046   float hy = hitPoint.y();
00047   float hz = hitPoint.z();
00048   float hr = std::sqrt( pow(hx,2)+pow(hy,2) );
00049 
00050   // Check if hit happened in first HO layer or second.
00051 
00052   if ( (hr > 3.*m) && (hr < 3.830*m) ) return scintID=17;
00053   if (hr > 3.830*m)                    return scintID=18;
00054 
00055   // Compute the scintID in the HB.
00056 
00057   float hR = hitPoint.mag();//sqrt( pow(hx,2)+pow(hy,2)+pow(hz,2) );
00058   float htheta =  (hR == 0. ? 0. : acos(max(min(hz/hR,float(1.)),float(-1.))));
00059   float hsintheta = sin(htheta);
00060   float hphi = (hR*hsintheta == 0. ? 0. :acos( max(min(hx/(hR*hsintheta),float(1.)),float(-1.)) ) );
00061   float heta = ( fabs(hsintheta) == 1.? 0. : -log(fabs(tan(htheta/2.))) );
00062   int eta = int(heta/0.087);
00063   int phi = int(hphi/(5.*degree));
00064 
00065   G4VPhysicalVolume*  thePV = preStepPoint->GetPhysicalVolume();
00066   int ilayer = ((thePV->GetCopyNo())/10)%100;
00067   LogDebug("HcalTBSim") << "HcalTB02HcalNumberingScheme:: Layer " 
00068                         << thePV->GetName() << " found at phi = " << phi
00069                         << " eta = " << eta << " lay = " << thePV->GetCopyNo()
00070                         << " " << ilayer;
00071 
00072   scintID = phiScale*phi + etaScale*eta + ilayer;
00073   if (hy<0.) scintID = -scintID;
00074 
00075   return scintID;
00076 }
00077 
00078 int HcalTB02HcalNumberingScheme::getlayerID(int sID) const {
00079 
00080   sID = abs(sID);
00081   int layerID = sID;
00082   if ( (layerID != 17) && (layerID != 18) )
00083     layerID = sID - int(float(sID)/float(etaScale))*etaScale;
00084 
00085   LogDebug("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID 
00086                         << " layer = " << layerID;
00087   return layerID;
00088 }      
00089 
00090 int HcalTB02HcalNumberingScheme::getphiID(int sID) const {
00091 
00092   float IDsign = 1.;
00093   if (sID<0) IDsign = -1;
00094   sID = abs(sID);
00095   int phiID = int(float(sID)/float(phiScale));
00096   LogDebug("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID 
00097                         << " phi = " << phiID;
00098   if (IDsign>0) {
00099     phiID += 4;
00100   } else {
00101     phiID = abs(phiID-3);
00102   }
00103   return phiID;
00104 }      
00105 
00106 int HcalTB02HcalNumberingScheme::getetaID(int sID) const {
00107 
00108   sID = abs(sID);
00109   int aux = sID - int(float(sID)/float(phiScale))*phiScale;
00110   int etaID = int(float(aux)/float(etaScale));
00111 
00112   LogDebug("HcalTBSim") << "HcalTB02HcalNumberingScheme:: scintID " << sID
00113                         << " eta = " << etaID;
00114   return etaID;
00115 
00116 }