CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/SimGeneral/GFlash/src/GflashHadronShowerProfile.cc

Go to the documentation of this file.
00001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00002 
00003 #include "SimGeneral/GFlash/interface/GflashHadronShowerProfile.h"
00004 #include "SimGeneral/GFlash/interface/GflashTrajectoryPoint.h"
00005 #include "SimGeneral/GFlash/interface/GflashHit.h"
00006 
00007 #include <CLHEP/GenericFunctions/IncompleteGamma.hh>
00008 #include <CLHEP/GenericFunctions/LogGamma.hh>
00009 #include <CLHEP/Units/PhysicalConstants.h>
00010 #include <CLHEP/Units/SystemOfUnits.h>
00011 #include <CLHEP/Random/Randomize.h>
00012 #include <CLHEP/Random/RandGaussQ.h>
00013 #include <CLHEP/Random/RandPoissonQ.h>
00014 
00015 #include <math.h>
00016 
00017 GflashHadronShowerProfile::GflashHadronShowerProfile(edm::ParameterSet parSet) : theParSet(parSet)
00018 {
00019   theBField = parSet.getParameter<double>("bField");
00020   theGflashHcalOuter = parSet.getParameter<bool>("GflashHcalOuter");
00021 
00022   theShowino = new GflashShowino();
00023   theHisto = GflashHistogram::instance();
00024 }
00025 
00026 GflashHadronShowerProfile::~GflashHadronShowerProfile() 
00027 {
00028   if(theShowino) delete theShowino;
00029 }
00030 
00031 void GflashHadronShowerProfile::initialize(int showerType, double energy, double globalTime,double charge, 
00032                                            Gflash3Vector &position,Gflash3Vector &momentum) { 
00033 
00034   //initialize GflashShowino for this track
00035   theShowino->initialize(showerType, energy, globalTime, charge, 
00036                          position, momentum, theBField);
00037 
00038 }
00039 
00040 void GflashHadronShowerProfile::hadronicParameterization()
00041 {
00042   // The skeleton of this method is based on the fortran code gfshow.F originally written  
00043   // by S. Peters and G. Grindhammer (also see NIM A290 (1990) 469-488), but longitudinal
00044   // parameterizations of hadron showers are significantly modified for the CMS calorimeter  
00045 
00046   // unit convention: energy in [GeV] and length in [cm]
00047   // intrinsic properties of hadronic showers (lateral shower profile)
00048   
00049   // The step size of showino along the helix trajectory in cm unit
00050   double showerDepthR50 = 0.0;
00051   bool firstHcalHit = true;
00052 
00053   //initial valuses that will be changed as the shower developes
00054   double stepLengthLeft = theShowino->getStepLengthToOut();
00055 
00056   double deltaStep = 0.0;
00057   double showerDepth = 0.0;
00058 
00059   Gflash::CalorimeterNumber whichCalor = Gflash::kNULL;
00060 
00061   theGflashHitList.clear();
00062 
00063   GflashHit aHit;
00064 
00065   while(stepLengthLeft > 0.0) {
00066 
00067     // update shower depth and stepLengthLeft
00068     if ( stepLengthLeft < Gflash::divisionStep ) {
00069       deltaStep = stepLengthLeft;
00070       stepLengthLeft  = 0.0;
00071     }
00072     else {
00073       deltaStep = Gflash::divisionStep;
00074       stepLengthLeft -= deltaStep;
00075     }
00076 
00077     showerDepth += deltaStep;
00078     showerDepthR50 += deltaStep;
00079 
00080     //update GflashShowino
00081     theShowino->updateShowino(deltaStep);
00082 
00083     //evaluate energy in this deltaStep along the longitudinal shower profile
00084     double heightProfile = 0.; 
00085     double deltaEnergy = 0.;
00086 
00087     whichCalor = Gflash::getCalorimeterNumber(theShowino->getPosition());
00088 
00089     //skip if Showino is outside envelopes
00090     if(whichCalor == Gflash::kNULL ) continue;
00091 
00092     heightProfile = longitudinalProfile();
00093 
00094     //skip if the delta energy for this step will be very small
00095     if(heightProfile < 1.00e-08 ) continue;
00096 
00097     //get energy deposition for this step 
00098     deltaEnergy =  heightProfile*Gflash::divisionStep*energyScale[whichCalor];    
00099     theShowino->addEnergyDeposited(deltaEnergy);
00100 
00101     //apply sampling fluctuation if showino is inside the sampling calorimeter
00102     double hadronicFraction = 1.0;
00103     double fluctuatedEnergy = deltaEnergy;
00104 
00105     int nSpotsInStep = std::max(1,static_cast<int>(getNumberOfSpots(whichCalor)*(deltaEnergy/energyScale[whichCalor]) ));
00106     double sampleSpotEnergy = hadronicFraction*fluctuatedEnergy/nSpotsInStep;
00107 
00108     // Lateral shape and fluctuations
00109 
00110     if((whichCalor==Gflash::kHB || whichCalor==Gflash::kHE) && firstHcalHit) {
00111       firstHcalHit = false;
00112       //reset the showerDepth used in the lateral parameterization inside Hcal
00113       showerDepthR50 = Gflash::divisionStep;
00114     }
00115 
00116     //evaluate the fluctuated median of the lateral distribution, R50
00117     double R50 = medianLateralArm(showerDepthR50,whichCalor);
00118 
00119     double hitEnergy = sampleSpotEnergy*CLHEP::GeV;
00120     double hitTime = theShowino->getGlobalTime()*CLHEP::nanosecond;
00121 
00122     Gflash::CalorimeterNumber hitCalor = Gflash::kNULL;
00123 
00124     for (int ispot = 0 ;  ispot < nSpotsInStep ; ispot++) {
00125         
00126       // Compute global position of generated spots with taking into account magnetic field
00127       // Divide deltaStep into nSpotsInStep and give a spot a global position
00128       double incrementPath = theShowino->getPathLength() 
00129         + (deltaStep/nSpotsInStep)*(ispot+0.5 - 0.5*nSpotsInStep);
00130         
00131       // trajectoryPoint give a spot an imaginary point along the shower development
00132       GflashTrajectoryPoint trajectoryPoint;
00133       theShowino->getHelix()->getGflashTrajectoryPoint(trajectoryPoint,incrementPath);
00134         
00135       Gflash3Vector hitPosition = locateHitPosition(trajectoryPoint,R50);   
00136         
00137       hitCalor = Gflash::getCalorimeterNumber(hitPosition);       
00138 
00139       if( hitCalor == Gflash::kNULL) continue;
00140 
00141       hitPosition *= CLHEP::cm;
00142 
00143       aHit.setTime(hitTime);
00144       aHit.setEnergy(hitEnergy);
00145       aHit.setPosition(hitPosition);
00146       theGflashHitList.push_back(aHit);
00147         
00148     } // end of for spot iteration
00149   } // end of while for longitudinal integration
00150 
00151   //HO parameterization
00152 
00153   if(theShowino->getEnergy()> Gflash::MinEnergyCutOffForHO &&
00154      fabs(theShowino->getPositionAtShower().pseudoRapidity()) < Gflash::EtaMax[Gflash::kHO] && 
00155      theGflashHcalOuter) {
00156     
00157     //non zero ho fraction to simulate based on geant4
00158     double nonzeroProb = 0.7*fTanh(theShowino->getEnergy(),Gflash::ho_nonzero);
00159     double r0 = CLHEP::HepUniformRand();
00160     double leftoverE =  theShowino->getEnergy() - theShowino->getEnergyDeposited();
00161 
00162     //@@@ nonzeroProb is not random - need further correlation for non-zero HO energy
00163 
00164     if( r0 < nonzeroProb && leftoverE > 0.0) {
00165             
00166       //starting path Length and stepLengthLeft
00167       double pathLength = theShowino->getHelix()->getPathLengthAtRhoEquals(Gflash::Rmin[Gflash::kHO]-10);
00168       stepLengthLeft = theShowino->getHelix()->getPathLengthAtRhoEquals(Gflash::Rmax[Gflash::kHO]+10) - pathLength;
00169       showerDepth = pathLength - theShowino->getPathLengthAtShower();
00170       
00171       theShowino->setPathLength(pathLength);
00172       
00173       double pathLengthx = theShowino->getHelix()->getPathLengthAtRhoEquals(Gflash::Rmax[Gflash::kHO]);
00174       double pathLengthy = theShowino->getHelix()->getPathLengthAtRhoEquals(Gflash::Rmax[Gflash::kHB]);
00175       
00176       while(stepLengthLeft > 0.0) {
00177         
00178         // update shower depth and stepLengthLeft
00179         if ( stepLengthLeft < Gflash::divisionStep ) {
00180           deltaStep = stepLengthLeft;
00181           stepLengthLeft  = 0.0;
00182         }
00183         else {
00184           deltaStep = Gflash::divisionStep;
00185           stepLengthLeft -= deltaStep;
00186         }
00187         
00188         showerDepth += deltaStep;
00189         
00190         //update GflashShowino
00191         theShowino->updateShowino(deltaStep);
00192         
00193         //evaluate energy in this deltaStep along the longitudinal shower profile
00194         double heightProfile = 0.; 
00195         double deltaEnergy = 0.;
00196         
00197         double hoScale  = leftoverE*(pathLengthx-pathLengthy)/(pathLengthx- theShowino->getPathLengthAtShower());
00198         double refDepth = theShowino->getPathLength() 
00199           - theShowino->getHelix()->getPathLengthAtRhoEquals(Gflash::Rmax[Gflash::kHB]);
00200         
00201         if( refDepth > 0) {
00202           heightProfile = hoProfile(theShowino->getPathLength(),refDepth);
00203           deltaEnergy = heightProfile*Gflash::divisionStep*hoScale;
00204         }
00205         
00206         int nSpotsInStep = std::max(50,static_cast<int>((160.+40* CLHEP::RandGaussQ::shoot())*std::log(theShowino->getEnergy())+50.));
00207         
00208         double hoFraction = 1.00;
00209         double poissonProb = CLHEP::RandPoissonQ::shoot(1.0);
00210         
00211         double fluctuatedEnergy = deltaEnergy*poissonProb;
00212         double sampleSpotEnergy = hoFraction*fluctuatedEnergy/nSpotsInStep;
00213         
00214         // Lateral shape and fluctuations
00215         
00216         //evaluate the fluctuated median of the lateral distribution, R50
00217         double R50 = medianLateralArm(showerDepth,Gflash::kHB);
00218 
00219         double hitEnergy = sampleSpotEnergy*CLHEP::GeV;
00220         double hitTime = theShowino->getGlobalTime()*CLHEP::nanosecond;
00221         
00222         for (int ispot = 0 ;  ispot < nSpotsInStep ; ispot++) {
00223           
00224           double incrementPath = theShowino->getPathLength() 
00225             + (deltaStep/nSpotsInStep)*(ispot+0.5 - 0.5*nSpotsInStep);
00226           
00227           // trajectoryPoint give a spot an imaginary point along the shower development
00228           GflashTrajectoryPoint trajectoryPoint;
00229           theShowino->getHelix()->getGflashTrajectoryPoint(trajectoryPoint,incrementPath);
00230           
00231           Gflash3Vector hitPosition = locateHitPosition(trajectoryPoint,R50);   
00232           hitPosition *= CLHEP::cm;
00233 
00234           if(std::fabs(hitPosition.getZ()/CLHEP::cm) > Gflash::Zmax[Gflash::kHO] ) continue;
00235 
00236           aHit.setTime(hitTime);
00237           aHit.setEnergy(hitEnergy);
00238           aHit.setPosition(hitPosition);
00239           theGflashHitList.push_back(aHit);
00240           
00241         } // end of for HO spot iteration
00242       } // end of while for HO longitudinal integration
00243     }
00244   }
00245 
00246   //  delete theGflashNavigator;
00247 
00248 }
00249 
00250 void GflashHadronShowerProfile::loadParameters()
00251 {
00252   edm::LogInfo("SimGeneralGFlash") << "GflashHadronShowerProfile::loadParameters() "
00253                                    << "should be implimented for each particle type";
00254 }
00255 
00256 double GflashHadronShowerProfile::medianLateralArm(double showerDepthR50, Gflash::CalorimeterNumber kCalor) 
00257 {
00258   double lateralArm = 0.0;
00259   if(kCalor != Gflash::kNULL) {
00260     
00261     double showerDepthR50X = std::min(showerDepthR50/22.4, Gflash::maxShowerDepthforR50);
00262     double R50          = lateralPar[kCalor][0] + std::max(0.0,lateralPar[kCalor][1]) * showerDepthR50X;
00263     double varinanceR50 = std::pow((lateralPar[kCalor][2] + lateralPar[kCalor][3] * showerDepthR50X) * R50, 2);
00264     
00265     // Simulation of lognormal distribution
00266     
00267     if(R50>0) {
00268       double sigmaSq  = std::log(varinanceR50/(R50*R50)+1.0);
00269       double sigmaR50 = std::sqrt(sigmaSq);
00270       double meanR50  = std::log(R50) - (sigmaSq/2.);
00271       
00272       lateralArm = std::exp(meanR50 + sigmaR50* CLHEP::RandGaussQ::shoot());
00273     }
00274   }
00275   return lateralArm;
00276 }
00277 
00278 Gflash3Vector GflashHadronShowerProfile::locateHitPosition(GflashTrajectoryPoint& point, double lateralArm) 
00279 {
00280   // Smearing in r according to f(r)= 2.*r*R50**2/(r**2+R50**2)**2
00281   double rnunif = CLHEP::HepUniformRand();
00282   double rxPDF = std::sqrt(rnunif/(1.-rnunif));
00283   double rShower = lateralArm*rxPDF;
00284   
00285   //rShower within maxLateralArmforR50
00286   rShower = std::min(Gflash::maxLateralArmforR50,rShower);
00287 
00288   // Uniform smearing in phi
00289   double azimuthalAngle = CLHEP::twopi*CLHEP::HepUniformRand(); 
00290 
00291   // actual spot position by adding a radial vector to a trajectoryPoint
00292   Gflash3Vector position = point.getPosition() +
00293     rShower*std::cos(azimuthalAngle)*point.getOrthogonalUnitVector() +
00294     rShower*std::sin(azimuthalAngle)*point.getCrossUnitVector();
00295   
00296   //@@@debugging histograms
00297   if(theHisto->getStoreFlag()) {
00298     theHisto->rshower->Fill(rShower);
00299     theHisto->lateralx->Fill(rShower*std::cos(azimuthalAngle));
00300     theHisto->lateraly->Fill(rShower*std::sin(azimuthalAngle));
00301   }
00302   return position;
00303 }
00304 
00305 //double GflashHadronShowerProfile::longitudinalProfile(double showerDepth, double pathLength) {
00306 double GflashHadronShowerProfile::longitudinalProfile() {
00307 
00308   double heightProfile = 0.0;
00309 
00310   Gflash3Vector pos = theShowino->getPosition();
00311   int showerType = theShowino->getShowerType();
00312   double showerDepth = theShowino->getDepth();
00313   double transDepth = theShowino->getStepLengthToHcal();
00314 
00315   // Energy in a delta step (dz) = (energy to deposite)*[Gamma(z+dz)-Gamma(z)]*dz
00316   // where the incomplete Gamma function gives an intergrate probability of the longitudinal 
00317   // shower up to the shower depth (z).
00318   // Instead, we use approximated energy; energy in dz = (energy to deposite)*gamma(z)*dz
00319   // where gamma is the Gamma-distributed probability function
00320 
00321   Gflash::CalorimeterNumber whichCalor = Gflash::getCalorimeterNumber(pos);
00322 
00323   if(showerType == 0 || showerType == 4 ) {
00324     double shiftDepth = theShowino->getPathLengthOnEcal() - theShowino->getPathLengthAtShower();
00325     if(shiftDepth > 0 ) {
00326       heightProfile = twoGammaProfile(longEcal,showerDepth-shiftDepth,whichCalor);
00327     }
00328     else  {
00329       heightProfile = 0.;
00330       //      std::cout << "negative shiftDepth for showerType 0 " << shiftDepth << std::endl;
00331     }
00332   }  
00333   else if(showerType == 1 || showerType == 5 ) {
00334     if(whichCalor == Gflash::kESPM || whichCalor == Gflash::kENCA ) {
00335       heightProfile = twoGammaProfile(longEcal,showerDepth,whichCalor);
00336     }
00337     else if(whichCalor == Gflash::kHB || whichCalor == Gflash::kHE) {
00338       heightProfile = twoGammaProfile(longHcal,showerDepth-transDepth,whichCalor);
00339     }
00340     else  heightProfile = 0.;
00341   }  
00342   else if (showerType == 2 || showerType == 6 ) {
00343     //two gammas between crystal and Hcal
00344     if((showerDepth - transDepth) > 0.0) {
00345       heightProfile = twoGammaProfile(longHcal,showerDepth-transDepth,Gflash::kHB);
00346     }
00347     else heightProfile = 0.;
00348   }
00349   else if (showerType == 3 || showerType == 7 ) {
00350     //two gammas inside Hcal
00351     heightProfile = twoGammaProfile(longHcal,showerDepth,Gflash::kHB);
00352   }
00353 
00354   return heightProfile;
00355 }
00356 
00357 double GflashHadronShowerProfile::hoProfile(double pathLength, double refDepth) {
00358 
00359   double heightProfile = 0;
00360 
00361   GflashTrajectoryPoint tempPoint;
00362   theShowino->getHelix()->getGflashTrajectoryPoint(tempPoint,pathLength);
00363 
00364   double dint = 1.4*Gflash::intLength[Gflash::kHO]*std::sin(tempPoint.getPosition().getTheta());
00365   heightProfile = std::exp(-1.0*refDepth/dint);
00366 
00367   return heightProfile;
00368 }
00369 
00370 void GflashHadronShowerProfile::getFluctuationVector(double *lowTriangle, double *correlationVector) {
00371 
00372     const int dim = Gflash::NPar;
00373 
00374     double **xr   = new double *[dim];
00375     double **xrho = new double *[dim];
00376     
00377     for(int j=0;j<dim;j++) {
00378       xr[j]   = new double [dim];
00379       xrho[j] = new double [dim];
00380     }
00381     
00382     for(int i = 0; i < dim; i++) {
00383       for(int j = 0; j < i+1 ; j++) {
00384         if(j==i) xrho[i][j] = 1.0;
00385         else {
00386           xrho[i][j] = lowTriangle[i*(i-1)/2 + j];
00387           xrho[j][i] = xrho[i][j];
00388         }
00389       }
00390     }
00391 
00392     doCholeskyReduction(xrho,xr,dim);
00393 
00394     for(int i = 0 ; i < dim ; i++) {
00395       for (int j = 0 ; j < i+1 ; j++){
00396         correlationVector[i*(i+1)/2 + j] = xr[i][j];
00397       }
00398     }
00399 
00400     for(int j=0;j<dim;j++) delete [] xr[j];
00401     delete [] xr;
00402     for(int j=0;j<dim;j++) delete [] xrho[j];
00403     delete [] xrho;
00404 }
00405 
00406 void GflashHadronShowerProfile::doCholeskyReduction(double **vv, double **cc, const int ndim) {
00407 
00408   double sumCjkSquare;
00409   double vjjLess;
00410   double sumCikjk;
00411 
00412   cc[0][0] = std::sqrt(vv[0][0]);
00413 
00414   for(int j=1 ; j < ndim ; j++) {
00415     cc[j][0] = vv[j][0]/cc[0][0];
00416   }
00417 
00418   for(int j=1 ; j < ndim ; j++) {
00419 
00420     sumCjkSquare = 0.0;
00421     for (int k=0 ; k < j ; k++) sumCjkSquare += cc[j][k]*cc[j][k];
00422 
00423     vjjLess =  vv[j][j] - sumCjkSquare;
00424 
00425     //check for the case that vjjLess is negative
00426     cc[j][j] = std::sqrt(std::fabs(vjjLess));
00427 
00428     for (int i=j+1 ; i < ndim ; i++) {
00429       sumCikjk = 0.;
00430       for(int k=0 ; k < j ; k++) sumCikjk += cc[i][k]*cc[j][k];
00431       cc[i][j] = (vv[i][j] - sumCikjk)/cc[j][j];
00432     }
00433   }
00434 }
00435 
00436 int GflashHadronShowerProfile::getNumberOfSpots(Gflash::CalorimeterNumber kCalor) {
00437   //generator number of spots: energy dependent Gamma distribution of Nspots based on Geant4
00438   //replacing old parameterization of H1,
00439   //int numberOfSpots = std::max( 50, static_cast<int>(80.*std::log(einc)+50.));
00440 
00441   double einc = theShowino->getEnergy();
00442   int showerType = theShowino->getShowerType();
00443 
00444   int numberOfSpots = 0;
00445   double nmean  = 0.0;
00446   double nsigma = 0.0;
00447 
00448   if(showerType == 0 || showerType == 1 || showerType == 4 || showerType == 5 ) {
00449     if(kCalor == Gflash::kESPM || kCalor == Gflash::kENCA) {
00450       nmean = 10000 + 5000*log(einc);
00451       nsigma = 1000;
00452     }
00453     if(kCalor == Gflash::kHB || kCalor == Gflash::kHE) {
00454       nmean =  5000 + 2500*log(einc);
00455       nsigma =  500;
00456     }
00457   }
00458   else if (showerType == 2 || showerType == 3 || showerType == 6 || showerType == 7 ) {
00459     if(kCalor == Gflash::kHB || kCalor == Gflash::kHE) {
00460       nmean =  5000 + 2500*log(einc);
00461       nsigma =  500;
00462     }
00463     else {
00464       nmean = 10000;
00465       nsigma = 1000;
00466     }
00467   }
00468   //@@@need correlation and individual fluctuation on alphaNspots and betaNspots here:
00469   //evaluating covariance should be straight forward since the distribution is 'one' Gamma
00470 
00471   numberOfSpots = std::max(500,static_cast<int> (nmean+nsigma* CLHEP::RandGaussQ::shoot()));
00472 
00473   //until we optimize the reduction scale in the number of Nspots
00474       
00475   if( kCalor == Gflash::kESPM ||  kCalor == Gflash::kENCA) {
00476     numberOfSpots = static_cast<int>(numberOfSpots/100);
00477   }
00478   else {
00479     numberOfSpots = static_cast<int>(numberOfSpots/3.0);
00480   }
00481 
00482   return numberOfSpots;
00483 }
00484 
00485 double GflashHadronShowerProfile::fTanh(double einc, const double *par) {
00486   double func = 0.0;
00487   if(einc>0.0) func = par[0]+par[1]*std::tanh(par[2]*(std::log(einc)-par[3])) + par[4]*std::log(einc);
00488   return func;
00489 }
00490 
00491 double GflashHadronShowerProfile::fLnE1(double einc, const double *par) {
00492   double func = 0.0;
00493   if(einc>0.0) func = par[0]+par[1]*std::log(einc);
00494   return func;
00495 }
00496 
00497 double GflashHadronShowerProfile::depthScale(double ssp, double ssp0, double length) {
00498   double func = 0.0;
00499   if(length>0.0) func = std::pow((ssp-ssp0)/length,2.0);
00500   return func;
00501 }
00502 
00503 double GflashHadronShowerProfile::gammaProfile(double alpha, double beta, double showerDepth, double lengthUnit) {
00504   double gamma = 0.0;
00505   //  if(alpha > 0 && beta > 0 && lengthUnit > 0) {
00506   if(showerDepth>0.0) {
00507     Genfun::LogGamma lgam;
00508     double x = showerDepth*(beta/lengthUnit);
00509     gamma = (beta/lengthUnit)*std::pow(x,alpha-1.0)*std::exp(-x)/std::exp(lgam(alpha));
00510   }
00511   return gamma;
00512 }
00513  
00514 double GflashHadronShowerProfile::twoGammaProfile(double *longPar, double depth, Gflash::CalorimeterNumber kIndex) {
00515   double twoGamma = 0.0;
00516 
00517   longPar[0] = std::min(1.0,longPar[0]);
00518   longPar[0] = std::max(0.0,longPar[0]);
00519 
00520   if(longPar[3] > 4.0 || longPar[4] > 4.0) {
00521     double rfactor = 2.0/std::max(longPar[3],longPar[4]);
00522     longPar[3] = rfactor*(longPar[3]+1.0);  
00523     longPar[4] *= rfactor;  
00524   }
00525 
00526   twoGamma  = longPar[0]*gammaProfile(exp(longPar[1]),exp(longPar[2]),depth,Gflash::radLength[kIndex])
00527          +(1-longPar[0])*gammaProfile(exp(longPar[3]),exp(longPar[4]),depth,Gflash::intLength[kIndex]);
00528   return twoGamma;
00529 }