CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/FastSimulation/ShowerDevelopment/src/HDShower.cc

Go to the documentation of this file.
00001 //FastSimulation Headers
00002 #include "FastSimulation/ShowerDevelopment/interface/HDShower.h"
00003 //#include "FastSimulation/Utilities/interface/Histos.h"
00004 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00005 
00006 
00008 // What's this?
00009 //#include "FastSimulation/FamosCalorimeters/interface/FASTCalorimeter.h"
00010 
00011 #include "FastSimulation/CaloHitMakers/interface/EcalHitMaker.h"
00012 #include "FastSimulation/CaloHitMakers/interface/HcalHitMaker.h"
00013 
00014 // CMSSW headers
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 
00018 // And This???? Doesn't seem to be needed
00019 // #include "Calorimetry/CaloDetector/interface/CellGeometry.h"
00020 
00021 #include <cmath>
00022 
00023 // number attempts for transverse distribution if exit on a spec. condition
00024 #define infinity 10000
00025 // debugging flag ( 0, 1, 2, 3)
00026 #define debug 0
00027 
00028 using namespace edm;
00029 
00030 HDShower::HDShower(const RandomEngine* engine,
00031                    HDShowerParametrization* myParam, 
00032                    EcalHitMaker* myGrid,
00033                    HcalHitMaker* myHcalHitMaker,
00034                    int onECAL,
00035                    double epart)
00036   : theParam(myParam), 
00037     theGrid(myGrid),
00038     theHcalHitMaker(myHcalHitMaker),
00039     onEcal(onECAL),
00040     e(epart),
00041     random(engine)
00042 { 
00043   // To get an access to constants read in FASTCalorimeter
00044   //  FASTCalorimeter * myCalorimeter= FASTCalorimeter::instance();
00045 
00046   // Values taken from FamosGeneric/FamosCalorimeter/src/FASTCalorimeter.cc
00047   lossesOpt       = myParam->hsParameters()->getHDlossesOpt();
00048   nDepthSteps     = myParam->hsParameters()->getHDnDepthSteps();
00049   nTRsteps        = myParam->hsParameters()->getHDnTRsteps();
00050   transParam      = myParam->hsParameters()->getHDtransParam();
00051   eSpotSize       = myParam->hsParameters()->getHDeSpotSize();
00052   depthStep       = myParam->hsParameters()->getHDdepthStep();
00053   criticalEnergy  = myParam->hsParameters()->getHDcriticalEnergy();
00054   maxTRfactor     = myParam->hsParameters()->getHDmaxTRfactor();
00055   balanceEH       = myParam->hsParameters()->getHDbalanceEH();
00056   hcalDepthFactor = myParam->hsParameters()->getHDhcalDepthFactor();
00057 
00058   // Special tr.size fluctuations 
00059   transParam *= (1. + random->flatShoot()); 
00060 
00061   // Special long. fluctuations
00062   hcalDepthFactor +=  0.05 * (2.* random->flatShoot() - 1.);
00063 
00064   transFactor = 1.;   // normally 1, in HF - might be smaller 
00065                       // to take into account
00066                       // a narrowness of the HF shower (Cherenkov light) 
00067 
00068   // simple protection ...
00069   if(e < 0) e = 0.;
00070 
00071   // Get the Famos Histos pointer
00072   //  myHistos = FamosHistos::instance();
00073   //  std::cout << " Hello FamosShower " << std::endl;
00074   
00075   theECALproperties = theParam->ecalProperties();
00076   theHCALproperties = theParam->hcalProperties();
00077 
00078   double emax = theParam->emax(); 
00079   double emid = theParam->emid(); 
00080   double emin = theParam->emin(); 
00081   double effective = e;
00082     
00083 
00084   if( e < emid ) {
00085     theParam->setCase(1);
00086     // avoid "underflow" below Emin (for parameters calculation only)
00087     if(e < emin) effective = emin; 
00088   } 
00089   else 
00090     theParam->setCase(2);
00091  
00092   // Avoid "overflow" beyond Emax (for parameters)
00093   if(effective > 0.5 * emax) {
00094     eSpotSize *= 2.5;
00095     if(effective > emax) {
00096       effective = emax; 
00097       eSpotSize *= 4.; 
00098       depthStep *= 2.;
00099       if(effective > 1000.)
00100       eSpotSize *= 2.; 
00101     }
00102   }
00103 
00104 
00105   if(debug == 2 )
00106     LogDebug("FastCalorimetry") <<  " HDShower : " << std::endl 
00107          << "       Energy   "  <<               e << std::endl     
00108          << "      lossesOpt "  <<       lossesOpt << std::endl  
00109          << "    nDepthSteps "  <<     nDepthSteps << std::endl
00110          << "       nTRsteps "  <<        nTRsteps << std::endl
00111          << "     transParam "  <<      transParam << std::endl
00112          << "      eSpotSize "  <<       eSpotSize << std::endl
00113          << " criticalEnergy "  <<  criticalEnergy << std::endl
00114          << "    maxTRfactor "  <<     maxTRfactor << std::endl
00115          << "      balanceEH "  <<       balanceEH << std::endl
00116          << "hcalDepthFactor "  << hcalDepthFactor << std::endl;
00117 
00118 
00119   double alpEM1 = theParam->alpe1();
00120   double alpEM2 = theParam->alpe2();
00121 
00122   double betEM1 = theParam->bete1();
00123   double betEM2 = theParam->bete2();
00124 
00125   double alpHD1 = theParam->alph1();
00126   double alpHD2 = theParam->alph2();
00127 
00128   double betHD1 = theParam->beth1();
00129   double betHD2 = theParam->beth2();
00130 
00131   double part1 = theParam->part1();
00132   double part2 = theParam->part2();
00133 
00134   aloge = std::log(effective);
00135  
00136   double edpar = (theParam->e1() + aloge * theParam->e2()) * effective;
00137   double aedep = std::log(edpar);
00138 
00139   if(debug == 2)
00140     LogDebug("FastCalorimetry") << " HDShower : " << std::endl
00141          << "     edpar " <<   edpar << "   aedep " << aedep << std::endl 
00142          << "    alpEM1 " <<  alpEM1 << std::endl  
00143          << "    alpEM2 " <<  alpEM2 << std::endl  
00144          << "    betEM1 " <<  betEM1 << std::endl  
00145          << "    betEM2 " <<  betEM2 << std::endl  
00146          << "    alpHD1 " <<  alpHD1 << std::endl  
00147          << "    alpHD2 " <<  alpHD2 << std::endl  
00148          << "    betHD1 " <<  betHD1 << std::endl  
00149          << "    betHD2 " <<  betHD2 << std::endl  
00150          << "     part1 " <<   part1 << std::endl  
00151          << "     part2 " <<   part2 << std::endl; 
00152 
00153   // private members to set
00154   theR1  = theParam->r1();
00155   theR2  = theParam->r2();
00156   theR3  = theParam->r3();
00157 
00158   alpEM  = alpEM1 + alpEM2 * aedep;
00159   tgamEM = tgamma(alpEM);
00160   betEM  = betEM1 - betEM2 * aedep;
00161   alpHD  = alpHD1 + alpHD2 * aedep;
00162   tgamHD = tgamma(alpHD);
00163   betHD  = betHD1 - betHD2 * aedep;
00164   part   = part1  -  part2 * aedep;
00165   if(part > 1.) part = 1.;          // protection - just in case of 
00166 
00167   if(debug  == 2 )
00168     LogDebug("FastCalorimetry") << " HDShower : " << std::endl 
00169          << "    alpEM " <<  alpEM << std::endl  
00170          << "   tgamEM " << tgamEM << std::endl
00171          << "    betEM " <<  betEM << std::endl
00172          << "    alpHD " <<  alpHD << std::endl
00173          << "   tgamHD " << tgamHD << std::endl
00174          << "    betHD " <<  betHD << std::endl
00175          << "     part " <<   part << std::endl;
00176        
00177 
00178   if(onECAL){
00179     lambdaEM = theParam->ecalProperties()->interactionLength();
00180     x0EM     = theParam->ecalProperties()->radLenIncm();
00181   } 
00182   else {
00183     lambdaEM = 0.;
00184     x0EM     = 0.;
00185   }
00186   lambdaHD = theParam->hcalProperties()->interactionLength();
00187   x0HD     = theParam->hcalProperties()->radLenIncm();
00188 
00189   if(debug == 2)
00190     LogDebug("FastCalorimetry") << " HDShower e " << e        << std::endl
00191          << "          x0EM = " << x0EM     << std::endl 
00192          << "          x0HD = " << x0HD     << std::endl 
00193          << "         lamEM = " << lambdaEM << std::endl
00194          << "         lamHD = " << lambdaHD << std::endl;
00195 
00196 
00197   // Starting point of the shower
00198   // try first with ECAL lambda  
00199 
00200   double sum1   = 0.;  // lambda path from the ECAL/HF entrance;  
00201   double sum2   = 0.;  // lambda path from the interaction point;
00202   double sum3   = 0.;  // x0     path from the interaction point;
00203   int  nsteps   = 0;   // full number of longitudinal steps (counter);
00204 
00205   int nmoresteps;      // how many longitudinal steps in addition to 
00206                        // one (if interaction happens there) in ECAL
00207 
00208   mip = 1;             // just to initiate particle as MIP in ECAL    
00209 
00210   if(e < criticalEnergy ) nmoresteps = 1;  
00211   else                    nmoresteps = nDepthSteps;
00212   
00213   double depthECAL  = 0.;
00214   double depthGAP   = 0.;
00215   double depthGAPx0 = 0.; 
00216   if(onECAL ) {
00217     depthECAL  = theGrid->ecalTotalL0();         // ECAL depth segment
00218     depthGAP   = theGrid->ecalHcalGapTotalL0();  // GAP  depth segment
00219     depthGAPx0 = theGrid->ecalHcalGapTotalX0();  // GAP  depth x0
00220   }
00221   
00222   double depthHCAL   = theGrid->hcalTotalL0();    // HCAL depth segment
00223   double depthToHCAL = depthECAL + depthGAP;    
00224 
00225   //---------------------------------------------------------------------------
00226   // Depth simulation & various protections, among them
00227   // if too deep - get flat random in the allowed region
00228   // if no HCAL material behind - force to deposit in ECAL
00229   double maxDepth    = depthToHCAL + depthHCAL - 1.1 * depthStep;
00230   double depthStart  = std::log(1./random->flatShoot()); // starting point lambda unts
00231 
00232   if(e < emin) {
00233     if(debug)
00234       LogDebug("FastCalorimetry") << " FamosHDShower : e <emin ->  depthStart = 0" << std::endl; 
00235     depthStart = 0.;
00236   }
00237  
00238   if(depthStart > maxDepth) {
00239     if(debug) LogDebug("FastCalorimetry") << " FamosHDShower : depthStart too big ...   = " 
00240                    << depthStart << std::endl; 
00241     
00242     depthStart = maxDepth *  random->flatShoot();
00243     if( depthStart < 0.) depthStart = 0.;
00244     if(debug) LogDebug("FastCalorimetry") << " FamosHDShower : depthStart re-calculated = " 
00245                    << depthStart << std::endl; 
00246   }
00247   
00248   if( onECAL && e < emid ) {
00249     if(depthECAL > depthStep && (depthECAL - depthStart)/depthECAL > 0.2) {
00250       
00251       depthStart = 0.5 * depthECAL * random->flatShoot();
00252       if(debug) 
00253         LogDebug("FastCalorimetry") << " FamosHDShower : small energy, "
00254              << " depthStart reduced to = " << depthStart << std::endl; 
00255       
00256     }
00257   }
00258 
00259   if( depthHCAL < depthStep) {
00260     if(debug) LogDebug("FastCalorimetry") << " FamosHDShower : depthHCAL  too small ... = " 
00261                    << depthHCAL << " depthStart -> forced to 0 !!!" 
00262                    << std::endl;
00263     depthStart = 0.;    
00264     nmoresteps = 0;
00265     
00266     if(depthECAL < depthStep) {
00267       nsteps = -1;
00268       LogInfo("FastCalorimetry") << " FamosHDShower : too small ECAL and HCAL depths - " 
00269            << " particle is lost !!! " << std::endl; 
00270     }
00271   }
00272 
00273 
00274 
00275   if(debug)
00276     LogDebug("FastCalorimetry") << " FamosHDShower  depths(lam) - "  << std::endl 
00277          << "          ECAL = " << depthECAL  << std::endl
00278          << "           GAP = " << depthGAP   << std::endl
00279          << "          HCAL = " << depthHCAL  << std::endl
00280          << "  starting point = " << depthStart << std::endl; 
00281 
00282   if( onEcal ) {
00283     if(debug) LogDebug("FastCalorimetry") << " FamosHDShower : onECAL" << std::endl;
00284     if(depthStart < depthECAL) {
00285       if(debug) LogDebug("FastCalorimetry") << " FamosHDShower : depthStart < depthECAL" << std::endl;
00286       if(depthECAL > depthStep && (depthECAL - depthStart)/depthECAL > 0.1) {
00287         if(debug) LogDebug("FastCalorimetry") << " FamosHDShower : enough space to make ECAL step"
00288                        << std::endl;
00289         //  ECAL - one step
00290         nsteps++; 
00291         sum1   += depthECAL;             // at the end of step
00292         sum2   += depthECAL-depthStart; 
00293         sum3   += sum2 * lambdaEM / x0EM;
00294         lamtotal.push_back(sum1);
00295         lamdepth.push_back(sum2);
00296         lamcurr.push_back(lambdaEM);
00297         lamstep.push_back(depthECAL-depthStart);
00298         x0depth.push_back(sum3);
00299         x0curr.push_back(x0EM);
00300         detector.push_back(1);
00301         mip = 0;        
00302 
00303         if(debug) LogDebug("FastCalorimetry") << " FamosHDShower : " << " in ECAL sum1, sum2 "
00304                        << sum1 << " " << sum2 << std::endl;
00305         
00306         //                           // Gap - no additional step after ECAL
00307         //                           // just move further to HCAL over the gap
00308         sum1   += depthGAP;          
00309         sum2   += depthGAP; 
00310         sum3   += depthGAPx0;
00311       }
00312       // Just shift starting point to HCAL
00313       else { 
00314         //      cout << " FamosHDShower : not enough space to make ECAL step" << std::endl;
00315         if(debug)  LogDebug("FastCalorimetry") << " FamosHDShower : goto HCAL" << std::endl;
00316 
00317         depthStart = depthToHCAL;
00318         sum1 += depthStart;     
00319       }
00320     }
00321     else { // GAP or HCAL   
00322       
00323       if(depthStart >= depthECAL &&  depthStart < depthToHCAL ) {
00324       depthStart = depthToHCAL; // just a shift to HCAL for simplicity
00325       }
00326       sum1 += depthStart;
00327 
00328       if(debug) LogDebug("FastCalorimetry") << " FamosHDShower : goto HCAL" << std::endl;
00329     }
00330   }
00331   else {   // Forward 
00332     if(debug)  LogDebug("FastCalorimetry") << " FamosHDShower : forward" << std::endl;
00333     sum1 += depthStart;
00334     transFactor = 0.5;   // makes narower tresverse size of shower     
00335   }
00336  
00337   for (int i = 0; i < nmoresteps ; i++) {
00338     sum1 += depthStep;
00339     if (sum1 > (depthECAL + depthGAP + depthHCAL)) break; 
00340     sum2 += depthStep;
00341     sum3 += sum2 * lambdaHD / x0HD;
00342     lamtotal.push_back(sum1);
00343     lamdepth.push_back(sum2);
00344     lamcurr.push_back(lambdaHD);
00345     lamstep.push_back(depthStep);
00346     x0depth.push_back(sum3);
00347     x0curr.push_back(x0HD);
00348     detector.push_back(3);
00349     nsteps++;
00350   }
00351   
00352 
00353   // Make fractions of energy and transverse radii at each step 
00354 
00355   if(nsteps > 0) { makeSteps(nsteps); }
00356 
00357 }
00358 
00359 void HDShower::makeSteps(int nsteps) {
00360 
00361   double sumes = 0.;
00362   double sum   = 0.;
00363   std::vector<double> temp;
00364 
00365 
00366   if(debug)
00367     LogDebug("FastCalorimetry") << " FamosHDShower::makeSteps - " 
00368        << " nsteps required : " << nsteps << std::endl;
00369 
00370   int count = 0;
00371   for (int i = 0; i < nsteps; i++) {    
00372     
00373     double deplam = lamdepth[i] - 0.5 * lamstep[i];
00374     double depx0  = x0depth[i]  - 0.5 * lamstep[i] / x0curr[i]; 
00375     double     x = betEM * depx0;
00376     double     y = betHD * deplam;
00377     
00378     if(debug == 2)
00379       LogDebug("FastCalorimetry") << " FamosHDShower::makeSteps " 
00380                                   << " - step " << i
00381                                   << "   depx0, x = " << depx0 << ", " << x 
00382                                   << "   deplam, y = " << deplam << ", "
00383                                   << y << std::endl;
00384     
00385     double est = (part * betEM * gam(x,alpEM) * lamcurr[i] /
00386                   (x0curr[i] * tgamEM) + 
00387                   (1.-part) * betHD * gam(y,alpHD) / tgamHD) * lamstep[i];
00388     
00389     // protection ...
00390     if(est < 0.) {
00391       LogDebug("FastCalorimetry") << "*** FamosHDShower::makeSteps " << " - negative step energy !!!" 
00392            << std::endl;
00393       est = 0.;
00394       break ; 
00395     }
00396 
00397     // for estimates only
00398     sum += est;
00399     int nPest = (int) (est * e / sum / eSpotSize) ;
00400 
00401     if(debug == 2)
00402       LogDebug("FastCalorimetry") << " FamosHDShower::makeSteps - nPoints estimate = " 
00403            <<  nPest << std::endl;
00404 
00405     if(nPest <= 1 && count !=0 ) break;
00406 
00407     // good step - to proceed
00408 
00409     temp.push_back(est);
00410     sumes += est;
00411     
00412     rlamStep.push_back(transParam * (theR1 + (theR2 - theR3 * aloge))
00413                        * deplam * transFactor); 
00414     count ++;
00415   }
00416 
00417   // fluctuations in ECAL and re-distribution of remaining energy in HCAL
00418   if(detector[0] == 1 && count > 1) {
00419     double oldECALenergy = temp[0];
00420     double oldHCALenergy = sumes - oldECALenergy ;
00421     double newECALenergy = 2. * sumes;
00422     for (int i = 0; newECALenergy > sumes && i < infinity; i++)
00423       newECALenergy = 2.* balanceEH * random->flatShoot() * oldECALenergy; 
00424      
00425     if(debug == 2)
00426       LogDebug("FastCalorimetry") << "*** FamosHDShower::makeSteps " << " ECAL fraction : old/new - "
00427            << oldECALenergy/sumes << "/" << newECALenergy/sumes << std::endl;
00428 
00429     temp[0] = newECALenergy;
00430     double newHCALenergy = sumes - newECALenergy;
00431     double newHCALreweight =  newHCALenergy / oldHCALenergy;
00432 
00433     for (int i = 1; i < count; i++) {
00434       temp[i] *= newHCALreweight;
00435     }
00436     
00437   }
00438 
00439   
00440   // final re-normalization of the energy fractions  
00441   for (int i = 0; i < count ; i++) {
00442     eStep.push_back(temp[i] * e / sumes );
00443     nspots.push_back((int)(eStep[i]/eSpotSize)+1);
00444 
00445    if(debug)
00446      LogDebug("FastCalorimetry") 
00447        << " step " << i
00448        << "  det: " << detector[i]   
00449        << "  xO and lamdepth at the end of step = " 
00450        << x0depth[i] << " "  
00451        << lamdepth[i] << "   Estep func = " <<  eStep[i]
00452        << "   Rstep = " << rlamStep[i] << "  Nspots = " <<  nspots[i]
00453        << "  espot = " <<  eStep[i] / (double)nspots[i]
00454        << std::endl; 
00455 
00456   }
00457 
00458 
00459   // The only step is in ECAL - let's make the size bigger ...  
00460   if(count == 1 and detector[0] == 1) rlamStep[0] *= 2.;
00461 
00462 
00463   if(debug) {
00464     if(eStep[0] > 0.95 * e && detector[0] == 1) 
00465       LogDebug("FastCalorimetry") << " FamosHDShower::makeSteps - " << "ECAL energy = " << eStep[0]
00466            << " out of total = " << e << std::endl;  
00467   }
00468 
00469 }
00470 
00471 bool HDShower::compute() {
00472   
00473   //  TimeMe theT("FamosHDShower::compute");
00474 
00475   bool status = false;
00476   int numLongit = eStep.size();
00477   if(debug)
00478     LogDebug("FastCalorimetry") << " FamosHDShower::compute - " 
00479             << " N_long.steps required : " << numLongit << std::endl;
00480 
00481   if(numLongit > 0) {
00482 
00483     status = true;    
00484     // Prepare the trsanverse probability function
00485     std::vector<double> Fhist;
00486     std::vector<double> rhist; 
00487     for (int j = 0; j < nTRsteps + 1; j++) {
00488       rhist.push_back(maxTRfactor * j / nTRsteps );  
00489       Fhist.push_back(transProb(maxTRfactor,1.,rhist[j]));
00490       if(debug == 3) 
00491         LogDebug("FastCalorimetry") << "indexFinder - i, Fhist[i] = " << j << " " << Fhist[j] << std::endl;
00492     }
00493     
00494     // Longitudinal steps
00495     for (int i = 0; i < numLongit ; i++) {
00496       
00497       double currentDepthL0 = lamtotal[i] - 0.5 * lamstep[i];
00498       // vary the longitudinal profile if needed
00499       if(detector[i] != 1) currentDepthL0 *= hcalDepthFactor;                     
00500       if(debug)
00501         LogDebug("FastCalorimetry") << " FamosHDShower::compute - detector = "
00502                                     << detector[i]
00503                                     << "  currentDepthL0 = " 
00504                                     << currentDepthL0 << std::endl;
00505       
00506       double maxTRsize   = maxTRfactor * rlamStep[i];     // in lambda units
00507       double rbinsize    = maxTRsize / nTRsteps; 
00508       double espot       = eStep[i] / (double)nspots[i];  // re-adjust espot
00509 
00510       if(espot > 2. || espot < 0. ) 
00511         LogDebug("FastCalorimetry") << " FamosHDShower::compute - unphysical espot = " 
00512              << espot << std::endl;
00513 
00514       int ecal = 0;
00515       if(detector[i] != 1) { 
00516         bool setHDdepth = theHcalHitMaker->setDepth(currentDepthL0);
00517         
00518         if(debug)
00519           LogDebug("FastCalorimetry") << " FamosHDShower::compute - status of " 
00520                << " theHcalHitMaker->setDepth(currentDepthL0) is " 
00521                << setHDdepth << std::endl;
00522         
00523         if(!setHDdepth) 
00524           {
00525             currentDepthL0 -= lamstep[i];
00526             setHDdepth =  theHcalHitMaker->setDepth(currentDepthL0);
00527           }
00528         if(!setHDdepth) continue;
00529         
00530         theHcalHitMaker->setSpotEnergy(espot);        
00531       }
00532       else {
00533         ecal = 1;
00534         bool status = theGrid->getPads(currentDepthL0);   
00535         
00536         if(debug)
00537           LogDebug("FastCalorimetry") << " FamosHDShower::compute - status of Grid = " 
00538                << status << std::endl;
00539         
00540         if(!status) continue; 
00541 
00542         int ntry = nspots[i] * 10;  
00543         if( ntry >= infinity ) {  // use max allowed in case of too many spots
00544           nspots[i] = 0.5 * infinity;  
00545           espot *= 0.1 * (double)ntry / double(nspots[i]);
00546         }     
00547         else {
00548           espot *= 0.1;  // fine-grain energy spots in ECAL
00549                          // to avoid false ECAL clustering 
00550           nspots[i] = ntry;
00551         }
00552 
00553 
00554         theGrid->setSpotEnergy(espot);
00555       }  
00556 
00557       
00558       // Transverse distribition
00559       int nok   = 0;                          // counter of OK  
00560       int count = 0;
00561       int inf   = infinity;
00562       if(lossesOpt) inf = nspots[i];          // if losses are enabled, otherwise
00563       // only OK points are counted ...
00564       if(nspots[i] > inf ){
00565         std::cout << " FamosHDShower::compute - at long.step " << i 
00566                   << "  too many spots required : "  <<  nspots[i] << " !!! " 
00567                   << std::endl;
00568       }
00569 
00570       for (int j = 0; j < inf; j++) {
00571         if(nok == nspots[i]) break;
00572         count ++;
00573         
00574         double prob   = random->flatShoot();
00575         int index     = indexFinder(prob,Fhist);
00576         double radius = rlamStep[i] * rhist[index] +
00577           random->flatShoot() * rbinsize; // in-bin  
00578         double phi = 2.*M_PI*random->flatShoot();
00579         
00580         if(debug == 2)
00581           LogDebug("FastCalorimetry") << std::endl << " FamosHDShower::compute " << " r = " << radius 
00582                << "    phi = " << phi << std::endl;
00583         
00584         bool result;
00585         if(ecal) {
00586           result = theGrid->addHit(radius,phi,0);
00587           
00588           if(debug == 2)
00589             LogDebug("FastCalorimetry") << " FamosHDShower::compute - " 
00590                  << " theGrid->addHit result = " 
00591                  << result << std::endl;
00592         }
00593         else {
00594           result = theHcalHitMaker->addHit(radius,phi,0); 
00595           
00596           if(debug == 2)
00597             LogDebug("FastCalorimetry") << " FamosHDShower::compute - " 
00598                  << " theHcalHitMaker->addHit result = " 
00599                  << result << std::endl;
00600         }    
00601         if(result) nok ++; 
00602         
00603       } // end of tranverse simulation
00604       if(count == infinity) { 
00605         if(debug)
00606           LogDebug("FastCalorimetry") 
00607             << " FamosHDShower::compute " << " maximum number of" 
00608             << " transverse points " << count << " is used !!!" << std::endl; 
00609       }
00610 
00611       if(debug)
00612         LogDebug("FastCalorimetry") 
00613           << " FamosHDShower::compute " << " long.step No." << i 
00614           << "   Ntry, Nok = " << count
00615           << " " << nok << std::endl; 
00616 
00617     } // end of longitudinal steps
00618   } // end of no steps
00619 
00620   return status;
00621 
00622 }
00623 
00624 int HDShower::indexFinder(double x, const std::vector<double> & Fhist) {
00625   // binary search in the vector of doubles
00626   int size = Fhist.size();
00627 
00628   int curr = size / 2;
00629   int step = size / 4;
00630   int iter;
00631   int prevdir = 0; 
00632   int actudir = 0; 
00633 
00634   for (iter = 0; iter < size ; iter++) {
00635 
00636     if( curr >= size || curr < 1 )
00637       LogWarning("FastCalorimetry") << " FamosHDShower::indexFinder - wrong current index = " 
00638            << curr << " !!!" << std::endl;
00639 
00640     if ((x <= Fhist[curr]) && (x > Fhist[curr-1])) break;
00641     prevdir = actudir;
00642     if(x > Fhist[curr]) {actudir =  1;}
00643     else                {actudir = -1;}
00644     if(prevdir * actudir < 0) { if(step > 1) step /= 2;}
00645     curr += actudir * step;
00646     if(curr > size) curr = size;
00647     else { if(curr < 1) {curr = 1;}}
00648  
00649     if(debug == 3)
00650       LogDebug("FastCalorimetry") << " indexFinder - end of iter." << iter 
00651            << " curr, F[curr-1], F[curr] = "
00652            << curr << " " << Fhist[curr-1] << " " << Fhist[curr] << std::endl;
00653     
00654   }
00655 
00656   if(debug == 3)
00657     LogDebug("FastCalorimetry") << " indexFinder x = " << x << "  found index = " << curr-1
00658          << std::endl;
00659 
00660 
00661   return curr-1;
00662 }