#include <TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h>
Public Member Functions | |
virtual TrajectoryStateClosestToPoint | operator() (const TSOS &originalTSOS, const GlobalPoint &referencePoint) const |
virtual TrajectoryStateClosestToPoint | operator() (const FTS &originalFTS, const GlobalPoint &referencePoint) const |
virtual | ~TSCPBuilderNoMaterial () |
Private Types | |
typedef Point3DBase< double, GlobalTag > | GlobalPointDouble |
typedef Vector3DBase< double, GlobalTag > | GlobalVectorDouble |
Private Member Functions | |
FreeTrajectoryState | createFTSatTransverseImpactPoint (const FTS &originalFTS, const GlobalPoint &referencePoint) const |
FreeTrajectoryState | createFTSatTransverseImpactPointCharged (const FTS &originalFTS, const GlobalPoint &referencePoint) const |
FreeTrajectoryState | createFTSatTransverseImpactPointNeutral (const FTS &originalFTS, const GlobalPoint &referencePoint) const |
This new state is then defined at the point of closest approach to the reference point.
Definition at line 15 of file TSCPBuilderNoMaterial.h.
typedef Point3DBase< double, GlobalTag> TSCPBuilderNoMaterial::GlobalPointDouble [private] |
Definition at line 30 of file TSCPBuilderNoMaterial.h.
typedef Vector3DBase< double, GlobalTag> TSCPBuilderNoMaterial::GlobalVectorDouble [private] |
Definition at line 31 of file TSCPBuilderNoMaterial.h.
virtual TSCPBuilderNoMaterial::~TSCPBuilderNoMaterial | ( | ) | [inline, virtual] |
FreeTrajectoryState TSCPBuilderNoMaterial::createFTSatTransverseImpactPoint | ( | const FTS & | originalFTS, | |
const GlobalPoint & | referencePoint | |||
) | const [private] |
Definition at line 41 of file TSCPBuilderNoMaterial.cc.
References createFTSatTransverseImpactPointCharged(), createFTSatTransverseImpactPointNeutral(), e, and FreeTrajectoryState::transverseCurvature().
Referenced by operator()().
00043 { 00044 // 00045 // Straight line approximation? |rho|<1.e-10 equivalent to ~ 1um 00046 // difference in transversal position at 10m. 00047 // 00048 if( fabs(originalFTS.transverseCurvature())<1.e-10 ) { 00049 return createFTSatTransverseImpactPointNeutral(originalFTS, referencePoint); 00050 } else { 00051 return createFTSatTransverseImpactPointCharged(originalFTS, referencePoint); 00052 } 00053 }
FreeTrajectoryState TSCPBuilderNoMaterial::createFTSatTransverseImpactPointCharged | ( | const FTS & | originalFTS, | |
const GlobalPoint & | referencePoint | |||
) | const [private] |
Definition at line 56 of file TSCPBuilderNoMaterial.cc.
References anyDirection, FreeTrajectoryState::charge(), FreeTrajectoryState::curvilinearError(), direction, HelixBarrelPlaneCrossingByCircle::direction(), FreeTrajectoryState::hasError(), MagneticField::inInverseGeV(), PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), CurvilinearTrajectoryError::matrix(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), HelixBarrelPlaneCrossingByCircle::pathLength(), HelixBarrelPlaneCrossingByCircle::position(), FreeTrajectoryState::position(), rot, s, FreeTrajectoryState::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), Basic3DVector< T >::x(), X, PV3DBase< T, PVType, FrameType >::y(), Basic3DVector< T >::y(), PV3DBase< T, PVType, FrameType >::z(), and Basic3DVector< T >::z().
Referenced by createFTSatTransverseImpactPoint().
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 // double fac = 1./originalFTS.charge()/MagneticField::inInverseGeV(referencePoint).z(); 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 // Need to define plane with orientation as the 00082 // ImpactPointSurface 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 // Using Teddy's HelixBarrelPlaneCrossingByCircle for general barrel planes. 00088 // A large variety of other, 00089 // direct solutions turned out to be not so stable numerically. 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 // direction (reconverted to GlobalVector, renormalised) 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 }
FreeTrajectoryState TSCPBuilderNoMaterial::createFTSatTransverseImpactPointNeutral | ( | const FTS & | originalFTS, | |
const GlobalPoint & | referencePoint | |||
) | const [private] |
Definition at line 126 of file TSCPBuilderNoMaterial.cc.
References FreeTrajectoryState::charge(), funct::cos(), FreeTrajectoryState::curvilinearError(), FreeTrajectoryState::hasError(), GlobalTrajectoryParameters::magneticField(), CurvilinearTrajectoryError::matrix(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), PV3DBase< T, PVType, FrameType >::phi(), phi, FreeTrajectoryState::position(), s, funct::sin(), funct::tan(), theta, PV3DBase< T, PVType, FrameType >::theta(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by createFTSatTransverseImpactPoint().
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 }
TrajectoryStateClosestToPoint TSCPBuilderNoMaterial::operator() | ( | const TSOS & | originalTSOS, | |
const GlobalPoint & | referencePoint | |||
) | const [virtual] |
Implements TrajectoryStateClosestToPointBuilder.
Definition at line 27 of file TSCPBuilderNoMaterial.cc.
References TrajectoryStateClosestToPointBuilder::constructTSCP(), createFTSatTransverseImpactPoint(), TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalPosition(), and TrajectoryStateClosestToPointBuilder::positionEqual().
00029 { 00030 if (positionEqual(referencePoint, originalTSOS.globalPosition())) 00031 return constructTSCP(*originalTSOS.freeState(), referencePoint); 00032 00033 // Now do the propagation 00034 00035 FreeTrajectoryState ftsAtTransverseImpactPoint = 00036 createFTSatTransverseImpactPoint(*originalTSOS.freeState(), referencePoint); 00037 return constructTSCP(ftsAtTransverseImpactPoint, referencePoint); 00038 }
TrajectoryStateClosestToPoint TSCPBuilderNoMaterial::operator() | ( | const FTS & | originalFTS, | |
const GlobalPoint & | referencePoint | |||
) | const [virtual] |
Implements TrajectoryStateClosestToPointBuilder.
Definition at line 12 of file TSCPBuilderNoMaterial.cc.
References TrajectoryStateClosestToPointBuilder::constructTSCP(), createFTSatTransverseImpactPoint(), FreeTrajectoryState::position(), and TrajectoryStateClosestToPointBuilder::positionEqual().
00014 { 00015 if (positionEqual(referencePoint, originalFTS.position())) 00016 return constructTSCP(originalFTS, referencePoint); 00017 00018 // Now do the propagation or whatever... 00019 00020 FreeTrajectoryState 00021 ftsAtTransverseImpactPoint = createFTSatTransverseImpactPoint(originalFTS, referencePoint); 00022 return constructTSCP(ftsAtTransverseImpactPoint, referencePoint); 00023 00024 }