CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/src/TrackingTools/MaterialEffects/src/VolumeEnergyLossEstimator.cc

Go to the documentation of this file.
00001 #include "TrackingTools/MaterialEffects/interface/VolumeEnergyLossEstimator.h"
00002 
00003 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00004 #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
00005 #include "TrackingTools/MaterialEffects/interface/VolumeMaterialEffectsEstimator.h"
00006 #include "TrackingTools/MaterialEffects/interface/VolumeMaterialEffectsEstimate.h"
00007 #include "TrackingTools/MaterialEffects/interface/VolumeMediumProperties.h"
00008 
00009 VolumeMaterialEffectsEstimate 
00010 VolumeEnergyLossEstimator::estimate (const TrajectoryStateOnSurface refTSOS,
00011                                      double pathLength,
00012                                      const VolumeMediumProperties& medium) const
00013 {
00014   //
00015   // Initialise dP and the update to the covariance matrix
00016   //
00017   double deltaP = 0.;
00018   double deltaCov00 = 0.;
00019   //
00020   // Bethe-Bloch
00021   //
00022   if ( mass()>0.001 )
00023     computeBetheBloch(refTSOS.localMomentum(),pathLength,medium,
00024                       deltaP,deltaCov00);
00025   //
00026   // Special treatment for electrons (currently rather crude
00027   // distinction using mass)
00028   //
00029   else
00030     computeElectrons(refTSOS.localMomentum(),pathLength,medium,
00031                      deltaP,deltaCov00);
00032 
00033   AlgebraicSymMatrix55 deltaCov;
00034   deltaCov(0,0) = deltaCov00;
00035   return VolumeMaterialEffectsEstimate(deltaP,deltaCov);
00036 }
00037 
00038 
00039 VolumeEnergyLossEstimator*
00040 VolumeEnergyLossEstimator::clone () const
00041 {
00042   return new VolumeEnergyLossEstimator(*this);
00043 }
00044 
00045 //
00046 // Computation of energy loss according to Bethe-Bloch
00047 //
00048 void
00049 VolumeEnergyLossEstimator::computeBetheBloch (const LocalVector& localP, double pathLength,
00050                                             const VolumeMediumProperties& medium,
00051                                             double& deltaP, double& deltaCov00) const {
00052   //
00053   // calculate absolute momentum and correction to path length from angle 
00054   // of incidence
00055   //
00056   double p = localP.mag();
00057   // constants
00058   const double m = mass();            // use mass hypothesis from constructor
00059   const double emass = 0.511e-3;
00060   // FIXME: replace constants for Si
00061   const double poti = 16.e-9 * 10.75; // = 16 eV * Z**0.9, for Si Z=14
00062   const double eplasma = 28.816e-9 * sqrt(2.33*0.498); // 28.816 eV * sqrt(rho*(Z/A)) for Si
00063   const double delta0 = 2*log(eplasma/poti) - 1.;
00064   // calculate general physics things
00065   double e     = sqrt(p*p + m*m);
00066   double beta  = p/e;
00067   double gamma = e/m;
00068   double eta2  = beta*gamma; eta2 *= eta2;
00069   double lnEta2 = log(eta2);
00070   double ratio = emass/m;
00071   double emax  = 2.*emass*eta2/(1. + 2.*ratio*gamma + ratio*ratio);
00072   double delta = delta0 + lnEta2;
00073   // calculate the mean and sigma of energy loss
00074   // xi = d[g/cm2] * 0.307075MeV/(g/cm2) * Z/A * 1/2
00075   double xi = pathLength * medium.xi(); xi /= (beta*beta);
00076 //   double dEdx = xi*(log(2.*emass*eta2*emax/(poti*poti)) - 2.*(beta*beta));
00077   double dEdx = xi*(log(2.*emass*emax/(poti*poti))+lnEta2 - 2.*(beta*beta) - delta);
00078   double dEdx2 = xi*emax*(1.-beta*beta/2.);
00079   double dP    = dEdx/beta;
00080   double sigp2 = dEdx2*e*e/(p*p*p*p*p*p);
00081   deltaP = -dP;
00082   deltaCov00 = sigp2;
00083 }
00084 //
00085 // Computation of energy loss for electrons
00086 //
00087 void 
00088 VolumeEnergyLossEstimator::computeElectrons (const LocalVector& localP, double pathLength,
00089                                              const VolumeMediumProperties& medium,
00090                                              double& deltaP, double& deltaCov00) const {
00091   //
00092   // calculate absolute momentum and correction to path length from angle 
00093   // of incidence
00094   //
00095   double p = localP.mag();
00096   double normalisedPath = pathLength / medium.x0();
00097   //
00098   // Energy loss and variance according to Bethe and Heitler, see also
00099   // Comp. Phys. Comm. 79 (1994) 157. 
00100   //
00101   double z = exp(-normalisedPath);
00102   double varz = (exp(-normalisedPath*log(3.)/log(2.))-
00103                  exp(-2*normalisedPath));
00104   // FIXME: need to know propagation direction at this point
00105 //   if ( propDir==oppositeToMomentum ) {
00106 //     //
00107 //     // for backward propagation: delta(1/p) is linear in z=p_outside/p_inside
00108 //     // convert to obtain equivalent delta(p). Sign of deltaP is corrected
00109 //     // in method compute -> deltaP<0 at this place!!!
00110 //     //
00111 //     deltaP = -p*(1/z-1);
00112 //     deltaCov00 = varz/p/p;
00113 //   }
00114 //   else {     
00115   //
00116   // for forward propagation: calculate in p (linear in 1/z=p_inside/p_outside),
00117   // then convert sig(p) to sig(1/p). 
00118   //
00119   deltaP = p*(z-1);
00120   //    double f = 1/p/z/z;
00121   // patch to ensure consistency between for- and backward propagation
00122   double f = 1./p/z;
00123   deltaCov00 += f*f*varz;
00124 }