ATLAS Offline Software
Loading...
Searching...
No Matches
ShiftingDerivCalcTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "TrkTrack/Track.h"
13
19
22
24
25#include "TGraph.h"
26#include "TF1.h"
27#include "TCanvas.h"
28#include "TPaveText.h"
29#include "TAxis.h"
30
31
32namespace Trk {
33
34 //________________________________________________________________________
36 const std::string& name,
37 const IInterface* parent)
38
39 : AthAlgTool(type,name,parent)
40 , m_unshiftedTrackChi2MeasType(std::make_unique<double[]>(TrackState::NumberOfMeasurementTypes))
41 {
42 declareInterface<IDerivCalcTool>(this);
43
44 m_logStream = nullptr;
45 }
46
47 //________________________________________________________________________
52
53 //________________________________________________________________________
55 {
56
57 msg(MSG::DEBUG) << "in ShiftingDerivCalcTool initialize()"<<endmsg;
58 ATH_CHECK(m_trackFitterTool.retrieve());
61 ATH_CHECK(m_alignModuleTool.retrieve());
62
64 msg(MSG::INFO) << "ParticleNumber: " << m_particleNumber << endmsg;
65 msg(MSG::INFO) << "ParticleHypothesis: " << m_particleHypothesis << endmsg;
66
67
68 if(!m_doFits){
69 m_nFits = 2;
70 }
71
72 m_traSize = 5.*m_traSize/(double)m_nFits;
73 m_rotSize = 5.*m_rotSize/(double)m_nFits;
74
75 msg(MSG::INFO) << "doFits: " << m_doFits << endmsg;
76 msg(MSG::INFO) << "nFits: " << m_nFits << endmsg;
77
78 return StatusCode::SUCCESS;
79 }
80
81 //________________________________________________________________________
83 {
84 ATH_MSG_INFO("number tracks processed: "<<m_ntracksProcessed<<
85 "\nnumber tracks passing initial scan: "<<m_ntracksPassInitScan<<
86 "\nnumber tracks passing setting unshifted residuals: "<< m_ntracksPassSetUnshiftedRes<<
87 "\nnumber tracks pass getting derivatives (1st pass): "<<m_ntracksPassGetDeriv<<
88 "\nnumber tracks pass getting derivatives (2nd pass): "<<m_ntracksPassGetDerivSecPass<<
89 "\nnumber tracks pass getting derivatives (3rd pass): "<<m_ntracksPassGetDerivLastPass<<
90 "\nnumber tracks pass setting derivatives: "<<m_ntracksPassDerivatives);
91 ATH_MSG_INFO("number tracks fail max iterations: "<<m_ntracksFailMaxIter<<
92 "\nnumber tracks fail track refit: "<<m_ntracksFailTrackRefit<<
93 "\nnumber tracks fail align param cut: "<<m_ntracksFailAlignParamCut<<
94 "\nnumber tracks fail final attempt: "<<m_ntracksFailFinalAttempt);
95
96 return StatusCode::SUCCESS;
97 }
98
99 //________________________________________________________________________
101 const std::vector<AlignModule*>& alignModules)
102 {
103 ATH_MSG_DEBUG("in scanShifts");
104
105 const Trk::Track* trackForRefit =
107 dynamic_cast<const Trk::Track*>(alignTrack);
108
109 // see whether straight track or not
110 m_fitter = alignTrack->isSLTrack() ?
112 ATH_MSG_DEBUG("refitting unshifted track with "<<m_fitter<<" (isSLTrack="
113 <<alignTrack->isSLTrack()<<")");
114
115 ATH_MSG_DEBUG("setting minNIterations to "<<m_nIterations);
116
117 // refit track
118
120 alignCache.m_minIterations = m_nIterations;
121 const Track* refittedTrack = m_fitter->alignmentFit(alignCache, *trackForRefit,
124 if (!refittedTrack) {
125 msg(MSG::WARNING) << "initial track refit failed" << endmsg;
126 return false;
127 }
128 else
129 ATH_MSG_DEBUG("initial track refit successful");
130
133 ATH_MSG_DEBUG("exceeded maximum number of iterations");
134 return false;
135 }
136 ATH_MSG_DEBUG("initial nIterations: "<<m_nIterations);
137
138 // loop over AlignModules
139 for (std::vector<AlignModule*>::const_iterator moduleIt=alignModules.begin();
140 moduleIt!=alignModules.end(); ++moduleIt) {
141
142 // loop over AlignPar
143 DataVector<AlignPar>* alignPars=m_alignModuleTool->getAlignPars(*moduleIt);
144 for (DataVector<AlignPar>::iterator alignParIt=alignPars->begin();
145 alignParIt!=alignPars->end(); ++alignParIt) {
146
147 for (int ishift=0;ishift<2;ishift++) {
148
149 double shiftsize = shiftSize(*alignParIt);
150 if (ishift>0) shiftsize*=-1.;
151 m_alignModuleTool->shiftModule(*moduleIt,alignTrack,(**alignParIt).paramType(),shiftsize);
152 refittedTrack = (m_fitter->fit(Gaudi::Hive::currentContext(),
153 *trackForRefit,m_runOutlierRemoval,
154 m_particleHypothesis)).release();
155 m_alignModuleTool->restoreModule(*moduleIt);
156 if (!refittedTrack) {
157 msg(MSG::WARNING) << "track refit failed!"<<endmsg;
159 return false;
160 }
161
162 int nIter=alignCache.m_iterationsOfLastFit;
163 ATH_MSG_DEBUG("nIter: "<<nIter);
164 if (nIter>m_maxIter) {
165 ATH_MSG_DEBUG("exceeded maximum number of iterations");
167 return false;
168 }
169
170 if (nIter>m_nIterations) m_nIterations=nIter;
171 }
172 } // loop over AlignPar
173 } // loop over AlignModules
174
175 ATH_MSG_DEBUG("done with scanShifts, m_nIterations="<<m_nIterations);
176 return true;
177 }
178
179 //________________________________________________________________________
181 {
182
183 // see whether straight track or not
184 m_fitter = alignTrack->isSLTrack() ?
186 ATH_MSG_DEBUG("refitting unshifted track with "<<m_fitter<<" (isSLTrack="
187 <<alignTrack->isSLTrack()<<")");
188
189 // refit track
190 ATH_MSG_DEBUG("\nsetting min number iterations to "<<m_nIterations);
192 alignCache.m_minIterations = m_nIterations;
193
194 const Trk::Track* trackForRefit =
196 dynamic_cast<const Trk::Track*>(alignTrack);
197 if (!trackForRefit) {
198 ATH_MSG_ERROR("no track for refit!");
199 return false;
200 }
201
202 const Track* refittedTrack = m_fitter->alignmentFit( alignCache,
203 *trackForRefit,
206
207 if (!refittedTrack) {
208 ATH_MSG_WARNING( "initial track refit failed" );
209 return false;
210 }
211 else
212 ATH_MSG_DEBUG("initial track refit successful");
213
214 // dump local track chi2 for debugging
215 double localChi2=m_residualCalculator->setResiduals(alignTrack,refittedTrack);
216 msg()<<MSG::DEBUG<<"local Chi2(unshifted) in setChi2VAlignParam="<<localChi2<<endmsg;
217 m_unshiftedTrackChi2 = localChi2;
218 for (int i=0;i<TrackState::NumberOfMeasurementTypes;i++) {
219 ATH_MSG_DEBUG("getting chi2 for measType "<<i);
221 }
222 ATH_MSG_DEBUG("done");
223
224 // create vector containing unshifted residuals and matrices containing errors
225 const int NMEAS=alignTrack->nAlignTSOSMeas();
226
227 // unshiftedResiduals owned by AlignTrack
229
230 // unshiftedResErrors owned by ShiftingDerivCalcTool
233
234 // loop over atsos and determine residuals and errors
235 int imeas=0;
236 AlignTSOSCollection::const_iterator atsosItr=alignTrack->firstAtsos();
237 for (; atsosItr != alignTrack->lastAtsos(); ++atsosItr) {
238 if (!(**atsosItr).isValid()) continue;
239 for (std::vector<Residual>::const_iterator itRes=(**atsosItr).firstResidual();
240 itRes!=(**atsosItr).lastResidual();++itRes,++imeas) {
241 double residual = itRes->residual();
242 double errSq = itRes->errSq();
243 (*m_unshiftedResiduals)[imeas]=residual;
244 (*m_unshiftedResErrors)[imeas]=std::sqrt(errSq);
245 //ATH_MSG_DEBUG("weight: "<<1./errSq<<", unshiftedRes["<<imeas<<"]="
246 // <<(*m_unshiftedResiduals)[imeas]
247 // <<", resNorm="<<itRes->residualNorm());
248 }
249 }
250 if (imeas!=NMEAS) {
251 msg(MSG::ERROR)<<"problem with nmeas, imeas="<<imeas<<", NMEAS="<<NMEAS<<endmsg;
252 throw std::runtime_error("Error in ShiftingDerivCalcTool::setUnshiftedResiduals");
253 }
255
256 delete refittedTrack; refittedTrack=nullptr;
257
258 return true;
259 }
260
261//________________________________________________________________________
263{
264 ATH_MSG_DEBUG("in ShiftingDerivCalcTool setDerivatives");
266
267 // loop over AlignTSOSCollection,
268 // find modules that are in the AlignModuleList,
269 std::vector<AlignModule*> alignModules;
270 for (AlignTSOSCollection::iterator atsosItr=alignTrack->firstAtsos();
271 atsosItr != alignTrack->lastAtsos(); ++atsosItr) {
272
273 ATH_MSG_VERBOSE("getting module");
274 AlignModule* module=(*atsosItr)->module();
275 if (module)
276 ATH_MSG_VERBOSE("have ATSOS for module "<<module->identify());
277 else
278 ATH_MSG_VERBOSE("no module!");
279
280 if (!(*atsosItr)->isValid() || !module) continue;
281 if (find(alignModules.begin(),alignModules.end(),module) == alignModules.end())
282 alignModules.push_back(module);
283 }
284
285 // find perigee of best track fit and use as starting perigee for all fits
287 if (m_setMinIterations && !scanShifts(alignTrack, alignModules)) {
288 return false;
289 };
290
292
293 // set unshifted residuals (this is done in AlignTrackDresser but redone here with track refit)
294 if (!setUnshiftedResiduals(alignTrack)) {
295 ATH_MSG_WARNING("problem with refitting track!");
296 return false;
297 };
298
300
301 // Determine derivatives from shifting these modules
302 std::vector<AlignModuleDerivatives> * derivatives = new std::vector<AlignModuleDerivatives>;
303 std::vector<AlignModuleDerivatives> * derivativeErr = new std::vector<AlignModuleDerivatives>;
304 std::vector<std::pair<AlignModule*, std::vector<double> > > * actualSecondDerivatives =
305 new std::vector<std::pair<AlignModule*, std::vector<double> > >;
307 for (auto *alignModule : alignModules) {
308
309 ATH_MSG_DEBUG("finding derivatives for module "<<(*alignModule).identify());
310
311 std::vector<Amg::VectorX> deriv_vec;
312 std::vector<Amg::VectorX> derivErr_vec;
313 std::vector<double> actualsecderiv_vec;
314
315 // get alignPars and create arrays to store chi2 vs. align pars
316 DataVector<AlignPar>* alignPars=m_alignModuleTool->getAlignPars(alignModule);
317 const int nAlignPar = alignPars->size();
318 m_tmpChi2VAlignParam = new double*[nAlignPar];
319 m_tmpChi2VAlignParamX = new double*[nAlignPar];
322 for (int i=0;i<TrackState::NumberOfMeasurementTypes;i++)
323 m_tmpChi2VAlignParamMeasType[i] = new double*[nAlignPar];
324 }
325
326
327 // get derivatives and arrays of chi2 vs. align params
328 bool resetIPar=false;
329 std::vector<Amg::VectorX> tmpderiv_vec;
330 std::vector<Amg::VectorX> tmpderivErr_vec;
331 std::vector<double> tmpactualsecderiv_vec;
332 m_secPass=false;
333
334 // first attempt with normal number of fitter iterations
335 bool success=getAllDerivatives(
336 alignTrack, alignModule,
337 tmpderiv_vec,tmpderivErr_vec,tmpactualsecderiv_vec,
338 resetIPar);
339 if (!success){
340 delete derivatives;
341 delete derivativeErr;
342 delete actualSecondDerivatives;
343 return false;
344 }
345
347
348 if (resetIPar) {
349 // second attempt with increased number of fitter iterations
350 m_secPass=true;
351 success=getAllDerivatives(alignTrack,alignModule,
352 tmpderiv_vec,tmpderivErr_vec,tmpactualsecderiv_vec,
353 resetIPar);
354 }
355
356 if (!success){
357 delete derivatives;
358 delete derivativeErr;
359 delete actualSecondDerivatives;
360 return false;
361 }
362
364
365 if (resetIPar) {
366 // third and last attempt with number of fitter iterations set to maximum
368 success=getAllDerivatives(alignTrack,alignModule,
369 tmpderiv_vec,tmpderivErr_vec,tmpactualsecderiv_vec,
370 resetIPar);
371 }
372
373 if (!success){
374 delete derivatives;
375 delete derivativeErr;
376 delete actualSecondDerivatives;
377 return false;
378 }
379
380
382
383 if (success && !resetIPar) {
384 for (int i=0;i<(int)tmpderiv_vec.size();i++) {
385 deriv_vec.push_back(tmpderiv_vec[i]);
386 derivErr_vec.push_back(tmpderivErr_vec[i]);
387 actualsecderiv_vec.push_back(tmpactualsecderiv_vec[i]);
388 }
389 }
390 else{
391 delete derivatives;
392 delete derivativeErr;
393 delete actualSecondDerivatives;
394 return false;
395 }
396 // set the chi2 vs. align param arrays
397 ATH_MSG_DEBUG("setting chi2 vs. align param arrays");
400 (*alignModule).setChi2VAlignParamArray (m_tmpChi2VAlignParam);
401 (*alignModule).setChi2VAlignParamXArray(m_tmpChi2VAlignParamX);
402
403 // arrays for measurement types
405 ATH_MSG_DEBUG("pushing back for measType");
407 for (int i=0;i<TrackState::NumberOfMeasurementTypes;i++)
408 (*alignModule).setChi2VAlignParamArrayMeasType(i,m_tmpChi2VAlignParamMeasType[i]);
409 }
410 ATH_MSG_DEBUG("done setting arrays");
411
412 derivatives->push_back(make_pair(alignModule,deriv_vec));
413 derivativeErr->push_back(make_pair(alignModule,derivErr_vec));
414 actualSecondDerivatives->push_back(make_pair(alignModule,actualsecderiv_vec));
415 }
416
418
419 alignTrack->setDerivatives(derivatives);
420 alignTrack->setDerivativeErr(derivativeErr);
421 alignTrack->setActualSecondDerivatives(actualSecondDerivatives);
422
423 // restore unshifted residuals in AlignTSOS
424 setUnshiftedResiduals(alignTrack);
425
426 return true;
427}
428
429//________________________________________________________________________
431 AlignTrack* alignTrack,
432 int ipar, AlignPar* alignPar,
433 Amg::VectorX& derivativeErr,
434 bool& resetIPar,
435 double& actualSecondDeriv)
436{
437 const Trk::Track* trackForRefit =
439 dynamic_cast<const Trk::Track*>(alignTrack);
440
441 ATH_MSG_DEBUG("m_nIterations: "<<m_nIterations);
442
443 // gets derivatives of residuals w.r.t. a specific alignment parameter given by alignPar
444 if (!m_fitter)
445 ATH_MSG_ERROR("set m_fitter before calling getDerivatives (by calling setUnshiftedResiduals)");
446
447 AlignModule* module=alignPar->alignModule();
448
449 // set derivatives for 2 shifts up and 2 shifts down
450 const int NFITS = m_nFits;
451 const int NMEAS = alignTrack->nAlignTSOSMeas();
452 module->setNChamberShifts(m_nFits);
453
454 ATH_MSG_DEBUG("NMEAS="<<NMEAS);
455 double** residuals=new double*[NFITS];
456 double** resErrors=new double*[NFITS];
457 double* chi2Array =new double[NFITS];
458 double* chi2ArrayX=new double[NFITS];
459
461 for (int i=0;i<TrackState::NumberOfMeasurementTypes;i++)
462 m_tmpChi2VAlignParamMeasType[i][ipar] =new double[NFITS];
463 }
464
465 for (int ifit=0;ifit<NFITS;ifit++) {
466 residuals[ifit]=new double[NMEAS];
467 resErrors[ifit]=new double[NMEAS];
468 }
469
470 // set the values for the unshifted track
471 const int unshiftedTrackIndex = m_doFits ? (m_nFits-1)/2 : 1;
472 chi2Array [unshiftedTrackIndex] = m_unshiftedTrackChi2;
473 ATH_MSG_DEBUG("chi2Array["<<unshiftedTrackIndex<<"]="<<chi2Array[unshiftedTrackIndex]);
474 chi2ArrayX[unshiftedTrackIndex] = 0.;
476 for (int i=0;i<TrackState::NumberOfMeasurementTypes;i++) {
477 m_tmpChi2VAlignParamMeasType[i][ipar][unshiftedTrackIndex]=
479 ATH_MSG_DEBUG("chi2ArrayMeasType["<<i<<"]["<<unshiftedTrackIndex<<"]="<<m_unshiftedTrackChi2MeasType[i]);
480 }
481 }
482
483
484 // get shift size
485 double shiftsize=shiftSize(alignPar);
486
488
489
490 ATH_MSG_VERBOSE("doing refits");
491 for (int ifit=0;ifit<NFITS;ifit++) {
492
493 ATH_MSG_VERBOSE("ifit="<<ifit);
494 int jfit=ifit;
495 if (ifit>unshiftedTrackIndex) {
496 jfit=NFITS-ifit+unshiftedTrackIndex;
497 }
498 if (m_doFits && ifit==unshiftedTrackIndex) {
499 for (int i=0;i<(int)m_unshiftedResiduals->rows();i++) {
500 residuals[ifit][i]=(*m_unshiftedResiduals)[i];
501 resErrors[ifit][i]=(*m_unshiftedResErrors)[i];
502 }
503 // change back in case it got changed on the other side of zero
504 shiftsize=shiftSize(alignPar);
505 continue;
506 }
507
508 // shift module and fit track
509 double currentshift = 0.;
510 if(m_doFits)
511 currentshift = shiftsize * (double)(jfit-unshiftedTrackIndex);
512 else
513 currentshift = (ifit==0) ? -1.*shiftsize : shiftsize;
514
515 ATH_MSG_DEBUG("current shift="<<currentshift<<" in getDerivatives");
516
517 m_alignModuleTool->shiftModule(module,alignTrack,
518 alignPar->paramType(),currentshift);
519
520
521 ATH_MSG_VERBOSE("fitting after shift");
522 const Track* refittedTrack=m_fitter->alignmentFit(alignCache,
523 *trackForRefit,
528 ATH_MSG_DEBUG("exceeded max number of iterations");
529 m_alignModuleTool->restoreModule(module);
530 resetIPar=false;
531 ATH_MSG_DEBUG("resetIPar set to false");
532 delete [] residuals; delete [] resErrors;
533 delete [] chi2Array; delete [] chi2ArrayX;
534 ATH_MSG_DEBUG("fail max iter");
536 Amg::VectorX derivatives(1);
537 return derivatives;
538 }
539 ATH_MSG_DEBUG("increasing m_nIterations to "<<m_nIterations<<" (not changing in fit yet)");
540 resetIPar=true;
541 ATH_MSG_DEBUG("resetIPar set to true");
542 }
543
544 // if resetIPar refit the rest of the tracks, but don't do anything with them until next pass
545 if (resetIPar) {
546 m_alignModuleTool->restoreModule(module);
547 continue;
548 }
549
550 if (!refittedTrack) {
551 msg(MSG::WARNING) << "track refit failed for jfit "<<jfit <<endmsg;
552 delete [] residuals; delete [] resErrors;
553 delete [] chi2Array; delete [] chi2ArrayX;
554 m_alignModuleTool->restoreModule(module);
555 if (!resetIPar || m_secPass) {
557 }
558 ATH_MSG_DEBUG("fail track refit, resetIPar "<<resetIPar<<", secPass "<<m_secPass);
559 Amg::VectorX derivatives(1);
560 return derivatives;
561 }
562 else
563 ATH_MSG_VERBOSE("track refit successful");
564
565 double chi2=refittedTrack->fitQuality()->chiSquared();
566
567 ATH_MSG_VERBOSE("jfit = "<<jfit);
568 double localChi2=m_residualCalculator->setResiduals(alignTrack,refittedTrack);
569 ATH_MSG_DEBUG("localChi2/fittedChi2="<<localChi2<<"/"<<chi2);
570
571 chi2ArrayX[jfit]= shiftsize * (double)(jfit-unshiftedTrackIndex);// / module->sigma(idof);
572 chi2Array[jfit]=localChi2;
573 ATH_MSG_DEBUG("chi2Array["<<jfit<<"]="<<chi2Array[jfit]);
575 for (int i=0;i<TrackState::NumberOfMeasurementTypes;i++) {
576 m_tmpChi2VAlignParamMeasType[i][ipar][jfit]= m_residualCalculator->chi2ForMeasType(i);
577 ATH_MSG_DEBUG("chi2ArrayMeasType["<<i<<"]["<<jfit<<"]="
578 <<m_tmpChi2VAlignParamMeasType[i][ipar][jfit]);
579 }
580 }
581
582 ATH_MSG_DEBUG("positions["<<jfit<<"]="<<chi2ArrayX[jfit]);
583
584 int imeas(0);
585 AlignTSOSCollection::const_iterator atsosItr=alignTrack->firstAtsos();
586 for (; atsosItr != alignTrack->lastAtsos(); ++atsosItr) {
587 if (!(*atsosItr)->isValid()) continue;
588 for (std::vector<Residual>::const_iterator itRes=(**atsosItr).firstResidual();
589 itRes!=(**atsosItr).lastResidual();++itRes,++imeas) {
590
591 if (refittedTrack) {
592 residuals[jfit][imeas]=itRes->residual();
593 resErrors[jfit][imeas]=std::sqrt(itRes->errSq());
594 }
595 else {
596 residuals[jfit][imeas]=resErrors[jfit][imeas]=0.;
597 }
598 ATH_MSG_DEBUG("residuals["<<jfit<<"]["<<imeas<<"]="<<residuals[jfit][imeas]);
599 ATH_MSG_DEBUG("resErrors["<<jfit<<"]["<<imeas<<"]="<<resErrors[jfit][imeas]);
600 }
601 }
602
603 delete refittedTrack; refittedTrack=nullptr;
604 ATH_MSG_VERBOSE("calling restoreModule");
605 m_alignModuleTool->restoreModule(module);
606 } // NFITS
607
608 int iimeas(0);
609 AlignTSOSCollection::const_iterator aatsosItr=alignTrack->firstAtsos();
610 for (; aatsosItr != alignTrack->lastAtsos(); ++aatsosItr) {
611 if (!(*aatsosItr)->isValid()) continue;
612 for (std::vector<Residual>::const_iterator itRes=(**aatsosItr).firstResidual();
613 itRes!=(**aatsosItr).lastResidual();++itRes,++iimeas) {
614 for (int ifit=0;ifit<NFITS;ifit++) {
615 ATH_MSG_DEBUG("["<<ifit<<"]["<<iimeas<<"] res="<<residuals[ifit][iimeas]<<
616 ", resErr="<<resErrors[ifit][iimeas]);
617 }
618 }
619 }
620
621 if (resetIPar) {
622 //resetIPar=false;
623 delete [] residuals; delete [] resErrors;
624 delete [] chi2Array; delete [] chi2ArrayX;
625 if (m_secPass) ATH_MSG_WARNING("failed second pass!");
626 ATH_MSG_DEBUG("returning to reset IPar");
627 Amg::VectorX derivatives;
628 return derivatives;
629 }
630
631 // check chi2 vs. chamber pos to see if discontinuous
632 TGraph* gr = new TGraph(m_nFits,chi2ArrayX,chi2Array);
633 gr->Fit("pol2","QF");
634 TF1* fit=gr->GetFunction("pol2");
635 double chi2 =fit->GetChisquare()/double(m_nFits-3);
636 double slope=fit->GetParameter(1);
637 actualSecondDeriv=fit->GetParameter(2);
638 delete gr;
639
640 ATH_MSG_DEBUG("discontinuity check: chi2="<<chi2);
641 alignTrack->setTrackAlignParamQuality(alignPar->paramType(),chi2);
642
643 // EventInfo
644 if (chi2>1.e-6 || std::fabs(slope)<1.e-10) {
645 const xAOD::EventInfo* eventInfo;
646 StatusCode sc=evtStore()->retrieve(eventInfo);
647 if (sc.isFailure())
648 ATH_MSG_ERROR("Couldn't retrieve event info");
649 int run=eventInfo->runNumber();
650 int evt=eventInfo->eventNumber();
651 ATH_MSG_DEBUG("discontinuity check: chi2="<<chi2<<", run/evt "<<run<<"/"<<evt);
652 }
653
654 //reset in case it got changed somewhere
655 shiftsize = shiftSize(alignPar);
656
657 //-----------------------------------------//
658 //-- get derivatives from residuals --//
659 //-----------------------------------------//
660 ATH_MSG_VERBOSE("calculating residuals");
661 Amg::VectorX derivatives(alignTrack->nAlignTSOSMeas(),0);
662 ATH_MSG_DEBUG("created derivatives with "<<derivatives.rows()<<" rows");
663
664 // if bad fit or first derivative close to zero, replace derivatives with zeros
665 if (chi2>m_trackAlignParamCut) {// || std::fabs(slope)<1.e-10 ) {
666 ATH_MSG_DEBUG("chi2/"<<m_nFits-3<<"="<<chi2);
667 delete [] residuals; delete [] resErrors;
668 delete [] chi2Array; delete [] chi2ArrayX;
669
672 ATH_MSG_DEBUG("exceeded max number of iterations");
673 resetIPar=false;
674 }
675 ATH_MSG_DEBUG("increasing m_nIterations to "<<m_nIterations<<" (not changing in fit yet)");
676 resetIPar=true;
677 ATH_MSG_INFO("fail align param cut, secPass "<<m_secPass);
678 if (m_secPass) {
680 }
681 Amg::VectorX emptyDerivatives;
682 return emptyDerivatives;
683 }
684
685 int imeas(0);
686 TCanvas* canv(nullptr);
687 std::vector<TGraph*> vecGraphs;
688 AlignTSOSCollection::const_iterator atsosItr=alignTrack->firstAtsos();
689 for (; atsosItr != alignTrack->lastAtsos(); ++atsosItr) {
690 if (!(*atsosItr)->isValid()) continue;
691 for (int idim=0;idim<(*atsosItr)->nResDim();idim++) {
692
693 double* gr_x = new double[NFITS];
694 double* gr_y = new double[NFITS]; // residuals only have float precision if determined from ESD
695 int ngoodfits=0;
696 for (int ifit=0;ifit<NFITS;ifit++) {
697 double residual=residuals[ifit][imeas];
698 double resError=resErrors[ifit][imeas];
699 if (residual>-999.) {
700 gr_x[ngoodfits] =chi2ArrayX[ifit];
701 gr_y[ngoodfits] =residual/resError;
702 ngoodfits++;
703 }
704 }
705
706 if (!m_doFits && ngoodfits==2) {
707 derivatives[imeas]=(residuals[1][imeas]-residuals[0][imeas])/(2.*shiftsize)*
708 resErrors[unshiftedTrackIndex][imeas];
709 }
710 else if (m_doFits && ngoodfits>3) {
711 TGraph* gr=new TGraph(ngoodfits,gr_x,gr_y);
712
714 gr->Fit("pol2","VF");
715 else
716 gr->Fit("pol2","QF");
717 TF1* fit=gr->GetFunction("pol2");
718
719 //double derivRatio=fit->GetParameter(2)/fit->GetParameter(1);
720 ATH_MSG_DEBUG("deriv["<<imeas<<"]="<<fit->GetParameter(1)<<" +/- "<<fit->GetParError(1)
721 <<", chi2="<<fit->GetChisquare());
722 derivatives[imeas]=fit->GetParameter(1)*resErrors[unshiftedTrackIndex][imeas]; // first derivative at x=0
723 derivativeErr[imeas]=fit->GetParError(1)*resErrors[unshiftedTrackIndex][imeas];
724
725
726 // plot residuals vs. chamber position
727 if (m_doResidualPlots) {
728 if (!canv) canv=new TCanvas("resPlots","resPlots");
729 canv->cd();
730 gr->SetMarkerStyle(20);
731 gr->Draw("AP");
732
733 gr->GetXaxis()->SetTitle("shift in chamber pos. from nominal (CLHEP::mm)");
734 gr->GetYaxis()->SetTitle("residual (CLHEP::mm)");
735
736 TPaveText* pave=new TPaveText(.4,.65,.97,.92,"NDC");
737 pave->SetFillColor(0);
738 pave->SetBorderSize(1);
739 std::stringstream measType; measType<<"meas type: ";
740 if ((*atsosItr)->measType()==TrackState::MDT) measType<<" MDT";
741 else if ((*atsosItr)->measType()==TrackState::TGC) measType<<" TGC";
742 else if ((*atsosItr)->measType()==TrackState::RPC) measType<<" RPC";
743 else measType<<" undefined";
744
745 pave->AddText(measType.str().c_str());
746
747 std::stringstream firstderivtxt,secndderivtxt,aptxt,chi2txt;
748 firstderivtxt<<fit->GetParameter(1)<<" +/- "<<fit->GetParError(1);
749 secndderivtxt<<fit->GetParameter(2)<<" +/- "<<fit->GetParError(2);
750 aptxt <<"alignPar "<<alignPar->paramType()<<", RIO in "<<(*atsosItr)->identify();
751 chi2txt<<"chi2="<<fit->GetChisquare();
752
753 pave->AddText(firstderivtxt.str().c_str());
754 pave->AddText(secndderivtxt.str().c_str());
755 pave->AddText(aptxt.str().c_str());
756 pave->AddText(chi2txt.str().c_str());
757 pave->Draw();
758
759 std::stringstream canvName;
760 canvName<<"resPlots_ap"<<alignPar->paramType()<<"_measType"
761 <<(*atsosItr)->measType()<<"_"<<imeas<<".eps";
762 canv->Print(canvName.str().c_str());
763 canv->Clear();
764
765 delete pave;
766 }
767 vecGraphs.push_back(gr);
768 }
769 else {
770 derivatives[imeas]=-999.;
771 derivativeErr[imeas]=-999.;
772 }
773
774 delete [] gr_y;
775 delete [] gr_x;
776
777 ++imeas;
778 }
779 }
780
781 // delete TGraphs and TCanvas
782 for (auto & vecGraph : vecGraphs)
783 delete vecGraph;
784 delete canv;
785
786 delete [] residuals;
787 delete [] resErrors;
788
789 // set chi2 v alignparam
790 for (int ifit=0;ifit<NFITS;ifit++) {
791 m_tmpChi2VAlignParamX[ipar]=chi2ArrayX;
792 m_tmpChi2VAlignParam [ipar]=chi2Array;
793 }
794
795 ATH_MSG_DEBUG("derivativeErr: "<<derivativeErr);
796 return derivatives;
797}
798
799//________________________________________________________________________
800double ShiftingDerivCalcTool::shiftSize(const AlignPar* alignPar) const {
801 bool rotation =
802 alignPar->paramType() == AlignModule::RotX ||
803 alignPar->paramType() == AlignModule::RotY ||
804 alignPar->paramType() == AlignModule::RotZ;
805
806 double shift = rotation ? m_rotSize : m_traSize;
807
808 //ok... this is kind of ugly.
809 double sigma=alignPar->sigma();
810 return shift * sigma;
811}
812
813//________________________________________________________________________
815{
816 Amg::MatrixX W(alignTrack->nAlignTSOSMeas(),alignTrack->nAlignTSOSMeas());
817
818 if (alignTrack->localErrorMatrixInv()) {
819 ATH_MSG_ERROR("Need to assign this matrix correctly: ShiftingDerivCalcTool.cxx:888");
820 W = *(alignTrack->localErrorMatrixInv());
821 //W.assign(*(alignTrack->localErrorMatrixInv()));
822 } else{
823 return false;
824 }
825 ATH_MSG_DEBUG("W: "<<W);
826
827 bool Wisvalid(true);
828 const double epsilon=1e-10;
829 for( int irow=0; irow<W.rows(); ++irow) {
830 Wisvalid = Wisvalid && W(irow,irow)>0;
831 if( !(W(irow,irow)>0) )
832 msg(MSG::WARNING) << "matrix invalid: " << W(irow,irow) << endmsg;
833
834 for(int icol=0; icol<=irow; ++icol) {
835
836 // this one must be true if everything else succeeded
837 double Wcorr = W(irow,icol)/sqrt(W(irow,irow)*W(icol,icol));
838 bool Wcorrisvalid = Wcorr+epsilon>=-1 && Wcorr-epsilon<=1;
839 Wisvalid = Wisvalid && Wcorrisvalid;
840 if( !Wcorrisvalid )
841 msg(MSG::WARNING) << "matrix corr invalid: " << Wcorr-1 << " " << Wcorr+1 << endmsg;
842 }
843 }
844
845 if (Wisvalid)
846 alignTrack->setWeightMatrix(new Amg::MatrixX(W));
847
848 alignTrack->setWeightMatrixFirstDeriv(new Amg::MatrixX(std::move(W)));
849
850 return true;
851}
852
853//________________________________________________________________________
855{
856 for (int i=0;i<(int)m_chi2VAlignParamVec.size();i++) {
857 delete [] m_chi2VAlignParamVec[i]; m_chi2VAlignParamVec[i]=nullptr;
858 delete [] m_chi2VAlignParamXVec[i]; m_chi2VAlignParamXVec[i]=nullptr;
859 }
860 m_chi2VAlignParamVec.clear();
861 m_chi2VAlignParamXVec.clear();
862
863 for (auto & i : m_chi2VAlignParamVecMeasType) {
864 delete [] i; i=nullptr;
865 }
867}
868
869//________________________________________________________________________
871 AlignTrack* alignTrack,
872 const AlignModule* alignModule,
873 std::vector<Amg::VectorX>& deriv_vec,
874 std::vector<Amg::VectorX>& derivErr_vec,
875 std::vector<double>& actualsecderiv_vec,
876 bool& resetIPar)
877{
878 resetIPar=false;
879
880 deriv_vec.clear();
881 derivErr_vec.clear();
882 actualsecderiv_vec.clear();
883
884 setUnshiftedResiduals(alignTrack); // this will set the min number of iterations to the new value
885
886 int ipar(0);
887 DataVector<AlignPar>* alignPars=m_alignModuleTool->getAlignPars(alignModule);
888 for (DataVector<AlignPar>::iterator it=alignPars->begin(); it!=alignPars->end(); ++it,ipar++) {
889 ATH_MSG_DEBUG("ipar: "<<ipar);
890 Amg::VectorX derivErr(alignTrack->nAlignTSOSMeas());
891 double actualSecondDeriv(0.);
892 const Amg::VectorX vec=getDerivatives(alignTrack,ipar,*it,derivErr,resetIPar,actualSecondDeriv);
893 ATH_MSG_DEBUG("vec size: "<<vec.rows());
894
895 ATH_MSG_DEBUG("resetIPar="<<resetIPar);
896 if (resetIPar) continue; // continue with derivatives to find max iteration
897
898 if (vec.rows()<1) return false; // derivatives won't be set for alignTrack because it's a bad track
899
900 deriv_vec.push_back(vec);
901 derivErr_vec.push_back(derivErr);
902 actualsecderiv_vec.push_back(actualSecondDeriv);
903
904 for (int i=0;i<m_nFits;i++) {
905 ATH_MSG_DEBUG("m_tmpChi2VAlignParam["<<ipar<<"]["
906 <<i<<"]="<<m_tmpChi2VAlignParam[ipar][i]);
907 }
908 }
909
910 return true;
911}
912
913} // end namespace
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
std::vector< size_t > vec
#define gr
static Double_t sc
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
ServiceHandle< StoreGateSvc > & evtStore()
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
Standard iterator.
Definition DataVector.h:842
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
double sigma() const
returns sigma
Definition AlignPar.h:65
AlignModule::TransformParameters paramType() const
returns the type of parameter (i.e.
Definition AlignPar.h:47
const AlignModule * alignModule() const
returns the AlignModule
Definition AlignPar.h:40
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
void setActualSecondDerivatives(std::vector< std::pair< AlignModule *, std::vector< double > > > *vec)
Definition AlignTrack.h:147
void setTrackAlignParamQuality(int i, double val)
Definition AlignTrack.h:182
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 Trk::Track * trackWithoutScattering() const
returns track with ScatteringAngle pointers all set to zero (used for refit by iPat)
void setDerivativeErr(std::vector< AlignModuleDerivatives > *vec)
Definition AlignTrack.h:141
void setResidualVector(Amg::VectorX *vec)
Definition AlignTrack.h:151
bool isSLTrack() const
method to determine whether a straight line track or not
double chiSquared() const
returns the of the overall track fit
Definition FitQuality.h:56
std::ostream * m_logStream
logfile output stream
Gaudi::Property< bool > m_removeScatteringBeforeRefit
Gaudi::Property< int > m_minIter
flag to remove scattering before refitting track
std::vector< double ** > m_chi2VAlignParamVec
track chi2[idof][ichambershift]
Gaudi::Property< int > m_particleNumber
ToolHandle< IGlobalTrackFitter > m_trackFitterTool
ToolHandle< IGlobalTrackFitter > m_fitter
int m_ntracksProcessed
number tracks processed
bool setUnshiftedResiduals(AlignTrack *alignTrack)
int m_ntracksPassSetUnshiftedRes
number tracks pass setting unshifted residuals
double shiftSize(const AlignPar *alignPar) const
Gaudi::Property< int > m_maxIter
set minimum number of iterations for first track fits
bool getAllDerivatives(AlignTrack *alignTrack, const AlignModule *alignModule, std::vector< Amg::VectorX > &deriv_vec, std::vector< Amg::VectorX > &derivErr_vec, std::vector< double > &actualsecderiv_vec, bool &resetIPar)
ShiftingDerivCalcTool(const std::string &type, const std::string &name, const IInterface *parent)
int m_ntracksPassInitScan
number tracks pass initial scan
std::unique_ptr< double[]> m_unshiftedTrackChi2MeasType
cut on value of track alignment parameter, determined from fit of chi2 vs.
Gaudi::Property< bool > m_doFits
int m_ntracksPassGetDeriv
number tracks pass getting derivatives
ToolHandle< IAlignResidualCalculator > m_residualCalculator
int m_ntracksPassGetDerivSecPass
number tracks pass 2nd pass of getting derivatives
Amg::VectorX getDerivatives(AlignTrack *alignTrack, int ipar, AlignPar *alignPar, Amg::VectorX &derivativeErr, bool &resetIPar, double &actualSecondDerivative)
int m_ntracksPassGetDerivLastPass
number tracks pass 2nd pass of getting derivatives
bool setDerivatives(AlignTrack *track)
sets derivatives of residuals w.r.t.
int m_ntracksPassDerivatives
number tracks pass setting derivatives
ToolHandle< IAlignModuleTool > m_alignModuleTool
ParticleHypothesis m_particleHypothesis
Gaudi::Property< double > m_rotSize
bool scanShifts(const AlignTrack *alignTrack, const std::vector< AlignModule * > &alignModules)
Gaudi::Property< bool > m_setMinIterations
reject track if exceed maximum number of iterations
Gaudi::Property< bool > m_doChi2VAlignParamMeasType
Gaudi::Property< bool > m_doResidualPlots
bool setResidualCovMatrix(AlignTrack *alignTrack) const
sets residual covariance matrix
Gaudi::Property< double > m_trackAlignParamCut
fit track with AlignModules shifted up and down in each extreme, find the number of iterations fitter...
ToolHandle< IGlobalTrackFitter > m_SLTrackFitterTool
std::vector< double *** > m_chi2VAlignParamVecMeasType
track chi2[idof][imeastype][ichambershift]
Gaudi::Property< double > m_traSize
Gaudi::Property< bool > m_runOutlierRemoval
std::vector< double ** > m_chi2VAlignParamXVec
chamber shift[idof][ichambershift]
const FitQuality * fitQuality() const
return a pointer to the fit quality const-overload
uint32_t runNumber() const
The current event's run number.
uint64_t eventNumber() const
The current event's event number.
double chi2(TH1 *h0, TH1 *h1)
std::string find(const std::string &s)
return a remapped string
Definition hcg.cxx:138
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
Ensure that the ATLAS eigen extensions are properly loaded.
STL namespace.
EventInfo_v1 EventInfo
Definition of the latest event info version.
int m_minIterations
sets the minimum number of iterations to be used in the track fit.
int m_iterationsOfLastFit
returns the number of iterations used by the last fit (count starts at 1 for a single-iteration fit)
int run(int argc, char *argv[])