CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/TrackingTools/MaterialEffects/src/MultipleScatteringUpdator.cc

Go to the documentation of this file.
00001 #include "TrackingTools/MaterialEffects/interface/MultipleScatteringUpdator.h"
00002 #include "DataFormats/GeometrySurface/interface/MediumProperties.h"
00003 
00004 
00005 void oldMUcompute (const TrajectoryStateOnSurface& TSoS, 
00006                    const PropagationDirection propDir, double mass, double ptmin);
00007 
00008 //#define DBG_MSU
00009 
00010 //
00011 // Computation of contribution of multiple scatterning to covariance matrix 
00012 //   of local parameters based on Highland formula for sigma(alpha) in plane.
00013 //
00014 void MultipleScatteringUpdator::compute (const TrajectoryStateOnSurface& TSoS, 
00015                                          const PropagationDirection propDir) const
00016 {
00017   //
00018   // Get surface
00019   //
00020   const Surface& surface = TSoS.surface();
00021   //
00022   // Initialise the update to the covariance matrix
00023   // (dP is constantly 0).
00024   //
00025   theDeltaCov(1,1) = 0.;
00026   theDeltaCov(1,2) = 0.;
00027   theDeltaCov(2,2) = 0.;
00028   //
00029   // Now get information on medium
00030   //
00031   if (surface.mediumProperties()) {
00032     // Momentum vector
00033     LocalVector d = TSoS.localMomentum();
00034     double p2 = d.mag2();
00035     d *= 1./sqrt(p2);
00036     // MediumProperties mp(0.02, .5e-4);
00037     const MediumProperties& mp = *surface.mediumProperties();
00038     double xf = 1./fabs(d.z());         // increase of path due to angle of incidence
00039     // calculate general physics things
00040     const double amscon = 1.8496e-4;    // (13.6MeV)**2
00041     const double m2 = mass()*mass();            // use mass hypothesis from constructor
00042     double e2     = p2 + m2;
00043     double beta2  = p2/e2;
00044     // calculate the multiple scattering angle
00045     double radLen = mp.radLen()*xf;     // effective rad. length
00046     double sigt2 = 0.;                  // sigma(alpha)**2
00047     if (radLen > 0) {
00048       // Calculated rms scattering angle squared.
00049       double fact = 1. + 0.038*log(radLen); fact *=fact;
00050       double a = fact/(beta2*p2);
00051       sigt2 = amscon*radLen*a;
00052       if (thePtMin > 0) {
00053 #ifdef DBG_MSU
00054         std::cout<<"Original rms scattering = "<<sqrt(sigt2);
00055 #endif
00056         // Inflate estimated rms scattering angle, to take into account 
00057         // that 1/p is not known precisely.
00058         AlgebraicSymMatrix55 const & covMatrix = TSoS.localError().matrix();
00059         double error2_QoverP = covMatrix(0,0);
00060         // Formula valid for ultra-relativistic particles.
00061         //      sigt2 *= (1. + p2 * error2_QoverP);
00062         // Exact formula
00063         sigt2 *= (1. + p2 * error2_QoverP *
00064                                (1. + 5.*m2/e2 + 3.*m2*beta2*error2_QoverP));
00065 #ifdef DBG_MSU
00066         std::cout<<" new = "<<sqrt(sigt2);
00067 #endif
00068         // Convert Pt constraint to P constraint, neglecting uncertainty in 
00069         // track angle.
00070         double pMin2 = thePtMin*thePtMin*(p2/TSoS.globalMomentum().perp2());       
00071         // Use P constraint to calculate rms maximum scattering angle.
00072         double betaMin2 = pMin2/(pMin2 + m2);
00073         double a_max = fact/(betaMin2 * pMin2);
00074         double sigt2_max = amscon*radLen*a_max;
00075         if (sigt2 > sigt2_max) sigt2 = sigt2_max;
00076 #ifdef DBG_MSU
00077         std::cout<<" after P constraint ("<<pMin<<") = "<<sqrt(sigt2);
00078         std::cout<<" for track with 1/p="<<1/p<<"+-"<<sqrt(error2_QoverP)<<std::endl;
00079 #endif
00080       }
00081     }
00082     double sl2 = d.perp2();
00083     double cl2 = (d.z()*d.z());
00084     double cf2 = (d.x()*d.x())/sl2;
00085     double sf2 = (d.y()*d.y())/sl2;
00086     // Create update (transformation of independant variations
00087     //   on angle in orthogonal planes to local parameters.
00088     double den = 1./(cl2*cl2);
00089     theDeltaCov(1,1) = (den*sigt2)*(sf2*cl2 + cf2);
00090     theDeltaCov(1,2) = (den*sigt2)*(d.x()*d.y()  );
00091     theDeltaCov(2,2) = (den*sigt2)*(cf2*cl2 + sf2);
00092 
00093     /*
00094     std::cout << "new " <<  theDeltaCov(1,1) << " " <<  theDeltaCov(1,2)  << " " <<  theDeltaCov(2,2) << std::endl;
00095     oldMUcompute(TSoS,propDir, mass(), thePtMin);
00096     */
00097   }
00098 
00099 }
00100 
00101 
00102 
00103 
00104 //
00105 // Computation of contribution of multiple scatterning to covariance matrix 
00106 //   of local parameters based on Highland formula for sigma(alpha) in plane.
00107 //
00108 void oldMUcompute (const TrajectoryStateOnSurface& TSoS, 
00109                    const PropagationDirection propDir, double mass, double thePtMin)
00110 {
00111   //
00112   // Get surface
00113   //
00114   const Surface& surface = TSoS.surface();
00115   //
00116   //
00117   // Now get information on medium
00118   //
00119   if (surface.mediumProperties()) {
00120     // Momentum vector
00121     LocalVector d = TSoS.localMomentum();
00122     double p = d.mag();
00123     d *= 1./p;
00124     // MediumProperties mp(0.02, .5e-4);
00125     const MediumProperties& mp = *surface.mediumProperties();
00126     double xf = 1./fabs(d.z());         // increase of path due to angle of incidence
00127     // calculate general physics things
00128     const double amscon = 1.8496e-4;    // (13.6MeV)**2
00129     const double m = mass;            // use mass hypothesis from constructor
00130     double e     = sqrt(p*p + m*m);
00131     double beta  = p/e;
00132     // calculate the multiple scattering angle
00133     double radLen = mp.radLen()*xf;     // effective rad. length
00134     double sigt2 = 0.;                  // sigma(alpha)**2
00135     if (radLen > 0) {
00136       // Calculated rms scattering angle squared.
00137       double a = (1. + 0.038*log(radLen))/(beta*p); a *= a;
00138       sigt2 = amscon*radLen*a;
00139       if (thePtMin > 0) {
00140 #ifdef DBG_MSU
00141         std::cout<<"Original rms scattering = "<<sqrt(sigt2);
00142 #endif
00143         // Inflate estimated rms scattering angle, to take into account 
00144         // that 1/p is not known precisely.
00145         AlgebraicSymMatrix55 covMatrix = TSoS.localError().matrix();
00146         double error2_QoverP = covMatrix(0,0);
00147         // Formula valid for ultra-relativistic particles.
00148 //      sigt2 *= (1. + (p*p) * error2_QoverP);
00149         // Exact formula
00150         sigt2 *= (1. + (p*p) * error2_QoverP *
00151                                (1. + 5*m*m/(e*e) + 3*m*m*beta*beta*error2_QoverP));
00152 #ifdef DBG_MSU
00153         std::cout<<" new = "<<sqrt(sigt2);
00154 #endif
00155         // Convert Pt constraint to P constraint, neglecting uncertainty in 
00156         // track angle.
00157         double pMin = thePtMin*(TSoS.globalMomentum().mag()/TSoS.globalMomentum().perp());       
00158         // Use P constraint to calculate rms maximum scattering angle.
00159         double betaMin = pMin/sqrt(pMin * pMin + m*m);
00160         double a_max = (1. + 0.038*log(radLen))/(betaMin * pMin); a_max *= a_max;
00161         double sigt2_max = amscon*radLen*a_max;
00162         if (sigt2 > sigt2_max) sigt2 = sigt2_max;
00163 #ifdef DBG_MSU
00164         std::cout<<" after P constraint ("<<pMin<<") = "<<sqrt(sigt2);
00165         std::cout<<" for track with 1/p="<<1/p<<"+-"<<sqrt(error2_QoverP)<<std::endl;
00166 #endif
00167       }
00168     }
00169     double sl = d.perp();
00170     double cl = d.z();
00171     double cf = d.x()/sl;
00172     double sf = d.y()/sl;
00173     // Create update (transformation of independant variations
00174     //   on angle in orthogonal planes to local parameters.
00175     std::cout << "old " << sigt2*(sf*sf*cl*cl + cf*cf)/(cl*cl*cl*cl)
00176               << " " << sigt2*(cf*sf*sl*sl        )/(cl*cl*cl*cl) 
00177               << " " << sigt2*(cf*cf*cl*cl + sf*sf)/(cl*cl*cl*cl) << std::endl;
00178   }
00179 }