]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibCalib.cxx
Macros moved to the TPC/scripts/OCDBscan directory
[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
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 = ESDfriend->GetTrack(i);
155  
156     const AliExternalTrackParam * trackIn  = track->GetInnerParam();
157     const AliExternalTrackParam * trackOut = track->GetOuterParam();
158     AliExternalTrackParam * tpcOut   = (AliExternalTrackParam *)friendTrack->GetTPCOut();
159     if (!trackIn) continue;
160     if (!trackOut) continue;
161     if (!tpcOut) continue;   
162     TObject *calibObject;
163     AliTPCseed *seed = 0;
164     for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
165       if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
166     }
167     if (!seed) continue;
168     RefitTrack(track, seed,event->GetMagneticField());
169     (*tpcOut)=*(track->GetOuterParam());  
170   }
171   return;
172 }
173
174 Bool_t  AliTPCcalibCalib::RefitTrack(AliESDtrack * track, AliTPCseed *seed, Float_t magesd){
175   //
176   // Refit track
177   // if magesd==0 forget the curvature
178
179   //
180   // 0 - Setup transform object
181   //
182   AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
183   AliTPCParam     *param     = AliTPCcalibDB::Instance()->GetParameters();
184   transform->SetCurrentRun(fRun);
185   transform->SetCurrentTimeStamp((UInt_t)fTime);
186   if(!fApplyExBCorrection) { // disable ExB correction in transform
187     if(transform->GetCurrentRecoParam())
188       transform->GetCurrentRecoParamNonConst()->SetUseExBCorrection(0);
189   }
190   if(!fApplyTOFCorrection) { // disable TOF correction in transform
191     if(transform->GetCurrentRecoParam())
192       transform->GetCurrentRecoParamNonConst()->SetUseTOFCorrection(kFALSE);
193   }
194
195   //
196   // First apply calibration
197   //
198   //  AliTPCPointCorrection * corr =  AliTPCPointCorrection::Instance();
199   for (Int_t irow=0;irow<159;irow++) {
200     AliTPCclusterMI *cluster=seed->GetClusterPointer(irow);
201     if (!cluster) continue; 
202     AliTPCclusterMI cl0(*cluster);
203     Double_t x[3]={cluster->GetRow(),cluster->GetPad(),cluster->GetTimeBin()};
204     Int_t i[1]={cluster->GetDetector()};
205
206     transform->Transform(x,i,0,1);
207     //
208     // get position correction
209     //
210     Int_t ipad=0;
211     if (cluster->GetDetector()>35) ipad=1;
212     Float_t dy =0;//AliTPCClusterParam::SPosCorrection(0,ipad,cluster->GetPad(),cluster->GetTimeBin(),cluster->GetZ(),cluster->GetSigmaY2(),cluster->GetSigmaZ2(),cluster->GetMax());
213     Float_t dz =0;//AliTPCClusterParam::SPosCorrection(1,ipad,cluster->GetPad(),cluster->GetTimeBin(),cluster->GetZ(),cluster->GetSigmaY2(),cluster->GetSigmaZ2(),cluster->GetMax());
214     // if(fApplyPositionCorrection) {
215     //x[1]-=dy;
216     //x[2]-=dz;
217     // }
218     //
219     // Apply sector alignment
220     //
221     Double_t dxq = AliTPCPointCorrection::SGetCorrectionSector(0,cluster->GetDetector()%36,cluster->GetX(),
222                                                                cluster->GetY(),cluster->GetZ());
223     Double_t dyq = AliTPCPointCorrection::SGetCorrectionSector(1,cluster->GetDetector()%36,cluster->GetX(),
224                                                                cluster->GetY(),cluster->GetZ());
225     Double_t dzq = AliTPCPointCorrection::SGetCorrectionSector(2,cluster->GetDetector()%36,cluster->GetX(),
226                                                                cluster->GetY(),cluster->GetZ());
227     if (0&fApplySectorAlignment){
228       x[0]-=dxq;
229       x[1]-=dyq;
230       x[2]-=dzq;
231     }
232 //     //
233 //     // Apply r-phi correction  - To be done on track level- knowing the track angle !!!
234 //     //
235 //     Double_t corrclY =  
236 //       corr->RPhiCOGCorrection(cluster->GetDetector(),cluster->GetRow(), cluster->GetPad(),
237 //                                cluster->GetY(),cluster->GetY(), cluster->GetZ(), 0., cluster->GetMax(),2.5);
238 //     // R correction
239 //     Double_t corrR   = corr->CorrectionOutR0(kFALSE,kFALSE,cluster->GetX(),cluster->GetY(),cluster->GetZ(),cluster->GetDetector());
240
241 //     if (0&fApplyRPhiCorrection){
242 //       if (cluster->GetY()>0) x[1]+=corrclY;  // rphi correction
243 //       if (cluster->GetY()<0) x[1]-=corrclY;  // rphi correction
244 //     }
245 //     if (0&fApplyRCorrection){      
246 //       x[0]+=corrR;                           // radial correction
247 //     }
248
249     //
250     //
251     //
252     cluster->SetX(x[0]);
253     cluster->SetY(x[1]);
254     cluster->SetZ(x[2]);
255     
256     //
257     // Apply alignemnt
258     //
259     if (1){
260       if (!param->IsGeoRead()) param->ReadGeoMatrices();
261       TGeoHMatrix  *mat = param->GetClusterMatrix(cluster->GetDetector());
262       //TGeoHMatrix  mat;
263       Double_t pos[3]= {cluster->GetX(),cluster->GetY(),cluster->GetZ()};
264       Double_t posC[3]={cluster->GetX(),cluster->GetY(),cluster->GetZ()};
265       if (mat) mat->LocalToMaster(pos,posC);
266       else{
267         // chack Loading of Geo matrices from GeoManager - TEMPORARY FIX
268       }
269       cluster->SetX(posC[0]);
270       cluster->SetY(posC[1]);
271       cluster->SetZ(posC[2]);
272     }
273
274
275
276     if (fStreamLevel>2){
277       TTreeSRedirector *cstream = GetDebugStreamer();
278       if (cstream){
279         (*cstream)<<"Clusters"<<
280           "run="<<fRun<<              //  run number
281           "event="<<fEvent<<          //  event number
282           "time="<<fTime<<            //  time stamp of event
283           "trigger="<<fTrigger<<      //  trigger
284           "triggerClass="<<&fTriggerClass<<      //  trigger      
285           "mag="<<fMagF<<             //  magnetic field
286           "cl0.="<<&cl0<<
287           "cl.="<<cluster<<
288           "cy="<<dy<<
289           "cz="<<dz<<
290           //      "cY="<<corrclY<<
291           //      "cR="<<corrR<<
292           "dxq="<<dxq<<
293           "dyq="<<dyq<<
294           "dzq="<<dzq<<
295           "\n";
296       }
297     }
298   }
299   //
300   //
301   //
302   Int_t ncl = seed->GetNumberOfClusters();
303   const Double_t kResetCov=4.;
304   const Double_t kSigma=5.;
305   Double_t covar[15];
306   for (Int_t i=0;i<15;i++) covar[i]=0;
307   covar[0]=kSigma*kSigma;
308   covar[2]=kSigma*kSigma;
309   covar[5]=kSigma*kSigma/Float_t(ncl*ncl);
310   covar[9]=kSigma*kSigma/Float_t(ncl*ncl);
311   covar[14]=0.2*0.2;
312   if (TMath::Abs(magesd)<0.05) {
313      covar[14]=0.025*0.025;
314   }
315   // 
316   // And now do refit
317   //
318   AliExternalTrackParam * trackInOld  = (AliExternalTrackParam*)track->GetInnerParam();
319   AliExternalTrackParam * trackOutOld = (AliExternalTrackParam*)track->GetOuterParam();
320
321
322   AliExternalTrackParam trackIn  = *trackOutOld;
323   trackIn.ResetCovariance(kResetCov);
324   trackIn.AddCovariance(covar);
325   if (TMath::Abs(magesd)<0.05) {
326     ((Double_t&)(trackIn.GetParameter()[4]))=0.000000001;
327     ((Double_t&)(trackIn.GetCovariance()[14]))=covar[14];  // fix the line
328   }  
329
330   Double_t xyz[3];
331   Int_t nclIn=0,nclOut=0;
332   //
333   // Refit in
334   //
335
336   for (Int_t irow=159; irow>0; irow--){
337     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
338     if (!cl) continue;
339     if (cl->GetX()<80) continue;
340     Int_t sector = cl->GetDetector();
341     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
342     if (TMath::Abs(dalpha)>0.01){
343       if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
344     }
345     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
346     Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation
347     AliTPCseed::GetError(cl, &trackIn,cov[0],cov[2]);
348     cov[0]*=cov[0];
349     cov[2]*=cov[2];
350     trackIn.GetXYZ(xyz);
351     Double_t bz = AliTracker::GetBz(xyz);
352
353     if (!trackIn.PropagateTo(r[0],bz)) continue;
354     if (RejectCluster(cl,&trackIn)) continue;
355     nclIn++;
356     trackIn.Update(&r[1],cov);    
357   }
358   //
359   AliExternalTrackParam trackOut = trackIn;
360   trackOut.ResetCovariance(kResetCov);
361   trackOut.AddCovariance(covar);
362   if (TMath::Abs(magesd)<0.05) {
363     ((Double_t&)(trackOut.GetParameter()[4]))=0.000000001;
364     ((Double_t&)(trackOut.GetCovariance()[14]))=covar[14];  // fix the line
365   }
366
367   //
368   // Refit out
369   //
370   //Bool_t lastEdge=kFALSE;
371   for (Int_t irow=0; irow<160; irow++){
372     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
373     if (!cl) continue;
374     if (cl->GetX()<80) continue;
375     Int_t sector = cl->GetDetector();
376     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
377
378     if (TMath::Abs(dalpha)>0.01){
379       if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
380     }
381     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
382
383     Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation    
384     AliTPCseed::GetError(cl, &trackOut,cov[0],cov[2]);
385     cov[0]*=cov[0];
386     cov[2]*=cov[2];
387     trackOut.GetXYZ(xyz);
388     Double_t bz = AliTracker::GetBz(xyz);
389     if (!trackOut.PropagateTo(r[0],bz)) continue;
390     if (RejectCluster(cl,&trackOut)) continue;
391     nclOut++;
392     trackOut.Update(&r[1],cov);         
393     //if (cl->GetType()<0) lastEdge=kTRUE;
394     //if (cl->GetType()>=0) lastEdge=kFALSE;    
395   }
396   //
397   //
398   //
399   nclIn=0;
400   trackIn  = trackOut;
401   trackIn.ResetCovariance(kResetCov);
402   if (TMath::Abs(magesd)<0.05) {
403     ((Double_t&)(trackIn.GetParameter()[4]))=0.000000001;
404     ((Double_t&)(trackIn.GetCovariance()[14]))=covar[14];  // fix the line
405   }
406   //
407   // Refit in one more time
408   //
409   for (Int_t irow=159; irow>0; irow--){
410     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
411     if (!cl) continue;
412     if (cl->GetX()<80) continue;
413     Int_t sector = cl->GetDetector();
414     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
415     if (TMath::Abs(dalpha)>0.01){
416       if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
417     }
418     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
419     Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation
420     AliTPCseed::GetError(cl, &trackIn,cov[0],cov[2]);
421     cov[0]*=cov[0];
422     cov[2]*=cov[2];
423     trackIn.GetXYZ(xyz);
424     Double_t bz = AliTracker::GetBz(xyz);
425
426     if (!trackIn.PropagateTo(r[0],bz)) continue;
427     if (RejectCluster(cl,&trackIn)) continue;
428     nclIn++;
429     trackIn.Update(&r[1],cov);    
430   }
431
432
433   trackIn.Rotate(trackInOld->GetAlpha());
434   trackOut.Rotate(trackOutOld->GetAlpha());
435   //
436   trackInOld->GetXYZ(xyz);
437   Double_t bz = AliTracker::GetBz(xyz);
438   trackIn.PropagateTo(trackInOld->GetX(),bz);
439   //
440   trackOutOld->GetXYZ(xyz);
441   bz = AliTracker::GetBz(xyz);  
442   trackOut.PropagateTo(trackOutOld->GetX(),bz);
443   
444
445   if (fStreamLevel>0){
446     TTreeSRedirector *cstream = GetDebugStreamer();
447     if (cstream){
448       (*cstream)<<"Tracks"<<
449         "run="<<fRun<<              //  run number
450         "event="<<fEvent<<          //  event number
451         "time="<<fTime<<            //  time stamp of event
452         "trigger="<<fTrigger<<      //  trigger
453         "triggerClass="<<&fTriggerClass<<      //  trigger
454         "mag="<<fMagF<<             //  magnetic field
455         "nclIn="<<nclIn<<
456         "nclOut="<<nclOut<<
457         "ncl="<<ncl<<
458         "TrIn0.="<<trackInOld<<
459         "TrOut0.="<<trackOutOld<<
460         "TrIn1.="<<&trackIn<<
461         "TrOut1.="<<&trackOut<<
462         "\n";
463     }
464   }
465   //
466   // And now rewrite ESDtrack and TPC seed
467   //
468  
469   (*trackInOld)  = trackIn;
470   (*trackOutOld) = trackOut;
471   AliExternalTrackParam *t = &trackIn;
472   track->Set(t->GetX(),t->GetAlpha(),t->GetParameter(),t->GetCovariance());
473   seed->Set(t->GetX(),t->GetAlpha(),t->GetParameter(),t->GetCovariance());
474   seed->SetNumberOfClusters((nclIn+nclOut)/2);
475   return kTRUE;
476 }
477
478
479
480 Bool_t AliTPCcalibCalib::RejectCluster(AliTPCclusterMI* cl, AliExternalTrackParam * param){
481   //
482   // check the acceptance of cluster
483   // Cut on edge effects
484   //
485   Float_t kEdgeCut=2.5;
486   Float_t kSigmaCut=6;
487
488   Bool_t isReject = kFALSE;
489   Float_t edgeY = cl->GetX()*TMath::Tan(TMath::Pi()/18);
490   Float_t dist  = edgeY - TMath::Abs(cl->GetY());
491   if (param)  dist  = TMath::Abs(edgeY - TMath::Abs(param->GetY()));
492   if (dist<kEdgeCut) isReject=kTRUE;
493
494   Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation    
495   AliTPCseed::GetError(cl, param,cov[0],cov[2]);
496   if (param->GetSigmaY2()<0 || param->GetSigmaZ2()<0){
497     AliError("Wrong parameters");
498     return kFALSE;
499   }
500   Double_t py = (cl->GetY()-param->GetY())/TMath::Sqrt(cov[0]*cov[0]+param->GetSigmaY2());
501   Double_t pz = (cl->GetZ()-param->GetZ())/TMath::Sqrt(cov[2]*cov[2]+param->GetSigmaZ2());
502   //
503   if ((py*py+pz*pz)>kSigmaCut*kSigmaCut) isReject=kTRUE;
504   
505   return isReject;
506 }
507
508
509