CMS 3D CMS Logo

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