CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/SimG4CMS/Calo/src/HFShowerLibrary.cc

Go to the documentation of this file.
00001 
00002 // File: HFShowerLibrary.cc
00003 // Description: Shower library for Very forward hadron calorimeter
00005 
00006 #include "SimG4CMS/Calo/interface/HFShowerLibrary.h"
00007 #include "SimDataFormats/CaloHit/interface/HFShowerLibraryEventInfo.h"
00008 #include "DetectorDescription/Core/interface/DDFilter.h"
00009 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00010 #include "DetectorDescription/Core/interface/DDValue.h"
00011 
00012 #include "FWCore/Utilities/interface/Exception.h"
00013 
00014 #include "G4VPhysicalVolume.hh"
00015 #include "G4NavigationHistory.hh"
00016 #include "G4Step.hh"
00017 #include "G4Track.hh"
00018 #include "Randomize.hh"
00019 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00020 
00021 //#define DebugLog
00022 
00023 HFShowerLibrary::HFShowerLibrary(std::string & name, const DDCompactView & cpv,
00024                                  edm::ParameterSet const & p) : fibre(0),hf(0),
00025                                                                 emBranch(0),
00026                                                                 hadBranch(0),
00027                                                                 npe(0) {
00028   
00029 
00030   edm::ParameterSet m_HF  = p.getParameter<edm::ParameterSet>("HFShower");
00031   probMax                 = m_HF.getParameter<double>("ProbMax");
00032 
00033   edm::ParameterSet m_HS= p.getParameter<edm::ParameterSet>("HFShowerLibrary");
00034   edm::FileInPath fp       = m_HS.getParameter<edm::FileInPath>("FileName");
00035   std::string pTreeName    = fp.fullPath();
00036   backProb                 = m_HS.getParameter<double>("BackProbability");  
00037   std::string emName       = m_HS.getParameter<std::string>("TreeEMID");
00038   std::string hadName      = m_HS.getParameter<std::string>("TreeHadID");
00039   std::string branchEvInfo = m_HS.getUntrackedParameter<std::string>("BranchEvt","HFShowerLibraryEventInfos_hfshowerlib_HFShowerLibraryEventInfo");
00040   std::string branchPre    = m_HS.getUntrackedParameter<std::string>("BranchPre","HFShowerPhotons_hfshowerlib_");
00041   std::string branchPost   = m_HS.getUntrackedParameter<std::string>("BranchPost","_R.obj");
00042   verbose                  = m_HS.getUntrackedParameter<bool>("Verbosity",false);
00043   applyFidCut              = m_HS.getParameter<bool>("ApplyFiducialCut");
00044 
00045   if (pTreeName.find(".") == 0) pTreeName.erase(0,2);
00046   const char* nTree = pTreeName.c_str();
00047   hf                = TFile::Open(nTree);
00048 
00049   if (!hf->IsOpen()) { 
00050     edm::LogError("HFShower") << "HFShowerLibrary: opening " << nTree 
00051                               << " failed";
00052     throw cms::Exception("Unknown", "HFShowerLibrary") 
00053       << "Opening of " << pTreeName << " fails\n";
00054   } else {
00055     edm::LogInfo("HFShower") << "HFShowerLibrary: opening " << nTree 
00056                              << " successfully"; 
00057   }
00058 
00059   TTree* event = (TTree *) hf ->Get("Events");
00060   if (event) {
00061     std::string info = branchEvInfo + branchPost;
00062     TBranch *evtInfo = event->GetBranch(info.c_str());
00063     if (evtInfo) {
00064       loadEventInfo(evtInfo);
00065     } else {
00066       edm::LogError("HFShower") << "HFShowerLibrary: HFShowerLibrayEventInfo"
00067                                 << " Branch does not exist in Event";
00068       throw cms::Exception("Unknown", "HFShowerLibrary")
00069         << "Event information absent\n";
00070     }
00071   } else {
00072     edm::LogError("HFShower") << "HFShowerLibrary: Events Tree does not "
00073                               << "exist";
00074     throw cms::Exception("Unknown", "HFShowerLibrary")
00075       << "Events tree absent\n";
00076   }
00077   
00078   edm::LogInfo("HFShower") << "HFShowerLibrary: Library " << libVers 
00079                            << " ListVersion "   << listVersion 
00080                            << " Events Total " << totEvents << " and "
00081                            << evtPerBin << " per bin";
00082   edm::LogInfo("HFShower") << "HFShowerLibrary: Energies (GeV) with " 
00083                            << nMomBin   << " bins";
00084   for (int i=0; i<nMomBin; i++)
00085     edm::LogInfo("HFShower") << "HFShowerLibrary: pmom[" << i << "] = "
00086                              << pmom[i]/GeV << " GeV";
00087 
00088   std::string nameBr = branchPre + emName + branchPost;
00089   emBranch         = event->GetBranch(nameBr.c_str());
00090   if (verbose) emBranch->Print();
00091   nameBr           = branchPre + hadName + branchPost;
00092   hadBranch        = event->GetBranch(nameBr.c_str());
00093   if (verbose) hadBranch->Print();
00094   edm::LogInfo("HFShower") << "HFShowerLibrary:Branch " << emName 
00095                            << " has " << emBranch->GetEntries() 
00096                            << " entries and Branch " << hadName 
00097                            << " has " << hadBranch->GetEntries() 
00098                            << " entries";
00099   edm::LogInfo("HFShower") << "HFShowerLibrary::No packing information -"
00100                            << " Assume x, y, z are not in packed form";
00101 
00102   edm::LogInfo("HFShower") << "HFShowerLibrary: Maximum probability cut off " 
00103                            << probMax << "  Back propagation of light prob. "
00104                            << backProb ;
00105   
00106   G4String attribute = "ReadOutName";
00107   G4String value     = name;
00108   DDSpecificsFilter filter;
00109   DDValue           ddv(attribute,value,0);
00110   filter.setCriteria(ddv,DDSpecificsFilter::equals);
00111   DDFilteredView fv(cpv);
00112   fv.addFilter(filter);
00113   bool dodet = fv.firstChild();
00114   if (dodet) {
00115     DDsvalues_type sv(fv.mergedSpecifics());
00116 
00117     //Radius (minimum and maximum)
00118     int nR     = -1;
00119     std::vector<double> rTable = getDDDArray("rTable",sv,nR);
00120     rMin = rTable[0];
00121     rMax = rTable[nR-1];
00122     edm::LogInfo("HFShower") << "HFShowerLibrary: rMIN " << rMin/cm 
00123                              << " cm and rMax " << rMax/cm;
00124 
00125     //Delta phi
00126     int nEta   = -1;
00127     std::vector<double> etaTable = getDDDArray("etaTable",sv,nEta);
00128     int nPhi   = nEta + nR - 2;
00129     std::vector<double> phibin   = getDDDArray("phibin",sv,nPhi);
00130     dphi       = phibin[nEta-1];
00131     edm::LogInfo("HFShower") << "HFShowerLibrary: (Half) Phi Width of wedge " 
00132                              << dphi/deg;
00133 
00134     //Special Geometry parameters
00135     int ngpar = 7;
00136     gpar      = getDDDArray("gparHF",sv,ngpar);
00137     edm::LogInfo("HFShower") << "HFShowerLibrary: " << ngpar << " gpar (cm)";
00138     for (int ig=0; ig<ngpar; ig++)
00139       edm::LogInfo("HFShower") << "HFShowerLibrary: gpar[" << ig << "] = "
00140                                << gpar[ig]/cm << " cm";
00141   } else {
00142     edm::LogError("HFShower") << "HFShowerLibrary: cannot get filtered "
00143                               << " view for " << attribute << " matching "
00144                               << name;
00145     throw cms::Exception("Unknown", "HFShowerLibrary")
00146       << "cannot match " << attribute << " to " << name <<"\n";
00147   }
00148   
00149   fibre = new HFFibre(name, cpv, p);
00150   emPDG = epPDG = gammaPDG = 0;
00151   pi0PDG = etaPDG = nuePDG = numuPDG = nutauPDG= 0;
00152   anuePDG= anumuPDG = anutauPDG = geantinoPDG = 0;
00153 }
00154 
00155 HFShowerLibrary::~HFShowerLibrary() {
00156   if (hf)     hf->Close();
00157   if (fibre)  delete   fibre;  fibre  = 0;
00158 }
00159 
00160 void HFShowerLibrary::initRun(G4ParticleTable * theParticleTable) {
00161 
00162   G4String parName;
00163   emPDG = theParticleTable->FindParticle(parName="e-")->GetPDGEncoding();
00164   epPDG = theParticleTable->FindParticle(parName="e+")->GetPDGEncoding();
00165   gammaPDG = theParticleTable->FindParticle(parName="gamma")->GetPDGEncoding();
00166   pi0PDG = theParticleTable->FindParticle(parName="pi0")->GetPDGEncoding();
00167   etaPDG = theParticleTable->FindParticle(parName="eta")->GetPDGEncoding();
00168   nuePDG = theParticleTable->FindParticle(parName="nu_e")->GetPDGEncoding();
00169   numuPDG = theParticleTable->FindParticle(parName="nu_mu")->GetPDGEncoding();
00170   nutauPDG= theParticleTable->FindParticle(parName="nu_tau")->GetPDGEncoding();
00171   anuePDG = theParticleTable->FindParticle(parName="anti_nu_e")->GetPDGEncoding();
00172   anumuPDG= theParticleTable->FindParticle(parName="anti_nu_mu")->GetPDGEncoding();
00173   anutauPDG= theParticleTable->FindParticle(parName="anti_nu_tau")->GetPDGEncoding();
00174   geantinoPDG= theParticleTable->FindParticle(parName="geantino")->GetPDGEncoding();
00175 #ifdef DebugLog
00176   edm::LogInfo("HFShower") << "HFShowerLibrary: Particle codes for e- = " 
00177                            << emPDG << ", e+ = " << epPDG << ", gamma = " 
00178                            << gammaPDG << ", pi0 = " << pi0PDG << ", eta = " 
00179                            << etaPDG << ", geantino = " << geantinoPDG 
00180                            << "\n        nu_e = " << nuePDG << ", nu_mu = " 
00181                            << numuPDG << ", nu_tau = " << nutauPDG 
00182                            << ", anti_nu_e = " << anuePDG << ", anti_nu_mu = " 
00183                            << anumuPDG << ", anti_nu_tau = " << anutauPDG;
00184 #endif
00185 }
00186 
00187 
00188 std::vector<HFShowerLibrary::Hit> HFShowerLibrary::getHits(G4Step * aStep,
00189                                                            bool & ok,
00190                                                            bool onlyLong) {
00191 
00192   G4StepPoint * preStepPoint  = aStep->GetPreStepPoint(); 
00193   G4StepPoint * postStepPoint = aStep->GetPostStepPoint(); 
00194   G4Track *     track    = aStep->GetTrack();
00195   // Get Z-direction 
00196   const G4DynamicParticle *aParticle = track->GetDynamicParticle();
00197   G4ThreeVector momDir = aParticle->GetMomentumDirection();
00198   //  double mom = aParticle->GetTotalMomentum();
00199 
00200   G4ThreeVector hitPoint = preStepPoint->GetPosition();   
00201   G4String      partType = track->GetDefinition()->GetParticleName();
00202   int           parCode  = track->GetDefinition()->GetPDGEncoding();
00203 
00204   std::vector<HFShowerLibrary::Hit> hit;
00205   ok = false;
00206   if (parCode == pi0PDG || parCode == etaPDG || parCode == nuePDG ||
00207       parCode == numuPDG || parCode == nutauPDG || parCode == anuePDG ||
00208       parCode == anumuPDG || parCode == anutauPDG || parCode == geantinoPDG) 
00209     return hit;
00210   ok = true;
00211 
00212   double tSlice = (postStepPoint->GetGlobalTime())/nanosecond;
00213   double pin    = preStepPoint->GetTotalEnergy();
00214   double pz     = momDir.z(); 
00215   double zint   = hitPoint.z(); 
00216 
00217   // if particle moves from interaction point or "backwards (halo)
00218   int backward = 0;
00219   if (pz * zint < 0.) backward = 1;
00220   
00221   double sphi   = sin(momDir.phi());
00222   double cphi   = cos(momDir.phi());
00223   double ctheta = cos(momDir.theta());
00224   double stheta = sin(momDir.theta());
00225 
00226 #ifdef DebugLog
00227   G4ThreeVector localPos = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
00228   double zoff   = localPos.z() + 0.5*gpar[1];
00229   //  if (zoff < 0) zoff = 0;
00230   edm::LogInfo("HFShower") << "HFShowerLibrary: getHits " << partType
00231                            << " of energy " << pin/GeV << " GeV"
00232                            << "  dir.orts " << momDir.x() << ", " <<momDir.y() 
00233                            << ", " << momDir.z() << "  Pos x,y,z = " 
00234                            << hitPoint.x() << "," << hitPoint.y() << "," 
00235                            << hitPoint.z() << " (" << zoff 
00236                            << ")   sphi,cphi,stheta,ctheta  = " << sphi 
00237                            << ","  << cphi << ", " << stheta << "," << ctheta; 
00238 #endif    
00239                        
00240   if (parCode == emPDG || parCode == epPDG || parCode == gammaPDG ) {
00241     if (pin<pmom[nMomBin-1]) {
00242       interpolate(0, pin);
00243     } else {
00244       extrapolate(0, pin);
00245     }
00246   } else {
00247     if (pin<pmom[nMomBin-1]) {
00248       interpolate(1, pin);
00249     } else {
00250       extrapolate(1, pin);
00251     }
00252   }
00253     
00254   int nHit = 0;
00255   HFShowerLibrary::Hit oneHit;
00256   for (int i = 0; i < npe; i++) {
00257     double zv = std::abs(pe[i].z()); // abs local z  
00258 #ifdef DebugLog
00259     LogDebug("HFShower") << "HFShowerLibrary: Hit " << i << " " << pe[i] << " zv " << zv;
00260 #endif
00261     if (zv <= gpar[1] && pe[i].lambda() > 0 && 
00262         (pe[i].z() >= 0 || (zv > gpar[0] && (!onlyLong)))) {
00263       int depth = 1;
00264       if (onlyLong) {
00265       } else if (backward == 0) {    // fully valid only for "front" particles
00266         if (pe[i].z() < 0) depth = 2;// with "front"-simulated shower lib.
00267       } else {                       // for "backward" particles - almost equal
00268         double r = G4UniformRand();  // share between L and S fibers
00269         if (r > 0.5) depth = 2;
00270       } 
00271       
00272 
00273       // Updated coordinate transformation from local
00274       //  back to global using two Euler angles: phi and theta
00275       double pex = pe[i].x();
00276       double pey = pe[i].y();
00277 
00278       double xx = pex*ctheta*cphi - pey*sphi + zv*stheta*cphi; 
00279       double yy = pex*ctheta*sphi + pey*cphi + zv*stheta*sphi;
00280       double zz = -pex*stheta + zv*ctheta;
00281 
00282       G4ThreeVector pos  = hitPoint + G4ThreeVector(xx,yy,zz);
00283       zv = std::abs(pos.z()) - gpar[4] - 0.5*gpar[1];
00284       G4ThreeVector lpos = G4ThreeVector(pos.x(),pos.y(),zv);
00285 
00286       zv = fibre->zShift(lpos,depth,0);     // distance to PMT !
00287 
00288       double r  = pos.perp();
00289       double p  = fibre->attLength(pe[i].lambda());
00290       double fi = pos.phi();
00291       if (fi < 0) fi += twopi;
00292       int    isect = int(fi/dphi) + 1;
00293       isect        = (isect + 1) / 2;
00294       double dfi   = ((isect*2-1)*dphi - fi);
00295       if (dfi < 0) dfi = -dfi;
00296       double dfir  = r * sin(dfi);
00297 #ifdef DebugLog
00298       LogDebug("HFShower") << "HFShowerLibrary: Position shift " << xx 
00299                            << ", " << yy << ", "  << zz << ": " << pos 
00300                            << " R " << r << " Phi " << fi << " Section " 
00301                            << isect << " R*Dfi " << dfir << " Dist " << zv;
00302 #endif
00303       zz           = std::abs(pos.z());
00304       double r1    = G4UniformRand();
00305       double r2    = G4UniformRand();
00306       double r3    = -9999.;
00307       if (backward)     r3    = G4UniformRand();
00308       if (!applyFidCut) dfir += gpar[5];
00309 
00310 #ifdef DebugLog
00311       LogDebug("HFShower") << "HFShowerLibrary: rLimits " << rInside(r)
00312                            << " attenuation " << r1 <<":" << exp(-p*zv) 
00313                            << " r2 " << r2 << " r3 " << r3 << " rDfi "  
00314                            << gpar[5] << " zz " 
00315                            << zz << " zLim " << gpar[4] << ":" 
00316                            << gpar[4]+gpar[1] << "\n"
00317                            << "  rInside(r) :" << rInside(r) 
00318                            << "  r1 <= exp(-p*zv) :" <<  (r1 <= exp(-p*zv))
00319                            << "  r2 <= probMax :"    <<  (r2 <= probMax)
00320                            << "  r3 <= backProb :"   <<  (r3 <= backProb) 
00321                            << "  dfir > gpar[5] :"   <<  (dfir > gpar[5])
00322                            << "  zz >= gpar[4] :"    <<  (zz >= gpar[4])
00323                            << "  zz <= gpar[4]+gpar[1] :" 
00324                            << (zz <= gpar[4]+gpar[1]);   
00325 #endif
00326       if (rInside(r) && r1 <= exp(-p*zv) && r2 <= probMax && dfir > gpar[5] &&
00327           zz >= gpar[4] && zz <= gpar[4]+gpar[1] && r3 <= backProb &&
00328           (depth != 2 || zz >= gpar[4]+gpar[0])) {
00329         oneHit.position = pos;
00330         oneHit.depth    = depth;
00331         oneHit.time     = (tSlice+(pe[i].t())+(fibre->tShift(lpos,depth,1)));
00332         hit.push_back(oneHit);
00333 #ifdef DebugLog
00334         LogDebug("HFShower") << "HFShowerLibrary: Final Hit " << nHit 
00335                              <<" position " << (hit[nHit].position) <<" Depth "
00336                              << (hit[nHit].depth) <<" Time " << tSlice << ":"
00337                              << pe[i].t() << ":" <<fibre->tShift(lpos,depth,1)
00338                              << ":" << (hit[nHit].time);
00339 #endif
00340         nHit++;
00341       }
00342 #ifdef DebugLog
00343       else  LogDebug("HFShower") << "HFShowerLibrary: REJECTED !!!";
00344 #endif
00345       if (onlyLong && zz >= gpar[4]+gpar[0] && zz <= gpar[4]+gpar[1]) {
00346         r1    = G4UniformRand();
00347         r2    = G4UniformRand();
00348         if (rInside(r) && r1 <= exp(-p*zv) && r2 <= probMax && dfir > gpar[5]){
00349           oneHit.position = pos;
00350           oneHit.depth    = 2;
00351           oneHit.time     = (tSlice+(pe[i].t())+(fibre->tShift(lpos,2,1)));
00352           hit.push_back(oneHit);
00353 #ifdef DebugLog
00354           LogDebug("HFShower") << "HFShowerLibrary: Final Hit " << nHit 
00355                                <<" position " << (hit[nHit].position) <<" Depth "
00356                                << (hit[nHit].depth) <<" Time " << tSlice << ":"
00357                                << pe[i].t() << ":" << fibre->tShift(lpos,2,1)
00358                                << ":" << (hit[nHit].time);
00359 #endif
00360           nHit++;
00361         }
00362       }
00363     }
00364   }
00365 
00366 #ifdef DebugLog
00367   LogDebug("HFShower") << "HFShowerLibrary: Total Hits " << nHit
00368                        << " out of " << npe << " PE";
00369 #endif
00370   if (nHit > npe && !onlyLong)
00371     edm::LogWarning("HFShower") << "HFShowerLibrary: Hit buffer " << npe 
00372                                 << " smaller than " << nHit << " Hits";
00373   return hit;
00374 
00375 }
00376 
00377 bool HFShowerLibrary::rInside(double r) {
00378 
00379   if (r >= rMin && r <= rMax) return true;
00380   else                        return false;
00381 }
00382 
00383 void HFShowerLibrary::getRecord(int type, int record) {
00384 
00385   int nrc     = record-1;
00386   photon.clear();
00387   if (type > 0) {
00388     hadBranch->SetAddress(&photon);
00389     hadBranch->GetEntry(nrc);
00390   } else {
00391     emBranch->SetAddress(&photon);
00392     emBranch->GetEntry(nrc);
00393   }
00394 #ifdef DebugLog
00395   int nPhoton = photon.size();
00396   LogDebug("HFShower") << "HFShowerLibrary::getRecord: Record " << record
00397                        << " of type " << type << " with " << nPhoton 
00398                        << " photons";
00399   for (int j = 0; j < nPhoton; j++) 
00400     LogDebug("HFShower") << "Photon " << j << " " << photon[j];
00401 #endif
00402 }
00403 
00404 void HFShowerLibrary::loadEventInfo(TBranch* branch) {
00405 
00406   std::vector<HFShowerLibraryEventInfo> eventInfoCollection;
00407   branch->SetAddress(&eventInfoCollection);
00408   branch->GetEntry(0);
00409   edm::LogInfo("HFShower") << "HFShowerLibrary::loadEventInfo loads "
00410                            << " EventInfo Collection of size "
00411                            << eventInfoCollection.size() << " records";
00412   totEvents   = eventInfoCollection[0].totalEvents();
00413   nMomBin     = eventInfoCollection[0].numberOfBins();
00414   evtPerBin   = eventInfoCollection[0].eventsPerBin();
00415   libVers     = eventInfoCollection[0].showerLibraryVersion();
00416   listVersion = eventInfoCollection[0].physListVersion();
00417   pmom        = eventInfoCollection[0].energyBins();
00418   for (unsigned int i=0; i<pmom.size(); i++) 
00419     pmom[i] *= GeV;
00420 }
00421 
00422 void HFShowerLibrary::interpolate(int type, double pin) {
00423 
00424 #ifdef DebugLog
00425   LogDebug("HFShower") << "HFShowerLibrary:: Interpolate for Energy " <<pin/GeV
00426                        << " GeV with " << nMomBin << " momentum bins and " 
00427                        << evtPerBin << " entries/bin -- total " << totEvents;
00428 #endif
00429   int irc[2];
00430   double w = 0.;
00431   double r = G4UniformRand();
00432 
00433   if (pin<pmom[0]) {
00434     w = pin/pmom[0];
00435     irc[1] = int(evtPerBin*r) + 1;
00436     irc[0] = 0;
00437   } else {
00438     for (int j=0; j<nMomBin-1; j++) {
00439       if (pin >= pmom[j] && pin < pmom[j+1]) {
00440         w = (pin-pmom[j])/(pmom[j+1]-pmom[j]);
00441         if (j == nMomBin-2) { 
00442           irc[1] = int(evtPerBin*0.5*r);
00443         } else {
00444           irc[1] = int(evtPerBin*r);
00445         }
00446         irc[1] += (j+1)*evtPerBin + 1;
00447         r = G4UniformRand();
00448         irc[0] = int(evtPerBin*r) + 1 + j*evtPerBin;
00449         if (irc[0]<0) {
00450           edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[0] = "
00451                                       << irc[0] << " now set to 0";
00452           irc[0] = 0;
00453         } else if (irc[0] > totEvents) {
00454           edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[0] = "
00455                                       << irc[0] << " now set to "<< totEvents;
00456           irc[0] = totEvents;
00457         }
00458       }
00459     }
00460   }
00461   if (irc[1]<1) {
00462     edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[1] = " 
00463                                 << irc[1] << " now set to 1";
00464     irc[1] = 1;
00465   } else if (irc[1] > totEvents) {
00466     edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[1] = " 
00467                                 << irc[1] << " now set to "<< totEvents;
00468     irc[1] = totEvents;
00469   }
00470 
00471 #ifdef DebugLog
00472   LogDebug("HFShower") << "HFShowerLibrary:: Select records " << irc[0] 
00473                        << " and " << irc[1] << " with weights " << 1-w 
00474                        << " and " << w;
00475 #endif
00476   pe.clear(); 
00477   npe       = 0;
00478   int npold = 0;
00479   for (int ir=0; ir < 2; ir++) {
00480     if (irc[ir]>0) {
00481       getRecord (type, irc[ir]);
00482       int nPhoton = photon.size();
00483       npold      += nPhoton;
00484       for (int j=0; j<nPhoton; j++) {
00485         r = G4UniformRand();
00486         if ((ir==0 && r > w) || (ir > 0 && r < w)) {
00487           storePhoton (j);
00488         }
00489       }
00490     }
00491   }
00492 
00493   if (npe > npold || (npold == 0 && irc[0] > 0)) 
00494     edm::LogWarning("HFShower") << "HFShowerLibrary: Interpolation Warning =="
00495                                 << " records " << irc[0] << " and " << irc[1]
00496                                 << " gives a buffer of " << npold 
00497                                 << " photons and fills " << npe << " *****";
00498 #ifdef DebugLog
00499   else
00500     LogDebug("HFShower") << "HFShowerLibrary: Interpolation == records " 
00501                          << irc[0] << " and " << irc[1] << " gives a "
00502                          << "buffer of " << npold << " photons and fills "
00503                          << npe << " PE";
00504   for (int j=0; j<npe; j++)
00505     LogDebug("HFShower") << "Photon " << j << " " << pe[j];
00506 #endif
00507 }
00508 
00509 void HFShowerLibrary::extrapolate(int type, double pin) {
00510 
00511   int nrec   = int(pin/pmom[nMomBin-1]);
00512   double w   = (pin - pmom[nMomBin-1]*nrec)/pmom[nMomBin-1];
00513   nrec++;
00514 #ifdef DebugLog
00515   LogDebug("HFShower") << "HFShowerLibrary:: Extrapolate for Energy " << pin 
00516                        << " GeV with " << nMomBin << " momentum bins and " 
00517                        << evtPerBin << " entries/bin -- total " << totEvents 
00518                        << " using " << nrec << " records";
00519 #endif
00520   std::vector<int> irc(nrec);
00521 
00522   for (int ir=0; ir<nrec; ir++) {
00523     double r = G4UniformRand();
00524     irc[ir] = int(evtPerBin*0.5*r) +(nMomBin-1)*evtPerBin + 1;
00525     if (irc[ir]<1) {
00526       edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[" << ir 
00527                                   << "] = " << irc[ir] << " now set to 1";
00528       irc[ir] = 1;
00529     } else if (irc[ir] > totEvents) {
00530       edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[" << ir 
00531                                   << "] = " << irc[ir] << " now set to "
00532                                   << totEvents;
00533       irc[ir] = totEvents;
00534 #ifdef DebugLog
00535     } else {
00536       LogDebug("HFShower") << "HFShowerLibrary::Extrapolation use irc[" 
00537                            << ir  << "] = " << irc[ir];
00538 #endif
00539     }
00540   }
00541 
00542   pe.clear(); 
00543   npe       = 0;
00544   int npold = 0;
00545   for (int ir=0; ir<nrec; ir++) {
00546     if (irc[ir]>0) {
00547       getRecord (type, irc[ir]);
00548       int nPhoton = photon.size();
00549       npold      += nPhoton;
00550       for (int j=0; j<nPhoton; j++) {
00551         double r = G4UniformRand();
00552         if (ir != nrec-1 || r < w) {
00553           storePhoton (j);
00554         }
00555       }
00556 #ifdef DebugLog
00557       LogDebug("HFShower") << "HFShowerLibrary: Record [" << ir << "] = " 
00558                            << irc[ir] << " npold = " << npold;
00559 #endif
00560     }
00561   }
00562 #ifdef DebugLog
00563   edm::LogInfo("HFShower") << "HFShowerLibrary:: uses " << npold << " photons";
00564 #endif
00565 
00566   if (npe > npold || npold == 0)
00567     edm::LogWarning("HFShower") << "HFShowerLibrary: Extrapolation Warning == "
00568                                 << nrec << " records " << irc[0] << ", " 
00569                                 << irc[1] << ", ... gives a buffer of " <<npold
00570                                 << " photons and fills " << npe 
00571                                 << " *****";
00572 #ifdef DebugLog
00573   else
00574     LogDebug("HFShower") << "HFShowerLibrary: Extrapolation == " << nrec
00575                          << " records " << irc[0] << ", " << irc[1] 
00576                          << ", ... gives a buffer of " << npold 
00577                          << " photons and fills " << npe << " PE";
00578   for (int j=0; j<npe; j++)
00579     LogDebug("HFShower") << "Photon " << j << " " << pe[j];
00580 #endif
00581 }
00582 
00583 void HFShowerLibrary::storePhoton(int j) {
00584 
00585   pe.push_back(photon[j]);
00586 #ifdef DebugLog
00587   LogDebug("HFShower") << "HFShowerLibrary: storePhoton " << j << " npe " 
00588                        << npe << " " << pe[npe];
00589 #endif
00590   npe++;
00591 }
00592 
00593 std::vector<double> HFShowerLibrary::getDDDArray(const std::string & str, 
00594                                                  const DDsvalues_type & sv, 
00595                                                  int & nmin) {
00596 
00597 #ifdef DebugLog
00598   LogDebug("HFShower") << "HFShowerLibrary:getDDDArray called for " << str 
00599                        << " with nMin " << nmin;
00600 #endif
00601   DDValue value(str);
00602   if (DDfetch(&sv,value)) {
00603 #ifdef DebugLog
00604     LogDebug("HFShower") << value;
00605 #endif
00606     const std::vector<double> & fvec = value.doubles();
00607     int nval = fvec.size();
00608     if (nmin > 0) {
00609       if (nval < nmin) {
00610         edm::LogError("HFShower") << "HFShowerLibrary : # of " << str 
00611                                   << " bins " << nval << " < " << nmin 
00612                                   << " ==> illegal";
00613         throw cms::Exception("Unknown", "HFShowerLibrary")
00614           << "nval < nmin for array " << str << "\n";
00615       }
00616     } else {
00617       if (nval < 2) {
00618         edm::LogError("HFShower") << "HFShowerLibrary : # of " << str 
00619                                   << " bins " << nval << " < 2 ==> illegal"
00620                                   << " (nmin=" << nmin << ")";
00621         throw cms::Exception("Unknown", "HFShowerLibrary")
00622           << "nval < 2 for array " << str << "\n";
00623       }
00624     }
00625     nmin = nval;
00626 
00627     return fvec;
00628   } else {
00629     edm::LogError("HFShower") << "HFShowerLibrary : cannot get array " << str;
00630     throw cms::Exception("Unknown", "HFShowerLibrary") 
00631       << "cannot get array " << str << "\n";
00632   }
00633 }