File indexing completed on 2024-07-24 04:45:12
0001 #ifndef PerigeeLinearizedTrackState_H
0002 #define PerigeeLinearizedTrackState_H
0003
0004 #include "RecoVertex/VertexPrimitives/interface/LinearizedTrackState.h"
0005 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
0006 #include "DataFormats/GeometrySurface/interface/ReferenceCounted.h"
0007 #include "RecoVertex/VertexPrimitives/interface/LinearizedTrackState.h"
0008 #include "Math/SMatrix.h"
0009 #include "DataFormats/CLHEP/interface/Migration.h"
0010 #include "FWCore/Utilities/interface/Likely.h"
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033 class PerigeeLinearizedTrackState final : public LinearizedTrackState<5> {
0034 public:
0035
0036
0037
0038 friend class LinearizedTrackStateFactory;
0039 typedef ReferenceCountingPointer<LinearizedTrackState<5> > RefCountedLinearizedTrackState;
0040
0041
0042
0043
0044
0045 RefCountedLinearizedTrackState stateWithNewLinearizationPoint(const GlobalPoint& newLP) const override;
0046
0047
0048
0049
0050 const GlobalPoint& linearizationPoint() const override { return theLinPoint; }
0051
0052 reco::TransientTrack track() const override { return theTrack; }
0053
0054 const TrajectoryStateOnSurface state() const { return theTSOS; }
0055
0056
0057
0058
0059 const AlgebraicVector5& constantTerm() const override;
0060
0061
0062
0063
0064 const AlgebraicMatrix53& positionJacobian() const override;
0065
0066
0067
0068
0069 const AlgebraicMatrix53& momentumJacobian() const override;
0070
0071
0072
0073 const AlgebraicVector5& parametersFromExpansion() const override;
0074
0075
0076
0077
0078
0079 const TrajectoryStateClosestToPoint& predictedState() const;
0080
0081
0082
0083
0084 AlgebraicVector5 predictedStateParameters() const override;
0085
0086
0087
0088
0089 AlgebraicVector3 predictedStateMomentumParameters() const override;
0090
0091
0092
0093
0094
0095 AlgebraicSymMatrix55 predictedStateWeight(int& error) const override;
0096
0097
0098
0099
0100 AlgebraicSymMatrix55 predictedStateError() const override;
0101
0102
0103
0104
0105 AlgebraicSymMatrix33 predictedStateMomentumError() const override;
0106
0107
0108
0109
0110
0111 TrackCharge charge() const override { return theCharge; }
0112
0113 bool hasError() const override;
0114
0115 bool operator==(const LinearizedTrackState<5>& other) const override;
0116
0117 bool operator==(const ReferenceCountingPointer<LinearizedTrackState<5> >& other) const;
0118
0119
0120
0121
0122 RefCountedRefittedTrackState createRefittedTrackState(const GlobalPoint& vertexPosition,
0123 const AlgebraicVector3& vectorParameters,
0124 const AlgebraicSymMatrix66& covarianceMatrix) const override;
0125
0126 AlgebraicVector5 refittedParamFromEquation(const RefCountedRefittedTrackState& theRefittedState) const override;
0127
0128 double weightInMixture() const override { return theTSOS.weight(); }
0129
0130 void inline checkParameters(AlgebraicVector5& parameters) const override;
0131
0132 std::vector<ReferenceCountingPointer<LinearizedTrackState<5> > > components() const override;
0133
0134 bool isValid() const override;
0135
0136 private:
0137
0138
0139
0140 PerigeeLinearizedTrackState(const GlobalPoint& linP,
0141 const reco::TransientTrack& track,
0142 const TrajectoryStateOnSurface& tsos)
0143 : theLinPoint(linP), theTrack(track), theTSOS(tsos), theCharge(theTrack.charge()), jacobiansAvailable(false) {}
0144
0145
0146
0147 void computeJacobians() const;
0148
0149
0150 void computeChargedJacobians() const;
0151
0152
0153 void computeNeutralJacobians() const;
0154
0155
0156 void compute3DImpactPoint() const;
0157
0158 GlobalPoint theLinPoint;
0159 reco::TransientTrack theTrack;
0160 const TrajectoryStateOnSurface theTSOS;
0161
0162 mutable AlgebraicMatrix53 thePositionJacobian, theMomentumJacobian;
0163 mutable TrajectoryStateClosestToPoint thePredState;
0164 mutable AlgebraicVector5 theConstantTerm;
0165 mutable AlgebraicVector5 theExpandedParams;
0166
0167 TSCPBuilderNoMaterial builder;
0168
0169 TrackCharge theCharge;
0170 mutable bool jacobiansAvailable;
0171 };
0172
0173
0174
0175
0176 inline const AlgebraicVector5& PerigeeLinearizedTrackState::constantTerm() const {
0177 if UNLIKELY (!jacobiansAvailable)
0178 computeJacobians();
0179 return theConstantTerm;
0180 }
0181
0182
0183
0184
0185 inline const AlgebraicMatrix53& PerigeeLinearizedTrackState::positionJacobian() const {
0186 if UNLIKELY (!jacobiansAvailable)
0187 computeJacobians();
0188 return thePositionJacobian;
0189 }
0190
0191
0192
0193
0194 inline const AlgebraicMatrix53& PerigeeLinearizedTrackState::momentumJacobian() const {
0195 if UNLIKELY (!jacobiansAvailable)
0196 computeJacobians();
0197 return theMomentumJacobian;
0198 }
0199
0200
0201
0202 inline const AlgebraicVector5& PerigeeLinearizedTrackState::parametersFromExpansion() const {
0203 if UNLIKELY (!jacobiansAvailable)
0204 computeJacobians();
0205 return theExpandedParams;
0206 }
0207
0208
0209
0210
0211
0212 inline const TrajectoryStateClosestToPoint& PerigeeLinearizedTrackState::predictedState() const {
0213 if UNLIKELY (!jacobiansAvailable)
0214 computeJacobians();
0215 return thePredState;
0216 }
0217
0218 inline bool PerigeeLinearizedTrackState::hasError() const {
0219 if UNLIKELY (!jacobiansAvailable)
0220 computeJacobians();
0221 return thePredState.hasError();
0222 }
0223
0224 inline AlgebraicVector5 PerigeeLinearizedTrackState::predictedStateParameters() const {
0225 if (!jacobiansAvailable)
0226 computeJacobians();
0227 return thePredState.perigeeParameters().vector();
0228 }
0229
0230 inline AlgebraicVector3 PerigeeLinearizedTrackState::predictedStateMomentumParameters() const {
0231 if UNLIKELY (!jacobiansAvailable)
0232 computeJacobians();
0233 auto v = thePredState.perigeeParameters().vector();
0234 return AlgebraicVector3(v[0], v[1], v[2]);
0235 }
0236
0237 inline AlgebraicSymMatrix55 PerigeeLinearizedTrackState::predictedStateWeight(int& error) const {
0238 if UNLIKELY (!jacobiansAvailable)
0239 computeJacobians();
0240 if (!thePredState.isValid()) {
0241 error = 1;
0242 return AlgebraicSymMatrix55();
0243 }
0244 return thePredState.perigeeError().weightMatrix(error);
0245 }
0246
0247 inline AlgebraicSymMatrix55 PerigeeLinearizedTrackState::predictedStateError() const {
0248 if UNLIKELY (!jacobiansAvailable)
0249 computeJacobians();
0250 return thePredState.perigeeError().covarianceMatrix();
0251 }
0252
0253 inline AlgebraicSymMatrix33 PerigeeLinearizedTrackState::predictedStateMomentumError() const {
0254 if UNLIKELY (!jacobiansAvailable)
0255 computeJacobians();
0256 return thePredState.perigeeError().covarianceMatrix().Sub<AlgebraicSymMatrix33>(0, 2);
0257 }
0258
0259 inline bool PerigeeLinearizedTrackState::isValid() const {
0260 if UNLIKELY (!theTSOS.isValid())
0261 return false;
0262 if UNLIKELY (!jacobiansAvailable)
0263 computeJacobians();
0264 return jacobiansAvailable;
0265 }
0266
0267 #endif