ATLAS Offline Software
Loading...
Searching...
No Matches
KalmanUpdatorSMatrix.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6// KalmanUpdatorSMatrix.cxx
7// Source file for class KalmanUpdatorSMatrix
8// uses SMatrix math library internally.
9// The main updating calculations can be found in
10// calculateFilterStep_1D() to _5D().
12// (c) ATLAS Detector software
14// Wolfgang Liebig <http://consult.cern.ch/xwho/people/54608>
16
18
20
24
25#include <cstring>
26
27// using namespace ROOT::Math;
28// using ROOT::Math::SMatrix;
29
31
32// constructor
33Trk::KalmanUpdatorSMatrix::KalmanUpdatorSMatrix(const std::string& t,const std::string& n,const IInterface* p) :
34 AthAlgTool (t,n,p),
35 m_cov_stdvec{250.,250.,0.25,0.25,0.000001}, // set defaults _before_ reading from job options
38{
39 // AlgTool stuff
40 declareProperty("InitialCovariances",m_cov_stdvec,"default covariance to be used at start of filter");
41 declareProperty("FastTrackStateCovCalculation",m_useFruehwirth8a=false,"toggles which formula to use for updated cov");
42 declareInterface<IUpdator>( this );
43}
44
45// destructor
47= default;
48
49// initialize
51{
52 if (m_cov_stdvec.size() < 5) {
53 ATH_MSG_WARNING( "Wrong-sized initial covariance given, so set to default: ");
54 m_cov_stdvec = {250.,250.,0.25, 0.25, 0.000001};
55 }
58 ATH_MSG_DEBUG( "Fast computation will be used for track state cov matrices (Fruehwirth-1987 eq. 8a)." );
59 } else{
60 ATH_MSG_DEBUG( "Track state cov matrix will be calculated according to Gelb-1975 p305." );
61 }
62
63 const SParVector5 IV(1.0, 1.0, 1.0, 1.0, 1.0);
64 m_unitMatrix.SetDiagonal(IV);
65
66 return StatusCode::SUCCESS;
67}
68
69// finalize
71{
72 return StatusCode::SUCCESS;
73}
74
75// updator #1 for Kalman Fitter - version with Amg::Vector2D (for example for PrepRawData)
76std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::addToState (const Trk::TrackParameters& inputTrkPar,
77 const Amg::Vector2D& measLocPos,
78 const Amg::MatrixX& measLocCov) const {
79 if (msgLvl(MSG::VERBOSE)) logStart("addToState(TP,LPOS,ERR)",inputTrkPar);
80 FitQualityOnSurface* fitQoS = nullptr;
81 const int updatingSign = 1;
82
83 SCovMatrix5 covTrk;
84 if (!getStartCov(covTrk,inputTrkPar,updatingSign)) return nullptr;
85
86
87 int nLocCoord = measLocCov.cols();
88 if (nLocCoord == 1) {
89 return calculateFilterStep_1D (inputTrkPar, // transmit the surface
90 SParVector5(&inputTrkPar.parameters()[0],5),covTrk,
91 measLocPos[0],1,measLocCov,
92 updatingSign,fitQoS,false);
93 } if (nLocCoord == 2) {
94 return calculateFilterStep_2D (inputTrkPar,// transmit the surface
95 SParVector5(&inputTrkPar.parameters()[0],5),covTrk,
96 SParVector2(measLocPos[0],measLocPos[1]),3,
97 measLocCov,
98 updatingSign,fitQoS,false);
99 }
100 ATH_MSG_WARNING( " number (" << nLocCoord << ") of local coordinates must be 1 or 2, can not update!" );
101 return nullptr;
102
103}
104
105// updator #2 for Kalman Fitter - version with LocalParameters (for example for RIO_OnTrack)
106std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::addToState (const Trk::TrackParameters& trkPar,
107 const LocalParameters& measmtPar,
108 const Amg::MatrixX& measmtCov) const {
109 if (msgLvl(MSG::VERBOSE)) logStart("addToState(TP,LPAR,ERR)",trkPar);
110 FitQualityOnSurface* fitQoS = nullptr;
111 return prepareFilterStep (trkPar, measmtPar, measmtCov, 1, fitQoS, false);
112}
113
114// updator #3 for Kalman Fitter - version with Amg::Vector2D (for example for PrepRawData)
115std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::addToState (const Trk::TrackParameters& inputTP,
116 const Amg::Vector2D& measLocPos,
117 const Amg::MatrixX& measLocCov,
118 FitQualityOnSurface*& fitQoS) const {
119 const int updatingSign = 1;
120 if (msgLvl(MSG::VERBOSE)) logStart("addToState(TP,LPOS,ERR,FQ)",inputTP);
121 if (fitQoS) {
122 ATH_MSG_WARNING( "expect nil FitQuality pointer, refuse operation to avoid mem leak!");
123 return nullptr;
124 }
125
126 SCovMatrix5 covTrk;
127 if (!getStartCov(covTrk,inputTP,updatingSign)) return nullptr;
128
129 int nLocCoord = measLocCov.cols();
130 if (nLocCoord == 1) {
131 return calculateFilterStep_1D (inputTP,
132 SParVector5(&inputTP.parameters()[0],5),covTrk,
133 measLocPos[0],1,measLocCov,
134 updatingSign,fitQoS,true);
135 } if (nLocCoord == 2) {
136 return calculateFilterStep_2D (inputTP,
137 SParVector5(&inputTP.parameters()[0],5),covTrk,
138 SParVector2(measLocPos[0],measLocPos[1]),3,
139 measLocCov,
140 updatingSign,fitQoS,true);
141 }
142 ATH_MSG_WARNING( " number (" << nLocCoord << ") of local coordinates must be 1 or 2, can not update!" );
143 return nullptr;
144
145
146}
147
148// updator #4 for Kalman Fitter - version with LocalParameters (for example for RIO_OnTrack)
149std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::addToState (const Trk::TrackParameters& trkPar,
150 const LocalParameters& measmtPar,
151 const Amg::MatrixX& measmtCov,
152 FitQualityOnSurface*& fitQoS) const {
153 if (msgLvl(MSG::VERBOSE)) logStart("addToState(TP,LPAR,ERR,FQ)",trkPar);
154 if (fitQoS) {
155 ATH_MSG_WARNING( "expect nil FitQuality pointer, refuse operation to"
156 << " avoid mem leak!" );
157 return nullptr;
158 }
159 return prepareFilterStep (trkPar, measmtPar, measmtCov, 1, fitQoS, true);
160
161}
162
163// inverse updator #1 for Kalman Fitter - version with Amg::Vector2D (for example for PrepRawData)
164std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::removeFromState (const Trk::TrackParameters& inputTP,
165 const Amg::Vector2D& measLocPos,
166 const Amg::MatrixX& measLocCov) const {
167 if (msgLvl(MSG::VERBOSE)) logStart("removeFromState(TP,LPOS,ERR)",inputTP);
168 FitQualityOnSurface* fitQoS = nullptr;
169 const int updatingSign = -1;
170 SCovMatrix5 covTrk;
171 if (!getStartCov(covTrk,inputTP,updatingSign)) return nullptr;
172
173 int nLocCoord = measLocCov.cols();
174 if (nLocCoord == 1) {
175 return calculateFilterStep_1D (inputTP,
176 SParVector5(&inputTP.parameters()[0],5),covTrk,
177 measLocPos[0],1,measLocCov,
178 updatingSign,fitQoS,false);
179 } if (nLocCoord == 2) {
180 return calculateFilterStep_2D (inputTP,
181 SParVector5(&inputTP.parameters()[0],5),covTrk,
182 SParVector2(measLocPos[0],measLocPos[1]),3,
183 measLocCov,
184 updatingSign,fitQoS,false);
185 }
186 ATH_MSG_WARNING( " number (" << nLocCoord << ") of local coordinates "
187 << "must be 1 or 2, can not un-update!" );
188 return nullptr;
189
190}
191
192// inverse updator #2 for Kalman Fitter - version with LocalParameters (for example for RIO_OnTrack)
193std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::removeFromState (const Trk::TrackParameters& trkPar,
194 const LocalParameters& measmtPar,
195 const Amg::MatrixX& measmtCov) const {
196 if (msgLvl(MSG::DEBUG)) logStart("removeFromState(TP,LPAR,ERR)",trkPar);
197 FitQualityOnSurface* fitQoS = nullptr;
198 return prepareFilterStep (trkPar, measmtPar, measmtCov,-1,fitQoS, false);
199}
200
201// inverse updator #3 for Kalman Fitter - version with Amg::Vector2D (for example for PrepRawData)
202std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::removeFromState (const Trk::TrackParameters& inputTP,
203 const Amg::Vector2D& measLocPos,
204 const Amg::MatrixX& measLocCov,
205 FitQualityOnSurface*& fitQoS) const {
206 const int updatingSign = -1;
207 if (msgLvl(MSG::VERBOSE)) logStart("removeFromState(TP,LPOS,ERR,FQ)",inputTP);
208 if (fitQoS) {
209 msg(MSG::WARNING) << "expect nil FitQuality pointer, refuse operation to"
210 << " avoid mem leak!" << endmsg;
211 return nullptr;
212 }
213 SCovMatrix5 covTrk;
214
215 if (!getStartCov(covTrk,inputTP,updatingSign)) return nullptr;
216
217 int nLocCoord = measLocCov.cols();
218 if (nLocCoord == 1) {
219 return calculateFilterStep_1D (inputTP,
220 SParVector5(&inputTP.parameters()[0],5),covTrk,
221 measLocPos[0],1,measLocCov,
222 updatingSign,fitQoS,true);
223 } if (nLocCoord == 2) {
224 return calculateFilterStep_2D (inputTP,
225 SParVector5(&inputTP.parameters()[0],5),covTrk,
226 SParVector2(measLocPos[0],measLocPos[1]),3,
227 measLocCov,
228 updatingSign,fitQoS,true);
229 }
230 ATH_MSG_WARNING( " number (" << nLocCoord << ") of local coordinates"
231 << " must be 1 or 2, can not un-update!" );
232 return nullptr;
233
234
235}
236
237// inverse updator #4 for Kalman Fitter - version with LocalParameters (for example for RIO_OnTrack)
238std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::removeFromState (const Trk::TrackParameters& trkPar,
239 const LocalParameters& measmtPar,
240 const Amg::MatrixX& measmtCov,
241 FitQualityOnSurface*& fitQoS) const {
242 if (msgLvl(MSG::VERBOSE)) logStart("removeFromState(TP,LPAR,ERR,FQ)",trkPar);
243 if (fitQoS) {
244 ATH_MSG_WARNING( "expect nil FitQuality pointer, refuse operation to avoid mem leak!" );
245 return nullptr;
246 }
247 return prepareFilterStep (trkPar, measmtPar, measmtCov, -1, fitQoS, true);
248
249}
250
251// state-to-state updator, trajectory combination - version without fitQuality
252std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::combineStates (const Trk::TrackParameters& one,
253 const Trk::TrackParameters& two) const {
254 // try if both Track Parameters are measured ones ?
255 const AmgSymMatrix(5)* covOne = one.covariance();
256 const AmgSymMatrix(5)* covTwo = two.covariance();
257
258 // remember, either one OR two might have no error, but not both !
259 if (!covOne && ! covTwo) {
260 ATH_MSG_WARNING( "both parameters have no errors, invalid use of Updator::combineStates()" );
261 return nullptr;
262 }
263 // if only one of two has an error, return that one
264 if (!covOne) {
265 if (msgLvl(MSG::VERBOSE)) logResult("combineStates(TP,TP)",two.parameters(),*covTwo);
266 return std::unique_ptr<Trk::TrackParameters>(two.clone());
267 }
268 if (!covTwo) {
269 if (msgLvl(MSG::VERBOSE)) logResult("combineStates(TP,TP)",one.parameters(),*covOne);
270 return std::unique_ptr<Trk::TrackParameters>(one.clone());
271 }
272
273 // ... FIXME - TRT is so difficult, need to check that both parameters are in the same frame
274 // otherwise go into frame of one !
275
276 // call 5dim updator to combine using gain matrix formalism
277 SCovMatrix5 covOneSM;
278 for (int i=0; i<5; ++i)
279 for (int j=0; j<=i; ++j) {
280 covOneSM(j,i) = (*covOne)(j,i);
281 }
282
283 FitQualityOnSurface* fitQoS = nullptr;
284 return calculateFilterStep_5D(one,SParVector5(&one.parameters()[0],5),covOneSM,
285 SParVector5(&two.parameters()[0],5),
286 *covTwo,
287 +1,fitQoS,false);
288}
289
290// state-to-state updator, trajectory combination - version with fitQuality
291std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::combineStates (const Trk::TrackParameters& one,
293 FitQualityOnSurface*& fitQoS) const {
294 // try if both Track Parameters are measured ones ?
295 const AmgSymMatrix(5)* covOne = one.covariance();
296 const AmgSymMatrix(5)* covTwo = two.covariance();
297
298 // remember, either one OR two might have no error, but not both !
299 if (!covOne && ! covTwo) {
300 ATH_MSG_WARNING( "both parameters have no errors, invalid use of Updator::combineStates()" );
301 return nullptr;
302 }
303
304 if (fitQoS) {
305 ATH_MSG_WARNING( "expect nil FitQuality pointer, refuse operation to avoid mem leak!");
306 return nullptr;
307 }
308
309 // if only one of two has an error, return that one
310 if (!covOne) {
311 if (msgLvl(MSG::VERBOSE)) logResult("combineStates(TP,TP)",two.parameters(),*covTwo);
312 return std::unique_ptr<Trk::TrackParameters>(two.clone());
313 }
314 if (!covTwo) {
315 if (msgLvl(MSG::VERBOSE)) logResult("combineStates(TP,TP)",one.parameters(),*covOne);
316 return std::unique_ptr<Trk::TrackParameters>(one.clone());
317 }
318
319 // call 5dim updator to combine using gain matrix formalism
320 SCovMatrix5 covOneSM;
321 for (int i=0; i<5; ++i)
322 for (int j=0; j<=i; ++j) {
323 covOneSM(j,i) = (*covOne)(j,i);
324 }
325
326 return calculateFilterStep_5D(one,SParVector5(&one.parameters()[0],5),covOneSM,
327 SParVector5(&two.parameters()[0],5),
328 (*covTwo),
329 +1,fitQoS,true);
330}
331
332
333// estimator for FitQuality on Surface (allows for setting of LR for straws)
336 const Amg::Vector2D& locPos,
337 const Amg::MatrixX& covRio) const {
338 ATH_MSG_DEBUG( "--> entered KalmanUpdatorSMatrix::fullStateFitQuality(TP,LPOS,ERR)" );
339 // try if Track Parameters are measured ones ?
340 if (!trkPar.covariance()) {
341 ATH_MSG_ERROR( "updated track state has no error matrix" );
342 return {};
343 }
344 // For the LocalPos. version, need to get # meas. coord. from covariance matrix.
345 int nLocCoord = covRio.cols();
346 if (nLocCoord == 1) {
347 return makeChi2_1D(SParVector5(&trkPar.parameters()[0],5),
348 (*trkPar.covariance()),
349 locPos[Trk::locX],covRio(0,0),
350 1,-1); // key=1, signForUpdatedChi2 = -1
351 } if (nLocCoord == 2) {
352 SCovMatrix2 SmeasCov;
353 SmeasCov(0,0) = covRio(0,0);
354 SmeasCov(1,0) = covRio(0,1);
355 SmeasCov(1,1) = covRio(1,1);
356 return makeChi2_2D(SParVector5(&trkPar.parameters()[0],5),
357 (*trkPar.covariance()),
358 SParVector2(locPos[Trk::locX],locPos[Trk::locY]),
359 SmeasCov, 2,-1);
360 }
361 ATH_MSG_WARNING( "Error in local position - must be 1D or 2D!" );
362 return {};
363
364}
365
366
367// estimator for FitQuality on Surface (allows for setting of LR for straws)
370 const Trk::LocalParameters& parRio,
371 const Amg::MatrixX& covRio) const {
372 ATH_MSG_VERBOSE( "--> entered KalmanUpdatorSMatrix::fullStateFitQuality(TP,LPAR,ERR)" );
373
374 // try if Track Parameters are measured ones ?
375 if (!trkPar.covariance()) {
376 ATH_MSG_ERROR( "updated track state has no error matrix" );
377 return {};
378 }
379 int nLocCoord = parRio.dimension();
380 if ( ! consistentParamDimensions(parRio,covRio.cols()) ) return {};
381 // local params can NOT be accessed like vector[i], therefore need some acrobatics:
382 ROOT::Math::SVector<int,5> intAccessor;
383 for (int i=0,k=0; i<5; ++i) { if (parRio.parameterKey() & (1<<i)) intAccessor(k++)=i; }
384 if (nLocCoord == 1) {
385 return makeChi2_1D(SParVector5(&trkPar.parameters()[0],5),
386 (*trkPar.covariance()),
387 parRio[Trk::ParamDefsAccessor::pardef[intAccessor(0)]],covRio(0,0),
388 parRio.parameterKey(),-1);
389 } if (nLocCoord == 2) {
390 SCovMatrix2 SmeasCov;
391 for (int i=0, irow=0; i<5; ++i) {
392 if (parRio.parameterKey() & (1<<i)) {
393 SmeasCov(0,irow) = covRio(0,irow);
394 SmeasCov(1,irow) = covRio(1,irow);
395 ++irow;
396 }
397 }
398 return makeChi2_2D(SParVector5(&trkPar.parameters()[0],5),
399 (*trkPar.covariance()),
400 SParVector2(parRio[Trk::ParamDefsAccessor::pardef[intAccessor(0)]],
401 parRio[Trk::ParamDefsAccessor::pardef[intAccessor(1)]]),
402 SmeasCov, parRio.parameterKey(),-1);
403 } if (nLocCoord == 5) {
404 return makeChi2_5D(SParVector5(&trkPar.parameters()[0],5),
405 (*trkPar.covariance()),
406 SParVector5(parRio[Trk::locX],parRio[Trk::locY],
407 parRio[Trk::phi0],parRio[Trk::theta],
408 parRio[Trk::qOverP]),
409 covRio,-1);
410 } // stay with Eigen for other cases
411
412 // State to measurement dimensional reduction Matrix ( n x m )
413 const Amg::MatrixX& H = parRio.expansionMatrix();
414
415 // residuals
417 if (parRio.parameterKey()==31) r = (parRio - trkPar.parameters());
418 else r = (parRio - H*trkPar.parameters());
419
420 // all the rest is outsourced to a common chi2 routine
421 return makeChi2Object(r,(*trkPar.covariance()),covRio,H,-1);
422
423}
424
425// estimator for FitQuality on Surface (allows for setting of LR for straws)
428 const Amg::Vector2D& rioLocPos,
429 const Amg::MatrixX& covRio) const {
430 ATH_MSG_VERBOSE( "--> entered KalmanUpdatorSMatrix::predictedStateFitQuality(TP,LPOS,ERR)" );
431 // try if Track Parameters are measured ones ?
432 if (!predPar.covariance()) {
433 ATH_MSG_WARNING( "input state has no error matrix in predictedStateFitQuality()" );
434 return {};
435 }
436 // For the LocalPos. version, need to get # meas. coord. from covariance matrix.
437 int nLocCoord = covRio.cols();
438 if (nLocCoord == 1) {
439 return makeChi2_1D(SParVector5(&predPar.parameters()[0],5),
440 (*predPar.covariance()),
441 rioLocPos[Trk::locX],covRio(0,0),
442 1,+1); // key=1, signForUpdatedChi2 = +1
443 } if (nLocCoord == 2) {
444 SCovMatrix2 SmeasCov;
445 SmeasCov(0,0) = covRio(0,0);
446 SmeasCov(1,0) = covRio(1,0);
447 SmeasCov(1,1) = covRio(1,1);
448 return makeChi2_2D(SParVector5(&predPar.parameters()[0],5),
449 (*predPar.covariance()),
450 SParVector2(rioLocPos[Trk::locX],rioLocPos[Trk::locY]),
451 SmeasCov, 2,+1);
452 }
453 ATH_MSG_WARNING( "Error in local position - must be 1D or 2D!" );
454 return {};
455
456}
457
458// estimator for FitQuality on Surface (allows for setting of LR for straws)
461 const Trk::LocalParameters& parRio,
462 const Amg::MatrixX& covRio) const {
463 ATH_MSG_VERBOSE( "--> entered KalmanUpdatorSMatrix::predictedStateFitQuality(TP,LPAR,ERR)" );
464
465 // try if Track Parameters are measured ones ?
466 if (!predPar.covariance()) {
467 ATH_MSG_WARNING( "input state has no error matrix in predictedStateFitQuality()" );
468 return {};
469 }
470 int nLocCoord = parRio.dimension();
471 if ( ! consistentParamDimensions(parRio,covRio.cols()) ) return {};
472
473 ROOT::Math::SVector<int,5> intAccessor;
474 for (int i=0,k=0; i<5; ++i) { if (parRio.parameterKey() & (1<<i)) intAccessor(k++)=i; }
475 if (nLocCoord == 1) {
476 return makeChi2_1D(SParVector5(&predPar.parameters()[0],5),
477 (*predPar.covariance()),
478 parRio[Trk::ParamDefsAccessor::pardef[intAccessor(0)]],covRio(0,0),
479 parRio.parameterKey(),+1);
480 } if (nLocCoord == 2) {
481 SCovMatrix2 SmeasCov;
482 for (int i=0, irow=0; i<5; ++i) {
483 if (parRio.parameterKey() & (1<<i)) {
484 SmeasCov(0,irow) = covRio(0,irow);
485 SmeasCov(1,irow) = covRio(1,irow);
486 ++irow;
487 }
488 }
489 return makeChi2_2D(SParVector5(&predPar.parameters()[0],5),
490 (*predPar.covariance()),
491 SParVector2(parRio[Trk::ParamDefsAccessor::pardef[intAccessor(0)]],
492 parRio[Trk::ParamDefsAccessor::pardef[intAccessor(1)]]),
493 SmeasCov, parRio.parameterKey(),+1);
494 } if (nLocCoord == 5 ) {
495 return makeChi2_5D(SParVector5(&predPar.parameters()[0],5),
496 (*predPar.covariance()),
497 SParVector5(parRio[Trk::locX],parRio[Trk::locY],
498 parRio[Trk::phi0],parRio[Trk::theta],
499 parRio[Trk::qOverP]),
500 covRio,+1);
501 } // stay with CLHEP for other cases
502
503 // State to measurement dimensional reduction Matrix ( n x m )
504 const Amg::MatrixX& H = parRio.expansionMatrix();
505
506 // residuals
508 if (parRio.parameterKey()==31) r = parRio - predPar.parameters();
509 else r = parRio - H * predPar.parameters();
510
511 // all the rest is outsourced to a common chi2 routine
512 return makeChi2Object(r,(*predPar.covariance()),
513 covRio,H,+1);
514
515}
516
517// estimator for FitQuality on Surface (allows for setting of LR for straws)
520 const Trk::TrackParameters& trkParTwo) const {
521 ATH_MSG_VERBOSE("--> entered KalmanUpdatorSMatrix::predictedStateFitQuality(TP,TP)");
522 // try if both Track Parameters are measured ones ?
523 const AmgSymMatrix(5)* covOne = trkParOne.covariance();
524 const AmgSymMatrix(5)* covTwo = trkParTwo.covariance();
525 // remember, either one OR two might have no error, but not both !
526 if (!covOne && ! covTwo) {
527 ATH_MSG_WARNING( "both parameters have no errors, invalid use of Updator::fitQuality()" );
528 return {};
529 }
530 // if only one of two has an error, place a message.
531 if (!covOne || ! covTwo) {
532 ATH_MSG_DEBUG( "One parameter does not have uncertainties, assume initial state and return chi2=0.0");
533 return {0.f, 5};
534 }
535 return makeChi2_5D(SParVector5(&trkParOne.parameters()[0],5),
536 *covOne,
537 SParVector5(&trkParTwo.parameters()[0],5),
538 *covTwo, +1);
539}
540
541std::vector<double> Trk::KalmanUpdatorSMatrix::initialErrors() const {
542 std::vector<double> E(5);
543 for (int i=0; i<5; ++i) E[i] = std::sqrt(m_cov0(i));
544 return E;
545}
546
547// analyse dimension of localParameters to call appropriate fast-access mathematics
548std::unique_ptr<Trk::TrackParameters> Trk::KalmanUpdatorSMatrix::prepareFilterStep (const Trk::TrackParameters& inputTrkPar,
549 const Trk::LocalParameters& parRio,
550 const Amg::MatrixX& covRio,
551 const int sign,
553 bool createFQoS ) const {
554
555 // try if Track Parameters are measured ones ?
556 SCovMatrix5 covTrk;
557 if (!getStartCov(covTrk,inputTrkPar,sign)) return nullptr;
558
559 int nLocCoord = covRio.cols();
560 if ( ! consistentParamDimensions(parRio,nLocCoord) ) return nullptr;
561 if (msgLvl(MSG::VERBOSE)) logInputCov(covTrk,parRio,covRio);
562
563 // local params can NOT be accessed like vector[i], therefore need some acrobatics:
564 ROOT::Math::SVector<int,5> intAccessor;
565 for (int i=0,k=0; i<5; ++i) { if (parRio.parameterKey() & (1<<i)) intAccessor(k++)=i; }
566
567 if (nLocCoord==1) return calculateFilterStep_1D (inputTrkPar,SParVector5(&inputTrkPar.parameters()[0],5),covTrk,
568 parRio[Trk::ParamDefsAccessor::pardef[intAccessor(0)]],
569 parRio.parameterKey(),covRio,
570 sign,fitQoS,createFQoS);
571 if (nLocCoord==2) return calculateFilterStep_2D (inputTrkPar,SParVector5(&inputTrkPar.parameters()[0],5),covTrk,
572 SParVector2(parRio[Trk::ParamDefsAccessor::pardef[intAccessor(0)]],parRio[Trk::ParamDefsAccessor::pardef[intAccessor(1)]]),
573 parRio.parameterKey(),covRio,
574 sign,fitQoS,createFQoS);
575
576 if (nLocCoord==3) return calculateFilterStep_3D (inputTrkPar,SParVector5(&inputTrkPar.parameters()[0],5),covTrk,
577 parRio,covRio,
578 sign,fitQoS,createFQoS);
579 if (nLocCoord==4) return calculateFilterStep_4D (inputTrkPar,SParVector5(&inputTrkPar.parameters()[0],5),covTrk,
580 parRio,covRio,
581 sign,fitQoS,createFQoS);
582 if (nLocCoord==5) return calculateFilterStep_5D (inputTrkPar,SParVector5(&inputTrkPar.parameters()[0],5),covTrk,
583 SParVector5(&(parRio[Trk::loc1]),5),covRio,
584 sign,fitQoS,createFQoS);
585 return nullptr;
586}
587
588// calculations for Kalman updator and inverse Kalman filter
589std::unique_ptr<Trk::TrackParameters>
591 const SParVector5& trkPar,
592 const SCovMatrix5& trkCov,
593 double measPar,
594 int paramKey,
595 const Amg::MatrixX& measCov,
596 const int sign,
598 bool createFQoS ) const {
599 // use measuring coordinate (variable "mk") instead of reduction matrix
600 int mk=0;
601 if (paramKey!=1) for (int i=0; i<5; ++i) if (paramKey & (1<<i)) mk=i;
602
603 double r = measPar - trkPar(mk);
604 double R = (sign * measCov(0,0)) + trkCov(mk,mk);
605 if (R == 0.0) {
606 ATH_MSG_DEBUG( "inversion of the error-on-the-residual failed.");
607 return nullptr;
608 } R = 1./R;
609
610 // --- compute Kalman gain matrix
611 ROOT::Math::SMatrix<double,5,1,ROOT::Math::MatRepStd<double, 5, 1> >
612 K = trkCov.Sub<ROOT::Math::SMatrix<double,5,1,ROOT::Math::MatRepStd<double, 5, 1> > >(0,mk)*R;
613 if (msgLvl(MSG::VERBOSE)) {
614 ATH_MSG_VERBOSE("-U- residual: r = " << r );
615 ATH_MSG_VERBOSE( "-U- inv. sigmaR = " << R);
616 ATH_MSG_VERBOSE( "-U- gain mtx K=("
617 <<std::setiosflags(std::ios::fixed | std::ios::showpoint | std::ios::right)
618 << std::setw(7) << std::setprecision(4) << K(0,0)<<", "
619 << std::setw(7) << std::setprecision(4) << K(1,0)<<", "
620 << std::setw(7) << std::setprecision(4) << K(2,0)<<", "
621 << std::setw(7) << std::setprecision(4) << K(3,0)<<", "
622 << std::setw(7) << std::setprecision(4) << K(4,0)<<")"
623 << std::resetiosflags(std::ios::fixed));
624 }
625
626 // --- compute local filtered state, here = TP+K*r = TP + TCov * H.T * R * r
627 SParVector5 newPar = trkPar + trkCov.Col(mk) * R * r;
628
629 if (!thetaWithinRange_5D(newPar)) {
630 if ( mk!=Trk::theta &&
631 (std::abs(R*r)>1.0 || trkCov(Trk::theta,Trk::theta) > 0.1*m_cov0(Trk::theta))) {
632 ATH_MSG_DEBUG( "calculateFS_1D(): decided to damp update of theta and re-calculate." );
633 SParVector5 dampedCov = trkCov.Col(mk);
634 dampedCov(Trk::theta) = dampedCov(Trk::theta)*m_thetaGainDampingValue;
635 newPar = trkPar + dampedCov * R * r;
637 ATH_MSG_DEBUG( "-U- damped gain K=("
638 <<std::setiosflags(std::ios::fixed | std::ios::showpoint | std::ios::right)
639 << std::setw(7) << std::setprecision(4) << K(0,0)<<", "
640 << std::setw(7) << std::setprecision(4) << K(1,0)<<", "
641 << std::setw(7) << std::setprecision(4) << K(2,0)<<", "
642 << std::setw(7) << std::setprecision(4) << K(3,0)<<", "
643 << std::setw(7) << std::setprecision(4) << K(4,0)<<")"
644 << std::resetiosflags(std::ios::fixed) );
645 } else {
646 ATH_MSG_DEBUG( "-U- theta out of range but can not damp this update.");
647 }
648 }
649
650 // --- compute covariance matrix of local filteres state
651 SGenMatrix5 KtimesH ; KtimesH.Place_at(K,0,mk);
652 SGenMatrix5 M(m_unitMatrix - KtimesH);
653 ATH_MSG_VERBOSE( "-U- matrix M: diag=("
654 << M(0,0)<<"," << M(1,1)<<","
655 << M(2,2)<<"," << M(3,3)<<","
656 << M(4,4) <<")" );
657 SCovMatrix5 newCov;
658 if (!m_useFruehwirth8a) {
659 // either: use formula C = M * trkCov * M.T() +/- K * covRio * K.T()
660 SCovMatrix1 measuredSMatrix1D(measCov(0,0));
661 newCov = ROOT::Math::Similarity(M,trkCov)
662 + sign*(ROOT::Math::Similarity(K,measuredSMatrix1D));
663 } else {
664 // or: original Fruehwirth eq. 8a is simpler, expression checked to be symm.
665 // C = (1 - K*H) * trkCov
666 /* SGenMatrix5 KtimesH ; KtimesH.Place_at(K,0,mk);
667 SCovMatrix5 M;
668 ROOT::Math::AssignSym::Evaluate(M, KtimesH * trkCov);
669 newCov -= M;*/
670 ROOT::Math::AssignSym::Evaluate(newCov, M * trkCov);
671 }
673 !correctThetaPhiRange_5D(newPar, newCov, Trk::absoluteCheck) : false ) {
674 ATH_MSG_WARNING("calculateFS_1D(): bad angles in filtered state!");
675 return nullptr;
676 }
677
678
679 if (createFQoS) {
680 double predictedResidual = (sign<0) ?r:measPar - newPar(mk);
681 SCovMatrix5 updatedCov = (sign<0) ?
682 trkCov : // when removing, the input are updated par
683 newCov ; // when adding, chi2 is made from upd. par
684
685 // for both signs (add/remove) the chi2 is now calculated like for updated states
686 double chiSquared = measCov(0,0) - updatedCov(mk,mk);
687 if (chiSquared == 0.0)
688 ATH_MSG_DEBUG( "division by zero in 1D chi2, set chi2 to 0.0 instead" );
689 else {
690 // get chi2 = r.T() * R^-1 * r
691 chiSquared = predictedResidual*predictedResidual/chiSquared;
692 ATH_MSG_VERBOSE( "-U- fitQuality of "<< (sign>0?" added ":"removed")
693 <<" state, chi2 :" << chiSquared << " / ndof= 1" );
694 }
695 fitQoS = new FitQualityOnSurface(chiSquared, 1);
696 }
697 return convertToClonedTrackPars(TP,newPar,newCov,sign,createFQoS,"1D");
698}
699
700// calculations for Kalman updator and inverse Kalman filter
701std::unique_ptr<Trk::TrackParameters>
703 const SParVector5& trkPar,
704 const SCovMatrix5& trkCov,
705 const SParVector2& SmeasPar,
706 int paramKey,
707 const Amg::MatrixX& measCov,
708 const int sign,
710 bool createFQoS ) const {
711 // make reduction matrix H from LocalParameters
712 ROOT::Math::SMatrix<double,2,5,ROOT::Math::MatRepStd<double, 2, 5> > H;
713 // and transform EDM to new math lib
714 SCovMatrix2 SmeasCov;
715 for (int i=0, irow=0; i<5; ++i) {
716 if (paramKey & (1<<i)) {
717 SParVector5 v; v(i)=1.0;
718 H.Place_in_row(v, irow, 0);
719 SmeasCov(0,irow) = measCov(0,irow);
720 SmeasCov(1,irow) = measCov(1,irow);
721 ++irow;
722 }
723 }
724 SParVector2 r = SmeasPar - H * trkPar;
725 // FIXME catchPiPi;
726 SCovMatrix2 R = sign * SmeasCov + projection_2D(trkCov,paramKey);
727 if (!R.Invert()) {
728 ATH_MSG_DEBUG( "inversion of residual error matrix (2D) failed.");
729 return nullptr;
730 }
731
732 // --- compute Kalman gain matrix
733 ROOT::Math::SMatrix<double,5,2,ROOT::Math::MatRepStd<double, 5, 2> >
734 K = trkCov * ROOT::Math::Transpose(H) * R;
735 if (msgLvl(MSG::VERBOSE)) {
736 SParVector5 trans_r; SCovMatrix5 trans_R; SGenMatrix5 trans_K;
737 logGainForm (2, trans_r.Place_at(r,0),
738 trans_R.Place_at(R,0,0),
739 trans_K.Place_at(K,0,0));
740 }
741
742 // --- compute local filtered state
743 SParVector5 newPar = trkPar + K * r;
744 if (!thetaWithinRange_5D(newPar)) {
745 if (H(0,Trk::theta) != 1.0 && H(1,Trk::theta) != 1.0 &&
746 ( std::abs(R(0,0)*r(0))>1.0 || std::abs(R(1,1)*r(1))>1.0 ||
747 trkCov(Trk::theta,Trk::theta) > 0.1*m_cov0(Trk::theta))) {
748 ATH_MSG_DEBUG( "calculateFS_2D(): decided to damp update of theta and re-calculate.");
751 newPar = trkPar + K * r;
752 if (msgLvl(MSG::DEBUG)) {
753 msg(MSG::DEBUG) << "-U- damped gain K0=("
754 <<std::setiosflags(std::ios::fixed | std::ios::showpoint | std::ios::right)
755 << std::setw(7) << std::setprecision(4) << K(0,0)<<", "
756 << std::setw(7) << std::setprecision(4) << K(1,0)<<", "
757 << std::setw(7) << std::setprecision(4) << K(2,0)<<", "
758 << std::setw(7) << std::setprecision(4) << K(3,0)<<", "
759 << std::setw(7) << std::setprecision(4) << K(4,0)<<")"
760 << std::resetiosflags(std::ios::fixed) << endmsg;
761 msg(MSG::DEBUG) << "-U- K1=("
762 <<std::setiosflags(std::ios::fixed | std::ios::showpoint | std::ios::right)
763 << std::setw(7) << std::setprecision(4) << K(0,1)<<", "
764 << std::setw(7) << std::setprecision(4) << K(1,1)<<", "
765 << std::setw(7) << std::setprecision(4) << K(2,1)<<", "
766 << std::setw(7) << std::setprecision(4) << K(3,1)<<", "
767 << std::setw(7) << std::setprecision(4) << K(4,1)<<")"
768 << std::resetiosflags(std::ios::fixed) << endmsg;
769 }
770 } else {
771 ATH_MSG_DEBUG( "-U- theta out of range but can not damp this update.");
772 }
773 }
774
775 // --- compute filtered covariance matrix
776 SGenMatrix5 M = m_unitMatrix - K * H;
777 SCovMatrix5 newCov;
778 if (!m_useFruehwirth8a) {
779 // compute covariance matrix of local filteres state
780 // C = M * covTrk * M.T() +/- K * covRio * K.T()
781 newCov = ROOT::Math::Similarity(M,trkCov)
782 + sign*(ROOT::Math::Similarity(K,SmeasCov));
783 } else {
784 // C = (1 - KH) trkCov = trkCov - K*H*trkCov
785 /* SGenMatrix5 KtimesH = K*H;
786 SCovMatrix5 M;
787 ROOT::Math::AssignSym::Evaluate(M, KtimesH * trkCov);
788 newCov -= M; */
789 ROOT::Math::AssignSym::Evaluate(newCov, M * trkCov);
790 }
792 !correctThetaPhiRange_5D(newPar, newCov, Trk::absoluteCheck) : false ) {
793 ATH_MSG_WARNING( "calculateFS_2D(): bad angles in filtered state!" );
794 return nullptr;
795 }
796
797 if (createFQoS) { // get chi2 = r.T() * R2^-1 * r
798 double chiSquared = (sign>0) ?
799 // when adding, the r and R are ready to for calculating the predicted chi2
800 ROOT::Math::Similarity(r,R) :
801 // when removing the r and -R are ready for calculating the updated chi2.
802 ROOT::Math::Similarity(r,-R);
803 ATH_MSG_DEBUG( "-U- fitQuality of "<< (sign>0?" added ":"removed")
804 <<" state, chi2 :" << chiSquared << " / ndof= 2" );
805 fitQoS = new FitQualityOnSurface(chiSquared, 2);
806 }
807 return convertToClonedTrackPars(TP,newPar,newCov,sign,createFQoS,"2D");
808}
809
810// calculations for Kalman updator and inverse Kalman filter
811std::unique_ptr<Trk::TrackParameters>
813 const SParVector5& trkPar,
814 const SCovMatrix5& trkCov,
815 const Trk::LocalParameters& measPar,
816 const Amg::MatrixX& measCov,
817 const int sign,
819 bool createFQoS ) const {
820 // transform EDM to new math lib
821 SParVector3 SmeasPar(measPar[Trk::locX],measPar[Trk::locY],measPar[Trk::phi]);
822 SCovMatrix3 SmeasCov;
823 // make reduction matrix H from LocalParameters
824 ROOT::Math::SMatrix<double,3,5,ROOT::Math::MatRepStd<double, 3, 5> > H;
825 for (int i=0, irow=0; i<5; ++i) {
826 if (measPar.parameterKey() & (1<<i)) {
827 SParVector5 v; v(i)=1.0;
828 H.Place_in_row(v, irow, 0);
829 SmeasCov(0,irow) = measCov(0,irow);
830 SmeasCov(1,irow) = measCov(1,irow);
831 SmeasCov(2,irow) = measCov(2,irow);
832 ++irow;
833 }
834 }
835 SParVector3 r = SmeasPar - H * trkPar;
836 // FIXME catchPiPi;
837 SCovMatrix3 R = sign * SmeasCov + projection_3D(trkCov,measPar.parameterKey());
838 if (!R.Invert()) {
839 ATH_MSG_DEBUG( "inversion of residual error matrix (3D) failed.");
840 return nullptr;
841 }
842
843 // compute Kalman gain matrix
844 ROOT::Math::SMatrix<double,5,3,ROOT::Math::MatRepStd<double, 5, 3> >
845 K = trkCov * ROOT::Math::Transpose(H) * R;
846 SGenMatrix5 M = m_unitMatrix - K * H;
847 if (msgLvl(MSG::VERBOSE)) {
848 SParVector5 trans_r; SCovMatrix5 trans_R; SGenMatrix5 trans_K;
849 logGainForm (3, trans_r.Place_at(r,0),
850 trans_R.Place_at(R,0,0), trans_K.Place_at(K,0,0));
851 }
852
853 // compute local filtered state
854 SParVector5 newPar = trkPar + K * r;
855
856 SCovMatrix5 newCov;
857 if (!m_useFruehwirth8a) {
858 // compute covariance matrix of local filteres state
859 // C = M * covTrk * M.T() +/- K * covRio * K.T()
860 newCov = ROOT::Math::Similarity(M,trkCov)
861 + sign*(ROOT::Math::Similarity(K,SmeasCov));
862 } else {
863 // C = (1 - KH) trkCov = trkCov - K*H*trkCov
864 /* SGenMatrix5 KtimesH = K*H;
865 SCovMatrix5 M;
866 ROOT::Math::AssignSym::Evaluate(M, KtimesH * trkCov);
867 newCov -= M; */
868 ROOT::Math::AssignSym::Evaluate(newCov, M * trkCov);
869 }
870
872 !correctThetaPhiRange_5D(newPar,newCov,Trk::absoluteCheck) : false ) {
873 ATH_MSG_WARNING( "calculateFS_3D(): bad angles in filtered state!" );
874 return nullptr;
875 }
876
877 if (createFQoS) {
878 SParVector3 predictedResidual = (sign<0) ?
879 r :
880 SmeasPar - H * newPar ;
881 SCovMatrix5 updatedCov = (sign<0) ?
882 trkCov : // when removing, the input are updated par
883 newCov ; // when adding, chi2 is made from upd. par
884
885 // for both signs (add/remove) the chi2 is now calculated like for updated states
886 SCovMatrix3 R2 = SmeasCov - projection_3D(updatedCov,measPar.parameterKey());
887 double chiSquared;
888 if ( !R2.Invert() ) {
889 ATH_MSG_DEBUG( "matrix (3D) inversion not possible, set chi2 to zero");
890 chiSquared = 0.f;
891 } else {
892 // get chi2 = r.T() * R2^-1 * r
893 chiSquared = ROOT::Math::Similarity(predictedResidual,R2);
894 ATH_MSG_VERBOSE( "-U- fitQuality of "<< (sign>0?" added ":"removed")
895 <<" state, chi2 :" << chiSquared << " / ndof= 3" );
896 }
897 fitQoS = new FitQualityOnSurface(chiSquared, 3);
898 }
899 return convertToClonedTrackPars(TP,newPar,newCov,sign,createFQoS,"3D");
900}
901
902// calculations for Kalman updator and inverse Kalman filter
903std::unique_ptr<Trk::TrackParameters>
905 const SParVector5& trkPar,
906 const SCovMatrix5& trkCov,
907 const Trk::LocalParameters& measPar,
908 const Amg::MatrixX& measCov,
909 const int sign,
911 bool createFQoS ) const {
912 // transform EDM to new math lib
913 SParVector4 SmeasPar(measPar[Trk::d0],measPar[Trk::z0],
914 measPar[Trk::phi0],measPar[Trk::theta]);
915 SCovMatrix4 SmeasCov;
916
917 // make reduction matrix H from LocalParameters
918 ROOT::Math::SMatrix<double,4,5,ROOT::Math::MatRepStd<double, 4, 5> > H;
919 for (int i=0, irow=0; i<5; ++i) {
920 if (measPar.parameterKey() & (1<<i)) {
921 SParVector5 v; v(i)=1.0;
922 H.Place_in_row(v, irow, 0);
923 SmeasCov(0,irow) = measCov(0,irow);
924 SmeasCov(1,irow) = measCov(1,irow);
925 SmeasCov(2,irow) = measCov(2,irow);
926 SmeasCov(3,irow) = measCov(3,irow);
927 ++irow;
928 }
929 }
930 SParVector4 r = SmeasPar - H * trkPar;
931 // FIXME catchPiPi;
932 SCovMatrix4 R = sign * SmeasCov + projection_4D(trkCov,measPar.parameterKey());
933 if (!R.Invert()) {
934 ATH_MSG_DEBUG( "inversion of residual error matrix (4D) failed." );
935 return nullptr;
936 }
937
938 // compute Kalman gain matrix
939 ROOT::Math::SMatrix<double,5,4,ROOT::Math::MatRepStd<double, 5, 4> >
940 K = trkCov * ROOT::Math::Transpose(H) * R;
941 SGenMatrix5 M = m_unitMatrix - K * H;
942 if (msgLvl(MSG::VERBOSE)) {
943 SParVector5 trans_r; SCovMatrix5 trans_R; SGenMatrix5 trans_K;
944 logGainForm (4, trans_r.Place_at(r,0),
945 trans_R.Place_at(R,0,0), trans_K.Place_at(K,0,0));
946 }
947
948 // compute local filtered state
949 SParVector5 newPar = trkPar + K * r;
950
951 SCovMatrix5 newCov;
952 if (!m_useFruehwirth8a) {
953 // compute covariance matrix of local filteres state
954 // C = M * covTrk * M.T() +/- K * covRio * K.T()
955 newCov = ROOT::Math::Similarity(M,trkCov)
956 + sign*(ROOT::Math::Similarity(K,SmeasCov));
957 } else {
958 // C = (1 - KH) trkCov = trkCov - K*H*trkCov
959 /* SGenMatrix5 KtimesH = K*H;
960 SCovMatrix5 M;
961 ROOT::Math::AssignSym::Evaluate(M, KtimesH * trkCov);
962 newCov -= M; */
963 ROOT::Math::AssignSym::Evaluate(newCov, M * trkCov);
964 }
965
967 !correctThetaPhiRange_5D(newPar,newCov,Trk::absoluteCheck) : false ) {
968 ATH_MSG_WARNING( "calculateFS_4D(): bad angles in filtered state!" );
969 return nullptr;
970 }
971
972 if (createFQoS) {
973 SParVector4 predictedResidual = (sign<0) ?r:SmeasPar - H * newPar;
974 SCovMatrix5 updatedCov = (sign<0) ?
975 trkCov : // when removing, the input are updated par
976 newCov ; // when adding, chi2 is made from upd. par
977
978 // for both signs (add/remove) the chi2 is now calculated like for updated states
979 SCovMatrix4 R2 = SmeasCov - projection_4D(updatedCov,measPar.parameterKey());
980 double chiSquared;
981 if ( !R2.Invert() ) {
982 ATH_MSG_DEBUG( "matrix (4D) inversion not possible, set chi2 to zero");
983 chiSquared = 0.f;
984 } else {
985 // get chi2 = r.T() * R2^-1 * r
986 chiSquared = ROOT::Math::Similarity(predictedResidual,R2);
987 ATH_MSG_VERBOSE( "-U- fitQuality of "<< (sign>0?" added ":"removed")
988 <<" state, chi2 :" << chiSquared << " / ndof= 4" );
989 }
990 fitQoS = new FitQualityOnSurface(chiSquared, 4);
991 }
992 return convertToClonedTrackPars(TP,newPar,newCov,sign,createFQoS,"4D");
993}
994
995// calculations for Kalman updator and inverse Kalman filter
996std::unique_ptr<Trk::TrackParameters>
998 const SParVector5& trkParOne,
999 const SCovMatrix5& trkCovOne,
1000 const SParVector5& trkParTwo,
1001 const Amg::MatrixX& measCov,
1002 const int sign,
1003 Trk::FitQualityOnSurface*& fitQoS,
1004 bool createFQoS ) const {
1005 SCovMatrix5 trkCovTwo; // transform EDM to new math lib
1006 for (int i=0; i<5; ++i)
1007 for (int j=0; j<=i; ++j) {
1008 trkCovTwo(j,i) = measCov(j,i);
1009 }
1010
1011 // this update is symmetric in state and measurement parameters, H == 1.
1012 SCovMatrix5 R = sign * trkCovTwo + trkCovOne;
1013 if (!R.Invert()) {
1014 ATH_MSG_DEBUG( "inversion of residual error matrix (5D) failed." );
1015 return nullptr;
1016 }
1017
1018 SParVector5 r = trkParTwo - trkParOne;
1021 ATH_MSG_WARNING( "calculateFS_5D(): bad angles in intermediate residual!" );
1022 return nullptr;
1023 }
1024 SGenMatrix5 K = trkCovOne * R;
1025 if (msgLvl(MSG::VERBOSE)) logGainForm (5, r, R, K);
1026
1027 // compute local filtered state
1028 SParVector5 newPar = trkParOne + K * r;
1029 SCovMatrix5 newCov;
1030 /* if (!m_useFruehwirth8a) {
1031 SGenMatrix5 M = m_unitMatrix - K;
1032 // compute covariance matrix of local filtered state
1033 // C = M * covTrk * M.T() +/- K * covRio * K.T()
1034 newCov = ROOT::Math::Similarity(M,trkCovOne)
1035 + sign*(ROOT::Math::Similarity(K,trkCovTwo));
1036 } else { */
1037 // original Fruehwirth eq. 8a is simpler, expression checked to be symm.
1038 ROOT::Math::AssignSym::Evaluate(newCov, K * trkCovTwo);
1039 // }
1040
1042 !correctThetaPhiRange_5D(newPar,newCov,Trk::absoluteCheck) : false ) {
1043 ATH_MSG_WARNING( "calculateFS_5D(): bad angles in filtered state!" );
1044 return nullptr;
1045 }
1046
1047 bool goodMatrix(true);
1048 for (int i=0; i<5; ++i) {
1049 if (newCov(i,i) < 0.0 && goodMatrix ) goodMatrix=false;
1050 }
1051 if (!goodMatrix) {
1052 SGenMatrix5 M = m_unitMatrix - K;
1053 // compute covariance matrix of local filtered state
1054 // C = M * covTrk * M.T() +/- K * covRio * K.T()
1055 newCov = ROOT::Math::Similarity(M,trkCovOne)
1056 + sign*(ROOT::Math::Similarity(K,trkCovTwo));
1057 goodMatrix = true;
1058 for (int i=0; i<5; ++i) {
1059 if (newCov(i,i) < 0.0 && goodMatrix ) goodMatrix=false;
1060 if (!goodMatrix) {
1061 ATH_MSG_DEBUG("calculateFS_5D(): unphysical cov!");
1062 return nullptr;
1063 }
1064 }
1065 }
1066
1067 if (createFQoS) { // get chi2 = r.T() * R2^-1 * r
1068 double chiSquared = (sign>0) ?
1069 // when adding, the r and R are ready to for calculating the predicted chi2
1070 ROOT::Math::Similarity(r,R) :
1071 // when removing the r and -R are ready for calculating the updated chi2.
1072 ROOT::Math::Similarity(r,-R);
1073 ATH_MSG_VERBOSE( "-U- fitQuality of "<< (sign>0?" added ":"removed")
1074 <<" state, chi2 :" << chiSquared << " / ndof= 5" );
1075 fitQoS = new FitQualityOnSurface(chiSquared, 5);
1076 }
1077 return convertToClonedTrackPars(TP,newPar,newCov,sign,createFQoS,"5D");
1078}
1079
1080
1082 const SParVector5& parTrk,
1083 const AmgSymMatrix(5)& covTrk,
1084 double valRio,
1085 double covRio,
1086 int key,
1087 int sign) const
1088{ // sign: -1 = updated, +1 = predicted parameters.
1089 int mk=0;
1090 if (key!=1) for (int i=0; i<5; ++i) if (key & (1<<i)) mk=i;
1091 double r = valRio - parTrk(mk);
1092 // if (mk==3) catchPiPi;
1093 double chiSquared = covRio + sign * covTrk(mk,mk);
1094 if (chiSquared == 0.0) {
1095 ATH_MSG_DEBUG( "inversion of the error-on-the-residual failed." );
1096 return {};
1097 }
1098 chiSquared = r*r/chiSquared;
1099
1100 return {chiSquared, 1};
1101}
1102
1104 const SParVector5& parTrk,
1105 const AmgSymMatrix(5)& covTrk,
1106 const SParVector2& parRio,
1107 const SCovMatrix2& covRio,
1108 int key,
1109 int sign) const
1110{ // sign: -1 = updated, +1 = predicted parameters.
1111
1112 ROOT::Math::SVector<int,2> index(0,1); // locX and Y if key==3
1113 if (key != 3) { // other localPar expansion vector if not key==3.
1114 for (int i=0, irow=0; i<5; ++i)
1115 if (key & (1<<i)) index(irow++) = i;
1116 }
1117 SParVector2 r(-parTrk(index(0)),-parTrk(index(1))); r += parRio;
1118 SCovMatrix2 R = sign*projection_2D(covTrk,key); R += covRio;
1119 double chiSquared = 0.0;
1120 if (!R.Invert()) {
1121 ATH_MSG_DEBUG( "matrix inversion not possible, set chi2 to zero" );
1122 } else {
1123 chiSquared = ROOT::Math::Similarity(r,R);
1124 ATH_MSG_VERBOSE( "-U- fitQuality of "<< (sign>0?"predicted":"updated")
1125 <<" state, chi2 :" << chiSquared << " / ndof= 2" );
1126 }
1127 return {chiSquared, 2};
1128}
1129
1131 const SParVector5& parOne,
1132 const AmgSymMatrix(5)& covOne,
1133 const SParVector5& parTwo,
1134 const AmgSymMatrix(5)& covTwo,
1135 int sign) const
1136{ // sign: -1 = updated, +1 = predicted parameters.
1137
1138 SCovMatrix5 ScovOne;
1139 SCovMatrix5 ScovTwo; // trafo EDM to new math lib
1140 for (int i=0; i<5; ++i)
1141 for (int j=0; j<=i; ++j) {
1142 ScovOne(i,j) = covOne(i,j);
1143 ScovTwo(i,j) = covTwo(i,j);
1144 }
1145 SParVector5 r = parTwo - parOne;
1146 SCovMatrix5 R = sign*ScovOne + ScovTwo;
1147 double chiSquared = 0.0;
1148 if (!R.Invert()) {
1149 ATH_MSG_DEBUG( "matrix inversion not possible, set chi2 to zero" );
1150 } else {
1151 chiSquared = ROOT::Math::Similarity(r,R);
1152 ATH_MSG_VERBOSE( "-U- fitQuality of "<< (sign>0?"predicted":"updated")
1153 <<" state, chi2 :" << chiSquared << " / ndof= 2" );
1154 }
1155 return {chiSquared, 5};
1156}
1157
1158std::unique_ptr<Trk::TrackParameters>
1160 const SParVector5& par,
1161 const SCovMatrix5& covpar,
1162 int sign,
1163 bool createFQoS,
1164 std::string_view ndtext) const {
1165 AmgSymMatrix(5) C;
1166 C.setZero();
1167 for (int i=0; i<5; ++i) {
1168 for (int j=0; j<=i; ++j) {
1169 C.fillSymmetric(i,j, covpar(i,j));
1170 }
1171 }
1172
1173 std::unique_ptr<Trk::TrackParameters> resultPar =
1174 TP.associatedSurface().createUniqueTrackParameters(par[0],par[1],par[2],par[3],par[4],C);
1175
1176 if (msgLvl(MSG::VERBOSE) && resultPar) {
1177 char reportCalledInterface[80];
1178 char ndtext2[5];
1179 memset(ndtext2, '\0', 5 ); ndtext.copy(ndtext2,2); // convert char to string
1180 if (sign>0)
1181 sprintf(reportCalledInterface,"%s-%s,%s)",
1182 (ndtext=="5D"?"combineStates(TP,TP":"addToState(TP,Meas"),
1183 ndtext2,(createFQoS?"Err,FQ":"Err"));
1184 else
1185 sprintf(reportCalledInterface,"%s,Meas-%s,%s)","removeFromState(TP,",
1186 ndtext2,(createFQoS?"Err,FQ":"Err"));
1187 logResult((std::string)reportCalledInterface, resultPar->parameters(),C);
1188 }
1189
1190 return resultPar;
1191}
1192
1194{
1195 if (key == 3) { // shortcut the most-used case
1196 SCovMatrix2 S = M.Sub<SCovMatrix2> (0,0);
1197 return S;
1198 }
1199 ROOT::Math::SVector<int,2> iv;
1200 for (int i=0,k=0; i<5; ++i) { if (key & (1<<i)) iv(k++)=i; }
1201 SCovMatrix2 covSubMatrix;
1202 for (int i=0; i<2; ++i) {
1203 for (int j=0; j<2; ++j) {
1204 covSubMatrix(i,j) = M(iv(i),iv(j));
1205 }
1206 }
1207 return covSubMatrix;
1208
1209}
1210
1212 int key)
1213{
1214 ROOT::Math::SVector<int,2> iv;
1215 for (int i=0,k=0; i<5; ++i) { if (key & (1<<i)) iv(k++)=i; }
1216 SCovMatrix2 covSubMatrix;
1217 for (int i=0; i<2; ++i) {
1218 for (int j=0; j<2; ++j) {
1219 covSubMatrix(i,j) = M(iv(i),iv(j));
1220 }
1221 }
1222 return covSubMatrix;
1223}
1224
1226{
1227 if (key == 7) { // shortcut the most-used case
1228 return M.Sub<SCovMatrix3> (0,0);
1229 }
1230 ROOT::Math::SVector<int,3> iv;
1231 for (int i=0,k=0; i<5; ++i) { if (key & (1<<i)) iv(k++)=i; }
1232 SCovMatrix3 covSubMatrix;
1233 for (int i=0; i<3; ++i) {
1234 for (int j=0; j<3; ++j) {
1235 covSubMatrix(i,j) = M(iv(i),iv(j));
1236 }
1237 }
1238 return covSubMatrix;
1239
1240}
1241
1243{
1244 if (key == 15) { // shortcut the most-used case
1245 return M.Sub<SCovMatrix4> (0,0);
1246 }
1247 ROOT::Math::SVector<int,4> iv;
1248 for (int i=0,k=0; i<5; ++i) { if (key & (1<<i)) iv(k++)=i; }
1249 SCovMatrix4 covSubMatrix;
1250 for (int i=0; i<4; ++i) {
1251 for (int j=0; j<4; ++j) {
1252 covSubMatrix(i,j) = M(iv(i),iv(j));
1253 }
1254 }
1255 return covSubMatrix;
1256
1257}
1258
1260 const Trk::TrackParameters& inputParameters,
1261 const int isign) const
1262{
1263
1264 const AmgSymMatrix(5)* covariance = inputParameters.covariance();
1265 if (!covariance) {
1266 if (isign<0) {
1267 ATH_MSG_WARNING ("MeasuredTrackParameters == Null, can not calculate updated parameter state.");
1268 return false;
1269 }
1270 // no error given - use a huge error matrix for the time
1271 ATH_MSG_VERBOSE ("-U- no covTrk at input - assign large error matrix for the time being.");
1272 M.SetDiagonal(m_cov0);
1273
1274 } else {
1275 // int k = measTrkPar->localAmg::MatrixX().covariance().cols(); is always 5
1276 for (int i=0; i<5; ++i) {
1277 for (int j=0; j<=i; ++j) {
1278 M(i,j) = (*covariance)(i,j); // can improve speed by apply() access
1279 }
1280 }
1281 }
1282 return true;
1283}
1284
1285
1287 const AmgSymMatrix(5)& covTrk,
1288 const Amg::MatrixX& covRio,
1289 const Amg::MatrixX& H,
1290 int sign) const
1291{ // sign: -1 = updated, +1 = predicted parameters.
1292 Amg::MatrixX R = covRio + sign * covTrk.similarity(H);
1293 double chiSquared = 0.;
1294 if (false) {
1295 ATH_MSG_DEBUG( "matrix inversion not possible, set chi2 to zero" );
1296 chiSquared = 0.f;
1297 } else {
1298 // get chi2 = r.T() * R^-1 * r
1299 chiSquared = Amg::chi2(R.inverse(), residual);
1300 ATH_MSG_VERBOSE( "-U- fitQuality of "<< (sign>0?"predicted":"updated")
1301 <<" state, chi2 :" << chiSquared << " / ndof= " << covRio.cols() );
1302 }
1303
1304 // number of degree of freedom added
1305 int numberDoF = covRio.cols();
1306
1307 return {chiSquared, numberDoF};
1308}
1309
1311 int dimCov) const {
1312 if (P.dimension() != dimCov ) {
1313 ATH_MSG_WARNING ("Inconsistency in dimension of local coord - problem with LocalParameters object?");
1314 ATH_MSG_WARNING ("dim of local parameters: "<<P.dimension()<<" vs. dim of error matrix: "<<dimCov);
1315 return false;
1316 } return true;
1317}
1318
1320 const RangeCheckDef rcd) const
1321{
1322 static const SParVector2 thetaMin(0.0, -M_PI);
1323
1324 // correct theta coordinate
1325 if ( V(Trk::theta)<thetaMin((int)rcd) || V(Trk::theta)> M_PI ) {
1326 ATH_MSG_DEBUG ("-U- theta angle out of range: theta= "<<V(Trk::theta)<<", phi= "<<V(Trk::phi));
1327 // absolute theta: repair if between -pi and +2pi.
1328 // differential theta: repair if between -pi and +pi
1329 if ( ( V(Trk::theta) < -M_PI ) ||
1330 ( V(Trk::theta) > (rcd==Trk::differentialCheck ? M_PI : 2*M_PI) )
1331 ) {
1332 ATH_MSG_WARNING ("-U- track theta = "<<V(Trk::theta)<<" -> this is too far from defined range, stop update.");
1333 return false;
1334 }
1335 if (V(Trk::theta) > M_PI) {
1336 V(Trk::theta) = 2*M_PI - V(Trk::theta);
1337 V(Trk::phi) += (V(Trk::phi)>0.0) ? -M_PI : M_PI;
1338 } else if (V(Trk::theta) < 0.0) {
1339 V(Trk::theta) = -V(Trk::theta);
1340 V(Trk::phi) += (V(Trk::phi)>0.0) ? -M_PI : M_PI;
1341 }
1342 // correct also cov matrix:
1343 M(0,3) = -M(0,3); // cov(polar,locX)
1344 M(1,3) = -M(1,3); // cov(polar,locY)
1345 M(2,3) = -M(2,3); // cov(polar,azimuthal)
1346 M(3,4) = -M(3,4); // cov(polar,qOverP)
1347 ATH_MSG_DEBUG ("-U- now use corrected value phi= "<<V(Trk::phi)<<" theta = "<<V(Trk::theta));
1348 }
1349
1350
1351 // correct phi coordinate if necessary
1352 if ( (V(Trk::phi) > M_PI) ) {
1353 if (msgLvl(MSG::DEBUG)) msg(MSG::DEBUG) << "-U- phi = " << V(Trk::phi);
1354 V(Trk::phi) = std::fmod(V(Trk::phi)+M_PI,2*M_PI)-M_PI;
1355 ATH_MSG_DEBUG( " out of range, now corrected to " << V(Trk::phi) );
1356 } else if( (V(Trk::phi)<-M_PI) ) {
1357 ATH_MSG_DEBUG( "-U- phi = " << V(Trk::phi));
1358 V(Trk::phi) = std::fmod(V(Trk::phi)-M_PI,2*M_PI)+M_PI;
1359 ATH_MSG_DEBUG( " out of range, now corrected to " << V(Trk::phi) );
1360 }
1361
1362 return true;
1363}
1364
1365void Trk::KalmanUpdatorSMatrix::logStart(const std::string& IDstring,
1366 const Trk::TrackParameters& tp) const
1367{
1368 ATH_MSG_DEBUG( "--> entered KalmanUpdatorSMatrix::" << IDstring );
1369 ATH_MSG_VERBOSE( "-U- TrkPar:" << std::setiosflags(std::ios::right)<<std::setprecision(4)
1370 << std::setw( 9)<<tp.parameters()[0]<< std::setw(10)<<tp.parameters()[1]<<std::setprecision(5)
1371 << std::setw(10)<<tp.parameters()[2]<< std::setw(10)<<tp.parameters()[3]<<std::setprecision(4)
1372 << std::setw(10)<<tp.parameters()[4]);
1373}
1374
1376 const Amg::VectorX& parRio, const Amg::MatrixX& covRio) const
1377{
1378 ATH_MSG_VERBOSE( "-U- cov "<<std::setiosflags(std::ios::right)<<std::setprecision(3)
1379 << std::setw(9)<<covTrk(0,0)<<" "<< std::setw(9)<<covTrk(0,1)<<" "
1380 << std::setw(9)<<covTrk(0,2)<<" "<< std::setw(9)<<covTrk(0,3)<<" "
1381 << std::setw(9)<<covTrk(0,4));
1382 ATH_MSG_VERBOSE( " " << " " << " "
1383 << std::setw(9)<<covTrk(1,1)<<" "<< std::setw(9)<<covTrk(1,2)<<" "
1384 << std::setw(9)<<covTrk(1,3)<<" "<< std::setw(9)<<covTrk(1,4));
1385 ATH_MSG_VERBOSE( " covariance matrix " << " " << " "
1386 << std::setw(9)<<covTrk(2,2)<<" "<< std::setw(9)<<covTrk(2,3)<<" "
1387 << std::setw(9)<<covTrk(2,4));
1388 ATH_MSG_VERBOSE( " of the PREDICTED track pars " << " " << " "
1389 << std::setw(9)<<covTrk(3,3)<<" "<< std::setw(9)<<covTrk(3,4));
1390 ATH_MSG_VERBOSE(" " << " " << " "
1391 << std::setw(9)<<covTrk(4,4)<<std::setprecision(6));
1392
1393 int nLocCoord = covRio.cols();
1394 ATH_MSG_VERBOSE( "-U- measurement locPar: ");
1395 for (int i=0; i<nLocCoord; i++) ATH_MSG_VERBOSE(parRio[i]);
1396 msg(MSG::VERBOSE) << "-U- measurement (err)^2: " <<std::setprecision(4)<<covRio(0,0);
1397 for (int i=1; i<nLocCoord; i++) msg(MSG::VERBOSE) << ", "<<covRio(i,i);
1398 msg(MSG::VERBOSE) << std::setprecision(6)<<"\n";
1399}
1400
1402 const SCovMatrix5& R,
1403 const SGenMatrix5& K) const
1404{
1405 // again some verbose debug output showing internals of updating
1406 msg(MSG::VERBOSE) << "-U- residual: r=("<<r(0);
1407 for (int i=1; i<nc; i++) msg(MSG::VERBOSE) <<","<<r(i);
1408 msg(MSG::VERBOSE) << ")\n";
1409 msg(MSG::VERBOSE) << "-U- inv. sigmaR=("<< R(0,0);
1410 for (int i=1; i<nc; i++) msg(MSG::VERBOSE) << "," << R(i,i);
1411 msg(MSG::VERBOSE) << ")\n";
1412 for (int i=0; i<nc; i++)
1413 msg(MSG::VERBOSE) // K is a row x col = 5 x nc matrix.
1414 << ( i==0 ? "-U- gain mtx K=(" : " (" )
1415 << std::setiosflags(std::ios::fixed | std::ios::showpoint | std::ios::right )
1416 << std::setw(7) << std::setprecision(4) << K(0,i)<<", "
1417 << std::setw(7) << std::setprecision(4) << K(1,i)<<", "
1418 << std::setw(7) << std::setprecision(4) << K(2,i)<<", "
1419 << std::setw(7) << std::setprecision(4) << K(3,i)<<", "
1420 << std::setw(7) << std::setprecision(4) << K(4,i)<<")"
1421 << std::resetiosflags(std::ios::fixed) << "\n";
1422}
1423
1424void Trk::KalmanUpdatorSMatrix::logResult(const std::string& methodName,
1425 const AmgVector(5)& par,
1426 const AmgSymMatrix(5)& covPar) const
1427{
1428 // again some verbose debug output
1429 msg(MSG::VERBOSE) << "-U- ==> result for KalmanUpdatorSMatrix::"<<methodName<<endmsg;
1430 msg(MSG::VERBOSE) << "-U- new par"<<std::setiosflags(std::ios::right)<<std::setprecision(4)
1431 << std::setw( 9)<<par[0]<< std::setw(10)<<par[1]<<std::setprecision(5)
1432 << std::setw(10)<<par[2]<< std::setw(10)<<par[3]<<std::setprecision(4)
1433 << std::setw(10)<<par[4] <<"\n";
1434 msg(MSG::VERBOSE) << "-U- new cov" <<std::setiosflags(std::ios::right)<<std::setprecision(3)
1435 << std::setw(9)<<(covPar)(0,0)<<" "<< std::setw(9)<<(covPar)(0,1)<<" "
1436 << std::setw(9)<<(covPar)(0,2)<<" "<< std::setw(9)<<(covPar)(0,3)
1437 << " " << std::setw(9)<<(covPar)(0,4)<< "\n";
1438 msg(MSG::VERBOSE) << " " << " " << " "
1439 << std::setw(9)<<(covPar)(1,1)<<" "<< std::setw(9)<<(covPar)(1,2)<<" "
1440 << std::setw(9)<<(covPar)(1,3)<<" "<< std::setw(9)<<(covPar)(1,4)<< "\n";
1441 msg(MSG::VERBOSE) << " covariance matrix " << " " << " "
1442 << std::setw(9)<<(covPar)(2,2) <<" "<< std::setw(9)<<(covPar)(2,3) <<" "
1443 << std::setw(9)<<(covPar)(2,4) << "\n";
1444 msg(MSG::VERBOSE) << " of the UPDATED track pars " << " "
1445 << " " <<std::setw(9)<<(covPar)(3,3) << " "
1446 << std::setw(9)<<(covPar)(3,4) << "\n";
1447 msg(MSG::VERBOSE) << " " << " "
1448 << " "
1449 << std::setw(9)<<(covPar)(4,4) <<std::setprecision(6)<< "\n";
1450}
#define M_PI
#define endmsg
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define AmgSymMatrix(dim)
#define AmgVector(rows)
ROOT::Math::SVector< double, 3 > SParVector3
ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepSym< double, 2 > > SCovMatrix2
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > SCovMatrix5
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > SGenMatrix5
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > SCovMatrix3
ROOT::Math::SMatrix< double, 1, 1, ROOT::Math::MatRepSym< double, 1 > > SCovMatrix1
ROOT::Math::SVector< double, 2 > SParVector2
ROOT::Math::SMatrix< double, 4, 4, ROOT::Math::MatRepSym< double, 4 > > SCovMatrix4
ROOT::Math::SVector< double, 5 > SParVector5
ROOT::Math::SVector< double, 4 > SParVector4
static Double_t P(Double_t *tt, Double_t *par)
#define H(x, y, z)
Definition MD5.cxx:114
int sign(int a)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
bool msgLvl(const MSG::Level lvl) const
std::unique_ptr< TrackParameters > calculateFilterStep_4D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, const LocalParameters &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
virtual FitQualityOnSurface fullStateFitQuality(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &) const override final
estimator for FitQuality on Surface from a full track state, that is a state which contains the curre...
bool thetaPhiWithinRange_5D(const SParVector5 &V, const RangeCheckDef rcd) const
Test if angles are inside boundaries.
std::unique_ptr< TrackParameters > calculateFilterStep_1D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, double, int, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
static const ParamDefsAccessor s_enumAccessor
std::unique_ptr< TrackParameters > calculateFilterStep_5D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, const SParVector5 &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
void logGainForm(int, const SParVector5 &, const SCovMatrix5 &, const SGenMatrix5 &) const
internal structuring: common logfile output during calculation
virtual FitQualityOnSurface predictedStateFitQuality(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &) const override final
estimator for FitQuality on Surface from a predicted track state, that is a state which contains the ...
SCovMatrix5 m_unitMatrix
avoid mem allocation at every call
bool getStartCov(SCovMatrix5 &, const TrackParameters &, const int) const
Helper method to transform Eigen cov matrix to SMatrix.
FitQualityOnSurface makeChi2_2D(const SParVector5 &, const AmgSymMatrix(5)&, const SParVector2 &, const SCovMatrix2 &, int, int) const
virtual StatusCode initialize() override
AlgTool initialisation.
FitQualityOnSurface makeChi2Object(const Amg::VectorX &, const AmgSymMatrix(5)&, const Amg::MatrixX &, const Amg::MatrixX &, int) const
std::unique_ptr< TrackParameters > convertToClonedTrackPars(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, int, bool, std::string_view) const
Helper method to convert internal results from SMatrix to Eigen. *‍/.
virtual std::vector< double > initialErrors() const override final
give back how updator is configured for inital covariances
std::unique_ptr< TrackParameters > calculateFilterStep_2D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, const SParVector2 &, int, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
virtual std::unique_ptr< TrackParameters > addToState(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &) const override final
measurement updator for the KalmanFitter getting the meas't coord' from Amg::Vector2D (use eg with PR...
static SCovMatrix3 projection_3D(const SCovMatrix5 &, int)
Avoid multiplications with sparse H matrices by cutting 3D rows&columns out of the full cov matrix.
virtual StatusCode finalize() override
AlgTool termination.
bool thetaWithinRange_5D(const SParVector5 &V) const
Test if theta angle is inside boundaries. No differential-check option.
virtual std::unique_ptr< TrackParameters > removeFromState(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &) const override final
reverse update eg for track property analysis (unbiased residuals) getting the measurement coordinate...
bool correctThetaPhiRange_5D(SParVector5 &, SCovMatrix5 &, const RangeCheckDef) const
method correcting the calculated angles back to their defined ranges phi (-pi, pi) and theta (0,...
virtual std::unique_ptr< TrackParameters > combineStates(const TrackParameters &, const TrackParameters &) const override final
trajectory state updator which combines two parts of a trajectory on a common surface.
std::unique_ptr< TrackParameters > prepareFilterStep(const TrackParameters &, const LocalParameters &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common code analysing the measurement's rank and calling the appropriate implementation for this rank...
KalmanUpdatorSMatrix(const std::string &, const std::string &, const IInterface *)
AlgTool standard constuctor.
SParVector5 m_cov0
initial cov values in SMatrix object
void logResult(const std::string &, const AmgVector(5) &, const AmgSymMatrix(5) &) const
internal structuring: common logfile output after calculation
static SCovMatrix2 projection_2D(const SCovMatrix5 &, int)
Avoid multiplications with sparse H matrices by cutting 2D rows&columns out of the full cov matrix.
bool m_useFruehwirth8a
job options controlling update formula for covariance matrix
void logInputCov(const SCovMatrix5 &, const Amg::VectorX &, const Amg::MatrixX &) const
internal structuring: common logfile output of the inputs
std::unique_ptr< TrackParameters > calculateFilterStep_3D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, const LocalParameters &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
std::vector< double > m_cov_stdvec
job options for initial cov values
static SCovMatrix4 projection_4D(const SCovMatrix5 &, int)
Avoid multiplications with sparse H matrices by cutting 4D rows&columns out of the full cov matrix.
FitQualityOnSurface makeChi2_5D(const SParVector5 &, const AmgSymMatrix(5)&, const SParVector5 &, const AmgSymMatrix(5)&, int) const
FitQualityOnSurface makeChi2_1D(const SParVector5 &, const AmgSymMatrix(5)&, double, double, int, int) const
also the chi2 calculation and FitQuality object creation is combined in an extra method.
bool consistentParamDimensions(const LocalParameters &, int) const
method testing correct use of LocalParameters
void logStart(const std::string &, const TrackParameters &) const
internal structuring: debugging output for start of method.
int parameterKey() const
Identifier key for matrix expansion/reduction.
const Amg::MatrixX & expansionMatrix() const
Expansion matrix from 5x5 to the [dimension()]x[dimension()].
int dimension() const
Dimension of this localParameters() vector.
virtual const Surface & associatedSurface() const override=0
Access to the Surface associated to the Parameters.
virtual ChargedTrackParametersUniquePtr createUniqueTrackParameters(double l1, double l2, double phi, double theat, double qop, std::optional< AmgSymMatrix(5)> cov=std::nullopt) const =0
Use the Surface as a ParametersBase constructor, from local parameters - charged.
int r
Definition globals.cxx:22
struct color C
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
double chi2(const T &precision, const U &residual, const int sign=1)
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.
@ locY
local cartesian
Definition ParamDefs.h:38
@ locX
Definition ParamDefs.h:37
@ v
Definition ParamDefs.h:78
@ phi0
Definition ParamDefs.h:65
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ phi
Definition ParamDefs.h:75
@ loc1
Definition ParamDefs.h:34
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParametersBase< TrackParametersDim, Charged > TrackParameters
Definition index.py:1
Simple struct to access the ParamDefs enum with ints.
Definition ParamDefs.h:92
static constexpr std::array< ParamDefs, 6 > pardef
Constructor.
Definition ParamDefs.h:94
MsgStream & msg
Definition testRead.cxx:32