]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibCalib.cxx
Fixes for Coverity defects
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibCalib.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
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //   Component for redoing the reconstruction from the clusters and tracks
20 //      
21 //   The new calibration data used 
22 //
23 //   In reality it overwrites the content of the ESD 
24 //    
25
26 /*
27
28   gSystem->Load("libANALYSIS");
29   gSystem->Load("libTPCcalib");
30   //
31   gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
32   gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
33   AliXRDPROOFtoolkit tool; 
34   TChain * chainCl = tool.MakeChain("calib.txt","Clusters",0,1);
35   chainCl->Lookup();
36   TChain * chainTr = tool.MakeChain("calib.txt","Tracks",0,1);
37   chainTr->Lookup();
38
39
40
41 */
42
43
44
45 //  marian.ivanov@cern.ch
46 // 
47 #include "AliTPCcalibCalib.h"
48 #include "TSystem.h"
49 #include "TFile.h"
50 #include "TTreeStream.h"
51 #include "AliLog.h"
52 #include "TTimeStamp.h"
53 #include "AliESDEvent.h"
54 #include "AliESDfriend.h"
55 #include "AliESDtrack.h"
56 #include "AliTracker.h"
57 #include "AliTPCClusterParam.h"
58 #include "AliTPCParam.h"
59
60 #include "AliTPCcalibDB.h"
61 #include "AliTPCTransform.h"
62 #include "AliTPCRecoParam.h"
63 #include "AliTPCclusterMI.h"
64 #include "AliTPCseed.h"
65 #include "AliTPCPointCorrection.h"
66 #include <TGeoManager.h>
67 #include <TGeoPhysicalNode.h>
68 #include "TDatabasePDG.h"
69 ClassImp(AliTPCcalibCalib)
70
71 AliTPCcalibCalib::AliTPCcalibCalib():
72 AliTPCcalibBase(),
73   fApplyExBCorrection(1),      // apply ExB correction
74   fApplyTOFCorrection(1),      // apply TOF correction
75   fApplyPositionCorrection(1), // apply position correction
76   fApplySectorAlignment(1),    // apply sector alignment
77   fApplyRPhiCorrection(1),     // apply R-Phi correction
78   fApplyRCorrection(1)         // apply Radial correction
79
80 {
81   //
82   // Constructor
83   //
84 }
85
86
87 AliTPCcalibCalib::AliTPCcalibCalib(const Text_t *name, const Text_t *title) 
88   :AliTPCcalibBase(),
89   fApplyExBCorrection(1),      // apply ExB correction
90   fApplyTOFCorrection(1),      // apply TOF correction
91   fApplyPositionCorrection(1), // apply position correction
92   fApplySectorAlignment(1),    // apply sector alignment
93   fApplyRPhiCorrection(1),     // apply R-Phi correction
94   fApplyRCorrection(1)         // apply Radial correction
95 {  
96   SetName(name);
97   SetTitle(title);
98 }
99
100
101 AliTPCcalibCalib::AliTPCcalibCalib(const AliTPCcalibCalib&calib):
102   AliTPCcalibBase(calib),
103   fApplyExBCorrection(calib.GetApplyExBCorrection()),
104   fApplyTOFCorrection(calib.GetApplyTOFCorrection()),
105   fApplyPositionCorrection(calib.GetApplyPositionCorrection()),
106   fApplySectorAlignment(calib.GetApplySectorAlignment()),
107   fApplyRPhiCorrection(calib.GetApplyRPhiCorrection()),
108   fApplyRCorrection(calib.GetApplyRCorrection())
109
110 {
111   //
112   // copy constructor
113   //
114 }
115
116 AliTPCcalibCalib &AliTPCcalibCalib::operator=(const AliTPCcalibCalib&calib){
117   //
118   //
119   //
120   ((AliTPCcalibBase *)this)->operator=(calib);
121   return *this;
122 }
123
124
125 AliTPCcalibCalib::~AliTPCcalibCalib() {
126   //
127   // destructor
128   //
129 }
130
131
132 void     AliTPCcalibCalib::Process(AliESDEvent *event){
133   //
134   // 
135   //
136   if (!event) {
137     return;
138   }  
139   AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
140   if (!ESDfriend) {
141    return;
142   }
143   
144   if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
145   Int_t ntracks=event->GetNumberOfTracks();   
146   //AliTPCcalibDB::Instance()->SetExBField(fMagF);
147
148   //
149   //
150   //
151
152   for (Int_t i=0;i<ntracks;++i) {
153     AliESDtrack *track = event->GetTrack(i);     
154     AliESDfriendTrack *friendTrack = (AliESDfriendTrack*) ESDfriend->GetTrack(i);
155     if (!friendTrack) continue;
156     //track->SetFriendTrack(friendTrack);
157     fCurrentFriendTrack=friendTrack;
158     const AliExternalTrackParam * trackIn  = track->GetInnerParam();
159     const AliExternalTrackParam * trackOut = track->GetOuterParam();
160     AliExternalTrackParam * tpcOut   = (AliExternalTrackParam *)friendTrack->GetTPCOut();
161     if (!trackIn) continue;
162     if (!trackOut) continue;
163     if (!tpcOut) continue;   
164     TObject *calibObject;
165     AliTPCseed *seed = 0;
166     for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
167       if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
168     }
169     if (!seed) continue;
170     RefitTrack(track, seed, event->GetMagneticField());
171     (*tpcOut)=*(track->GetOuterParam());  
172   }
173   return;
174 }
175
176 Bool_t  AliTPCcalibCalib::RefitTrack(AliESDtrack * track, AliTPCseed *seed, Float_t magesd){
177   //
178   // Refit track
179   // if magesd==0 forget the curvature
180
181   //
182   // 0 - Setup transform object
183   //
184   static Int_t streamCounter=0;
185   streamCounter++;
186   AliESDfriendTrack *friendTrack = fCurrentFriendTrack;
187
188   AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
189   AliTPCParam     *param     = AliTPCcalibDB::Instance()->GetParameters();
190   transform->SetCurrentRun(fRun);
191   transform->SetCurrentTimeStamp((UInt_t)fTime);
192   if(!fApplyExBCorrection) { // disable ExB correction in transform
193     if(transform->GetCurrentRecoParam())
194       transform->GetCurrentRecoParamNonConst()->SetUseExBCorrection(0);
195   }
196   if(!fApplyTOFCorrection) { // disable TOF correction in transform
197     if(transform->GetCurrentRecoParam())
198       transform->GetCurrentRecoParamNonConst()->SetUseTOFCorrection(kFALSE);
199   }
200
201   //
202   // First apply calibration
203   //
204   //  AliTPCPointCorrection * corr =  AliTPCPointCorrection::Instance();
205   for (Int_t irow=0;irow<159;irow++) {
206     AliTPCclusterMI *cluster=seed->GetClusterPointer(irow);
207     if (!cluster) continue; 
208     AliTPCclusterMI cl0(*cluster);
209     Double_t x[3]={cluster->GetRow(),cluster->GetPad(),cluster->GetTimeBin()};
210     Int_t i[1]={cluster->GetDetector()};
211
212     transform->Transform(x,i,0,1);
213     //
214     // get position correction
215     //
216     Int_t ipad=0;
217     if (cluster->GetDetector()>35) ipad=1;
218     Float_t dy =0;//AliTPCClusterParam::SPosCorrection(0,ipad,cluster->GetPad(),cluster->GetTimeBin(),cluster->GetZ(),cluster->GetSigmaY2(),cluster->GetSigmaZ2(),cluster->GetMax());
219     Float_t dz =0;//AliTPCClusterParam::SPosCorrection(1,ipad,cluster->GetPad(),cluster->GetTimeBin(),cluster->GetZ(),cluster->GetSigmaY2(),cluster->GetSigmaZ2(),cluster->GetMax());
220     //
221     cluster->SetX(x[0]);
222     cluster->SetY(x[1]);
223     cluster->SetZ(x[2]);
224     
225     //
226     // Apply alignemnt
227     //
228     if (transform->GetCurrentRecoParam()->GetUseSectorAlignment()){
229       if (!param->IsGeoRead()) param->ReadGeoMatrices();
230       TGeoHMatrix  *mat = param->GetClusterMatrix(cluster->GetDetector());
231       //TGeoHMatrix  mat;
232       Double_t pos[3]= {cluster->GetX(),cluster->GetY(),cluster->GetZ()};
233       Double_t posC[3]={cluster->GetX(),cluster->GetY(),cluster->GetZ()};
234       if (mat) mat->LocalToMaster(pos,posC);
235       else{
236         // chack Loading of Geo matrices from GeoManager - TEMPORARY FIX
237       }
238       cluster->SetX(posC[0]);
239       cluster->SetY(posC[1]);
240       cluster->SetZ(posC[2]);
241     }
242
243
244
245     if (fStreamLevel>2 && streamCounter<20*fStreamLevel ){
246       // dump debug info if required
247       TTreeSRedirector *cstream = GetDebugStreamer();
248       if (cstream){
249         (*cstream)<<"Clusters"<<
250           "run="<<fRun<<              //  run number
251           "event="<<fEvent<<          //  event number
252           "time="<<fTime<<            //  time stamp of event
253           "trigger="<<fTrigger<<      //  trigger
254           "triggerClass="<<&fTriggerClass<<      //  trigger      
255           "mag="<<fMagF<<             //  magnetic field
256           "cl0.="<<&cl0<<
257           "cl.="<<cluster<<
258           "cy="<<dy<<
259           "cz="<<dz<<
260           "\n";
261       }
262     }
263   }
264   //
265   //
266   //
267   Int_t ncl = seed->GetNumberOfClusters();
268   const Double_t kResetCov=4.;
269   const Double_t kSigma=5.;
270   Double_t covar[15];
271   for (Int_t i=0;i<15;i++) covar[i]=0;
272   covar[0]=kSigma*kSigma;
273   covar[2]=kSigma*kSigma;
274   covar[5]=kSigma*kSigma/Float_t(ncl*ncl);
275   covar[9]=kSigma*kSigma/Float_t(ncl*ncl);
276   covar[14]=0.2*0.2;
277   if (TMath::Abs(magesd)<0.05) {
278      covar[14]=0.025*0.025;
279   }
280   // 
281   // And now do refit
282   //
283   AliExternalTrackParam * trackInOld  = (AliExternalTrackParam*)track->GetInnerParam();
284   AliExternalTrackParam * trackOuter = (AliExternalTrackParam*)track->GetOuterParam();
285   AliExternalTrackParam * trackOutOld   = (AliExternalTrackParam *)friendTrack->GetTPCOut();
286   Double_t mass =    TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
287
288     
289
290   AliExternalTrackParam trackIn  = *trackOutOld;
291   trackIn.ResetCovariance(kResetCov);
292   trackIn.AddCovariance(covar);
293   if (TMath::Abs(magesd)<0.05) {
294     ((Double_t&)(trackIn.GetParameter()[4]))=0.000000001;
295     ((Double_t&)(trackIn.GetCovariance()[14]))=covar[14];  // fix the line
296   }  
297
298   Double_t xyz[3];
299   Int_t nclIn=0,nclOut=0;
300   //
301   // Refit in
302   //
303   for (Int_t irow=159; irow>0; irow--){
304     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
305     if (!cl) continue;
306     if (cl->GetX()<80) continue;
307     Int_t sector = cl->GetDetector();
308     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
309     if (TMath::Abs(dalpha)>0.01){
310       if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
311     }
312     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
313     Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation
314     AliTPCseed::GetError(cl, &trackIn,cov[0],cov[2]);
315     cov[0]*=cov[0];
316     cov[2]*=cov[2];
317     trackIn.GetXYZ(xyz);
318     //    Double_t bz = AliTracker::GetBz(xyz);
319
320     //    if (!trackIn.PropagateTo(r[0],bz)) continue;
321     if (!AliTracker::PropagateTrackToBxByBz(&trackIn, r[0],mass,1.,kFALSE)) continue;
322
323     if (RejectCluster(cl,&trackIn)) continue;
324     nclIn++;
325     trackIn.Update(&r[1],cov);    
326   }
327   //
328   AliExternalTrackParam trackOut = trackIn;
329   trackOut.ResetCovariance(kResetCov);
330   trackOut.AddCovariance(covar);
331   if (TMath::Abs(magesd)<0.05) {
332     ((Double_t&)(trackOut.GetParameter()[4]))=0.000000001;
333     ((Double_t&)(trackOut.GetCovariance()[14]))=covar[14];  // fix the line
334   }
335
336   //
337   // Refit out
338   //
339   //Bool_t lastEdge=kFALSE;
340   for (Int_t irow=0; irow<160; irow++){
341     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
342     if (!cl) continue;
343     if (cl->GetX()<80) continue;
344     Int_t sector = cl->GetDetector();
345     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
346
347     if (TMath::Abs(dalpha)>0.01){
348       if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
349     }
350     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
351
352     Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation    
353     AliTPCseed::GetError(cl, &trackOut,cov[0],cov[2]);
354     cov[0]*=cov[0];
355     cov[2]*=cov[2];
356     trackOut.GetXYZ(xyz);
357     //Double_t bz = AliTracker::GetBz(xyz);
358     //    if (!trackOut.PropagateTo(r[0],bz)) continue;
359     if (!AliTracker::PropagateTrackToBxByBz(&trackOut, r[0],mass,1.,kFALSE)) continue;
360
361     if (RejectCluster(cl,&trackOut)) continue;
362     nclOut++;
363     trackOut.Update(&r[1],cov);         
364     //if (cl->GetType()<0) lastEdge=kTRUE;
365     //if (cl->GetType()>=0) lastEdge=kFALSE;    
366   }
367   //
368   //
369   //
370   nclIn=0;
371   trackIn  = trackOut;
372   trackIn.ResetCovariance(kResetCov);
373   if (TMath::Abs(magesd)<0.05) {
374     ((Double_t&)(trackIn.GetParameter()[4]))=0.000000001;
375     ((Double_t&)(trackIn.GetCovariance()[14]))=covar[14];  // fix the line
376   }
377   //
378   // Refit in one more time
379   //
380   for (Int_t irow=159; irow>0; irow--){
381     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
382     if (!cl) continue;
383     if (cl->GetX()<80) continue;
384     Int_t sector = cl->GetDetector();
385     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
386     if (TMath::Abs(dalpha)>0.01){
387       if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
388     }
389     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
390     Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation
391     AliTPCseed::GetError(cl, &trackIn,cov[0],cov[2]);
392     cov[0]*=cov[0];
393     cov[2]*=cov[2];
394     trackIn.GetXYZ(xyz);
395     //Double_t bz = AliTracker::GetBz(xyz);
396
397     //    if (!trackIn.PropagateTo(r[0],bz)) continue;
398     if (!AliTracker::PropagateTrackToBxByBz(&trackIn, r[0],mass,1,kFALSE)) continue;
399
400     if (RejectCluster(cl,&trackIn)) continue;
401     nclIn++;
402     trackIn.Update(&r[1],cov);    
403   }
404
405
406   trackIn.Rotate(trackInOld->GetAlpha());
407   trackOut.Rotate(trackOutOld->GetAlpha());
408   //
409   trackInOld->GetXYZ(xyz);
410   Double_t bz = AliTracker::GetBz(xyz);
411   trackIn.PropagateTo(trackInOld->GetX(),bz);
412   //
413   trackOutOld->GetXYZ(xyz);
414   bz = AliTracker::GetBz(xyz);  
415   trackOut.PropagateTo(trackOutOld->GetX(),bz);
416   
417
418   if (fStreamLevel>0 && streamCounter<100*fStreamLevel){
419     TTreeSRedirector *cstream = GetDebugStreamer();
420     if (cstream){
421       (*cstream)<<"Tracks"<<
422         "run="<<fRun<<              //  run number
423         "event="<<fEvent<<          //  event number
424         "time="<<fTime<<            //  time stamp of event
425         "trigger="<<fTrigger<<      //  trigger
426         "triggerClass="<<&fTriggerClass<<      //  trigger
427         "mag="<<fMagF<<             //  magnetic field
428         "nclIn="<<nclIn<<
429         "nclOut="<<nclOut<<
430         "ncl="<<ncl<<
431         "TrIn0.="<<trackInOld<<
432         "TrOut0.="<<trackOutOld<<
433         "TrIn1.="<<&trackIn<<
434         "TrOut1.="<<&trackOut<<
435         "\n";
436     }
437   }
438   //
439   // And now rewrite ESDtrack and TPC seed
440   //
441  
442   (*trackInOld)  = trackIn;
443   (*trackOutOld) = trackOut;
444   (*trackOuter) = trackOut;
445   AliExternalTrackParam *t = &trackIn;
446   //track->Set(t->GetX(),t->GetAlpha(),t->GetParameter(),t->GetCovariance());
447   seed->Set(t->GetX(),t->GetAlpha(),t->GetParameter(),t->GetCovariance());
448   seed->SetNumberOfClusters((nclIn+nclOut)/2);
449   return kTRUE;
450 }
451
452
453
454 Bool_t AliTPCcalibCalib::RejectCluster(AliTPCclusterMI* cl, AliExternalTrackParam * param){
455   //
456   // check the acceptance of cluster
457   // Cut on edge effects
458   //
459   Float_t kEdgeCut=2.5;
460   Float_t kSigmaCut=6;
461
462   Bool_t isReject = kFALSE;
463   Float_t edgeY = cl->GetX()*TMath::Tan(TMath::Pi()/18);
464   Float_t dist  = edgeY - TMath::Abs(cl->GetY());
465   if (param)  dist  = TMath::Abs(edgeY - TMath::Abs(param->GetY()));
466   if (dist<kEdgeCut) isReject=kTRUE;
467
468   Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation    
469   AliTPCseed::GetError(cl, param,cov[0],cov[2]);
470   if (param->GetSigmaY2()<0 || param->GetSigmaZ2()<0){
471     AliError("Wrong parameters");
472     return kFALSE;
473   }
474   Double_t py = (cl->GetY()-param->GetY())/TMath::Sqrt(cov[0]*cov[0]+param->GetSigmaY2());
475   Double_t pz = (cl->GetZ()-param->GetZ())/TMath::Sqrt(cov[2]*cov[2]+param->GetSigmaZ2());
476   //
477   if ((py*py+pz*pz)>kSigmaCut*kSigmaCut) isReject=kTRUE;
478   
479   return isReject;
480 }
481
482
483