ATLAS Offline Software
Loading...
Searching...
No Matches
AnalyticalDerivCalcTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "TrkTrack/Track.h"
9
12
15
21
23
25
26#include <string>
27
28namespace Trk {
29
30 //________________________________________________________________________
32 const std::string & name,
33 const IInterface * parent)
34 : AthAlgTool(type,name,parent)
35 {
36 declareInterface<IDerivCalcTool>(this);
37 }
38
39 //________________________________________________________________________
41 {
42 if (m_alignModuleTool.retrieve().isFailure()) {
43 ATH_MSG_FATAL("Could not get " << m_alignModuleTool);
44 return StatusCode::FAILURE;
45 }
46 ATH_MSG_INFO("Retrieved " << m_alignModuleTool);
47
48 if (detStore()->retrieve(m_idHelper, "AtlasID").isFailure()) {
49 ATH_MSG_FATAL("Could not get AtlasDetectorID helper");
50 return StatusCode::FAILURE;
51 }
53
54 return StatusCode::SUCCESS;
55 }
56
57 //________________________________________________________________________
59 {
60 ATH_MSG_DEBUG("in AnalyticalDerivCalcTool::finalize()");
61 return StatusCode::SUCCESS;
62 }
63
64 //________________________________________________________________________
66 {
68 checkResidualType(alignTrack);
69
70 // create table of modules for checking whether the module
71 // is hit by the track or not
72 int nModules = m_alignModuleTool->alignModules1D()->size();
73 std::vector<bool> hitModules(nModules,false);
74
75 // loop over AlignTSOSCollection,
76 // find modules that are in the AlignModuleList,
77 std::vector<AlignModule *> alignModules;
78 AlignTSOSCollection::iterator atsosItr = alignTrack->firstAtsos();
79 for (; atsosItr != alignTrack->lastAtsos(); ++atsosItr) {
80 AlignModule * module=(*atsosItr)->module();
81 if (module)
82 ATH_MSG_DEBUG("have ATSOS for module "<<module->identify());
83 else
84 ATH_MSG_DEBUG("no module!");
85
86 if (!(*atsosItr)->isValid() || !module)
87 continue;
88
89 // if the module is not yet in the list for this track, add it
90 if(!hitModules[module->identifyHash()]) {
91 hitModules[module->identifyHash()] = true;
92 alignModules.push_back(module);
93 }
94 }
95
96 // Determine derivatives from shifting these modules
97 std::vector<AlignModuleDerivatives> * derivatives = new std::vector<AlignModuleDerivatives>;
98 std::vector<AlignModule *>::iterator moduleIt = alignModules.begin();
99 for ( ; moduleIt!=alignModules.end(); ++moduleIt) {
100 std::vector<Amg::VectorX> deriv_vec = getDerivatives(alignTrack,*moduleIt);
101 derivatives->push_back(make_pair(*moduleIt,deriv_vec));
102 }
103
104 // alignTrack takes care of deleting the derivatives
105 ATH_MSG_DEBUG("setting matrix derivatives");
106 alignTrack->setDerivatives(derivatives);
107
108 ATH_MSG_DEBUG("returning from setDerivatives");
109
110 return true;
111 }
112
113 //________________________________________________________________________
115 {
116 static std::once_flag flag;
117 std::call_once(flag, [&]() {
118 if(m_logStream) {
119 *m_logStream<<"*************************************************************"<<std::endl;
120 *m_logStream<<"*************************************************************"<<std::endl;
121 *m_logStream<<"*** *****"<<std::endl;
122
124 *m_logStream<<"*** Running full LOCAL Chi2 method *****"<<std::endl;
125 else
126 *m_logStream<<"*** Running full GLOBAL Chi2 method *****"<<std::endl;
127
128 *m_logStream<<"*** using analytical derivatives *****"<<std::endl;
129
131 std::string resType = "HitOnly ";
133 resType = "Unbiased";
134 *m_logStream<<"*** with residual type: "<<resType<<" *****"<<std::endl;
135 }
136
137 *m_logStream<<"*** *****"<<std::endl;
138 *m_logStream<<"*************************************************************"<<std::endl;
139 }
140 });
141
142 // get inverse local error matrix of the track
143 const Amg::MatrixX * Vinv = alignTrack->localErrorMatrixInv();
144 if(!checkValidity(*Vinv)) {
145 ATH_MSG_WARNING("Inverse local error matrix is invalid, skipping the track");
146 return false;
147 }
148 ATH_MSG_DEBUG("V inverse (diagonal only):");
149 if (msgLvl(MSG::DEBUG)) {
150 for (int i=0;i<Vinv->rows();i++) msg()<<(*Vinv)(i,i)<<" ";
151 msg()<<endmsg;
152 }
153
154 // ========================
155 // setup local chi2 method
156 // ========================
157 if (m_useLocalSetting) {
158 ATH_MSG_DEBUG("setting Residual covariance matrix for local method");
159
160 // for local method we only need the Vinv
161 const Amg::MatrixX * Vinv = alignTrack->localErrorMatrixInv();
162 Amg::MatrixX * W = new Amg::MatrixX(*Vinv);
163
164 // if we want to ignore the real errors of measurements for a particular subretector
165 // (we should add something similar for Muon detectors)
167
169 ATH_MSG_DEBUG("Ignoring measured errors for Pixel clusters, using intrinsic errors");
171 ATH_MSG_DEBUG("Ignoring measured errors for SCT clusters, using intrinsic errors");
173 ATH_MSG_DEBUG("Ignoring measured errors for TRT clusters, using intrinsic errors");
174
175 int index(0);
176 AlignTSOSCollection::const_iterator itAtsos = alignTrack->firstAtsos();
177 AlignTSOSCollection::const_iterator itAtsos_end = alignTrack->lastAtsos();
178 for ( ; itAtsos != itAtsos_end; ++itAtsos) {
179 const AlignTSOS * atsos = *itAtsos;
180 if (!atsos->isValid())
181 continue;
182
184 (*W)(index,index) = 1./(0.05*0.05/12);
185 (*W)(index+1,index+1) = 1./(0.4*0.4/12);
186 }
188 (*W)(index,index) = 1./(0.1*0.1/12);
190 (*W)(index,index) = 1./(0.5*0.5/12);
191
192 index += atsos->nResDim();
193 }
194 }
195
196 ATH_MSG_DEBUG("setting local weight matrix W ( "<<W->rows()<<" x "<<W->cols()<<" ) (diagonal only):");
197 if (msgLvl(MSG::DEBUG)) {
198 for (int i=0;i<W->rows();i++) msg()<<(*W)(i,i)<<" ";
199 msg()<<endmsg;
200 }
201 alignTrack->setWeightMatrix(W);
202
203 Amg::MatrixX * W1st = new Amg::MatrixX(*W);
204 alignTrack->setWeightMatrixFirstDeriv(W1st);
205
206 return true;
207 }
208
209 // ========================
210 // setup global chi2 method
211 // ========================
212
213 // get the ingredients for calculating R = V - Q
214 const Amg::MatrixX * V = alignTrack->localErrorMatrix();
215 ATH_MSG_DEBUG("V ( "<<V->rows()<<" x "<<V->cols()<<" ) (diagonal only):");
216 if (msgLvl(MSG::DEBUG)) {
217 for (int i=0;i<V->rows();i++) msg()<<(*V)(i,i)<<" ";
218 msg()<<endmsg;
219 }
220
221 const int outputdim = alignTrack->nAlignTSOSMeas();
222
223 Amg::MatrixX Q(outputdim, outputdim); //symmetric matrix
224 if (!getTrkParamCovMatrix(alignTrack, Q))
225 return false;
226 ATH_MSG_DEBUG("Q ( "<<Q.rows()<<" x "<<Q.cols()<<" )");
227 ATH_MSG_DEBUG("Q: "<<Q);
228
229 // calculate R
230 const Amg::MatrixX R = (*V) - Q;
231 if (!checkValidity(R)) {
232 ATH_MSG_WARNING("Matrix R = V - HCH is invalid, skipping the track");
233 return false;
234 }
235 ATH_MSG_DEBUG("R ( "<<R.rows()<<" x "<<R.cols()<<" )");
236 ATH_MSG_DEBUG("R: "<<R);
237
238 // calculate weight matrix and set it
239 Amg::MatrixX * W = new Amg::MatrixX(outputdim,outputdim); //symmetric matrix
240 *W = ((*Vinv) * R * (*Vinv));
241
242 if (!checkValidity(*W)) {
243 ATH_MSG_DEBUG("Weight matrix is invalid, skipping the track");
244 delete W;
245 return false;
246 }
247 ATH_MSG_DEBUG("setting weight: "<<(*W));
248 alignTrack->setWeightMatrix(W);
249
250 // for 1st derivatives the weight matrix is just the V^-1
251 // it has been chacked above so no need to do it again
252 Amg::MatrixX * W1st = new Amg::MatrixX(*Vinv); //symmetric matrix
253 ATH_MSG_DEBUG("setting weight for 1st derivatives (diagonal only): ");
254 if (msgLvl(MSG::DEBUG)) {
255 for (int i=0;i<W1st->rows();i++) msg()<<(*W1st)(i,i)<<" ";
256 msg()<<endmsg;
257 }
258 alignTrack->setWeightMatrixFirstDeriv(W1st);
259
260 return true;
261 }
262
263 //________________________________________________________________________
264 bool AnalyticalDerivCalcTool::getTrkParamCovMatrix(const AlignTrack * alignTrack, Amg::MatrixX & Q /*symmetric matrix*/) const
265 {
266
267 // get derivative matrices from fitter
268 const Amg::MatrixX * H0 = alignTrack->derivativeMatrix();
269 const Amg::MatrixX * C = alignTrack->fullCovarianceMatrix(); //symmetric matrix
270
271 // H0 is a q0 x p matrix,
272 // C is a p x p matrix,
273 // q0 = number of tsos measurements
274 // p = number perigee params+2*nscat+nbrem
275
276 if( H0==nullptr || C==nullptr) {
277 ATH_MSG_ERROR("no derivative matrix or cov matrix stored on AlignTrack!"
278 << "This should have been done in AlignTrackPreProcessor!"
279 << H0 << " " << C );
280 return false;
281 }
282
283 ATH_MSG_DEBUG("H0 ( "<<H0->rows()<<" x "<<H0->cols()<<" )");
284 ATH_MSG_DEBUG("H0: "<<(*H0));
285 ATH_MSG_DEBUG("C ( "<<C->rows()<<" x "<<C->cols()<<" )");
286 ATH_MSG_DEBUG("C: "<<(*C));
287
288
289 int Csize(C->rows());
290 Amg::MatrixX CC(*C); // take a copy of C //symmetric matrix
291 int ierr(0);
292 bool amendC( !(alignTrack->refitD0() && alignTrack->refitZ0() &&
293 alignTrack->refitPhi() && alignTrack->refitTheta() &&
294 alignTrack->refitQovP()) );
295
296 // amendC = true; // FUDGE!
297
298 if( amendC ) {
299 // test with AlSymMat:
300 // build AlSymMat instance from
301 AlSymMat* CA = new AlSymMat(Csize);
302
303 // take a copy of C:
304 for( int ii=0; ii<Csize; ++ii ) {
305 for( int jj=ii; jj<Csize; ++jj ) {
306 CA->elemr(ii,jj) = (*C)(ii,jj);
307 }
308 }
309
310 // first inversion:
311 ierr = CA->invert();
312 if( ierr ) {
313 ATH_MSG_ERROR("First inversion of matrix CA failed with LAPACK status flag " << ierr);
314 return false;
315 } else {
316 // disable selected track parametrs (remove from refit).
317 // It is your duty to assure that corresponding constraints have been imposed!
318 if( !(alignTrack->refitD0()) ) {
319 for(int ii=0; ii<(Csize); ++ii) CA->elemr(0,ii)=0.0; // should do the trick. It is a CLHEP::HepSymMatrix!
320 CA->elemr(0,0)=1.0;
321 }
322 if( !(alignTrack->refitZ0()) ) {
323 for(int ii=0; ii<(Csize); ++ii) CA->elemr(1,ii)=0.0; // should do the trick. It is a CLHEP::HepSymMatrix!
324 CA->elemr(1,1)=1.0;
325 }
326 if( !(alignTrack->refitPhi()) ) {
327 for(int ii=0; ii<(Csize); ++ii) CA->elemr(2,ii)=0.0; // should do the trick. It is a CLHEP::HepSymMatrix!
328 CA->elemr(2,2)=1.0;
329 }
330 if( !(alignTrack->refitTheta()) ) {
331 for(int ii=0; ii<(Csize); ++ii) CA->elemr(3,ii)=0.0; // should do the trick. It is a CLHEP::HepSymMatrix!
332 CA->elemr(3,3)=1.0;
333 }
334 if( !(alignTrack->refitQovP()) ) {
335 for(int ii=0; ii<(Csize); ++ii) CA->elemr(4,ii)=0.0; // should do the trick. It is a CLHEP::HepSymMatrix!
336 CA->elemr(4,4)=1.0;
337 }
338
339
340 // invert back:
341 ierr = CA->invert();
342 if( ierr ) {
343 ATH_MSG_ERROR("Second inversion of matrix CA failed with LAPACK status flag " << ierr);
344 return false;
345 }
346
347 // copy back to CC:
348 for( int ii=0; ii<Csize; ++ii ) {
349 for( int jj=ii; jj<Csize; ++jj ) {
350 CC(ii,jj) = CA->elemc(ii,jj);
351 }
352 }
353
354 // clear the disabled rows/collumns
355 if( !(alignTrack->refitD0()) ) for(int ii=0; ii<(Csize); ++ii){ CC(0,ii)=0.0; CC(ii,0)=0.0; }; // should do the trick.
356 if( !(alignTrack->refitZ0()) ) for(int ii=0; ii<(Csize); ++ii){ CC(1,ii)=0.0; CC(ii,1)=0.0; }; // should do the trick.
357 if( !(alignTrack->refitPhi()) ) for(int ii=0; ii<(Csize); ++ii){ CC(2,ii)=0.0; CC(ii,2)=0.0; }; // should do the trick.
358 if( !(alignTrack->refitTheta()) ) for(int ii=0; ii<(Csize); ++ii){ CC(3,ii)=0.0; CC(ii,3)=0.0; }; // should do the trick.
359 if( !(alignTrack->refitQovP()) ) for(int ii=0; ii<(Csize); ++ii){ CC(4,ii)=0.0; CC(ii,4)=0.0; }; // should do the trick.
360
361
362 }
363
364 //garbage collection:
365 delete CA;
366 }
367
368
369
370 int nMeas = H0->rows();
371 int nAtsos = alignTrack->nAlignTSOSMeas();
372 ATH_MSG_DEBUG("nMeas: "<<nMeas);
373
374 //CLHEP::HepSymMatrix HCH(nMeas,0);
375 Amg::MatrixX HCH(nAtsos,nAtsos);
376 HCH = CC.similarity( *H0 );
377
378 ATH_MSG_DEBUG("HCH ( "<<HCH.rows()<<" x "<<HCH.cols()<<" )");
379 ATH_MSG_DEBUG("HCH: "<<HCH);
380
381 //
382 // get indices of HCH matrix corresponding to alignTSOSs in alignTrack
383 const AlignTSOSCollection * alignTSOSCollection = alignTrack->alignTSOSCollection();
384
385 std::vector<int> matrixIndices(nMeas);
386
387 int imeas(-1);
388 for (const TrackStateOnSurface* tsos : *alignTrack->trackStateOnSurfaces()){
389
390 ATH_MSG_DEBUG("tsos: "<<tsos->dumpType());
391
392 // get tsos and make sure it is a RIO_OnTrack
393 if (tsos->type(TrackStateOnSurface::Outlier))
394 continue;
395
396 // RIO
397 if (!tsos->type(TrackStateOnSurface::Scatterer)) {
398 ATH_MSG_DEBUG("not scatterer, trying rio");
399
400 const MeasurementBase * mesb = tsos->measurementOnTrack();
401 const RIO_OnTrack * rio = dynamic_cast<const RIO_OnTrack *>(mesb);
402 const CompetingRIOsOnTrack * crio = dynamic_cast<const CompetingRIOsOnTrack *>(mesb);
403 if (!rio && crio)
404 rio=&crio->rioOnTrack(0);
405
406 if (rio==nullptr)
407 continue;
408
409 ++imeas;
410 matrixIndices[imeas]=-1;
411
412 Identifier tsosId = rio->identify();
413 ATH_MSG_DEBUG("have tsos with Id "<<tsosId);
414
415 // get matching alignTSOS and store track index in goodMatrixIndices
416 int iameas(0);
417 for (const AlignTSOS* atsos : *alignTSOSCollection) {
418
419 if (!atsos->isValid())
420 continue;
421
422 if (atsos->type(TrackStateOnSurface::Scatterer))
423 ATH_MSG_ERROR("can't use scatterers on AlignTrack yet for analytical derivatives!");
424
425 const RIO_OnTrack * atsos_rio = atsos->rio();
426 if (atsos_rio) {
427// ATH_MSG_DEBUG("tsosId / atsosId : "<<tsosId<<" / "<<atsos_rio->identify());
428
429 if (atsos_rio->identify()==tsosId) {
430 matrixIndices[imeas]=iameas;
431 ATH_MSG_DEBUG("matrixIndices["<<imeas<<"]="<<iameas);
432
433 // for Pixel we have two measurements
434 if (atsos->nResDim()>1)
435 {
436 imeas++;
437 iameas++;
438 matrixIndices[imeas]=iameas;
439 ATH_MSG_DEBUG("matrixIndices["<<imeas<<"]="<<iameas);
440 }
441 break;
442 }
443 // for Pixel we have two measurements
444 else if(atsos->nResDim()>1)
445 iameas++;
446 }
447 iameas++;
448 }
449
450 // even when the Pixel is not aligned we have to take into account
451 // that it has two measurements
452 if(matrixIndices[imeas]==-1 && m_measTypeIdHelper->defineType(mesb) == TrackState::Pixel) {
453 imeas++;
454 matrixIndices[imeas]=-1;
455 ATH_MSG_DEBUG("matrixIndices["<<imeas<<"]="<<matrixIndices[imeas]);
456 }
457
458 } // end tsos==rio
459
460 } // end loop over track tsos
461
462 // strip elements in the HCH matrix which don't correspond
463 // to AlignTSOS on alignTrack
464 ATH_MSG_DEBUG("Filling the Q matrix:");
465 for (int k=0;k<nMeas;k++) {
466
467 int iameas=matrixIndices[k];
468 if (iameas==-1)
469 continue;
470
471 for (int l=0;l<nMeas;l++) {
472 int jameas=matrixIndices[l];
473 if (jameas==-1)
474 continue;
475
476 Q(iameas,jameas) = HCH(k,l);
477
478 }
479 }
480
481 ATH_MSG_DEBUG("before check Q ( "<<Q.rows()<<" x "<<Q.cols()<<" )");
482 ATH_MSG_DEBUG("before check Q: "<<Q);
483
484 if(!checkValidity(Q)) {
485 ATH_MSG_DEBUG("Matrix Q = HCH is invalid, skipping the track");
486 return false;
487 }
488
489 return true;
490 }
491
492 //________________________________________________________________________
493 // is this method used used anywhere? maybe we should remove it
494 bool AnalyticalDerivCalcTool::getMeasErrorMatrix(const AlignTrack* alignTrack, Amg::MatrixX& V /*Symmetric Matrix*/) const
495 {
496 int index(0);
498 for (; itAtsos != alignTrack->lastAtsos(); ++itAtsos) {
499
500 std::vector<Residual>::const_iterator itRes=(**itAtsos).firstResidual();
501 for (; itRes!=(**itAtsos).lastResidual(); ++itRes,index++) {
502
503 V(index,index) = itRes->errSq();
504 }
505 }
506 return checkValidity(V);
507 }
508
509 //________________________________________________________________________
510 bool AnalyticalDerivCalcTool::checkValidity(const Amg::MatrixX& R /*symmetric matrix*/) const
511 {
512 // perform some sort of sanity check. depending on how many hits
513 // on the track we use, R can have zero determinant, but it
514 // should definitely not be negative. for now, we'll just check
515 // that all diagonal elements are positive and that all correlation
516 // coefficients are within bounds
517
518 bool Risvalid(true);
519 const double epsilon=1e-10;
520 for( int irow=0; irow<R.rows(); ++irow) {
521
522 Risvalid = Risvalid && R(irow,irow)>0;
523 if ( msgLvl(MSG::DEBUG) ) {
524 if( !(R(irow,irow)>0) )
525 msg(MSG::DEBUG) << "matrix invalid: (" << irow << "," << irow<<") = " << R(irow,irow) << endmsg;
526 }
527 else if (!Risvalid)
528 break;
529
530 for(int icol=0; icol<=irow; ++icol) {
531 // this one must be true if everything else succeeded
532 double Rcorr = R(irow,icol)/sqrt(R(irow,irow)*R(icol,icol));
533 if( Rcorr+epsilon<-1 || Rcorr-epsilon>1 )
534 {
535 Risvalid = false;
536 if (msgLvl(MSG::DEBUG))
537 ATH_MSG_DEBUG("matrix corr invalid for (" << irow << "," << icol << ") Rcorr = " << Rcorr);
538 else
539 break;
540 }
541 }
542 }
543
544 if( !Risvalid ) {
545 ATH_MSG_WARNING("Checked matrix is invalid.");
546 ATH_MSG_WARNING("R: \n"<<R);
547 }
548 return Risvalid;
549 }
550
551 //________________________________________________________________________
552 std::vector<Amg::VectorX> AnalyticalDerivCalcTool::getDerivatives(AlignTrack * alignTrack, const AlignModule * module)
553 {
554 // module-specific transforms
555 Amg::Transform3D globalFrameToAlignFrame = module->globalFrameToAlignFrame();
556 ATH_MSG_DEBUG("globalFrameToAlignFrame: ");
557 ATH_MSG_DEBUG(globalFrameToAlignFrame(0,0)<<" "<<
558 globalFrameToAlignFrame(0,1)<<" "<<
559 globalFrameToAlignFrame(0,2));
560 ATH_MSG_DEBUG(globalFrameToAlignFrame(1,0)<<" "<<
561 globalFrameToAlignFrame(1,1)<<" "<<
562 globalFrameToAlignFrame(1,2));
563 ATH_MSG_DEBUG(globalFrameToAlignFrame(2,0)<<" "<<
564 globalFrameToAlignFrame(2,1)<<" "<<
565 globalFrameToAlignFrame(2,2));
566
567 Amg::RotationMatrix3D globalToAlignFrameRotation = module->globalToAlignFrameRotation();
568 ATH_MSG_DEBUG("globalToAlignFrameRotation: ");
569 ATH_MSG_DEBUG(globalToAlignFrameRotation(0,0)<<" "<<
570 globalToAlignFrameRotation(0,1)<<" "<<
571 globalToAlignFrameRotation(0,2));
572 ATH_MSG_DEBUG(globalToAlignFrameRotation(1,0)<<" "<<
573 globalToAlignFrameRotation(1,1)<<" "<<
574 globalToAlignFrameRotation(1,2));
575 ATH_MSG_DEBUG(globalToAlignFrameRotation(2,0)<<" "<<
576 globalToAlignFrameRotation(2,1)<<" "<<
577 globalToAlignFrameRotation(2,2));
578
579 DataVector<AlignPar> * alignPars = m_alignModuleTool->getAlignPars(module);
580 const int nAlignPar = alignPars->size();
581
582 //Create derivatives storage vector and initialise them ro zero
583 std::vector<Amg::VectorX> derivatives( nAlignPar+3 , Amg::VectorX(alignTrack->nAlignTSOSMeas()));
584 for(int i(0); i<nAlignPar+3; ++i) derivatives[i].setZero();
585
586 int imeas(0);
587 AlignTSOSCollection::iterator iatsos = alignTrack->firstAtsos();
588 for (; iatsos != alignTrack->lastAtsos(); ++iatsos) {
589
590 AlignTSOS * alignTSOS = *iatsos;
591 if (!alignTSOS->isValid() || nullptr==alignTSOS->module())
592 continue;
593
594 // we only calculate the derivatives if the AlignTSOS belongs to the align module
595 int nResDim = alignTSOS->nResDim();
596 if (alignTSOS->module() != module) {
597 imeas += nResDim;
598 continue;
599 }
600
601 // derivatives to be stored on the AlignTSOS
602 std::vector<Amg::VectorX> * atsosDerivs = nullptr;
603 std::vector<Amg::VectorX> * atsosDerVtx = nullptr;
604 if (m_storeDerivatives) {
605 atsosDerivs = new std::vector<Amg::VectorX>(nResDim,Amg::VectorX(nAlignPar));
606 atsosDerVtx = new std::vector<Amg::VectorX>(nResDim,Amg::VectorX(3));
607 ATH_MSG_DEBUG("nResDim = "<<nResDim<<" vector size is "<<atsosDerivs->size());
608 ATH_MSG_DEBUG("nAlignPar = "<<nAlignPar<<" CLHEP::HepVector size is "<<atsosDerivs->at(0).rows());
609 }
610
611 // Get the rotation to the frame in which the residual is
612 // defined. In this frame the trkdistance is the x-coordinate.
613
614 const TrackParameters * mtp = alignTSOS->trackParameters();
615 if (!mtp || !(mtp->covariance()) ) continue;
616 Amg::RotationMatrix3D localToGlobalRotation = mtp->measurementFrame();
617
618
619 ATH_MSG_DEBUG( "localToGlobalRotation:");
620 ATH_MSG_DEBUG(localToGlobalRotation(0,0) << " " <<
621 localToGlobalRotation(0,1) << " " <<
622 localToGlobalRotation(0,2));
623 ATH_MSG_DEBUG(localToGlobalRotation(1,0) << " " <<
624 localToGlobalRotation(1,1) << " " <<
625 localToGlobalRotation(1,2));
626 ATH_MSG_DEBUG(localToGlobalRotation(2,0) << " " <<
627 localToGlobalRotation(2,1) << " " <<
628 localToGlobalRotation(2,2));
629
630 if(double alphastrip=alignTSOS->alphaStrip()) {
631 ATH_MSG_DEBUG( "applying fanout rotation : " << alphastrip );
632 localToGlobalRotation = localToGlobalRotation * Amg::AngleAxis3D(alphastrip, Amg::Vector3D(0.,0.,1.));
633 ATH_MSG_DEBUG( "localToGlobalRotation * fanout_rotation:");
634 ATH_MSG_DEBUG(localToGlobalRotation(0,0) << " " <<
635 localToGlobalRotation(0,1) << " " <<
636 localToGlobalRotation(0,2));
637 ATH_MSG_DEBUG(localToGlobalRotation(1,0) << " " <<
638 localToGlobalRotation(1,1) << " " <<
639 localToGlobalRotation(1,2));
640 ATH_MSG_DEBUG(localToGlobalRotation(2,0) << " " <<
641 localToGlobalRotation(2,1) << " " <<
642 localToGlobalRotation(2,2));
643 }
644
645 // get the position of the track in the alignmentframe.
646 Amg::Vector3D refPos = globalFrameToAlignFrame * alignTSOS->trackParameters()->position();
647 ATH_MSG_DEBUG("refPos: "<<refPos);
648
649
650 const Amg::RotationMatrix3D R = globalToAlignFrameRotation * localToGlobalRotation;
651 ATH_MSG_DEBUG("R:");
652 ATH_MSG_DEBUG(R(0,0) << " " << R(0,1) << " " << R(0,2));
653 ATH_MSG_DEBUG(R(1,0) << " " << R(1,1) << " " << R(1,2));
654 ATH_MSG_DEBUG(R(2,0) << " " << R(2,1) << " " << R(2,2));
655
656 // In the SCT measurement frame:
657 // x --> perpendicular to strips in wafer plane
658 // y --> along strips in wafer plane
659 // z --> perpendicular to wafer plane
660
661 // In the TRT measurement frame:
662 // x --> perpendicular to track and straw
663 // y --> along straw wire
664 // z --> perpendicular to x and y (but not parallel to track!)
665
666 // now 'correct' for the track angle in the measurement frame.
667 const TrackParameters * trkpars = nullptr;
668 if(m_residualType == HitOnly) {
669 ATH_MSG_DEBUG("using BIASED track parameters");
670 trkpars = alignTSOS->trackParameters();
671 }
672 else {
673 ATH_MSG_DEBUG("using UNBIASED track parameters");
674 trkpars = alignTSOS->unbiasedTrackPars();
675 }
676
677
678 Amg::Vector3D trackdir = localToGlobalRotation.inverse() * trkpars->momentum();
679 ATH_MSG_DEBUG( "trackdir " << trackdir[0] << " " << trackdir[1] << " " << trackdir[2]);
680
681 // for 1-dimensional measurements and Pixel-x
682 ATH_MSG_DEBUG( "trackdir.z(): " << trackdir.z() );
683 double cotphi_x = trackdir.x() / trackdir.z();
684
685 // some 1D measurements are in Y direction (e.g. in CSC)
686 // so we need the other angle
687 if (alignTSOS->measDir() == Trk::y)
688 cotphi_x = trackdir.y() / trackdir.z();
689
690
691 double Rxx = R(0,0) - cotphi_x * R(0,2);
692 double Ryx = R(1,0) - cotphi_x * R(1,2);
693 double Rzx = R(2,0) - cotphi_x * R(2,2);
694 ATH_MSG_DEBUG("Rxx/Ryx/Rzx: " << Rxx << "/" << Ryx << "/" << Rzx);
695
696 double projR[AlignModule::NTransformPar];
697 projR[AlignModule::TransX] = Rxx;
698 projR[AlignModule::TransY] = Ryx;
699 projR[AlignModule::TransZ] = Rzx;
700 projR[AlignModule::RotX] = -Ryx * refPos.z() + Rzx * refPos.y();
701 projR[AlignModule::RotY] = -Rzx * refPos.x() + Rxx * refPos.z();
702 projR[AlignModule::RotZ] = -Rxx * refPos.y() + Ryx * refPos.x();
703 projR[AlignModule::BowX] = 0;
704 projR[AlignModule::BowY] = 0;
705 projR[AlignModule::BowZ] = 0;
706
716
725
726
728 const double localz = alignTSOS->trackParameters()->position().z(); // - globalToAlignFrameTranslation().z(); // the last term to be doublechecked!
729 // stave length in the IBL -- we will see if there is a more generic way of doing this
730 const double z0z0 = 366.5*366.5;
731
732 projR[AlignModule::BowX] = ( localz*localz - z0z0) / z0z0; // this formula should work for both L11 ans L16, sign to be checked!
733
734
735 // prepare derivatives w.r.t. the vertex position:
736 Amg::Vector3D RxLoc(Rxx, Ryx, Rzx);
737 Amg::Vector3D RxGlob=-1.0 * (globalToAlignFrameRotation.inverse() * RxLoc); // to be double checked!!!
738
739 for (int ipar=0; ipar<nAlignPar; ipar++) {
740 const AlignPar * alignPar = (*alignPars)[ipar];
741 int paramType = alignPar->paramType();
742 //double sigma = alignPar->sigma();
743 ATH_MSG_DEBUG("ipar="<<ipar<<", paramType="<<paramType);
744 derivatives[ipar][imeas] = projR[paramType];//*sigma;
746 (*atsosDerivs)[0][ipar] = projR[paramType];//*sigma;
747 }
748 // the dr/db bit:
749 for (int ipar=0; ipar<3; ipar++) {
750 derivatives[nAlignPar+ipar][imeas] = RxGlob[ipar];
751 if (m_storeDerivatives) (*atsosDerVtx)[0][ipar] = RxGlob[ipar];
752 }
753
754 for (int i=0;i<nAlignPar+3;i++)
755 ATH_MSG_DEBUG("derivatives["<<i<<"]["<<imeas<<"]="<<derivatives[i][imeas]);
756
757 imeas++;
758
759 if (nResDim>1) {
760 // for Pixel the second measurement has to be corrected
761 // for the second angle
762 double cotphi_y = trackdir.y() / trackdir.z() ;
763 double Rxy = R(0,1) - cotphi_y * R(0,2) ;
764 double Ryy = R(1,1) - cotphi_y * R(1,2) ;
765 double Rzy = R(2,1) - cotphi_y * R(2,2) ;
766 ATH_MSG_DEBUG("Rxy/Ryy/Rzy: " << Rxy << "/" << Ryy << "/" << Rzy);
767
768 projR[AlignModule::TransX] = Rxy;
769 projR[AlignModule::TransY] = Ryy;
770 projR[AlignModule::TransZ] = Rzy;
771 projR[AlignModule::RotX] = -Ryy * refPos.z() + Rzy * refPos.y();
772 projR[AlignModule::RotY] = -Rzy * refPos.x() + Rxy * refPos.z();
773 projR[AlignModule::RotZ] = -Rxy * refPos.y() + Ryy * refPos.x();
774
775 //Possibly could add the bowing correction -- very weakly coupled to y residuals
776 projR[AlignModule::BowX] = 0;
777 projR[AlignModule::BowY] = 0;
778 projR[AlignModule::BowZ] = 0;
779
780
781 // prepare derivatives w.r.t. the vertex position:
782 Amg::Vector3D RyLoc(Rxy, Ryy, Rzy);
783 Amg::Vector3D RyGlob=-1.0 * (globalToAlignFrameRotation.inverse() * RyLoc); // to be double checked!!!
784
785 for (int ipar=0; ipar<nAlignPar; ipar++) {
786 const AlignPar * alignPar = (*alignPars)[ipar];
787 int paramType = alignPar->paramType();
788 ATH_MSG_DEBUG("2nd dim, ipar="<<ipar<<", paramType="<<paramType);
789 //double sigma=alignPar->sigma();
790 derivatives[ipar][imeas] = projR[paramType];//*sigma;
792 (*atsosDerivs)[1][ipar] = projR[paramType];//*sigma;
793 }
794
795 // the dr/db bit:
796 for (int ipar=0; ipar<3; ipar++) {
797 derivatives[nAlignPar+ipar][imeas] = RyGlob[ipar];
798 if (m_storeDerivatives) (*atsosDerVtx)[1][ipar] = RyGlob[ipar];
799 }
800
801 for (int i=0;i<nAlignPar+3;i++)
802 ATH_MSG_DEBUG("2nd dim: derivatives["<<i<<"]["<<imeas<<"]="<<derivatives[i][imeas]);
803
804 imeas++;
805 }
806
807 alignTSOS->setDerivatives(atsosDerivs);
808 alignTSOS->setDerivativesVtx(atsosDerVtx);
809 }
810 ATH_MSG_DEBUG("returning derivatives");
811 return derivatives;
812 }
813
814 //________________________________________________________________________
816 {
817 // get first AlignTSOS of the AlignTrack
818 // this assumes that for unbiased or DCA residuals the scatterers
819 // and energy deposits are not included in the AlignTSOSSollection
820 const AlignTSOS * atsos = *(alignTrack->firstAtsos());
821
822 // get residual type of the first residual
823 m_residualType = atsos->firstResidual()->residualType();
824 ATH_MSG_DEBUG("setting residualType to "<<m_residualType);
825
826 m_residualTypeSet = true;
827 }
828
829} // end namespace
#define endmsg
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
const ServiceHandle< StoreGateSvc > & detStore() const
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
Derived DataVector<T>.
Definition DataVector.h:795
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
DataModel_detail::iterator< DataVector > iterator
Definition DataVector.h:842
size_type size() const noexcept
Returns the number of elements in the collection.
contains the implementation for handling symmetric matrices in triangular representation
Definition AlSymMat.h:26
AlignModule::TransformParameters paramType() const
returns the type of parameter (i.e.
Definition AlignPar.h:47
double alphaStrip() const
returns strip angle for fan-out structured modules (SCT endcap)
Definition AlignTSOS.h:122
const TrackParameters * unbiasedTrackPars() const
returns pointer to unbiased track parameters if present
Definition AlignTSOS.h:134
Trk::ParamDefs measDir() const
retrieve the measurement direction
Definition AlignTSOS.h:128
int nResDim() const
returns number of measurement residual + scatterer residual dimensions
Definition AlignTSOS.h:77
void setDerivativesVtx(std::vector< Amg::VectorX > *derivs)
setter for the derivatives w.r.t.
Definition AlignTSOS.h:149
std::vector< Residual >::const_iterator firstResidual() const
returns first Residual iterator
Definition AlignTSOS.h:104
const AlignModule * module() const
accessor method for AlignModule to which tsos belongs
Definition AlignTSOS.h:69
void setDerivatives(std::vector< Amg::VectorX > *derivs)
setter for the derivatives
Definition AlignTSOS.h:146
bool isValid() const
Definition AlignTSOS.h:74
TrackState::MeasurementType measType() const
returns measurement type enum
Definition AlignTSOS.h:80
const Amg::SymMatrixX * fullCovarianceMatrix() const
set and get full covariance matrix
Definition AlignTrack.h:177
bool refitZ0() const
Definition AlignTrack.h:194
void setWeightMatrix(Amg::SymMatrixX *mat)
Definition AlignTrack.h:157
AlignTSOSCollection::const_iterator lastAtsos() const
returns iterator pointer to last element in collection
Definition AlignTrack.h:280
const Amg::SymMatrixX * localErrorMatrixInv() const
inverse local error matrix, calculated by AlignTrack by calling atsos->hitDistanceVar()
Definition AlignTrack.h:125
int nAlignTSOSMeas() const
number of alignTSOS (including scatterers if included on AlignTrack
Definition AlignTrack.h:128
const Amg::SymMatrixX * localErrorMatrix() const
local error matrix, calculated by AlignTrack by calling atsos->hitDistanceVar()
Definition AlignTrack.h:122
bool refitQovP() const
Definition AlignTrack.h:197
void setDerivatives(std::vector< AlignModuleDerivatives > *vec)
Definition AlignTrack.h:134
AlignTSOSCollection::const_iterator firstAtsos() const
retrieve iterator pointer to first element in collection
Definition AlignTrack.h:279
void setWeightMatrixFirstDeriv(Amg::SymMatrixX *mat)
Definition AlignTrack.h:162
const AlignTSOSCollection * alignTSOSCollection() const
returns collection of alignTSOS
Definition AlignTrack.h:267
const Amg::MatrixX * derivativeMatrix() const
set and get derivative matrix
Definition AlignTrack.h:173
bool refitPhi() const
Definition AlignTrack.h:195
bool refitTheta() const
Definition AlignTrack.h:196
bool refitD0() const
get refit flags
Definition AlignTrack.h:193
PublicToolHandle< IAlignModuleTool > m_alignModuleTool
bool checkValidity(const Amg::MatrixX &R) const
bool getMeasErrorMatrix(const AlignTrack *alignTrack, Amg::MatrixX &V) const
AnalyticalDerivCalcTool(const std::string &type, const std::string &name, const IInterface *parent)
bool getTrkParamCovMatrix(const AlignTrack *alignTrack, Amg::MatrixX &HCH) const
int m_residualType
residual type to be used in the calculations
bool setResidualCovMatrix(AlignTrack *alignTrack) const override
sets residual covariance matrix
bool m_residualTypeSet
do we have the residual type set?
bool setDerivatives(AlignTrack *alignTrack) override
sets analytical partial derivatives of residuals w.r.t alignment parameters for TSOS on alignTrack.
void checkResidualType(const AlignTrack *alignTrack)
std::vector< Amg::VectorX > getDerivatives(AlignTrack *alignTrack, const AlignModule *module)
Base class for all CompetingRIOsOnTack implementations, extends the common MeasurementBase.
virtual const RIO_OnTrack & rioOnTrack(unsigned int) const =0
returns the RIO_OnTrack (also known as ROT) objects depending on the integer.
std::ostream * m_logStream
logfile output stream
This class is the pure abstract base class for all fittable tracking measurements.
classifies a MeasurementBase into one of the known inherited flavours or one of the detector types fo...
const Amg::Vector3D & momentum() const
Access method for the momentum.
virtual Amg::RotationMatrix3D measurementFrame() const override=0
Return the measurement frame - this is needed for alignment, in particular for StraightLine and Perig...
const Amg::Vector3D & position() const
Access method for the position.
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
Identifier identify() const
return the identifier -extends MeasurementBase
represents the track state (measurement, material, fit parameters and quality) at a surface.
const TrackParameters * trackParameters() const
return ptr to trackparameters const overload
@ Outlier
This TSoS contains an outlier, that is, it contains a MeasurementBase/RIO_OnTrack which was not used ...
@ Scatterer
This represents a scattering point on the track, and so will contain TrackParameters and MaterialEffe...
const Trk::TrackStates * trackStateOnSurfaces() const
return a pointer to a const DataVector of const TrackStateOnSurfaces.
struct color C
Eigen::AngleAxisd AngleAxis3D
Eigen::Matrix< double, 3, 3 > RotationMatrix3D
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.
Ensure that the ATLAS eigen extensions are properly loaded.
DataVector< AlignTSOS > AlignTSOSCollection
Definition AlignTrack.h:37
@ y
Definition ParamDefs.h:56
ParametersBase< TrackParametersDim, Charged > TrackParameters
Definition index.py:1