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