CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions
TSCPBuilderNoMaterial Class Reference

#include <TSCPBuilderNoMaterial.h>

Inheritance diagram for TSCPBuilderNoMaterial:
TrajectoryStateClosestToPointBuilder

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 ()
 
- Public Member Functions inherited from TrajectoryStateClosestToPointBuilder
bool positionEqual (const GlobalPoint &ptB, const GlobalPoint &ptA) const
 
virtual ~TrajectoryStateClosestToPointBuilder ()
 

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
 

Additional Inherited Members

- Public Types inherited from TrajectoryStateClosestToPointBuilder
typedef FreeTrajectoryState FTS
 
typedef TrajectoryStateOnSurface TSOS
 
- Protected Member Functions inherited from TrajectoryStateClosestToPointBuilder
TrajectoryStateClosestToPoint constructTSCP (const FTS &originalFTS, const GlobalPoint &referencePoint) const
 

Detailed Description

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.

Member Typedef Documentation

Definition at line 31 of file TSCPBuilderNoMaterial.h.

Definition at line 32 of file TSCPBuilderNoMaterial.h.

typedef std::pair<bool, FreeTrajectoryState> TSCPBuilderNoMaterial::PairBoolFTS
private

Definition at line 33 of file TSCPBuilderNoMaterial.h.

Constructor & Destructor Documentation

virtual TSCPBuilderNoMaterial::~TSCPBuilderNoMaterial ( )
inlinevirtual

Definition at line 21 of file TSCPBuilderNoMaterial.h.

21 {}

Member Function Documentation

TSCPBuilderNoMaterial::PairBoolFTS TSCPBuilderNoMaterial::createFTSatTransverseImpactPoint ( const FTS originalFTS,
const GlobalPoint referencePoint 
) const
private

Definition at line 49 of file TSCPBuilderNoMaterial.cc.

References createFTSatTransverseImpactPointCharged(), createFTSatTransverseImpactPointNeutral(), ExpressReco_HICollisions_FallBack::e, and FreeTrajectoryState::transverseCurvature().

Referenced by operator()().

51 {
52  //
53  // Straight line approximation? |rho|<1.e-10 equivalent to ~ 1um
54  // difference in transversal position at 10m.
55  //
56  if( fabs(originalFTS.transverseCurvature())<1.e-10 ) {
57  return createFTSatTransverseImpactPointNeutral(originalFTS, referencePoint);
58  } else {
59  return createFTSatTransverseImpactPointCharged(originalFTS, referencePoint);
60  }
61 }
PairBoolFTS createFTSatTransverseImpactPointCharged(const FTS &originalFTS, const GlobalPoint &referencePoint) const
double transverseCurvature() const
PairBoolFTS createFTSatTransverseImpactPointNeutral(const FTS &originalFTS, const GlobalPoint &referencePoint) const
TSCPBuilderNoMaterial::PairBoolFTS TSCPBuilderNoMaterial::createFTSatTransverseImpactPointCharged ( const FTS originalFTS,
const GlobalPoint referencePoint 
) const
private

Definition at line 64 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(), asciidump::s, FreeTrajectoryState::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), X, PV3DBase< T, PVType, FrameType >::x(), Basic3DVector< T >::x(), PV3DBase< T, PVType, FrameType >::y(), Basic3DVector< T >::y(), PV3DBase< T, PVType, FrameType >::z(), and Basic3DVector< T >::z().

Referenced by createFTSatTransverseImpactPoint().

66 {
67 
68  GlobalVector pvecOrig = originalFTS.momentum();
69  GlobalPoint xvecOrig = originalFTS.position();
70  double kappa = originalFTS.transverseCurvature();
71  double pxOrig = pvecOrig.x();
72  double pyOrig = pvecOrig.y();
73  double pzOrig = pvecOrig.z();
74  double xOrig = xvecOrig.x();
75  double yOrig = xvecOrig.y();
76  double zOrig = xvecOrig.z();
77 
78 // double fac = 1./originalFTS.charge()/MagneticField::inInverseGeV(referencePoint).z();
79  double fac = 1./originalFTS.charge()/
80  (originalFTS.parameters().magneticField().inInverseGeV(referencePoint).z());
81  GlobalVectorDouble xOrig2Centre = GlobalVectorDouble(fac * pyOrig, -fac * pxOrig, 0.);
82  GlobalVectorDouble xOrigProj = GlobalVectorDouble(xOrig, yOrig, 0.);
83  GlobalVectorDouble xRefProj = GlobalVectorDouble(referencePoint.x(), referencePoint.y(), 0.);
84  GlobalVectorDouble deltax = xRefProj-xOrigProj-xOrig2Centre;
85  GlobalVectorDouble ndeltax = deltax.unit();
86 
88  Surface::PositionType pos(referencePoint);
89  // Need to define plane with orientation as the
90  // ImpactPointSurface
91  GlobalVector X(ndeltax.x(), ndeltax.y(), ndeltax.z());
92  GlobalVector Y(0.,0.,1.);
93  Surface::RotationType rot(X,Y);
94  BoundPlane* plane = new BoundPlane(pos,rot);
95  // Using Teddy's HelixBarrelPlaneCrossingByCircle for general barrel planes.
96  // A large variety of other,
97  // direct solutions turned out to be not so stable numerically.
99  planeCrossing(HelixPlaneCrossing::PositionType(xOrig, yOrig, zOrig),
100  HelixPlaneCrossing::DirectionType(pxOrig, pyOrig, pzOrig),
101  kappa, direction);
102  std::pair<bool,double> propResult = planeCrossing.pathLength(*plane);
103  if ( !propResult.first ) {
104  edm::LogWarning ("TSCPBuilderNoMaterial") << "Propagation to perigee plane failed!";
105  return PairBoolFTS(false, FreeTrajectoryState() );
106  }
107  double s = propResult.second;
108  HelixPlaneCrossing::PositionType xGen = planeCrossing.position(s);
109  GlobalPoint xPerigee = GlobalPoint(xGen.x(),xGen.y(),xGen.z());
110  // direction (reconverted to GlobalVector, renormalised)
111  HelixPlaneCrossing::DirectionType pGen = planeCrossing.direction(s);
112  pGen *= pvecOrig.mag()/pGen.mag();
113  GlobalVector pPerigee = GlobalVector(pGen.x(),pGen.y(),pGen.z());
114  delete plane;
115 
116  if (originalFTS.hasError()) {
117  const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix();
118  AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee,
119  pPerigee, s);
120  const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian();
121  CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) );
122 
123  return PairBoolFTS(true,
124  FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
125  &(originalFTS.parameters().magneticField())), cte) );
126  }
127  else {
128  return PairBoolFTS(true,
129  FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
130  &(originalFTS.parameters().magneticField()))) );
131  }
132 
133 }
T y() const
Cartesian y coordinate.
T x() const
Cartesian x coordinate.
const GlobalTrajectoryParameters & parameters() const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:57
#define X(str)
Definition: MuonsGrabber.cc:49
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
PropagationDirection
TrackCharge charge() const
std::pair< bool, FreeTrajectoryState > PairBoolFTS
const CurvilinearTrajectoryError & curvilinearError() const
virtual GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
T mag() const
Definition: PV3DBase.h:61
T z() const
Cartesian z coordinate.
T z() const
Definition: PV3DBase.h:58
Vector3DBase< double, GlobalTag > GlobalVectorDouble
virtual const AlgebraicMatrix55 & jacobian() const
GlobalVector momentum() const
GlobalPoint position() const
double transverseCurvature() const
const AlgebraicSymMatrix55 & matrix() const
const MagneticField & magneticField() const
string s
Definition: asciidump.py:422
T x() const
Definition: PV3DBase.h:56
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
Global3DVector GlobalVector
Definition: GlobalVector.h:10
TSCPBuilderNoMaterial::PairBoolFTS TSCPBuilderNoMaterial::createFTSatTransverseImpactPointNeutral ( const FTS originalFTS,
const GlobalPoint referencePoint 
) const
private

Definition at line 137 of file TSCPBuilderNoMaterial.cc.

References FreeTrajectoryState::charge(), funct::cos(), FreeTrajectoryState::curvilinearError(), FreeTrajectoryState::hasError(), AnalyticalCurvilinearJacobian::jacobian(), GlobalTrajectoryParameters::magneticField(), CurvilinearTrajectoryError::matrix(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), phi, PV3DBase< T, PVType, FrameType >::phi(), FreeTrajectoryState::position(), asciidump::s, funct::sin(), funct::tan(), PV3DBase< T, PVType, FrameType >::theta(), theta(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by createFTSatTransverseImpactPoint().

139 {
140 
141  GlobalPoint xvecOrig = originalFTS.position();
142  double phi = originalFTS.momentum().phi();
143  double theta = originalFTS.momentum().theta();
144  double xOrig = xvecOrig.x();
145  double yOrig = xvecOrig.y();
146  double zOrig = xvecOrig.z();
147  double xR = referencePoint.x();
148  double yR = referencePoint.y();
149 
150  double s2D = (xR - xOrig) * cos(phi) + (yR - yOrig) * sin(phi);
151  double s = s2D / sin(theta);
152  double xGen = xOrig + s2D*cos(phi);
153  double yGen = yOrig + s2D*sin(phi);
154  double zGen = zOrig + s2D/tan(theta);
155  GlobalPoint xPerigee = GlobalPoint(xGen, yGen, zGen);
156 
157  GlobalVector pPerigee = originalFTS.momentum();
158 
159  if (originalFTS.hasError()) {
160  const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix();
161  AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee,
162  pPerigee, s);
163  const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian();
164  CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) );
165 
166  return PairBoolFTS(true,
167  FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
168  &(originalFTS.parameters().magneticField())), cte));
169  }
170  else {
171  return PairBoolFTS(true,
172  FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
173  &(originalFTS.parameters().magneticField()))) );
174  }
175 
176 }
const GlobalTrajectoryParameters & parameters() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:63
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:57
Geom::Theta< T > theta() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
TrackCharge charge() const
std::pair< bool, FreeTrajectoryState > PairBoolFTS
const CurvilinearTrajectoryError & curvilinearError() const
Geom::Theta< T > theta() const
Definition: PV3DBase.h:69
T z() const
Definition: PV3DBase.h:58
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
virtual const AlgebraicMatrix55 & jacobian() const
GlobalVector momentum() const
GlobalPoint position() const
const AlgebraicSymMatrix55 & matrix() const
const MagneticField & magneticField() const
string s
Definition: asciidump.py:422
T x() const
Definition: PV3DBase.h:56
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
Definition: DDAxes.h:10
TrajectoryStateClosestToPoint TSCPBuilderNoMaterial::operator() ( const FTS originalFTS,
const GlobalPoint referencePoint 
) const
virtual

Implements TrajectoryStateClosestToPointBuilder.

Definition at line 13 of file TSCPBuilderNoMaterial.cc.

References TrajectoryStateClosestToPointBuilder::constructTSCP(), createFTSatTransverseImpactPoint(), FreeTrajectoryState::position(), and TrajectoryStateClosestToPointBuilder::positionEqual().

15 {
16  if (positionEqual(referencePoint, originalFTS.position()))
17  return constructTSCP(originalFTS, referencePoint);
18 
19  // Now do the propagation or whatever...
20 
21  PairBoolFTS newStatePair =
22  createFTSatTransverseImpactPoint(originalFTS, referencePoint);
23  if (newStatePair.first) {
24  return constructTSCP(newStatePair.second, referencePoint);
25  } else {
27  }
28 }
TrajectoryStateClosestToPoint constructTSCP(const FTS &originalFTS, const GlobalPoint &referencePoint) const
std::pair< bool, FreeTrajectoryState > PairBoolFTS
bool positionEqual(const GlobalPoint &ptB, const GlobalPoint &ptA) const
GlobalPoint position() const
PairBoolFTS createFTSatTransverseImpactPoint(const FTS &originalFTS, const GlobalPoint &referencePoint) const
TrajectoryStateClosestToPoint TSCPBuilderNoMaterial::operator() ( const TSOS originalTSOS,
const GlobalPoint referencePoint 
) const
virtual

Implements TrajectoryStateClosestToPointBuilder.

Definition at line 31 of file TSCPBuilderNoMaterial.cc.

References TrajectoryStateClosestToPointBuilder::constructTSCP(), createFTSatTransverseImpactPoint(), TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalPosition(), and TrajectoryStateClosestToPointBuilder::positionEqual().

33 {
34  if (positionEqual(referencePoint, originalTSOS.globalPosition()))
35  return constructTSCP(*originalTSOS.freeState(), referencePoint);
36 
37  // Now do the propagation
38 
39  PairBoolFTS newStatePair =
40  createFTSatTransverseImpactPoint(*originalTSOS.freeState(), referencePoint);
41  if (newStatePair.first) {
42  return constructTSCP(newStatePair.second, referencePoint);
43  } else {
45  }
46 }
TrajectoryStateClosestToPoint constructTSCP(const FTS &originalFTS, const GlobalPoint &referencePoint) const
GlobalPoint globalPosition() const
std::pair< bool, FreeTrajectoryState > PairBoolFTS
FreeTrajectoryState * freeState(bool withErrors=true) const
bool positionEqual(const GlobalPoint &ptB, const GlobalPoint &ptA) const
PairBoolFTS createFTSatTransverseImpactPoint(const FTS &originalFTS, const GlobalPoint &referencePoint) const