CMS 3D CMS Logo

TwoTrackMinimumDistanceLineLine.cc

Go to the documentation of this file.
00001 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistanceLineLine.h"
00002 #include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
00003 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "MagneticField/Engine/interface/MagneticField.h"
00006 
00007 
00008 bool TwoTrackMinimumDistanceLineLine::calculate(
00009     const GlobalTrajectoryParameters & theG,
00010     const GlobalTrajectoryParameters & theH)
00011 {
00012   GlobalPoint gOrig = theG.position();
00013   GlobalPoint hOrig = theH.position();
00014   if ( ( ( theH.charge() != 0. || theG.charge() != 0. ) ) && 
00015     ((theG.magneticField().inTesla(gOrig).z() != 0.)|| 
00016         (theH.magneticField().inTesla(hOrig).z() != 0.)) )
00017   {
00018     edm::LogWarning ("TwoTrackMinimumDistanceLineLine")
00019       << "charge of input track is not zero or field non zero";
00020     return true;
00021   };
00022 
00023   GlobalVector gVec = theG.momentum();
00024   const double gMag= theG.momentum().mag(); double gMag2 = gMag*gMag;
00025   GlobalVector hVec = theH.momentum();
00026   const double hMag= theH.momentum().mag(); double hMag2 = hMag*hMag;
00027 
00028   if ( gMag == 0. || hMag == 0. )
00029   {
00030     edm::LogWarning ("TwoTrackMinimumDistanceLineLine")
00031       << "momentum of input trajectory is zero.";
00032     return true;
00033   };
00034 
00035   phiG = theG.momentum().phi();
00036   phiH = theH.momentum().phi();
00037 
00038   double gVec_Dot_hVec = gVec.dot(hVec);
00039   double norm  = gVec_Dot_hVec*gVec_Dot_hVec - gMag2*hMag2;
00040 
00041   if ( norm == 0 )
00042   {
00043     edm::LogWarning ("TwoTrackMinimumDistanceLineLine")
00044       << "Tracks are parallel.";
00045     return true;
00046   }
00047 
00048   GlobalVector posDiff = gOrig - hOrig;
00049 
00050   double tG = (posDiff.dot(gVec) * hMag2 - gVec_Dot_hVec * posDiff.dot(hVec))/ norm;
00051   double tH = (gVec_Dot_hVec * posDiff.dot(gVec) - posDiff.dot(hVec) * gMag2)/ norm;
00052 
00053   gPos = GlobalPoint ( gOrig + tG * gVec);
00054   hPos = GlobalPoint ( hOrig + tH * hVec);
00055 //   cout << "TwoTrackMinimumDistanceLineLine "<<gPos<<hPos<<endl;
00056 
00057   pathG = tG * gMag;
00058   pathH = tH * hMag;
00059   return false;
00060 }
00061 
00062 std::pair<GlobalPoint, GlobalPoint> TwoTrackMinimumDistanceLineLine::points() const
00063 {
00064   return std::pair<GlobalPoint, GlobalPoint > ( gPos, hPos );
00065 }
00066 
00067 std::pair<double, double> TwoTrackMinimumDistanceLineLine::pathLength() const
00068 {
00069   return std::pair<double, double> ( pathG, pathH);
00070 }

Generated on Tue Jun 9 17:48:27 2009 for CMSSW by  doxygen 1.5.4