CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/SimG4CMS/Calo/src/HFShowerPMT.cc

Go to the documentation of this file.
00001 
00002 // File: HFShowerPMT.cc
00003 // Description: Parametrized version of HF hits
00005 
00006 #include "SimG4CMS/Calo/interface/HFShowerPMT.h"
00007 #include "DetectorDescription/Core/interface/DDFilter.h"
00008 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00009 #include "DetectorDescription/Core/interface/DDValue.h"
00010 
00011 #include "FWCore/Utilities/interface/Exception.h"
00012 
00013 #include "G4NavigationHistory.hh"
00014 #include "G4VPhysicalVolume.hh"
00015 #include "G4Step.hh"
00016 #include "G4Track.hh"
00017 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00018 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00019 
00020 //#define DebugLog
00021 
00022 HFShowerPMT::HFShowerPMT(std::string & name, const DDCompactView & cpv,
00023                          edm::ParameterSet const & p) : cherenkov(0) {
00024 
00025   edm::ParameterSet m_HF  = p.getParameter<edm::ParameterSet>("HFShowerPMT");
00026   pePerGeV                = m_HF.getParameter<double>("PEPerGeVPMT");
00027   
00028   G4String attribute = "ReadOutName";
00029   G4String value     = name;
00030   DDSpecificsFilter filter0;
00031   DDValue           ddv0(attribute,value,0);
00032   filter0.setCriteria(ddv0,DDSpecificsFilter::equals);
00033   DDFilteredView fv0(cpv);
00034   fv0.addFilter(filter0);
00035   if (fv0.firstChild()) {
00036     DDsvalues_type sv0(fv0.mergedSpecifics());
00037 
00038     //Special Geometry parameters
00039     rTable   = getDDDArray("rTable",sv0);
00040     edm::LogInfo("HFShower") << "HFShowerPMT: " << rTable.size() 
00041                              << " rTable (cm)";
00042     for (unsigned int ig=0; ig<rTable.size(); ig++)
00043       edm::LogInfo("HFShower") << "HFShowerPMT: rTable[" << ig << "] = "
00044                                << rTable[ig]/cm << " cm";
00045   } else {
00046     edm::LogError("HFShower") << "HFShowerPMT: cannot get filtered "
00047                               << " view for " << attribute << " matching "
00048                               << value;
00049     throw cms::Exception("Unknown", "HFShowerPMT")
00050       << "cannot match " << attribute << " to " << name <<"\n";
00051   }
00052 
00053   attribute = "Volume";
00054   value     = "HFPMT";
00055   DDSpecificsFilter filter1;
00056   DDValue           ddv1(attribute,value,0);
00057   filter1.setCriteria(ddv1,DDSpecificsFilter::equals);
00058   DDFilteredView fv1(cpv);
00059   fv1.addFilter(filter1);
00060   if (fv1.firstChild()) {
00061     DDsvalues_type sv1(fv1.mergedSpecifics());
00062     std::vector<double> neta;
00063     neta = getDDDArray("indexPMTR",sv1);
00064     for (unsigned int ii=0; ii<neta.size(); ii++) {
00065       int index = static_cast<int>(neta[ii]);
00066       int ir=-1, ifib=-1;
00067       if (index >= 0) {
00068         ir   = index/10; ifib = index%10;
00069       }
00070       pmtR1.push_back(ir);
00071       pmtFib1.push_back(ifib);
00072     }
00073     neta = getDDDArray("indexPMTL",sv1);
00074     for (unsigned int ii=0; ii<neta.size(); ii++) {
00075       int index = static_cast<int>(neta[ii]);
00076       int ir=-1, ifib=-1;
00077       if (index >= 0) {
00078         ir   = index/10; ifib = index%10;
00079       }
00080       pmtR2.push_back(ir);
00081       pmtFib2.push_back(ifib);
00082     }
00083     edm::LogInfo("HFShower") << "HFShowerPMT: gets the Index matches for "
00084                              << neta.size() << " PMTs";
00085     for (unsigned int ii=0; ii<neta.size(); ii++) 
00086       edm::LogInfo("HFShower") << "HFShowerPMT: rIndexR[" << ii << "] = "
00087                                << pmtR1[ii] << " fibreR[" << ii << "] = "
00088                                << pmtFib1[ii] << " rIndexL[" << ii << "] = "
00089                                << pmtR2[ii] << " fibreL[" << ii << "] = "
00090                                << pmtFib2[ii];
00091   } else {
00092     edm::LogWarning("HFShower") << "HFShowerPMT: cannot get filtered "
00093                                 << " view for " << attribute << " matching "
00094                                 << value;
00095   }
00096 
00097   cherenkov = new HFCherenkov(m_HF);
00098 }
00099 
00100 HFShowerPMT::~HFShowerPMT() {
00101   if (cherenkov) delete cherenkov;
00102 }
00103 
00104 double HFShowerPMT::getHits(G4Step * aStep) {
00105 
00106   indexR = indexF = -1;
00107 
00108   G4StepPoint * preStepPoint  = aStep->GetPreStepPoint(); 
00109   const G4VTouchable* touch   = preStepPoint->GetTouchable();
00110   int                 boxNo   = touch->GetReplicaNumber(2);
00111   int                 pmtNo   = touch->GetReplicaNumber(1);
00112   if (boxNo <= 1) {
00113     indexR = pmtR1[pmtNo-1];
00114     indexF = pmtFib1[pmtNo-1];
00115   } else {
00116     indexR = pmtR2[pmtNo-1];
00117     indexF = pmtFib2[pmtNo-1];
00118   }
00119 
00120 #ifdef DebugLog
00121   double edep = aStep->GetTotalEnergyDeposit();
00122   LogDebug("HFShower") << "HFShowerPMT: Box " << boxNo << " PMT "
00123                        << pmtNo << " Mapped Indices " << indexR << ", "
00124                        << indexF << " Edeposit " << edep/MeV << " MeV; PE "
00125                        << edep*pePerGeV/GeV;
00126 #endif
00127 
00128   double photons = 0;
00129   if (indexR >= 0 && indexF > 0) {
00130     G4Track *aTrack = aStep->GetTrack();
00131     G4ParticleDefinition *particleDef = aTrack->GetDefinition();
00132     double stepl = aStep->GetStepLength();
00133     double beta  = preStepPoint->GetBeta();
00134     G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
00135     G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->
00136       GetTopTransform().TransformAxis(pDir);
00137     photons = cherenkov->computeNPEinPMT(particleDef, beta, localMom.x(),
00138                                          localMom.y(), localMom.z(), stepl);
00139 #ifdef DebugLog
00140   LogDebug("HFShower") << "HFShowerPMT::getHits: for particle " 
00141                        << particleDef->GetParticleName() << " Step " << stepl
00142                        << " Beta " << beta << " Direction " << pDir
00143                        << " Local " << localMom << " p.e. " << photons;
00144 #endif 
00145 
00146   }
00147   return photons;
00148 }
00149  
00150 double HFShowerPMT::getRadius() {
00151    
00152   double r = 0.;
00153   if (indexR >= 0 && indexR+1 < (int)(rTable.size()))
00154     r = 0.5*(rTable[indexR]+rTable[indexR+1]);
00155 #ifdef DebugLog
00156   else
00157     LogDebug("HFShower") << "HFShowerPMT::getRadius: R " << indexR
00158                          << " F " << indexF;
00159 #endif
00160   if (indexF == 2)  r =-r;
00161 #ifdef DebugLog
00162   LogDebug("HFShower") << "HFShowerPMT: Radius (" << indexR << "/" << indexF 
00163                        << ") " << r;
00164 #endif
00165   return r;
00166 }
00167 
00168 std::vector<double> HFShowerPMT::getDDDArray(const std::string & str, 
00169                                              const DDsvalues_type & sv) {
00170 
00171 #ifdef DebugLog
00172   LogDebug("HFShower") << "HFShowerPMT:getDDDArray called for " << str;
00173 #endif
00174   DDValue value(str);
00175   if (DDfetch(&sv,value)) {
00176 #ifdef DebugLog
00177     LogDebug("HFShower") << value;
00178 #endif
00179     const std::vector<double> & fvec = value.doubles();
00180     int nval = fvec.size();
00181     if (nval < 2) {
00182       edm::LogError("HFShower") << "HFShowerPMT: # of " << str 
00183                                 << " bins " << nval << " < 2 ==> illegal";
00184       throw cms::Exception("Unknown", "HFShowerPMT")
00185         << "nval < 2 for array " << str << "\n";
00186     }
00187 
00188     return fvec;
00189   } else {
00190     edm::LogError("HFShower") << "HFShowerPMT: cannot get array " << str;
00191     throw cms::Exception("Unknown", "HFShowerPMT") 
00192       << "cannot get array " << str << "\n";
00193   }
00194 }