CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/TrackingTools/PatternTools/src/TwoTrackMinimumDistanceHelixLine.cc

Go to the documentation of this file.
00001 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistanceHelixLine.h"
00002 #include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
00003 #include "MagneticField/Engine/interface/MagneticField.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include <iostream>
00006 #include <iomanip>
00007 
00008 using namespace std;
00009 
00010 bool TwoTrackMinimumDistanceHelixLine::updateCoeffs()
00011 {
00012 
00013   if (firstGTP->charge() == 0. && secondGTP->charge() != 0.) {
00014     theL= firstGTP;
00015     theH= secondGTP;
00016   } else if (firstGTP->charge() != 0. && secondGTP->charge() == 0.) {
00017     theH= firstGTP;
00018     theL= secondGTP;
00019   } else {
00020     edm::LogWarning ("TwoTrackMinimumDistanceHelixLine")
00021       << "Error in track charge: "
00022       << "One of the tracks has to be charged, and the other not." << endl
00023       << "Track Charges: "<<firstGTP->charge() << " and " <<secondGTP->charge();
00024     return true;
00025   }
00026 
00027   Hn = theH->momentum().mag();
00028   Ln = theL->momentum().mag();
00029 
00030   if ( Hn == 0. || Ln == 0. )
00031   {
00032     edm::LogWarning ("TwoTrackMinimumDistanceHelixLine")
00033       << "Momentum of input trajectory is zero.";
00034     return true;
00035   };
00036 
00037   GlobalPoint lOrig = theL->position();
00038   GlobalPoint hOrig = theH->position();
00039   posDiff = GlobalVector((lOrig - hOrig).basicVector());
00040   X = posDiff.x();
00041   Y = posDiff.y();
00042   Z = posDiff.z();
00043   theLp = theL->momentum();
00044   px = theLp.x(); px2 = px*px;
00045   py = theLp.y(); py2 = py*py;
00046   pz = theLp.z(); pz2 = pz*pz;
00047   
00048   const double Bc2kH = theH->magneticField().inTesla(hOrig).z() * 2.99792458e-3;
00049 //   MagneticField::inInverseGeV ( hOrig ).z();
00050 
00051   if ( Bc2kH == 0. )
00052   {
00053     edm::LogWarning ("TwoTrackMinimumDistanceHelixLine")
00054       << "Magnetic field at point " << hOrig << " is zero.";
00055     return true;
00056   };
00057 
00058   theh= - Hn / (theH->charge() * Bc2kH ) *
00059      sqrt( 1 - ( ( (theH->momentum().z()*theH->momentum().z()) / (Hn*Hn) )));
00060 
00061   thetanlambdaH = - theH->momentum().z() / ( theH->charge() * Bc2kH * theh);
00062 
00063   thePhiH0 = theH->momentum().phi();
00064   thesinPhiH0= sin(thePhiH0);
00065   thecosPhiH0= cos(thePhiH0);
00066 
00067   aa = (X + theh*thesinPhiH0)*(py2 + pz2) - px*(py*Y + pz*Z);
00068   bb = (Y - theh*thecosPhiH0)*(px2 + pz2) - py*(px*X + pz*Z);
00069   cc = pz*theh*thetanlambdaH;
00070   dd = theh* px *py;
00071   ee = theh*(px2 - py2);
00072   ff = (px2 + py2)*theh*thetanlambdaH*thetanlambdaH;
00073 
00074   baseFct = thetanlambdaH * (Z*(px2+py2) - pz*(px*X + py*Y));
00075   baseDer = - ff;
00076   return false;
00077 }
00078 
00079 bool TwoTrackMinimumDistanceHelixLine::oneIteration(
00080     double & thePhiH, double & fct, double & derivative ) const
00081 {
00082 
00083   double thesinPhiH = sin(thePhiH);
00084   double thecosPhiH = cos(thePhiH);
00085 
00086         // Fonction of which the root is to be found:
00087 
00088   fct = baseFct;
00089   fct -= ff*(thePhiH - thePhiH0);
00090   fct += thecosPhiH * aa;
00091   fct += thesinPhiH * bb;
00092   fct += cc*(thePhiH - thePhiH0)*(px * thecosPhiH + py * thesinPhiH);
00093   fct += cc * (px * (thesinPhiH - thesinPhiH0) - py * (thecosPhiH - thecosPhiH0));
00094   fct += dd * (thesinPhiH* (thesinPhiH - thesinPhiH0) - 
00095                 thecosPhiH*(thecosPhiH - thecosPhiH0));
00096   fct += ee * thecosPhiH * thesinPhiH;
00097 
00098           // Its derivative:
00099 
00100   derivative = baseDer;
00101   derivative += - thesinPhiH * aa;
00102   derivative += thecosPhiH * bb;
00103   derivative += cc*(thePhiH - thePhiH0)*(py * thecosPhiH - px * thesinPhiH);
00104   derivative += 2* cc*(px * thecosPhiH + py * thesinPhiH);
00105   derivative += dd *(4 * thecosPhiH * thesinPhiH - thecosPhiH * thesinPhiH0 - 
00106                         thesinPhiH * thecosPhiH0);
00107   derivative += ee * (thecosPhiH*thecosPhiH-thesinPhiH*thesinPhiH);
00108 
00109   return false;
00110 }
00111 
00112 bool TwoTrackMinimumDistanceHelixLine::calculate(
00113     const GlobalTrajectoryParameters & theFirstGTP,
00114     const GlobalTrajectoryParameters & theSecondGTP, const float qual )
00115 {
00116   pointsUpdated = false;
00117   firstGTP  = (GlobalTrajectoryParameters *) &theFirstGTP;
00118   secondGTP = (GlobalTrajectoryParameters *) &theSecondGTP;
00119 
00120   if ( updateCoeffs () )
00121   {
00122     return true;
00123   };
00124 
00125   double fctVal, derVal, dPhiH;
00126   thePhiH = thePhiH0;
00127   
00128   double x1=thePhiH0-M_PI, x2=thePhiH0+M_PI;
00129   for (int j=1; j<=themaxiter; ++j) { 
00130     oneIteration(thePhiH, fctVal, derVal);
00131     dPhiH=fctVal/derVal;
00132     thePhiH -= dPhiH;
00133     if ((x1-thePhiH)*(thePhiH-x2) < 0.0) {
00134       LogDebug ("TwoTrackMinimumDistanceHelixLine")
00135         << "Jumped out of brackets in root finding. Will be moved closer.";
00136       thePhiH += (dPhiH*0.8);
00137     }
00138     if (fabs(dPhiH) < qual) {return false;}
00139   }
00140   LogDebug ("TwoTrackMinimumDistanceHelixLine")
00141     <<"Number of steps exceeded. Has not converged.";
00142   return true;
00143 }
00144 
00145 double TwoTrackMinimumDistanceHelixLine::firstAngle()  const
00146 {
00147   if (firstGTP==theL) return theL->momentum().phi();
00148   else return thePhiH;
00149 }
00150 
00151 double TwoTrackMinimumDistanceHelixLine::secondAngle() const
00152 {
00153   if (secondGTP==theL) return theL->momentum().phi();
00154   else return thePhiH;
00155 }
00156 
00157 pair <GlobalPoint, GlobalPoint> TwoTrackMinimumDistanceHelixLine::points()
00158     const 
00159 {
00160   if (!pointsUpdated)finalPoints();
00161   if (firstGTP==theL) 
00162     return pair<GlobalPoint, GlobalPoint> (linePoint, helixPoint);
00163   else return pair<GlobalPoint, GlobalPoint> (helixPoint, linePoint);
00164 }
00165 
00166 pair <double, double> TwoTrackMinimumDistanceHelixLine::pathLength() const
00167 {
00168   if (!pointsUpdated)finalPoints();
00169   if (firstGTP==theL) 
00170     return pair<double, double> (linePath, helixPath);
00171   else return pair<double, double> (helixPath, linePath);
00172 }
00173 
00174 void TwoTrackMinimumDistanceHelixLine::finalPoints() const
00175 {
00176   helixPoint = GlobalPoint (
00177       theH->position().x() + theh * ( sin ( thePhiH) - thesinPhiH0 ),
00178       theH->position().y() + theh * ( - cos ( thePhiH) + thecosPhiH0 ),
00179       theH->position().z() + theh * ( thetanlambdaH * ( thePhiH- thePhiH0 ))
00180   );
00181   helixPath = ( thePhiH- thePhiH0 ) * (theh*sqrt(1+thetanlambdaH*thetanlambdaH)) ;
00182 
00183   GlobalVector diff((theL->position() -helixPoint).basicVector());
00184   tL = ( - diff.dot(theLp)) / (Ln*Ln);
00185   linePoint = GlobalPoint (
00186       theL->position().x() + tL * px,
00187       theL->position().y() + tL * py,
00188       theL->position().z() + tL * pz );
00189   linePath = tL * theLp.mag();
00190   pointsUpdated = true;
00191 }