00001 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
00002 #include "DataFormats/GeometrySurface/interface/Surface.h"
00003 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00004 #include "DataFormats/CLHEP/interface/Migration.h"
00005 #include "TrackingTools/AnalyticalJacobians/interface/AnalyticalCurvilinearJacobian.h"
00006 #include "DataFormats/GeometrySurface/interface/BoundPlane.h"
00007 #include "TrackingTools/GeomPropagators/interface/HelixBarrelPlaneCrossingByCircle.h"
00008 #include "TrackingTools/TrajectoryParametrization/interface/TrajectoryStateExceptions.h"
00009 #include "MagneticField/Engine/interface/MagneticField.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
00019
00020 FreeTrajectoryState
00021 ftsAtTransverseImpactPoint = createFTSatTransverseImpactPoint(originalFTS, referencePoint);
00022 return constructTSCP(ftsAtTransverseImpactPoint, referencePoint);
00023
00024 }
00025
00026 TrajectoryStateClosestToPoint
00027 TSCPBuilderNoMaterial::operator() (const TSOS& originalTSOS,
00028 const GlobalPoint& referencePoint) const
00029 {
00030 if (positionEqual(referencePoint, originalTSOS.globalPosition()))
00031 return constructTSCP(*originalTSOS.freeState(), referencePoint);
00032
00033
00034
00035 FreeTrajectoryState ftsAtTransverseImpactPoint =
00036 createFTSatTransverseImpactPoint(*originalTSOS.freeState(), referencePoint);
00037 return constructTSCP(ftsAtTransverseImpactPoint, referencePoint);
00038 }
00039
00040 FreeTrajectoryState
00041 TSCPBuilderNoMaterial::createFTSatTransverseImpactPoint(
00042 const FTS& originalFTS, const GlobalPoint& referencePoint) const
00043 {
00044
00045
00046
00047
00048 if( fabs(originalFTS.transverseCurvature())<1.e-10 ) {
00049 return createFTSatTransverseImpactPointNeutral(originalFTS, referencePoint);
00050 } else {
00051 return createFTSatTransverseImpactPointCharged(originalFTS, referencePoint);
00052 }
00053 }
00054
00055 FreeTrajectoryState
00056 TSCPBuilderNoMaterial::createFTSatTransverseImpactPointCharged(
00057 const FTS& originalFTS, const GlobalPoint& referencePoint) const
00058 {
00059
00060 GlobalVector pvecOrig = originalFTS.momentum();
00061 GlobalPoint xvecOrig = originalFTS.position();
00062 double kappa = originalFTS.transverseCurvature();
00063 double pxOrig = pvecOrig.x();
00064 double pyOrig = pvecOrig.y();
00065 double pzOrig = pvecOrig.z();
00066 double xOrig = xvecOrig.x();
00067 double yOrig = xvecOrig.y();
00068 double zOrig = xvecOrig.z();
00069
00070
00071 double fac = 1./originalFTS.charge()/
00072 (originalFTS.parameters().magneticField().inInverseGeV(referencePoint).z());
00073 GlobalVectorDouble xOrig2Centre = GlobalVectorDouble(fac * pyOrig, -fac * pxOrig, 0.);
00074 GlobalVectorDouble xOrigProj = GlobalVectorDouble(xOrig, yOrig, 0.);
00075 GlobalVectorDouble xRefProj = GlobalVectorDouble(referencePoint.x(), referencePoint.y(), 0.);
00076 GlobalVectorDouble deltax = xRefProj-xOrigProj-xOrig2Centre;
00077 GlobalVectorDouble ndeltax = deltax.unit();
00078
00079 PropagationDirection direction = anyDirection;
00080 Surface::PositionType pos(referencePoint);
00081
00082
00083 GlobalVector X(ndeltax.x(), ndeltax.y(), ndeltax.z());
00084 GlobalVector Y(0.,0.,1.);
00085 Surface::RotationType rot(X,Y);
00086 BoundPlane* plane = new BoundPlane(pos,rot);
00087
00088
00089
00090 HelixBarrelPlaneCrossingByCircle
00091 planeCrossing(HelixPlaneCrossing::PositionType(xOrig, yOrig, zOrig),
00092 HelixPlaneCrossing::DirectionType(pxOrig, pyOrig, pzOrig),
00093 kappa, direction);
00094 std::pair<bool,double> propResult = planeCrossing.pathLength(*plane);
00095 if ( !propResult.first )
00096 throw TrajectoryStateException("Propagation to perigee plane failed!");
00097 double s = propResult.second;
00098 HelixPlaneCrossing::PositionType xGen = planeCrossing.position(s);
00099 GlobalPoint xPerigee = GlobalPoint(xGen.x(),xGen.y(),xGen.z());
00100
00101 HelixPlaneCrossing::DirectionType pGen = planeCrossing.direction(s);
00102 pGen *= pvecOrig.mag()/pGen.mag();
00103 GlobalVector pPerigee = GlobalVector(pGen.x(),pGen.y(),pGen.z());
00104 delete plane;
00105
00106 if (originalFTS.hasError()) {
00107 const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix();
00108 AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee,
00109 pPerigee, s);
00110 const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian();
00111 CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) );
00112
00113 return FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
00114 &(originalFTS.parameters().magneticField())),
00115 cte);
00116 }
00117 else {
00118 return FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
00119 &(originalFTS.parameters().magneticField())));
00120 }
00121
00122 }
00123
00124
00125 FreeTrajectoryState
00126 TSCPBuilderNoMaterial::createFTSatTransverseImpactPointNeutral(const FTS& originalFTS,
00127 const GlobalPoint& referencePoint) const
00128 {
00129
00130 GlobalPoint xvecOrig = originalFTS.position();
00131 double phi = originalFTS.momentum().phi();
00132 double theta = originalFTS.momentum().theta();
00133 double xOrig = xvecOrig.x();
00134 double yOrig = xvecOrig.y();
00135 double zOrig = xvecOrig.z();
00136 double xR = referencePoint.x();
00137 double yR = referencePoint.y();
00138
00139 double s2D = (xR - xOrig) * cos(phi) + (yR - yOrig) * sin(phi);
00140 double s = s2D / sin(theta);
00141 double xGen = xOrig + s2D*cos(phi);
00142 double yGen = yOrig + s2D*sin(phi);
00143 double zGen = zOrig + s2D/tan(theta);
00144 GlobalPoint xPerigee = GlobalPoint(xGen, yGen, zGen);
00145
00146 GlobalVector pPerigee = originalFTS.momentum();
00147
00148 if (originalFTS.hasError()) {
00149 const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix();
00150 AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee,
00151 pPerigee, s);
00152 const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian();
00153 CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) );
00154
00155 return FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
00156 &(originalFTS.parameters().magneticField())),
00157 cte);
00158 }
00159 else {
00160 return FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
00161 &(originalFTS.parameters().magneticField())));
00162 }
00163
00164 }