]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCcalibCalib.cxx
fixing warnings (Federico C.)
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibCalib.cxx
CommitLineData
9dcfce73 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
3bf5c7a6 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;
11c7a854 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();
3bf5c7a6 38
39
40
41*/
42
43
44
9dcfce73 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"
11c7a854 57#include "AliTPCClusterParam.h"
9dcfce73 58
59#include "AliTPCcalibDB.h"
60#include "AliTPCTransform.h"
7af539c6 61#include "AliTPCRecoParam.h"
9dcfce73 62#include "AliTPCclusterMI.h"
63#include "AliTPCseed.h"
15e48021 64#include "AliTPCPointCorrection.h"
9dcfce73 65
66ClassImp(AliTPCcalibCalib)
67
68AliTPCcalibCalib::AliTPCcalibCalib():
7af539c6 69AliTPCcalibBase(),
70 fApplyExBCorrection(1), // apply ExB correction
71 fApplyTOFCorrection(1), // apply TOF correction
72 fApplyPositionCorrection(1), // apply position correction
73 fApplySectorAlignment(1), // apply sector alignment
74 fApplyRPhiCorrection(1), // apply R-Phi correction
75 fApplyRCorrection(1) // apply Radial correction
76
9dcfce73 77{
78 //
79 // Constructor
80 //
81}
82
83
84AliTPCcalibCalib::AliTPCcalibCalib(const Text_t *name, const Text_t *title)
7af539c6 85 :AliTPCcalibBase(),
86 fApplyExBCorrection(1), // apply ExB correction
87 fApplyTOFCorrection(1), // apply TOF correction
88 fApplyPositionCorrection(1), // apply position correction
89 fApplySectorAlignment(1), // apply sector alignment
90 fApplyRPhiCorrection(1), // apply R-Phi correction
91 fApplyRCorrection(1) // apply Radial correction
9dcfce73 92{
93 SetName(name);
94 SetTitle(title);
95}
96
97
98AliTPCcalibCalib::AliTPCcalibCalib(const AliTPCcalibCalib&calib):
7af539c6 99 AliTPCcalibBase(calib),
100 fApplyExBCorrection(calib.GetApplyExBCorrection()),
101 fApplyTOFCorrection(calib.GetApplyTOFCorrection()),
102 fApplyPositionCorrection(calib.GetApplyPositionCorrection()),
103 fApplySectorAlignment(calib.GetApplySectorAlignment()),
104 fApplyRPhiCorrection(calib.GetApplyRPhiCorrection()),
105 fApplyRCorrection(calib.GetApplyRCorrection())
106
9dcfce73 107{
108 //
109 // copy constructor
110 //
111}
112
113AliTPCcalibCalib &AliTPCcalibCalib::operator=(const AliTPCcalibCalib&calib){
114 //
115 //
116 //
117 ((AliTPCcalibBase *)this)->operator=(calib);
118 return *this;
119}
120
121
122AliTPCcalibCalib::~AliTPCcalibCalib() {
123 //
124 // destructor
125 //
126}
127
128
129void AliTPCcalibCalib::Process(AliESDEvent *event){
130 //
131 //
132 //
133 if (!event) {
134 return;
135 }
136 AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
137 if (!ESDfriend) {
138 return;
139 }
140
141 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
142 Int_t ntracks=event->GetNumberOfTracks();
91350f5b 143 //AliTPCcalibDB::Instance()->SetExBField(fMagF);
5c3e0d17 144
9dcfce73 145 //
146 //
147 //
148
149 for (Int_t i=0;i<ntracks;++i) {
150 AliESDtrack *track = event->GetTrack(i);
151 const AliExternalTrackParam * trackIn = track->GetInnerParam();
152 const AliExternalTrackParam * trackOut = track->GetOuterParam();
153 if (!trackIn) continue;
154 if (!trackOut) continue;
155
156 AliESDfriendTrack *friendTrack = ESDfriend->GetTrack(i);
157 TObject *calibObject;
158 AliTPCseed *seed = 0;
159 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
160 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
161 }
162 if (!seed) continue;
38b1a1ca 163 RefitTrack(track, seed,event->GetMagneticField());
9dcfce73 164 }
165 return;
166}
167
38b1a1ca 168Bool_t AliTPCcalibCalib::RefitTrack(AliESDtrack * track, AliTPCseed *seed, Float_t magesd){
9dcfce73 169 //
170 // Refit track
38b1a1ca 171 // if magesd==0 forget the curvature
9dcfce73 172
15e48021 173 //
174 // 0 - Setup transform object
175 //
176 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
177 transform->SetCurrentRun(fRun);
178 transform->SetCurrentTimeStamp((UInt_t)fTime);
7af539c6 179 if(!fApplyExBCorrection) { // disable ExB correction in transform
180 if(transform->GetCurrentRecoParam())
181 transform->GetCurrentRecoParamNonConst()->SetUseExBCorrection(0);
182 }
183 if(!fApplyTOFCorrection) { // disable TOF correction in transform
184 if(transform->GetCurrentRecoParam())
185 transform->GetCurrentRecoParamNonConst()->SetUseTOFCorrection(kFALSE);
186 }
187
9dcfce73 188 //
189 // First apply calibration
190 //
15e48021 191 AliTPCPointCorrection * corr = AliTPCPointCorrection::Instance();
9dcfce73 192 for (Int_t irow=0;irow<159;irow++) {
193 AliTPCclusterMI *cluster=seed->GetClusterPointer(irow);
194 if (!cluster) continue;
195 AliTPCclusterMI cl0(*cluster);
196 Double_t x[3]={cluster->GetRow(),cluster->GetPad(),cluster->GetTimeBin()};
197 Int_t i[1]={cluster->GetDetector()};
7af539c6 198
9dcfce73 199 transform->Transform(x,i,0,1);
11c7a854 200 //
201 // get position correction
202 //
203 Int_t ipad=0;
204 if (cluster->GetDetector()>35) ipad=1;
205 Float_t dy =AliTPCClusterParam::SPosCorrection(0,ipad,cluster->GetPad(),cluster->GetTimeBin(),cluster->GetZ(),cluster->GetSigmaY2(),cluster->GetSigmaZ2(),cluster->GetMax());
206 Float_t dz =AliTPCClusterParam::SPosCorrection(1,ipad,cluster->GetPad(),cluster->GetTimeBin(),cluster->GetZ(),cluster->GetSigmaY2(),cluster->GetSigmaZ2(),cluster->GetMax());
7af539c6 207 // if(fApplyPositionCorrection) {
11c7a854 208 //x[1]-=dy;
209 //x[2]-=dz;
7af539c6 210 // }
15e48021 211 //
212 // Apply sector alignment
213 //
214 Double_t dxq = AliTPCPointCorrection::SGetCorrectionSector(0,cluster->GetDetector()%36,cluster->GetX(),
215 cluster->GetY(),cluster->GetZ());
216 Double_t dyq = AliTPCPointCorrection::SGetCorrectionSector(1,cluster->GetDetector()%36,cluster->GetX(),
217 cluster->GetY(),cluster->GetZ());
218 Double_t dzq = AliTPCPointCorrection::SGetCorrectionSector(2,cluster->GetDetector()%36,cluster->GetX(),
219 cluster->GetY(),cluster->GetZ());
7af539c6 220 if (fApplySectorAlignment){
15e48021 221 x[0]-=dxq;
222 x[1]-=dyq;
223 x[2]-=dzq;
224 }
225 //
226 // Apply r-phi correction - To be done on track level- knowing the track angle !!!
227 //
228 Double_t corrclY =
229 corr->RPhiCOGCorrection(cluster->GetDetector(),cluster->GetRow(), cluster->GetPad(),
230 cluster->GetY(),cluster->GetY(), cluster->GetZ(), 0., cluster->GetMax(),2.5);
231 // R correction
232 Double_t corrR = corr->CorrectionOutR0(kFALSE,kFALSE,cluster->GetX(),cluster->GetY(),cluster->GetZ(),cluster->GetDetector());
233
7af539c6 234 if (fApplyRPhiCorrection){
15e48021 235 if (cluster->GetY()>0) x[1]+=corrclY; // rphi correction
236 if (cluster->GetY()<0) x[1]-=corrclY; // rphi correction
7af539c6 237 }
238 if (fApplyRCorrection){
15e48021 239 x[0]+=corrR; // radial correction
240 }
11c7a854 241
15e48021 242 //
243 //
9dcfce73 244 //
245 cluster->SetX(x[0]);
246 cluster->SetY(x[1]);
247 cluster->SetZ(x[2]);
248 if (fStreamLevel>2){
249 TTreeSRedirector *cstream = GetDebugStreamer();
250 if (cstream){
251 (*cstream)<<"Clusters"<<
108953e9 252 "run="<<fRun<< // run number
253 "event="<<fEvent<< // event number
254 "time="<<fTime<< // time stamp of event
255 "trigger="<<fTrigger<< // trigger
15e48021 256 "triggerClass="<<&fTriggerClass<< // trigger
108953e9 257 "mag="<<fMagF<< // magnetic field
9dcfce73 258 "cl0.="<<&cl0<<
259 "cl.="<<cluster<<
11c7a854 260 "cy="<<dy<<
261 "cz="<<dz<<
fe8454d9 262 "cY="<<corrclY<<
15e48021 263 "cR="<<corrR<<
264 "dxq="<<dxq<<
265 "dyq="<<dyq<<
266 "dzq="<<dzq<<
9dcfce73 267 "\n";
268 }
269 }
270 }
15e48021 271 //
272 //
273 //
11c7a854 274 Int_t ncl = seed->GetNumberOfClusters();
275 Double_t covar[15];
276 for (Int_t i=0;i<15;i++) covar[i]=0;
277 covar[0]=10.*10.;
278 covar[2]=10.*10.;
279 covar[5]=10.*10./(64.*64.);
280 covar[9]=10.*10./(64.*64.);
38b1a1ca 281 covar[14]=0.3*0.3;
282 if (TMath::Abs(magesd)<0.05) {
283 covar[14]=0.025*0.025;
284 }
9dcfce73 285 //
286 // And now do refit
287 //
11c7a854 288 AliExternalTrackParam * trackInOld = (AliExternalTrackParam*)track->GetInnerParam();
9dcfce73 289 AliExternalTrackParam * trackOutOld = (AliExternalTrackParam*)track->GetOuterParam();
290
291 AliExternalTrackParam trackIn = *trackOutOld;
38b1a1ca 292 if (TMath::Abs(magesd)<0.05) {
293 ((Double_t&)(trackIn.GetParameter()[4]))=0.000000001;
294 }
295 trackIn.ResetCovariance(20.);
11c7a854 296 trackIn.AddCovariance(covar);
9dcfce73 297 Double_t xyz[3];
298 Int_t nclIn=0,nclOut=0;
299 //
15e48021 300 // Refit in
301 //
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();
91f8fa1d 309 if (TMath::Abs(dalpha)>0.01){
310 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
311 }
15e48021 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 (RejectCluster(cl,&trackIn)) continue;
322 nclIn++;
323 trackIn.Update(&r[1],cov);
324 }
325 //
326 AliExternalTrackParam trackOut = trackIn;
38b1a1ca 327 if (TMath::Abs(magesd)<0.05) {
328 ((Double_t&)(trackOut.GetParameter()[4]))=0.000000001;
329 }
330 trackOut.ResetCovariance(20.);
15e48021 331 trackOut.AddCovariance(covar);
332 //
9dcfce73 333 // Refit out
334 //
15e48021 335 //Bool_t lastEdge=kFALSE;
9dcfce73 336 for (Int_t irow=0; irow<160; 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.)-trackOut.GetAlpha();
342
91f8fa1d 343 if (TMath::Abs(dalpha)>0.01){
344 if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
345 }
9dcfce73 346 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
91350f5b 347
348 Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation
349 AliTPCseed::GetError(cl, &trackOut,cov[0],cov[2]);
350 cov[0]*=cov[0];
351 cov[2]*=cov[2];
9dcfce73 352 trackOut.GetXYZ(xyz);
353 Double_t bz = AliTracker::GetBz(xyz);
15e48021 354 if (!trackOut.PropagateTo(r[0],bz)) continue;
42b40d07 355 if (RejectCluster(cl,&trackOut)) continue;
15e48021 356 nclOut++;
357 trackOut.Update(&r[1],cov);
358 //if (cl->GetType()<0) lastEdge=kTRUE;
359 //if (cl->GetType()>=0) lastEdge=kFALSE;
9dcfce73 360 }
361 //
9dcfce73 362 //
15e48021 363 //
364 nclIn=0;
365 trackIn = trackOut;
366 trackIn.ResetCovariance(10.);
367 //
368 // Refit in one more time
369 //
9dcfce73 370 for (Int_t irow=159; irow>0; irow--){
371 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
372 if (!cl) continue;
373 if (cl->GetX()<80) continue;
374 Int_t sector = cl->GetDetector();
375 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
91f8fa1d 376 if (TMath::Abs(dalpha)>0.01){
377 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
378 }
9dcfce73 379 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
380 Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation
91350f5b 381 AliTPCseed::GetError(cl, &trackIn,cov[0],cov[2]);
382 cov[0]*=cov[0];
383 cov[2]*=cov[2];
384 trackIn.GetXYZ(xyz);
9dcfce73 385 Double_t bz = AliTracker::GetBz(xyz);
386
15e48021 387 if (!trackIn.PropagateTo(r[0],bz)) continue;
42b40d07 388 if (RejectCluster(cl,&trackIn)) continue;
15e48021 389 nclIn++;
9dcfce73 390 trackIn.Update(&r[1],cov);
391 }
15e48021 392
393
3bf5c7a6 394 trackIn.Rotate(trackInOld->GetAlpha());
395 trackOut.Rotate(trackOutOld->GetAlpha());
396 //
397 trackInOld->GetXYZ(xyz);
398 Double_t bz = AliTracker::GetBz(xyz);
399 trackIn.PropagateTo(trackInOld->GetX(),bz);
400 //
401 trackOutOld->GetXYZ(xyz);
402 bz = AliTracker::GetBz(xyz);
403 trackOut.PropagateTo(trackOutOld->GetX(),bz);
404
9dcfce73 405
406 if (fStreamLevel>0){
407 TTreeSRedirector *cstream = GetDebugStreamer();
408 if (cstream){
409 (*cstream)<<"Tracks"<<
108953e9 410 "run="<<fRun<< // run number
411 "event="<<fEvent<< // event number
412 "time="<<fTime<< // time stamp of event
413 "trigger="<<fTrigger<< // trigger
15e48021 414 "triggerClass="<<&fTriggerClass<< // trigger
108953e9 415 "mag="<<fMagF<< // magnetic field
9dcfce73 416 "nclIn="<<nclIn<<
417 "nclOut="<<nclOut<<
418 "ncl="<<ncl<<
419 "TrIn0.="<<trackInOld<<
420 "TrOut0.="<<trackOutOld<<
421 "TrIn1.="<<&trackIn<<
422 "TrOut1.="<<&trackOut<<
423 "\n";
424 }
425 }
3bf5c7a6 426 //
981e9de5 427 // And now rewrite ESDtrack and TPC seed
3bf5c7a6 428 //
429
430 (*trackInOld) = trackIn;
431 (*trackOutOld) = trackOut;
432 AliExternalTrackParam *t = &trackIn;
433 track->Set(t->GetX(),t->GetAlpha(),t->GetParameter(),t->GetCovariance());
981e9de5 434 seed->Set(t->GetX(),t->GetAlpha(),t->GetParameter(),t->GetCovariance());
fe8454d9 435 seed->SetNumberOfClusters((nclIn+nclOut)/2);
2cfb8d90 436 return kTRUE;
9dcfce73 437}
438
439
440
441Bool_t AliTPCcalibCalib::RejectCluster(AliTPCclusterMI* cl, AliExternalTrackParam * param){
442 //
443 // check the acceptance of cluster
444 // Cut on edge effects
445 //
15e48021 446 Float_t kEdgeCut=2.5;
447 Float_t kSigmaCut=6;
448
9dcfce73 449 Bool_t isReject = kFALSE;
450 Float_t edgeY = cl->GetX()*TMath::Tan(TMath::Pi()/18);
451 Float_t dist = edgeY - TMath::Abs(cl->GetY());
452 if (param) dist = TMath::Abs(edgeY - TMath::Abs(param->GetY()));
15e48021 453 if (dist<kEdgeCut) isReject=kTRUE;
454
455 Double_t cov[3]={0.01,0.,0.01}; //TODO: correct error parametrisation
456 AliTPCseed::GetError(cl, param,cov[0],cov[2]);
91f8fa1d 457 if (param->GetSigmaY2()<0 || param->GetSigmaZ2()<0){
d3ce44cb 458 AliError("Wrong parameters");
459 return kFALSE;
460 }
15e48021 461 Double_t py = (cl->GetY()-param->GetY())/TMath::Sqrt(cov[0]*cov[0]+param->GetSigmaY2());
462 Double_t pz = (cl->GetZ()-param->GetZ())/TMath::Sqrt(cov[2]*cov[2]+param->GetSigmaZ2());
463 //
464 if ((py*py+pz*pz)>kSigmaCut*kSigmaCut) isReject=kTRUE;
465
9dcfce73 466 return isReject;
467}
468
469
91350f5b 470