CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/TrackingTools/PatternTools/src/TSCPBuilderNoMaterial.cc

Go to the documentation of this file.
00001 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
00002 #include "DataFormats/GeometrySurface/interface/Surface.h" 
00003 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h" 
00004 #include "TrackingTools/AnalyticalJacobians/interface/AnalyticalCurvilinearJacobian.h"
00005 #include "DataFormats/GeometrySurface/interface/BoundPlane.h"
00006 #include "TrackingTools/GeomPropagators/interface/HelixBarrelPlaneCrossingByCircle.h"
00007 #include "TrackingTools/TrajectoryParametrization/interface/TrajectoryStateExceptions.h"
00008 #include "MagneticField/Engine/interface/MagneticField.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 
00011 TrajectoryStateClosestToPoint 
00012 TSCPBuilderNoMaterial::operator() (const FTS& originalFTS, 
00013     const GlobalPoint& referencePoint) const
00014 {
00015   if (positionEqual(referencePoint, originalFTS.position()))
00016     return constructTSCP(originalFTS, referencePoint);
00017 
00018   // Now do the propagation or whatever...
00019 
00020   PairBoolFTS newStatePair =
00021     createFTSatTransverseImpactPoint(originalFTS, referencePoint);
00022   if (newStatePair.first) {
00023     return constructTSCP(newStatePair.second, referencePoint);
00024   } else {
00025     return TrajectoryStateClosestToPoint();
00026   }  
00027 }
00028 
00029 TrajectoryStateClosestToPoint 
00030 TSCPBuilderNoMaterial::operator() (const TSOS& originalTSOS, 
00031     const GlobalPoint& referencePoint) const
00032 {
00033   if (positionEqual(referencePoint, originalTSOS.globalPosition()))
00034     return constructTSCP(*originalTSOS.freeState(), referencePoint);
00035 
00036   // Now do the propagation
00037   
00038   PairBoolFTS newStatePair =
00039     createFTSatTransverseImpactPoint(*originalTSOS.freeState(), referencePoint);
00040   if (newStatePair.first) {
00041     return constructTSCP(newStatePair.second, referencePoint);
00042   } else {
00043     return TrajectoryStateClosestToPoint();
00044   }
00045 }
00046 
00047 TSCPBuilderNoMaterial::PairBoolFTS 
00048 TSCPBuilderNoMaterial::createFTSatTransverseImpactPoint(
00049         const FTS& originalFTS, const GlobalPoint& referencePoint) const 
00050 {
00051   //
00052   // Straight line approximation? |rho|<1.e-10 equivalent to ~ 1um 
00053   // difference in transversal position at 10m.
00054   //
00055   if( fabs(originalFTS.transverseCurvature())<1.e-10 ) {
00056     return createFTSatTransverseImpactPointNeutral(originalFTS, referencePoint);
00057   } else {
00058     return createFTSatTransverseImpactPointCharged(originalFTS, referencePoint);
00059   }
00060 }
00061 
00062 TSCPBuilderNoMaterial::PairBoolFTS 
00063 TSCPBuilderNoMaterial::createFTSatTransverseImpactPointCharged(
00064         const FTS& originalFTS, const GlobalPoint& referencePoint) const 
00065 {
00066 
00067   GlobalVector pvecOrig = originalFTS.momentum();
00068   GlobalPoint xvecOrig = originalFTS.position();
00069   double kappa = originalFTS.transverseCurvature();
00070   double pxOrig = pvecOrig.x();
00071   double pyOrig = pvecOrig.y();
00072   double pzOrig = pvecOrig.z();
00073   double xOrig = xvecOrig.x();
00074   double yOrig = xvecOrig.y();
00075   double zOrig = xvecOrig.z();
00076 
00077 //  double fac = 1./originalFTS.charge()/MagneticField::inInverseGeV(referencePoint).z();
00078   double fac = 1./originalFTS.charge()/
00079     (originalFTS.parameters().magneticField().inInverseGeV(referencePoint).z());
00080   GlobalVectorDouble xOrig2Centre = GlobalVectorDouble(fac * pyOrig, -fac * pxOrig, 0.);
00081   GlobalVectorDouble xOrigProj = GlobalVectorDouble(xOrig, yOrig, 0.);
00082   GlobalVectorDouble xRefProj = GlobalVectorDouble(referencePoint.x(), referencePoint.y(), 0.);
00083   GlobalVectorDouble deltax = xRefProj-xOrigProj-xOrig2Centre;
00084   GlobalVectorDouble ndeltax = deltax.unit();
00085 
00086   PropagationDirection direction = anyDirection;
00087   Surface::PositionType pos(referencePoint);
00088   // Need to define plane with orientation as the
00089   // ImpactPointSurface
00090   GlobalVector X(ndeltax.x(), ndeltax.y(), ndeltax.z());
00091   GlobalVector Y(0.,0.,1.);
00092   Surface::RotationType rot(X,Y);
00093   BoundPlane* plane = new BoundPlane(pos,rot);
00094   // Using Teddy's HelixBarrelPlaneCrossingByCircle for general barrel planes. 
00095   // A large variety of other,
00096   // direct solutions turned out to be not so stable numerically.
00097   HelixBarrelPlaneCrossingByCircle 
00098     planeCrossing(HelixPlaneCrossing::PositionType(xOrig, yOrig, zOrig),
00099                   HelixPlaneCrossing::DirectionType(pxOrig, pyOrig, pzOrig), 
00100                   kappa, direction);
00101   std::pair<bool,double> propResult = planeCrossing.pathLength(*plane);
00102   if ( !propResult.first ) {
00103     edm::LogWarning ("TSCPBuilderNoMaterial") << "Propagation to perigee plane failed!";
00104     return PairBoolFTS(false, FreeTrajectoryState() );
00105   }
00106   double s = propResult.second;
00107   HelixPlaneCrossing::PositionType xGen = planeCrossing.position(s);
00108   GlobalPoint xPerigee = GlobalPoint(xGen.x(),xGen.y(),xGen.z());
00109   // direction (reconverted to GlobalVector, renormalised)
00110   HelixPlaneCrossing::DirectionType pGen = planeCrossing.direction(s);
00111   pGen *= pvecOrig.mag()/pGen.mag();
00112   GlobalVector pPerigee = GlobalVector(pGen.x(),pGen.y(),pGen.z());
00113   delete plane;
00114                   
00115   if (originalFTS.hasError()) {
00116     const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix();
00117     AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee,
00118                                                    pPerigee, s);
00119     const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian();
00120     CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) );
00121   
00122     return PairBoolFTS(true,
00123         FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(), 
00124                         &(originalFTS.parameters().magneticField())), cte) );
00125   } 
00126   else {
00127     return PairBoolFTS(true,
00128         FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(), 
00129                         &(originalFTS.parameters().magneticField()))) );
00130   }
00131   
00132 }
00133 
00134 
00135 TSCPBuilderNoMaterial::PairBoolFTS 
00136 TSCPBuilderNoMaterial::createFTSatTransverseImpactPointNeutral(const FTS& originalFTS, 
00137     const GlobalPoint& referencePoint) const 
00138 {
00139 
00140   GlobalPoint xvecOrig = originalFTS.position();
00141   double phi = originalFTS.momentum().phi();
00142   double theta = originalFTS.momentum().theta();
00143   double xOrig = xvecOrig.x();
00144   double yOrig = xvecOrig.y();
00145   double zOrig = xvecOrig.z();
00146   double xR = referencePoint.x();
00147   double yR = referencePoint.y();
00148 
00149   double s2D = (xR - xOrig)  * cos(phi) + (yR - yOrig)  * sin(phi);
00150   double s = s2D / sin(theta);
00151   double xGen = xOrig + s2D*cos(phi);
00152   double yGen = yOrig + s2D*sin(phi);
00153   double zGen = zOrig + s2D/tan(theta);
00154   GlobalPoint xPerigee = GlobalPoint(xGen, yGen, zGen);
00155 
00156   GlobalVector pPerigee = originalFTS.momentum();
00157                   
00158   if (originalFTS.hasError()) {
00159     const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix();
00160     AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee,
00161                                                    pPerigee, s);
00162     const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian();
00163     CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) );
00164   
00165     return PairBoolFTS(true,
00166         FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(), 
00167                         &(originalFTS.parameters().magneticField())), cte));
00168   } 
00169   else {
00170     return PairBoolFTS(true,
00171         FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
00172                         &(originalFTS.parameters().magneticField()))) );
00173   }
00174   
00175 }