]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliTrackFitterRieman.cxx
Adding includes required by ROOT
[u/mrichter/AliRoot.git] / STEER / AliTrackFitterRieman.cxx
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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //
20 // Class to the track points on the Riemann sphere. Inputs are
21 // the set of id's (volids) of the volumes in which residuals are
22 // calculated to construct a chi2 function to be minimized during 
23 // the alignment procedures. For the moment the track extrapolation is 
24 // taken at the space-point reference plane. The reference plane is
25 // found using the covariance matrix of the point
26 // (assuming sigma(x)=0 at the reference coordinate system.
27 //
28 // Internal usage of AliRieman class for minimization
29 //
30 //////////////////////////////////////////////////////////////////////////////
31
32 #include <TArrayI.h>
33 #include <TLinearFitter.h>
34 #include <TMath.h>
35 #include <TMatrixD.h>
36 #include <TMatrixDSym.h>
37 #include <TRandom.h>
38 #include <TTreeStream.h>
39
40 #include "AliLog.h"
41 #include "AliLog.h"
42 #include "AliRieman.h"
43 #include "AliTrackFitterRieman.h"
44
45 ClassImp(AliTrackFitterRieman)
46
47
48 AliTrackFitterRieman::AliTrackFitterRieman():
49   AliTrackFitter(),
50   fBCorrection(kFALSE),
51   fAlpha(0.),
52   fNUsed(0),
53   fConv(kFALSE),
54   fMaxDelta(3),
55   fRieman(new AliRieman(10000)),  // allocate rieman
56   fDebugStream(new TTreeSRedirector("RiemanAlignDebug.root"))
57 {
58   //
59   // default constructor
60   //
61 }
62
63
64 AliTrackFitterRieman::AliTrackFitterRieman(AliTrackPointArray *array, Bool_t owner):
65   AliTrackFitter(array,owner),
66   fBCorrection(kFALSE),
67   fAlpha(0.),
68   fNUsed(0),
69   fConv(kFALSE),
70   fMaxDelta(3),
71   fRieman(new AliRieman(10000)),  //allocate rieman
72   fDebugStream(0)
73 {
74   //
75   // Constructor
76   //
77   if (AliLog::GetDebugLevel("","AliTrackFitterRieman")) fDebugStream = new TTreeSRedirector("RiemanAlignDebug.root");
78 }
79
80 AliTrackFitterRieman::AliTrackFitterRieman(const AliTrackFitterRieman &rieman):
81   AliTrackFitter(rieman),
82   fBCorrection(rieman.fBCorrection),
83   fAlpha(rieman.fAlpha),
84   fNUsed(rieman.fNUsed),
85   fConv(rieman.fConv),
86   fMaxDelta(rieman.fMaxDelta),
87   fRieman(new AliRieman(*(rieman.fRieman))),
88   fDebugStream(0)
89 {
90   //
91   // copy constructor
92   //
93   if (AliLog::GetDebugLevel("","AliTrackFitterRieman")) fDebugStream = new TTreeSRedirector("RiemanAlignDebug.root");
94 }
95
96 //_____________________________________________________________________________
97 AliTrackFitterRieman &AliTrackFitterRieman::operator =(const AliTrackFitterRieman& rieman)
98 {
99   //
100   // Assignment operator
101   //
102   if(this==&rieman) return *this;
103   ((AliTrackFitter *)this)->operator=(rieman);
104
105   fBCorrection = rieman.fBCorrection;
106   fAlpha  = rieman.fAlpha;
107   fNUsed  = rieman.fNUsed;
108   fConv   = rieman.fConv;
109   fMaxDelta = rieman.fMaxDelta;
110   fRieman = new AliRieman(*(rieman.fRieman));
111   fDebugStream = 0;
112   if (AliLog::GetDebugLevel("","AliTrackFitterRieman")) fDebugStream = new TTreeSRedirector("RiemanAlignDebug.root");
113   return *this;
114 }
115
116
117 AliTrackFitterRieman::~AliTrackFitterRieman(){
118   //
119   //
120   //
121   delete fRieman;
122   delete fDebugStream;
123 }
124
125 void AliTrackFitterRieman::Reset()
126 {
127   // Reset the track parameters and
128   // rieman sums
129   //
130   AliTrackFitter::Reset();
131   fRieman->Reset();
132   fAlpha = 0.;
133   fNUsed = 0;
134   fConv =kFALSE;
135 }
136
137 Bool_t AliTrackFitterRieman::Fit(const TArrayI *volIds,const TArrayI *volIdsFit,
138                                  AliAlignObj::ELayerID layerRangeMin,
139                          AliAlignObj::ELayerID layerRangeMax)
140 {
141   // Fit the track points. The method takes as an input
142   // the set of id's (volids) of the volumes in which
143   // one wants to calculate the residuals.
144   // The following parameters are used to define the
145   // range of volumes to be used in the fitting
146   // As a result two AliTrackPointArray's obects are filled.
147   // The first one contains the space points with
148   // volume id's from volids list. The second array of points represents
149   // the track extrapolations corresponding to the space points
150   // in the first array. The two arrays can be used to find
151   // the residuals in the volids and consequently construct a
152   // chi2 function to be minimized during the alignment
153   // procedures. For the moment the track extrapolation is taken
154   // at the space-point reference plane. The reference plane is
155   // found using the covariance matrix of the point
156   // (assuming sigma(x)=0 at the reference coordinate system.
157   Int_t debugLevel = AliLog::GetDebugLevel("","AliTrackFitterRieman");
158   Float_t debugRatio = 1./(1.+debugLevel);
159
160   const Int_t kMinPoints =1;
161   Int_t npoints = fPoints->GetNPoints();
162   if ( npoints<fMinNPoints) return kFALSE;
163   //
164   // fast count points
165   Int_t countFit = 0;  
166   for (Int_t ifit=0; ifit<volIdsFit->GetSize(); ifit++){
167     Int_t volIdFit  = volIdsFit->At(ifit);
168     for (Int_t ipoint = 0; ipoint < npoints; ipoint++)
169       if (volIdFit==fPoints->GetVolumeID()[ipoint]) countFit++;    
170   }
171   if (countFit<fMinNPoints) return kFALSE;
172   //
173   Int_t countPoint = 0;  
174   for (Int_t jpoint=0; jpoint<volIds->GetSize(); jpoint++){
175     Int_t volIdPoint  = volIds->At(jpoint);
176     for (Int_t ipoint = 0; ipoint < npoints; ipoint++)
177       if (volIdPoint==fPoints->GetVolumeID()[ipoint]) countPoint++;    
178   }
179   if (countPoint<kMinPoints) return kFALSE;
180   //
181   //
182
183   Reset();
184
185   if (fPoints && AliLog::GetDebugLevel("","AliTrackFitterRieman")>1&& gRandom->Rndm()<debugRatio){
186     Int_t nVol    = volIds->GetSize();
187     Int_t nVolFit = volIdsFit->GetSize();
188     Int_t volId  = volIds->At(0);    
189     (*fDebugStream)<<"PInput"<<
190       "NPoints="<<npoints<<    // number of points
191       "VolId="<<volId<<        // first vol ID
192       "NVol="<<nVol<<          // number of volumes 
193       "NvolFit="<<nVolFit<<    // number of volumes to fit
194       "fPoints.="<<fPoints<<   // input points
195       "\n";
196   }
197
198   Bool_t isAlphaCalc = kFALSE;
199   AliTrackPoint p,plocal;
200 //   fPoints->GetPoint(p,0);
201 //   fAlpha = TMath::ATan2(p.GetY(),p.GetX());
202
203   Int_t npVolId = 0;
204   fNUsed = 0;
205   Int_t *pindex = new Int_t[npoints];
206   for (Int_t ipoint = 0; ipoint < npoints; ipoint++)
207     {
208       fPoints->GetPoint(p,ipoint);
209       UShort_t iVolId = p.GetVolumeID();
210       if (FindVolId(volIds,iVolId)) {
211         pindex[npVolId] = ipoint;
212         npVolId++;
213       }
214       if (volIdsFit != 0x0) {
215         if (!FindVolId(volIdsFit,iVolId)) continue;
216       }
217       else {
218         if (iVolId < AliAlignObj::LayerToVolUID(layerRangeMin,0) ||
219             iVolId > AliAlignObj::LayerToVolUID(layerRangeMax,
220                                                 AliAlignObj::LayerSize(layerRangeMax))) continue;
221       }
222       if (!isAlphaCalc) {
223         fAlpha = p.GetAngle();
224         isAlphaCalc = kTRUE;
225       }
226       plocal = p.Rotate(fAlpha);
227       if (TMath::Abs(plocal.GetX())>500 || TMath::Abs(plocal.GetX())<2 || plocal.GetCov()[3]<=0 ||plocal.GetCov()[5]<=0 ){
228         printf("<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<</n");
229         p.Dump();
230         plocal.Dump();
231         printf("Problematic point\n");  
232         printf("<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<</n");
233       }else{
234         AddPoint(plocal.GetX(),plocal.GetY(),plocal.GetZ(),
235                  TMath::Sqrt(plocal.GetCov()[3]),TMath::Sqrt(plocal.GetCov()[5]));
236       }
237       // fNUsed++;  AddPoint should be responsible
238     }
239
240   if (npVolId == 0 || fNUsed < fMinNPoints) {
241     delete [] pindex;
242     return kFALSE;
243   }
244
245   Update();
246  
247   if (!fConv) {
248     delete [] pindex;
249     return kFALSE;
250   }
251
252   if ((fParams[0] == 0) ||
253       ((-fParams[2]*fParams[0]+fParams[1]*fParams[1]+1) <= 0)) {
254     delete [] pindex;
255     return kFALSE;
256   }
257
258
259   if (fNUsed < fMinNPoints) {
260     delete [] pindex;
261     return kFALSE;
262   }
263
264   fPVolId = new AliTrackPointArray(npVolId);
265   fPTrack = new AliTrackPointArray(npVolId);
266   AliTrackPoint p2;
267   for (Int_t ipoint = 0; ipoint < npVolId; ipoint++)
268     {
269       Int_t index = pindex[ipoint];
270       fPoints->GetPoint(p,index);
271       if (GetPCA(p,p2) && (
272           TMath::Abs(p.GetX()-p2.GetX())<fMaxDelta && 
273           TMath::Abs(p.GetY()-p2.GetY())<fMaxDelta &&
274           TMath::Abs(p.GetZ()-p2.GetZ())<fMaxDelta
275           )) {
276         Float_t xyz[3],xyz2[3];
277         p.GetXYZ(xyz); p2.GetXYZ(xyz2);
278         //      printf("residuals %f %d %d %f %f %f %f %f %f\n",fChi2,fNUsed,fConv,xyz[0],xyz[1],xyz[2],xyz2[0]-xyz[0],xyz2[1]-xyz[1],xyz2[2]-xyz[2]);  
279         fPVolId->AddPoint(ipoint,&p);
280         fPTrack->AddPoint(ipoint,&p2);
281       }else{
282         // what should be default bahavior -  
283         delete [] pindex;
284         delete fPVolId;
285         delete fPTrack;
286         fPVolId =0;
287         fPTrack =0;
288         return kFALSE;
289       }
290     }  
291   
292   if (AliLog::GetDebugLevel("","AliTrackFitterRieman")>0 && gRandom->Rndm()<debugRatio){
293     //
294     // Debug Info
295     //
296     AliTrackPointArray *lPVolId = new AliTrackPointArray(npVolId);
297     AliTrackPointArray *lPTrack = new AliTrackPointArray(npVolId);
298     AliTrackPointArray *lPTrackE = new AliTrackPointArray(npVolId);
299     AliRieman * residual = fRieman->MakeResiduals();
300     for (Int_t ipoint = 0; ipoint < npVolId; ipoint++){
301       AliTrackPoint p0, p0local;      
302       AliTrackPoint pFit, pFitlocal, pFitLocalE;      
303       fPVolId->GetPoint(p0,ipoint);
304       Float_t lAngle = p0.GetAngle();
305       p0local= p0.MasterToLocal();
306       fPTrack->GetPoint(pFit,ipoint);
307       pFitlocal= pFit.Rotate(lAngle);
308       //
309       Float_t xyz[3], cov[3];
310       xyz[0] = pFitlocal.GetX();
311       xyz[1] = pFitlocal.GetY();
312       xyz[2] = pFitlocal.GetZ();
313       for (Int_t icov=0; icov<6; icov++) cov[icov]=0;
314       cov[3] = GetErrY2at(xyz[0]);
315       cov[5] = GetErrZ2at(xyz[0]);      
316       pFitLocalE.SetXYZ(xyz,cov);
317       //
318       lPVolId->AddPoint(ipoint,&p0local);
319       lPTrack->AddPoint(ipoint,&pFitlocal);
320       lPTrackE->AddPoint(ipoint,&pFitLocalE);
321     }      
322     //
323     // debug info
324     //
325     Int_t nVol    = volIds->GetSize();
326     Int_t nVolFit = volIdsFit->GetSize();
327     Int_t volId   = volIds->At(0);   
328     Int_t modId   =0;
329     Int_t layer   =  AliAlignObj::VolUIDToLayer(volId,modId);
330     Int_t volIdFit   = volIdsFit->At(0);   
331     Int_t modIdFit   =0;
332     Int_t layerFit   =  AliAlignObj::VolUIDToLayer(volIdFit,modIdFit);
333     
334     (*fDebugStream)<<"Fit"<<
335       "VolId="<<volId<<        // volume ID
336       "Layer="<<layer<<        // layer ID
337       "Module="<<modId<<       // module ID
338       "LayerFit="<<layerFit<<        // layer ID fit
339       "ModuleFit="<<modIdFit<<       // module ID fit
340       "NVol="<<nVol<<          // number of volumes 
341       "NvolFit="<<nVolFit<<    // number of volumes to fit
342       "Points0.="<<fPVolId<<    // original points
343       "Points1.="<<fPTrack<<    // fitted points 
344       "LPoints0.="<<lPVolId<<   // original points - local frame
345       "LPoints1.="<<lPTrack<<   // fitted points   - local frame      
346       "LPointsE.="<<lPTrackE<<   // fitted points with ext error - local frame      
347       "Rieman.="<<this<<        // original rieman fit
348       "Res.="<<residual<<       // residuals of rieman fit  
349       "\n";
350     delete lPVolId;
351     delete lPTrack;
352     delete residual;
353   }
354
355   delete [] pindex;
356   return kTRUE;
357 }
358
359
360 void AliTrackFitterRieman::AddPoint(Float_t x, Float_t y, Float_t z, Float_t sy, Float_t sz)
361 {
362   //
363   // add point to rieman fitter
364   //
365   fRieman->AddPoint(x,y,z,sy,sz);
366   fNUsed = fRieman->GetN();
367 }
368
369
370
371 void AliTrackFitterRieman::Update(){
372   //
373   // 
374   //
375   fRieman->Update();
376   fConv = kFALSE;
377   if (fRieman->IsValid()){
378     for (Int_t ipar=0; ipar<6; ipar++){
379       fParams[ipar] = fRieman->GetParam()[ipar];
380     }
381     fChi2  =  fRieman->GetChi2();
382     fNdf   =  fRieman->GetN()- 2;
383     fNUsed =  fRieman->GetN();
384     fConv = kTRUE;
385   }
386   //
387   //
388   //
389   TLinearFitter fitY(3,"pol2");
390   TLinearFitter fitZ(3,"pol2");
391   for (Int_t ip=0; ip<fRieman->GetN();ip++){
392     Double_t x = fRieman->GetX()[ip]; 
393     fitY.AddPoint(&x,fRieman->GetY()[ip]-fRieman->GetYat(x),1);
394     fitZ.AddPoint(&x,fRieman->GetZ()[ip]-fRieman->GetZat(x),1);    
395   }
396   fitY.Eval();
397   fitZ.Eval();
398   for (Int_t iparam=0; iparam<3; iparam++){
399     fCorrY[iparam]=fitY.GetParameter(iparam);
400     fCorrZ[iparam]=fitZ.GetParameter(iparam);
401   }
402   fCorrY[3]=fitY.GetChisquare()/Float_t(fRieman->GetN()-3);   
403   fCorrZ[3]=fitZ.GetChisquare()/Float_t(fRieman->GetN()-3);
404 }
405
406
407
408 //_____________________________________________________________________________
409 Bool_t AliTrackFitterRieman::GetPCA(const AliTrackPoint &p, AliTrackPoint &p2) const
410 {
411   //
412   // Get the closest to a given spacepoint track trajectory point
413   // Look for details in the description of the Fit() method
414   //
415   if (!fConv) return kFALSE;
416
417   // First X and Y coordinates
418   Double_t sin = TMath::Sin(fAlpha);
419   Double_t cos = TMath::Cos(fAlpha);
420   //  fParam[0]  = 1/y0
421   //  fParam[1]  = -x0/y0
422   //  fParam[2]  = - (R^2 - x0^2 - y0^2)/y0
423   if (fParams[0] == 0) return kFALSE;
424   // Track parameters in the global coordinate system
425   Double_t x0 = -fParams[1]/fParams[0]*cos -         1./fParams[0]*sin;
426   Double_t y0 =          1./fParams[0]*cos - fParams[1]/fParams[0]*sin;
427   if ((-fParams[2]*fParams[0]+fParams[1]*fParams[1]+1) <= 0) return kFALSE;
428   Double_t r  = TMath::Sqrt(-fParams[2]*fParams[0]+fParams[1]*fParams[1]+1)/
429                 fParams[0];
430
431   // Define space-point refence plane
432   Double_t alphap = p.GetAngle();
433   Double_t sinp = TMath::Sin(alphap);
434   Double_t cosp = TMath::Cos(alphap);
435   Double_t x  = p.GetX()*cosp + p.GetY()*sinp;
436   Double_t y  = p.GetY()*cosp - p.GetX()*sinp;
437   Double_t x0p= x0*cosp + y0*sinp;
438   Double_t y0p= y0*cosp - x0*sinp;
439   if ((r*r - (x-x0p)*(x-x0p))<0) {
440     AliWarning(Form("Track extrapolation failed ! (Track radius = %f, track circle x = %f, space-point x = %f, reference plane angle = %f\n",r,x0p,x,alphap));
441     return kFALSE;
442   }
443   Double_t temp = TMath::Sqrt(r*r - (x-x0p)*(x-x0p));
444   Double_t y1 = y0p + temp;
445   Double_t y2 = y0p - temp;
446   Double_t yprime = y1;
447   if(TMath::Abs(y2-y) < TMath::Abs(y1-y)) yprime = y2;
448
449   // Back to the global coordinate system
450   Double_t xsecond = x*cosp - yprime*sinp;
451   Double_t ysecond = yprime*cosp + x*sinp;
452
453   // Now Z coordinate and track angles
454   Double_t x2 = xsecond*cos + ysecond*sin;
455   Double_t zsecond = GetZat(x2);
456   Double_t dydx = GetDYat(x2);
457   Double_t dzdx = GetDZat(x2);
458
459   // Fill the cov matrix of the track extrapolation point
460   Double_t cov[6] = {0,0,0,0,0,0};
461   Double_t sigmax = 100*100.;
462   cov[0] = sigmax;           cov[1] = sigmax*dydx;      cov[2] = sigmax*dzdx;
463   cov[3] = sigmax*dydx*dydx; cov[4] = sigmax*dydx*dzdx;
464   cov[5] = sigmax*dzdx*dzdx;
465
466   Double_t sigmay2 = GetErrY2at(x2);
467   Double_t sigmaz2 = GetErrZ2at(x2);
468   cov[3] += sigmay2;
469   cov[5] += sigmaz2;
470    
471
472   Float_t  newcov[6];
473   newcov[0] = cov[0]*cos*cos-
474             2*cov[1]*sin*cos+
475               cov[3]*sin*sin;
476   newcov[1] = cov[1]*(cos*cos-sin*sin)-
477              (cov[3]-cov[0])*sin*cos;
478   newcov[2] = cov[2]*cos-
479               cov[4]*sin;
480   newcov[3] = cov[0]*sin*sin+
481             2*cov[1]*sin*cos+
482               cov[3]*cos*cos;
483   newcov[4] = cov[4]*cos+
484               cov[2]*sin;
485   newcov[5] = cov[5];
486
487   p2.SetXYZ(xsecond,ysecond,zsecond,newcov);
488   Int_t debugLevel = AliLog::GetDebugLevel("","AliTrackFitterRieman");
489   Float_t debugRatio = 1./(1.+debugLevel);
490   if (AliLog::GetDebugLevel("","AliTrackFitterRieman")>0 && gRandom->Rndm()<debugRatio){
491     AliTrackPoint lp0(p);
492     AliTrackPoint lp2(p2);
493     (*fDebugStream)<<"PCA"<<
494       "P0.="<<&lp0<<
495       "P2.="<<&lp2<<
496       "\n";
497   }
498   return kTRUE;
499 }
500
501 Double_t AliTrackFitterRieman::GetYat(Double_t x) const  {
502   //
503   // get y position at given point 
504   //
505   Double_t correction=0;
506   if (fBCorrection){ // systematic effect correction
507     correction  =  fCorrY[0]+fCorrY[1]*x +fCorrY[2]*x*x;
508   }
509   return fRieman->GetYat(x)+correction;
510 }
511
512 Double_t AliTrackFitterRieman::GetZat(Double_t x) const  {
513   //
514   // get z position at given point 
515   //
516   Double_t correction=0;
517   if (fBCorrection){ // systematic effect correction
518     correction  =  fCorrZ[0]+fCorrZ[1]*x +fCorrZ[2]*x*x;
519   }
520   return fRieman->GetZat(x)+correction;
521 }
522
523 Double_t AliTrackFitterRieman::GetErrY2at(Double_t x) const {
524   //
525   // get estimate of extrapolation error
526   //
527   Double_t error = fRieman->GetErrY(x);
528   Double_t correction=0;
529   if (fBCorrection){ // everestimate error due systematic effect
530     error      *=fCorrY[3];
531     correction  =  fCorrY[0]+fCorrY[1]*x +fCorrY[2]*x*x;
532     correction *=correction;
533   }
534   return TMath::Sqrt(error+correction);
535 }
536
537 Double_t AliTrackFitterRieman::GetErrZ2at(Double_t x) const {
538   //
539   // get estimate of extrapolation error
540   //
541   Double_t error = fRieman->GetErrZ(x)*fCorrZ[3];
542   Double_t correction=0;
543   if (fBCorrection){  
544      error    *=  fCorrZ[3];
545     correction =  fCorrZ[0]+fCorrZ[1]*x +fCorrZ[2]*x*x;
546     correction*=  correction;
547   }
548   return TMath::Sqrt(error+correction);
549 }
550