CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/SimG4CMS/Calo/src/HFShowerFibreBundle.cc

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