CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes | Friends
ParticleKinematicLinearizedTrackState Class Reference

#include <ParticleKinematicLinearizedTrackState.h>

Inheritance diagram for ParticleKinematicLinearizedTrackState:
LinearizedTrackState< 6 > ReferenceCounted

Public Types

typedef ReferenceCountingPointer< LinearizedTrackState< 6 > > RefCountedLinearizedTrackState
 
- Public Types inherited from LinearizedTrackState< 6 >
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, NAlgebraicVectorN
 
typedef ReferenceCountingPointer< RefittedTrackState< N > > RefCountedRefittedTrackState
 

Public Member Functions

TrackCharge charge () const override
 
void checkParameters (AlgebraicVectorN &parameters) const override
 
std::vector< ReferenceCountingPointer< LinearizedTrackState< 6 > > > components () const override
 
const AlgebraicVector6constantTerm () const override
 
RefCountedRefittedTrackState createRefittedTrackState (const GlobalPoint &vertexPosition, const AlgebraicVectorM &vectorParameters, const AlgebraicSymMatrix77 &covarianceMatrix) const override
 
bool hasError () const override
 
const GlobalPointlinearizationPoint () const override
 
const AlgebraicMatrix64momentumJacobian () const override
 
bool operator== (LinearizedTrackState< 6 > &other) const override
 
const AlgebraicVector6parametersFromExpansion () const override
 
RefCountedKinematicParticle particle () const
 
 ParticleKinematicLinearizedTrackState ()
 
const AlgebraicMatrix63positionJacobian () const override
 
AlgebraicSymMatrix66 predictedStateError () const override
 
AlgebraicSymMatrix44 predictedStateMomentumError () const override
 
AlgebraicVectorM predictedStateMomentumParameters () const override
 
AlgebraicVector6 predictedStateParameters () const override
 
AlgebraicSymMatrix66 predictedStateWeight (int &error) const override
 
AlgebraicVectorN refittedParamFromEquation (const RefCountedRefittedTrackState &theRefittedState) const override
 
ReferenceCountingPointer< LinearizedTrackState< 6 > > stateWithNewLinearizationPoint (const GlobalPoint &newLP) const override
 
reco::TransientTrack track () const override
 
double weightInMixture () const override
 
- Public Member Functions inherited from LinearizedTrackState< 6 >
virtual void checkParameters (AlgebraicVectorN &parameters) const
 
virtual RefCountedRefittedTrackState createRefittedTrackState (const GlobalPoint &vertexPosition, const AlgebraicVectorM &vectorParameters, const AlgebraicSymMatrixOO &covarianceMatrix) const=0
 
virtual bool isValid () const
 
virtual bool operator== (LinearizedTrackState< N > &other) const=0
 
virtual AlgebraicVectorN refittedParamFromEquation (const RefCountedRefittedTrackState &theRefittedState) const=0
 
 ~LinearizedTrackState () override
 

Private Member Functions

void computeChargedJacobians () const
 
void computeJacobians () const
 
void computeNeutralJacobians () const
 
 ParticleKinematicLinearizedTrackState (const GlobalPoint &linP, RefCountedKinematicParticle &prt)
 

Private Attributes

TransientTrackKinematicStateBuilder builder
 
bool errorAvailable
 
bool impactPointAvailable
 
bool jacobiansAvailable
 
RefCountedKinematicParticle part
 
TrackCharge theCharge
 
AlgebraicVector6 theConstantTerm
 
AlgebraicVector6 theExpandedParams
 
GlobalPoint theLinPoint
 
AlgebraicMatrix64 theMomentumJacobian
 
AlgebraicMatrix63 thePositionJacobian
 
PerigeeKinematicState thePredState
 

Friends

class ParticleKinematicLinearizedTrackStateFactory
 

Detailed Description

Definition at line 13 of file ParticleKinematicLinearizedTrackState.h.

Member Typedef Documentation

◆ RefCountedLinearizedTrackState

Definition at line 16 of file ParticleKinematicLinearizedTrackState.h.

Constructor & Destructor Documentation

◆ ParticleKinematicLinearizedTrackState() [1/2]

ParticleKinematicLinearizedTrackState::ParticleKinematicLinearizedTrackState ( )
inline

◆ ParticleKinematicLinearizedTrackState() [2/2]

ParticleKinematicLinearizedTrackState::ParticleKinematicLinearizedTrackState ( const GlobalPoint linP,
RefCountedKinematicParticle prt 
)
inlineprivate

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

Definition at line 115 of file ParticleKinematicLinearizedTrackState.h.

Member Function Documentation

◆ charge()

TrackCharge ParticleKinematicLinearizedTrackState::charge ( void  ) const
overridevirtual

Method returning the impact point measurement

Implements LinearizedTrackState< 6 >.

Definition at line 54 of file ParticleKinematicLinearizedTrackState.cc.

Referenced by createRefittedTrackState().

54 { return part->initialState().particleCharge(); }
part
Definition: HCALResponse.h:20

◆ checkParameters()

void ParticleKinematicLinearizedTrackState::checkParameters ( AlgebraicVectorN parameters) const
inlineoverride

Definition at line 172 of file ParticleKinematicLinearizedTrackState.cc.

References M_PI.

172  {
173  if (parameters(2) > M_PI)
174  parameters(2) -= 2 * M_PI;
175  if (parameters(2) < -M_PI)
176  parameters(2) += 2 * M_PI;
177 }
#define M_PI

◆ components()

std::vector< ReferenceCountingPointer< LinearizedTrackState< 6 > > > ParticleKinematicLinearizedTrackState::components ( ) const
overridevirtual

Implements LinearizedTrackState< 6 >.

Definition at line 140 of file ParticleKinematicLinearizedTrackState.cc.

141  {
142  std::vector<ReferenceCountingPointer<LinearizedTrackState<6> > > res;
143  res.reserve(1);
144  res.push_back(RefCountedLinearizedTrackState(const_cast<ParticleKinematicLinearizedTrackState*>(this)));
145  return res;
146 }
Definition: Electron.h:6
ReferenceCountingPointer< LinearizedTrackState< 6 > > RefCountedLinearizedTrackState

◆ computeChargedJacobians()

void ParticleKinematicLinearizedTrackState::computeChargedJacobians ( ) const
private

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

Definition at line 179 of file ParticleKinematicLinearizedTrackState.cc.

References funct::abs(), funct::cos(), KinematicState::globalMomentum(), KinematicState::globalPosition(), M_PI, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), funct::sin(), mathSSE::sqrt(), funct::tan(), theConstantTerm, theExpandedParams, theLinPoint, theMomentumJacobian, thePositionJacobian, thePredState, PerigeeKinematicState::theState(), PV3DBase< T, PVType, FrameType >::theta(), X, PV3DBase< T, PVType, FrameType >::x(), BeamSpotPI::Y, PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by computeJacobians().

179  {
180  GlobalPoint paramPt(theLinPoint);
181 
182  double field = part->currentState().magneticField()->inInverseGeV(thePredState.theState().globalPosition()).z();
183  double signTC = -part->currentState().particleCharge();
184 
185  double thetaAtEP = thePredState.theState().globalMomentum().theta();
186  double phiAtEP = thePredState.theState().globalMomentum().phi();
187  double ptAtEP = thePredState.theState().globalMomentum().perp();
188  double transverseCurvatureAtEP = field / ptAtEP * signTC;
189 
190  // Fix calculation for case where magnetic field swaps sign between previous state and current state
191  if (field * part->currentState().magneticField()->inInverseGeV(part->currentState().globalPosition()).z() < 0.) {
192  signTC = -signTC;
193  }
194 
195  double x_v = thePredState.theState().globalPosition().x();
196  double y_v = thePredState.theState().globalPosition().y();
197  double z_v = thePredState.theState().globalPosition().z();
198  double X = x_v - paramPt.x() - sin(phiAtEP) / transverseCurvatureAtEP;
199  double Y = y_v - paramPt.y() + cos(phiAtEP) / transverseCurvatureAtEP;
200  double SS = X * X + Y * Y;
201  double S = sqrt(SS);
202 
203  // The track parameters at the expansion point
204  theExpandedParams[0] = transverseCurvatureAtEP;
205  theExpandedParams[1] = thetaAtEP;
206  theExpandedParams[3] = 1 / transverseCurvatureAtEP - signTC * S;
207 
208  theExpandedParams[5] = part->currentState().mass();
209 
210  double phiFEP;
211  if (std::abs(X) > std::abs(Y)) {
212  double signX = (X > 0.0 ? +1.0 : -1.0);
213  phiFEP = -signTC * signX * acos(signTC * Y / S);
214  } else {
215  phiFEP = asin(-signTC * X / S);
216  if ((signTC * Y) < 0.0)
217  phiFEP = M_PI - phiFEP;
218  }
219  if (phiFEP > M_PI)
220  phiFEP -= 2 * M_PI;
221  theExpandedParams[2] = phiFEP;
222  theExpandedParams[4] =
223  z_v - paramPt.z() - (phiAtEP - theExpandedParams[2]) / tan(thetaAtEP) / transverseCurvatureAtEP;
224 
225  thePositionJacobian(2, 0) = -Y / (SS);
226  thePositionJacobian(2, 1) = X / (SS);
227  thePositionJacobian(3, 0) = -signTC * X / S;
228  thePositionJacobian(3, 1) = -signTC * Y / S;
229  thePositionJacobian(4, 0) = thePositionJacobian(2, 0) / tan(thetaAtEP) / transverseCurvatureAtEP;
230  thePositionJacobian(4, 1) = thePositionJacobian(2, 1) / tan(thetaAtEP) / transverseCurvatureAtEP;
231  thePositionJacobian(4, 2) = 1;
232 
233  //debug code - to be removed later
234  // cout<<"parameters for momentum jacobian: X "<<X<<endl;
235  // cout<<"parameters for momentum jacobian: Y "<<Y<<endl;
236  // cout<<"parameters for momentum jacobian: SS "<<SS<<endl;
237  // cout<<"parameters for momentum jacobian: PhiAtEP "<<phiAtEP<<endl;
238  // cout<<"parameters for momentum jacobian: curv "<<transverseCurvatureAtEP<<endl;
239  // cout<<"sin phi Atep "<<sin(phiAtEP)<<endl;
240  // cout<<"cos phi at EP "<<cos(phiAtEP)<<endl;
241  // cout<<"upper part is "<<X*cos(phiAtEP) + Y*sin(phiAtEP) <<endl;
242  // cout<<"lower part is"<<SS*transverseCurvatureAtEP*transverseCurvatureAtEP<<endl;
243 
244  theMomentumJacobian(0, 0) = 1;
245  theMomentumJacobian(1, 1) = 1;
246 
247  theMomentumJacobian(2, 0) =
248  -(X * cos(phiAtEP) + Y * sin(phiAtEP)) / (SS * transverseCurvatureAtEP * transverseCurvatureAtEP);
249 
250  theMomentumJacobian(2, 2) = (Y * cos(phiAtEP) - X * sin(phiAtEP)) / (SS * transverseCurvatureAtEP);
251 
252  theMomentumJacobian(3, 0) =
253  (signTC * (Y * cos(phiAtEP) - X * sin(phiAtEP)) / S - 1) / (transverseCurvatureAtEP * transverseCurvatureAtEP);
254 
255  theMomentumJacobian(3, 2) = signTC * (X * cos(phiAtEP) + Y * sin(phiAtEP)) / (S * transverseCurvatureAtEP);
256 
257  theMomentumJacobian(4, 0) =
258  (phiAtEP - theExpandedParams(2)) / tan(thetaAtEP) / (transverseCurvatureAtEP * transverseCurvatureAtEP) +
259  theMomentumJacobian(2, 0) / tan(thetaAtEP) / transverseCurvatureAtEP;
260 
261  theMomentumJacobian(4, 1) =
262  (phiAtEP - theExpandedParams(2)) * (1 + 1 / (tan(thetaAtEP) * tan(thetaAtEP))) / transverseCurvatureAtEP;
263 
264  theMomentumJacobian(4, 2) = (theMomentumJacobian(2, 2) - 1) / tan(thetaAtEP) / transverseCurvatureAtEP;
265 
266  theMomentumJacobian(5, 3) = 1;
267 
268  // And finally the residuals:
269  AlgebraicVector3 expansionPoint;
270  expansionPoint[0] = thePredState.theState().globalPosition().x();
271  expansionPoint[1] = thePredState.theState().globalPosition().y();
272  expansionPoint[2] = thePredState.theState().globalPosition().z();
273  AlgebraicVector4 momentumAtExpansionPoint;
274  momentumAtExpansionPoint[0] = transverseCurvatureAtEP; // Transverse Curv
275  momentumAtExpansionPoint[1] = thetaAtEP;
276  momentumAtExpansionPoint[2] = phiAtEP;
277  momentumAtExpansionPoint[3] = theExpandedParams[5];
278 
280  theMomentumJacobian * momentumAtExpansionPoint);
281 }
T perp() const
Definition: PV3DBase.h:69
T z() const
Definition: PV3DBase.h:61
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define X(str)
Definition: MuonsGrabber.cc:38
GlobalPoint globalPosition() const
const KinematicState & theState() const
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
GlobalVector globalMomentum() const
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ROOT::Math::SVector< double, 4 > AlgebraicVector4
#define M_PI
part
Definition: HCALResponse.h:20
ROOT::Math::SVector< double, 3 > AlgebraicVector3
ROOT::Math::SVector< double, 6 > AlgebraicVector6
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72

◆ computeJacobians()

void ParticleKinematicLinearizedTrackState::computeJacobians ( ) const
private

Method calculating the track parameters and the Jacobians.

Definition at line 79 of file ParticleKinematicLinearizedTrackState.cc.

References funct::abs(), builder, computeChargedJacobians(), computeNeutralJacobians(), MillePedeFileConverter_cfg::e, jacobiansAvailable, theCharge, theLinPoint, and thePredState.

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

79  {
80  GlobalPoint paramPt(theLinPoint);
81  thePredState = builder(part->currentState(), paramPt);
82  // bool valid = thePredState.isValid();
83  // if (!valid) std::cout <<"Help!!!!!!!!! State is invalid\n";
84  // if (!valid) return;
85  if (std::abs(theCharge) < 1e-5) {
86  //neutral track
88  } else {
89  //charged track
91  }
92  jacobiansAvailable = true;
93 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
part
Definition: HCALResponse.h:20

◆ computeNeutralJacobians()

void ParticleKinematicLinearizedTrackState::computeNeutralJacobians ( ) const
private

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

Definition at line 283 of file ParticleKinematicLinearizedTrackState.cc.

References funct::cos(), KinematicState::globalMomentum(), KinematicState::globalPosition(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), funct::sin(), funct::tan(), theConstantTerm, theExpandedParams, theLinPoint, theMomentumJacobian, thePositionJacobian, thePredState, PerigeeKinematicState::theState(), PV3DBase< T, PVType, FrameType >::theta(), X, PV3DBase< T, PVType, FrameType >::x(), BeamSpotPI::Y, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by computeJacobians().

283  {
284  GlobalPoint paramPt(theLinPoint);
285  double thetaAtEP = thePredState.theState().globalMomentum().theta();
286  double phiAtEP = thePredState.theState().globalMomentum().phi();
287  double ptAtEP = thePredState.theState().globalMomentum().perp();
288 
289  double x_v = thePredState.theState().globalPosition().x();
290  double y_v = thePredState.theState().globalPosition().y();
291  double z_v = thePredState.theState().globalPosition().z();
292  double X = x_v - paramPt.x();
293  double Y = y_v - paramPt.y();
294 
295  // The track parameters at the expansion point
296  theExpandedParams[0] = 1 / ptAtEP;
297  theExpandedParams[1] = thetaAtEP;
298  theExpandedParams[2] = phiAtEP;
299  theExpandedParams[3] = X * sin(phiAtEP) - Y * cos(phiAtEP);
300  theExpandedParams[4] = z_v - paramPt.z() - (X * cos(phiAtEP) + Y * sin(phiAtEP)) / tan(thetaAtEP);
301  theExpandedParams[5] = part->currentState().mass();
302 
303  // The Jacobian: (all at the expansion point)
304  // [i,j]
305  // i = 0: rho = 1/pt , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
306  // j = 0: x_v, 1: y_v, 2: z_v
307  thePositionJacobian(3, 0) = sin(phiAtEP);
308  thePositionJacobian(3, 1) = -cos(phiAtEP);
309  thePositionJacobian(4, 0) = -cos(phiAtEP) / tan(thetaAtEP);
310  thePositionJacobian(4, 1) = -sin(phiAtEP) / tan(thetaAtEP);
311  thePositionJacobian(4, 2) = 1;
312 
313  // [i,j]
314  // i = 0: rho = 1/pt , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
315  // j = 0: rho = 1/pt , 1: theta, 2: phi_v
316  theMomentumJacobian(0, 0) = 1;
317  theMomentumJacobian(1, 1) = 1;
318  theMomentumJacobian(2, 2) = 1;
319 
320  theMomentumJacobian(3, 2) = X * cos(phiAtEP) + Y * sin(phiAtEP);
321  theMomentumJacobian(4, 1) = theMomentumJacobian(3, 2) * (1 + 1 / (tan(thetaAtEP) * tan(thetaAtEP)));
322 
323  theMomentumJacobian(4, 2) = (X * sin(phiAtEP) - Y * cos(phiAtEP)) / tan(thetaAtEP);
324  theMomentumJacobian(5, 3) = 1;
325 
326  // And finally the residuals:
327  AlgebraicVector3 expansionPoint;
328  expansionPoint[0] = thePredState.theState().globalPosition().x();
329  expansionPoint[1] = thePredState.theState().globalPosition().y();
330  expansionPoint[2] = thePredState.theState().globalPosition().z();
331  AlgebraicVector4 momentumAtExpansionPoint;
332  momentumAtExpansionPoint[0] = 1 / ptAtEP;
333  momentumAtExpansionPoint[1] = thetaAtEP;
334  momentumAtExpansionPoint[2] = phiAtEP;
335  momentumAtExpansionPoint[3] = theExpandedParams[5];
336 
338  theMomentumJacobian * momentumAtExpansionPoint);
339 }
T perp() const
Definition: PV3DBase.h:69
T z() const
Definition: PV3DBase.h:61
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define X(str)
Definition: MuonsGrabber.cc:38
GlobalPoint globalPosition() const
const KinematicState & theState() const
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
GlobalVector globalMomentum() const
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
ROOT::Math::SVector< double, 4 > AlgebraicVector4
part
Definition: HCALResponse.h:20
ROOT::Math::SVector< double, 3 > AlgebraicVector3
ROOT::Math::SVector< double, 6 > AlgebraicVector6
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72

◆ constantTerm()

const AlgebraicVector6 & ParticleKinematicLinearizedTrackState::constantTerm ( ) const
overridevirtual

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

Implements LinearizedTrackState< 6 >.

Definition at line 5 of file ParticleKinematicLinearizedTrackState.cc.

References computeJacobians(), jacobiansAvailable, and theConstantTerm.

Referenced by refittedParamFromEquation().

◆ createRefittedTrackState()

ParticleKinematicLinearizedTrackState::RefCountedRefittedTrackState ParticleKinematicLinearizedTrackState::createRefittedTrackState ( const GlobalPoint vertexPosition,
const AlgebraicVectorM vectorParameters,
const AlgebraicSymMatrix77 covarianceMatrix 
) const
override

Definition at line 102 of file ParticleKinematicLinearizedTrackState.cc.

References charge(), and pwdgSkimBPark_cfi::conversions.

104  {
106  KinematicState lst = conversions.kinematicState(
107  vectorParameters, vertexPosition, charge(), covarianceMatrix, part->currentState().magneticField());
109  RefCountedRefittedTrackState(new KinematicRefittedTrackState(lst, vectorParameters));
110  return rst;
111 }
ReferenceCountingPointer< RefittedTrackState< N > > RefCountedRefittedTrackState
part
Definition: HCALResponse.h:20

◆ hasError()

bool ParticleKinematicLinearizedTrackState::hasError ( void  ) const
overridevirtual

◆ linearizationPoint()

const GlobalPoint& ParticleKinematicLinearizedTrackState::linearizationPoint ( ) const
inlineoverridevirtual

The point at which the track state has been linearized

Implements LinearizedTrackState< 6 >.

Definition at line 30 of file ParticleKinematicLinearizedTrackState.h.

References theLinPoint.

◆ momentumJacobian()

const AlgebraicMatrix64 & ParticleKinematicLinearizedTrackState::momentumJacobian ( ) const
overridevirtual

◆ operator==()

bool ParticleKinematicLinearizedTrackState::operator== ( LinearizedTrackState< 6 > &  other) const
override

Definition at line 58 of file ParticleKinematicLinearizedTrackState.cc.

References trackingPlots::other, and particle().

58  {
60  dynamic_cast<const ParticleKinematicLinearizedTrackState*>(&other);
61  if (otherP == nullptr) {
62  throw VertexException(
63  " ParticleKinematicLinearizedTrackState:: don't know how to compare myself to non-kinematic track state");
64  }
65  return (*(otherP->particle()) == *part);
66 }
Common base class.
part
Definition: HCALResponse.h:20

◆ parametersFromExpansion()

const AlgebraicVector6 & ParticleKinematicLinearizedTrackState::parametersFromExpansion ( ) const
overridevirtual

◆ particle()

RefCountedKinematicParticle ParticleKinematicLinearizedTrackState::particle ( ) const

Definition at line 56 of file ParticleKinematicLinearizedTrackState.cc.

References part.

Referenced by operator==().

56 { return part; }

◆ positionJacobian()

const AlgebraicMatrix63 & ParticleKinematicLinearizedTrackState::positionJacobian ( ) const
overridevirtual

◆ predictedStateError()

AlgebraicSymMatrix66 ParticleKinematicLinearizedTrackState::predictedStateError ( ) const
overridevirtual

◆ predictedStateMomentumError()

AlgebraicSymMatrix44 ParticleKinematicLinearizedTrackState::predictedStateMomentumError ( ) const
overridevirtual

4x4 error matrix ofe xtended perigee mometum components

Implements LinearizedTrackState< 6 >.

Definition at line 124 of file ParticleKinematicLinearizedTrackState.cc.

References computeJacobians(), relativeConstraints::error, jacobiansAvailable, PerigeeKinematicState::perigeeError(), thePredState, and ExtendedPerigeeTrajectoryError::weightMatrix().

124  {
125  int error;
126  if (!jacobiansAvailable)
130  res.Place_at(m3, 0, 0);
135  return res;
136 }
const AlgebraicSymMatrix66 & weightMatrix(int &error) const
Definition: Electron.h:6
const ExtendedPerigeeTrajectoryError & perigeeError() const
ROOT::Math::SMatrix< double, 4, 4, ROOT::Math::MatRepSym< double, 4 > > AlgebraicSymMatrix44
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33

◆ predictedStateMomentumParameters()

AlgebraicVector4 ParticleKinematicLinearizedTrackState::predictedStateMomentumParameters ( ) const
overridevirtual

returns predicted 4-momentum in extended perigee parametrization

Implements LinearizedTrackState< 6 >.

Definition at line 113 of file ParticleKinematicLinearizedTrackState.cc.

References computeJacobians(), jacobiansAvailable, PerigeeKinematicState::perigeeParameters(), thePredState, and ExtendedPerigeeTrajectoryParameters::vector().

Referenced by refittedParamFromEquation().

113  {
114  if (!jacobiansAvailable)
121  return res;
122 }
const ExtendedPerigeeTrajectoryParameters & perigeeParameters() const
Definition: Electron.h:6
ROOT::Math::SVector< double, 4 > AlgebraicVector4

◆ predictedStateParameters()

AlgebraicVector6 ParticleKinematicLinearizedTrackState::predictedStateParameters ( ) const
overridevirtual

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).
extended perigee predicted parameters

Implements LinearizedTrackState< 6 >.

Definition at line 29 of file ParticleKinematicLinearizedTrackState.cc.

References computeJacobians(), jacobiansAvailable, PerigeeKinematicState::perigeeParameters(), thePredState, and ExtendedPerigeeTrajectoryParameters::vector().

◆ predictedStateWeight()

AlgebraicSymMatrix66 ParticleKinematicLinearizedTrackState::predictedStateWeight ( int &  error) const
overridevirtual

Full predicted weight matrix

Implements LinearizedTrackState< 6 >.

Definition at line 35 of file ParticleKinematicLinearizedTrackState.cc.

References computeJacobians(), relativeConstraints::error, mps_fire::i, jacobiansAvailable, PerigeeKinematicState::perigeeError(), thePredState, ExtendedPerigeeTrajectoryError::weightMatrix(), and z.

35  {
36  if (!jacobiansAvailable)
38  int i = 0;
40  error = i;
41  return z;
42  // return thePredState.perigeeError().weightMatrix(error);
43 }
const AlgebraicSymMatrix66 & weightMatrix(int &error) const
const ExtendedPerigeeTrajectoryError & perigeeError() const
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepSym< double, 6 > > AlgebraicSymMatrix66

◆ refittedParamFromEquation()

AlgebraicVector6 ParticleKinematicLinearizedTrackState::refittedParamFromEquation ( const RefCountedRefittedTrackState theRefittedState) const
override

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

Definition at line 148 of file ParticleKinematicLinearizedTrackState.cc.

References constantTerm(), M_PI, momentumJacobian(), positionJacobian(), and predictedStateMomentumParameters().

149  {
150  AlgebraicVectorM momentum = theRefittedState->momentumVector();
151  if ((momentum(2) * predictedStateMomentumParameters()(2) < 0) && (std::fabs(momentum(2)) > M_PI / 2)) {
152  if (predictedStateMomentumParameters()(2) < 0.)
153  momentum(2) -= 2 * M_PI;
154  if (predictedStateMomentumParameters()(2) > 0.)
155  momentum(2) += 2 * M_PI;
156  }
157 
158  AlgebraicVector3 vertexPosition;
159  vertexPosition(0) = theRefittedState->position().x();
160  vertexPosition(1) = theRefittedState->position().y();
161  vertexPosition(2) = theRefittedState->position().z();
162  AlgebraicVector6 param = constantTerm() + positionJacobian() * vertexPosition + momentumJacobian() * momentum;
163 
164  if (param(2) > M_PI)
165  param(2) -= 2 * M_PI;
166  if (param(2) < -M_PI)
167  param(2) += 2 * M_PI;
168 
169  return param;
170 }
const AlgebraicMatrix63 & positionJacobian() const override
const AlgebraicMatrix64 & momentumJacobian() const override
ROOT::Math::SVector< double, N - 2 > AlgebraicVectorM
AlgebraicVectorM predictedStateMomentumParameters() const override
const AlgebraicVector6 & constantTerm() const override
#define M_PI
ROOT::Math::SVector< double, 3 > AlgebraicVector3
ROOT::Math::SVector< double, 6 > AlgebraicVector6

◆ stateWithNewLinearizationPoint()

ReferenceCountingPointer< LinearizedTrackState< 6 > > ParticleKinematicLinearizedTrackState::stateWithNewLinearizationPoint ( const GlobalPoint newLP) const
overridevirtual

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< 6 >.

Definition at line 96 of file ParticleKinematicLinearizedTrackState.cc.

References part, and ParticleKinematicLinearizedTrackState().

◆ track()

reco::TransientTrack ParticleKinematicLinearizedTrackState::track ( ) const
overridevirtual

Implements LinearizedTrackState< 6 >.

Definition at line 341 of file ParticleKinematicLinearizedTrackState.cc.

341  {
342  throw VertexException(" ParticleKinematicLinearizedTrackState:: no TransientTrack to return");
343 }
Common base class.

◆ weightInMixture()

double ParticleKinematicLinearizedTrackState::weightInMixture ( ) const
overridevirtual

Implements LinearizedTrackState< 6 >.

Definition at line 138 of file ParticleKinematicLinearizedTrackState.cc.

138 { return 1.; }

Friends And Related Function Documentation

◆ ParticleKinematicLinearizedTrackStateFactory

Definition at line 15 of file ParticleKinematicLinearizedTrackState.h.

Member Data Documentation

◆ builder

TransientTrackKinematicStateBuilder ParticleKinematicLinearizedTrackState::builder
private

Definition at line 140 of file ParticleKinematicLinearizedTrackState.h.

Referenced by computeJacobians().

◆ errorAvailable

bool ParticleKinematicLinearizedTrackState::errorAvailable
mutableprivate

Definition at line 142 of file ParticleKinematicLinearizedTrackState.h.

◆ impactPointAvailable

bool ParticleKinematicLinearizedTrackState::impactPointAvailable
mutableprivate

Definition at line 152 of file ParticleKinematicLinearizedTrackState.h.

◆ jacobiansAvailable

bool ParticleKinematicLinearizedTrackState::jacobiansAvailable
mutableprivate

◆ part

RefCountedKinematicParticle ParticleKinematicLinearizedTrackState::part
private

◆ theCharge

TrackCharge ParticleKinematicLinearizedTrackState::theCharge
private

Definition at line 150 of file ParticleKinematicLinearizedTrackState.h.

Referenced by computeJacobians().

◆ theConstantTerm

AlgebraicVector6 ParticleKinematicLinearizedTrackState::theConstantTerm
mutableprivate

◆ theExpandedParams

AlgebraicVector6 ParticleKinematicLinearizedTrackState::theExpandedParams
mutableprivate

◆ theLinPoint

GlobalPoint ParticleKinematicLinearizedTrackState::theLinPoint
private

◆ theMomentumJacobian

AlgebraicMatrix64 ParticleKinematicLinearizedTrackState::theMomentumJacobian
mutableprivate

◆ thePositionJacobian

AlgebraicMatrix63 ParticleKinematicLinearizedTrackState::thePositionJacobian
mutableprivate

◆ thePredState

PerigeeKinematicState ParticleKinematicLinearizedTrackState::thePredState
mutableprivate