CMS 3D CMS Logo

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

#include <PerigeeLinearizedTrackState.h>

Inheritance diagram for PerigeeLinearizedTrackState:
LinearizedTrackState< 5 > ReferenceCounted

Public Types

typedef
ReferenceCountingPointer
< LinearizedTrackState< 5 > > 
RefCountedLinearizedTrackState
 
- Public Types inherited from LinearizedTrackState< 5 >
typedef ROOT::Math::SMatrix
< double, N-2,
3, ROOT::Math::MatRepStd
< double, N-2, 3 > > 
AlgebraicMatrixM3
 
typedef ROOT::Math::SMatrix
< double, N,
3, ROOT::Math::MatRepStd
< double, N, 3 > > 
AlgebraicMatrixN3
 
typedef ROOT::Math::SMatrix
< double, N, N-2,
ROOT::Math::MatRepStd< double,
N, N-2 > > 
AlgebraicMatrixNM
 
typedef ROOT::Math::SMatrix
< double, N-2, N-2,
ROOT::Math::MatRepSym< double,
N-2 > > 
AlgebraicSymMatrixMM
 
typedef ROOT::Math::SMatrix
< double, N, N,
ROOT::Math::MatRepSym< double,
N > > 
AlgebraicSymMatrixNN
 
typedef ROOT::Math::SMatrix
< double, N+1, N+1,
ROOT::Math::MatRepSym< double,
N+1 > > 
AlgebraicSymMatrixOO
 
typedef ROOT::Math::SVector
< double, N-2 > 
AlgebraicVectorM
 
typedef ROOT::Math::SVector
< double, N
AlgebraicVectorN
 
typedef
ReferenceCountingPointer
< RefittedTrackState< N > > 
RefCountedRefittedTrackState
 

Public Member Functions

TrackCharge charge () const
 
virtual void checkParameters (AlgebraicVector5 &parameters) const
 
virtual std::vector
< ReferenceCountingPointer
< LinearizedTrackState< 5 > > > 
components () const
 
const AlgebraicVector5constantTerm () const
 
virtual
RefCountedRefittedTrackState 
createRefittedTrackState (const GlobalPoint &vertexPosition, const AlgebraicVector3 &vectorParameters, const AlgebraicSymMatrix66 &covarianceMatrix) const
 
bool hasError () const
 
virtual bool isValid () const
 
const GlobalPointlinearizationPoint () const
 
const AlgebraicMatrix53momentumJacobian () const
 
bool operator== (LinearizedTrackState< 5 > &other) const
 
bool operator== (ReferenceCountingPointer< LinearizedTrackState< 5 > > &other) const
 
const AlgebraicVector5parametersFromExpansion () const
 
const AlgebraicMatrix53positionJacobian () const
 
const
TrajectoryStateClosestToPoint
predictedState () const
 
AlgebraicSymMatrix55 predictedStateError () const
 
AlgebraicSymMatrix33 predictedStateMomentumError () const
 
virtual AlgebraicVector3 predictedStateMomentumParameters () const
 
AlgebraicVector5 predictedStateParameters () const
 
AlgebraicSymMatrix55 predictedStateWeight (int &error) const
 
virtual AlgebraicVector5 refittedParamFromEquation (const RefCountedRefittedTrackState &theRefittedState) const
 
const TrajectoryStateOnSurface state () const
 
virtual
RefCountedLinearizedTrackState 
stateWithNewLinearizationPoint (const GlobalPoint &newLP) const
 
virtual reco::TransientTrack track () const
 
virtual double weightInMixture () const
 
- Public Member Functions inherited from LinearizedTrackState< 5 >
virtual void checkParameters (AlgebraicVectorN &parameters) const
 
virtual
RefCountedRefittedTrackState 
createRefittedTrackState (const GlobalPoint &vertexPosition, const AlgebraicVectorM &vectorParameters, const AlgebraicSymMatrixOO &covarianceMatrix) const =0
 
virtual bool operator== (LinearizedTrackState< N > &other) const =0
 
virtual ~LinearizedTrackState ()
 

Private Member Functions

void compute3DImpactPoint () const
 
void computeChargedJacobians () const
 
void computeJacobians () const
 
void computeNeutralJacobians () const
 
 PerigeeLinearizedTrackState (const GlobalPoint &linP, const reco::TransientTrack &track, const TrajectoryStateOnSurface &tsos)
 

Private Attributes

TSCPBuilderNoMaterial builder
 
bool jacobiansAvailable
 
TrackCharge theCharge
 
AlgebraicVector5 theConstantTerm
 
AlgebraicVector5 theExpandedParams
 
GlobalPoint theLinPoint
 
AlgebraicMatrix53 theMomentumJacobian
 
AlgebraicMatrix53 thePositionJacobian
 
TrajectoryStateClosestToPoint thePredState
 
reco::TransientTrack theTrack
 
const TrajectoryStateOnSurface theTSOS
 

Friends

class LinearizedTrackStateFactory
 

Detailed Description

Calculates and stores the ImpactPointMeasurement of the impact point (point of closest approach in 3D) to the given linearization point. (see V.Karimaki, HIP-1997-77 / EXP)

Computes the parameters of the trajectory state at the transverse point of closest approach (in the global transverse plane) to the linearization point, and the jacobiam matrices. (see R.Fruehwirth et al. Data Analysis Techniques in HEP Experiments Second Edition, Cambridge University Press 2000, or R.Fruehwirth et al. Vertex reconstruction and track bundling at the LEP collider using robust algorithms. Computer Physics Communications 96 (1996) 189-208).

Both are done `on-demand' to improve CPU performance.

This particular implementation works with "perigee" helix parametrization: see Billoir et al. NIM in PR A311(1992) 139-150

Definition at line 35 of file PerigeeLinearizedTrackState.h.

Member Typedef Documentation

Definition at line 44 of file PerigeeLinearizedTrackState.h.

Constructor & Destructor Documentation

PerigeeLinearizedTrackState::PerigeeLinearizedTrackState ( const GlobalPoint linP,
const reco::TransientTrack track,
const TrajectoryStateOnSurface tsos 
)
inlineprivate

Constructor with the linearization point and the track. Private, can only be used by LinearizedTrackFactory.

Definition at line 151 of file PerigeeLinearizedTrackState.h.

Member Function Documentation

TrackCharge PerigeeLinearizedTrackState::charge ( void  ) const
inlinevirtual

Method returning the impact point measurement

Implements LinearizedTrackState< 5 >.

Definition at line 118 of file PerigeeLinearizedTrackState.h.

References theCharge.

Referenced by createRefittedTrackState().

void PerigeeLinearizedTrackState::checkParameters ( AlgebraicVector5 parameters) const
inlinevirtual

Definition at line 104 of file PerigeeLinearizedTrackState.cc.

References M_PI, and Parameters::parameters.

105 {
106  if (parameters(2) > M_PI) parameters(2)-= 2*M_PI;
107  if (parameters(2) < -M_PI) parameters(2)+= 2*M_PI;
108 }
dictionary parameters
Definition: Parameters.py:2
#define M_PI
std::vector< PerigeeLinearizedTrackState::RefCountedLinearizedTrackState > PerigeeLinearizedTrackState::components ( ) const
virtual

Implements LinearizedTrackState< 5 >.

Definition at line 75 of file PerigeeLinearizedTrackState.cc.

References query::result.

76 {
77  std::vector<RefCountedLinearizedTrackState> result; result.reserve(1);
78  result.push_back(RefCountedLinearizedTrackState(
79  const_cast<PerigeeLinearizedTrackState*>(this)));
80  return result;
81 }
ReferenceCountingPointer< LinearizedTrackState< 5 > > RefCountedLinearizedTrackState
tuple result
Definition: query.py:137
void PerigeeLinearizedTrackState::compute3DImpactPoint ( ) const
private

Method calculating the 3D impact point measurement

void PerigeeLinearizedTrackState::computeChargedJacobians ( ) const
private

Method calculating the track parameters and the Jacobians for charged particles.

Definition at line 110 of file PerigeeLinearizedTrackState.cc.

References funct::abs(), funct::cos(), reco::TransientTrack::field(), MagneticField::inInverseGeV(), M_PI, AlCaHLTBitMon_ParallelJobs::p, TrajectoryStateClosestToPoint::perigeeParameters(), PerigeeTrajectoryParameters::phi(), FreeTrajectoryState::position(), TrajectoryStateClosestToPoint::pt(), S(), funct::sin(), mathSSE::sqrt(), funct::tan(), theCharge, theConstantTerm, theExpandedParams, theLinPoint, theMomentumJacobian, thePositionJacobian, thePredState, TrajectoryStateClosestToPoint::theState(), PerigeeTrajectoryParameters::theta(), theTrack, X, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by computeJacobians().

111 {
112  GlobalPoint paramPt(theLinPoint);
113  //tarjectory parameters
114  double field = theTrack.field()->inInverseGeV(thePredState.theState().position()).z();
115  double signTC = -theCharge;
116 
117  double thetaAtEP = thePredState.perigeeParameters().theta();
118  double phiAtEP = thePredState.perigeeParameters().phi();
119  double ptAtEP = thePredState.pt();
120  double transverseCurvatureAtEP = field / ptAtEP*signTC;
121 
122  double x_v = thePredState.theState().position().x();
123  double y_v = thePredState.theState().position().y();
124  double z_v = thePredState.theState().position().z();
125  double X = x_v - paramPt.x() - sin(phiAtEP) / transverseCurvatureAtEP;
126  double Y = y_v - paramPt.y() + cos(phiAtEP) / transverseCurvatureAtEP;
127  double SS = X*X + Y*Y;
128  double S = sqrt(SS);
129 
130  // The track parameters at the expansion point
131 
132  theExpandedParams[0] = transverseCurvatureAtEP;
133  theExpandedParams[1] = thetaAtEP;
134  theExpandedParams[3] = 1/transverseCurvatureAtEP - signTC * S;
135  double phiFEP;
136  if (std::abs(X)>std::abs(Y)) {
137  double signX = (X>0.0? +1.0:-1.0);
138  phiFEP = -signTC * signX*acos(signTC*Y/S);
139  } else {
140  phiFEP = asin(-signTC*X/S);
141  if ((signTC*Y)<0.0)
142  phiFEP = M_PI - phiFEP;
143  }
144  if (phiFEP>M_PI) phiFEP-= 2*M_PI;
145  theExpandedParams[2] = phiFEP;
146  theExpandedParams[4] = z_v - paramPt.z() -
147  (phiAtEP - theExpandedParams[2]) / tan(thetaAtEP)/transverseCurvatureAtEP;
148 
149  // The Jacobian: (all at the expansion point)
150  // [i,j]
151  // i = 0: rho , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
152  // j = 0: x_v, 1: y_v, 2: z_v
153 
154  thePositionJacobian(2,0) = - Y / (SS);
155  thePositionJacobian(2,1) = X / (SS);
156  thePositionJacobian(3,0) = - signTC * X / S;
157  thePositionJacobian(3,1) = - signTC * Y / S;
158  thePositionJacobian(4,0) = thePositionJacobian(2,0)/tan(thetaAtEP)/transverseCurvatureAtEP;
159  thePositionJacobian(4,1) = thePositionJacobian(2,1)/tan(thetaAtEP)/transverseCurvatureAtEP;
160  thePositionJacobian(4,2) = 1;
161 
162  // [i,j]
163  // i = 0: rho , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
164  // j = 0: rho, 1: theta, 2: phi_v
165  theMomentumJacobian(0,0) = 1;
166  theMomentumJacobian(1,1) = 1;
167 
168  theMomentumJacobian(2,0) = -
169  (X*cos(phiAtEP) + Y*sin(phiAtEP))/
170  (SS*transverseCurvatureAtEP*transverseCurvatureAtEP);
171 
172  theMomentumJacobian(2,2) = (Y*cos(phiAtEP) - X*sin(phiAtEP)) /
173  (SS*transverseCurvatureAtEP);
174 
175  theMomentumJacobian(3,0) =
176  (signTC * (Y*cos(phiAtEP) - X*sin(phiAtEP)) / S - 1)/
177  (transverseCurvatureAtEP*transverseCurvatureAtEP);
178 
179  theMomentumJacobian(3,2) = signTC *(X*cos(phiAtEP) + Y*sin(phiAtEP))/
180  (S*transverseCurvatureAtEP);
181 
182  theMomentumJacobian(4,0) = (phiAtEP - theExpandedParams[2]) /
183  tan(thetaAtEP)/(transverseCurvatureAtEP*transverseCurvatureAtEP)+
184  theMomentumJacobian(2,0) / tan(thetaAtEP)/transverseCurvatureAtEP;
185 
186  theMomentumJacobian(4,1) = (phiAtEP - theExpandedParams[2]) *
187  (1 + 1/(tan(thetaAtEP)*tan(thetaAtEP)))/transverseCurvatureAtEP;
188 
189  theMomentumJacobian(4,2) = (theMomentumJacobian(2,2) - 1) /
190  tan(thetaAtEP)/transverseCurvatureAtEP;
191 
192  // And finally the residuals:
193 
194  auto p = thePredState.theState().position();
195  AlgebraicVector3 expansionPoint(p.x(),p.y(),p.z());
196  AlgebraicVector3 momentumAtExpansionPoint( transverseCurvatureAtEP,thetaAtEP,phiAtEP);
197 
199  thePositionJacobian * expansionPoint -
200  theMomentumJacobian * momentumAtExpansionPoint );
201 
202 }
const FreeTrajectoryState & theState() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T y() const
Definition: PV3DBase.h:63
#define X(str)
Definition: MuonsGrabber.cc:48
const MagneticField * field() const
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:39
TrajectoryStateClosestToPoint thePredState
const PerigeeTrajectoryParameters & perigeeParameters() const
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
ROOT::Math::SVector< double, 3 > AlgebraicVector3
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define M_PI
GlobalPoint position() const
ROOT::Math::SVector< double, 5 > AlgebraicVector5
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
T x() const
Definition: PV3DBase.h:62
void PerigeeLinearizedTrackState::computeJacobians ( ) const
private

Method calculating the track parameters and the Jacobians.

Definition at line 10 of file PerigeeLinearizedTrackState.cc.

References funct::abs(), builder, computeChargedJacobians(), computeNeutralJacobians(), alignCSCRings::e, reco::TransientTrack::field(), MagneticField::inInverseGeV(), TrajectoryStateClosestToPoint::isValid(), jacobiansAvailable, FreeTrajectoryState::position(), theCharge, theLinPoint, thePredState, TrajectoryStateClosestToPoint::theState(), theTrack, theTSOS, unlikely, and z.

Referenced by constantTerm(), hasError(), isValid(), momentumJacobian(), parametersFromExpansion(), positionJacobian(), predictedState(), predictedStateError(), predictedStateMomentumError(), predictedStateMomentumParameters(), predictedStateParameters(), and predictedStateWeight().

11 {
12  GlobalPoint paramPt(theLinPoint);
13 
14  thePredState = builder(theTSOS, paramPt);
15  if unlikely(!thePredState.isValid()) return;
16 
17  double field = theTrack.field()->inInverseGeV(thePredState.theState().position()).z();
18 
19  if ((std::abs(theCharge)<1e-5)||(fabs(field)<1.e-10)){
20  //neutral track
22  } else {
23  //charged track
25  }
26 
27  jacobiansAvailable = true;
28 }
const FreeTrajectoryState & theState() const
const MagneticField * field() const
#define unlikely(x)
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:39
TrajectoryStateClosestToPoint thePredState
const TrajectoryStateOnSurface theTSOS
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GlobalPoint position() const
void PerigeeLinearizedTrackState::computeNeutralJacobians ( ) const
private

Method calculating the track parameters and the Jacobians for neutral particles.

Definition at line 208 of file PerigeeLinearizedTrackState.cc.

References funct::cos(), FreeTrajectoryState::momentum(), AlCaHLTBitMon_ParallelJobs::p, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), FreeTrajectoryState::position(), funct::sin(), funct::tan(), theConstantTerm, theExpandedParams, theLinPoint, theMomentumJacobian, thePositionJacobian, thePredState, TrajectoryStateClosestToPoint::theState(), PV3DBase< T, PVType, FrameType >::theta(), X, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by computeJacobians().

209 {
210  GlobalPoint paramPt(theLinPoint);
211 
212  //tarjectory parameters
213  double thetaAtEP = thePredState.theState().momentum().theta();
214  double phiAtEP = thePredState.theState().momentum().phi();
215  double ptAtEP = thePredState.theState().momentum().perp();
216 
217  double x_v = thePredState.theState().position().x();
218  double y_v = thePredState.theState().position().y();
219  double z_v = thePredState.theState().position().z();
220  double X = x_v - paramPt.x();
221  double Y = y_v - paramPt.y();
222 
223  // The track parameters at the expansion point
224 
225  theExpandedParams(0) = 1 / ptAtEP;
226  theExpandedParams(1) = thetaAtEP;
227  theExpandedParams(2) = phiAtEP;
228  theExpandedParams(3) = X*sin(phiAtEP) - Y*cos(phiAtEP);
229  theExpandedParams(4) = z_v - paramPt.z() -
230  (X*cos(phiAtEP) + Y*sin(phiAtEP)) / tan(thetaAtEP);
231 
232  // The Jacobian: (all at the expansion point)
233  // [i,j]
234  // i = 0: rho = 1/pt , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
235  // j = 0: x_v, 1: y_v, 2: z_v
236 
237  thePositionJacobian(3,0) = sin(phiAtEP);
238  thePositionJacobian(3,1) = - cos(phiAtEP);
239  thePositionJacobian(4,0) = - cos(phiAtEP)/tan(thetaAtEP);
240  thePositionJacobian(4,1) = - sin(phiAtEP)/tan(thetaAtEP);
241  thePositionJacobian(4,2) = 1;
242 
243  // [i,j]
244  // i = 0: rho = 1/pt , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
245  // j = 0: rho = 1/pt , 1: theta, 2: phi_v
246 
247  theMomentumJacobian(0,0) = 1;
248  theMomentumJacobian(1,1) = 1;
249  theMomentumJacobian(2,2) = 1;
250 
251  theMomentumJacobian(3,2) = X*cos(phiAtEP) + Y*sin(phiAtEP);
252 
254  (1 + 1/(tan(thetaAtEP)*tan(thetaAtEP)));
255 
256  theMomentumJacobian(4,2) = (X*sin(phiAtEP) - Y*cos(phiAtEP))/tan(thetaAtEP);
257 
258  // And finally the residuals:
259 
260  auto p = thePredState.theState().position();
261  AlgebraicVector3 expansionPoint(p.x(),p.y(),p.z());
262  AlgebraicVector3 momentumAtExpansionPoint(1./ptAtEP,thetaAtEP,phiAtEP);
263 
265  thePositionJacobian * expansionPoint -
266  theMomentumJacobian * momentumAtExpansionPoint );
267 
268 }
T perp() const
Definition: PV3DBase.h:72
const FreeTrajectoryState & theState() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
#define X(str)
Definition: MuonsGrabber.cc:48
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
TrajectoryStateClosestToPoint thePredState
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
ROOT::Math::SVector< double, 3 > AlgebraicVector3
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
GlobalVector momentum() const
GlobalPoint position() const
ROOT::Math::SVector< double, 5 > AlgebraicVector5
T x() const
Definition: PV3DBase.h:62
const AlgebraicVector5 & PerigeeLinearizedTrackState::constantTerm ( ) const
inlinevirtual

Method returning the constant term of the Taylor expansion of the measurement equation

Implements LinearizedTrackState< 5 >.

Definition at line 192 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), jacobiansAvailable, theConstantTerm, and unlikely.

Referenced by TwoBodyDecayEstimator::constructMatrices(), and refittedParamFromEquation().

PerigeeLinearizedTrackState::RefCountedRefittedTrackState PerigeeLinearizedTrackState::createRefittedTrackState ( const GlobalPoint vertexPosition,
const AlgebraicVector3 vectorParameters,
const AlgebraicSymMatrix66 covarianceMatrix 
) const
virtual

Creates the correct refitted state according to the results of the track refit.

Definition at line 63 of file PerigeeLinearizedTrackState.cc.

References charge(), reco::TransientTrack::field(), theTrack, and PerigeeConversions::trajectoryStateClosestToPoint().

67 {
68  TrajectoryStateClosestToPoint refittedTSCP =
70  vectorParameters, vertexPosition, charge(), covarianceMatrix, theTrack.field());
71  return RefCountedRefittedTrackState(new PerigeeRefittedTrackState(refittedTSCP, vectorParameters));
72 }
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const AlgebraicVector3 &momentum, const GlobalPoint &referencePoint, const TrackCharge &charge, const AlgebraicSymMatrix66 &theCovarianceMatrix, const MagneticField *field)
const MagneticField * field() const
ReferenceCountingPointer< RefittedTrackState< N > > RefCountedRefittedTrackState
bool PerigeeLinearizedTrackState::hasError ( void  ) const
inlinevirtual
bool PerigeeLinearizedTrackState::isValid ( void  ) const
inlinevirtual
const GlobalPoint& PerigeeLinearizedTrackState::linearizationPoint ( ) const
inlinevirtual

The point at which the track state has been linearized

Implements LinearizedTrackState< 5 >.

Definition at line 57 of file PerigeeLinearizedTrackState.h.

References theLinPoint.

Referenced by TwoBodyDecayEstimator::constructMatrices().

const AlgebraicMatrix53 & PerigeeLinearizedTrackState::momentumJacobian ( ) const
inlinevirtual

Method returning the Momentum Jacobian from the Taylor expansion (Matrix B)

Method returning the Momentum Jacobian (Matrix B)

Implements LinearizedTrackState< 5 >.

Definition at line 212 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), jacobiansAvailable, theMomentumJacobian, and unlikely.

Referenced by TwoBodyDecayEstimator::constructMatrices(), and refittedParamFromEquation().

bool PerigeeLinearizedTrackState::operator== ( LinearizedTrackState< 5 > &  other) const

Definition at line 32 of file PerigeeLinearizedTrackState.cc.

References theTrack, and track().

33 {
34  const PerigeeLinearizedTrackState* otherP =
35  dynamic_cast<const PerigeeLinearizedTrackState*>(&other);
36  if (otherP == 0) {
37  throw VertexException("PerigeeLinearizedTrackState: don't know how to compare myself to non-perigee track state");
38  }
39  return (otherP->track() == theTrack);
40 }
Common base class.
virtual reco::TransientTrack track() const
bool PerigeeLinearizedTrackState::operator== ( ReferenceCountingPointer< LinearizedTrackState< 5 > > &  other) const

Definition at line 43 of file PerigeeLinearizedTrackState.cc.

References theTrack, and track().

44 {
45  const PerigeeLinearizedTrackState* otherP =
46  dynamic_cast<const PerigeeLinearizedTrackState*>(other.get());
47  if (otherP == 0) {
48  throw VertexException("PerigeeLinearizedTrackState: don't know how to compare myself to non-perigee track state");
49  }
50  return (otherP->track() == theTrack);
51 }
Common base class.
virtual reco::TransientTrack track() const
const AlgebraicVector5 & PerigeeLinearizedTrackState::parametersFromExpansion ( ) const
inlinevirtual

Method returning the parameters of the Taylor expansion

Implements LinearizedTrackState< 5 >.

Definition at line 221 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), jacobiansAvailable, theExpandedParams, and unlikely.

const AlgebraicMatrix53 & PerigeeLinearizedTrackState::positionJacobian ( ) const
inlinevirtual

Method returning the Position Jacobian from the Taylor expansion (Matrix A)

Method returning the Position Jacobian (Matrix A)

Implements LinearizedTrackState< 5 >.

Definition at line 202 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), jacobiansAvailable, thePositionJacobian, and unlikely.

Referenced by TwoBodyDecayEstimator::constructMatrices(), and refittedParamFromEquation().

const TrajectoryStateClosestToPoint & PerigeeLinearizedTrackState::predictedState ( ) const
inline

Method returning the track state at the point of closest approach to the linearization point, in the transverse plane (a.k.a. transverse impact point).

Method returning the TrajectoryStateClosestToPoint at the point of closest approch to the z-axis (a.k.a. transverse impact point)

Definition at line 232 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), jacobiansAvailable, thePredState, and unlikely.

Referenced by TwoBodyDecayLinearizationPointFinder::getLinearizationPoint(), and PerigeeMultiLTS::predictedState().

233 {
236 }
return((rh^lh)&mask)
#define unlikely(x)
TrajectoryStateClosestToPoint thePredState
AlgebraicSymMatrix55 PerigeeLinearizedTrackState::predictedStateError ( ) const
inlinevirtual

Method returning the covariance matrix of the track state at the transverse impact point.

Implements LinearizedTrackState< 5 >.

Definition at line 274 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), PerigeeTrajectoryError::covarianceMatrix(), jacobiansAvailable, TrajectoryStateClosestToPoint::perigeeError(), thePredState, and unlikely.

Referenced by TwoBodyDecayEstimator::constructMatrices().

275 {
277  return thePredState.perigeeError().covarianceMatrix();
278 }
return((rh^lh)&mask)
#define unlikely(x)
TrajectoryStateClosestToPoint thePredState
AlgebraicSymMatrix33 PerigeeLinearizedTrackState::predictedStateMomentumError ( ) const
inlinevirtual

Method returning the momentum covariance matrix of the track state at the transverse impact point.

Implements LinearizedTrackState< 5 >.

Definition at line 281 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), PerigeeTrajectoryError::covarianceMatrix(), jacobiansAvailable, TrajectoryStateClosestToPoint::perigeeError(), thePredState, and unlikely.

282 {
284  return thePredState.perigeeError().covarianceMatrix().Sub<AlgebraicSymMatrix33>(0,2);
285 }
return((rh^lh)&mask)
#define unlikely(x)
TrajectoryStateClosestToPoint thePredState
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
AlgebraicVector3 PerigeeLinearizedTrackState::predictedStateMomentumParameters ( ) const
inlinevirtual

Method returning the momentum part of the parameters of the track state at the linearization point.

Implements LinearizedTrackState< 5 >.

Definition at line 254 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), jacobiansAvailable, TrajectoryStateClosestToPoint::perigeeParameters(), thePredState, unlikely, findQualityFiles::v, and PerigeeTrajectoryParameters::vector().

Referenced by TwoBodyDecayEstimator::constructMatrices(), and refittedParamFromEquation().

255 {
257  auto v= thePredState.perigeeParameters().vector();
258  return AlgebraicVector3(v[0],v[1],v[2]);
259 
260 }
return((rh^lh)&mask)
#define unlikely(x)
TrajectoryStateClosestToPoint thePredState
ROOT::Math::SVector< double, 3 > AlgebraicVector3
AlgebraicVector5 PerigeeLinearizedTrackState::predictedStateParameters ( ) const
inlinevirtual

Method returning the parameters of the track state at the transverse impact point.

Implements LinearizedTrackState< 5 >.

Definition at line 247 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), jacobiansAvailable, TrajectoryStateClosestToPoint::perigeeParameters(), thePredState, and PerigeeTrajectoryParameters::vector().

Referenced by TwoBodyDecayEstimator::constructMatrices().

248 {
251 }
const AlgebraicVector5 & vector() const
TrajectoryStateClosestToPoint thePredState
const PerigeeTrajectoryParameters & perigeeParameters() const
AlgebraicSymMatrix55 PerigeeLinearizedTrackState::predictedStateWeight ( int &  error) const
inlinevirtual

Method returning the weight matrix of the track state at the transverse impact point. The error variable is 0 in case of success.

Implements LinearizedTrackState< 5 >.

Definition at line 263 of file PerigeeLinearizedTrackState.h.

References computeJacobians(), TrajectoryStateClosestToPoint::isValid(), jacobiansAvailable, TrajectoryStateClosestToPoint::perigeeError(), thePredState, unlikely, and PerigeeTrajectoryError::weightMatrix().

264 {
266  if (!thePredState.isValid()) {
267  error = 1;
268  return AlgebraicSymMatrix55();
269  }
271 }
const PerigeeTrajectoryError & perigeeError() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
#define unlikely(x)
TrajectoryStateClosestToPoint thePredState
if(c.getParameter< edm::InputTag >("puppiValueMap").label().size()!=0)
const AlgebraicSymMatrix55 & weightMatrix(int &error) const
AlgebraicVector5 PerigeeLinearizedTrackState::refittedParamFromEquation ( const RefCountedRefittedTrackState theRefittedState) const
virtual

Method returning the parameters of the Taylor expansion evaluated with the refitted state.

Implements LinearizedTrackState< 5 >.

Definition at line 84 of file PerigeeLinearizedTrackState.cc.

References constantTerm(), M_PI, momentumJacobian(), AlCaHLTBitMon_ParallelJobs::p, positionJacobian(), and predictedStateMomentumParameters().

86 {
87  auto p = theRefittedState->position();
88  AlgebraicVector3 vertexPosition(p.x(),p.y(),p.z());
89  AlgebraicVector3 momentum = theRefittedState->momentumVector();
90  if ((momentum(2)*predictedStateMomentumParameters()(2) < 0)&&(fabs(momentum(2))>M_PI/2) ) {
91  if (predictedStateMomentumParameters()(2) < 0.) momentum(2)-= 2*M_PI;
92  if (predictedStateMomentumParameters()(2) > 0.) momentum(2)+= 2*M_PI;
93  }
94  AlgebraicVectorN param = constantTerm() +
95  positionJacobian() * vertexPosition +
96  momentumJacobian() * momentum;
97  if (param(2) > M_PI) param(2)-= 2*M_PI;
98  if (param(2) < -M_PI) param(2)+= 2*M_PI;
99 
100  return param;
101 }
virtual AlgebraicVector3 predictedStateMomentumParameters() const
const AlgebraicMatrix53 & momentumJacobian() const
ROOT::Math::SVector< double, N > AlgebraicVectorN
const AlgebraicMatrix53 & positionJacobian() const
ROOT::Math::SVector< double, 3 > AlgebraicVector3
#define M_PI
const AlgebraicVector5 & constantTerm() const
const TrajectoryStateOnSurface PerigeeLinearizedTrackState::state ( ) const
inline

Definition at line 61 of file PerigeeLinearizedTrackState.h.

References theTSOS.

61 { return theTSOS; }
const TrajectoryStateOnSurface theTSOS
PerigeeLinearizedTrackState::RefCountedLinearizedTrackState PerigeeLinearizedTrackState::stateWithNewLinearizationPoint ( const GlobalPoint newLP) const
virtual

Returns a new linearized state with respect to a new linearization point. A new object of the same type is returned, without change to the existing one.

Implements LinearizedTrackState< 5 >.

Definition at line 56 of file PerigeeLinearizedTrackState.cc.

57 {
60 }
virtual reco::TransientTrack track() const
ReferenceCountingPointer< LinearizedTrackState< 5 > > RefCountedLinearizedTrackState
const TrajectoryStateOnSurface theTSOS
PerigeeLinearizedTrackState(const GlobalPoint &linP, const reco::TransientTrack &track, const TrajectoryStateOnSurface &tsos)
virtual reco::TransientTrack PerigeeLinearizedTrackState::track ( ) const
inlinevirtual
virtual double PerigeeLinearizedTrackState::weightInMixture ( ) const
inlinevirtual

Implements LinearizedTrackState< 5 >.

Definition at line 138 of file PerigeeLinearizedTrackState.h.

References theTSOS, and TrajectoryStateOnSurface::weight().

138 {return theTSOS.weight();}
const TrajectoryStateOnSurface theTSOS

Friends And Related Function Documentation

friend class LinearizedTrackStateFactory
friend

Friend class properly dealing with creation of reference-counted pointers to LinearizedTrack objects

Definition at line 43 of file PerigeeLinearizedTrackState.h.

Member Data Documentation

TSCPBuilderNoMaterial PerigeeLinearizedTrackState::builder
private

Definition at line 178 of file PerigeeLinearizedTrackState.h.

Referenced by computeJacobians().

bool PerigeeLinearizedTrackState::jacobiansAvailable
mutableprivate
TrackCharge PerigeeLinearizedTrackState::theCharge
private
AlgebraicVector5 PerigeeLinearizedTrackState::theConstantTerm
mutableprivate
AlgebraicVector5 PerigeeLinearizedTrackState::theExpandedParams
mutableprivate
GlobalPoint PerigeeLinearizedTrackState::theLinPoint
private
AlgebraicMatrix53 PerigeeLinearizedTrackState::theMomentumJacobian
mutableprivate
AlgebraicMatrix53 PerigeeLinearizedTrackState::thePositionJacobian
mutableprivate
TrajectoryStateClosestToPoint PerigeeLinearizedTrackState::thePredState
mutableprivate
reco::TransientTrack PerigeeLinearizedTrackState::theTrack
private
const TrajectoryStateOnSurface PerigeeLinearizedTrackState::theTSOS
private

Definition at line 171 of file PerigeeLinearizedTrackState.h.

Referenced by computeJacobians(), isValid(), state(), and weightInMixture().