CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/FastSimulation/ShowerDevelopment/src/HFShower.cc

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