CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/TrackingTools/PatternTools/src/TwoTrackMinimumDistanceHelixHelix.cc

Go to the documentation of this file.
00001 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistanceHelixHelix.h"
00002 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00003 #include "MagneticField/Engine/interface/MagneticField.h"
00004 #include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 
00007 using namespace std;
00008 
00009 namespace {
00010   inline GlobalPoint operator - ( const GlobalPoint & a, const GlobalPoint & b ){
00011     return GlobalPoint ( a.basicVector() - b.basicVector() );
00012   }
00013 
00014   inline GlobalPoint operator + ( const GlobalPoint & a, const GlobalPoint & b ){
00015     return GlobalPoint ( a.basicVector() + b.basicVector() );
00016   }
00017 
00018   inline GlobalPoint operator / ( const GlobalPoint & a, const double b ){
00019     return GlobalPoint ( a.basicVector() / b );
00020   }
00021 
00022   inline GlobalPoint operator * ( const GlobalPoint & a, const float b ){
00023     return GlobalPoint ( a.basicVector() * b );
00024   }
00025 
00026   inline GlobalPoint operator * ( const float b , const GlobalPoint & a ){
00027     return GlobalPoint ( a.basicVector() * b );
00028   }
00029 
00030   inline double square ( const double s ) { return s*s; }
00031 }
00032 
00033 TwoTrackMinimumDistanceHelixHelix::TwoTrackMinimumDistanceHelixHelix():
00034 theH(0), theG(0), pointsUpdated(false), themaxjump(20),thesingjacI(1./0.1), themaxiter(4)
00035 { }
00036 
00037 TwoTrackMinimumDistanceHelixHelix::~TwoTrackMinimumDistanceHelixHelix() {}
00038 
00039 bool TwoTrackMinimumDistanceHelixHelix::updateCoeffs(
00040     const GlobalPoint & gpH, const GlobalPoint & gpG ) {
00041 
00042   const double Bc2kH = theH->magneticField().inInverseGeV(gpH).z();
00043   const double Bc2kG = theG->magneticField().inInverseGeV(gpG).z();
00044   const double Ht= theH->momentum().perp();
00045   const double Gt= theG->momentum().perp();
00046   //  thelambdaH=asin ( theH->momentum().z() / Hn );
00047   
00048   if ( Ht == 0. || Gt == 0. ) {
00049     edm::LogWarning ("TwoTrackMinimumDistanceHelixHelix")
00050       << "transverse momentum of input trajectory is zero.";
00051     return true;
00052   };
00053   
00054   if ( theH->charge() == 0. || theG->charge() == 0. ) {
00055     edm::LogWarning ("TwoTrackMinimumDistanceHelixHelix")
00056       << "charge of input track is zero.";
00057     return true;
00058   };
00059   
00060   if ( Bc2kG == 0.) {
00061     edm::LogWarning ("TwoTrackMinimumDistanceHelixHelix")
00062       << "magnetic field at point " << gpG << " is zero.";
00063     return true;
00064   };
00065   
00066   if ( Bc2kH == 0. ) {
00067     edm::LogWarning ("TwoTrackMinimumDistanceHelixHelix")
00068       << "magnetic field at point " << gpH << " is zero.";
00069     return true;
00070   };
00071   
00072   theh= Ht / (theH->charge() * Bc2kH );
00073   thesinpH0= - theH->momentum().y() / ( Ht );
00074   thecospH0= - theH->momentum().x() / ( Ht );
00075   thetanlambdaH = - theH->momentum().z() / ( Ht);
00076   thepH0 = atan2 ( thesinpH0 , thecospH0 );
00077   
00078   // thelambdaG=asin ( theG->momentum().z()/( Gn ) );
00079   
00080   theg= Gt / (theG->charge() * Bc2kG );
00081   thesinpG0= - theG->momentum().y() / ( Gt);
00082   thecospG0= - theG->momentum().x() / (Gt);
00083   thetanlambdaG = - theG->momentum().z() / ( Gt);
00084   thepG0 = atan2 ( thesinpG0 , thecospG0 );
00085   
00086   thea=theH->position().x() - theG->position().x() + theg * thesinpG0 -
00087     theh * thesinpH0;
00088   theb=theH->position().y() - theG->position().y() - theg * thecospG0 +
00089     theh * thecospH0;
00090   thec1=  theh * thetanlambdaH * thetanlambdaH;
00091   thec2= -theg * thetanlambdaG * thetanlambdaG;
00092   thed1= -theg * thetanlambdaG * thetanlambdaH;
00093   thed2=  theh * thetanlambdaG * thetanlambdaH;
00094   thee1= thetanlambdaH * ( theH->position().z() - theG->position().z() -
00095                            theh * thepH0 * thetanlambdaH + theg * thepG0 * thetanlambdaG );
00096   thee2= thetanlambdaG * ( theH->position().z() - theG->position().z() -
00097                            theh * thepH0 * thetanlambdaH + theg * thepG0 * thetanlambdaG );
00098   return false;
00099 }
00100 
00101 bool TwoTrackMinimumDistanceHelixHelix::oneIteration(double & dH, double & dG ) const {
00102   thesinpH=sin(thepH);
00103   thecospH=cos(thepH);
00104   thesinpG=sin(thepG);
00105   thecospG=cos(thepG);
00106   
00107   const double A11= theh * ( - thesinpH * ( thea - theg * thesinpG ) +
00108                              thecospH * ( theb + theg * thecospG ) + thec1);
00109   if (A11 < 0) { return true; };
00110   const double A22= -theg * (- thesinpG * ( thea + theh * thesinpH ) +
00111                              thecospG*( theb - theh * thecospH ) + thec2);
00112   if (A22 < 0) { return true; };
00113   const double A12= theh * (-theg * thecospG * thecospH -
00114                             theg * thesinpH * thesinpG + thed1);
00115   const double A21= -theg * (theh * thecospG * thecospH +
00116                              theh * thesinpH * thesinpG + thed2);
00117   const double detaI = 1./(A11 * A22 - A12 * A21);
00118   const double z1=theh * ( thecospH * ( thea - theg * thesinpG ) + thesinpH *
00119                            ( theb + theg*thecospG ) + thec1 * thepH + thed1 * thepG + thee1);
00120   const double z2=-theg * (thecospG * ( thea + theh * thesinpH ) + thesinpG *
00121                            ( theb - theh*thecospH ) + thec2 * thepG + thed2 * thepH + thee2);
00122   
00123   dH=( z1 * A22 - z2 * A12 ) * detaI;
00124   dG=( z2 * A11 - z1 * A21 ) * detaI;
00125   if ( fabs(detaI) > thesingjacI ) { return true; };
00126   
00127   thepH -= dH;
00128   thepG -= dG;
00129   return false;
00130 }
00131 
00132 /*
00133 bool TwoTrackMinimumDistanceHelixHelix::parallelTracks() const {
00134   return (fabs(theH->momentum().x() - theG->momentum().x()) < .00000001 )
00135     && (fabs(theH->momentum().y() - theG->momentum().y()) < .00000001 )
00136     && (fabs(theH->momentum().z() - theG->momentum().z()) < .00000001 )
00137     && (theH->charge()==theG->charge()) 
00138     ;
00139 }
00140 */
00141 
00142 bool TwoTrackMinimumDistanceHelixHelix::calculate(
00143                                                   const GlobalTrajectoryParameters & G,
00144                                                   const GlobalTrajectoryParameters & H, const float qual ) {
00145   pointsUpdated = false;
00146   theG= &G;
00147   theH= &H;
00148   bool retval=false;
00149   
00150   if ( updateCoeffs ( theG->position(), theH->position() ) ) return true;
00151   
00152   thepG = thepG0;
00153   thepH = thepH0;
00154   
00155   int counter=0;
00156   double pH=0; double pG=0;
00157   do {
00158     retval=oneIteration ( pG, pH );
00159     if ( std::isinf(pG) || std::isinf(pH) ) retval=true;
00160     if ( counter++>themaxiter ) retval=true;
00161   } while ( (!retval) && ( fabs(pG) > qual || fabs(pH) > qual ));
00162   if ( fabs ( theg * ( thepG - thepG0 ) ) > themaxjump ) retval=true;
00163   if ( fabs ( theh * ( thepH - thepH0 ) ) > themaxjump ) retval=true;
00164   return retval;
00165 }
00166 
00167 
00168 void TwoTrackMinimumDistanceHelixHelix::finalPoints() const {
00169   GlobalVector tmpG( sin(thepG) - thesinpG0,
00170                    - cos(thepG) + thecospG0,
00171                    thetanlambdaG * ( thepG- thepG0 ) 
00172                    );
00173   pointG = theG->position() + theg * tmpG;
00174   pathG = ( thepG- thepG0 ) * (theg*sqrt(1+thetanlambdaG*thetanlambdaG)) ;
00175 
00176   GlobalVector tmpH( sin(thepH) - thesinpH0,
00177                    - cos(thepH) + thecospH0,
00178                    thetanlambdaH * ( thepH- thepH0 ) 
00179                    );
00180   pointH = theH->position() + theh * tmpH;
00181   pathH = ( thepH- thepH0 ) * (theh*sqrt(1+thetanlambdaH*thetanlambdaH)) ;
00182 
00183   pointsUpdated = true;
00184 }
00185