]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliRelAlignerKalman.cxx
Updates from Mikolaj, requested by TPC
[u/mrichter/AliRoot.git] / STEER / AliRelAlignerKalman.cxx
CommitLineData
043badeb 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16///////////////////////////////////////////////////////////////////////////////
17//
18// Kalman filter based aligner:
044eb03e 19// Finds alignement constants for two tracking volumes (by default ITS
20// and TPC)
55d5da9e 21// Determines the inverse transformation of the second volume (TPC)
22// with respect to the first (ITS) (how to realign TPC to ITS)
23// by measuring the residual between the 2 tracks.
044eb03e 24// Additionally calculates some callibration parameters for TPC
043badeb 25// Fit parameters are:
26// - 3 shifts, x,y,z
27// - 3 Cardan angles, psi, theta, phi (see definition in alignment docs),
c3b5bfc1 28// - TPC drift velocity correction,
55d5da9e 29// - TPC time offset correction.
043badeb 30//
31// Basic usage:
3ebfbf52 32// When aligning two volumes, at any given time a single instance of
043badeb 33// the class should be active. The fit of the parameters is updated
55d5da9e 34// by adding new data using one of the Add.... methods:
043badeb 35//
3ebfbf52 36// In collision events add an ESD event to update the fit (adds all tracks):
55d5da9e 37//
3ebfbf52 38// Bool_t AddESDevent( AliESDevent* pTrack );
39//
40// or add each individual track
41//
42// AddESDtrack( AliESDtrack* pTrack );
55d5da9e 43//
044eb03e 44// For cosmic data, the assumption is that the tracking is done twice:
55d5da9e 45// once global and once only ITS and the tracklets are saved inside
044eb03e 46// one AliESDEvent. The method
55d5da9e 47//
48// Bool_t AddCosmicEvent( AliESDEvent* pEvent );
49//
044eb03e 50// then searches the event for matching tracklets and upon succes it updates.
51// One cosmic ideally triggers two updates: for the upper and lower half of
52// the cosmic (upper ITS tracklet+upper TPC tracklet, idem dito for lower)
55d5da9e 53//
3ebfbf52 54// by default give misalignment parameters for TPC as they appear to be.
55// TPC calibration parameters are always given as correction to values used in reco.
56//
044eb03e 57// _________________________________________________________________________
58// Expert options:
3ebfbf52 59// look at AddESDevent() and AddCosmicEvent() to get the idea of how the
60// aligner works, it's safe to repeat the needed steps outside of the class,
61// only public methods are used.
043badeb 62//
63// Origin: Mikolaj Krzewicki, Nikhef, Mikolaj.Krzewicki@cern.ch
64//
65//////////////////////////////////////////////////////////////////////////////
66
55d5da9e 67#include <iostream>
68#include <TObject.h>
69#include <TMath.h>
70#include <TMatrix.h>
71#include <TVector.h>
72#include <TVector3.h>
73#include <TDecompLU.h>
74#include <TArrayI.h>
75#include <TH1D.h>
76#include <TF1.h>
77
78#include "AliESDtrack.h"
55d5da9e 79#include "AliESDEvent.h"
55d5da9e 80#include "AliExternalTrackParam.h"
81
043badeb 82#include "AliRelAlignerKalman.h"
83
84ClassImp(AliRelAlignerKalman)
85
044eb03e 86//______________________________________________________________________________
043badeb 87AliRelAlignerKalman::AliRelAlignerKalman():
3ebfbf52 88 TObject(),
044eb03e 89 fAlpha(0.),
90 fLocalX(80.),
3ebfbf52 91 fPTrackParamArr1(new AliExternalTrackParam[fgkNTracksPerMeasurement]),
92 fPTrackParamArr2(new AliExternalTrackParam[fgkNTracksPerMeasurement]),
93 fMagField(0.),
16cbbc2c 94 fNMeasurementParams(4),
782e5230 95 fPX(new TVectorD( fgkNSystemParams )),
96 fPXcov(new TMatrixDSym( fgkNSystemParams )),
16cbbc2c 97 fPH(new TMatrixD( fNMeasurementParams, fgkNSystemParams )),
55d5da9e 98 fQ(1.e-15),
16cbbc2c 99 fPMeasurement(new TVectorD( fNMeasurementParams )),
100 fPMeasurementCov(new TMatrixDSym( fNMeasurementParams )),
101 fPMeasurementPrediction(new TVectorD( fNMeasurementParams )),
044eb03e 102 fOutRejSigmas(1.),
329b5744 103 //fDelta(new Double_t[fgkNSystemParams]),
16cbbc2c 104 fYZOnly(kFALSE),
3ebfbf52 105 fNumericalParanoia(kTRUE),
044eb03e 106 fRejectOutliers(kTRUE),
c3b5bfc1 107 fRequireMatchInTPC(kFALSE),
043badeb 108 fCuts(kFALSE),
c3b5bfc1 109 fMinPointsVol1(3),
044eb03e 110 fMinPointsVol2(50),
043badeb 111 fMinMom(0.),
55d5da9e 112 fMaxMom(1.e100),
c3b5bfc1 113 fMaxMatchingAngle(0.1),
55d5da9e 114 fMaxMatchingDistance(10.), //in cm
115 fCorrectionMode(kFALSE),
044eb03e 116 fNTracks(0),
117 fNUpdates(0),
118 fNOutliers(0),
119 fNMatchedCosmics(0),
c3b5bfc1 120 fNMatchedTPCtracklets(0),
044eb03e 121 fNProcessedEvents(0),
3ebfbf52 122 fTrackInBuffer(0),
a80268df 123 fTimeStamp(0),
329b5744 124 fRunNumber(0),
55d5da9e 125 fTPCvd(2.64),
3ebfbf52 126 fTPCZLengthA(2.4972500e02),
127 fTPCZLengthC(2.4969799e02)
043badeb 128{
55d5da9e 129 //Default constructor
130
131 //default seed: zero, reset errors to large default
132 Reset();
55d5da9e 133}
134
135//______________________________________________________________________________
136AliRelAlignerKalman::AliRelAlignerKalman(const AliRelAlignerKalman& a):
3ebfbf52 137 TObject(static_cast<TObject>(a)),
55d5da9e 138 fAlpha(a.fAlpha),
139 fLocalX(a.fLocalX),
3ebfbf52 140 fPTrackParamArr1(new AliExternalTrackParam[fgkNTracksPerMeasurement]),
141 fPTrackParamArr2(new AliExternalTrackParam[fgkNTracksPerMeasurement]),
142 fMagField(a.fMagField),
16cbbc2c 143 fNMeasurementParams(a.fNMeasurementParams),
55d5da9e 144 fPX(new TVectorD( *a.fPX )),
145 fPXcov(new TMatrixDSym( *a.fPXcov )),
146 fPH(new TMatrixD( *a.fPH )),
147 fQ(a.fQ),
148 fPMeasurement(new TVectorD( *a.fPMeasurement )),
149 fPMeasurementCov(new TMatrixDSym( *a.fPMeasurementCov )),
a80268df 150 fPMeasurementPrediction(new TVectorD( *a.fPMeasurement )),
55d5da9e 151 fOutRejSigmas(a.fOutRejSigmas),
329b5744 152 //fDelta(new Double_t[fgkNSystemParams]),
16cbbc2c 153 fYZOnly(a.fYZOnly),
3ebfbf52 154 fNumericalParanoia(a.fNumericalParanoia),
55d5da9e 155 fRejectOutliers(a.fRejectOutliers),
55d5da9e 156 fRequireMatchInTPC(a.fRequireMatchInTPC),
55d5da9e 157 fCuts(a.fCuts),
158 fMinPointsVol1(a.fMinPointsVol1),
159 fMinPointsVol2(a.fMinPointsVol2),
160 fMinMom(a.fMinMom),
161 fMaxMom(a.fMaxMom),
162 fMaxMatchingAngle(a.fMaxMatchingAngle),
163 fMaxMatchingDistance(a.fMaxMatchingDistance), //in cm
164 fCorrectionMode(a.fCorrectionMode),
165 fNTracks(a.fNTracks),
166 fNUpdates(a.fNUpdates),
167 fNOutliers(a.fNOutliers),
168 fNMatchedCosmics(a.fNMatchedCosmics),
169 fNMatchedTPCtracklets(a.fNMatchedTPCtracklets),
170 fNProcessedEvents(a.fNProcessedEvents),
3ebfbf52 171 fTrackInBuffer(a.fTrackInBuffer),
a80268df 172 fTimeStamp(a.fTimeStamp),
329b5744 173 fRunNumber(a.fRunNumber),
55d5da9e 174 fTPCvd(a.fTPCvd),
175 fTPCZLengthA(a.fTPCZLengthA),
176 fTPCZLengthC(a.fTPCZLengthC)
3ebfbf52 177 //fApplyCovarianceCorrection(a.fApplyCovarianceCorrection),
178 //fCalibrationMode(a.fCalibrationMode),
179 //fFillHistograms(a.fFillHistograms),
180 //fNHistogramBins(a.fNHistogramBins),
181 //fPMes0Hist(new TH1D(*a.fPMes0Hist)),
182 //fPMes1Hist(new TH1D(*a.fPMes1Hist)),
183 //fPMes2Hist(new TH1D(*a.fPMes2Hist)),
184 //fPMes3Hist(new TH1D(*a.fPMes3Hist)),
185 //fPMesErr0Hist(new TH1D(*a.fPMesErr0Hist)),
186 //fPMesErr1Hist(new TH1D(*a.fPMesErr1Hist)),
187 //fPMesErr2Hist(new TH1D(*a.fPMesErr2Hist)),
188 //fPMesErr3Hist(new TH1D(*a.fPMesErr3Hist)),
189 //fPMeasurementCovCorr(new TMatrixDSym(*a.fPMeasurementCovCorr)),
55d5da9e 190{
191 //copy constructor
192 memcpy(fDelta,a.fDelta,fgkNSystemParams*sizeof(Double_t));
193}
194
195//______________________________________________________________________________
196AliRelAlignerKalman& AliRelAlignerKalman::operator=(const AliRelAlignerKalman& a)
197{
198 //assignment operator
3ebfbf52 199 //fAlpha=a.fAlpha;
200 //fLocalX=a.fLocalX;
201 //memcpy(fPTrackParamArr1,a.fPTrackParamArr1,fgkNTracksPerMeasurement*sizeof(AliExternalTrackParam));
202 //memcpy(fPTrackParamArr2,a.fPTrackParamArr2,fgkNTracksPerMeasurement*sizeof(AliExternalTrackParam));
203 fMagField=a.fMagField,
55d5da9e 204 *fPX = *a.fPX;
205 *fPXcov = *a.fPXcov;
3ebfbf52 206 //*fPH = *a.fPH;
55d5da9e 207 fQ=a.fQ;
3ebfbf52 208 //*fPMeasurement=*a.fPMeasurement;
209 //*fPMeasurementCov=*a.fPMeasurementCov;
55d5da9e 210 fOutRejSigmas=a.fOutRejSigmas;
211 memcpy(fDelta,a.fDelta,fgkNSystemParams*sizeof(Double_t));
16cbbc2c 212 fYZOnly=a.fYZOnly;
3ebfbf52 213 fNumericalParanoia=a.fNumericalParanoia;
55d5da9e 214 fRejectOutliers=a.fRejectOutliers;
55d5da9e 215 fRequireMatchInTPC=a.fRequireMatchInTPC;
55d5da9e 216 fCuts=a.fCuts;
217 fMinPointsVol1=a.fMinPointsVol1;
218 fMinPointsVol2=a.fMinPointsVol2;
219 fMinMom=a.fMinMom;
220 fMaxMom=a.fMaxMom;
221 fMaxMatchingAngle=a.fMaxMatchingAngle;
222 fMaxMatchingDistance=a.fMaxMatchingDistance; //in c;
223 fCorrectionMode=a.fCorrectionMode;
224 fNTracks=a.fNTracks;
225 fNUpdates=a.fNUpdates;
226 fNOutliers=a.fNOutliers;
227 fNMatchedCosmics=a.fNMatchedCosmics;
228 fNMatchedTPCtracklets=a.fNMatchedTPCtracklets;
229 fNProcessedEvents=a.fNProcessedEvents;
3ebfbf52 230 fTrackInBuffer=a.fTrackInBuffer;
a80268df 231 fTimeStamp=a.fTimeStamp;
329b5744 232 fRunNumber=a.fRunNumber;
3ebfbf52 233 //fApplyCovarianceCorrection=a.fApplyCovarianceCorrection;
234 //fCalibrationMode=a.fCalibrationMode;
235 //fFillHistograms=a.fFillHistograms;
236 //fNHistogramBins=a.fNHistogramBins;
237 //*fPMes0Hist=*(a.fPMes0Hist);
238 //*fPMes1Hist=*(a.fPMes1Hist);
239 //*fPMes2Hist=*(a.fPMes2Hist);
240 //*fPMes3Hist=*(a.fPMes3Hist);
241 //*fPMesErr0Hist=*(a.fPMesErr0Hist);
242 //*fPMesErr1Hist=*(a.fPMesErr1Hist);
243 //*fPMesErr2Hist=*(a.fPMesErr2Hist);
244 //*fPMesErr3Hist=*(a.fPMesErr3Hist);
245 //*fPMeasurementCovCorr=*(a.fPMeasurementCovCorr);
55d5da9e 246 fTPCvd=a.fTPCvd;
247 fTPCZLengthA=a.fTPCZLengthA;
248 fTPCZLengthC=a.fTPCZLengthC;
249 return *this;
043badeb 250}
251
044eb03e 252//______________________________________________________________________________
253AliRelAlignerKalman::~AliRelAlignerKalman()
043badeb 254{
55d5da9e 255 //destructor
3ebfbf52 256 if (fPTrackParamArr1) delete [] fPTrackParamArr1;
257 if (fPTrackParamArr2) delete [] fPTrackParamArr2;
258 if (fPX) delete fPX;
259 if (fPXcov) delete fPXcov;
260 if (fPH) delete fPH;
261 if (fPMeasurement) delete fPMeasurement;
262 if (fPMeasurementCov) delete fPMeasurementCov;
329b5744 263 //if (fDelta) delete [] fDelta;
3ebfbf52 264 //if (fPMes0Hist) delete fPMes0Hist;
265 //if (fPMes1Hist) delete fPMes1Hist;
266 //if (fPMes2Hist) delete fPMes2Hist;
267 //if (fPMes3Hist) delete fPMes3Hist;
268 //if (fPMesErr0Hist) delete fPMesErr0Hist;
269 //if (fPMesErr1Hist) delete fPMesErr1Hist;
270 //if (fPMesErr2Hist) delete fPMesErr2Hist;
271 //if (fPMesErr3Hist) delete fPMesErr3Hist;
272 //if (fPMeasurementCovCorr) delete fPMeasurementCovCorr;
273}
274
275//______________________________________________________________________________
276Bool_t AliRelAlignerKalman::AddESDevent( const AliESDEvent* pEvent )
277{
278 //Add all tracks in an ESD event
279
280 fNProcessedEvents++; //update the counter
281
282 Bool_t success=kFALSE;
283 SetMagField( pEvent->GetMagneticField() );
284 AliESDtrack* track;
285
286 for (Int_t i=0; i<pEvent->GetNumberOfTracks(); i++)
287 {
288 track = pEvent->GetTrack(i);
289 if (!track) continue;
290 if ( ((track->GetStatus()&AliESDtrack::kTPCin)>0)&&
291 ((track->GetStatus()&AliESDtrack::kITSout)>0)&&
292 (track->GetNcls(0)>=fMinPointsVol1)&&
293 (track->GetNcls(1)>=fMinPointsVol2) )
294 {
295 success = ( AddESDtrack( track ) || success );
296 }
297 }
329b5744 298 if (success)
299 {
300 fTimeStamp = pEvent->GetTimeStamp();
301 fRunNumber = pEvent->GetRunNumber();
302 }
3ebfbf52 303 return success;
043badeb 304}
305
044eb03e 306//______________________________________________________________________________
3ebfbf52 307Bool_t AliRelAlignerKalman::AddESDtrack( const AliESDtrack* pTrack )
043badeb 308{
a80268df 309 //Adds a full track, returns true if results in a new estimate
310 // gets the inner TPC parameters from AliESDTrack::GetInnerParam()
311 // gets the outer ITS parameters from AliESDTrack::GetOuterParam()
312
313 const AliExternalTrackParam* pconstparamsITS = pTrack->GetOuterParam();
314 if (!pconstparamsITS) return kFALSE;
315 const AliExternalTrackParam* pconstparamsTPC = pTrack->GetInnerParam();
316 if (!pconstparamsTPC) return kFALSE;
3ebfbf52 317
a80268df 318 //TPC part
319 AliExternalTrackParam paramsTPC = (*pconstparamsTPC);
320 paramsTPC.Rotate(pconstparamsITS->GetAlpha());
321 paramsTPC.PropagateTo(pconstparamsITS->GetX(), fMagField);
322
323 return (AddTrackParams(pconstparamsITS, &paramsTPC));
043badeb 324}
325
18f65700 326//______________________________________________________________________________
327Bool_t AliRelAlignerKalman::AddTrackParams( const AliExternalTrackParam* p1, const AliExternalTrackParam* p2 )
328{
329 //Update the estimate using new matching tracklets
330
331 if (!SetTrackParams(p1, p2)) return kFALSE;
332 return Update();
333}
334
044eb03e 335//______________________________________________________________________________
55d5da9e 336Bool_t AliRelAlignerKalman::AddCosmicEvent( const AliESDEvent* pEvent )
043badeb 337{
55d5da9e 338 //Add an cosmic with separately tracked ITS and TPC parts, do trackmatching
339
3ebfbf52 340 fNProcessedEvents++; //update the counter
341
55d5da9e 342 Bool_t success=kFALSE;
343 TArrayI trackTArrITS(1);
344 TArrayI trackTArrTPC(1);
345 if (!FindCosmicTrackletNumbersInEvent( trackTArrITS, trackTArrTPC, pEvent )) return kFALSE;
3ebfbf52 346 SetMagField( pEvent->GetMagneticField() );
55d5da9e 347 AliESDtrack* ptrack;
3ebfbf52 348 const AliExternalTrackParam* pconstparams1;
349 const AliExternalTrackParam* pconstparams2;
350 AliExternalTrackParam params1;
351 AliExternalTrackParam params2;
352
55d5da9e 353 ////////////////////////////////
354 for (Int_t i=0;i<trackTArrITS.GetSize();i++)
355 {
356 //ITS track
357 ptrack = pEvent->GetTrack(trackTArrITS[i]);
3ebfbf52 358 pconstparams1 = ptrack->GetOuterParam();
359 if (!pconstparams1) continue;
360 params1 = *pconstparams1; //make copy to be safe
044eb03e 361
55d5da9e 362 //TPC track
363 ptrack = pEvent->GetTrack(trackTArrTPC[i]);
3ebfbf52 364 pconstparams2 = ptrack->GetInnerParam();
365 if (!pconstparams2) continue;
366 params2 = *pconstparams2; //make copy
367 params2.Rotate(params1.GetAlpha());
a80268df 368 params2.PropagateTo( params1.GetX(), fMagField );
3ebfbf52 369
370 if (!SetTrackParams( &params1, &params2 )) continue;
55d5da9e 371
372 //do some accounting and update
55d5da9e 373 if (Update())
374 success = kTRUE;
375 else
376 continue;
377 }
329b5744 378 fTimeStamp=pEvent->GetTimeStamp(); //always update timestamp even when no update performed
379 fRunNumber=pEvent->GetRunNumber();
55d5da9e 380 return success;
043badeb 381}
382
16cbbc2c 383//______________________________________________________________________________
384void AliRelAlignerKalman::SetPoint2Track( Bool_t set )
385{
386 fNMeasurementParams = (set)?2:4;
387 if (fPH) delete fPH;
388 fPH = new TMatrixD( fNMeasurementParams, fgkNSystemParams );
389 if (fPMeasurement) delete fPMeasurement;
390 fPMeasurement = new TVectorD( fNMeasurementParams );
391 if (fPMeasurementCov) delete fPMeasurementCov;
392 fPMeasurementCov = new TMatrixDSym( fNMeasurementParams );
393 if (fPMeasurementPrediction) delete fPMeasurementPrediction;
394 fPMeasurementPrediction = new TVectorD( fNMeasurementParams );
395 fYZOnly = set;
396}
397
044eb03e 398//______________________________________________________________________________
399void AliRelAlignerKalman::Print(Option_t*) const
043badeb 400{
55d5da9e 401 //Print some useful info
402 Double_t rad2deg = 180./TMath::Pi();
3ebfbf52 403 printf("\nAliRelAlignerKalman\n");
404 if (fCorrectionMode) printf("(Correction mode)\n");
55d5da9e 405 printf(" %i pairs, %i updates, %i outliers,\n", fNTracks, fNUpdates, fNOutliers );
3ebfbf52 406 printf(" %i TPC matches, %i ITS-TPC matches in %i events\n", fNMatchedTPCtracklets, fNMatchedCosmics, fNProcessedEvents );
55d5da9e 407 printf(" psi(x): % .3f ± (%.2f) mrad | % .3f ± (%.2f) deg\n",1e3*(*fPX)(0), 1e3*TMath::Sqrt((*fPXcov)(0,0)),(*fPX)(0)*rad2deg,TMath::Sqrt((*fPXcov)(0,0))*rad2deg);
408 printf(" theta(y): % .3f ± (%.2f) mrad | % .3f ± (%.2f) deg\n",1e3*(*fPX)(1), 1e3*TMath::Sqrt((*fPXcov)(1,1)),(*fPX)(1)*rad2deg,TMath::Sqrt((*fPXcov)(1,1))*rad2deg);
409 printf(" phi(z): % .3f ± (%.2f) mrad | % .3f ± (%.2f) deg\n",1e3*(*fPX)(2), 1e3*TMath::Sqrt((*fPXcov)(2,2)),(*fPX)(2)*rad2deg,TMath::Sqrt((*fPXcov)(2,2))*rad2deg);
410 printf(" x: % .3f ± (%.2f) micron\n", 1e4*(*fPX)(3), 1e4*TMath::Sqrt((*fPXcov)(3,3)));
411 printf(" y: % .3f ± (%.2f) micron\n", 1e4*(*fPX)(4), 1e4*TMath::Sqrt((*fPXcov)(4,4)));
412 printf(" z: % .3f ± (%.2f) micron\n", 1e4*(*fPX)(5), 1e4*TMath::Sqrt((*fPXcov)(5,5)));
3ebfbf52 413 if (fgkNSystemParams>6) printf(" vd corr % .5g ± (%.2g) [ vd should be %.4g (was %.4g in reco) ]\n", (*fPX)(6), TMath::Sqrt((*fPXcov)(6,6)), (*fPX)(6)*fTPCvd, fTPCvd);
414 if (fgkNSystemParams>7) printf(" t0 % .5g ± (%.2g) us | %.4g ± (%.2g) cm [ t0_real = t0_rec+t0 ]\n",(*fPX)(7), TMath::Sqrt((*fPXcov)(7,7)), fTPCvd*(*fPX)(7), fTPCvd*TMath::Sqrt((*fPXcov)(7,7)));
415 if (fgkNSystemParams>8) printf(" vd/dy % .5f ± (%.2f) (cm/us)/m\n", (*fPX)(8), TMath::Sqrt((*fPXcov)(8,8)));
329b5744 416 printf(" run: %i, timestamp: %i\n", fRunNumber, fTimeStamp);
55d5da9e 417 printf("\n");
418 return;
043badeb 419}
420
044eb03e 421//______________________________________________________________________________
422void AliRelAlignerKalman::PrintSystemMatrix()
043badeb 423{
55d5da9e 424 //Print the system matrix for this measurement
425 printf("Kalman system matrix:\n");
16cbbc2c 426 for ( Int_t i=0; i<fNMeasurementParams; i++ )
55d5da9e 427 {
428 for ( Int_t j=0; j<fgkNSystemParams; j++ )
044eb03e 429 {
55d5da9e 430 printf("% -2.2f ", (*fPH)(i,j) );
431 }//for i
044eb03e 432 printf("\n");
55d5da9e 433 }//for j
434 printf("\n");
435 return;
043badeb 436}
437
044eb03e 438//______________________________________________________________________________
3ebfbf52 439Bool_t AliRelAlignerKalman::SetTrackParams( const AliExternalTrackParam* exparam1, const AliExternalTrackParam* exparam2 )
043badeb 440{
3ebfbf52 441 //Set the parameters, exparam1 will normally be ITS and exparam 2 tht TPC
442
443 fPTrackParamArr1[fTrackInBuffer] = *exparam1;
444 fPTrackParamArr2[fTrackInBuffer] = *exparam2;
445
446 fTrackInBuffer++;
043badeb 447
3ebfbf52 448 if (fTrackInBuffer == fgkNTracksPerMeasurement)
449 {
450 fTrackInBuffer = 0;
451 return kTRUE;
452 }
453 return kFALSE;
043badeb 454}
455
044eb03e 456//______________________________________________________________________________
3ebfbf52 457Bool_t AliRelAlignerKalman::Update()
043badeb 458{
3ebfbf52 459 //perform the update
460
461 //if (fCalibrationMode) return UpdateCalibration();
462 //if (fFillHistograms)
463 //{
464 // if (!UpdateEstimateKalman()) return kFALSE;
465 // return UpdateCalibration(); //Update histograms only when update ok.
466 //}
467 //else return UpdateEstimateKalman();
55d5da9e 468 if (!PrepareMeasurement()) return kFALSE;
469 if (!PrepareSystemMatrix()) return kFALSE;
a80268df 470 if (!PreparePrediction()) return kFALSE;
3ebfbf52 471 return UpdateEstimateKalman();
043badeb 472}
473
044eb03e 474//______________________________________________________________________________
043badeb 475void AliRelAlignerKalman::RotMat( TMatrixD &R, const TVectorD& angles )
476{
55d5da9e 477 //Get Rotation matrix R given the Cardan angles psi, theta, phi (around x, y, z).
478 Double_t sinpsi = TMath::Sin(angles(0));
479 Double_t sintheta = TMath::Sin(angles(1));
480 Double_t sinphi = TMath::Sin(angles(2));
481 Double_t cospsi = TMath::Cos(angles(0));
482 Double_t costheta = TMath::Cos(angles(1));
483 Double_t cosphi = TMath::Cos(angles(2));
484
485 R(0,0) = costheta*cosphi;
486 R(0,1) = -costheta*sinphi;
487 R(0,2) = sintheta;
488 R(1,0) = sinpsi*sintheta*cosphi + cospsi*sinphi;
489 R(1,1) = -sinpsi*sintheta*sinphi + cospsi*cosphi;
490 R(1,2) = -costheta*sinpsi;
491 R(2,0) = -cospsi*sintheta*cosphi + sinpsi*sinphi;
492 R(2,1) = cospsi*sintheta*sinphi + sinpsi*cosphi;
493 R(2,2) = costheta*cospsi;
043badeb 494}
495
044eb03e 496//______________________________________________________________________________
497Bool_t AliRelAlignerKalman::PrepareMeasurement()
043badeb 498{
55d5da9e 499 //Calculate the residuals and their covariance matrix
3ebfbf52 500
501 for (Int_t i=0;i<fgkNTracksPerMeasurement;i++)
502 {
503 const Double_t* pararr1 = fPTrackParamArr1[i].GetParameter();
504 const Double_t* pararr2 = fPTrackParamArr2[i].GetParameter();
505
506 //Take the track parameters and calculate the input to the Kalman filter
507 Int_t x = i*4;
508 (*fPMeasurement)(x+0) = pararr2[0]-pararr1[0];
509 (*fPMeasurement)(x+1) = pararr2[1]-pararr1[1];
16cbbc2c 510 if (!fYZOnly)
511 {
512 (*fPMeasurement)(x+2) = pararr2[2]-pararr1[2];
513 (*fPMeasurement)(x+3) = pararr2[3]-pararr1[3];
514 }
3ebfbf52 515
516 //the covariance
517 const Double_t* parcovarr1 = fPTrackParamArr1[i].GetCovariance();
518 const Double_t* parcovarr2 = fPTrackParamArr2[i].GetCovariance();
519 (*fPMeasurementCov)(x+0,x+0)=parcovarr1[0];
520 (*fPMeasurementCov)(x+0,x+1)=parcovarr1[1];
3ebfbf52 521 (*fPMeasurementCov)(x+1,x+0)=parcovarr1[1];
522 (*fPMeasurementCov)(x+1,x+1)=parcovarr1[2];
3ebfbf52 523 (*fPMeasurementCov)(x+0,x+0)+=parcovarr2[0];
524 (*fPMeasurementCov)(x+0,x+1)+=parcovarr2[1];
3ebfbf52 525 (*fPMeasurementCov)(x+1,x+0)+=parcovarr2[1];
526 (*fPMeasurementCov)(x+1,x+1)+=parcovarr2[2];
16cbbc2c 527 if (!fYZOnly)
528 {
529 (*fPMeasurementCov)(x+0,x+2)=parcovarr1[3];
530 (*fPMeasurementCov)(x+0,x+3)=parcovarr1[6];
531 (*fPMeasurementCov)(x+1,x+2)=parcovarr1[4];
532 (*fPMeasurementCov)(x+1,x+3)=parcovarr1[7];
533 (*fPMeasurementCov)(x+2,x+0)=parcovarr1[3];
534 (*fPMeasurementCov)(x+2,x+1)=parcovarr1[4];
535 (*fPMeasurementCov)(x+2,x+2)=parcovarr1[5];
536 (*fPMeasurementCov)(x+2,x+3)=parcovarr1[8];
537 (*fPMeasurementCov)(x+3,x+0)=parcovarr1[6];
538 (*fPMeasurementCov)(x+3,x+1)=parcovarr1[7];
539 (*fPMeasurementCov)(x+3,x+2)=parcovarr1[8];
540 (*fPMeasurementCov)(x+3,x+3)=parcovarr1[9];
541 (*fPMeasurementCov)(x+0,x+2)+=parcovarr2[3];
542 (*fPMeasurementCov)(x+0,x+3)+=parcovarr2[6];
543 (*fPMeasurementCov)(x+1,x+2)+=parcovarr2[4];
544 (*fPMeasurementCov)(x+1,x+3)+=parcovarr2[7];
545 (*fPMeasurementCov)(x+2,x+0)+=parcovarr2[3];
546 (*fPMeasurementCov)(x+2,x+1)+=parcovarr2[4];
547 (*fPMeasurementCov)(x+2,x+2)+=parcovarr2[5];
548 (*fPMeasurementCov)(x+2,x+3)+=parcovarr2[8];
549 (*fPMeasurementCov)(x+3,x+0)+=parcovarr2[6];
550 (*fPMeasurementCov)(x+3,x+1)+=parcovarr2[7];
551 (*fPMeasurementCov)(x+3,x+2)+=parcovarr2[8];
552 (*fPMeasurementCov)(x+3,x+3)+=parcovarr2[9];
553 }
3ebfbf52 554
555 fNTracks++; //count added track sets
556 }
557 //if (fApplyCovarianceCorrection)
558 // *fPMeasurementCov += *fPMeasurementCovCorr;
55d5da9e 559 return kTRUE;
043badeb 560}
561
044eb03e 562//______________________________________________________________________________
563Bool_t AliRelAlignerKalman::PrepareSystemMatrix()
043badeb 564{
55d5da9e 565 //Calculate the system matrix for the Kalman filter
566 //approximate the system using as reference the track in the first volume
567
16cbbc2c 568 TVectorD z1( fNMeasurementParams );
569 TVectorD z2( fNMeasurementParams );
55d5da9e 570 TVectorD x1( fgkNSystemParams );
571 TVectorD x2( fgkNSystemParams );
55d5da9e 572 //get the derivatives
573 for ( Int_t i=0; i<fgkNSystemParams; i++ )
574 {
575 x1 = *fPX;
576 x2 = *fPX;
3ebfbf52 577 x1(i) = x1(i) - fDelta[i]/(2.0);
578 x2(i) = x2(i) + fDelta[i]/(2.0);
55d5da9e 579 if (!PredictMeasurement( z1, x1 )) return kFALSE;
580 if (!PredictMeasurement( z2, x2 )) return kFALSE;
16cbbc2c 581 for (Int_t j=0; j<fNMeasurementParams; j++ )
3ebfbf52 582 {
583 (*fPH)(j,i) = ( z2(j)-z1(j) ) / fDelta[i];
584 }
55d5da9e 585 }
586 return kTRUE;
043badeb 587}
588
a80268df 589//______________________________________________________________________________
590Bool_t AliRelAlignerKalman::PreparePrediction()
591{
592 //Prepare the prediction of the measurement using state vector
593 return PredictMeasurement( (*fPMeasurementPrediction), (*fPX) );
594}
595
044eb03e 596//______________________________________________________________________________
043badeb 597Bool_t AliRelAlignerKalman::PredictMeasurement( TVectorD& pred, const TVectorD& state )
598{
55d5da9e 599 // Implements a system model for the Kalman fit
600 // pred is [dy,dz,dsinphi,dtgl]
601 // state is [psi,theta,phi,x,y,z,driftTPC,offsetTPC]
602 // note: the measurement is in a local frame, so the prediction also has to be
603 // note: state is the misalignment in global reference system
604
55d5da9e 605 if (fCorrectionMode)
606 {
3ebfbf52 607 for (Int_t i=0;i<fgkNTracksPerMeasurement;i++)
608 {
609 AliExternalTrackParam track(fPTrackParamArr2[i]); //make a copy track
610 if (!CorrectTrack( &track, state )) return kFALSE; //predict what the ideal track would be by applying correction
611
612 const Double_t* oldparam = fPTrackParamArr2[i].GetParameter();
613 const Double_t* newparam = track.GetParameter();
614
615 Int_t x = 4*i;
616 //calculate the predicted residual
617 pred(x+0) = oldparam[0] - newparam[0];
618 pred(x+1) = oldparam[1] - newparam[1];
16cbbc2c 619 if (!fYZOnly)
620 {
621 pred(x+2) = oldparam[2] - newparam[2];
622 pred(x+3) = oldparam[3] - newparam[3];
623 }
3ebfbf52 624 return kTRUE;
625 }
55d5da9e 626 }
627 else
628 {
3ebfbf52 629 for (Int_t i=0;i<fgkNTracksPerMeasurement;i++)
630 {
631 AliExternalTrackParam track(fPTrackParamArr1[i]); //make a copy track
632 if (!MisalignTrack( &track, state )) return kFALSE; //predict what the measured track would be by applying misalignment
633
634 const Double_t* oldparam = fPTrackParamArr1[i].GetParameter();
635 const Double_t* newparam = track.GetParameter();
636
637 Int_t x = 4*i;
638 //calculate the predicted residual
639 pred(x+0) = newparam[0] - oldparam[0];
640 pred(x+1) = newparam[1] - oldparam[1];
16cbbc2c 641 if (!fYZOnly)
642 {
643 pred(x+2) = newparam[2] - oldparam[2];
644 pred(x+3) = newparam[3] - oldparam[3];
645 }
3ebfbf52 646 return kTRUE;
647 }
55d5da9e 648 }
649 return kFALSE;
043badeb 650}
651
044eb03e 652//______________________________________________________________________________
043badeb 653Bool_t AliRelAlignerKalman::UpdateEstimateKalman()
654{
55d5da9e 655 //Kalman estimation of noisy constants: in the model A=1
656 //The arguments are (following the usual convention):
657 // fPX - the state vector (parameters)
658 // fPXcov - the state covariance matrix (parameter errors)
659 // fPMeasurement - measurement vector
660 // fPMeasurementCov - measurement covariance matrix
661 // fPH - measurement model matrix ( fPMeasurement = Hx + v ) v being measurement noise (error fR)
662
663 TMatrixDSym identity(TMatrixDSym::kUnit, (*fPXcov)); //unit matrix
664
665 //predict the state
666 //(*fPXcov) = (*fPXcov) + fQ*identity; //add some process noise (diagonal)
667
668 // update prediction with measurement
669 // calculate Kalman gain
670 // K = PH/(HPH+fPMeasurementCov)
671 TMatrixD pht( (*fPXcov), TMatrixD::kMultTranspose, (*fPH) ); //common factor (used twice)
672 TMatrixD hpht( (*fPH), TMatrixD::kMult, pht );
673 hpht += (*fPMeasurementCov);
3ebfbf52 674
55d5da9e 675 //shit happens so protect yourself!
3ebfbf52 676// if (fNumericalParanoia)
677// {
678// TDecompLU lu(hpht);
679// if (lu.Condition() > 1e12) return kFALSE;
680// lu.Invert(hpht);
681// }
682// else
683// {
684 Double_t det;
685 hpht.InvertFast(&det); //since the matrix is small...
686 if (det < 2e-55) return kFALSE; //we need some sort of protection even in this case....
687// }
688 //printf("KalmanUpdate: det(hpht): %.4g\n",det);
689
690 TMatrixD k(pht, TMatrixD::kMult, hpht ); //compute K (hpht is already inverted)
55d5da9e 691
692 // update the state and its covariance matrix
693 TVectorD xupdate(fgkNSystemParams);
a80268df 694 xupdate = k*((*fPMeasurement)-(*fPMeasurementPrediction));
55d5da9e 695
696 //SIMPLE OUTLIER REJECTION
697 if ( IsOutlier( xupdate, (*fPXcov) ) && fRejectOutliers )
698 {
699 fNOutliers++;
700 return kFALSE;
701 }
702
55d5da9e 703 TMatrixD kh( k, TMatrixD::kMult, (*fPH) );
704 TMatrixD ikh(fgkNSystemParams,fgkNSystemParams); //this is because for some reason TMatrixD::kAdd didn't work
705 ikh = identity - kh;
706 TMatrixD ikhp( ikh, TMatrixD::kMult, (*fPXcov) ); // (identity-KH)fPXcov
3ebfbf52 707 if (!IsPositiveDefinite(ikhp)) return kFALSE;
708
709 (*fPX) += xupdate;
710 TMatrixDSymFromTMatrixD( (*fPXcov), ikhp ); //make the matrix completely symetrical
711
55d5da9e 712 fNUpdates++;
3ebfbf52 713
55d5da9e 714 return kTRUE;
043badeb 715}
716
044eb03e 717//______________________________________________________________________________
718Bool_t AliRelAlignerKalman::IsOutlier( const TVectorD& update, const TMatrixDSym& covmatrix )
719{
55d5da9e 720 //check whether an update is an outlier given the covariance matrix of the fit
721
722 Bool_t is=kFALSE;
723 for (Int_t i=0;i<fgkNSystemParams;i++)
724 {
3ebfbf52 725 if (covmatrix(i,i)<0.) return kTRUE; //if cov matrix has neg diagonals something went wrong
55d5da9e 726 is = (is) || (TMath::Abs(update(i)) > fOutRejSigmas*TMath::Sqrt((covmatrix)(i,i)));
727 }
728 return is;
044eb03e 729}
730
3ebfbf52 731//______________________________________________________________________________
732Bool_t AliRelAlignerKalman::IsPositiveDefinite( const TMatrixD& mat ) const
733{
18f65700 734 //check for positive definiteness
735
3ebfbf52 736 for (Int_t i=0; i<mat.GetNcols(); i++)
737 {
738 if (mat(i,i)<=0.) return kFALSE;
739 }
740
741 if (!fNumericalParanoia) return kTRUE;
742
743 TDecompLU lu(mat);
744 return (lu.Decompose());
745}
746
044eb03e 747//______________________________________________________________________________
55d5da9e 748void AliRelAlignerKalman::TMatrixDSymFromTMatrixD( TMatrixDSym& matsym, const TMatrixD& mat )
043badeb 749{
55d5da9e 750 //Produce a valid symmetric matrix out of an almost symmetric TMatrixD
043badeb 751
55d5da9e 752 for (Int_t i=0; i<mat.GetNcols(); i++)
753 {
754 matsym(i,i) = mat(i,i); //copy diagonal
755 for (Int_t j=i+1; j<mat.GetNcols(); j++)
043badeb 756 {
55d5da9e 757 //copy the rest
758 Double_t average = (mat(i,j)+mat(j,i))/2.;
759 matsym(i,j)=average;
760 matsym(j,i)=average;
043badeb 761 }
55d5da9e 762 }
3ebfbf52 763 matsym.MakeValid();
55d5da9e 764 return;
043badeb 765}
766
044eb03e 767//______________________________________________________________________________
043badeb 768void AliRelAlignerKalman::Angles( TVectorD &angles, const TMatrixD &rotmat )
769{
55d5da9e 770 //Calculate the Cardan angles (psi,theta,phi) from rotation matrix
771 //b = R*a
772 angles(0) = TMath::ATan2( -rotmat(1,2), rotmat(2,2) );
773 angles(1) = TMath::ASin( rotmat(0,2) );
774 angles(2) = TMath::ATan2( -rotmat(0,1), rotmat(0,0) );
775 return;
043badeb 776}
777
044eb03e 778//______________________________________________________________________________
043badeb 779void AliRelAlignerKalman::PrintCorrelationMatrix()
780{
55d5da9e 781 //Print the correlation matrix for the fitted parameters
782 printf("Correlation matrix for system parameters:\n");
783 for ( Int_t i=0; i<fgkNSystemParams; i++ )
784 {
785 for ( Int_t j=0; j<i+1; j++ )
043badeb 786 {
3ebfbf52 787 if ((*fPXcov)(i,i)==0. || (*fPXcov)(j,j)==0.) printf(" NaN ");
788 else
789 printf("% -1.3f ", (*fPXcov)(i,j)/TMath::Sqrt( (*fPXcov)(i,i) * (*fPXcov)(j,j) ) );
55d5da9e 790 }//for j
043badeb 791 printf("\n");
55d5da9e 792 }//for i
793 printf("\n");
794 return;
043badeb 795}
796
044eb03e 797//______________________________________________________________________________
c3b5bfc1 798Bool_t AliRelAlignerKalman::FindCosmicTrackletNumbersInEvent( TArrayI& outITSindexTArr, TArrayI& outTPCindexTArr, const AliESDEvent* pEvent )
043badeb 799{
55d5da9e 800 //Find matching track segments in an event with tracks in TPC and ITS(standalone)
043badeb 801
55d5da9e 802 //Sanity cuts on tracks + check which tracks are ITS which are TPC
803 Int_t ntracks = pEvent->GetNumberOfTracks(); ////printf("number of tracks in event: %i\n", ntracks);
3ebfbf52 804 fMagField = pEvent->GetMagneticField();
55d5da9e 805 if (ntracks<2)
806 {
807 //printf("TrackFinder: less than 2 tracks!\n");
808 return kFALSE;
809 }
810 Float_t* phiArr = new Float_t[ntracks];
811 Float_t* thetaArr = new Float_t[ntracks];
812 Int_t* goodtracksArr = new Int_t[ntracks];
813 Int_t* candidateTPCtracksArr = new Int_t[ntracks];
814 Int_t* matchedITStracksArr = new Int_t[ntracks];
815 Int_t* matchedTPCtracksArr = new Int_t[ntracks];
816 Int_t* tracksArrITS = new Int_t[ntracks];
817 Int_t* tracksArrTPC = new Int_t[ntracks];
818 Int_t* nPointsArr = new Int_t[ntracks];
819 Int_t nITStracks = 0;
820 Int_t nTPCtracks = 0;
821 Int_t nGoodTracks = 0;
822 Int_t nCandidateTPCtracks = 0;
823 Int_t nMatchedITStracks = 0;
824 AliESDtrack* pTrack = NULL;
825 Bool_t foundMatchTPC = kFALSE;
826
827 //select and clasify tracks
828 for (Int_t itrack=0; itrack < ntracks; itrack++)
829 {
830 pTrack = pEvent->GetTrack(itrack);
831 if (!pTrack)
043badeb 832 {
3ebfbf52 833 //std::cout<<"no track!"<<std::endl;
55d5da9e 834 continue;
043badeb 835 }
55d5da9e 836 if (fCuts)
c3b5bfc1 837 {
55d5da9e 838 if (pTrack->GetP()<fMinMom || pTrack->GetP()>fMaxMom) continue;
043badeb 839 }
55d5da9e 840 goodtracksArr[nGoodTracks]=itrack;
841 Float_t phi = pTrack->GetAlpha()+TMath::ASin(pTrack->GetSnp());
842 Float_t theta = 0.5*TMath::Pi()-TMath::ATan(pTrack->GetTgl());
843 phiArr[nGoodTracks]=phi;
844 thetaArr[nGoodTracks]=theta;
845
846 //check if track is ITS
847 Int_t nClsITS = pTrack->GetNcls(0);
848 Int_t nClsTPC = pTrack->GetNcls(1);
849 if ( ((pTrack->GetStatus()&AliESDtrack::kITSout)>0)&&
850 !((pTrack->GetStatus()&AliESDtrack::kTPCin)>0)&&
851 !(nClsITS<fMinPointsVol1) ) //enough points
043badeb 852 {
55d5da9e 853 tracksArrITS[nITStracks] = nGoodTracks;
854 nITStracks++;
855 nGoodTracks++;
856 continue;
043badeb 857 }
858
55d5da9e 859 //check if track is TPC
860 if ( ((pTrack->GetStatus()&AliESDtrack::kTPCin)>0)&&
861 !(nClsTPC<fMinPointsVol2) ) //enough points
043badeb 862 {
55d5da9e 863 tracksArrTPC[nTPCtracks] = nGoodTracks;
864 nTPCtracks++;
865 nGoodTracks++;
866 //printf("tracksArrTPC[%d]=%d, goodtracksArr[%d]=%d\n",nTPCtracks-1,tracksArrTPC[nTPCtracks-1],nGoodTracks-1,goodtracksArr[nGoodTracks-1]);
867 continue;
044eb03e 868 }
55d5da9e 869 }//for itrack - selection fo tracks
044eb03e 870
55d5da9e 871 //printf("TrackFinder: %d ITS | %d TPC out of %d tracks in event\n", nITStracks,nTPCtracks,ntracks);
872
873 if ( nITStracks < 1 || nTPCtracks < 1 )
874 {
875 delete [] phiArr;
876 delete [] thetaArr;
877 delete [] goodtracksArr;
878 delete [] matchedITStracksArr;
879 delete [] candidateTPCtracksArr;
880 delete [] matchedTPCtracksArr;
881 delete [] tracksArrITS;
882 delete [] tracksArrTPC;
883 delete [] nPointsArr;
884 return kFALSE;
885 }
886
887 //find matching in TPC
888 if (nTPCtracks>1) //if there is something to be matched, try and match it
889 {
890 Float_t min = 10000000.;
891 for (Int_t itr1=0; itr1<nTPCtracks; itr1++)
044eb03e 892 {
55d5da9e 893 for (Int_t itr2=itr1+1; itr2<nTPCtracks; itr2++)
894 {
895 Float_t deltatheta = TMath::Abs(TMath::Pi()-thetaArr[tracksArrTPC[itr1]]-thetaArr[tracksArrTPC[itr2]]);
896 if (deltatheta > fMaxMatchingAngle) continue;
897 Float_t deltaphi = TMath::Abs(TMath::Abs(phiArr[tracksArrTPC[itr1]]-phiArr[tracksArrTPC[itr2]])-TMath::Pi());
898 if (deltaphi > fMaxMatchingAngle) continue;
899 if (deltatheta+deltaphi<min) //only the best matching pair
043badeb 900 {
55d5da9e 901 min=deltatheta+deltaphi;
902 candidateTPCtracksArr[0] = tracksArrTPC[itr1]; //store the index of track in goodtracksArr[]
903 candidateTPCtracksArr[1] = tracksArrTPC[itr2];
904 nCandidateTPCtracks = 2;
905 foundMatchTPC = kTRUE;
906 //printf("TrackFinder: Matching TPC tracks candidates:\n");
907 //printf("TrackFinder: candidateTPCtracksArr[0]=%d\n",tracksArrTPC[itr1]);
908 //printf("TrackFinder: candidateTPCtracksArr[1]=%d\n",tracksArrTPC[itr2]);
043badeb 909 }
55d5da9e 910 }
043badeb 911 }
55d5da9e 912 }//if nTPCtracks>1
913 else //if nTPCtracks==1 - if nothing to match, take the only one we've got
914 {
915 candidateTPCtracksArr[0] = tracksArrTPC[0];
916 nCandidateTPCtracks = 1;
917 foundMatchTPC = kFALSE;
918 }
919 if (foundMatchTPC) fNMatchedTPCtracklets++;
920 //if no match but the requirement is set return kFALSE
921 if (fRequireMatchInTPC && !foundMatchTPC)
922 {
923 delete [] phiArr;
924 delete [] thetaArr;
925 delete [] goodtracksArr;
926 delete [] candidateTPCtracksArr;
927 delete [] matchedITStracksArr;
928 delete [] matchedTPCtracksArr;
929 delete [] tracksArrITS;
930 delete [] tracksArrTPC;
931 delete [] nPointsArr;
932 //printf("TrackFinder: no match in TPC && required\n");
933 return kFALSE;
934 }
043badeb 935
55d5da9e 936 //if no match and more than one track take all TPC tracks
937 if (!fRequireMatchInTPC && !foundMatchTPC)
938 {
939 for (Int_t i=0;i<nTPCtracks;i++)
c3b5bfc1 940 {
55d5da9e 941 candidateTPCtracksArr[i] = tracksArrTPC[i];
c3b5bfc1 942 }
55d5da9e 943 nCandidateTPCtracks = nTPCtracks;
944 }
945 //printf("TrackFinder: nCandidateTPCtracks: %i\n", nCandidateTPCtracks);
946
947 Double_t* minDifferenceArr = new Double_t[nCandidateTPCtracks];
948
949 //find ITS matches for good TPC tracks
950 Bool_t matchedITStracks=kFALSE;
951 for (Int_t itpc=0;itpc<nCandidateTPCtracks;itpc++)
952 {
953 minDifferenceArr[nMatchedITStracks] = 10000000.;
954 matchedITStracks=kFALSE;
955 for (Int_t iits=0; iits<nITStracks; iits++)
c3b5bfc1 956 {
55d5da9e 957 AliESDtrack* itstrack = pEvent->GetTrack(goodtracksArr[tracksArrITS[iits]]);
958 const AliExternalTrackParam* parits = itstrack->GetOuterParam();
959 AliESDtrack* tpctrack = pEvent->GetTrack(goodtracksArr[candidateTPCtracksArr[itpc]]);
960 const AliExternalTrackParam* tmp = tpctrack->GetInnerParam();
961 AliExternalTrackParam partpc(*tmp); //make a copy to avoid tampering with original params
962 partpc.Rotate(parits->GetAlpha());
3ebfbf52 963 partpc.PropagateTo(parits->GetX(),fMagField);
55d5da9e 964 Float_t dtgl = TMath::Abs(partpc.GetTgl()-parits->GetTgl());
965 if (dtgl > fMaxMatchingAngle) continue;
966 Float_t dsnp = TMath::Abs(partpc.GetSnp()-parits->GetSnp());
967 if (dsnp > fMaxMatchingAngle) continue;
968 Float_t dy = TMath::Abs(partpc.GetY()-parits->GetY());
969 Float_t dz = TMath::Abs(partpc.GetZ()-parits->GetZ());
970 if (TMath::Sqrt(dy*dy+dz*dz) > fMaxMatchingDistance) continue;
971 if (dtgl+dsnp<minDifferenceArr[nMatchedITStracks]) //only the best matching pair
972 {
973 minDifferenceArr[nMatchedITStracks]=dtgl+dsnp;
974 matchedITStracksArr[nMatchedITStracks] = tracksArrITS[iits];
975 matchedTPCtracksArr[nMatchedITStracks] = candidateTPCtracksArr[itpc]; //this tells us minDifferenceArrwhich TPC track this ITS track belongs to
976 //printf("TrackFinder: Matching ITS to TPC:\n");
977 //printf("TrackFinder: minDifferenceArr[%i]=%.2f\n",nMatchedITStracks,minDifferenceArr[nMatchedITStracks]);
978 //printf("TrackFinder: matchedITStracksArr[%i]=%i\n",nMatchedITStracks,matchedITStracksArr[nMatchedITStracks]);
979 //printf("TrackFinder: matchedTPCtracksArr[%i]=%i\n",nMatchedITStracks,matchedTPCtracksArr[nMatchedITStracks]);
980 matchedITStracks=kTRUE;;
981 }
c3b5bfc1 982 }
55d5da9e 983 if (matchedITStracks) nMatchedITStracks++;
984 }
044eb03e 985
55d5da9e 986 if (nMatchedITStracks==0) //no match in ITS
987 {
c3b5bfc1 988 delete [] phiArr;
989 delete [] thetaArr;
990 delete [] minDifferenceArr;
044eb03e 991 delete [] goodtracksArr;
c3b5bfc1 992 delete [] matchedITStracksArr;
55d5da9e 993 delete [] candidateTPCtracksArr;
c3b5bfc1 994 delete [] matchedTPCtracksArr;
55d5da9e 995 delete [] tracksArrITS;
996 delete [] tracksArrTPC;
044eb03e 997 delete [] nPointsArr;
55d5da9e 998 //printf("TrackFinder: No match in ITS\n");
999 return kFALSE;
1000 }
1001
1002 //printf("TrackFinder: nMatchedITStracks: %i\n",nMatchedITStracks);
1003 //we found a cosmic
1004 fNMatchedCosmics++;
1005
1006 //Now we may have ended up with more matches than we want in the case there was
1007 //no TPC match and there were many TPC tracks
1008 //a cosmic in a scenario like this will only have produced 1 pair, the rest is garbage
1009 //so take only the best pair.
1010 //same applies when there are more matches than ITS tracks - means that one ITS track
1011 //matches more TPC tracks.
1012 if ((nMatchedITStracks>2 && !foundMatchTPC) || nMatchedITStracks>nITStracks)
1013 {
1014 Int_t imin = TMath::LocMin(nMatchedITStracks,minDifferenceArr);
1015 matchedITStracksArr[0] = matchedITStracksArr[imin];
1016 matchedTPCtracksArr[0] = matchedTPCtracksArr[imin];
1017 nMatchedITStracks = 1;
1018 //printf("TrackFinder: too many matches - take only the best one\n");
1019 //printf("TrackFinder: LocMin in matched its tracks: %d\n",imin);
1020 //printf("TrackFinder: matchedITStracksArr[0]=%d\n",matchedITStracksArr[0]);
1021 //printf("TrackFinder: matchedTPCtracksArr[0]=%d\n",matchedTPCtracksArr[0]);
1022 }
1023
1024 ///////////////////////////////////////////////////////////////////////////
1025 outITSindexTArr.Set(nMatchedITStracks);
1026 outTPCindexTArr.Set(nMatchedITStracks);
1027 for (Int_t i=0;i<nMatchedITStracks;i++)
1028 {
1029 outITSindexTArr.AddAt( goodtracksArr[matchedITStracksArr[i]], i );
1030 outTPCindexTArr.AddAt( goodtracksArr[matchedTPCtracksArr[i]], i );
1031 //printf("TrackFinder: Fill the output\n");
1032 //printf("TrackFinder: matchedITStracksArr[%d]=%d\n",i,matchedITStracksArr[i]);
1033 //printf("TrackFinder: matchedTPCtracksArr[%d]=%d\n",i,matchedTPCtracksArr[i]);
1034 }
1035 //printf("TrackFinder: Size of outputarrays: %d, %d\n", outITSindexTArr.GetSize(), outTPCindexTArr.GetSize());
1036 ///////////////////////////////////////////////////////////////////////////
1037
1038 delete [] phiArr;
1039 delete [] thetaArr;
1040 delete [] minDifferenceArr;
1041 delete [] goodtracksArr;
1042 delete [] candidateTPCtracksArr;
1043 delete [] matchedITStracksArr;
1044 delete [] matchedTPCtracksArr;
1045 delete [] tracksArrITS;
1046 delete [] tracksArrTPC;
1047 delete [] nPointsArr;
1048 return kTRUE;
043badeb 1049}
55d5da9e 1050
044eb03e 1051//______________________________________________________________________________
3ebfbf52 1052Bool_t AliRelAlignerKalman::CorrectTrack( AliExternalTrackParam* tr, const TVectorD& misal ) const
043badeb 1053{
55d5da9e 1054 //implements the system model -
1055 //applies correction for misalignment and calibration to track
3ebfbf52 1056 //track needs to be already propagated to the global reference plane
55d5da9e 1057
1058 Double_t x = tr->GetX();
1059 Double_t alpha = tr->GetAlpha();
1060 Double_t point[3],dir[3];
1061 tr->GetXYZ(point);
1062 tr->GetDirection(dir);
1063 TVector3 Point(point);
1064 TVector3 Dir(dir);
1065
1066 //Apply corrections to track
1067
1068 //Shift
3ebfbf52 1069 Point(0) -= misal(3); //add shift in x
1070 Point(1) -= misal(4); //add shift in y
1071 Point(2) -= misal(5); //add shift in z
55d5da9e 1072 //Rotation
1073 TMatrixD rotmat(3,3);
1074 RotMat( rotmat, misal );
3ebfbf52 1075 Point = rotmat.T() * Point;
55d5da9e 1076 Dir = rotmat * Dir;
1077
1078 //TPC vdrift and T0 corrections
3ebfbf52 1079 TVector3 Point2; //second point of the track
1080 Point2 = Point + Dir;
1081 Double_t vdCorr = 1./misal(6);
1082 Double_t t0 = misal(7);
1083 Double_t vdY = 0.0;
1084 if (fgkNSystemParams>8) vdY = misal(8)/100.; //change over 100cm.
1085
1086 //my model
55d5da9e 1087 if (Point(2)>0)
1088 {
1089 //A-Side
1090 Point(2) = Point(2) - (fTPCZLengthA-Point(2)) * (vdCorr-1.+vdY*Point(1)/fTPCvd) - (fTPCvd*vdCorr+vdY*Point(1))*t0;
1091 Point2(2) = Point2(2) - (fTPCZLengthA-Point2(2)) * (vdCorr-1.+vdY*Point2(1)/fTPCvd) - (fTPCvd*vdCorr+vdY*Point2(1))*t0;
1092 }
1093 else
1094 {
1095 //C-side
1096 Point(2) = Point(2) - (fTPCZLengthC+Point(2)) * (1.-vdCorr-vdY*Point(1)/fTPCvd) + (fTPCvd*vdCorr+vdY*Point(1))*t0;
1097 Point2(2) = Point2(2) - (fTPCZLengthC+Point2(2)) * (1.-vdCorr-vdY*Point2(1)/fTPCvd) + (fTPCvd*vdCorr+vdY*Point2(1))*t0;
1098 }
3ebfbf52 1099
1100 //Stefan's model
1101 //if (Point(2)>0)
1102 //{
1103 // //A-Side
1104 // Point(2) = Point(2) - (fTPCZLengthA-Point(2)) * (1.-vdCorr+vdY*Point(1)/fTPCvd) - (fTPCvd*vdCorr+vdY*Point(1))*t0;
1105 // Point2(2) = Point2(2) - (fTPCZLengthA-Point2(2)) * (1.-vdCorr+vdY*Point2(1)/fTPCvd) - (fTPCvd*vdCorr+vdY*Point2(1))*t0;
1106 //}
1107 //else
1108 //{
1109 // //C-side
1110 // Point(2) = Point(2) + (fTPCZLengthC+Point(2)) * (1.-vdCorr+vdY*Point(1)/fTPCvd) + (fTPCvd*vdCorr+vdY*Point(1))*t0;
1111 // Point2(2) = Point2(2) + (fTPCZLengthC+Point2(2)) * (1.-vdCorr+vdY*Point2(1)/fTPCvd) + (fTPCvd*vdCorr+vdY*Point2(1))*t0;
1112 //}
1113
55d5da9e 1114 Dir = Point2-Point;
1115 Dir=Dir.Unit(); //keep unit length
1116
1117 //Turn back to local system
1118 Dir.GetXYZ(dir);
1119 Point.GetXYZ(point);
1120 tr->Global2LocalPosition(point,alpha);
1121 tr->Global2LocalPosition(dir,alpha);
1122
1123 //Calculate new intersection point with ref plane
1124 Double_t p[5],pcov[15];
1125 if (dir[0]==0.0) return kFALSE;
1126 Double_t s=(x-point[0])/dir[0];
1127 p[0] = point[1]+s*dir[1];
1128 p[1] = point[2]+s*dir[2];
1129 Double_t pt = TMath::Sqrt(dir[0]*dir[0]+dir[1]*dir[1]);
1130 if (pt==0.0) return kFALSE;
1131 p[2] = dir[1]/pt;
1132 p[3] = dir[2]/pt;
55d5da9e 1133 //insert everything back into track
1134 const Double_t* pcovtmp = tr->GetCovariance();
3ebfbf52 1135 p[4] = tr->GetSigned1Pt(); //copy the momentum
55d5da9e 1136 memcpy(pcov,pcovtmp,15*sizeof(Double_t));
1137 tr->Set(x,alpha,p,pcov);
55d5da9e 1138 return kTRUE;
3ebfbf52 1139
1140 ////put params back into track and propagate to ref
1141 //Double_t p[5],pcov[15];
1142 //p[0] = point[1];
1143 //p[1] = point[2];
1144 //Double_t xnew = point[0];
1145 //Double_t pt = TMath::Sqrt(dir[0]*dir[0]+dir[1]*dir[1]);
1146 //if (pt==0.0) return kFALSE;
1147 //p[2] = dir[1]/pt;
1148 //p[3] = dir[2]/pt;
1149 //p[4] = tr->GetSigned1Pt(); //copy the momentum
1150 //const Double_t* pcovtmp = tr->GetCovariance();
1151 //memcpy(pcov,pcovtmp,15*sizeof(Double_t));
1152 //tr->Set(xnew,alpha,p,pcov);
1153 //return tr->PropagateTo(x,fMagField);
043badeb 1154}
1155
044eb03e 1156//______________________________________________________________________________
3ebfbf52 1157Bool_t AliRelAlignerKalman::MisalignTrack( AliExternalTrackParam* tr, const TVectorD& misal ) const
043badeb 1158{
55d5da9e 1159 //implements the system model -
1160 //applies misalignment and miscalibration to reference track
3ebfbf52 1161 //trackparams have to be at the global reference plane
55d5da9e 1162
1163 Double_t x = tr->GetX();
1164 Double_t alpha = tr->GetAlpha();
1165 Double_t point[3],dir[3];
1166 tr->GetXYZ(point);
1167 tr->GetDirection(dir);
1168 TVector3 Point(point);
1169 TVector3 Dir(dir);
1170
1171 //Apply misalignment to track
1172
1173 //TPC vdrift and T0 corrections
3ebfbf52 1174 TVector3 Point2; //second point of the track
1175 Point2 = Point + Dir;
1176 Double_t vdCorr = 1./misal(6);
1177 Double_t t0 = misal(7);
1178 Double_t vdY = 0.0;
1179 if (fgkNSystemParams>8) vdY = misal(8)/100.; //change over 100cm.
1180
55d5da9e 1181 if (Point(2)>0)
1182 {
1183 //A-Side
1184 Point(2) = Point(2) + ((fTPCZLengthA-Point(2))/(vdCorr*fTPCvd+vdY*Point(1)))
1185 * (fTPCvd*(vdCorr-1.)+vdY*Point(1)) + fTPCvd*t0;
1186 Point2(2) = Point2(2) + ((fTPCZLengthA-Point2(2))/(vdCorr*fTPCvd+vdY*Point2(1)))
1187 * (fTPCvd*(vdCorr-1.)+vdY*Point2(1)) + fTPCvd*t0;
1188 }
1189 else
1190 {
1191 //C-side
1192 Point(2) = Point(2) + (fTPCZLengthC+Point(2))/(vdCorr*fTPCvd+vdY*Point(1))
1193 * (fTPCvd*(1.-vdCorr)-vdY*Point(1)) - fTPCvd*t0;
1194 Point2(2) = Point2(2) + (fTPCZLengthC+Point2(2))/(vdCorr*fTPCvd+vdY*Point2(1))
1195 * (fTPCvd*(1.-vdCorr)-vdY*Point2(1)) - fTPCvd*t0;
1196 }
1197 Dir = Point2-Point;
1198 Dir=Dir.Unit(); //keep unit length
1199
1200 //Rotation
1201 TMatrixD rotmat(3,3);
1202 RotMat( rotmat, misal );
1203 Point = rotmat * Point;
1204 Dir = rotmat * Dir;
1205 //Shift
1206 Point(0) += misal(3); //add shift in x
1207 Point(1) += misal(4); //add shift in y
1208 Point(2) += misal(5); //add shift in z
1209
1210 //Turn back to local system
1211 Dir.GetXYZ(dir);
1212 Point.GetXYZ(point);
1213 tr->Global2LocalPosition(point,alpha);
1214 tr->Global2LocalPosition(dir,alpha);
1215
1216 //Calculate new intersection point with ref plane
1217 Double_t p[5],pcov[15];
3ebfbf52 1218 if (dir[0]==0.0) return kFALSE;
55d5da9e 1219 Double_t s=(x-point[0])/dir[0];
1220 p[0] = point[1]+s*dir[1];
1221 p[1] = point[2]+s*dir[2];
1222 Double_t pt = TMath::Sqrt(dir[0]*dir[0]+dir[1]*dir[1]);
3ebfbf52 1223 if (pt==0.0) return kFALSE;
55d5da9e 1224 p[2] = dir[1]/pt;
1225 p[3] = dir[2]/pt;
55d5da9e 1226 //insert everything back into track
1227 const Double_t* pcovtmp = tr->GetCovariance();
3ebfbf52 1228 p[4] = tr->GetSigned1Pt(); //copy the momentum
55d5da9e 1229 memcpy(pcov,pcovtmp,15*sizeof(Double_t));
1230 tr->Set(x,alpha,p,pcov);
55d5da9e 1231 return kTRUE;
3ebfbf52 1232
1233 ////put params back into track and propagate to ref
1234 //Double_t p[5];
1235 //Double_t pcov[15];
1236 //p[0] = point[1];
1237 //p[1] = point[2];
1238 //Double_t xnew = point[0];
1239 //Double_t pt = TMath::Sqrt(dir[0]*dir[0]+dir[1]*dir[1]);
1240 //if (pt==0.0) return kFALSE;
1241 //p[2] = dir[1]/pt;
1242 //p[3] = dir[2]/pt;
1243 //p[4] = tr->GetSigned1Pt(); //copy the momentum
1244 //const Double_t* pcovtmp = tr->GetCovariance();
1245 //memcpy(pcov,pcovtmp,15*sizeof(Double_t));
1246 //printf("x before: %.5f, after: %.5f\n",x, xnew);
1247 //printf("before: %.4f %.4f %.4f %.4f %.4f \n",tr->GetParameter()[0],tr->GetParameter()[1],tr->GetParameter()[2],tr->GetParameter()[3],tr->GetParameter()[4]);
1248 //printf("after: %.4f %.4f %.4f %.4f %.4f \n",p[0],p[1],p[2],p[3],p[4]);
1249 //tr->Set(xnew,alpha,p,pcov);
1250 //return tr->PropagateTo(x,fMagField);
043badeb 1251}
1252
044eb03e 1253//______________________________________________________________________________
55d5da9e 1254void AliRelAlignerKalman::Reset()
043badeb 1255{
3ebfbf52 1256 //full reset to defaults
55d5da9e 1257 fPX->Zero();
1258 (*fPX)(6)=1.;
1259 ResetCovariance();
3ebfbf52 1260
1261 //initialize the differentials per parameter
1262 for (Int_t i=0;i<fgkNSystemParams;i++) fDelta[i] = 1.e-6;
1263
1264 fNMatchedCosmics=0;
1265 fNMatchedTPCtracklets=0;
1266 fNUpdates=0;
1267 fNOutliers=0;
1268 fNTracks=0;
1269 fNProcessedEvents=0;
329b5744 1270 fRunNumber=0;
1271 fTimeStamp=0;
043badeb 1272}
1273
044eb03e 1274//______________________________________________________________________________
c3b5bfc1 1275void AliRelAlignerKalman::ResetCovariance( const Double_t number )
043badeb 1276{
55d5da9e 1277 //Resets the covariance to the default if arg=0 or resets the off diagonals
1278 //to zero and releases the diagonals by factor arg.
1279 if (number!=0.)
1280 {
3ebfbf52 1281 for (Int_t z=0;z<6;z++)
c3b5bfc1 1282 {
3ebfbf52 1283 for (Int_t zz=0;zz<6;zz++)
55d5da9e 1284 {
1285 if (zz==z) continue; //don't touch diagonals
1286 (*fPXcov)(zz,z) = 0.;
1287 (*fPXcov)(z,zz) = 0.;
1288 }
1289 (*fPXcov)(z,z) = (*fPXcov)(z,z) * number;
c3b5bfc1 1290 }
55d5da9e 1291 }
1292 else
1293 {
1294 //Resets the covariance of the fit to a default value
1295 fPXcov->Zero();
3ebfbf52 1296 (*fPXcov)(0,0) = .08*.08; //psi (rad)
1297 (*fPXcov)(1,1) = .08*.08; //theta (rad
1298 (*fPXcov)(2,2) = .08*.08; //phi (rad)
1299 (*fPXcov)(3,3) = .3*.3; //x (cm)
1300 (*fPXcov)(4,4) = .3*.3; //y (cm)
1301 (*fPXcov)(5,5) = .3*.3; //z (cm)
55d5da9e 1302 }
3ebfbf52 1303 ResetTPCparamsCovariance(number);
55d5da9e 1304}
1305
1306//______________________________________________________________________________
1307void AliRelAlignerKalman::ResetTPCparamsCovariance( const Double_t number )
1308{
1309 //Resets the covariance to the default if arg=0 or resets the off diagonals
1310 //to zero and releases the diagonals by factor arg.
1311
1312 //release diagonals
1313 if (number==0.)
1314 {
3ebfbf52 1315 if (fgkNSystemParams>6) (*fPXcov)(6,6) = .1*.1;
1316 if (fgkNSystemParams>7) (*fPXcov)(7,7) = 1.*1.;
a80268df 1317 if (fgkNSystemParams>8) (*fPXcov)(8,8) = .1*.1;
55d5da9e 1318 }
1319 else
1320 {
3ebfbf52 1321 if (fgkNSystemParams>6) (*fPXcov)(6,6) = number * (*fPXcov)(6,6);
1322 if (fgkNSystemParams>7) (*fPXcov)(7,7) = number * (*fPXcov)(7,7);
1323 if (fgkNSystemParams>8) (*fPXcov)(8,8) = number * (*fPXcov)(8,8);
55d5da9e 1324 }
1325
1326 //set crossterms to zero
1327 for (Int_t i=0;i<fgkNSystemParams;i++)
1328 {
3ebfbf52 1329 for (Int_t j=6;j<fgkNSystemParams;j++) //TPC params
c3b5bfc1 1330 {
55d5da9e 1331 if (i==j) continue; //don't touch diagonals
1332 (*fPXcov)(i,j) = 0.;
1333 (*fPXcov)(j,i) = 0.;
c3b5bfc1 1334 }
55d5da9e 1335 }
043badeb 1336}
55d5da9e 1337
329b5744 1338//______________________________________________________________________________
a80268df 1339Bool_t AliRelAlignerKalman::Merge( const AliRelAlignerKalman* al )
1340{
1341 //Merge two aligners
1342
1343 if (!al) return kFALSE;
329b5744 1344 if (al==this) return kTRUE;
a80268df 1345 if (al->fgkNSystemParams != fgkNSystemParams) return kFALSE;
329b5744 1346 if (fRunNumber != al->fRunNumber) return kFALSE;
1347 if (al->fNUpdates == 0) return kTRUE; //no point in merging with an empty one
a80268df 1348
1349 //store the pointers to current stuff
1350 TVectorD* pmes = fPMeasurement;
1351 TMatrixDSym* pmescov = fPMeasurementCov;
1352 TVectorD* pmespred = fPMeasurementPrediction;
1353 TMatrixD* ph = fPH;
1354
1355 //make a unity system matrix
1356 TMatrixD tmp(fgkNSystemParams,fgkNSystemParams);
1357 fPH = new TMatrixD(TMatrixD::kUnit, tmp);
1358
1359 //mesurement is the state of the new aligner
1360 fPMeasurement = al->fPX;
1361 fPMeasurementCov = al->fPXcov;
1362
1363 //the mesurement prediction is the state
1364 fPMeasurementPrediction = fPX; //this is safe as fPX doesn't change until end
1365
1366 //do the merging
1367 Bool_t success = UpdateEstimateKalman();
1368
1369 //restore pointers to old stuff
1370 fPMeasurement = pmes;
1371 fPMeasurementCov = pmescov;
1372 fPMeasurementPrediction = pmespred;
1373 delete fPH;
1374 fPH = ph;
1375
1376 //merge stats
9b7d9478 1377 if (!success) return kFALSE; //no point in merging stats if merge not succesful
a80268df 1378 fNProcessedEvents += al->fNProcessedEvents;
1379 fNUpdates += al->fNUpdates;
1380 fNOutliers += al->fNOutliers;
1381 fNTracks += al->fNTracks;
1382 fNMatchedTPCtracklets += al->fNMatchedTPCtracklets;
1383 fNMatchedCosmics += al->fNMatchedCosmics;
329b5744 1384 if (fTimeStamp < al->fTimeStamp) fTimeStamp = al->fTimeStamp; //take the older one
a80268df 1385
1386 return success;
1387}
1388
329b5744 1389//______________________________________________________________________________
1390Int_t AliRelAlignerKalman::Compare(const TObject *obj) const
1391{
1392 if (this == obj) return 0;
1393 const AliRelAlignerKalman* aobj = dynamic_cast<const AliRelAlignerKalman*>(obj);
1394 if (!aobj) return 0;
1395 if (fTimeStamp < aobj->fTimeStamp) return -1;
1396 else if (fTimeStamp > aobj->fTimeStamp) return 1;
1397 else return 0;
1398}
1399
1400//______________________________________________________________________________
1401//Int_t AliRelAlignerKalman::Compare(const AliRelAlignerKalman *al) const
1402//{
1403// if (this == al) return 0;
1404// if (fTimeStamp > al->fTimeStamp) return -1;
1405// else if (fTimeStamp < al->fTimeStamp) return 1;
1406// else return 0;
1407//}
1408
3ebfbf52 1409//______________________________________________________________________________
1410//void AliRelAlignerKalman::PrintCovarianceCorrection()
1411//{
1412// //Print the measurement covariance correction matrix
1413// printf("Covariance correction matrix:\n");
16cbbc2c 1414// for ( Int_t i=0; i<fNMeasurementParams; i++ )
3ebfbf52 1415// {
1416// for ( Int_t j=0; j<i+1; j++ )
1417// {
1418// printf("% -2.2f ", (*fPMeasurementCovCorr)(i,j) );
1419// }//for i
1420// printf("\n");
1421// }//for j
1422// printf("\n");
1423// return;
1424//}
1425
1426//_______________________________________________________________________________
1427//Bool_t AliRelAlignerKalman::UpdateCalibration()
1428//{
1429// //Update the calibration with new data (in calibration mode)
1430//
1431// fPMes0Hist->Fill( (*fPMeasurement)(0) );
1432// fPMes1Hist->Fill( (*fPMeasurement)(1) );
1433// fPMes2Hist->Fill( (*fPMeasurement)(2) );
1434// fPMes3Hist->Fill( (*fPMeasurement)(3) );
1435// fPMesErr0Hist->Fill( TMath::Sqrt((*fPMeasurementCov)(0,0)) );
1436// fPMesErr1Hist->Fill( TMath::Sqrt((*fPMeasurementCov)(1,1)) );
1437// fPMesErr2Hist->Fill( TMath::Sqrt((*fPMeasurementCov)(2,2)) );
1438// fPMesErr3Hist->Fill( TMath::Sqrt((*fPMeasurementCov)(3,3)) );
1439// return kTRUE;
1440//}
1441
1442//______________________________________________________________________________
1443//Bool_t AliRelAlignerKalman::SetCalibrationMode( const Bool_t cp )
1444//{
1445// //sets the calibration mode
1446// if (cp)
1447// {
1448// fCalibrationMode=kTRUE;
1449// return kTRUE;
1450// }//if (cp)
1451// else
1452// {
1453// if (fCalibrationMode) // do it only after the calibration pass
1454// {
1455// CalculateCovarianceCorrection();
1456// SetApplyCovarianceCorrection();
1457// fCalibrationMode=kFALSE;
1458// return kTRUE;
1459// }//if (fCalibrationMode)
1460// }//else (cp)
1461// return kFALSE;
1462//}
1463
1464//______________________________________________________________________________
1465//Bool_t AliRelAlignerKalman::CalculateCovarianceCorrection()
1466//{
1467// //Calculates the correction to the measurement covariance
1468// //using the calibration histograms
1469//
1470// fPMeasurementCovCorr->Zero(); //reset the correction
1471//
1472// Double_t s,m,c; //sigma,meansigma,correction
1473//
1474// //TF1* fitformula;
1475// //fPMes0Hist->Fit("gaus");
1476// //fitformula = fPMes0Hist->GetFunction("gaus");
1477// //s = fitformula->GetParameter(2); //spread of the measurement
1478// //fPMesErr0Hist->Fit("gaus");
1479// //fitformula = fPMesErr0Hist->GetFunction("gaus"); //average error from cov matrices
1480// //m = fitformula->GetParameter(1);
1481// s = fPMes0Hist->GetRMS();
1482// m = fPMesErr0Hist->GetMean();
1483// c = s-m; //the difference between the average error and real spread of the data
1484// if (c>0) //only correct is spread bigger than average error
1485// (*fPMeasurementCovCorr)(0,0) = c*c;
1486//
1487// //fPMes1Hist->Fit("gaus");
1488// //fitformula = fPMes1Hist->GetFunction("gaus");
1489// //s = fitformula->GetParameter(2);
1490// //fPMesErr1Hist->Fit("gaus");
1491// //fitformula = fPMesErr1Hist->GetFunction("gaus");
1492// //m = fitformula->GetParameter(1);
1493// s = fPMes1Hist->GetRMS();
1494// m = fPMesErr1Hist->GetMean();
1495// c = s-m;
1496// if (c>0) //only correct is spread bigger than average error
1497// (*fPMeasurementCovCorr)(1,1) = c*c;
1498//
1499// //fPMes2Hist->Fit("gaus");
1500// //fitformula = fPMes2Hist->GetFunction("gaus");
1501// //s = fitformula->GetParameter(2);
1502// //fPMesErr2Hist->Fit("gaus");
1503// //fitformula = fPMesErr2Hist->GetFunction("gaus");
1504// //m = fitformula->GetParameter(1);
1505// s = fPMes2Hist->GetRMS();
1506// m = fPMesErr2Hist->GetMean();
1507// c = s-m;
1508// if (c>0) //only correct is spread bigger than average error
1509// (*fPMeasurementCovCorr)(2,2) = c*c;
1510//
1511// //fPMes3Hist->Fit("gaus");
1512// //fitformula = fPMes3Hist->GetFunction("gaus");
1513// //s = fitformula->GetParameter(2);
1514// //fPMesErr3Hist->Fit("gaus");
1515// //fitformula = fPMesErr3Hist->GetFunction("gaus");
1516// //m = fitformula->GetParameter(1);
1517// s = fPMes3Hist->GetRMS();
1518// m = fPMesErr3Hist->GetMean();
1519// c = s-m;
1520// if (c>0) //only correct is spread bigger than average error
1521// (*fPMeasurementCovCorr)(3,3) = c*c;
1522//
1523// return kTRUE;
1524//}
1525
16cbbc2c 1526