#include <TSCPBuilderNoMaterial.h>
Public Member Functions | |
virtual TrajectoryStateClosestToPoint | operator() (const FTS &originalFTS, const GlobalPoint &referencePoint) const |
virtual TrajectoryStateClosestToPoint | operator() (const TSOS &originalTSOS, const GlobalPoint &referencePoint) const |
virtual | ~TSCPBuilderNoMaterial () |
Private Types | |
typedef Point3DBase< double, GlobalTag > | GlobalPointDouble |
typedef Vector3DBase< double, GlobalTag > | GlobalVectorDouble |
typedef std::pair< bool, FreeTrajectoryState > | PairBoolFTS |
Private Member Functions | |
PairBoolFTS | createFTSatTransverseImpactPoint (const FTS &originalFTS, const GlobalPoint &referencePoint) const |
PairBoolFTS | createFTSatTransverseImpactPointCharged (const FTS &originalFTS, const GlobalPoint &referencePoint) const |
PairBoolFTS | createFTSatTransverseImpactPointNeutral (const FTS &originalFTS, const GlobalPoint &referencePoint) const |
This class builds a TrajectoryStateClosestToPoint given an original TrajectoryStateOnSurface or FreeTrajectoryState. This new state is then defined at the point of closest approach to the reference point. In case the propagation was not successful, this state can be invalid.
Definition at line 16 of file TSCPBuilderNoMaterial.h.
typedef Point3DBase< double, GlobalTag> TSCPBuilderNoMaterial::GlobalPointDouble [private] |
Definition at line 31 of file TSCPBuilderNoMaterial.h.
typedef Vector3DBase< double, GlobalTag> TSCPBuilderNoMaterial::GlobalVectorDouble [private] |
Definition at line 32 of file TSCPBuilderNoMaterial.h.
typedef std::pair<bool, FreeTrajectoryState> TSCPBuilderNoMaterial::PairBoolFTS [private] |
Definition at line 33 of file TSCPBuilderNoMaterial.h.
virtual TSCPBuilderNoMaterial::~TSCPBuilderNoMaterial | ( | ) | [inline, virtual] |
Definition at line 21 of file TSCPBuilderNoMaterial.h.
{}
TSCPBuilderNoMaterial::PairBoolFTS TSCPBuilderNoMaterial::createFTSatTransverseImpactPoint | ( | const FTS & | originalFTS, |
const GlobalPoint & | referencePoint | ||
) | const [private] |
Definition at line 48 of file TSCPBuilderNoMaterial.cc.
References createFTSatTransverseImpactPointCharged(), createFTSatTransverseImpactPointNeutral(), alignCSCRings::e, and FreeTrajectoryState::transverseCurvature().
Referenced by operator()().
{ // // Straight line approximation? |rho|<1.e-10 equivalent to ~ 1um // difference in transversal position at 10m. // if( fabs(originalFTS.transverseCurvature())<1.e-10 ) { return createFTSatTransverseImpactPointNeutral(originalFTS, referencePoint); } else { return createFTSatTransverseImpactPointCharged(originalFTS, referencePoint); } }
TSCPBuilderNoMaterial::PairBoolFTS TSCPBuilderNoMaterial::createFTSatTransverseImpactPointCharged | ( | const FTS & | originalFTS, |
const GlobalPoint & | referencePoint | ||
) | const [private] |
Definition at line 63 of file TSCPBuilderNoMaterial.cc.
References anyDirection, FreeTrajectoryState::charge(), FreeTrajectoryState::curvilinearError(), HelixBarrelPlaneCrossingByCircle::direction(), FreeTrajectoryState::hasError(), MagneticField::inInverseGeV(), AnalyticalCurvilinearJacobian::jacobian(), PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), CurvilinearTrajectoryError::matrix(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), HelixBarrelPlaneCrossingByCircle::pathLength(), pos, HelixBarrelPlaneCrossingByCircle::position(), FreeTrajectoryState::position(), makeMuonMisalignmentScenario::rot, alignCSCRings::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().
{ GlobalVector pvecOrig = originalFTS.momentum(); GlobalPoint xvecOrig = originalFTS.position(); double kappa = originalFTS.transverseCurvature(); double pxOrig = pvecOrig.x(); double pyOrig = pvecOrig.y(); double pzOrig = pvecOrig.z(); double xOrig = xvecOrig.x(); double yOrig = xvecOrig.y(); double zOrig = xvecOrig.z(); // double fac = 1./originalFTS.charge()/MagneticField::inInverseGeV(referencePoint).z(); double fac = 1./originalFTS.charge()/ (originalFTS.parameters().magneticField().inInverseGeV(referencePoint).z()); GlobalVectorDouble xOrig2Centre = GlobalVectorDouble(fac * pyOrig, -fac * pxOrig, 0.); GlobalVectorDouble xOrigProj = GlobalVectorDouble(xOrig, yOrig, 0.); GlobalVectorDouble xRefProj = GlobalVectorDouble(referencePoint.x(), referencePoint.y(), 0.); GlobalVectorDouble deltax = xRefProj-xOrigProj-xOrig2Centre; GlobalVectorDouble ndeltax = deltax.unit(); PropagationDirection direction = anyDirection; Surface::PositionType pos(referencePoint); // Need to define plane with orientation as the // ImpactPointSurface GlobalVector X(ndeltax.x(), ndeltax.y(), ndeltax.z()); GlobalVector Y(0.,0.,1.); Surface::RotationType rot(X,Y); BoundPlane* plane = new BoundPlane(pos,rot); // Using Teddy's HelixBarrelPlaneCrossingByCircle for general barrel planes. // A large variety of other, // direct solutions turned out to be not so stable numerically. HelixBarrelPlaneCrossingByCircle planeCrossing(HelixPlaneCrossing::PositionType(xOrig, yOrig, zOrig), HelixPlaneCrossing::DirectionType(pxOrig, pyOrig, pzOrig), kappa, direction); std::pair<bool,double> propResult = planeCrossing.pathLength(*plane); if ( !propResult.first ) { edm::LogWarning ("TSCPBuilderNoMaterial") << "Propagation to perigee plane failed!"; return PairBoolFTS(false, FreeTrajectoryState() ); } double s = propResult.second; HelixPlaneCrossing::PositionType xGen = planeCrossing.position(s); GlobalPoint xPerigee = GlobalPoint(xGen.x(),xGen.y(),xGen.z()); // direction (reconverted to GlobalVector, renormalised) HelixPlaneCrossing::DirectionType pGen = planeCrossing.direction(s); pGen *= pvecOrig.mag()/pGen.mag(); GlobalVector pPerigee = GlobalVector(pGen.x(),pGen.y(),pGen.z()); delete plane; if (originalFTS.hasError()) { const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix(); AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee, pPerigee, s); const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian(); CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) ); return PairBoolFTS(true, FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(), &(originalFTS.parameters().magneticField())), cte) ); } else { return PairBoolFTS(true, FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(), &(originalFTS.parameters().magneticField()))) ); } }
TSCPBuilderNoMaterial::PairBoolFTS TSCPBuilderNoMaterial::createFTSatTransverseImpactPointNeutral | ( | const FTS & | originalFTS, |
const GlobalPoint & | referencePoint | ||
) | const [private] |
Definition at line 136 of file TSCPBuilderNoMaterial.cc.
References FreeTrajectoryState::charge(), funct::cos(), FreeTrajectoryState::curvilinearError(), FreeTrajectoryState::hasError(), AnalyticalCurvilinearJacobian::jacobian(), GlobalTrajectoryParameters::magneticField(), CurvilinearTrajectoryError::matrix(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), PV3DBase< T, PVType, FrameType >::phi(), phi, FreeTrajectoryState::position(), alignCSCRings::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().
{ GlobalPoint xvecOrig = originalFTS.position(); double phi = originalFTS.momentum().phi(); double theta = originalFTS.momentum().theta(); double xOrig = xvecOrig.x(); double yOrig = xvecOrig.y(); double zOrig = xvecOrig.z(); double xR = referencePoint.x(); double yR = referencePoint.y(); double s2D = (xR - xOrig) * cos(phi) + (yR - yOrig) * sin(phi); double s = s2D / sin(theta); double xGen = xOrig + s2D*cos(phi); double yGen = yOrig + s2D*sin(phi); double zGen = zOrig + s2D/tan(theta); GlobalPoint xPerigee = GlobalPoint(xGen, yGen, zGen); GlobalVector pPerigee = originalFTS.momentum(); if (originalFTS.hasError()) { const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix(); AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee, pPerigee, s); const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian(); CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) ); return PairBoolFTS(true, FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(), &(originalFTS.parameters().magneticField())), cte)); } else { return PairBoolFTS(true, FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(), &(originalFTS.parameters().magneticField()))) ); } }
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().
{ if (positionEqual(referencePoint, originalFTS.position())) return constructTSCP(originalFTS, referencePoint); // Now do the propagation or whatever... PairBoolFTS newStatePair = createFTSatTransverseImpactPoint(originalFTS, referencePoint); if (newStatePair.first) { return constructTSCP(newStatePair.second, referencePoint); } else { return TrajectoryStateClosestToPoint(); } }
TrajectoryStateClosestToPoint TSCPBuilderNoMaterial::operator() | ( | const TSOS & | originalTSOS, |
const GlobalPoint & | referencePoint | ||
) | const [virtual] |
Implements TrajectoryStateClosestToPointBuilder.
Definition at line 30 of file TSCPBuilderNoMaterial.cc.
References TrajectoryStateClosestToPointBuilder::constructTSCP(), createFTSatTransverseImpactPoint(), TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalPosition(), and TrajectoryStateClosestToPointBuilder::positionEqual().
{ if (positionEqual(referencePoint, originalTSOS.globalPosition())) return constructTSCP(*originalTSOS.freeState(), referencePoint); // Now do the propagation PairBoolFTS newStatePair = createFTSatTransverseImpactPoint(*originalTSOS.freeState(), referencePoint); if (newStatePair.first) { return constructTSCP(newStatePair.second, referencePoint); } else { return TrajectoryStateClosestToPoint(); } }